Actual source code: plexfem.c

  1: #include <petsc/private/dmpleximpl.h>
  2: #include <petscsf.h>

  4: #include <petscblaslapack.h>
  5: #include <petsc/private/hashsetij.h>
  6: #include <petsc/private/petscfeimpl.h>
  7: #include <petsc/private/petscfvimpl.h>

  9: PetscBool  Clementcite       = PETSC_FALSE;
 10: const char ClementCitation[] = "@article{clement1975approximation,\n"
 11:                                "  title   = {Approximation by finite element functions using local regularization},\n"
 12:                                "  author  = {Philippe Cl{\\'e}ment},\n"
 13:                                "  journal = {Revue fran{\\c{c}}aise d'automatique, informatique, recherche op{\\'e}rationnelle. Analyse num{\\'e}rique},\n"
 14:                                "  volume  = {9},\n"
 15:                                "  number  = {R2},\n"
 16:                                "  pages   = {77--84},\n"
 17:                                "  year    = {1975}\n}\n";

 19: static PetscErrorCode DMPlexConvertPlex(DM dm, DM *plex, PetscBool copy)
 20: {
 21:   PetscBool isPlex;

 23:   PetscFunctionBegin;
 24:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
 25:   if (isPlex) {
 26:     *plex = dm;
 27:     PetscCall(PetscObjectReference((PetscObject)dm));
 28:   } else {
 29:     PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex));
 30:     if (!*plex) {
 31:       PetscCall(DMConvert(dm, DMPLEX, plex));
 32:       PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex));
 33:     } else {
 34:       PetscCall(PetscObjectReference((PetscObject)*plex));
 35:     }
 36:     if (copy) {
 37:       DMSubDomainHookLink link;

 39:       PetscCall(DMCopyDS(dm, *plex));
 40:       PetscCall(DMCopyAuxiliaryVec(dm, *plex));
 41:       /* Run the subdomain hook (this will copy the DMSNES/DMTS) */
 42:       for (link = dm->subdomainhook; link; link = link->next) {
 43:         if (link->ddhook) PetscCall((*link->ddhook)(dm, *plex, link->ctx));
 44:       }
 45:     }
 46:   }
 47:   PetscFunctionReturn(PETSC_SUCCESS);
 48: }

 50: static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx)
 51: {
 52:   PetscFEGeom *geom = (PetscFEGeom *)ctx;

 54:   PetscFunctionBegin;
 55:   PetscCall(PetscFEGeomDestroy(&geom));
 56:   PetscFunctionReturn(PETSC_SUCCESS);
 57: }

 59: static PetscErrorCode DMPlexGetFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
 60: {
 61:   char           composeStr[33] = {0};
 62:   PetscObjectId  id;
 63:   PetscContainer container;

 65:   PetscFunctionBegin;
 66:   PetscCall(PetscObjectGetId((PetscObject)quad, &id));
 67:   PetscCall(PetscSNPrintf(composeStr, 32, "DMPlexGetFEGeom_%" PetscInt64_FMT "\n", id));
 68:   PetscCall(PetscObjectQuery((PetscObject)pointIS, composeStr, (PetscObject *)&container));
 69:   if (container) {
 70:     PetscCall(PetscContainerGetPointer(container, (void **)geom));
 71:   } else {
 72:     PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, faceData, geom));
 73:     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
 74:     PetscCall(PetscContainerSetPointer(container, (void *)*geom));
 75:     PetscCall(PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom));
 76:     PetscCall(PetscObjectCompose((PetscObject)pointIS, composeStr, (PetscObject)container));
 77:     PetscCall(PetscContainerDestroy(&container));
 78:   }
 79:   PetscFunctionReturn(PETSC_SUCCESS);
 80: }

 82: static PetscErrorCode DMPlexRestoreFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
 83: {
 84:   PetscFunctionBegin;
 85:   *geom = NULL;
 86:   PetscFunctionReturn(PETSC_SUCCESS);
 87: }

 89: /*@
 90:   DMPlexGetScale - Get the scale for the specified fundamental unit

 92:   Not Collective

 94:   Input Parameters:
 95: + dm   - the `DM`
 96: - unit - The SI unit

 98:   Output Parameter:
 99: . scale - The value used to scale all quantities with this unit

101:   Level: advanced

103: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetScale()`, `PetscUnit`
104: @*/
105: PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
106: {
107:   DM_Plex *mesh = (DM_Plex *)dm->data;

109:   PetscFunctionBegin;
111:   PetscAssertPointer(scale, 3);
112:   *scale = mesh->scale[unit];
113:   PetscFunctionReturn(PETSC_SUCCESS);
114: }

116: /*@
117:   DMPlexSetScale - Set the scale for the specified fundamental unit

119:   Not Collective

121:   Input Parameters:
122: + dm    - the `DM`
123: . unit  - The SI unit
124: - scale - The value used to scale all quantities with this unit

126:   Level: advanced

128: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetScale()`, `PetscUnit`
129: @*/
130: PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
131: {
132:   DM_Plex *mesh = (DM_Plex *)dm->data;

134:   PetscFunctionBegin;
136:   mesh->scale[unit] = scale;
137:   PetscFunctionReturn(PETSC_SUCCESS);
138: }

140: PetscErrorCode DMPlexGetUseCeed_Plex(DM dm, PetscBool *useCeed)
141: {
142:   DM_Plex *mesh = (DM_Plex *)dm->data;

144:   PetscFunctionBegin;
145:   *useCeed = mesh->useCeed;
146:   PetscFunctionReturn(PETSC_SUCCESS);
147: }
148: PetscErrorCode DMPlexSetUseCeed_Plex(DM dm, PetscBool useCeed)
149: {
150:   DM_Plex *mesh = (DM_Plex *)dm->data;

152:   PetscFunctionBegin;
153:   mesh->useCeed = useCeed;
154:   PetscFunctionReturn(PETSC_SUCCESS);
155: }

157: /*@
158:   DMPlexGetUseCeed - Get flag for using the LibCEED backend

160:   Not collective

162:   Input Parameter:
163: . dm - The `DM`

165:   Output Parameter:
166: . useCeed - The flag

168:   Level: intermediate

170: .seealso: `DMPlexSetUseCeed()`
171: @*/
172: PetscErrorCode DMPlexGetUseCeed(DM dm, PetscBool *useCeed)
173: {
174:   PetscFunctionBegin;
176:   PetscAssertPointer(useCeed, 2);
177:   *useCeed = PETSC_FALSE;
178:   PetscTryMethod(dm, "DMPlexGetUseCeed_C", (DM, PetscBool *), (dm, useCeed));
179:   PetscFunctionReturn(PETSC_SUCCESS);
180: }

182: /*@
183:   DMPlexSetUseCeed - Set flag for using the LibCEED backend

185:   Not collective

187:   Input Parameters:
188: + dm      - The `DM`
189: - useCeed - The flag

191:   Level: intermediate

193: .seealso: `DMPlexGetUseCeed()`
194: @*/
195: PetscErrorCode DMPlexSetUseCeed(DM dm, PetscBool useCeed)
196: {
197:   PetscFunctionBegin;
200:   PetscUseMethod(dm, "DMPlexSetUseCeed_C", (DM, PetscBool), (dm, useCeed));
201:   PetscFunctionReturn(PETSC_SUCCESS);
202: }

204: /*@
205:   DMPlexGetUseMatClosurePermutation - Get flag for using a closure permutation for matrix insertion

207:   Not collective

209:   Input Parameter:
210: . dm - The `DM`

212:   Output Parameter:
213: . useClPerm - The flag

215:   Level: intermediate

217: .seealso: `DMPlexSetUseMatClosurePermutation()`
218: @*/
219: PetscErrorCode DMPlexGetUseMatClosurePermutation(DM dm, PetscBool *useClPerm)
220: {
221:   DM_Plex *mesh = (DM_Plex *)dm->data;

223:   PetscFunctionBegin;
225:   PetscAssertPointer(useClPerm, 2);
226:   *useClPerm = mesh->useMatClPerm;
227:   PetscFunctionReturn(PETSC_SUCCESS);
228: }

230: /*@
231:   DMPlexSetUseMatClosurePermutation - Set flag for using a closure permutation for matrix insertion

233:   Not collective

235:   Input Parameters:
236: + dm        - The `DM`
237: - useClPerm - The flag

239:   Level: intermediate

241: .seealso: `DMPlexGetUseMatClosurePermutation()`
242: @*/
243: PetscErrorCode DMPlexSetUseMatClosurePermutation(DM dm, PetscBool useClPerm)
244: {
245:   DM_Plex *mesh = (DM_Plex *)dm->data;

247:   PetscFunctionBegin;
250:   mesh->useMatClPerm = useClPerm;
251:   PetscFunctionReturn(PETSC_SUCCESS);
252: }

254: static PetscErrorCode DMPlexProjectRigidBody_Private(PetscInt dim, PetscReal t, const PetscReal X[], PetscInt Nc, PetscScalar *mode, void *ctx)
255: {
256:   const PetscInt eps[3][3][3] = {
257:     {{0, 0, 0},  {0, 0, 1},  {0, -1, 0}},
258:     {{0, 0, -1}, {0, 0, 0},  {1, 0, 0} },
259:     {{0, 1, 0},  {-1, 0, 0}, {0, 0, 0} }
260:   };
261:   PetscInt *ctxInt = (PetscInt *)ctx;
262:   PetscInt  dim2   = ctxInt[0];
263:   PetscInt  d      = ctxInt[1];
264:   PetscInt  i, j, k = dim > 2 ? d - dim : d;

266:   PetscFunctionBegin;
267:   PetscCheck(dim == dim2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input dimension %" PetscInt_FMT " does not match context dimension %" PetscInt_FMT, dim, dim2);
268:   for (i = 0; i < dim; i++) mode[i] = 0.;
269:   if (d < dim) {
270:     mode[d] = 1.; /* Translation along axis d */
271:   } else {
272:     for (i = 0; i < dim; i++) {
273:       for (j = 0; j < dim; j++) { mode[j] += eps[i][j][k] * X[i]; /* Rotation about axis d */ }
274:     }
275:   }
276:   PetscFunctionReturn(PETSC_SUCCESS);
277: }

279: /*@
280:   DMPlexCreateRigidBody - For the default global section, create rigid body modes by function space interpolation

282:   Collective

284:   Input Parameters:
285: + dm    - the `DM`
286: - field - The field number for the rigid body space, or 0 for the default

288:   Output Parameter:
289: . sp - the null space

291:   Level: advanced

293:   Note:
294:   This is necessary to provide a suitable coarse space for algebraic multigrid

296: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `MatNullSpaceCreate()`, `PCGAMG`
297: @*/
298: PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscInt field, MatNullSpace *sp)
299: {
300:   PetscErrorCode (**func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *);
301:   MPI_Comm     comm;
302:   Vec          mode[6];
303:   PetscSection section, globalSection;
304:   PetscInt     dim, dimEmbed, Nf, n, m, mmin, d, i, j;
305:   void       **ctxs;

307:   PetscFunctionBegin;
308:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
309:   PetscCall(DMGetDimension(dm, &dim));
310:   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
311:   PetscCall(DMGetNumFields(dm, &Nf));
312:   PetscCheck(!Nf || !(field < 0 || field >= Nf), comm, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", field, Nf);
313:   if (dim == 1 && Nf < 2) {
314:     PetscCall(MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp));
315:     PetscFunctionReturn(PETSC_SUCCESS);
316:   }
317:   PetscCall(DMGetLocalSection(dm, &section));
318:   PetscCall(DMGetGlobalSection(dm, &globalSection));
319:   PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &n));
320:   PetscCall(PetscCalloc2(Nf, &func, Nf, &ctxs));
321:   m = (dim * (dim + 1)) / 2;
322:   PetscCall(VecCreate(comm, &mode[0]));
323:   PetscCall(VecSetType(mode[0], dm->vectype));
324:   PetscCall(VecSetSizes(mode[0], n, PETSC_DETERMINE));
325:   PetscCall(VecSetUp(mode[0]));
326:   PetscCall(VecGetSize(mode[0], &n));
327:   mmin        = PetscMin(m, n);
328:   func[field] = DMPlexProjectRigidBody_Private;
329:   for (i = 1; i < m; ++i) PetscCall(VecDuplicate(mode[0], &mode[i]));
330:   for (d = 0; d < m; d++) {
331:     PetscInt ctx[2];

333:     ctxs[field] = (void *)(&ctx[0]);
334:     ctx[0]      = dimEmbed;
335:     ctx[1]      = d;
336:     PetscCall(DMProjectFunction(dm, 0.0, func, ctxs, INSERT_VALUES, mode[d]));
337:   }
338:   /* Orthonormalize system */
339:   for (i = 0; i < mmin; ++i) {
340:     PetscScalar dots[6];

342:     PetscCall(VecNormalize(mode[i], NULL));
343:     PetscCall(VecMDot(mode[i], mmin - i - 1, mode + i + 1, dots + i + 1));
344:     for (j = i + 1; j < mmin; ++j) {
345:       dots[j] *= -1.0;
346:       PetscCall(VecAXPY(mode[j], dots[j], mode[i]));
347:     }
348:   }
349:   PetscCall(MatNullSpaceCreate(comm, PETSC_FALSE, mmin, mode, sp));
350:   for (i = 0; i < m; ++i) PetscCall(VecDestroy(&mode[i]));
351:   PetscCall(PetscFree2(func, ctxs));
352:   PetscFunctionReturn(PETSC_SUCCESS);
353: }

355: /*@
356:   DMPlexCreateRigidBodies - For the default global section, create rigid body modes by function space interpolation

358:   Collective

360:   Input Parameters:
361: + dm    - the `DM`
362: . nb    - The number of bodies
363: . label - The `DMLabel` marking each domain
364: . nids  - The number of ids per body
365: - ids   - An array of the label ids in sequence for each domain

367:   Output Parameter:
368: . sp - the null space

370:   Level: advanced

372:   Note:
373:   This is necessary to provide a suitable coarse space for algebraic multigrid

375: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `MatNullSpaceCreate()`
376: @*/
377: PetscErrorCode DMPlexCreateRigidBodies(DM dm, PetscInt nb, DMLabel label, const PetscInt nids[], const PetscInt ids[], MatNullSpace *sp)
378: {
379:   MPI_Comm     comm;
380:   PetscSection section, globalSection;
381:   Vec         *mode;
382:   PetscScalar *dots;
383:   PetscInt     dim, dimEmbed, n, m, b, d, i, j, off;

385:   PetscFunctionBegin;
386:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
387:   PetscCall(DMGetDimension(dm, &dim));
388:   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
389:   PetscCall(DMGetLocalSection(dm, &section));
390:   PetscCall(DMGetGlobalSection(dm, &globalSection));
391:   PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &n));
392:   m = nb * (dim * (dim + 1)) / 2;
393:   PetscCall(PetscMalloc2(m, &mode, m, &dots));
394:   PetscCall(VecCreate(comm, &mode[0]));
395:   PetscCall(VecSetSizes(mode[0], n, PETSC_DETERMINE));
396:   PetscCall(VecSetUp(mode[0]));
397:   for (i = 1; i < m; ++i) PetscCall(VecDuplicate(mode[0], &mode[i]));
398:   for (b = 0, off = 0; b < nb; ++b) {
399:     for (d = 0; d < m / nb; ++d) {
400:       PetscInt ctx[2];
401:       PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody_Private;
402:       void *voidctx                                                                                   = (void *)(&ctx[0]);

404:       ctx[0] = dimEmbed;
405:       ctx[1] = d;
406:       PetscCall(DMProjectFunctionLabel(dm, 0.0, label, nids[b], &ids[off], 0, NULL, &func, &voidctx, INSERT_VALUES, mode[d]));
407:       off += nids[b];
408:     }
409:   }
410:   /* Orthonormalize system */
411:   for (i = 0; i < m; ++i) {
412:     PetscScalar dots[6];

414:     PetscCall(VecNormalize(mode[i], NULL));
415:     PetscCall(VecMDot(mode[i], m - i - 1, mode + i + 1, dots + i + 1));
416:     for (j = i + 1; j < m; ++j) {
417:       dots[j] *= -1.0;
418:       PetscCall(VecAXPY(mode[j], dots[j], mode[i]));
419:     }
420:   }
421:   PetscCall(MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp));
422:   for (i = 0; i < m; ++i) PetscCall(VecDestroy(&mode[i]));
423:   PetscCall(PetscFree2(mode, dots));
424:   PetscFunctionReturn(PETSC_SUCCESS);
425: }

427: /*@
428:   DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs
429:   are computed by associating the basis function with one of the mesh points in its transitively-closed support, and
430:   evaluating the dual space basis of that point.

432:   Input Parameters:
433: + dm     - the `DMPLEX` object
434: - height - the maximum projection height >= 0

436:   Level: advanced

438:   Notes:
439:   A basis function is associated with the point in its transitively-closed support whose mesh
440:   height is highest (w.r.t. DAG height), but not greater than the maximum projection height,
441:   which is set with this function.  By default, the maximum projection height is zero, which
442:   means that only mesh cells are used to project basis functions.  A height of one, for
443:   example, evaluates a cell-interior basis functions using its cells dual space basis, but all
444:   other basis functions with the dual space basis of a face.

446: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetMaxProjectionHeight()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabelLocal()`
447: @*/
448: PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
449: {
450:   DM_Plex *plex = (DM_Plex *)dm->data;

452:   PetscFunctionBegin;
454:   plex->maxProjectionHeight = height;
455:   PetscFunctionReturn(PETSC_SUCCESS);
456: }

458: /*@
459:   DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in
460:   DMPlexProjectXXXLocal() functions.

462:   Input Parameter:
463: . dm - the `DMPLEX` object

465:   Output Parameter:
466: . height - the maximum projection height

468:   Level: intermediate

470: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetMaxProjectionHeight()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabelLocal()`
471: @*/
472: PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
473: {
474:   DM_Plex *plex = (DM_Plex *)dm->data;

476:   PetscFunctionBegin;
478:   *height = plex->maxProjectionHeight;
479:   PetscFunctionReturn(PETSC_SUCCESS);
480: }

482: typedef struct {
483:   PetscReal    alpha; /* The first Euler angle, and in 2D the only one */
484:   PetscReal    beta;  /* The second Euler angle */
485:   PetscReal    gamma; /* The third Euler angle */
486:   PetscInt     dim;   /* The dimension of R */
487:   PetscScalar *R;     /* The rotation matrix, transforming a vector in the local basis to the global basis */
488:   PetscScalar *RT;    /* The transposed rotation matrix, transforming a vector in the global basis to the local basis */
489: } RotCtx;

491: /*
492:   Note: Following https://en.wikipedia.org/wiki/Euler_angles, we will specify Euler angles by extrinsic rotations, meaning that
493:   we rotate with respect to a fixed initial coordinate system, the local basis (x-y-z). The global basis (X-Y-Z) is reached as follows:
494:   $ The XYZ system rotates about the z axis by alpha. The X axis is now at angle alpha with respect to the x axis.
495:   $ The XYZ system rotates again about the x axis by beta. The Z axis is now at angle beta with respect to the z axis.
496:   $ The XYZ system rotates a third time about the z axis by gamma.
497: */
498: static PetscErrorCode DMPlexBasisTransformSetUp_Rotation_Internal(DM dm, void *ctx)
499: {
500:   RotCtx   *rc  = (RotCtx *)ctx;
501:   PetscInt  dim = rc->dim;
502:   PetscReal c1, s1, c2, s2, c3, s3;

504:   PetscFunctionBegin;
505:   PetscCall(PetscMalloc2(PetscSqr(dim), &rc->R, PetscSqr(dim), &rc->RT));
506:   switch (dim) {
507:   case 2:
508:     c1       = PetscCosReal(rc->alpha);
509:     s1       = PetscSinReal(rc->alpha);
510:     rc->R[0] = c1;
511:     rc->R[1] = s1;
512:     rc->R[2] = -s1;
513:     rc->R[3] = c1;
514:     PetscCall(PetscArraycpy(rc->RT, rc->R, PetscSqr(dim)));
515:     DMPlex_Transpose2D_Internal(rc->RT);
516:     break;
517:   case 3:
518:     c1       = PetscCosReal(rc->alpha);
519:     s1       = PetscSinReal(rc->alpha);
520:     c2       = PetscCosReal(rc->beta);
521:     s2       = PetscSinReal(rc->beta);
522:     c3       = PetscCosReal(rc->gamma);
523:     s3       = PetscSinReal(rc->gamma);
524:     rc->R[0] = c1 * c3 - c2 * s1 * s3;
525:     rc->R[1] = c3 * s1 + c1 * c2 * s3;
526:     rc->R[2] = s2 * s3;
527:     rc->R[3] = -c1 * s3 - c2 * c3 * s1;
528:     rc->R[4] = c1 * c2 * c3 - s1 * s3;
529:     rc->R[5] = c3 * s2;
530:     rc->R[6] = s1 * s2;
531:     rc->R[7] = -c1 * s2;
532:     rc->R[8] = c2;
533:     PetscCall(PetscArraycpy(rc->RT, rc->R, PetscSqr(dim)));
534:     DMPlex_Transpose3D_Internal(rc->RT);
535:     break;
536:   default:
537:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported", dim);
538:   }
539:   PetscFunctionReturn(PETSC_SUCCESS);
540: }

542: static PetscErrorCode DMPlexBasisTransformDestroy_Rotation_Internal(DM dm, void *ctx)
543: {
544:   RotCtx *rc = (RotCtx *)ctx;

546:   PetscFunctionBegin;
547:   PetscCall(PetscFree2(rc->R, rc->RT));
548:   PetscCall(PetscFree(rc));
549:   PetscFunctionReturn(PETSC_SUCCESS);
550: }

552: static PetscErrorCode DMPlexBasisTransformGetMatrix_Rotation_Internal(DM dm, const PetscReal x[], PetscBool l2g, const PetscScalar **A, void *ctx)
553: {
554:   RotCtx *rc = (RotCtx *)ctx;

556:   PetscFunctionBeginHot;
557:   PetscAssertPointer(ctx, 5);
558:   if (l2g) {
559:     *A = rc->R;
560:   } else {
561:     *A = rc->RT;
562:   }
563:   PetscFunctionReturn(PETSC_SUCCESS);
564: }

566: PetscErrorCode DMPlexBasisTransformApplyReal_Internal(DM dm, const PetscReal x[], PetscBool l2g, PetscInt dim, const PetscReal *y, PetscReal *z, void *ctx)
567: {
568:   PetscFunctionBegin;
569: #if defined(PETSC_USE_COMPLEX)
570:   switch (dim) {
571:   case 2: {
572:     PetscScalar yt[2] = {y[0], y[1]}, zt[2] = {0.0, 0.0};

574:     PetscCall(DMPlexBasisTransformApply_Internal(dm, x, l2g, dim, yt, zt, ctx));
575:     z[0] = PetscRealPart(zt[0]);
576:     z[1] = PetscRealPart(zt[1]);
577:   } break;
578:   case 3: {
579:     PetscScalar yt[3] = {y[0], y[1], y[2]}, zt[3] = {0.0, 0.0, 0.0};

581:     PetscCall(DMPlexBasisTransformApply_Internal(dm, x, l2g, dim, yt, zt, ctx));
582:     z[0] = PetscRealPart(zt[0]);
583:     z[1] = PetscRealPart(zt[1]);
584:     z[2] = PetscRealPart(zt[2]);
585:   } break;
586:   }
587: #else
588:   PetscCall(DMPlexBasisTransformApply_Internal(dm, x, l2g, dim, y, z, ctx));
589: #endif
590:   PetscFunctionReturn(PETSC_SUCCESS);
591: }

593: PetscErrorCode DMPlexBasisTransformApply_Internal(DM dm, const PetscReal x[], PetscBool l2g, PetscInt dim, const PetscScalar *y, PetscScalar *z, void *ctx)
594: {
595:   const PetscScalar *A;

597:   PetscFunctionBeginHot;
598:   PetscCall((*dm->transformGetMatrix)(dm, x, l2g, &A, ctx));
599:   switch (dim) {
600:   case 2:
601:     DMPlex_Mult2D_Internal(A, 1, y, z);
602:     break;
603:   case 3:
604:     DMPlex_Mult3D_Internal(A, 1, y, z);
605:     break;
606:   }
607:   PetscFunctionReturn(PETSC_SUCCESS);
608: }

610: static PetscErrorCode DMPlexBasisTransformField_Internal(DM dm, DM tdm, Vec tv, PetscInt p, PetscInt f, PetscBool l2g, PetscScalar *a)
611: {
612:   PetscSection       ts;
613:   const PetscScalar *ta, *tva;
614:   PetscInt           dof;

616:   PetscFunctionBeginHot;
617:   PetscCall(DMGetLocalSection(tdm, &ts));
618:   PetscCall(PetscSectionGetFieldDof(ts, p, f, &dof));
619:   PetscCall(VecGetArrayRead(tv, &ta));
620:   PetscCall(DMPlexPointLocalFieldRead(tdm, p, f, ta, &tva));
621:   if (l2g) {
622:     switch (dof) {
623:     case 4:
624:       DMPlex_Mult2D_Internal(tva, 1, a, a);
625:       break;
626:     case 9:
627:       DMPlex_Mult3D_Internal(tva, 1, a, a);
628:       break;
629:     }
630:   } else {
631:     switch (dof) {
632:     case 4:
633:       DMPlex_MultTranspose2D_Internal(tva, 1, a, a);
634:       break;
635:     case 9:
636:       DMPlex_MultTranspose3D_Internal(tva, 1, a, a);
637:       break;
638:     }
639:   }
640:   PetscCall(VecRestoreArrayRead(tv, &ta));
641:   PetscFunctionReturn(PETSC_SUCCESS);
642: }

644: static PetscErrorCode DMPlexBasisTransformFieldTensor_Internal(DM dm, DM tdm, Vec tv, PetscInt pf, PetscInt f, PetscInt pg, PetscInt g, PetscBool l2g, PetscInt lda, PetscScalar *a)
645: {
646:   PetscSection       s, ts;
647:   const PetscScalar *ta, *tvaf, *tvag;
648:   PetscInt           fdof, gdof, fpdof, gpdof;

650:   PetscFunctionBeginHot;
651:   PetscCall(DMGetLocalSection(dm, &s));
652:   PetscCall(DMGetLocalSection(tdm, &ts));
653:   PetscCall(PetscSectionGetFieldDof(s, pf, f, &fpdof));
654:   PetscCall(PetscSectionGetFieldDof(s, pg, g, &gpdof));
655:   PetscCall(PetscSectionGetFieldDof(ts, pf, f, &fdof));
656:   PetscCall(PetscSectionGetFieldDof(ts, pg, g, &gdof));
657:   PetscCall(VecGetArrayRead(tv, &ta));
658:   PetscCall(DMPlexPointLocalFieldRead(tdm, pf, f, ta, &tvaf));
659:   PetscCall(DMPlexPointLocalFieldRead(tdm, pg, g, ta, &tvag));
660:   if (l2g) {
661:     switch (fdof) {
662:     case 4:
663:       DMPlex_MatMult2D_Internal(tvaf, gpdof, lda, a, a);
664:       break;
665:     case 9:
666:       DMPlex_MatMult3D_Internal(tvaf, gpdof, lda, a, a);
667:       break;
668:     }
669:     switch (gdof) {
670:     case 4:
671:       DMPlex_MatMultTransposeLeft2D_Internal(tvag, fpdof, lda, a, a);
672:       break;
673:     case 9:
674:       DMPlex_MatMultTransposeLeft3D_Internal(tvag, fpdof, lda, a, a);
675:       break;
676:     }
677:   } else {
678:     switch (fdof) {
679:     case 4:
680:       DMPlex_MatMultTranspose2D_Internal(tvaf, gpdof, lda, a, a);
681:       break;
682:     case 9:
683:       DMPlex_MatMultTranspose3D_Internal(tvaf, gpdof, lda, a, a);
684:       break;
685:     }
686:     switch (gdof) {
687:     case 4:
688:       DMPlex_MatMultLeft2D_Internal(tvag, fpdof, lda, a, a);
689:       break;
690:     case 9:
691:       DMPlex_MatMultLeft3D_Internal(tvag, fpdof, lda, a, a);
692:       break;
693:     }
694:   }
695:   PetscCall(VecRestoreArrayRead(tv, &ta));
696:   PetscFunctionReturn(PETSC_SUCCESS);
697: }

699: PetscErrorCode DMPlexBasisTransformPoint_Internal(DM dm, DM tdm, Vec tv, PetscInt p, PetscBool fieldActive[], PetscBool l2g, PetscScalar *a)
700: {
701:   PetscSection    s;
702:   PetscSection    clSection;
703:   IS              clPoints;
704:   const PetscInt *clp;
705:   PetscInt       *points = NULL;
706:   PetscInt        Nf, f, Np, cp, dof, d = 0;

708:   PetscFunctionBegin;
709:   PetscCall(DMGetLocalSection(dm, &s));
710:   PetscCall(PetscSectionGetNumFields(s, &Nf));
711:   PetscCall(DMPlexGetCompressedClosure(dm, s, p, 0, &Np, &points, &clSection, &clPoints, &clp));
712:   for (f = 0; f < Nf; ++f) {
713:     for (cp = 0; cp < Np * 2; cp += 2) {
714:       PetscCall(PetscSectionGetFieldDof(s, points[cp], f, &dof));
715:       if (!dof) continue;
716:       if (fieldActive[f]) PetscCall(DMPlexBasisTransformField_Internal(dm, tdm, tv, points[cp], f, l2g, &a[d]));
717:       d += dof;
718:     }
719:   }
720:   PetscCall(DMPlexRestoreCompressedClosure(dm, s, p, &Np, &points, &clSection, &clPoints, &clp));
721:   PetscFunctionReturn(PETSC_SUCCESS);
722: }

724: PetscErrorCode DMPlexBasisTransformPointTensor_Internal(DM dm, DM tdm, Vec tv, PetscInt p, PetscBool l2g, PetscInt lda, PetscScalar *a)
725: {
726:   PetscSection    s;
727:   PetscSection    clSection;
728:   IS              clPoints;
729:   const PetscInt *clp;
730:   PetscInt       *points = NULL;
731:   PetscInt        Nf, f, g, Np, cpf, cpg, fdof, gdof, r, c = 0;

733:   PetscFunctionBegin;
734:   PetscCall(DMGetLocalSection(dm, &s));
735:   PetscCall(PetscSectionGetNumFields(s, &Nf));
736:   PetscCall(DMPlexGetCompressedClosure(dm, s, p, 0, &Np, &points, &clSection, &clPoints, &clp));
737:   for (f = 0, r = 0; f < Nf; ++f) {
738:     for (cpf = 0; cpf < Np * 2; cpf += 2) {
739:       PetscCall(PetscSectionGetFieldDof(s, points[cpf], f, &fdof));
740:       for (g = 0, c = 0; g < Nf; ++g) {
741:         for (cpg = 0; cpg < Np * 2; cpg += 2) {
742:           PetscCall(PetscSectionGetFieldDof(s, points[cpg], g, &gdof));
743:           PetscCall(DMPlexBasisTransformFieldTensor_Internal(dm, tdm, tv, points[cpf], f, points[cpg], g, l2g, lda, &a[r * lda + c]));
744:           c += gdof;
745:         }
746:       }
747:       PetscCheck(c == lda, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of columns %" PetscInt_FMT " should be %" PetscInt_FMT, c, lda);
748:       r += fdof;
749:     }
750:   }
751:   PetscCheck(r == lda, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of rows %" PetscInt_FMT " should be %" PetscInt_FMT, c, lda);
752:   PetscCall(DMPlexRestoreCompressedClosure(dm, s, p, &Np, &points, &clSection, &clPoints, &clp));
753:   PetscFunctionReturn(PETSC_SUCCESS);
754: }

756: static PetscErrorCode DMPlexBasisTransform_Internal(DM dm, Vec lv, PetscBool l2g)
757: {
758:   DM                 tdm;
759:   Vec                tv;
760:   PetscSection       ts, s;
761:   const PetscScalar *ta;
762:   PetscScalar       *a, *va;
763:   PetscInt           pStart, pEnd, p, Nf, f;

765:   PetscFunctionBegin;
766:   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
767:   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
768:   PetscCall(DMGetLocalSection(tdm, &ts));
769:   PetscCall(DMGetLocalSection(dm, &s));
770:   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
771:   PetscCall(PetscSectionGetNumFields(s, &Nf));
772:   PetscCall(VecGetArray(lv, &a));
773:   PetscCall(VecGetArrayRead(tv, &ta));
774:   for (p = pStart; p < pEnd; ++p) {
775:     for (f = 0; f < Nf; ++f) {
776:       PetscCall(DMPlexPointLocalFieldRef(dm, p, f, a, &va));
777:       PetscCall(DMPlexBasisTransformField_Internal(dm, tdm, tv, p, f, l2g, va));
778:     }
779:   }
780:   PetscCall(VecRestoreArray(lv, &a));
781:   PetscCall(VecRestoreArrayRead(tv, &ta));
782:   PetscFunctionReturn(PETSC_SUCCESS);
783: }

785: /*@
786:   DMPlexGlobalToLocalBasis - Transform the values in the given local vector from the global basis to the local basis

788:   Input Parameters:
789: + dm - The `DM`
790: - lv - A local vector with values in the global basis

792:   Output Parameter:
793: . lv - A local vector with values in the local basis

795:   Level: developer

797:   Note:
798:   This method is only intended to be called inside `DMGlobalToLocal()`. It is unlikely that a user will have a local vector full of coefficients for the global basis unless they are reimplementing GlobalToLocal.

800: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexLocalToGlobalBasis()`, `DMGetLocalSection()`, `DMPlexCreateBasisRotation()`
801: @*/
802: PetscErrorCode DMPlexGlobalToLocalBasis(DM dm, Vec lv)
803: {
804:   PetscFunctionBegin;
807:   PetscCall(DMPlexBasisTransform_Internal(dm, lv, PETSC_FALSE));
808:   PetscFunctionReturn(PETSC_SUCCESS);
809: }

811: /*@
812:   DMPlexLocalToGlobalBasis - Transform the values in the given local vector from the local basis to the global basis

814:   Input Parameters:
815: + dm - The `DM`
816: - lv - A local vector with values in the local basis

818:   Output Parameter:
819: . lv - A local vector with values in the global basis

821:   Level: developer

823:   Note:
824:   This method is only intended to be called inside `DMGlobalToLocal()`. It is unlikely that a user would want a local vector full of coefficients for the global basis unless they are reimplementing GlobalToLocal.

826: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGlobalToLocalBasis()`, `DMGetLocalSection()`, `DMPlexCreateBasisRotation()`
827: @*/
828: PetscErrorCode DMPlexLocalToGlobalBasis(DM dm, Vec lv)
829: {
830:   PetscFunctionBegin;
833:   PetscCall(DMPlexBasisTransform_Internal(dm, lv, PETSC_TRUE));
834:   PetscFunctionReturn(PETSC_SUCCESS);
835: }

837: /*@
838:   DMPlexCreateBasisRotation - Create an internal transformation from the global basis, used to specify boundary conditions
839:   and global solutions, to a local basis, appropriate for discretization integrals and assembly.

841:   Input Parameters:
842: + dm    - The `DM`
843: . alpha - The first Euler angle, and in 2D the only one
844: . beta  - The second Euler angle
845: - gamma - The third Euler angle

847:   Level: developer

849:   Note:
850:   Following https://en.wikipedia.org/wiki/Euler_angles, we will specify Euler angles by extrinsic rotations, meaning that
851:   we rotate with respect to a fixed initial coordinate system, the local basis (x-y-z). The global basis (X-Y-Z) is reached as follows
852: .vb
853:    The XYZ system rotates about the z axis by alpha. The X axis is now at angle alpha with respect to the x axis.
854:    The XYZ system rotates again about the x axis by beta. The Z axis is now at angle beta with respect to the z axis.
855:    The XYZ system rotates a third time about the z axis by gamma.
856: .ve

858: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGlobalToLocalBasis()`, `DMPlexLocalToGlobalBasis()`
859: @*/
860: PetscErrorCode DMPlexCreateBasisRotation(DM dm, PetscReal alpha, PetscReal beta, PetscReal gamma)
861: {
862:   RotCtx  *rc;
863:   PetscInt cdim;

865:   PetscFunctionBegin;
866:   PetscCall(DMGetCoordinateDim(dm, &cdim));
867:   PetscCall(PetscMalloc1(1, &rc));
868:   dm->transformCtx       = rc;
869:   dm->transformSetUp     = DMPlexBasisTransformSetUp_Rotation_Internal;
870:   dm->transformDestroy   = DMPlexBasisTransformDestroy_Rotation_Internal;
871:   dm->transformGetMatrix = DMPlexBasisTransformGetMatrix_Rotation_Internal;
872:   rc->dim                = cdim;
873:   rc->alpha              = alpha;
874:   rc->beta               = beta;
875:   rc->gamma              = gamma;
876:   PetscCall((*dm->transformSetUp)(dm, dm->transformCtx));
877:   PetscCall(DMConstructBasisTransform_Internal(dm));
878:   PetscFunctionReturn(PETSC_SUCCESS);
879: }

881: /*@C
882:   DMPlexInsertBoundaryValuesEssential - Insert boundary values into a local vector using a function of the coordinates

884:   Input Parameters:
885: + dm     - The `DM`, with a `PetscDS` that matches the problem being constrained
886: . time   - The time
887: . field  - The field to constrain
888: . Nc     - The number of constrained field components, or 0 for all components
889: . comps  - An array of constrained component numbers, or `NULL` for all components
890: . label  - The `DMLabel` defining constrained points
891: . numids - The number of `DMLabel` ids for constrained points
892: . ids    - An array of ids for constrained points
893: . func   - A pointwise function giving boundary values
894: - ctx    - An optional user context for bcFunc

896:   Output Parameter:
897: . locX - A local vector to receives the boundary values

899:   Level: developer

901: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLabel`, `DMPlexInsertBoundaryValuesEssentialField()`, `DMPlexInsertBoundaryValuesEssentialBdField()`, `DMAddBoundary()`
902: @*/
903: PetscErrorCode DMPlexInsertBoundaryValuesEssential(DM dm, PetscReal time, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, Vec locX)
904: {
905:   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
906:   void   **ctxs;
907:   PetscInt numFields;

909:   PetscFunctionBegin;
910:   PetscCall(DMGetNumFields(dm, &numFields));
911:   PetscCall(PetscCalloc2(numFields, &funcs, numFields, &ctxs));
912:   funcs[field] = func;
913:   ctxs[field]  = ctx;
914:   PetscCall(DMProjectFunctionLabelLocal(dm, time, label, numids, ids, Nc, comps, funcs, ctxs, INSERT_BC_VALUES, locX));
915:   PetscCall(PetscFree2(funcs, ctxs));
916:   PetscFunctionReturn(PETSC_SUCCESS);
917: }

919: /*@C
920:   DMPlexInsertBoundaryValuesEssentialField - Insert boundary values into a local vector using a function of the coordinates and field data

922:   Input Parameters:
923: + dm     - The `DM`, with a `PetscDS` that matches the problem being constrained
924: . time   - The time
925: . locU   - A local vector with the input solution values
926: . field  - The field to constrain
927: . Nc     - The number of constrained field components, or 0 for all components
928: . comps  - An array of constrained component numbers, or `NULL` for all components
929: . label  - The `DMLabel` defining constrained points
930: . numids - The number of `DMLabel` ids for constrained points
931: . ids    - An array of ids for constrained points
932: . func   - A pointwise function giving boundary values
933: - ctx    - An optional user context for bcFunc

935:   Output Parameter:
936: . locX - A local vector to receives the boundary values

938:   Level: developer

940: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexInsertBoundaryValuesEssential()`, `DMPlexInsertBoundaryValuesEssentialBdField()`, `DMAddBoundary()`
941: @*/
942: PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM dm, PetscReal time, Vec locU, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], void (*func)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void *ctx, Vec locX)
943: {
944:   void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
945:   void   **ctxs;
946:   PetscInt numFields;

948:   PetscFunctionBegin;
949:   PetscCall(DMGetNumFields(dm, &numFields));
950:   PetscCall(PetscCalloc2(numFields, &funcs, numFields, &ctxs));
951:   funcs[field] = func;
952:   ctxs[field]  = ctx;
953:   PetscCall(DMProjectFieldLabelLocal(dm, time, label, numids, ids, Nc, comps, locU, funcs, INSERT_BC_VALUES, locX));
954:   PetscCall(PetscFree2(funcs, ctxs));
955:   PetscFunctionReturn(PETSC_SUCCESS);
956: }

958: /*@C
959:   DMPlexInsertBoundaryValuesEssentialBdField - Insert boundary values into a local vector using a function of the coordinates and boundary field data

961:   Collective

963:   Input Parameters:
964: + dm     - The `DM`, with a `PetscDS` that matches the problem being constrained
965: . time   - The time
966: . locU   - A local vector with the input solution values
967: . field  - The field to constrain
968: . Nc     - The number of constrained field components, or 0 for all components
969: . comps  - An array of constrained component numbers, or `NULL` for all components
970: . label  - The `DMLabel` defining constrained points
971: . numids - The number of `DMLabel` ids for constrained points
972: . ids    - An array of ids for constrained points
973: . func   - A pointwise function giving boundary values, the calling sequence is given in `DMProjectBdFieldLabelLocal()`
974: - ctx    - An optional user context for `func`

976:   Output Parameter:
977: . locX - A local vector to receive the boundary values

979:   Level: developer

981: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectBdFieldLabelLocal()`, `DMPlexInsertBoundaryValuesEssential()`, `DMPlexInsertBoundaryValuesEssentialField()`, `DMAddBoundary()`
982: @*/
983: PetscErrorCode DMPlexInsertBoundaryValuesEssentialBdField(DM dm, PetscReal time, Vec locU, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], void (*func)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void *ctx, Vec locX)
984: {
985:   void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
986:   void   **ctxs;
987:   PetscInt numFields;

989:   PetscFunctionBegin;
990:   PetscCall(DMGetNumFields(dm, &numFields));
991:   PetscCall(PetscCalloc2(numFields, &funcs, numFields, &ctxs));
992:   funcs[field] = func;
993:   ctxs[field]  = ctx;
994:   PetscCall(DMProjectBdFieldLabelLocal(dm, time, label, numids, ids, Nc, comps, locU, funcs, INSERT_BC_VALUES, locX));
995:   PetscCall(PetscFree2(funcs, ctxs));
996:   PetscFunctionReturn(PETSC_SUCCESS);
997: }

999: /*@C
1000:   DMPlexInsertBoundaryValuesRiemann - Insert boundary values into a local vector

1002:   Input Parameters:
1003: + dm           - The `DM`, with a `PetscDS` that matches the problem being constrained
1004: . time         - The time
1005: . faceGeometry - A vector with the FVM face geometry information
1006: . cellGeometry - A vector with the FVM cell geometry information
1007: . Grad         - A vector with the FVM cell gradient information
1008: . field        - The field to constrain
1009: . Nc           - The number of constrained field components, or 0 for all components
1010: . comps        - An array of constrained component numbers, or `NULL` for all components
1011: . label        - The `DMLabel` defining constrained points
1012: . numids       - The number of `DMLabel` ids for constrained points
1013: . ids          - An array of ids for constrained points
1014: . func         - A pointwise function giving boundary values
1015: - ctx          - An optional user context for bcFunc

1017:   Output Parameter:
1018: . locX - A local vector to receives the boundary values

1020:   Level: developer

1022:   Note:
1023:   This implementation currently ignores the numcomps/comps argument from `DMAddBoundary()`

1025: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexInsertBoundaryValuesEssential()`, `DMPlexInsertBoundaryValuesEssentialField()`, `DMAddBoundary()`
1026: @*/
1027: PetscErrorCode DMPlexInsertBoundaryValuesRiemann(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscReal, const PetscReal *, const PetscReal *, const PetscScalar *, PetscScalar *, void *), void *ctx, Vec locX)
1028: {
1029:   PetscDS            prob;
1030:   PetscSF            sf;
1031:   DM                 dmFace, dmCell, dmGrad;
1032:   const PetscScalar *facegeom, *cellgeom = NULL, *grad;
1033:   const PetscInt    *leaves;
1034:   PetscScalar       *x, *fx;
1035:   PetscInt           dim, nleaves, loc, fStart, fEnd, pdim, i;
1036:   PetscErrorCode     ierru = PETSC_SUCCESS;

1038:   PetscFunctionBegin;
1039:   PetscCall(DMGetPointSF(dm, &sf));
1040:   PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL));
1041:   nleaves = PetscMax(0, nleaves);
1042:   PetscCall(DMGetDimension(dm, &dim));
1043:   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
1044:   PetscCall(DMGetDS(dm, &prob));
1045:   PetscCall(VecGetDM(faceGeometry, &dmFace));
1046:   PetscCall(VecGetArrayRead(faceGeometry, &facegeom));
1047:   if (cellGeometry) {
1048:     PetscCall(VecGetDM(cellGeometry, &dmCell));
1049:     PetscCall(VecGetArrayRead(cellGeometry, &cellgeom));
1050:   }
1051:   if (Grad) {
1052:     PetscFV fv;

1054:     PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fv));
1055:     PetscCall(VecGetDM(Grad, &dmGrad));
1056:     PetscCall(VecGetArrayRead(Grad, &grad));
1057:     PetscCall(PetscFVGetNumComponents(fv, &pdim));
1058:     PetscCall(DMGetWorkArray(dm, pdim, MPIU_SCALAR, &fx));
1059:   }
1060:   PetscCall(VecGetArray(locX, &x));
1061:   for (i = 0; i < numids; ++i) {
1062:     IS              faceIS;
1063:     const PetscInt *faces;
1064:     PetscInt        numFaces, f;

1066:     PetscCall(DMLabelGetStratumIS(label, ids[i], &faceIS));
1067:     if (!faceIS) continue; /* No points with that id on this process */
1068:     PetscCall(ISGetLocalSize(faceIS, &numFaces));
1069:     PetscCall(ISGetIndices(faceIS, &faces));
1070:     for (f = 0; f < numFaces; ++f) {
1071:       const PetscInt   face = faces[f], *cells;
1072:       PetscFVFaceGeom *fg;

1074:       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
1075:       PetscCall(PetscFindInt(face, nleaves, (PetscInt *)leaves, &loc));
1076:       if (loc >= 0) continue;
1077:       PetscCall(DMPlexPointLocalRead(dmFace, face, facegeom, &fg));
1078:       PetscCall(DMPlexGetSupport(dm, face, &cells));
1079:       if (Grad) {
1080:         PetscFVCellGeom *cg;
1081:         PetscScalar     *cx, *cgrad;
1082:         PetscScalar     *xG;
1083:         PetscReal        dx[3];
1084:         PetscInt         d;

1086:         PetscCall(DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg));
1087:         PetscCall(DMPlexPointLocalRead(dm, cells[0], x, &cx));
1088:         PetscCall(DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad));
1089:         PetscCall(DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG));
1090:         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
1091:         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d * dim], dx);
1092:         PetscCall((*func)(time, fg->centroid, fg->normal, fx, xG, ctx));
1093:       } else {
1094:         PetscScalar *xI;
1095:         PetscScalar *xG;

1097:         PetscCall(DMPlexPointLocalRead(dm, cells[0], x, &xI));
1098:         PetscCall(DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG));
1099:         ierru = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);
1100:         if (ierru) {
1101:           PetscCall(ISRestoreIndices(faceIS, &faces));
1102:           PetscCall(ISDestroy(&faceIS));
1103:           goto cleanup;
1104:         }
1105:       }
1106:     }
1107:     PetscCall(ISRestoreIndices(faceIS, &faces));
1108:     PetscCall(ISDestroy(&faceIS));
1109:   }
1110: cleanup:
1111:   PetscCall(VecRestoreArray(locX, &x));
1112:   if (Grad) {
1113:     PetscCall(DMRestoreWorkArray(dm, pdim, MPIU_SCALAR, &fx));
1114:     PetscCall(VecRestoreArrayRead(Grad, &grad));
1115:   }
1116:   if (cellGeometry) PetscCall(VecRestoreArrayRead(cellGeometry, &cellgeom));
1117:   PetscCall(VecRestoreArrayRead(faceGeometry, &facegeom));
1118:   PetscCall(ierru);
1119:   PetscFunctionReturn(PETSC_SUCCESS);
1120: }

1122: static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
1123: {
1124:   PetscInt c;
1125:   for (c = 0; c < Nc; ++c) u[c] = 0.0;
1126:   return PETSC_SUCCESS;
1127: }

1129: PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
1130: {
1131:   PetscObject isZero;
1132:   PetscDS     prob;
1133:   PetscInt    numBd, b;

1135:   PetscFunctionBegin;
1136:   PetscCall(DMGetDS(dm, &prob));
1137:   PetscCall(PetscDSGetNumBoundary(prob, &numBd));
1138:   PetscCall(PetscObjectQuery((PetscObject)locX, "__Vec_bc_zero__", &isZero));
1139:   PetscCall(PetscDSUpdateBoundaryLabels(prob, dm));
1140:   for (b = 0; b < numBd; ++b) {
1141:     PetscWeakForm           wf;
1142:     DMBoundaryConditionType type;
1143:     const char             *name;
1144:     DMLabel                 label;
1145:     PetscInt                field, Nc;
1146:     const PetscInt         *comps;
1147:     PetscObject             obj;
1148:     PetscClassId            id;
1149:     void (*bvfunc)(void);
1150:     PetscInt        numids;
1151:     const PetscInt *ids;
1152:     void           *ctx;

1154:     PetscCall(PetscDSGetBoundary(prob, b, &wf, &type, &name, &label, &numids, &ids, &field, &Nc, &comps, &bvfunc, NULL, &ctx));
1155:     if (insertEssential != (type & DM_BC_ESSENTIAL)) continue;
1156:     PetscCall(DMGetField(dm, field, NULL, &obj));
1157:     PetscCall(PetscObjectGetClassId(obj, &id));
1158:     if (id == PETSCFE_CLASSID) {
1159:       switch (type) {
1160:         /* for FEM, there is no insertion to be done for non-essential boundary conditions */
1161:       case DM_BC_ESSENTIAL: {
1162:         PetscSimplePointFn *func = (PetscSimplePointFn *)bvfunc;

1164:         if (isZero) func = zero;
1165:         PetscCall(DMPlexLabelAddCells(dm, label));
1166:         PetscCall(DMPlexInsertBoundaryValuesEssential(dm, time, field, Nc, comps, label, numids, ids, func, ctx, locX));
1167:         PetscCall(DMPlexLabelClearCells(dm, label));
1168:       } break;
1169:       case DM_BC_ESSENTIAL_FIELD: {
1170:         PetscPointFunc func = (PetscPointFunc)bvfunc;

1172:         PetscCall(DMPlexLabelAddCells(dm, label));
1173:         PetscCall(DMPlexInsertBoundaryValuesEssentialField(dm, time, locX, field, Nc, comps, label, numids, ids, func, ctx, locX));
1174:         PetscCall(DMPlexLabelClearCells(dm, label));
1175:       } break;
1176:       default:
1177:         break;
1178:       }
1179:     } else if (id == PETSCFV_CLASSID) {
1180:       {
1181:         PetscErrorCode (*func)(PetscReal, const PetscReal *, const PetscReal *, const PetscScalar *, PetscScalar *, void *) = (PetscErrorCode(*)(PetscReal, const PetscReal *, const PetscReal *, const PetscScalar *, PetscScalar *, void *))bvfunc;

1183:         if (!faceGeomFVM) continue;
1184:         PetscCall(DMPlexInsertBoundaryValuesRiemann(dm, time, faceGeomFVM, cellGeomFVM, gradFVM, field, Nc, comps, label, numids, ids, func, ctx, locX));
1185:       }
1186:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1187:   }
1188:   PetscFunctionReturn(PETSC_SUCCESS);
1189: }

1191: PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues_Plex(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
1192: {
1193:   PetscObject isZero;
1194:   PetscDS     prob;
1195:   PetscInt    numBd, b;

1197:   PetscFunctionBegin;
1198:   if (!locX) PetscFunctionReturn(PETSC_SUCCESS);
1199:   PetscCall(DMGetDS(dm, &prob));
1200:   PetscCall(PetscDSGetNumBoundary(prob, &numBd));
1201:   PetscCall(PetscObjectQuery((PetscObject)locX, "__Vec_bc_zero__", &isZero));
1202:   for (b = 0; b < numBd; ++b) {
1203:     PetscWeakForm           wf;
1204:     DMBoundaryConditionType type;
1205:     const char             *name;
1206:     DMLabel                 label;
1207:     PetscInt                field, Nc;
1208:     const PetscInt         *comps;
1209:     PetscObject             obj;
1210:     PetscClassId            id;
1211:     PetscInt                numids;
1212:     const PetscInt         *ids;
1213:     void (*bvfunc)(void);
1214:     void *ctx;

1216:     PetscCall(PetscDSGetBoundary(prob, b, &wf, &type, &name, &label, &numids, &ids, &field, &Nc, &comps, NULL, &bvfunc, &ctx));
1217:     if (insertEssential != (type & DM_BC_ESSENTIAL)) continue;
1218:     PetscCall(DMGetField(dm, field, NULL, &obj));
1219:     PetscCall(PetscObjectGetClassId(obj, &id));
1220:     if (id == PETSCFE_CLASSID) {
1221:       switch (type) {
1222:         /* for FEM, there is no insertion to be done for non-essential boundary conditions */
1223:       case DM_BC_ESSENTIAL: {
1224:         PetscSimplePointFn *func_t = (PetscSimplePointFn *)bvfunc;

1226:         if (isZero) func_t = zero;
1227:         PetscCall(DMPlexLabelAddCells(dm, label));
1228:         PetscCall(DMPlexInsertBoundaryValuesEssential(dm, time, field, Nc, comps, label, numids, ids, func_t, ctx, locX));
1229:         PetscCall(DMPlexLabelClearCells(dm, label));
1230:       } break;
1231:       case DM_BC_ESSENTIAL_FIELD: {
1232:         PetscPointFunc func_t = (PetscPointFunc)bvfunc;

1234:         PetscCall(DMPlexLabelAddCells(dm, label));
1235:         PetscCall(DMPlexInsertBoundaryValuesEssentialField(dm, time, locX, field, Nc, comps, label, numids, ids, func_t, ctx, locX));
1236:         PetscCall(DMPlexLabelClearCells(dm, label));
1237:       } break;
1238:       default:
1239:         break;
1240:       }
1241:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1242:   }
1243:   PetscFunctionReturn(PETSC_SUCCESS);
1244: }

1246: /*@
1247:   DMPlexInsertBoundaryValues - Puts coefficients which represent boundary values into the local solution vector

1249:   Not Collective

1251:   Input Parameters:
1252: + dm              - The `DM`
1253: . insertEssential - Should I insert essential (e.g. Dirichlet) or inessential (e.g. Neumann) boundary conditions
1254: . time            - The time
1255: . faceGeomFVM     - Face geometry data for FV discretizations
1256: . cellGeomFVM     - Cell geometry data for FV discretizations
1257: - gradFVM         - Gradient reconstruction data for FV discretizations

1259:   Output Parameter:
1260: . locX - Solution updated with boundary values

1262:   Level: intermediate

1264: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunctionLabelLocal()`, `DMAddBoundary()`
1265: @*/
1266: PetscErrorCode DMPlexInsertBoundaryValues(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
1267: {
1268:   PetscFunctionBegin;
1274:   PetscTryMethod(dm, "DMPlexInsertBoundaryValues_C", (DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec), (dm, insertEssential, locX, time, faceGeomFVM, cellGeomFVM, gradFVM));
1275:   PetscFunctionReturn(PETSC_SUCCESS);
1276: }

1278: /*@
1279:   DMPlexInsertTimeDerivativeBoundaryValues - Puts coefficients which represent boundary values of the time derivative into the local solution vector

1281:   Input Parameters:
1282: + dm              - The `DM`
1283: . insertEssential - Should I insert essential (e.g. Dirichlet) or inessential (e.g. Neumann) boundary conditions
1284: . time            - The time
1285: . faceGeomFVM     - Face geometry data for FV discretizations
1286: . cellGeomFVM     - Cell geometry data for FV discretizations
1287: - gradFVM         - Gradient reconstruction data for FV discretizations

1289:   Output Parameter:
1290: . locX_t - Solution updated with boundary values

1292:   Level: developer

1294: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunctionLabelLocal()`
1295: @*/
1296: PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues(DM dm, PetscBool insertEssential, Vec locX_t, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
1297: {
1298:   PetscFunctionBegin;
1304:   PetscTryMethod(dm, "DMPlexInsertTimeDerivativeBoundaryValues_C", (DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec), (dm, insertEssential, locX_t, time, faceGeomFVM, cellGeomFVM, gradFVM));
1305:   PetscFunctionReturn(PETSC_SUCCESS);
1306: }

1308: // Handle non-essential (e.g. outflow) boundary values
1309: PetscErrorCode DMPlexInsertBoundaryValuesFVM(DM dm, PetscFV fv, Vec locX, PetscReal time, Vec *locGradient)
1310: {
1311:   DM  dmGrad;
1312:   Vec cellGeometryFVM, faceGeometryFVM, locGrad = NULL;

1314:   PetscFunctionBegin;
1318:   if (locGradient) {
1319:     PetscAssertPointer(locGradient, 5);
1320:     *locGradient = NULL;
1321:   }
1322:   PetscCall(DMPlexGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL));
1323:   /* Reconstruct and limit cell gradients */
1324:   PetscCall(DMPlexGetGradientDM(dm, fv, &dmGrad));
1325:   if (dmGrad) {
1326:     Vec      grad;
1327:     PetscInt fStart, fEnd;

1329:     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
1330:     PetscCall(DMGetGlobalVector(dmGrad, &grad));
1331:     PetscCall(DMPlexReconstructGradients_Internal(dm, fv, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad));
1332:     /* Communicate gradient values */
1333:     PetscCall(DMGetLocalVector(dmGrad, &locGrad));
1334:     PetscCall(DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad));
1335:     PetscCall(DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad));
1336:     PetscCall(DMRestoreGlobalVector(dmGrad, &grad));
1337:   }
1338:   PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, time, faceGeometryFVM, cellGeometryFVM, locGrad));
1339:   if (locGradient) *locGradient = locGrad;
1340:   else if (locGrad) PetscCall(DMRestoreLocalVector(dmGrad, &locGrad));
1341:   PetscFunctionReturn(PETSC_SUCCESS);
1342: }

1344: PetscErrorCode DMComputeL2Diff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
1345: {
1346:   Vec localX;

1348:   PetscFunctionBegin;
1349:   PetscCall(DMGetLocalVector(dm, &localX));
1350:   PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localX, time, NULL, NULL, NULL));
1351:   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX));
1352:   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX));
1353:   PetscCall(DMPlexComputeL2DiffLocal(dm, time, funcs, ctxs, localX, diff));
1354:   PetscCall(DMRestoreLocalVector(dm, &localX));
1355:   PetscFunctionReturn(PETSC_SUCCESS);
1356: }

1358: /*@C
1359:   DMPlexComputeL2DiffLocal - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.

1361:   Collective

1363:   Input Parameters:
1364: + dm     - The `DM`
1365: . time   - The time
1366: . funcs  - The functions to evaluate for each field component
1367: . ctxs   - Optional array of contexts to pass to each function, or `NULL`.
1368: - localX - The coefficient vector u_h, a local vector

1370:   Output Parameter:
1371: . diff - The diff ||u - u_h||_2

1373:   Level: developer

1375: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunction()`, `DMComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
1376: @*/
1377: PetscErrorCode DMPlexComputeL2DiffLocal(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec localX, PetscReal *diff)
1378: {
1379:   const PetscInt   debug = ((DM_Plex *)dm->data)->printL2;
1380:   DM               tdm;
1381:   Vec              tv;
1382:   PetscSection     section;
1383:   PetscQuadrature  quad;
1384:   PetscFEGeom      fegeom;
1385:   PetscScalar     *funcVal, *interpolant;
1386:   PetscReal       *coords, *gcoords;
1387:   PetscReal        localDiff = 0.0;
1388:   const PetscReal *quadWeights;
1389:   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cellHeight, cStart, cEnd, c, field, fieldOffset;
1390:   PetscBool        transform;

1392:   PetscFunctionBegin;
1393:   PetscCall(DMGetDimension(dm, &dim));
1394:   PetscCall(DMGetCoordinateDim(dm, &coordDim));
1395:   fegeom.dimEmbed = coordDim;
1396:   PetscCall(DMGetLocalSection(dm, &section));
1397:   PetscCall(PetscSectionGetNumFields(section, &numFields));
1398:   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
1399:   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
1400:   PetscCall(DMHasBasisTransform(dm, &transform));
1401:   PetscCheck(numFields, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fields is zero!");
1402:   for (field = 0; field < numFields; ++field) {
1403:     PetscObject  obj;
1404:     PetscClassId id;
1405:     PetscInt     Nc;

1407:     PetscCall(DMGetField(dm, field, NULL, &obj));
1408:     PetscCall(PetscObjectGetClassId(obj, &id));
1409:     if (id == PETSCFE_CLASSID) {
1410:       PetscFE fe = (PetscFE)obj;

1412:       PetscCall(PetscFEGetQuadrature(fe, &quad));
1413:       PetscCall(PetscFEGetNumComponents(fe, &Nc));
1414:     } else if (id == PETSCFV_CLASSID) {
1415:       PetscFV fv = (PetscFV)obj;

1417:       PetscCall(PetscFVGetQuadrature(fv, &quad));
1418:       PetscCall(PetscFVGetNumComponents(fv, &Nc));
1419:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1420:     numComponents += Nc;
1421:   }
1422:   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights));
1423:   PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents);
1424:   PetscCall(PetscMalloc6(numComponents, &funcVal, numComponents, &interpolant, coordDim * (Nq + 1), &coords, Nq, &fegeom.detJ, coordDim * coordDim * Nq, &fegeom.J, coordDim * coordDim * Nq, &fegeom.invJ));
1425:   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
1426:   PetscCall(DMPlexGetSimplexOrBoxCells(dm, cellHeight, &cStart, &cEnd));
1427:   for (c = cStart; c < cEnd; ++c) {
1428:     PetscScalar *x        = NULL;
1429:     PetscReal    elemDiff = 0.0;
1430:     PetscInt     qc       = 0;

1432:     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
1433:     PetscCall(DMPlexVecGetOrientedClosure_Internal(dm, NULL, PETSC_FALSE, localX, c, 0, NULL, &x));

1435:     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1436:       PetscObject  obj;
1437:       PetscClassId id;
1438:       void *const  ctx = ctxs ? ctxs[field] : NULL;
1439:       PetscInt     Nb, Nc, q, fc;

1441:       PetscCall(DMGetField(dm, field, NULL, &obj));
1442:       PetscCall(PetscObjectGetClassId(obj, &id));
1443:       if (id == PETSCFE_CLASSID) {
1444:         PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
1445:         PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
1446:       } else if (id == PETSCFV_CLASSID) {
1447:         PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
1448:         Nb = 1;
1449:       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1450:       if (debug) {
1451:         char title[1024];
1452:         PetscCall(PetscSNPrintf(title, 1023, "Solution for Field %" PetscInt_FMT, field));
1453:         PetscCall(DMPrintCellVector(c, title, Nb, &x[fieldOffset]));
1454:       }
1455:       for (q = 0; q < Nq; ++q) {
1456:         PetscFEGeom    qgeom;
1457:         PetscErrorCode ierr;

1459:         qgeom.dimEmbed = fegeom.dimEmbed;
1460:         qgeom.J        = &fegeom.J[q * coordDim * coordDim];
1461:         qgeom.invJ     = &fegeom.invJ[q * coordDim * coordDim];
1462:         qgeom.detJ     = &fegeom.detJ[q];
1463:         PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT ", point %" PetscInt_FMT, (double)fegeom.detJ[q], c, q);
1464:         if (transform) {
1465:           gcoords = &coords[coordDim * Nq];
1466:           PetscCall(DMPlexBasisTransformApplyReal_Internal(dm, &coords[coordDim * q], PETSC_TRUE, coordDim, &coords[coordDim * q], gcoords, dm->transformCtx));
1467:         } else {
1468:           gcoords = &coords[coordDim * q];
1469:         }
1470:         PetscCall(PetscArrayzero(funcVal, Nc));
1471:         ierr = (*funcs[field])(coordDim, time, gcoords, Nc, funcVal, ctx);
1472:         if (ierr) {
1473:           PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x));
1474:           PetscCall(DMRestoreLocalVector(dm, &localX));
1475:           PetscCall(PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
1476:         }
1477:         if (transform) PetscCall(DMPlexBasisTransformApply_Internal(dm, &coords[coordDim * q], PETSC_FALSE, Nc, funcVal, funcVal, dm->transformCtx));
1478:         if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolate_Static((PetscFE)obj, &x[fieldOffset], &qgeom, q, interpolant));
1479:         else if (id == PETSCFV_CLASSID) PetscCall(PetscFVInterpolate_Static((PetscFV)obj, &x[fieldOffset], q, interpolant));
1480:         else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1481:         for (fc = 0; fc < Nc; ++fc) {
1482:           const PetscReal wt = quadWeights[q * qNc + (qNc == 1 ? 0 : qc + fc)];
1483:           if (debug)
1484:             PetscCall(PetscPrintf(PETSC_COMM_SELF, "    elem %" PetscInt_FMT " field %" PetscInt_FMT ",%" PetscInt_FMT " point %g %g %g diff %g (%g, %g)\n", c, field, fc, (double)(coordDim > 0 ? coords[coordDim * q] : 0.), (double)(coordDim > 1 ? coords[coordDim * q + 1] : 0.), (double)(coordDim > 2 ? coords[coordDim * q + 2] : 0.),
1485:                                   (double)(PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q]), (double)PetscRealPart(interpolant[fc]), (double)PetscRealPart(funcVal[fc])));
1486:           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q];
1487:         }
1488:       }
1489:       fieldOffset += Nb;
1490:       qc += Nc;
1491:     }
1492:     PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x));
1493:     if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  elem %" PetscInt_FMT " diff %g\n", c, (double)elemDiff));
1494:     localDiff += elemDiff;
1495:   }
1496:   PetscCall(PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
1497:   PetscCall(MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm)));
1498:   *diff = PetscSqrtReal(*diff);
1499:   PetscFunctionReturn(PETSC_SUCCESS);
1500: }

1502: PetscErrorCode DMComputeL2GradientDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
1503: {
1504:   const PetscInt   debug = ((DM_Plex *)dm->data)->printL2;
1505:   DM               tdm;
1506:   PetscSection     section;
1507:   PetscQuadrature  quad;
1508:   Vec              localX, tv;
1509:   PetscScalar     *funcVal, *interpolant;
1510:   const PetscReal *quadWeights;
1511:   PetscFEGeom      fegeom;
1512:   PetscReal       *coords, *gcoords;
1513:   PetscReal        localDiff = 0.0;
1514:   PetscInt         dim, coordDim, qNc = 0, Nq = 0, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset;
1515:   PetscBool        transform;

1517:   PetscFunctionBegin;
1518:   PetscCall(DMGetDimension(dm, &dim));
1519:   PetscCall(DMGetCoordinateDim(dm, &coordDim));
1520:   fegeom.dimEmbed = coordDim;
1521:   PetscCall(DMGetLocalSection(dm, &section));
1522:   PetscCall(PetscSectionGetNumFields(section, &numFields));
1523:   PetscCall(DMGetLocalVector(dm, &localX));
1524:   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX));
1525:   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX));
1526:   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
1527:   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
1528:   PetscCall(DMHasBasisTransform(dm, &transform));
1529:   for (field = 0; field < numFields; ++field) {
1530:     PetscFE  fe;
1531:     PetscInt Nc;

1533:     PetscCall(DMGetField(dm, field, NULL, (PetscObject *)&fe));
1534:     PetscCall(PetscFEGetQuadrature(fe, &quad));
1535:     PetscCall(PetscFEGetNumComponents(fe, &Nc));
1536:     numComponents += Nc;
1537:   }
1538:   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights));
1539:   PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents);
1540:   /* PetscCall(DMProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX)); */
1541:   PetscCall(PetscMalloc6(numComponents, &funcVal, coordDim * (Nq + 1), &coords, coordDim * coordDim * Nq, &fegeom.J, coordDim * coordDim * Nq, &fegeom.invJ, numComponents * coordDim, &interpolant, Nq, &fegeom.detJ));
1542:   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
1543:   for (c = cStart; c < cEnd; ++c) {
1544:     PetscScalar *x        = NULL;
1545:     PetscReal    elemDiff = 0.0;
1546:     PetscInt     qc       = 0;

1548:     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
1549:     PetscCall(DMPlexVecGetOrientedClosure_Internal(dm, NULL, PETSC_FALSE, localX, c, 0, NULL, &x));

1551:     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1552:       PetscFE     fe;
1553:       void *const ctx = ctxs ? ctxs[field] : NULL;
1554:       PetscInt    Nb, Nc, q, fc;

1556:       PetscCall(DMGetField(dm, field, NULL, (PetscObject *)&fe));
1557:       PetscCall(PetscFEGetDimension(fe, &Nb));
1558:       PetscCall(PetscFEGetNumComponents(fe, &Nc));
1559:       if (debug) {
1560:         char title[1024];
1561:         PetscCall(PetscSNPrintf(title, 1023, "Solution for Field %" PetscInt_FMT, field));
1562:         PetscCall(DMPrintCellVector(c, title, Nb, &x[fieldOffset]));
1563:       }
1564:       for (q = 0; q < Nq; ++q) {
1565:         PetscFEGeom    qgeom;
1566:         PetscErrorCode ierr;

1568:         qgeom.dimEmbed = fegeom.dimEmbed;
1569:         qgeom.J        = &fegeom.J[q * coordDim * coordDim];
1570:         qgeom.invJ     = &fegeom.invJ[q * coordDim * coordDim];
1571:         qgeom.detJ     = &fegeom.detJ[q];
1572:         PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT ", quadrature points %" PetscInt_FMT, (double)fegeom.detJ[q], c, q);
1573:         if (transform) {
1574:           gcoords = &coords[coordDim * Nq];
1575:           PetscCall(DMPlexBasisTransformApplyReal_Internal(dm, &coords[coordDim * q], PETSC_TRUE, coordDim, &coords[coordDim * q], gcoords, dm->transformCtx));
1576:         } else {
1577:           gcoords = &coords[coordDim * q];
1578:         }
1579:         PetscCall(PetscArrayzero(funcVal, Nc));
1580:         ierr = (*funcs[field])(coordDim, time, gcoords, n, Nc, funcVal, ctx);
1581:         if (ierr) {
1582:           PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x));
1583:           PetscCall(DMRestoreLocalVector(dm, &localX));
1584:           PetscCall(PetscFree6(funcVal, coords, fegeom.J, fegeom.invJ, interpolant, fegeom.detJ));
1585:         }
1586:         if (transform) PetscCall(DMPlexBasisTransformApply_Internal(dm, &coords[coordDim * q], PETSC_FALSE, Nc, funcVal, funcVal, dm->transformCtx));
1587:         PetscCall(PetscFEInterpolateGradient_Static(fe, 1, &x[fieldOffset], &qgeom, q, interpolant));
1588:         /* Overwrite with the dot product if the normal is given */
1589:         if (n) {
1590:           for (fc = 0; fc < Nc; ++fc) {
1591:             PetscScalar sum = 0.0;
1592:             PetscInt    d;
1593:             for (d = 0; d < dim; ++d) sum += interpolant[fc * dim + d] * n[d];
1594:             interpolant[fc] = sum;
1595:           }
1596:         }
1597:         for (fc = 0; fc < Nc; ++fc) {
1598:           const PetscReal wt = quadWeights[q * qNc + (qNc == 1 ? 0 : qc + fc)];
1599:           if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    elem %" PetscInt_FMT " fieldDer %" PetscInt_FMT ",%" PetscInt_FMT " diff %g\n", c, field, fc, (double)(PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q])));
1600:           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q];
1601:         }
1602:       }
1603:       fieldOffset += Nb;
1604:       qc += Nc;
1605:     }
1606:     PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x));
1607:     if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  elem %" PetscInt_FMT " diff %g\n", c, (double)elemDiff));
1608:     localDiff += elemDiff;
1609:   }
1610:   PetscCall(PetscFree6(funcVal, coords, fegeom.J, fegeom.invJ, interpolant, fegeom.detJ));
1611:   PetscCall(DMRestoreLocalVector(dm, &localX));
1612:   PetscCall(MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm)));
1613:   *diff = PetscSqrtReal(*diff);
1614:   PetscFunctionReturn(PETSC_SUCCESS);
1615: }

1617: PetscErrorCode DMComputeL2FieldDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
1618: {
1619:   const PetscInt debug = ((DM_Plex *)dm->data)->printL2;
1620:   DM             tdm;
1621:   DMLabel        depthLabel;
1622:   PetscSection   section;
1623:   Vec            localX, tv;
1624:   PetscReal     *localDiff;
1625:   PetscInt       dim, depth, dE, Nf, f, Nds, s;
1626:   PetscBool      transform;

1628:   PetscFunctionBegin;
1629:   PetscCall(DMGetDimension(dm, &dim));
1630:   PetscCall(DMGetCoordinateDim(dm, &dE));
1631:   PetscCall(DMGetLocalSection(dm, &section));
1632:   PetscCall(DMGetLocalVector(dm, &localX));
1633:   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
1634:   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
1635:   PetscCall(DMHasBasisTransform(dm, &transform));
1636:   PetscCall(DMGetNumFields(dm, &Nf));
1637:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
1638:   PetscCall(DMLabelGetNumValues(depthLabel, &depth));

1640:   PetscCall(VecSet(localX, 0.0));
1641:   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX));
1642:   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX));
1643:   PetscCall(DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX));
1644:   PetscCall(DMGetNumDS(dm, &Nds));
1645:   PetscCall(PetscCalloc1(Nf, &localDiff));
1646:   for (s = 0; s < Nds; ++s) {
1647:     PetscDS          ds;
1648:     DMLabel          label;
1649:     IS               fieldIS, pointIS;
1650:     const PetscInt  *fields, *points = NULL;
1651:     PetscQuadrature  quad;
1652:     const PetscReal *quadPoints, *quadWeights;
1653:     PetscFEGeom      fegeom;
1654:     PetscReal       *coords, *gcoords;
1655:     PetscScalar     *funcVal, *interpolant;
1656:     PetscBool        isCohesive;
1657:     PetscInt         qNc, Nq, totNc, cStart = 0, cEnd, c, dsNf;

1659:     PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
1660:     PetscCall(ISGetIndices(fieldIS, &fields));
1661:     PetscCall(PetscDSIsCohesive(ds, &isCohesive));
1662:     PetscCall(PetscDSGetNumFields(ds, &dsNf));
1663:     PetscCall(PetscDSGetTotalComponents(ds, &totNc));
1664:     PetscCall(PetscDSGetQuadrature(ds, &quad));
1665:     PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
1666:     PetscCheck(!(qNc != 1) || !(qNc != totNc), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, totNc);
1667:     PetscCall(PetscCalloc6(totNc, &funcVal, totNc, &interpolant, dE * (Nq + 1), &coords, Nq, &fegeom.detJ, dE * dE * Nq, &fegeom.J, dE * dE * Nq, &fegeom.invJ));
1668:     if (!label) {
1669:       PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
1670:     } else {
1671:       PetscCall(DMLabelGetStratumIS(label, 1, &pointIS));
1672:       PetscCall(ISGetLocalSize(pointIS, &cEnd));
1673:       PetscCall(ISGetIndices(pointIS, &points));
1674:     }
1675:     for (c = cStart; c < cEnd; ++c) {
1676:       const PetscInt  cell = points ? points[c] : c;
1677:       PetscScalar    *x    = NULL;
1678:       const PetscInt *cone;
1679:       PetscInt        qc = 0, fOff = 0, dep;

1681:       PetscCall(DMLabelGetValue(depthLabel, cell, &dep));
1682:       if (dep != depth - 1) continue;
1683:       if (isCohesive) {
1684:         PetscCall(DMPlexGetCone(dm, cell, &cone));
1685:         PetscCall(DMPlexComputeCellGeometryFEM(dm, cone[0], quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
1686:       } else {
1687:         PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
1688:       }
1689:       PetscCall(DMPlexVecGetOrientedClosure_Internal(dm, NULL, PETSC_FALSE, localX, cell, 0, NULL, &x));
1690:       for (f = 0; f < dsNf; ++f) {
1691:         PetscObject  obj;
1692:         PetscClassId id;
1693:         void *const  ctx = ctxs ? ctxs[fields[f]] : NULL;
1694:         PetscInt     Nb, Nc, q, fc;
1695:         PetscReal    elemDiff = 0.0;
1696:         PetscBool    cohesive;

1698:         PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
1699:         PetscCall(PetscDSGetDiscretization(ds, f, &obj));
1700:         PetscCall(PetscObjectGetClassId(obj, &id));
1701:         if (id == PETSCFE_CLASSID) {
1702:           PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
1703:           PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
1704:         } else if (id == PETSCFV_CLASSID) {
1705:           PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
1706:           Nb = 1;
1707:         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, fields[f]);
1708:         if (isCohesive && !cohesive) {
1709:           fOff += Nb * 2;
1710:           qc += Nc;
1711:           continue;
1712:         }
1713:         if (debug) {
1714:           char title[1024];
1715:           PetscCall(PetscSNPrintf(title, 1023, "Solution for Field %" PetscInt_FMT, fields[f]));
1716:           PetscCall(DMPrintCellVector(cell, title, Nb, &x[fOff]));
1717:         }
1718:         for (q = 0; q < Nq; ++q) {
1719:           PetscFEGeom    qgeom;
1720:           PetscErrorCode ierr;

1722:           qgeom.dimEmbed = fegeom.dimEmbed;
1723:           qgeom.J        = &fegeom.J[q * dE * dE];
1724:           qgeom.invJ     = &fegeom.invJ[q * dE * dE];
1725:           qgeom.detJ     = &fegeom.detJ[q];
1726:           PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for cell %" PetscInt_FMT ", quadrature point %" PetscInt_FMT, (double)fegeom.detJ[q], cell, q);
1727:           if (transform) {
1728:             gcoords = &coords[dE * Nq];
1729:             PetscCall(DMPlexBasisTransformApplyReal_Internal(dm, &coords[dE * q], PETSC_TRUE, dE, &coords[dE * q], gcoords, dm->transformCtx));
1730:           } else {
1731:             gcoords = &coords[dE * q];
1732:           }
1733:           for (fc = 0; fc < Nc; ++fc) funcVal[fc] = 0.;
1734:           ierr = (*funcs[fields[f]])(dE, time, gcoords, Nc, funcVal, ctx);
1735:           if (ierr) {
1736:             PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, cell, NULL, &x));
1737:             PetscCall(DMRestoreLocalVector(dm, &localX));
1738:             PetscCall(PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
1739:           }
1740:           if (transform) PetscCall(DMPlexBasisTransformApply_Internal(dm, &coords[dE * q], PETSC_FALSE, Nc, funcVal, funcVal, dm->transformCtx));
1741:           /* Call once for each face, except for lagrange field */
1742:           if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolate_Static((PetscFE)obj, &x[fOff], &qgeom, q, interpolant));
1743:           else if (id == PETSCFV_CLASSID) PetscCall(PetscFVInterpolate_Static((PetscFV)obj, &x[fOff], q, interpolant));
1744:           else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, fields[f]);
1745:           for (fc = 0; fc < Nc; ++fc) {
1746:             const PetscReal wt = quadWeights[q * qNc + (qNc == 1 ? 0 : qc + fc)];
1747:             if (debug)
1748:               PetscCall(PetscPrintf(PETSC_COMM_SELF, "    cell %" PetscInt_FMT " field %" PetscInt_FMT ",%" PetscInt_FMT " point %g %g %g diff %g\n", cell, fields[f], fc, (double)(dE > 0 ? coords[dE * q] : 0.), (double)(dE > 1 ? coords[dE * q + 1] : 0.), (double)(dE > 2 ? coords[dE * q + 2] : 0.),
1749:                                     (double)(PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q])));
1750:             elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q];
1751:           }
1752:         }
1753:         fOff += Nb;
1754:         qc += Nc;
1755:         localDiff[fields[f]] += elemDiff;
1756:         if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  cell %" PetscInt_FMT " field %" PetscInt_FMT " cum diff %g\n", cell, fields[f], (double)localDiff[fields[f]]));
1757:       }
1758:       PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, cell, NULL, &x));
1759:     }
1760:     if (label) {
1761:       PetscCall(ISRestoreIndices(pointIS, &points));
1762:       PetscCall(ISDestroy(&pointIS));
1763:     }
1764:     PetscCall(ISRestoreIndices(fieldIS, &fields));
1765:     PetscCall(PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
1766:   }
1767:   PetscCall(DMRestoreLocalVector(dm, &localX));
1768:   PetscCall(MPIU_Allreduce(localDiff, diff, Nf, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm)));
1769:   PetscCall(PetscFree(localDiff));
1770:   for (f = 0; f < Nf; ++f) diff[f] = PetscSqrtReal(diff[f]);
1771:   PetscFunctionReturn(PETSC_SUCCESS);
1772: }

1774: /*@C
1775:   DMPlexComputeL2DiffVec - This function computes the cellwise L_2 difference between a function u and an FEM interpolant solution u_h, and stores it in a Vec.

1777:   Collective

1779:   Input Parameters:
1780: + dm    - The `DM`
1781: . time  - The time
1782: . funcs - The functions to evaluate for each field component: `NULL` means that component does not contribute to error calculation
1783: . ctxs  - Optional array of contexts to pass to each function, or `NULL`.
1784: - X     - The coefficient vector u_h

1786:   Output Parameter:
1787: . D - A `Vec` which holds the difference ||u - u_h||_2 for each cell

1789:   Level: developer

1791: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMPlexComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
1792: @*/
1793: PetscErrorCode DMPlexComputeL2DiffVec(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, Vec D)
1794: {
1795:   PetscSection     section;
1796:   PetscQuadrature  quad;
1797:   Vec              localX;
1798:   PetscFEGeom      fegeom;
1799:   PetscScalar     *funcVal, *interpolant;
1800:   PetscReal       *coords;
1801:   const PetscReal *quadPoints, *quadWeights;
1802:   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, c, field, fieldOffset;

1804:   PetscFunctionBegin;
1805:   PetscCall(VecSet(D, 0.0));
1806:   PetscCall(DMGetDimension(dm, &dim));
1807:   PetscCall(DMGetCoordinateDim(dm, &coordDim));
1808:   PetscCall(DMGetLocalSection(dm, &section));
1809:   PetscCall(PetscSectionGetNumFields(section, &numFields));
1810:   PetscCall(DMGetLocalVector(dm, &localX));
1811:   PetscCall(DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX));
1812:   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX));
1813:   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX));
1814:   for (field = 0; field < numFields; ++field) {
1815:     PetscObject  obj;
1816:     PetscClassId id;
1817:     PetscInt     Nc;

1819:     PetscCall(DMGetField(dm, field, NULL, &obj));
1820:     PetscCall(PetscObjectGetClassId(obj, &id));
1821:     if (id == PETSCFE_CLASSID) {
1822:       PetscFE fe = (PetscFE)obj;

1824:       PetscCall(PetscFEGetQuadrature(fe, &quad));
1825:       PetscCall(PetscFEGetNumComponents(fe, &Nc));
1826:     } else if (id == PETSCFV_CLASSID) {
1827:       PetscFV fv = (PetscFV)obj;

1829:       PetscCall(PetscFVGetQuadrature(fv, &quad));
1830:       PetscCall(PetscFVGetNumComponents(fv, &Nc));
1831:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1832:     numComponents += Nc;
1833:   }
1834:   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
1835:   PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents);
1836:   PetscCall(PetscMalloc6(numComponents, &funcVal, numComponents, &interpolant, coordDim * Nq, &coords, Nq, &fegeom.detJ, coordDim * coordDim * Nq, &fegeom.J, coordDim * coordDim * Nq, &fegeom.invJ));
1837:   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
1838:   for (c = cStart; c < cEnd; ++c) {
1839:     PetscScalar *x        = NULL;
1840:     PetscScalar  elemDiff = 0.0;
1841:     PetscInt     qc       = 0;

1843:     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
1844:     PetscCall(DMPlexVecGetOrientedClosure_Internal(dm, NULL, PETSC_FALSE, localX, c, 0, NULL, &x));

1846:     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1847:       PetscObject  obj;
1848:       PetscClassId id;
1849:       void *const  ctx = ctxs ? ctxs[field] : NULL;
1850:       PetscInt     Nb, Nc, q, fc;

1852:       PetscCall(DMGetField(dm, field, NULL, &obj));
1853:       PetscCall(PetscObjectGetClassId(obj, &id));
1854:       if (id == PETSCFE_CLASSID) {
1855:         PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
1856:         PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
1857:       } else if (id == PETSCFV_CLASSID) {
1858:         PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
1859:         Nb = 1;
1860:       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1861:       if (funcs[field]) {
1862:         for (q = 0; q < Nq; ++q) {
1863:           PetscFEGeom qgeom;

1865:           qgeom.dimEmbed = fegeom.dimEmbed;
1866:           qgeom.J        = &fegeom.J[q * coordDim * coordDim];
1867:           qgeom.invJ     = &fegeom.invJ[q * coordDim * coordDim];
1868:           qgeom.detJ     = &fegeom.detJ[q];
1869:           PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT ", quadrature points %" PetscInt_FMT, (double)fegeom.detJ[q], c, q);
1870:           PetscCall((*funcs[field])(coordDim, time, &coords[q * coordDim], Nc, funcVal, ctx));
1871: #if defined(needs_fix_with_return_code_argument)
1872:           if (ierr) {
1873:             PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x));
1874:             PetscCall(PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
1875:             PetscCall(DMRestoreLocalVector(dm, &localX));
1876:           }
1877: #endif
1878:           if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolate_Static((PetscFE)obj, &x[fieldOffset], &qgeom, q, interpolant));
1879:           else if (id == PETSCFV_CLASSID) PetscCall(PetscFVInterpolate_Static((PetscFV)obj, &x[fieldOffset], q, interpolant));
1880:           else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1881:           for (fc = 0; fc < Nc; ++fc) {
1882:             const PetscReal wt = quadWeights[q * qNc + (qNc == 1 ? 0 : qc + fc)];
1883:             elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q];
1884:           }
1885:         }
1886:       }
1887:       fieldOffset += Nb;
1888:       qc += Nc;
1889:     }
1890:     PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x));
1891:     PetscCall(VecSetValue(D, c - cStart, elemDiff, INSERT_VALUES));
1892:   }
1893:   PetscCall(PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
1894:   PetscCall(DMRestoreLocalVector(dm, &localX));
1895:   PetscCall(VecSqrtAbs(D));
1896:   PetscFunctionReturn(PETSC_SUCCESS);
1897: }

1899: /*@
1900:   DMPlexComputeClementInterpolant - This function computes the L2 projection of the cellwise values of a function u onto P1

1902:   Collective

1904:   Input Parameters:
1905: + dm   - The `DM`
1906: - locX - The coefficient vector u_h

1908:   Output Parameter:
1909: . locC - A `Vec` which holds the Clement interpolant of the function

1911:   Level: developer

1913:   Note:
1914:   $ u_h(v_i) = \sum_{T_i \in support(v_i)} |T_i| u_h(T_i) / \sum_{T_i \in support(v_i)} |T_i| $ where $ |T_i| $ is the cell volume

1916: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMPlexComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
1917: @*/
1918: PetscErrorCode DMPlexComputeClementInterpolant(DM dm, Vec locX, Vec locC)
1919: {
1920:   PetscInt         debug = ((DM_Plex *)dm->data)->printFEM;
1921:   DM               dmc;
1922:   PetscQuadrature  quad;
1923:   PetscScalar     *interpolant, *valsum;
1924:   PetscFEGeom      fegeom;
1925:   PetscReal       *coords;
1926:   const PetscReal *quadPoints, *quadWeights;
1927:   PetscInt         dim, cdim, Nf, f, Nc = 0, Nq, qNc, cStart, cEnd, vStart, vEnd, v;

1929:   PetscFunctionBegin;
1930:   PetscCall(PetscCitationsRegister(ClementCitation, &Clementcite));
1931:   PetscCall(VecGetDM(locC, &dmc));
1932:   PetscCall(VecSet(locC, 0.0));
1933:   PetscCall(DMGetDimension(dm, &dim));
1934:   PetscCall(DMGetCoordinateDim(dm, &cdim));
1935:   fegeom.dimEmbed = cdim;
1936:   PetscCall(DMGetNumFields(dm, &Nf));
1937:   PetscCheck(Nf > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fields is zero!");
1938:   for (f = 0; f < Nf; ++f) {
1939:     PetscObject  obj;
1940:     PetscClassId id;
1941:     PetscInt     fNc;

1943:     PetscCall(DMGetField(dm, f, NULL, &obj));
1944:     PetscCall(PetscObjectGetClassId(obj, &id));
1945:     if (id == PETSCFE_CLASSID) {
1946:       PetscFE fe = (PetscFE)obj;

1948:       PetscCall(PetscFEGetQuadrature(fe, &quad));
1949:       PetscCall(PetscFEGetNumComponents(fe, &fNc));
1950:     } else if (id == PETSCFV_CLASSID) {
1951:       PetscFV fv = (PetscFV)obj;

1953:       PetscCall(PetscFVGetQuadrature(fv, &quad));
1954:       PetscCall(PetscFVGetNumComponents(fv, &fNc));
1955:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
1956:     Nc += fNc;
1957:   }
1958:   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
1959:   PetscCheck(qNc == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " > 1", qNc);
1960:   PetscCall(PetscMalloc6(Nc * 2, &valsum, Nc, &interpolant, cdim * Nq, &coords, Nq, &fegeom.detJ, cdim * cdim * Nq, &fegeom.J, cdim * cdim * Nq, &fegeom.invJ));
1961:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1962:   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
1963:   for (v = vStart; v < vEnd; ++v) {
1964:     PetscScalar volsum = 0.0;
1965:     PetscInt   *star   = NULL;
1966:     PetscInt    starSize, st, fc;

1968:     PetscCall(PetscArrayzero(valsum, Nc));
1969:     PetscCall(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
1970:     for (st = 0; st < starSize * 2; st += 2) {
1971:       const PetscInt cell = star[st];
1972:       PetscScalar   *val  = &valsum[Nc];
1973:       PetscScalar   *x    = NULL;
1974:       PetscReal      vol  = 0.0;
1975:       PetscInt       foff = 0;

1977:       if ((cell < cStart) || (cell >= cEnd)) continue;
1978:       PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
1979:       PetscCall(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x));
1980:       for (f = 0; f < Nf; ++f) {
1981:         PetscObject  obj;
1982:         PetscClassId id;
1983:         PetscInt     Nb, fNc, q;

1985:         PetscCall(PetscArrayzero(val, Nc));
1986:         PetscCall(DMGetField(dm, f, NULL, &obj));
1987:         PetscCall(PetscObjectGetClassId(obj, &id));
1988:         if (id == PETSCFE_CLASSID) {
1989:           PetscCall(PetscFEGetNumComponents((PetscFE)obj, &fNc));
1990:           PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
1991:         } else if (id == PETSCFV_CLASSID) {
1992:           PetscCall(PetscFVGetNumComponents((PetscFV)obj, &fNc));
1993:           Nb = 1;
1994:         } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
1995:         for (q = 0; q < Nq; ++q) {
1996:           const PetscReal wt = quadWeights[q] * fegeom.detJ[q];
1997:           PetscFEGeom     qgeom;

1999:           qgeom.dimEmbed = fegeom.dimEmbed;
2000:           qgeom.J        = &fegeom.J[q * cdim * cdim];
2001:           qgeom.invJ     = &fegeom.invJ[q * cdim * cdim];
2002:           qgeom.detJ     = &fegeom.detJ[q];
2003:           PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT ", quadrature points %" PetscInt_FMT, (double)fegeom.detJ[q], cell, q);
2004:           if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolate_Static((PetscFE)obj, &x[foff], &qgeom, q, interpolant));
2005:           else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
2006:           for (fc = 0; fc < fNc; ++fc) val[foff + fc] += interpolant[fc] * wt;
2007:           vol += wt;
2008:         }
2009:         foff += Nb;
2010:       }
2011:       PetscCall(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x));
2012:       for (fc = 0; fc < Nc; ++fc) valsum[fc] += val[fc];
2013:       volsum += vol;
2014:       if (debug) {
2015:         PetscCall(PetscPrintf(PETSC_COMM_SELF, "Vertex %" PetscInt_FMT " Cell %" PetscInt_FMT " value: [", v, cell));
2016:         for (fc = 0; fc < Nc; ++fc) {
2017:           if (fc) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
2018:           PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(val[fc])));
2019:         }
2020:         PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
2021:       }
2022:     }
2023:     for (fc = 0; fc < Nc; ++fc) valsum[fc] /= volsum;
2024:     PetscCall(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
2025:     PetscCall(DMPlexVecSetClosure(dmc, NULL, locC, v, valsum, INSERT_VALUES));
2026:   }
2027:   PetscCall(PetscFree6(valsum, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
2028:   PetscFunctionReturn(PETSC_SUCCESS);
2029: }

2031: /*@
2032:   DMPlexComputeGradientClementInterpolant - This function computes the L2 projection of the cellwise gradient of a function u onto P1

2034:   Collective

2036:   Input Parameters:
2037: + dm   - The `DM`
2038: - locX - The coefficient vector u_h

2040:   Output Parameter:
2041: . locC - A `Vec` which holds the Clement interpolant of the gradient

2043:   Level: developer

2045:   Note:
2046:   $\nabla u_h(v_i) = \sum_{T_i \in support(v_i)} |T_i| \nabla u_h(T_i) / \sum_{T_i \in support(v_i)} |T_i| $ where $ |T_i| $ is the cell volume

2048: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMPlexComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
2049: @*/
2050: PetscErrorCode DMPlexComputeGradientClementInterpolant(DM dm, Vec locX, Vec locC)
2051: {
2052:   DM_Plex         *mesh  = (DM_Plex *)dm->data;
2053:   PetscInt         debug = mesh->printFEM;
2054:   DM               dmC;
2055:   PetscQuadrature  quad;
2056:   PetscScalar     *interpolant, *gradsum;
2057:   PetscFEGeom      fegeom;
2058:   PetscReal       *coords;
2059:   const PetscReal *quadPoints, *quadWeights;
2060:   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset;

2062:   PetscFunctionBegin;
2063:   PetscCall(PetscCitationsRegister(ClementCitation, &Clementcite));
2064:   PetscCall(VecGetDM(locC, &dmC));
2065:   PetscCall(VecSet(locC, 0.0));
2066:   PetscCall(DMGetDimension(dm, &dim));
2067:   PetscCall(DMGetCoordinateDim(dm, &coordDim));
2068:   fegeom.dimEmbed = coordDim;
2069:   PetscCall(DMGetNumFields(dm, &numFields));
2070:   PetscCheck(numFields, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fields is zero!");
2071:   for (field = 0; field < numFields; ++field) {
2072:     PetscObject  obj;
2073:     PetscClassId id;
2074:     PetscInt     Nc;

2076:     PetscCall(DMGetField(dm, field, NULL, &obj));
2077:     PetscCall(PetscObjectGetClassId(obj, &id));
2078:     if (id == PETSCFE_CLASSID) {
2079:       PetscFE fe = (PetscFE)obj;

2081:       PetscCall(PetscFEGetQuadrature(fe, &quad));
2082:       PetscCall(PetscFEGetNumComponents(fe, &Nc));
2083:     } else if (id == PETSCFV_CLASSID) {
2084:       PetscFV fv = (PetscFV)obj;

2086:       PetscCall(PetscFVGetQuadrature(fv, &quad));
2087:       PetscCall(PetscFVGetNumComponents(fv, &Nc));
2088:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
2089:     numComponents += Nc;
2090:   }
2091:   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
2092:   PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents);
2093:   PetscCall(PetscMalloc6(coordDim * numComponents * 2, &gradsum, coordDim * numComponents, &interpolant, coordDim * Nq, &coords, Nq, &fegeom.detJ, coordDim * coordDim * Nq, &fegeom.J, coordDim * coordDim * Nq, &fegeom.invJ));
2094:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2095:   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
2096:   for (v = vStart; v < vEnd; ++v) {
2097:     PetscScalar volsum = 0.0;
2098:     PetscInt   *star   = NULL;
2099:     PetscInt    starSize, st, d, fc;

2101:     PetscCall(PetscArrayzero(gradsum, coordDim * numComponents));
2102:     PetscCall(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
2103:     for (st = 0; st < starSize * 2; st += 2) {
2104:       const PetscInt cell = star[st];
2105:       PetscScalar   *grad = &gradsum[coordDim * numComponents];
2106:       PetscScalar   *x    = NULL;
2107:       PetscReal      vol  = 0.0;

2109:       if ((cell < cStart) || (cell >= cEnd)) continue;
2110:       PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
2111:       PetscCall(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x));
2112:       for (field = 0, fieldOffset = 0; field < numFields; ++field) {
2113:         PetscObject  obj;
2114:         PetscClassId id;
2115:         PetscInt     Nb, Nc, q, qc = 0;

2117:         PetscCall(PetscArrayzero(grad, coordDim * numComponents));
2118:         PetscCall(DMGetField(dm, field, NULL, &obj));
2119:         PetscCall(PetscObjectGetClassId(obj, &id));
2120:         if (id == PETSCFE_CLASSID) {
2121:           PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
2122:           PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
2123:         } else if (id == PETSCFV_CLASSID) {
2124:           PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
2125:           Nb = 1;
2126:         } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
2127:         for (q = 0; q < Nq; ++q) {
2128:           PetscFEGeom qgeom;

2130:           qgeom.dimEmbed = fegeom.dimEmbed;
2131:           qgeom.J        = &fegeom.J[q * coordDim * coordDim];
2132:           qgeom.invJ     = &fegeom.invJ[q * coordDim * coordDim];
2133:           qgeom.detJ     = &fegeom.detJ[q];
2134:           PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT ", quadrature points %" PetscInt_FMT, (double)fegeom.detJ[q], cell, q);
2135:           if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolateGradient_Static((PetscFE)obj, 1, &x[fieldOffset], &qgeom, q, interpolant));
2136:           else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
2137:           for (fc = 0; fc < Nc; ++fc) {
2138:             const PetscReal wt = quadWeights[q * qNc + qc];

2140:             for (d = 0; d < coordDim; ++d) grad[fc * coordDim + d] += interpolant[fc * dim + d] * wt * fegeom.detJ[q];
2141:           }
2142:           vol += quadWeights[q * qNc] * fegeom.detJ[q];
2143:         }
2144:         fieldOffset += Nb;
2145:         qc += Nc;
2146:       }
2147:       PetscCall(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x));
2148:       for (fc = 0; fc < numComponents; ++fc) {
2149:         for (d = 0; d < coordDim; ++d) gradsum[fc * coordDim + d] += grad[fc * coordDim + d];
2150:       }
2151:       volsum += vol;
2152:       if (debug) {
2153:         PetscCall(PetscPrintf(PETSC_COMM_SELF, "Vertex %" PetscInt_FMT " Cell %" PetscInt_FMT " gradient: [", v, cell));
2154:         for (fc = 0; fc < numComponents; ++fc) {
2155:           for (d = 0; d < coordDim; ++d) {
2156:             if (fc || d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
2157:             PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc * coordDim + d])));
2158:           }
2159:         }
2160:         PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
2161:       }
2162:     }
2163:     for (fc = 0; fc < numComponents; ++fc) {
2164:       for (d = 0; d < coordDim; ++d) gradsum[fc * coordDim + d] /= volsum;
2165:     }
2166:     PetscCall(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
2167:     PetscCall(DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES));
2168:   }
2169:   PetscCall(PetscFree6(gradsum, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
2170:   PetscFunctionReturn(PETSC_SUCCESS);
2171: }

2173: PetscErrorCode DMPlexComputeIntegral_Internal(DM dm, Vec locX, PetscInt cStart, PetscInt cEnd, PetscScalar *cintegral, void *user)
2174: {
2175:   DM           dmAux = NULL, plexA = NULL;
2176:   PetscDS      prob, probAux       = NULL;
2177:   PetscSection section, sectionAux;
2178:   Vec          locA;
2179:   PetscInt     dim, numCells = cEnd - cStart, c, f;
2180:   PetscBool    useFVM = PETSC_FALSE;
2181:   /* DS */
2182:   PetscInt           Nf, totDim, *uOff, *uOff_x, numConstants;
2183:   PetscInt           NfAux, totDimAux, *aOff;
2184:   PetscScalar       *u, *a = NULL;
2185:   const PetscScalar *constants;
2186:   /* Geometry */
2187:   PetscFEGeom       *cgeomFEM;
2188:   DM                 dmGrad;
2189:   PetscQuadrature    affineQuad      = NULL;
2190:   Vec                cellGeometryFVM = NULL, faceGeometryFVM = NULL, locGrad = NULL;
2191:   PetscFVCellGeom   *cgeomFVM;
2192:   const PetscScalar *lgrad;
2193:   PetscInt           maxDegree;
2194:   DMField            coordField;
2195:   IS                 cellIS;

2197:   PetscFunctionBegin;
2198:   PetscCall(DMGetDS(dm, &prob));
2199:   PetscCall(DMGetDimension(dm, &dim));
2200:   PetscCall(DMGetLocalSection(dm, &section));
2201:   PetscCall(DMGetNumFields(dm, &Nf));
2202:   /* Determine which discretizations we have */
2203:   for (f = 0; f < Nf; ++f) {
2204:     PetscObject  obj;
2205:     PetscClassId id;

2207:     PetscCall(PetscDSGetDiscretization(prob, f, &obj));
2208:     PetscCall(PetscObjectGetClassId(obj, &id));
2209:     if (id == PETSCFV_CLASSID) useFVM = PETSC_TRUE;
2210:   }
2211:   /* Read DS information */
2212:   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
2213:   PetscCall(PetscDSGetComponentOffsets(prob, &uOff));
2214:   PetscCall(PetscDSGetComponentDerivativeOffsets(prob, &uOff_x));
2215:   PetscCall(ISCreateStride(PETSC_COMM_SELF, numCells, cStart, 1, &cellIS));
2216:   PetscCall(PetscDSGetConstants(prob, &numConstants, &constants));
2217:   /* Read Auxiliary DS information */
2218:   PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &locA));
2219:   if (locA) {
2220:     PetscCall(VecGetDM(locA, &dmAux));
2221:     PetscCall(DMConvert(dmAux, DMPLEX, &plexA));
2222:     PetscCall(DMGetDS(dmAux, &probAux));
2223:     PetscCall(PetscDSGetNumFields(probAux, &NfAux));
2224:     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
2225:     PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
2226:     PetscCall(PetscDSGetComponentOffsets(probAux, &aOff));
2227:   }
2228:   /* Allocate data  arrays */
2229:   PetscCall(PetscCalloc1(numCells * totDim, &u));
2230:   if (dmAux) PetscCall(PetscMalloc1(numCells * totDimAux, &a));
2231:   /* Read out geometry */
2232:   PetscCall(DMGetCoordinateField(dm, &coordField));
2233:   PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
2234:   if (maxDegree <= 1) {
2235:     PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &affineQuad));
2236:     if (affineQuad) PetscCall(DMFieldCreateFEGeom(coordField, cellIS, affineQuad, PETSC_FALSE, &cgeomFEM));
2237:   }
2238:   if (useFVM) {
2239:     PetscFV   fv = NULL;
2240:     Vec       grad;
2241:     PetscInt  fStart, fEnd;
2242:     PetscBool compGrad;

2244:     for (f = 0; f < Nf; ++f) {
2245:       PetscObject  obj;
2246:       PetscClassId id;

2248:       PetscCall(PetscDSGetDiscretization(prob, f, &obj));
2249:       PetscCall(PetscObjectGetClassId(obj, &id));
2250:       if (id == PETSCFV_CLASSID) {
2251:         fv = (PetscFV)obj;
2252:         break;
2253:       }
2254:     }
2255:     PetscCall(PetscFVGetComputeGradients(fv, &compGrad));
2256:     PetscCall(PetscFVSetComputeGradients(fv, PETSC_TRUE));
2257:     PetscCall(DMPlexComputeGeometryFVM(dm, &cellGeometryFVM, &faceGeometryFVM));
2258:     PetscCall(DMPlexComputeGradientFVM(dm, fv, faceGeometryFVM, cellGeometryFVM, &dmGrad));
2259:     PetscCall(PetscFVSetComputeGradients(fv, compGrad));
2260:     PetscCall(VecGetArrayRead(cellGeometryFVM, (const PetscScalar **)&cgeomFVM));
2261:     /* Reconstruct and limit cell gradients */
2262:     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
2263:     PetscCall(DMGetGlobalVector(dmGrad, &grad));
2264:     PetscCall(DMPlexReconstructGradients_Internal(dm, fv, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad));
2265:     /* Communicate gradient values */
2266:     PetscCall(DMGetLocalVector(dmGrad, &locGrad));
2267:     PetscCall(DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad));
2268:     PetscCall(DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad));
2269:     PetscCall(DMRestoreGlobalVector(dmGrad, &grad));
2270:     /* Handle non-essential (e.g. outflow) boundary values */
2271:     PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, 0.0, faceGeometryFVM, cellGeometryFVM, locGrad));
2272:     PetscCall(VecGetArrayRead(locGrad, &lgrad));
2273:   }
2274:   /* Read out data from inputs */
2275:   for (c = cStart; c < cEnd; ++c) {
2276:     PetscScalar *x = NULL;
2277:     PetscInt     i;

2279:     PetscCall(DMPlexVecGetClosure(dm, section, locX, c, NULL, &x));
2280:     for (i = 0; i < totDim; ++i) u[c * totDim + i] = x[i];
2281:     PetscCall(DMPlexVecRestoreClosure(dm, section, locX, c, NULL, &x));
2282:     if (dmAux) {
2283:       PetscCall(DMPlexVecGetClosure(plexA, sectionAux, locA, c, NULL, &x));
2284:       for (i = 0; i < totDimAux; ++i) a[c * totDimAux + i] = x[i];
2285:       PetscCall(DMPlexVecRestoreClosure(plexA, sectionAux, locA, c, NULL, &x));
2286:     }
2287:   }
2288:   /* Do integration for each field */
2289:   for (f = 0; f < Nf; ++f) {
2290:     PetscObject  obj;
2291:     PetscClassId id;
2292:     PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;

2294:     PetscCall(PetscDSGetDiscretization(prob, f, &obj));
2295:     PetscCall(PetscObjectGetClassId(obj, &id));
2296:     if (id == PETSCFE_CLASSID) {
2297:       PetscFE         fe = (PetscFE)obj;
2298:       PetscQuadrature q;
2299:       PetscFEGeom    *chunkGeom = NULL;
2300:       PetscInt        Nq, Nb;

2302:       PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
2303:       PetscCall(PetscFEGetQuadrature(fe, &q));
2304:       PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL));
2305:       PetscCall(PetscFEGetDimension(fe, &Nb));
2306:       blockSize = Nb * Nq;
2307:       batchSize = numBlocks * blockSize;
2308:       PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
2309:       numChunks = numCells / (numBatches * batchSize);
2310:       Ne        = numChunks * numBatches * batchSize;
2311:       Nr        = numCells % (numBatches * batchSize);
2312:       offset    = numCells - Nr;
2313:       if (!affineQuad) PetscCall(DMFieldCreateFEGeom(coordField, cellIS, q, PETSC_FALSE, &cgeomFEM));
2314:       PetscCall(PetscFEGeomGetChunk(cgeomFEM, 0, offset, &chunkGeom));
2315:       PetscCall(PetscFEIntegrate(prob, f, Ne, chunkGeom, u, probAux, a, cintegral));
2316:       PetscCall(PetscFEGeomGetChunk(cgeomFEM, offset, numCells, &chunkGeom));
2317:       PetscCall(PetscFEIntegrate(prob, f, Nr, chunkGeom, &u[offset * totDim], probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), &cintegral[offset * Nf]));
2318:       PetscCall(PetscFEGeomRestoreChunk(cgeomFEM, offset, numCells, &chunkGeom));
2319:       if (!affineQuad) PetscCall(PetscFEGeomDestroy(&cgeomFEM));
2320:     } else if (id == PETSCFV_CLASSID) {
2321:       PetscInt       foff;
2322:       PetscPointFunc obj_func;

2324:       PetscCall(PetscDSGetObjective(prob, f, &obj_func));
2325:       PetscCall(PetscDSGetFieldOffset(prob, f, &foff));
2326:       if (obj_func) {
2327:         for (c = 0; c < numCells; ++c) {
2328:           PetscScalar *u_x;
2329:           PetscScalar  lint = 0.;

2331:           PetscCall(DMPlexPointLocalRead(dmGrad, c, lgrad, &u_x));
2332:           obj_func(dim, Nf, NfAux, uOff, uOff_x, &u[totDim * c + foff], NULL, u_x, aOff, NULL, &a[totDimAux * c], NULL, NULL, 0.0, cgeomFVM[c].centroid, numConstants, constants, &lint);
2333:           cintegral[c * Nf + f] += PetscRealPart(lint) * cgeomFVM[c].volume;
2334:         }
2335:       }
2336:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
2337:   }
2338:   /* Cleanup data arrays */
2339:   if (useFVM) {
2340:     PetscCall(VecRestoreArrayRead(locGrad, &lgrad));
2341:     PetscCall(VecRestoreArrayRead(cellGeometryFVM, (const PetscScalar **)&cgeomFVM));
2342:     PetscCall(DMRestoreLocalVector(dmGrad, &locGrad));
2343:     PetscCall(VecDestroy(&faceGeometryFVM));
2344:     PetscCall(VecDestroy(&cellGeometryFVM));
2345:     PetscCall(DMDestroy(&dmGrad));
2346:   }
2347:   if (dmAux) PetscCall(PetscFree(a));
2348:   PetscCall(DMDestroy(&plexA));
2349:   PetscCall(PetscFree(u));
2350:   /* Cleanup */
2351:   if (affineQuad) PetscCall(PetscFEGeomDestroy(&cgeomFEM));
2352:   PetscCall(PetscQuadratureDestroy(&affineQuad));
2353:   PetscCall(ISDestroy(&cellIS));
2354:   PetscFunctionReturn(PETSC_SUCCESS);
2355: }

2357: /*@
2358:   DMPlexComputeIntegralFEM - Form the integral over the domain from the global input X using pointwise functions specified by the user

2360:   Input Parameters:
2361: + dm   - The mesh
2362: . X    - Global input vector
2363: - user - The user context

2365:   Output Parameter:
2366: . integral - Integral for each field

2368:   Level: developer

2370: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSNESComputeResidualFEM()`
2371: @*/
2372: PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscScalar *integral, void *user)
2373: {
2374:   PetscInt     printFEM;
2375:   PetscScalar *cintegral, *lintegral;
2376:   PetscInt     Nf, f, cellHeight, cStart, cEnd, cell;
2377:   Vec          locX;

2379:   PetscFunctionBegin;
2382:   PetscAssertPointer(integral, 3);
2383:   PetscCall(PetscLogEventBegin(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2384:   PetscCall(DMPlexConvertPlex(dm, &dm, PETSC_TRUE));
2385:   PetscCall(DMGetNumFields(dm, &Nf));
2386:   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
2387:   PetscCall(DMPlexGetSimplexOrBoxCells(dm, cellHeight, &cStart, &cEnd));
2388:   /* TODO Introduce a loop over large chunks (right now this is a single chunk) */
2389:   PetscCall(PetscCalloc2(Nf, &lintegral, (cEnd - cStart) * Nf, &cintegral));
2390:   /* Get local solution with boundary values */
2391:   PetscCall(DMGetLocalVector(dm, &locX));
2392:   PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locX, 0.0, NULL, NULL, NULL));
2393:   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
2394:   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
2395:   PetscCall(DMPlexComputeIntegral_Internal(dm, locX, cStart, cEnd, cintegral, user));
2396:   PetscCall(DMRestoreLocalVector(dm, &locX));
2397:   printFEM = ((DM_Plex *)dm->data)->printFEM;
2398:   /* Sum up values */
2399:   for (cell = cStart; cell < cEnd; ++cell) {
2400:     const PetscInt c = cell - cStart;

2402:     if (printFEM > 1) PetscCall(DMPrintCellVector(cell, "Cell Integral", Nf, &cintegral[c * Nf]));
2403:     for (f = 0; f < Nf; ++f) lintegral[f] += cintegral[c * Nf + f];
2404:   }
2405:   PetscCall(MPIU_Allreduce(lintegral, integral, Nf, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)dm)));
2406:   if (printFEM) {
2407:     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Integral:"));
2408:     for (f = 0; f < Nf; ++f) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), " %g", (double)PetscRealPart(integral[f])));
2409:     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "\n"));
2410:   }
2411:   PetscCall(PetscFree2(lintegral, cintegral));
2412:   PetscCall(PetscLogEventEnd(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2413:   PetscCall(DMDestroy(&dm));
2414:   PetscFunctionReturn(PETSC_SUCCESS);
2415: }

2417: /*@
2418:   DMPlexComputeCellwiseIntegralFEM - Form the vector of cellwise integrals F from the global input X using pointwise functions specified by the user

2420:   Input Parameters:
2421: + dm   - The mesh
2422: . X    - Global input vector
2423: - user - The user context

2425:   Output Parameter:
2426: . F - Cellwise integrals for each field

2428:   Level: developer

2430: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSNESComputeResidualFEM()`
2431: @*/
2432: PetscErrorCode DMPlexComputeCellwiseIntegralFEM(DM dm, Vec X, Vec F, void *user)
2433: {
2434:   PetscInt     printFEM;
2435:   DM           dmF;
2436:   PetscSection sectionF = NULL;
2437:   PetscScalar *cintegral, *af;
2438:   PetscInt     Nf, f, cellHeight, cStart, cEnd, cell, n;
2439:   Vec          locX;

2441:   PetscFunctionBegin;
2445:   PetscCall(PetscLogEventBegin(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2446:   PetscCall(DMPlexConvertPlex(dm, &dm, PETSC_TRUE));
2447:   PetscCall(DMGetNumFields(dm, &Nf));
2448:   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
2449:   PetscCall(DMPlexGetSimplexOrBoxCells(dm, cellHeight, &cStart, &cEnd));
2450:   /* TODO Introduce a loop over large chunks (right now this is a single chunk) */
2451:   PetscCall(PetscCalloc1((cEnd - cStart) * Nf, &cintegral));
2452:   /* Get local solution with boundary values */
2453:   PetscCall(DMGetLocalVector(dm, &locX));
2454:   PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locX, 0.0, NULL, NULL, NULL));
2455:   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
2456:   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
2457:   PetscCall(DMPlexComputeIntegral_Internal(dm, locX, cStart, cEnd, cintegral, user));
2458:   PetscCall(DMRestoreLocalVector(dm, &locX));
2459:   /* Put values in F */
2460:   PetscCall(VecGetArray(F, &af));
2461:   PetscCall(VecGetDM(F, &dmF));
2462:   if (dmF) PetscCall(DMGetLocalSection(dmF, &sectionF));
2463:   PetscCall(VecGetLocalSize(F, &n));
2464:   PetscCheck(n >= (cEnd - cStart) * Nf, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Vector size %" PetscInt_FMT " < %" PetscInt_FMT, n, (cEnd - cStart) * Nf);
2465:   printFEM = ((DM_Plex *)dm->data)->printFEM;
2466:   for (cell = cStart; cell < cEnd; ++cell) {
2467:     const PetscInt c   = cell - cStart;
2468:     PetscInt       dof = Nf, off = c * Nf;

2470:     if (printFEM > 1) PetscCall(DMPrintCellVector(cell, "Cell Integral", Nf, &cintegral[c * Nf]));
2471:     if (sectionF) {
2472:       PetscCall(PetscSectionGetDof(sectionF, cell, &dof));
2473:       PetscCall(PetscSectionGetOffset(sectionF, cell, &off));
2474:     }
2475:     PetscCheck(dof == Nf, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of cell dofs %" PetscInt_FMT " != %" PetscInt_FMT, dof, Nf);
2476:     for (f = 0; f < Nf; ++f) af[off + f] = cintegral[c * Nf + f];
2477:   }
2478:   PetscCall(VecRestoreArray(F, &af));
2479:   PetscCall(PetscFree(cintegral));
2480:   PetscCall(PetscLogEventEnd(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2481:   PetscCall(DMDestroy(&dm));
2482:   PetscFunctionReturn(PETSC_SUCCESS);
2483: }

2485: static PetscErrorCode DMPlexComputeBdIntegral_Internal(DM dm, Vec locX, IS pointIS, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscScalar *fintegral, void *user)
2486: {
2487:   DM                 plex = NULL, plexA = NULL;
2488:   DMEnclosureType    encAux;
2489:   PetscDS            prob, probAux       = NULL;
2490:   PetscSection       section, sectionAux = NULL;
2491:   Vec                locA = NULL;
2492:   DMField            coordField;
2493:   PetscInt           Nf, totDim, *uOff, *uOff_x;
2494:   PetscInt           NfAux = 0, totDimAux = 0, *aOff = NULL;
2495:   PetscScalar       *u, *a = NULL;
2496:   const PetscScalar *constants;
2497:   PetscInt           numConstants, f;

2499:   PetscFunctionBegin;
2500:   PetscCall(DMGetCoordinateField(dm, &coordField));
2501:   PetscCall(DMConvert(dm, DMPLEX, &plex));
2502:   PetscCall(DMGetDS(dm, &prob));
2503:   PetscCall(DMGetLocalSection(dm, &section));
2504:   PetscCall(PetscSectionGetNumFields(section, &Nf));
2505:   /* Determine which discretizations we have */
2506:   for (f = 0; f < Nf; ++f) {
2507:     PetscObject  obj;
2508:     PetscClassId id;

2510:     PetscCall(PetscDSGetDiscretization(prob, f, &obj));
2511:     PetscCall(PetscObjectGetClassId(obj, &id));
2512:     PetscCheck(id != PETSCFV_CLASSID, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not supported for FVM (field %" PetscInt_FMT ")", f);
2513:   }
2514:   /* Read DS information */
2515:   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
2516:   PetscCall(PetscDSGetComponentOffsets(prob, &uOff));
2517:   PetscCall(PetscDSGetComponentDerivativeOffsets(prob, &uOff_x));
2518:   PetscCall(PetscDSGetConstants(prob, &numConstants, &constants));
2519:   /* Read Auxiliary DS information */
2520:   PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &locA));
2521:   if (locA) {
2522:     DM dmAux;

2524:     PetscCall(VecGetDM(locA, &dmAux));
2525:     PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
2526:     PetscCall(DMConvert(dmAux, DMPLEX, &plexA));
2527:     PetscCall(DMGetDS(dmAux, &probAux));
2528:     PetscCall(PetscDSGetNumFields(probAux, &NfAux));
2529:     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
2530:     PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
2531:     PetscCall(PetscDSGetComponentOffsets(probAux, &aOff));
2532:   }
2533:   /* Integrate over points */
2534:   {
2535:     PetscFEGeom    *fgeom, *chunkGeom = NULL;
2536:     PetscInt        maxDegree;
2537:     PetscQuadrature qGeom = NULL;
2538:     const PetscInt *points;
2539:     PetscInt        numFaces, face, Nq, field;
2540:     PetscInt        numChunks, chunkSize, chunk, Nr, offset;

2542:     PetscCall(ISGetLocalSize(pointIS, &numFaces));
2543:     PetscCall(ISGetIndices(pointIS, &points));
2544:     PetscCall(PetscCalloc2(numFaces * totDim, &u, locA ? numFaces * totDimAux : 0, &a));
2545:     PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
2546:     for (face = 0; face < numFaces; ++face) {
2547:       const PetscInt point = points[face], *support;
2548:       PetscScalar   *x     = NULL;

2550:       PetscCall(DMPlexGetSupport(dm, point, &support));
2551:       PetscCall(DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x));
2552:       for (PetscInt i = 0; i < totDim; ++i) u[face * totDim + i] = x[i];
2553:       PetscCall(DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x));
2554:       if (locA) {
2555:         PetscInt subp;
2556:         PetscCall(DMGetEnclosurePoint(plexA, dm, encAux, support[0], &subp));
2557:         PetscCall(DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x));
2558:         for (PetscInt i = 0; i < totDimAux; ++i) a[f * totDimAux + i] = x[i];
2559:         PetscCall(DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x));
2560:       }
2561:     }
2562:     for (field = 0; field < Nf; ++field) {
2563:       PetscFE fe;

2565:       PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe));
2566:       if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &qGeom));
2567:       if (!qGeom) {
2568:         PetscCall(PetscFEGetFaceQuadrature(fe, &qGeom));
2569:         PetscCall(PetscObjectReference((PetscObject)qGeom));
2570:       }
2571:       PetscCall(PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL));
2572:       PetscCall(DMPlexGetFEGeom(coordField, pointIS, qGeom, PETSC_TRUE, &fgeom));
2573:       /* Get blocking */
2574:       {
2575:         PetscQuadrature q;
2576:         PetscInt        numBatches, batchSize, numBlocks, blockSize;
2577:         PetscInt        Nq, Nb;

2579:         PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
2580:         PetscCall(PetscFEGetQuadrature(fe, &q));
2581:         PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL));
2582:         PetscCall(PetscFEGetDimension(fe, &Nb));
2583:         blockSize = Nb * Nq;
2584:         batchSize = numBlocks * blockSize;
2585:         chunkSize = numBatches * batchSize;
2586:         PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
2587:         numChunks = numFaces / chunkSize;
2588:         Nr        = numFaces % chunkSize;
2589:         offset    = numFaces - Nr;
2590:       }
2591:       /* Do integration for each field */
2592:       for (chunk = 0; chunk < numChunks; ++chunk) {
2593:         PetscCall(PetscFEGeomGetChunk(fgeom, chunk * chunkSize, (chunk + 1) * chunkSize, &chunkGeom));
2594:         PetscCall(PetscFEIntegrateBd(prob, field, funcs[field], chunkSize, chunkGeom, &u[chunk * chunkSize * totDim], probAux, PetscSafePointerPlusOffset(a, chunk * chunkSize * totDimAux), &fintegral[chunk * chunkSize * Nf]));
2595:         PetscCall(PetscFEGeomRestoreChunk(fgeom, 0, offset, &chunkGeom));
2596:       }
2597:       PetscCall(PetscFEGeomGetChunk(fgeom, offset, numFaces, &chunkGeom));
2598:       PetscCall(PetscFEIntegrateBd(prob, field, funcs[field], Nr, chunkGeom, &u[offset * totDim], probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), &fintegral[offset * Nf]));
2599:       PetscCall(PetscFEGeomRestoreChunk(fgeom, offset, numFaces, &chunkGeom));
2600:       /* Cleanup data arrays */
2601:       PetscCall(DMPlexRestoreFEGeom(coordField, pointIS, qGeom, PETSC_TRUE, &fgeom));
2602:       PetscCall(PetscQuadratureDestroy(&qGeom));
2603:     }
2604:     PetscCall(PetscFree2(u, a));
2605:     PetscCall(ISRestoreIndices(pointIS, &points));
2606:   }
2607:   if (plex) PetscCall(DMDestroy(&plex));
2608:   if (plexA) PetscCall(DMDestroy(&plexA));
2609:   PetscFunctionReturn(PETSC_SUCCESS);
2610: }

2612: /*@C
2613:   DMPlexComputeBdIntegral - Form the integral over the specified boundary from the global input X using pointwise functions specified by the user

2615:   Input Parameters:
2616: + dm      - The mesh
2617: . X       - Global input vector
2618: . label   - The boundary `DMLabel`
2619: . numVals - The number of label values to use, or `PETSC_DETERMINE` for all values
2620: . vals    - The label values to use, or NULL for all values
2621: . funcs   - The functions to integrate along the boundary for each field
2622: - user    - The user context

2624:   Output Parameter:
2625: . integral - Integral for each field

2627:   Level: developer

2629: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexComputeIntegralFEM()`, `DMPlexComputeBdResidualFEM()`
2630: @*/
2631: PetscErrorCode DMPlexComputeBdIntegral(DM dm, Vec X, DMLabel label, PetscInt numVals, const PetscInt vals[], void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscScalar *integral, void *user)
2632: {
2633:   Vec          locX;
2634:   PetscSection section;
2635:   DMLabel      depthLabel;
2636:   IS           facetIS;
2637:   PetscInt     dim, Nf, f, v;

2639:   PetscFunctionBegin;
2643:   if (vals) PetscAssertPointer(vals, 5);
2644:   PetscAssertPointer(integral, 7);
2645:   PetscCall(PetscLogEventBegin(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2646:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
2647:   PetscCall(DMGetDimension(dm, &dim));
2648:   PetscCall(DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS));
2649:   PetscCall(DMGetLocalSection(dm, &section));
2650:   PetscCall(PetscSectionGetNumFields(section, &Nf));
2651:   /* Get local solution with boundary values */
2652:   PetscCall(DMGetLocalVector(dm, &locX));
2653:   PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locX, 0.0, NULL, NULL, NULL));
2654:   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
2655:   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
2656:   /* Loop over label values */
2657:   PetscCall(PetscArrayzero(integral, Nf));
2658:   for (v = 0; v < numVals; ++v) {
2659:     IS           pointIS;
2660:     PetscInt     numFaces, face;
2661:     PetscScalar *fintegral;

2663:     PetscCall(DMLabelGetStratumIS(label, vals[v], &pointIS));
2664:     if (!pointIS) continue; /* No points with that id on this process */
2665:     {
2666:       IS isectIS;

2668:       /* TODO: Special cases of ISIntersect where it is quick to check a priori if one is a superset of the other */
2669:       PetscCall(ISIntersect_Caching_Internal(facetIS, pointIS, &isectIS));
2670:       PetscCall(ISDestroy(&pointIS));
2671:       pointIS = isectIS;
2672:     }
2673:     PetscCall(ISGetLocalSize(pointIS, &numFaces));
2674:     PetscCall(PetscCalloc1(numFaces * Nf, &fintegral));
2675:     PetscCall(DMPlexComputeBdIntegral_Internal(dm, locX, pointIS, funcs, fintegral, user));
2676:     /* Sum point contributions into integral */
2677:     for (f = 0; f < Nf; ++f)
2678:       for (face = 0; face < numFaces; ++face) integral[f] += fintegral[face * Nf + f];
2679:     PetscCall(PetscFree(fintegral));
2680:     PetscCall(ISDestroy(&pointIS));
2681:   }
2682:   PetscCall(DMRestoreLocalVector(dm, &locX));
2683:   PetscCall(ISDestroy(&facetIS));
2684:   PetscCall(PetscLogEventEnd(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2685:   PetscFunctionReturn(PETSC_SUCCESS);
2686: }

2688: /*@
2689:   DMPlexComputeInterpolatorNested - Form the local portion of the interpolation matrix from the coarse `DM` to a uniformly refined `DM`.

2691:   Input Parameters:
2692: + dmc       - The coarse mesh
2693: . dmf       - The fine mesh
2694: . isRefined - Flag indicating regular refinement, rather than the same topology
2695: - user      - The user context

2697:   Output Parameter:
2698: . In - The interpolation matrix

2700:   Level: developer

2702: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexComputeInterpolatorGeneral()`, `DMPlexComputeJacobianFEM()`
2703: @*/
2704: PetscErrorCode DMPlexComputeInterpolatorNested(DM dmc, DM dmf, PetscBool isRefined, Mat In, void *user)
2705: {
2706:   DM_Plex     *mesh = (DM_Plex *)dmc->data;
2707:   const char  *name = "Interpolator";
2708:   PetscFE     *feRef;
2709:   PetscFV     *fvRef;
2710:   PetscSection fsection, fglobalSection;
2711:   PetscSection csection, cglobalSection;
2712:   PetscScalar *elemMat;
2713:   PetscInt     dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
2714:   PetscInt     cTotDim = 0, rTotDim = 0;

2716:   PetscFunctionBegin;
2717:   PetscCall(PetscLogEventBegin(DMPLEX_InterpolatorFEM, dmc, dmf, 0, 0));
2718:   PetscCall(DMGetDimension(dmf, &dim));
2719:   PetscCall(DMGetLocalSection(dmf, &fsection));
2720:   PetscCall(DMGetGlobalSection(dmf, &fglobalSection));
2721:   PetscCall(DMGetLocalSection(dmc, &csection));
2722:   PetscCall(DMGetGlobalSection(dmc, &cglobalSection));
2723:   PetscCall(PetscSectionGetNumFields(fsection, &Nf));
2724:   PetscCall(DMPlexGetSimplexOrBoxCells(dmc, 0, &cStart, &cEnd));
2725:   PetscCall(PetscCalloc2(Nf, &feRef, Nf, &fvRef));
2726:   for (f = 0; f < Nf; ++f) {
2727:     PetscObject  obj, objc;
2728:     PetscClassId id, idc;
2729:     PetscInt     rNb = 0, Nc = 0, cNb = 0;

2731:     PetscCall(DMGetField(dmf, f, NULL, &obj));
2732:     PetscCall(PetscObjectGetClassId(obj, &id));
2733:     if (id == PETSCFE_CLASSID) {
2734:       PetscFE fe = (PetscFE)obj;

2736:       if (isRefined) {
2737:         PetscCall(PetscFERefine(fe, &feRef[f]));
2738:       } else {
2739:         PetscCall(PetscObjectReference((PetscObject)fe));
2740:         feRef[f] = fe;
2741:       }
2742:       PetscCall(PetscFEGetDimension(feRef[f], &rNb));
2743:       PetscCall(PetscFEGetNumComponents(fe, &Nc));
2744:     } else if (id == PETSCFV_CLASSID) {
2745:       PetscFV        fv = (PetscFV)obj;
2746:       PetscDualSpace Q;

2748:       if (isRefined) {
2749:         PetscCall(PetscFVRefine(fv, &fvRef[f]));
2750:       } else {
2751:         PetscCall(PetscObjectReference((PetscObject)fv));
2752:         fvRef[f] = fv;
2753:       }
2754:       PetscCall(PetscFVGetDualSpace(fvRef[f], &Q));
2755:       PetscCall(PetscDualSpaceGetDimension(Q, &rNb));
2756:       PetscCall(PetscFVGetDualSpace(fv, &Q));
2757:       PetscCall(PetscFVGetNumComponents(fv, &Nc));
2758:     }
2759:     PetscCall(DMGetField(dmc, f, NULL, &objc));
2760:     PetscCall(PetscObjectGetClassId(objc, &idc));
2761:     if (idc == PETSCFE_CLASSID) {
2762:       PetscFE fe = (PetscFE)objc;

2764:       PetscCall(PetscFEGetDimension(fe, &cNb));
2765:     } else if (id == PETSCFV_CLASSID) {
2766:       PetscFV        fv = (PetscFV)obj;
2767:       PetscDualSpace Q;

2769:       PetscCall(PetscFVGetDualSpace(fv, &Q));
2770:       PetscCall(PetscDualSpaceGetDimension(Q, &cNb));
2771:     }
2772:     rTotDim += rNb;
2773:     cTotDim += cNb;
2774:   }
2775:   PetscCall(PetscMalloc1(rTotDim * cTotDim, &elemMat));
2776:   PetscCall(PetscArrayzero(elemMat, rTotDim * cTotDim));
2777:   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
2778:     PetscDualSpace   Qref;
2779:     PetscQuadrature  f;
2780:     const PetscReal *qpoints, *qweights;
2781:     PetscReal       *points;
2782:     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;

2784:     /* Compose points from all dual basis functionals */
2785:     if (feRef[fieldI]) {
2786:       PetscCall(PetscFEGetDualSpace(feRef[fieldI], &Qref));
2787:       PetscCall(PetscFEGetNumComponents(feRef[fieldI], &Nc));
2788:     } else {
2789:       PetscCall(PetscFVGetDualSpace(fvRef[fieldI], &Qref));
2790:       PetscCall(PetscFVGetNumComponents(fvRef[fieldI], &Nc));
2791:     }
2792:     PetscCall(PetscDualSpaceGetDimension(Qref, &fpdim));
2793:     for (i = 0; i < fpdim; ++i) {
2794:       PetscCall(PetscDualSpaceGetFunctional(Qref, i, &f));
2795:       PetscCall(PetscQuadratureGetData(f, NULL, NULL, &Np, NULL, NULL));
2796:       npoints += Np;
2797:     }
2798:     PetscCall(PetscMalloc1(npoints * dim, &points));
2799:     for (i = 0, k = 0; i < fpdim; ++i) {
2800:       PetscCall(PetscDualSpaceGetFunctional(Qref, i, &f));
2801:       PetscCall(PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL));
2802:       for (p = 0; p < Np; ++p, ++k)
2803:         for (d = 0; d < dim; ++d) points[k * dim + d] = qpoints[p * dim + d];
2804:     }

2806:     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
2807:       PetscObject  obj;
2808:       PetscClassId id;
2809:       PetscInt     NcJ = 0, cpdim = 0, j, qNc;

2811:       PetscCall(DMGetField(dmc, fieldJ, NULL, &obj));
2812:       PetscCall(PetscObjectGetClassId(obj, &id));
2813:       if (id == PETSCFE_CLASSID) {
2814:         PetscFE         fe = (PetscFE)obj;
2815:         PetscTabulation T  = NULL;

2817:         /* Evaluate basis at points */
2818:         PetscCall(PetscFEGetNumComponents(fe, &NcJ));
2819:         PetscCall(PetscFEGetDimension(fe, &cpdim));
2820:         /* For now, fields only interpolate themselves */
2821:         if (fieldI == fieldJ) {
2822:           PetscCheck(Nc == NcJ, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %" PetscInt_FMT " does not match coarse field %" PetscInt_FMT, Nc, NcJ);
2823:           PetscCall(PetscFECreateTabulation(fe, 1, npoints, points, 0, &T));
2824:           for (i = 0, k = 0; i < fpdim; ++i) {
2825:             PetscCall(PetscDualSpaceGetFunctional(Qref, i, &f));
2826:             PetscCall(PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights));
2827:             PetscCheck(qNc == NcJ, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %" PetscInt_FMT " does not match coarse field %" PetscInt_FMT, qNc, NcJ);
2828:             for (p = 0; p < Np; ++p, ++k) {
2829:               for (j = 0; j < cpdim; ++j) {
2830:                 /*
2831:                    cTotDim:            Total columns in element interpolation matrix, sum of number of dual basis functionals in each field
2832:                    offsetI, offsetJ:   Offsets into the larger element interpolation matrix for different fields
2833:                    fpdim, i, cpdim, j: Dofs for fine and coarse grids, correspond to dual space basis functionals
2834:                    qNC, Nc, Ncj, c:    Number of components in this field
2835:                    Np, p:              Number of quad points in the fine grid functional i
2836:                    k:                  i*Np + p, overall point number for the interpolation
2837:                 */
2838:                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i) * cTotDim + offsetJ + j] += T->T[0][k * cpdim * NcJ + j * Nc + c] * qweights[p * qNc + c];
2839:               }
2840:             }
2841:           }
2842:           PetscCall(PetscTabulationDestroy(&T));
2843:         }
2844:       } else if (id == PETSCFV_CLASSID) {
2845:         PetscFV fv = (PetscFV)obj;

2847:         /* Evaluate constant function at points */
2848:         PetscCall(PetscFVGetNumComponents(fv, &NcJ));
2849:         cpdim = 1;
2850:         /* For now, fields only interpolate themselves */
2851:         if (fieldI == fieldJ) {
2852:           PetscCheck(Nc == NcJ, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %" PetscInt_FMT " does not match coarse field %" PetscInt_FMT, Nc, NcJ);
2853:           for (i = 0, k = 0; i < fpdim; ++i) {
2854:             PetscCall(PetscDualSpaceGetFunctional(Qref, i, &f));
2855:             PetscCall(PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights));
2856:             PetscCheck(qNc == NcJ, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %" PetscInt_FMT " does not match coarse field %" PetscInt_FMT, qNc, NcJ);
2857:             for (p = 0; p < Np; ++p, ++k) {
2858:               for (j = 0; j < cpdim; ++j) {
2859:                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i) * cTotDim + offsetJ + j] += 1.0 * qweights[p * qNc + c];
2860:               }
2861:             }
2862:           }
2863:         }
2864:       }
2865:       offsetJ += cpdim;
2866:     }
2867:     offsetI += fpdim;
2868:     PetscCall(PetscFree(points));
2869:   }
2870:   if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat));
2871:   /* Preallocate matrix */
2872:   {
2873:     Mat          preallocator;
2874:     PetscScalar *vals;
2875:     PetscInt    *cellCIndices, *cellFIndices;
2876:     PetscInt     locRows, locCols, cell;

2878:     PetscCall(MatGetLocalSize(In, &locRows, &locCols));
2879:     PetscCall(MatCreate(PetscObjectComm((PetscObject)In), &preallocator));
2880:     PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
2881:     PetscCall(MatSetSizes(preallocator, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE));
2882:     PetscCall(MatSetUp(preallocator));
2883:     PetscCall(PetscCalloc3(rTotDim * cTotDim, &vals, cTotDim, &cellCIndices, rTotDim, &cellFIndices));
2884:     for (cell = cStart; cell < cEnd; ++cell) {
2885:       if (isRefined) {
2886:         PetscCall(DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices));
2887:         PetscCall(MatSetValues(preallocator, rTotDim, cellFIndices, cTotDim, cellCIndices, vals, INSERT_VALUES));
2888:       } else {
2889:         PetscCall(DMPlexMatSetClosureGeneral(dmf, fsection, fglobalSection, PETSC_FALSE, dmc, csection, cglobalSection, PETSC_FALSE, preallocator, cell, vals, INSERT_VALUES));
2890:       }
2891:     }
2892:     PetscCall(PetscFree3(vals, cellCIndices, cellFIndices));
2893:     PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
2894:     PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
2895:     PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, In));
2896:     PetscCall(MatDestroy(&preallocator));
2897:   }
2898:   /* Fill matrix */
2899:   PetscCall(MatZeroEntries(In));
2900:   for (c = cStart; c < cEnd; ++c) {
2901:     if (isRefined) {
2902:       PetscCall(DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES));
2903:     } else {
2904:       PetscCall(DMPlexMatSetClosureGeneral(dmf, fsection, fglobalSection, PETSC_FALSE, dmc, csection, cglobalSection, PETSC_FALSE, In, c, elemMat, INSERT_VALUES));
2905:     }
2906:   }
2907:   for (f = 0; f < Nf; ++f) PetscCall(PetscFEDestroy(&feRef[f]));
2908:   PetscCall(PetscFree2(feRef, fvRef));
2909:   PetscCall(PetscFree(elemMat));
2910:   PetscCall(MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY));
2911:   PetscCall(MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY));
2912:   if (mesh->printFEM > 1) {
2913:     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)In), "%s:\n", name));
2914:     PetscCall(MatFilter(In, 1.0e-10, PETSC_FALSE, PETSC_FALSE));
2915:     PetscCall(MatView(In, NULL));
2916:   }
2917:   PetscCall(PetscLogEventEnd(DMPLEX_InterpolatorFEM, dmc, dmf, 0, 0));
2918:   PetscFunctionReturn(PETSC_SUCCESS);
2919: }

2921: PetscErrorCode DMPlexComputeMassMatrixNested(DM dmc, DM dmf, Mat mass, void *user)
2922: {
2923:   SETERRQ(PetscObjectComm((PetscObject)dmc), PETSC_ERR_SUP, "Laziness");
2924: }

2926: /*@
2927:   DMPlexComputeInterpolatorGeneral - Form the local portion of the interpolation matrix from the coarse `DM` to a non-nested fine `DM`.

2929:   Input Parameters:
2930: + dmf  - The fine mesh
2931: . dmc  - The coarse mesh
2932: - user - The user context

2934:   Output Parameter:
2935: . In - The interpolation matrix

2937:   Level: developer

2939: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexComputeInterpolatorNested()`, `DMPlexComputeJacobianFEM()`
2940: @*/
2941: PetscErrorCode DMPlexComputeInterpolatorGeneral(DM dmc, DM dmf, Mat In, void *user)
2942: {
2943:   DM_Plex     *mesh = (DM_Plex *)dmf->data;
2944:   const char  *name = "Interpolator";
2945:   PetscDS      prob;
2946:   Mat          interp;
2947:   PetscSection fsection, globalFSection;
2948:   PetscSection csection, globalCSection;
2949:   PetscInt     locRows, locCols;
2950:   PetscReal   *x, *v0, *J, *invJ, detJ;
2951:   PetscReal   *v0c, *Jc, *invJc, detJc;
2952:   PetscScalar *elemMat;
2953:   PetscInt     dim, Nf, field, totDim, cStart, cEnd, cell, ccell, s;

2955:   PetscFunctionBegin;
2956:   PetscCall(PetscLogEventBegin(DMPLEX_InterpolatorFEM, dmc, dmf, 0, 0));
2957:   PetscCall(DMGetCoordinateDim(dmc, &dim));
2958:   PetscCall(DMGetDS(dmc, &prob));
2959:   PetscCall(PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL));
2960:   PetscCall(PetscDSGetNumFields(prob, &Nf));
2961:   PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ));
2962:   PetscCall(PetscMalloc3(dim, &v0c, dim * dim, &Jc, dim * dim, &invJc));
2963:   PetscCall(DMGetLocalSection(dmf, &fsection));
2964:   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
2965:   PetscCall(DMGetLocalSection(dmc, &csection));
2966:   PetscCall(DMGetGlobalSection(dmc, &globalCSection));
2967:   PetscCall(DMPlexGetSimplexOrBoxCells(dmf, 0, &cStart, &cEnd));
2968:   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
2969:   PetscCall(PetscMalloc1(totDim, &elemMat));

2971:   PetscCall(MatGetLocalSize(In, &locRows, &locCols));
2972:   PetscCall(MatCreate(PetscObjectComm((PetscObject)In), &interp));
2973:   PetscCall(MatSetType(interp, MATPREALLOCATOR));
2974:   PetscCall(MatSetSizes(interp, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE));
2975:   PetscCall(MatSetUp(interp));
2976:   for (s = 0; s < 2; ++s) {
2977:     for (field = 0; field < Nf; ++field) {
2978:       PetscObject      obj;
2979:       PetscClassId     id;
2980:       PetscDualSpace   Q = NULL;
2981:       PetscTabulation  T = NULL;
2982:       PetscQuadrature  f;
2983:       const PetscReal *qpoints, *qweights;
2984:       PetscInt         Nc, qNc, Np, fpdim, off, i, d;

2986:       PetscCall(PetscDSGetFieldOffset(prob, field, &off));
2987:       PetscCall(PetscDSGetDiscretization(prob, field, &obj));
2988:       PetscCall(PetscObjectGetClassId(obj, &id));
2989:       if (id == PETSCFE_CLASSID) {
2990:         PetscFE fe = (PetscFE)obj;

2992:         PetscCall(PetscFEGetDualSpace(fe, &Q));
2993:         PetscCall(PetscFEGetNumComponents(fe, &Nc));
2994:         if (s) PetscCall(PetscFECreateTabulation(fe, 1, 1, x, 0, &T));
2995:       } else if (id == PETSCFV_CLASSID) {
2996:         PetscFV fv = (PetscFV)obj;

2998:         PetscCall(PetscFVGetDualSpace(fv, &Q));
2999:         Nc = 1;
3000:       } else SETERRQ(PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
3001:       PetscCall(PetscDualSpaceGetDimension(Q, &fpdim));
3002:       /* For each fine grid cell */
3003:       for (cell = cStart; cell < cEnd; ++cell) {
3004:         PetscInt *findices, *cindices;
3005:         PetscInt  numFIndices, numCIndices;

3007:         PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3008:         PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
3009:         PetscCheck(numFIndices == totDim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %" PetscInt_FMT " != %" PetscInt_FMT " dual basis vecs", numFIndices, totDim);
3010:         for (i = 0; i < fpdim; ++i) {
3011:           Vec                pointVec;
3012:           PetscScalar       *pV;
3013:           PetscSF            coarseCellSF = NULL;
3014:           const PetscSFNode *coarseCells;
3015:           PetscInt           numCoarseCells, cpdim, row = findices[i + off], q, c, j;

3017:           /* Get points from the dual basis functional quadrature */
3018:           PetscCall(PetscDualSpaceGetFunctional(Q, i, &f));
3019:           PetscCall(PetscQuadratureGetData(f, NULL, &qNc, &Np, &qpoints, &qweights));
3020:           PetscCheck(qNc == Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %" PetscInt_FMT " does not match coarse field %" PetscInt_FMT, qNc, Nc);
3021:           PetscCall(VecCreateSeq(PETSC_COMM_SELF, Np * dim, &pointVec));
3022:           PetscCall(VecSetBlockSize(pointVec, dim));
3023:           PetscCall(VecGetArray(pointVec, &pV));
3024:           for (q = 0; q < Np; ++q) {
3025:             const PetscReal xi0[3] = {-1., -1., -1.};

3027:             /* Transform point to real space */
3028:             CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q * dim], x);
3029:             for (d = 0; d < dim; ++d) pV[q * dim + d] = x[d];
3030:           }
3031:           PetscCall(VecRestoreArray(pointVec, &pV));
3032:           /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
3033:           /* OPT: Read this out from preallocation information */
3034:           PetscCall(DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF));
3035:           /* Update preallocation info */
3036:           PetscCall(PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells));
3037:           PetscCheck(numCoarseCells == Np, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not all closure points located");
3038:           PetscCall(VecGetArray(pointVec, &pV));
3039:           for (ccell = 0; ccell < numCoarseCells; ++ccell) {
3040:             PetscReal       pVReal[3];
3041:             const PetscReal xi0[3] = {-1., -1., -1.};

3043:             PetscCall(DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL));
3044:             if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetDimension((PetscFE)obj, &cpdim));
3045:             else cpdim = 1;

3047:             if (s) {
3048:               /* Transform points from real space to coarse reference space */
3049:               PetscCall(DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc));
3050:               for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell * dim + d]);
3051:               CoordinatesRealToRef(dim, dim, xi0, v0c, invJc, pVReal, x);

3053:               if (id == PETSCFE_CLASSID) {
3054:                 /* Evaluate coarse basis on contained point */
3055:                 PetscCall(PetscFEComputeTabulation((PetscFE)obj, 1, x, 0, T));
3056:                 PetscCall(PetscArrayzero(elemMat, cpdim));
3057:                 /* Get elemMat entries by multiplying by weight */
3058:                 for (j = 0; j < cpdim; ++j) {
3059:                   for (c = 0; c < Nc; ++c) elemMat[j] += T->T[0][j * Nc + c] * qweights[ccell * qNc + c];
3060:                 }
3061:               } else {
3062:                 for (j = 0; j < cpdim; ++j) {
3063:                   for (c = 0; c < Nc; ++c) elemMat[j] += 1.0 * qweights[ccell * qNc + c];
3064:                 }
3065:               }
3066:               if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
3067:             }
3068:             /* Update interpolator */
3069:             PetscCheck(numCIndices == totDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %" PetscInt_FMT " != %" PetscInt_FMT, numCIndices, totDim);
3070:             PetscCall(MatSetValues(interp, 1, &row, cpdim, &cindices[off], elemMat, INSERT_VALUES));
3071:             PetscCall(DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL));
3072:           }
3073:           PetscCall(VecRestoreArray(pointVec, &pV));
3074:           PetscCall(PetscSFDestroy(&coarseCellSF));
3075:           PetscCall(VecDestroy(&pointVec));
3076:         }
3077:         PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3078:       }
3079:       if (s && id == PETSCFE_CLASSID) PetscCall(PetscTabulationDestroy(&T));
3080:     }
3081:     if (!s) {
3082:       PetscCall(MatAssemblyBegin(interp, MAT_FINAL_ASSEMBLY));
3083:       PetscCall(MatAssemblyEnd(interp, MAT_FINAL_ASSEMBLY));
3084:       PetscCall(MatPreallocatorPreallocate(interp, PETSC_TRUE, In));
3085:       PetscCall(MatDestroy(&interp));
3086:       interp = In;
3087:     }
3088:   }
3089:   PetscCall(PetscFree3(v0, J, invJ));
3090:   PetscCall(PetscFree3(v0c, Jc, invJc));
3091:   PetscCall(PetscFree(elemMat));
3092:   PetscCall(MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY));
3093:   PetscCall(MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY));
3094:   PetscCall(PetscLogEventEnd(DMPLEX_InterpolatorFEM, dmc, dmf, 0, 0));
3095:   PetscFunctionReturn(PETSC_SUCCESS);
3096: }

3098: /*@
3099:   DMPlexComputeMassMatrixGeneral - Form the local portion of the mass matrix from the coarse `DM` to a non-nested fine `DM`.

3101:   Input Parameters:
3102: + dmf  - The fine mesh
3103: . dmc  - The coarse mesh
3104: - user - The user context

3106:   Output Parameter:
3107: . mass - The mass matrix

3109:   Level: developer

3111: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexComputeMassMatrixNested()`, `DMPlexComputeInterpolatorNested()`, `DMPlexComputeInterpolatorGeneral()`, `DMPlexComputeJacobianFEM()`
3112: @*/
3113: PetscErrorCode DMPlexComputeMassMatrixGeneral(DM dmc, DM dmf, Mat mass, void *user)
3114: {
3115:   DM_Plex     *mesh = (DM_Plex *)dmf->data;
3116:   const char  *name = "Mass Matrix";
3117:   PetscDS      prob;
3118:   PetscSection fsection, csection, globalFSection, globalCSection;
3119:   PetscHSetIJ  ht;
3120:   PetscLayout  rLayout;
3121:   PetscInt    *dnz, *onz;
3122:   PetscInt     locRows, rStart, rEnd;
3123:   PetscReal   *x, *v0, *J, *invJ, detJ;
3124:   PetscReal   *v0c, *Jc, *invJc, detJc;
3125:   PetscScalar *elemMat;
3126:   PetscInt     dim, Nf, field, totDim, cStart, cEnd, cell, ccell;

3128:   PetscFunctionBegin;
3129:   PetscCall(DMGetCoordinateDim(dmc, &dim));
3130:   PetscCall(DMGetDS(dmc, &prob));
3131:   PetscCall(PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL));
3132:   PetscCall(PetscDSGetNumFields(prob, &Nf));
3133:   PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ));
3134:   PetscCall(PetscMalloc3(dim, &v0c, dim * dim, &Jc, dim * dim, &invJc));
3135:   PetscCall(DMGetLocalSection(dmf, &fsection));
3136:   PetscCall(DMGetGlobalSection(dmf, &globalFSection));
3137:   PetscCall(DMGetLocalSection(dmc, &csection));
3138:   PetscCall(DMGetGlobalSection(dmc, &globalCSection));
3139:   PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
3140:   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
3141:   PetscCall(PetscMalloc1(totDim, &elemMat));

3143:   PetscCall(MatGetLocalSize(mass, &locRows, NULL));
3144:   PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)mass), &rLayout));
3145:   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
3146:   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
3147:   PetscCall(PetscLayoutSetUp(rLayout));
3148:   PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
3149:   PetscCall(PetscLayoutDestroy(&rLayout));
3150:   PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
3151:   PetscCall(PetscHSetIJCreate(&ht));
3152:   for (field = 0; field < Nf; ++field) {
3153:     PetscObject      obj;
3154:     PetscClassId     id;
3155:     PetscQuadrature  quad;
3156:     const PetscReal *qpoints;
3157:     PetscInt         Nq, Nc, i, d;

3159:     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
3160:     PetscCall(PetscObjectGetClassId(obj, &id));
3161:     if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetQuadrature((PetscFE)obj, &quad));
3162:     else PetscCall(PetscFVGetQuadrature((PetscFV)obj, &quad));
3163:     PetscCall(PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, NULL));
3164:     /* For each fine grid cell */
3165:     for (cell = cStart; cell < cEnd; ++cell) {
3166:       Vec                pointVec;
3167:       PetscScalar       *pV;
3168:       PetscSF            coarseCellSF = NULL;
3169:       const PetscSFNode *coarseCells;
3170:       PetscInt           numCoarseCells, q, c;
3171:       PetscInt          *findices, *cindices;
3172:       PetscInt           numFIndices, numCIndices;

3174:       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3175:       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
3176:       /* Get points from the quadrature */
3177:       PetscCall(VecCreateSeq(PETSC_COMM_SELF, Nq * dim, &pointVec));
3178:       PetscCall(VecSetBlockSize(pointVec, dim));
3179:       PetscCall(VecGetArray(pointVec, &pV));
3180:       for (q = 0; q < Nq; ++q) {
3181:         const PetscReal xi0[3] = {-1., -1., -1.};

3183:         /* Transform point to real space */
3184:         CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q * dim], x);
3185:         for (d = 0; d < dim; ++d) pV[q * dim + d] = x[d];
3186:       }
3187:       PetscCall(VecRestoreArray(pointVec, &pV));
3188:       /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
3189:       PetscCall(DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF));
3190:       PetscCall(PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view"));
3191:       /* Update preallocation info */
3192:       PetscCall(PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells));
3193:       PetscCheck(numCoarseCells == Nq, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not all closure points located");
3194:       {
3195:         PetscHashIJKey key;
3196:         PetscBool      missing;

3198:         for (i = 0; i < numFIndices; ++i) {
3199:           key.i = findices[i];
3200:           if (key.i >= 0) {
3201:             /* Get indices for coarse elements */
3202:             for (ccell = 0; ccell < numCoarseCells; ++ccell) {
3203:               PetscCall(DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL));
3204:               for (c = 0; c < numCIndices; ++c) {
3205:                 key.j = cindices[c];
3206:                 if (key.j < 0) continue;
3207:                 PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
3208:                 if (missing) {
3209:                   if ((key.j >= rStart) && (key.j < rEnd)) ++dnz[key.i - rStart];
3210:                   else ++onz[key.i - rStart];
3211:                 }
3212:               }
3213:               PetscCall(DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL));
3214:             }
3215:           }
3216:         }
3217:       }
3218:       PetscCall(PetscSFDestroy(&coarseCellSF));
3219:       PetscCall(VecDestroy(&pointVec));
3220:       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3221:     }
3222:   }
3223:   PetscCall(PetscHSetIJDestroy(&ht));
3224:   PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
3225:   PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
3226:   PetscCall(PetscFree2(dnz, onz));
3227:   for (field = 0; field < Nf; ++field) {
3228:     PetscObject      obj;
3229:     PetscClassId     id;
3230:     PetscTabulation  T, Tfine;
3231:     PetscQuadrature  quad;
3232:     const PetscReal *qpoints, *qweights;
3233:     PetscInt         Nq, Nc, i, d;

3235:     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
3236:     PetscCall(PetscObjectGetClassId(obj, &id));
3237:     if (id == PETSCFE_CLASSID) {
3238:       PetscCall(PetscFEGetQuadrature((PetscFE)obj, &quad));
3239:       PetscCall(PetscFEGetCellTabulation((PetscFE)obj, 1, &Tfine));
3240:       PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, 1, x, 0, &T));
3241:     } else {
3242:       PetscCall(PetscFVGetQuadrature((PetscFV)obj, &quad));
3243:     }
3244:     PetscCall(PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, &qweights));
3245:     /* For each fine grid cell */
3246:     for (cell = cStart; cell < cEnd; ++cell) {
3247:       Vec                pointVec;
3248:       PetscScalar       *pV;
3249:       PetscSF            coarseCellSF = NULL;
3250:       const PetscSFNode *coarseCells;
3251:       PetscInt           numCoarseCells, cpdim, q, c, j;
3252:       PetscInt          *findices, *cindices;
3253:       PetscInt           numFIndices, numCIndices;

3255:       PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3256:       PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
3257:       /* Get points from the quadrature */
3258:       PetscCall(VecCreateSeq(PETSC_COMM_SELF, Nq * dim, &pointVec));
3259:       PetscCall(VecSetBlockSize(pointVec, dim));
3260:       PetscCall(VecGetArray(pointVec, &pV));
3261:       for (q = 0; q < Nq; ++q) {
3262:         const PetscReal xi0[3] = {-1., -1., -1.};

3264:         /* Transform point to real space */
3265:         CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q * dim], x);
3266:         for (d = 0; d < dim; ++d) pV[q * dim + d] = x[d];
3267:       }
3268:       PetscCall(VecRestoreArray(pointVec, &pV));
3269:       /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
3270:       PetscCall(DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF));
3271:       /* Update matrix */
3272:       PetscCall(PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells));
3273:       PetscCheck(numCoarseCells == Nq, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not all closure points located");
3274:       PetscCall(VecGetArray(pointVec, &pV));
3275:       for (ccell = 0; ccell < numCoarseCells; ++ccell) {
3276:         PetscReal       pVReal[3];
3277:         const PetscReal xi0[3] = {-1., -1., -1.};

3279:         PetscCall(DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL));
3280:         /* Transform points from real space to coarse reference space */
3281:         PetscCall(DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc));
3282:         for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell * dim + d]);
3283:         CoordinatesRealToRef(dim, dim, xi0, v0c, invJc, pVReal, x);

3285:         if (id == PETSCFE_CLASSID) {
3286:           PetscFE fe = (PetscFE)obj;

3288:           /* Evaluate coarse basis on contained point */
3289:           PetscCall(PetscFEGetDimension(fe, &cpdim));
3290:           PetscCall(PetscFEComputeTabulation(fe, 1, x, 0, T));
3291:           /* Get elemMat entries by multiplying by weight */
3292:           for (i = 0; i < numFIndices; ++i) {
3293:             PetscCall(PetscArrayzero(elemMat, cpdim));
3294:             for (j = 0; j < cpdim; ++j) {
3295:               for (c = 0; c < Nc; ++c) elemMat[j] += T->T[0][j * Nc + c] * Tfine->T[0][(ccell * numFIndices + i) * Nc + c] * qweights[ccell * Nc + c] * detJ;
3296:             }
3297:             /* Update interpolator */
3298:             if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
3299:             PetscCheck(numCIndices == cpdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %" PetscInt_FMT " != %" PetscInt_FMT, numCIndices, cpdim);
3300:             PetscCall(MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES));
3301:           }
3302:         } else {
3303:           cpdim = 1;
3304:           for (i = 0; i < numFIndices; ++i) {
3305:             PetscCall(PetscArrayzero(elemMat, cpdim));
3306:             for (j = 0; j < cpdim; ++j) {
3307:               for (c = 0; c < Nc; ++c) elemMat[j] += 1.0 * 1.0 * qweights[ccell * Nc + c] * detJ;
3308:             }
3309:             /* Update interpolator */
3310:             if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
3311:             PetscCall(PetscPrintf(PETSC_COMM_SELF, "Nq: %" PetscInt_FMT " %" PetscInt_FMT " Nf: %" PetscInt_FMT " %" PetscInt_FMT " Nc: %" PetscInt_FMT " %" PetscInt_FMT "\n", ccell, Nq, i, numFIndices, j, numCIndices));
3312:             PetscCheck(numCIndices == cpdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %" PetscInt_FMT " != %" PetscInt_FMT, numCIndices, cpdim);
3313:             PetscCall(MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES));
3314:           }
3315:         }
3316:         PetscCall(DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL));
3317:       }
3318:       PetscCall(VecRestoreArray(pointVec, &pV));
3319:       PetscCall(PetscSFDestroy(&coarseCellSF));
3320:       PetscCall(VecDestroy(&pointVec));
3321:       PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3322:     }
3323:     if (id == PETSCFE_CLASSID) PetscCall(PetscTabulationDestroy(&T));
3324:   }
3325:   PetscCall(PetscFree3(v0, J, invJ));
3326:   PetscCall(PetscFree3(v0c, Jc, invJc));
3327:   PetscCall(PetscFree(elemMat));
3328:   PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
3329:   PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
3330:   PetscFunctionReturn(PETSC_SUCCESS);
3331: }

3333: /*@
3334:   DMPlexComputeInjectorFEM - Compute a mapping from coarse unknowns to fine unknowns

3336:   Input Parameters:
3337: + dmc  - The coarse mesh
3338: . dmf  - The fine mesh
3339: - user - The user context

3341:   Output Parameter:
3342: . sc - The mapping

3344:   Level: developer

3346: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexComputeInterpolatorNested()`, `DMPlexComputeJacobianFEM()`
3347: @*/
3348: PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
3349: {
3350:   PetscDS      prob;
3351:   PetscFE     *feRef;
3352:   PetscFV     *fvRef;
3353:   Vec          fv, cv;
3354:   IS           fis, cis;
3355:   PetscSection fsection, fglobalSection, csection, cglobalSection;
3356:   PetscInt    *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
3357:   PetscInt     cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, c, dim, d, startC, endC, offsetC, offsetF, m;
3358:   PetscBool   *needAvg;

3360:   PetscFunctionBegin;
3361:   PetscCall(PetscLogEventBegin(DMPLEX_InjectorFEM, dmc, dmf, 0, 0));
3362:   PetscCall(DMGetDimension(dmf, &dim));
3363:   PetscCall(DMGetLocalSection(dmf, &fsection));
3364:   PetscCall(DMGetGlobalSection(dmf, &fglobalSection));
3365:   PetscCall(DMGetLocalSection(dmc, &csection));
3366:   PetscCall(DMGetGlobalSection(dmc, &cglobalSection));
3367:   PetscCall(PetscSectionGetNumFields(fsection, &Nf));
3368:   PetscCall(DMPlexGetSimplexOrBoxCells(dmc, 0, &cStart, &cEnd));
3369:   PetscCall(DMGetDS(dmc, &prob));
3370:   PetscCall(PetscCalloc3(Nf, &feRef, Nf, &fvRef, Nf, &needAvg));
3371:   for (f = 0; f < Nf; ++f) {
3372:     PetscObject  obj;
3373:     PetscClassId id;
3374:     PetscInt     fNb = 0, Nc = 0;

3376:     PetscCall(PetscDSGetDiscretization(prob, f, &obj));
3377:     PetscCall(PetscObjectGetClassId(obj, &id));
3378:     if (id == PETSCFE_CLASSID) {
3379:       PetscFE    fe = (PetscFE)obj;
3380:       PetscSpace sp;
3381:       PetscInt   maxDegree;

3383:       PetscCall(PetscFERefine(fe, &feRef[f]));
3384:       PetscCall(PetscFEGetDimension(feRef[f], &fNb));
3385:       PetscCall(PetscFEGetNumComponents(fe, &Nc));
3386:       PetscCall(PetscFEGetBasisSpace(fe, &sp));
3387:       PetscCall(PetscSpaceGetDegree(sp, NULL, &maxDegree));
3388:       if (!maxDegree) needAvg[f] = PETSC_TRUE;
3389:     } else if (id == PETSCFV_CLASSID) {
3390:       PetscFV        fv = (PetscFV)obj;
3391:       PetscDualSpace Q;

3393:       PetscCall(PetscFVRefine(fv, &fvRef[f]));
3394:       PetscCall(PetscFVGetDualSpace(fvRef[f], &Q));
3395:       PetscCall(PetscDualSpaceGetDimension(Q, &fNb));
3396:       PetscCall(PetscFVGetNumComponents(fv, &Nc));
3397:       needAvg[f] = PETSC_TRUE;
3398:     }
3399:     fTotDim += fNb;
3400:   }
3401:   PetscCall(PetscDSGetTotalDimension(prob, &cTotDim));
3402:   PetscCall(PetscMalloc1(cTotDim, &cmap));
3403:   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
3404:     PetscFE        feC;
3405:     PetscFV        fvC;
3406:     PetscDualSpace QF, QC;
3407:     PetscInt       order = -1, NcF, NcC, fpdim, cpdim;

3409:     if (feRef[field]) {
3410:       PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&feC));
3411:       PetscCall(PetscFEGetNumComponents(feC, &NcC));
3412:       PetscCall(PetscFEGetNumComponents(feRef[field], &NcF));
3413:       PetscCall(PetscFEGetDualSpace(feRef[field], &QF));
3414:       PetscCall(PetscDualSpaceGetOrder(QF, &order));
3415:       PetscCall(PetscDualSpaceGetDimension(QF, &fpdim));
3416:       PetscCall(PetscFEGetDualSpace(feC, &QC));
3417:       PetscCall(PetscDualSpaceGetDimension(QC, &cpdim));
3418:     } else {
3419:       PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fvC));
3420:       PetscCall(PetscFVGetNumComponents(fvC, &NcC));
3421:       PetscCall(PetscFVGetNumComponents(fvRef[field], &NcF));
3422:       PetscCall(PetscFVGetDualSpace(fvRef[field], &QF));
3423:       PetscCall(PetscDualSpaceGetDimension(QF, &fpdim));
3424:       PetscCall(PetscFVGetDualSpace(fvC, &QC));
3425:       PetscCall(PetscDualSpaceGetDimension(QC, &cpdim));
3426:     }
3427:     PetscCheck(NcF == NcC, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %" PetscInt_FMT " does not match coarse field %" PetscInt_FMT, NcF, NcC);
3428:     for (c = 0; c < cpdim; ++c) {
3429:       PetscQuadrature  cfunc;
3430:       const PetscReal *cqpoints, *cqweights;
3431:       PetscInt         NqcC, NpC;
3432:       PetscBool        found = PETSC_FALSE;

3434:       PetscCall(PetscDualSpaceGetFunctional(QC, c, &cfunc));
3435:       PetscCall(PetscQuadratureGetData(cfunc, NULL, &NqcC, &NpC, &cqpoints, &cqweights));
3436:       PetscCheck(NqcC == NcC, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of quadrature components %" PetscInt_FMT " must match number of field components %" PetscInt_FMT, NqcC, NcC);
3437:       PetscCheck(NpC == 1 || !feRef[field], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
3438:       for (f = 0; f < fpdim; ++f) {
3439:         PetscQuadrature  ffunc;
3440:         const PetscReal *fqpoints, *fqweights;
3441:         PetscReal        sum = 0.0;
3442:         PetscInt         NqcF, NpF;

3444:         PetscCall(PetscDualSpaceGetFunctional(QF, f, &ffunc));
3445:         PetscCall(PetscQuadratureGetData(ffunc, NULL, &NqcF, &NpF, &fqpoints, &fqweights));
3446:         PetscCheck(NqcF == NcF, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of quadrature components %" PetscInt_FMT " must match number of field components %" PetscInt_FMT, NqcF, NcF);
3447:         if (NpC != NpF) continue;
3448:         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
3449:         if (sum > 1.0e-9) continue;
3450:         for (d = 0; d < NcC; ++d) sum += PetscAbsReal(cqweights[d] * fqweights[d]);
3451:         if (sum < 1.0e-9) continue;
3452:         cmap[offsetC + c] = offsetF + f;
3453:         found             = PETSC_TRUE;
3454:         break;
3455:       }
3456:       if (!found) {
3457:         /* TODO We really want the average here, but some asshole put VecScatter in the interface */
3458:         if (fvRef[field] || (feRef[field] && order == 0)) {
3459:           cmap[offsetC + c] = offsetF + 0;
3460:         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection");
3461:       }
3462:     }
3463:     offsetC += cpdim;
3464:     offsetF += fpdim;
3465:   }
3466:   for (f = 0; f < Nf; ++f) {
3467:     PetscCall(PetscFEDestroy(&feRef[f]));
3468:     PetscCall(PetscFVDestroy(&fvRef[f]));
3469:   }
3470:   PetscCall(PetscFree3(feRef, fvRef, needAvg));

3472:   PetscCall(DMGetGlobalVector(dmf, &fv));
3473:   PetscCall(DMGetGlobalVector(dmc, &cv));
3474:   PetscCall(VecGetOwnershipRange(cv, &startC, &endC));
3475:   PetscCall(PetscSectionGetConstrainedStorageSize(cglobalSection, &m));
3476:   PetscCall(PetscMalloc2(cTotDim, &cellCIndices, fTotDim, &cellFIndices));
3477:   PetscCall(PetscMalloc1(m, &cindices));
3478:   PetscCall(PetscMalloc1(m, &findices));
3479:   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
3480:   for (c = cStart; c < cEnd; ++c) {
3481:     PetscCall(DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices));
3482:     for (d = 0; d < cTotDim; ++d) {
3483:       if ((cellCIndices[d] < startC) || (cellCIndices[d] >= endC)) continue;
3484:       PetscCheck(!(findices[cellCIndices[d] - startC] >= 0) || !(findices[cellCIndices[d] - startC] != cellFIndices[cmap[d]]), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %" PetscInt_FMT " Coarse dof %" PetscInt_FMT " maps to both %" PetscInt_FMT " and %" PetscInt_FMT, c, cindices[cellCIndices[d] - startC], findices[cellCIndices[d] - startC], cellFIndices[cmap[d]]);
3485:       cindices[cellCIndices[d] - startC] = cellCIndices[d];
3486:       findices[cellCIndices[d] - startC] = cellFIndices[cmap[d]];
3487:     }
3488:   }
3489:   PetscCall(PetscFree(cmap));
3490:   PetscCall(PetscFree2(cellCIndices, cellFIndices));

3492:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis));
3493:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis));
3494:   PetscCall(VecScatterCreate(cv, cis, fv, fis, sc));
3495:   PetscCall(ISDestroy(&cis));
3496:   PetscCall(ISDestroy(&fis));
3497:   PetscCall(DMRestoreGlobalVector(dmf, &fv));
3498:   PetscCall(DMRestoreGlobalVector(dmc, &cv));
3499:   PetscCall(PetscLogEventEnd(DMPLEX_InjectorFEM, dmc, dmf, 0, 0));
3500:   PetscFunctionReturn(PETSC_SUCCESS);
3501: }

3503: /*@C
3504:   DMPlexGetCellFields - Retrieve the field values values for a chunk of cells

3506:   Input Parameters:
3507: + dm     - The `DM`
3508: . cellIS - The cells to include
3509: . locX   - A local vector with the solution fields
3510: . locX_t - A local vector with solution field time derivatives, or NULL
3511: - locA   - A local vector with auxiliary fields, or NULL

3513:   Output Parameters:
3514: + u   - The field coefficients
3515: . u_t - The fields derivative coefficients
3516: - a   - The auxiliary field coefficients

3518:   Level: developer

3520: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetFaceFields()`
3521: @*/
3522: PetscErrorCode DMPlexGetCellFields(DM dm, IS cellIS, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a)
3523: {
3524:   DM              plex, plexA = NULL;
3525:   DMEnclosureType encAux;
3526:   PetscSection    section, sectionAux;
3527:   PetscDS         prob;
3528:   const PetscInt *cells;
3529:   PetscInt        cStart, cEnd, numCells, totDim, totDimAux, c;

3531:   PetscFunctionBegin;
3536:   PetscAssertPointer(u, 6);
3537:   PetscAssertPointer(u_t, 7);
3538:   PetscAssertPointer(a, 8);
3539:   PetscCall(DMPlexConvertPlex(dm, &plex, PETSC_FALSE));
3540:   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
3541:   PetscCall(DMGetLocalSection(dm, &section));
3542:   PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob, NULL));
3543:   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
3544:   if (locA) {
3545:     DM      dmAux;
3546:     PetscDS probAux;

3548:     PetscCall(VecGetDM(locA, &dmAux));
3549:     PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
3550:     PetscCall(DMPlexConvertPlex(dmAux, &plexA, PETSC_FALSE));
3551:     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
3552:     PetscCall(DMGetDS(dmAux, &probAux));
3553:     PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
3554:   }
3555:   numCells = cEnd - cStart;
3556:   PetscCall(DMGetWorkArray(dm, numCells * totDim, MPIU_SCALAR, u));
3557:   if (locX_t) PetscCall(DMGetWorkArray(dm, numCells * totDim, MPIU_SCALAR, u_t));
3558:   else *u_t = NULL;
3559:   if (locA) PetscCall(DMGetWorkArray(dm, numCells * totDimAux, MPIU_SCALAR, a));
3560:   else *a = NULL;
3561:   for (c = cStart; c < cEnd; ++c) {
3562:     const PetscInt cell = cells ? cells[c] : c;
3563:     const PetscInt cind = c - cStart;
3564:     PetscScalar   *x = NULL, *x_t = NULL, *ul = *u, *ul_t = *u_t, *al = *a;
3565:     PetscInt       i;

3567:     PetscCall(DMPlexVecGetClosure(plex, section, locX, cell, NULL, &x));
3568:     for (i = 0; i < totDim; ++i) ul[cind * totDim + i] = x[i];
3569:     PetscCall(DMPlexVecRestoreClosure(plex, section, locX, cell, NULL, &x));
3570:     if (locX_t) {
3571:       PetscCall(DMPlexVecGetClosure(plex, section, locX_t, cell, NULL, &x_t));
3572:       for (i = 0; i < totDim; ++i) ul_t[cind * totDim + i] = x_t[i];
3573:       PetscCall(DMPlexVecRestoreClosure(plex, section, locX_t, cell, NULL, &x_t));
3574:     }
3575:     if (locA) {
3576:       PetscInt subcell;
3577:       PetscCall(DMGetEnclosurePoint(plexA, dm, encAux, cell, &subcell));
3578:       PetscCall(DMPlexVecGetClosure(plexA, sectionAux, locA, subcell, NULL, &x));
3579:       for (i = 0; i < totDimAux; ++i) al[cind * totDimAux + i] = x[i];
3580:       PetscCall(DMPlexVecRestoreClosure(plexA, sectionAux, locA, subcell, NULL, &x));
3581:     }
3582:   }
3583:   PetscCall(DMDestroy(&plex));
3584:   if (locA) PetscCall(DMDestroy(&plexA));
3585:   PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
3586:   PetscFunctionReturn(PETSC_SUCCESS);
3587: }

3589: /*@C
3590:   DMPlexRestoreCellFields - Restore the field values values for a chunk of cells

3592:   Input Parameters:
3593: + dm     - The `DM`
3594: . cellIS - The cells to include
3595: . locX   - A local vector with the solution fields
3596: . locX_t - A local vector with solution field time derivatives, or NULL
3597: - locA   - A local vector with auxiliary fields, or NULL

3599:   Output Parameters:
3600: + u   - The field coefficients
3601: . u_t - The fields derivative coefficients
3602: - a   - The auxiliary field coefficients

3604:   Level: developer

3606: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetFaceFields()`
3607: @*/
3608: PetscErrorCode DMPlexRestoreCellFields(DM dm, IS cellIS, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a)
3609: {
3610:   PetscFunctionBegin;
3611:   PetscCall(DMRestoreWorkArray(dm, 0, MPIU_SCALAR, u));
3612:   if (locX_t) PetscCall(DMRestoreWorkArray(dm, 0, MPIU_SCALAR, u_t));
3613:   if (locA) PetscCall(DMRestoreWorkArray(dm, 0, MPIU_SCALAR, a));
3614:   PetscFunctionReturn(PETSC_SUCCESS);
3615: }

3617: static PetscErrorCode DMPlexGetHybridCellFields(DM dm, IS cellIS, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a)
3618: {
3619:   DM              plex, plexA = NULL;
3620:   DMEnclosureType encAux;
3621:   PetscSection    section, sectionAux;
3622:   PetscDS         ds, dsIn;
3623:   const PetscInt *cells;
3624:   PetscInt        cStart, cEnd, numCells, c, totDim, totDimAux, Nf, f;

3626:   PetscFunctionBegin;
3632:   PetscAssertPointer(u, 6);
3633:   PetscAssertPointer(u_t, 7);
3634:   PetscAssertPointer(a, 8);
3635:   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
3636:   numCells = cEnd - cStart;
3637:   PetscCall(DMPlexConvertPlex(dm, &plex, PETSC_FALSE));
3638:   PetscCall(DMGetLocalSection(dm, &section));
3639:   PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &ds, &dsIn));
3640:   PetscCall(PetscDSGetNumFields(dsIn, &Nf));
3641:   PetscCall(PetscDSGetTotalDimension(dsIn, &totDim));
3642:   if (locA) {
3643:     DM      dmAux;
3644:     PetscDS probAux;

3646:     PetscCall(VecGetDM(locA, &dmAux));
3647:     PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
3648:     PetscCall(DMPlexConvertPlex(dmAux, &plexA, PETSC_FALSE));
3649:     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
3650:     PetscCall(DMGetDS(dmAux, &probAux));
3651:     PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
3652:   }
3653:   PetscCall(DMGetWorkArray(dm, numCells * totDim, MPIU_SCALAR, u));
3654:   if (locX_t) PetscCall(DMGetWorkArray(dm, numCells * totDim, MPIU_SCALAR, u_t));
3655:   else {
3656:     *u_t = NULL;
3657:   }
3658:   if (locA) PetscCall(DMGetWorkArray(dm, numCells * totDimAux, MPIU_SCALAR, a));
3659:   else {
3660:     *a = NULL;
3661:   }
3662:   // Loop over cohesive cells
3663:   for (c = cStart; c < cEnd; ++c) {
3664:     const PetscInt  cell = cells ? cells[c] : c;
3665:     const PetscInt  cind = c - cStart;
3666:     PetscScalar    *xf = NULL, *xc = NULL, *x = NULL, *xf_t = NULL, *xc_t = NULL;
3667:     PetscScalar    *ul = &(*u)[cind * totDim], *ul_t = PetscSafePointerPlusOffset(*u_t, cind * totDim);
3668:     const PetscInt *cone, *ornt;
3669:     PetscInt        Nx = 0, Nxf, s;

3671:     PetscCall(DMPlexGetCone(dm, cell, &cone));
3672:     PetscCall(DMPlexGetConeOrientation(dm, cell, &ornt));
3673:     // Put in cohesive unknowns
3674:     PetscCall(DMPlexVecGetClosure(plex, section, locX, cell, &Nxf, &xf));
3675:     if (locX_t) PetscCall(DMPlexVecGetClosure(plex, section, locX_t, cell, NULL, &xf_t));
3676:     for (f = 0; f < Nf; ++f) {
3677:       PetscInt  fdofIn, foff, foffIn;
3678:       PetscBool cohesive;

3680:       PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
3681:       if (!cohesive) continue;
3682:       PetscCall(PetscDSGetFieldSize(dsIn, f, &fdofIn));
3683:       PetscCall(PetscDSGetFieldOffsetCohesive(ds, f, &foff));
3684:       PetscCall(PetscDSGetFieldOffsetCohesive(dsIn, f, &foffIn));
3685:       for (PetscInt i = 0; i < fdofIn; ++i) ul[foffIn + i] = xf[foff + i];
3686:       if (locX_t)
3687:         for (PetscInt i = 0; i < fdofIn; ++i) ul_t[foffIn + i] = xf_t[foff + i];
3688:       Nx += fdofIn;
3689:     }
3690:     PetscCall(DMPlexVecRestoreClosure(plex, section, locX, cell, &Nxf, &xf));
3691:     if (locX_t) PetscCall(DMPlexVecRestoreClosure(plex, section, locX_t, cell, NULL, &xf_t));
3692:     // Loop over sides of surface
3693:     for (s = 0; s < 2; ++s) {
3694:       const PetscInt *support;
3695:       const PetscInt  face = cone[s];
3696:       PetscInt        ssize, ncell, Nxc;

3698:       // I don't think I need the face to have 0 orientation in the hybrid cell
3699:       //PetscCheck(!ornt[s], PETSC_COMM_SELF, PETSC_ERR_SUP, "Face %" PetscInt_FMT " in hybrid cell %" PetscInt_FMT " has orientation %" PetscInt_FMT " != 0", face, cell, ornt[s]);
3700:       PetscCall(DMPlexGetSupport(dm, face, &support));
3701:       PetscCall(DMPlexGetSupportSize(dm, face, &ssize));
3702:       if (support[0] == cell) ncell = support[1];
3703:       else if (support[1] == cell) ncell = support[0];
3704:       else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " does not have cell %" PetscInt_FMT " in its support", face, cell);
3705:       // Get closure of both face and cell, stick in cell for normal fields and face for cohesive fields
3706:       PetscCall(DMPlexVecGetClosure(plex, section, locX, ncell, &Nxc, &xc));
3707:       if (locX_t) PetscCall(DMPlexVecGetClosure(plex, section, locX_t, ncell, NULL, &xc_t));
3708:       for (f = 0; f < Nf; ++f) {
3709:         PetscInt  fdofIn, foffIn;
3710:         PetscBool cohesive;

3712:         PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
3713:         if (cohesive) continue;
3714:         PetscCall(PetscDSGetFieldSize(dsIn, f, &fdofIn));
3715:         PetscCall(PetscDSGetFieldOffsetCohesive(dsIn, f, &foffIn));
3716:         for (PetscInt i = 0; i < fdofIn; ++i) ul[foffIn + s * fdofIn + i] = xc[foffIn + i];
3717:         if (locX_t)
3718:           for (PetscInt i = 0; i < fdofIn; ++i) ul_t[foffIn + s * fdofIn + i] = xc_t[foffIn + i];
3719:         Nx += fdofIn;
3720:       }
3721:       PetscCall(DMPlexVecRestoreClosure(plex, section, locX, ncell, &Nxc, &xc));
3722:       if (locX_t) PetscCall(DMPlexVecRestoreClosure(plex, section, locX_t, ncell, NULL, &xc_t));
3723:     }
3724:     PetscCheck(Nx == totDim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Closure size %" PetscInt_FMT " for cell %" PetscInt_FMT " does not match DS size %" PetscInt_FMT, Nx, cell, totDim);

3726:     if (locA) {
3727:       PetscScalar *al = &(*a)[cind * totDimAux];
3728:       PetscInt     subcell;

3730:       PetscCall(DMGetEnclosurePoint(plexA, dm, encAux, cell, &subcell));
3731:       PetscCall(DMPlexVecGetClosure(plexA, sectionAux, locA, subcell, &Nx, &x));
3732:       PetscCheck(Nx == totDimAux, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Closure size %" PetscInt_FMT " for subcell %" PetscInt_FMT "does not match DS size %" PetscInt_FMT, Nx, subcell, totDimAux);
3733:       for (PetscInt i = 0; i < totDimAux; ++i) al[i] = x[i];
3734:       PetscCall(DMPlexVecRestoreClosure(plexA, sectionAux, locA, subcell, &Nx, &x));
3735:     }
3736:   }
3737:   PetscCall(DMDestroy(&plex));
3738:   PetscCall(DMDestroy(&plexA));
3739:   PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
3740:   PetscFunctionReturn(PETSC_SUCCESS);
3741: }

3743: /*
3744:   DMPlexGetHybridFields - Get the field values for the negative side (s = 0) and positive side (s = 1) of the interface

3746:   Input Parameters:
3747: + dm      - The full domain DM
3748: . dmX     - An array of DM for the field, say an auxiliary DM, indexed by s
3749: . dsX     - An array of PetscDS for the field, indexed by s
3750: . cellIS  - The interface cells for which we want values
3751: . locX    - An array of local vectors with the field values, indexed by s
3752: - useCell - Flag to have values come from neighboring cell rather than endcap face

3754:   Output Parameter:
3755: . x       - An array of field values, indexed by s

3757:   Note:
3758:   The arrays in `x` will be allocated using `DMGetWorkArray()`, and must be returned using `DMPlexRestoreHybridFields()`.

3760:   Level: advanced

3762: .seealso: `DMPlexRestoreHybridFields()`, `DMGetWorkArray()`
3763: */
3764: static PetscErrorCode DMPlexGetHybridFields(DM dm, DM dmX[], PetscDS dsX[], IS cellIS, Vec locX[], PetscBool useCell, PetscScalar *x[])
3765: {
3766:   DM              plexX[2];
3767:   DMEnclosureType encX[2];
3768:   PetscSection    sectionX[2];
3769:   const PetscInt *cells;
3770:   PetscInt        cStart, cEnd, numCells, c, s, totDimX[2];

3772:   PetscFunctionBegin;
3773:   PetscAssertPointer(locX, 5);
3774:   if (!locX[0] || !locX[1]) PetscFunctionReturn(PETSC_SUCCESS);
3775:   PetscAssertPointer(dmX, 2);
3776:   PetscAssertPointer(dsX, 3);
3778:   PetscAssertPointer(x, 7);
3779:   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
3780:   numCells = cEnd - cStart;
3781:   for (s = 0; s < 2; ++s) {
3785:     PetscCall(DMPlexConvertPlex(dmX[s], &plexX[s], PETSC_FALSE));
3786:     PetscCall(DMGetEnclosureRelation(dmX[s], dm, &encX[s]));
3787:     PetscCall(DMGetLocalSection(dmX[s], &sectionX[s]));
3788:     PetscCall(PetscDSGetTotalDimension(dsX[s], &totDimX[s]));
3789:     PetscCall(DMGetWorkArray(dmX[s], numCells * totDimX[s], MPIU_SCALAR, &x[s]));
3790:   }
3791:   for (c = cStart; c < cEnd; ++c) {
3792:     const PetscInt  cell = cells ? cells[c] : c;
3793:     const PetscInt  cind = c - cStart;
3794:     const PetscInt *cone, *ornt;

3796:     PetscCall(DMPlexGetCone(dm, cell, &cone));
3797:     PetscCall(DMPlexGetConeOrientation(dm, cell, &ornt));
3798:     //PetscCheck(!ornt[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "Face %" PetscInt_FMT " in hybrid cell %" PetscInt_FMT " has orientation %" PetscInt_FMT " != 0", cone[0], cell, ornt[0]);
3799:     for (s = 0; s < 2; ++s) {
3800:       const PetscInt tdX     = totDimX[s];
3801:       PetscScalar   *closure = NULL, *xl = &x[s][cind * tdX];
3802:       PetscInt       face = cone[s], point = face, subpoint, Nx, i;

3804:       if (useCell) {
3805:         const PetscInt *support;
3806:         PetscInt        ssize;

3808:         PetscCall(DMPlexGetSupport(dm, face, &support));
3809:         PetscCall(DMPlexGetSupportSize(dm, face, &ssize));
3810:         PetscCheck(ssize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " from cell %" PetscInt_FMT " has support size %" PetscInt_FMT " != 2", face, cell, ssize);
3811:         if (support[0] == cell) point = support[1];
3812:         else if (support[1] == cell) point = support[0];
3813:         else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " does not have cell %" PetscInt_FMT " in its support", face, cell);
3814:       }
3815:       PetscCall(DMGetEnclosurePoint(plexX[s], dm, encX[s], point, &subpoint));
3816:       PetscCall(DMPlexVecGetOrientedClosure_Internal(plexX[s], sectionX[s], PETSC_FALSE, locX[s], subpoint, ornt[s], &Nx, &closure));
3817:       PetscCheck(Nx == tdX, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Closure size %" PetscInt_FMT " for subpoint %" PetscInt_FMT " does not match DS size %" PetscInt_FMT, Nx, subpoint, tdX);
3818:       for (i = 0; i < Nx; ++i) xl[i] = closure[i];
3819:       PetscCall(DMPlexVecRestoreClosure(plexX[s], sectionX[s], locX[s], subpoint, &Nx, &closure));
3820:     }
3821:   }
3822:   for (s = 0; s < 2; ++s) PetscCall(DMDestroy(&plexX[s]));
3823:   PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
3824:   PetscFunctionReturn(PETSC_SUCCESS);
3825: }

3827: static PetscErrorCode DMPlexRestoreHybridFields(DM dm, DM dmX[], PetscDS dsX[], IS cellIS, Vec locX[], PetscBool useCell, PetscScalar *x[])
3828: {
3829:   PetscFunctionBegin;
3830:   if (!locX[0] || !locX[1]) PetscFunctionReturn(PETSC_SUCCESS);
3831:   PetscCall(DMRestoreWorkArray(dmX[0], 0, MPIU_SCALAR, &x[0]));
3832:   PetscCall(DMRestoreWorkArray(dmX[1], 0, MPIU_SCALAR, &x[1]));
3833:   PetscFunctionReturn(PETSC_SUCCESS);
3834: }

3836: /*@C
3837:   DMPlexGetFaceFields - Retrieve the field values values for a chunk of faces

3839:   Input Parameters:
3840: + dm           - The `DM`
3841: . fStart       - The first face to include
3842: . fEnd         - The first face to exclude
3843: . locX         - A local vector with the solution fields
3844: . locX_t       - A local vector with solution field time derivatives, or NULL
3845: . faceGeometry - A local vector with face geometry
3846: . cellGeometry - A local vector with cell geometry
3847: - locGrad      - A local vector with field gradients, or NULL

3849:   Output Parameters:
3850: + Nface - The number of faces with field values
3851: . uL    - The field values at the left side of the face
3852: - uR    - The field values at the right side of the face

3854:   Level: developer

3856: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCellFields()`
3857: @*/
3858: PetscErrorCode DMPlexGetFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscInt *Nface, PetscScalar **uL, PetscScalar **uR)
3859: {
3860:   DM                 dmFace, dmCell, dmGrad = NULL;
3861:   PetscSection       section;
3862:   PetscDS            prob;
3863:   DMLabel            ghostLabel;
3864:   const PetscScalar *facegeom, *cellgeom, *x, *lgrad;
3865:   PetscBool         *isFE;
3866:   PetscInt           dim, Nf, f, Nc, numFaces = fEnd - fStart, iface, face;

3868:   PetscFunctionBegin;
3875:   PetscAssertPointer(uL, 10);
3876:   PetscAssertPointer(uR, 11);
3877:   PetscCall(DMGetDimension(dm, &dim));
3878:   PetscCall(DMGetDS(dm, &prob));
3879:   PetscCall(DMGetLocalSection(dm, &section));
3880:   PetscCall(PetscDSGetNumFields(prob, &Nf));
3881:   PetscCall(PetscDSGetTotalComponents(prob, &Nc));
3882:   PetscCall(PetscMalloc1(Nf, &isFE));
3883:   for (f = 0; f < Nf; ++f) {
3884:     PetscObject  obj;
3885:     PetscClassId id;

3887:     PetscCall(PetscDSGetDiscretization(prob, f, &obj));
3888:     PetscCall(PetscObjectGetClassId(obj, &id));
3889:     if (id == PETSCFE_CLASSID) {
3890:       isFE[f] = PETSC_TRUE;
3891:     } else if (id == PETSCFV_CLASSID) {
3892:       isFE[f] = PETSC_FALSE;
3893:     } else {
3894:       isFE[f] = PETSC_FALSE;
3895:     }
3896:   }
3897:   PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
3898:   PetscCall(VecGetArrayRead(locX, &x));
3899:   PetscCall(VecGetDM(faceGeometry, &dmFace));
3900:   PetscCall(VecGetArrayRead(faceGeometry, &facegeom));
3901:   PetscCall(VecGetDM(cellGeometry, &dmCell));
3902:   PetscCall(VecGetArrayRead(cellGeometry, &cellgeom));
3903:   if (locGrad) {
3904:     PetscCall(VecGetDM(locGrad, &dmGrad));
3905:     PetscCall(VecGetArrayRead(locGrad, &lgrad));
3906:   }
3907:   PetscCall(DMGetWorkArray(dm, numFaces * Nc, MPIU_SCALAR, uL));
3908:   PetscCall(DMGetWorkArray(dm, numFaces * Nc, MPIU_SCALAR, uR));
3909:   /* Right now just eat the extra work for FE (could make a cell loop) */
3910:   for (face = fStart, iface = 0; face < fEnd; ++face) {
3911:     const PetscInt  *cells;
3912:     PetscFVFaceGeom *fg;
3913:     PetscFVCellGeom *cgL, *cgR;
3914:     PetscScalar     *xL, *xR, *gL, *gR;
3915:     PetscScalar     *uLl = *uL, *uRl = *uR;
3916:     PetscInt         ghost, nsupp, nchild;

3918:     PetscCall(DMLabelGetValue(ghostLabel, face, &ghost));
3919:     PetscCall(DMPlexGetSupportSize(dm, face, &nsupp));
3920:     PetscCall(DMPlexGetTreeChildren(dm, face, &nchild, NULL));
3921:     if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
3922:     PetscCall(DMPlexPointLocalRead(dmFace, face, facegeom, &fg));
3923:     PetscCall(DMPlexGetSupport(dm, face, &cells));
3924:     PetscCall(DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL));
3925:     PetscCall(DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR));
3926:     for (f = 0; f < Nf; ++f) {
3927:       PetscInt off;

3929:       PetscCall(PetscDSGetComponentOffset(prob, f, &off));
3930:       if (isFE[f]) {
3931:         const PetscInt *cone;
3932:         PetscInt        comp, coneSizeL, coneSizeR, faceLocL, faceLocR, ldof, rdof, d;

3934:         xL = xR = NULL;
3935:         PetscCall(PetscSectionGetFieldComponents(section, f, &comp));
3936:         PetscCall(DMPlexVecGetClosure(dm, section, locX, cells[0], &ldof, (PetscScalar **)&xL));
3937:         PetscCall(DMPlexVecGetClosure(dm, section, locX, cells[1], &rdof, (PetscScalar **)&xR));
3938:         PetscCall(DMPlexGetCone(dm, cells[0], &cone));
3939:         PetscCall(DMPlexGetConeSize(dm, cells[0], &coneSizeL));
3940:         for (faceLocL = 0; faceLocL < coneSizeL; ++faceLocL)
3941:           if (cone[faceLocL] == face) break;
3942:         PetscCall(DMPlexGetCone(dm, cells[1], &cone));
3943:         PetscCall(DMPlexGetConeSize(dm, cells[1], &coneSizeR));
3944:         for (faceLocR = 0; faceLocR < coneSizeR; ++faceLocR)
3945:           if (cone[faceLocR] == face) break;
3946:         PetscCheck(faceLocL != coneSizeL || faceLocR != coneSizeR, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %" PetscInt_FMT " in cone of cell %" PetscInt_FMT " or cell %" PetscInt_FMT, face, cells[0], cells[1]);
3947:         /* Check that FEM field has values in the right cell (sometimes its an FV ghost cell) */
3948:         /* TODO: this is a hack that might not be right for nonconforming */
3949:         if (faceLocL < coneSizeL) {
3950:           PetscCall(PetscFEEvaluateFaceFields_Internal(prob, f, faceLocL, xL, &uLl[iface * Nc + off]));
3951:           if (rdof == ldof && faceLocR < coneSizeR) PetscCall(PetscFEEvaluateFaceFields_Internal(prob, f, faceLocR, xR, &uRl[iface * Nc + off]));
3952:           else {
3953:             for (d = 0; d < comp; ++d) uRl[iface * Nc + off + d] = uLl[iface * Nc + off + d];
3954:           }
3955:         } else {
3956:           PetscCall(PetscFEEvaluateFaceFields_Internal(prob, f, faceLocR, xR, &uRl[iface * Nc + off]));
3957:           PetscCall(PetscSectionGetFieldComponents(section, f, &comp));
3958:           for (d = 0; d < comp; ++d) uLl[iface * Nc + off + d] = uRl[iface * Nc + off + d];
3959:         }
3960:         PetscCall(DMPlexVecRestoreClosure(dm, section, locX, cells[0], &ldof, (PetscScalar **)&xL));
3961:         PetscCall(DMPlexVecRestoreClosure(dm, section, locX, cells[1], &rdof, (PetscScalar **)&xR));
3962:       } else {
3963:         PetscFV  fv;
3964:         PetscInt numComp, c;

3966:         PetscCall(PetscDSGetDiscretization(prob, f, (PetscObject *)&fv));
3967:         PetscCall(PetscFVGetNumComponents(fv, &numComp));
3968:         PetscCall(DMPlexPointLocalFieldRead(dm, cells[0], f, x, &xL));
3969:         PetscCall(DMPlexPointLocalFieldRead(dm, cells[1], f, x, &xR));
3970:         if (dmGrad) {
3971:           PetscReal dxL[3], dxR[3];

3973:           PetscCall(DMPlexPointLocalRead(dmGrad, cells[0], lgrad, &gL));
3974:           PetscCall(DMPlexPointLocalRead(dmGrad, cells[1], lgrad, &gR));
3975:           DMPlex_WaxpyD_Internal(dim, -1, cgL->centroid, fg->centroid, dxL);
3976:           DMPlex_WaxpyD_Internal(dim, -1, cgR->centroid, fg->centroid, dxR);
3977:           for (c = 0; c < numComp; ++c) {
3978:             uLl[iface * Nc + off + c] = xL[c] + DMPlex_DotD_Internal(dim, &gL[c * dim], dxL);
3979:             uRl[iface * Nc + off + c] = xR[c] + DMPlex_DotD_Internal(dim, &gR[c * dim], dxR);
3980:           }
3981:         } else {
3982:           for (c = 0; c < numComp; ++c) {
3983:             uLl[iface * Nc + off + c] = xL[c];
3984:             uRl[iface * Nc + off + c] = xR[c];
3985:           }
3986:         }
3987:       }
3988:     }
3989:     ++iface;
3990:   }
3991:   *Nface = iface;
3992:   PetscCall(VecRestoreArrayRead(locX, &x));
3993:   PetscCall(VecRestoreArrayRead(faceGeometry, &facegeom));
3994:   PetscCall(VecRestoreArrayRead(cellGeometry, &cellgeom));
3995:   if (locGrad) PetscCall(VecRestoreArrayRead(locGrad, &lgrad));
3996:   PetscCall(PetscFree(isFE));
3997:   PetscFunctionReturn(PETSC_SUCCESS);
3998: }

4000: /*@C
4001:   DMPlexRestoreFaceFields - Restore the field values values for a chunk of faces

4003:   Input Parameters:
4004: + dm           - The `DM`
4005: . fStart       - The first face to include
4006: . fEnd         - The first face to exclude
4007: . locX         - A local vector with the solution fields
4008: . locX_t       - A local vector with solution field time derivatives, or NULL
4009: . faceGeometry - A local vector with face geometry
4010: . cellGeometry - A local vector with cell geometry
4011: - locGrad      - A local vector with field gradients, or NULL

4013:   Output Parameters:
4014: + Nface - The number of faces with field values
4015: . uL    - The field values at the left side of the face
4016: - uR    - The field values at the right side of the face

4018:   Level: developer

4020: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetFaceFields()`
4021: @*/
4022: PetscErrorCode DMPlexRestoreFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscInt *Nface, PetscScalar **uL, PetscScalar **uR)
4023: {
4024:   PetscFunctionBegin;
4025:   PetscCall(DMRestoreWorkArray(dm, 0, MPIU_SCALAR, uL));
4026:   PetscCall(DMRestoreWorkArray(dm, 0, MPIU_SCALAR, uR));
4027:   PetscFunctionReturn(PETSC_SUCCESS);
4028: }

4030: /*@C
4031:   DMPlexGetFaceGeometry - Retrieve the geometric values for a chunk of faces

4033:   Input Parameters:
4034: + dm           - The `DM`
4035: . fStart       - The first face to include
4036: . fEnd         - The first face to exclude
4037: . faceGeometry - A local vector with face geometry
4038: - cellGeometry - A local vector with cell geometry

4040:   Output Parameters:
4041: + Nface - The number of faces with field values
4042: . fgeom - The extract the face centroid and normal
4043: - vol   - The cell volume

4045:   Level: developer

4047: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCellFields()`
4048: @*/
4049: PetscErrorCode DMPlexGetFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscInt *Nface, PetscFVFaceGeom **fgeom, PetscReal **vol)
4050: {
4051:   DM                 dmFace, dmCell;
4052:   DMLabel            ghostLabel;
4053:   const PetscScalar *facegeom, *cellgeom;
4054:   PetscInt           dim, numFaces = fEnd - fStart, iface, face;

4056:   PetscFunctionBegin;
4060:   PetscAssertPointer(fgeom, 7);
4061:   PetscAssertPointer(vol, 8);
4062:   PetscCall(DMGetDimension(dm, &dim));
4063:   PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
4064:   PetscCall(VecGetDM(faceGeometry, &dmFace));
4065:   PetscCall(VecGetArrayRead(faceGeometry, &facegeom));
4066:   PetscCall(VecGetDM(cellGeometry, &dmCell));
4067:   PetscCall(VecGetArrayRead(cellGeometry, &cellgeom));
4068:   PetscCall(PetscMalloc1(numFaces, fgeom));
4069:   PetscCall(DMGetWorkArray(dm, numFaces * 2, MPIU_SCALAR, vol));
4070:   for (face = fStart, iface = 0; face < fEnd; ++face) {
4071:     const PetscInt  *cells;
4072:     PetscFVFaceGeom *fg;
4073:     PetscFVCellGeom *cgL, *cgR;
4074:     PetscFVFaceGeom *fgeoml = *fgeom;
4075:     PetscReal       *voll   = *vol;
4076:     PetscInt         ghost, d, nchild, nsupp;

4078:     PetscCall(DMLabelGetValue(ghostLabel, face, &ghost));
4079:     PetscCall(DMPlexGetSupportSize(dm, face, &nsupp));
4080:     PetscCall(DMPlexGetTreeChildren(dm, face, &nchild, NULL));
4081:     if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
4082:     PetscCall(DMPlexPointLocalRead(dmFace, face, facegeom, &fg));
4083:     PetscCall(DMPlexGetSupport(dm, face, &cells));
4084:     PetscCall(DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL));
4085:     PetscCall(DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR));
4086:     for (d = 0; d < dim; ++d) {
4087:       fgeoml[iface].centroid[d] = fg->centroid[d];
4088:       fgeoml[iface].normal[d]   = fg->normal[d];
4089:     }
4090:     voll[iface * 2 + 0] = cgL->volume;
4091:     voll[iface * 2 + 1] = cgR->volume;
4092:     ++iface;
4093:   }
4094:   *Nface = iface;
4095:   PetscCall(VecRestoreArrayRead(faceGeometry, &facegeom));
4096:   PetscCall(VecRestoreArrayRead(cellGeometry, &cellgeom));
4097:   PetscFunctionReturn(PETSC_SUCCESS);
4098: }

4100: /*@C
4101:   DMPlexRestoreFaceGeometry - Restore the field values values for a chunk of faces

4103:   Input Parameters:
4104: + dm           - The `DM`
4105: . fStart       - The first face to include
4106: . fEnd         - The first face to exclude
4107: . faceGeometry - A local vector with face geometry
4108: - cellGeometry - A local vector with cell geometry

4110:   Output Parameters:
4111: + Nface - The number of faces with field values
4112: . fgeom - The extract the face centroid and normal
4113: - vol   - The cell volume

4115:   Level: developer

4117: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetFaceFields()`
4118: @*/
4119: PetscErrorCode DMPlexRestoreFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscInt *Nface, PetscFVFaceGeom **fgeom, PetscReal **vol)
4120: {
4121:   PetscFunctionBegin;
4122:   PetscCall(PetscFree(*fgeom));
4123:   PetscCall(DMRestoreWorkArray(dm, 0, MPIU_REAL, vol));
4124:   PetscFunctionReturn(PETSC_SUCCESS);
4125: }

4127: PetscErrorCode DMSNESGetFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
4128: {
4129:   char           composeStr[33] = {0};
4130:   PetscObjectId  id;
4131:   PetscContainer container;

4133:   PetscFunctionBegin;
4134:   PetscCall(PetscObjectGetId((PetscObject)quad, &id));
4135:   PetscCall(PetscSNPrintf(composeStr, 32, "DMSNESGetFEGeom_%" PetscInt64_FMT "\n", id));
4136:   PetscCall(PetscObjectQuery((PetscObject)pointIS, composeStr, (PetscObject *)&container));
4137:   if (container) {
4138:     PetscCall(PetscContainerGetPointer(container, (void **)geom));
4139:   } else {
4140:     PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, faceData, geom));
4141:     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
4142:     PetscCall(PetscContainerSetPointer(container, (void *)*geom));
4143:     PetscCall(PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom));
4144:     PetscCall(PetscObjectCompose((PetscObject)pointIS, composeStr, (PetscObject)container));
4145:     PetscCall(PetscContainerDestroy(&container));
4146:   }
4147:   PetscFunctionReturn(PETSC_SUCCESS);
4148: }

4150: PetscErrorCode DMSNESRestoreFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
4151: {
4152:   PetscFunctionBegin;
4153:   *geom = NULL;
4154:   PetscFunctionReturn(PETSC_SUCCESS);
4155: }

4157: PetscErrorCode DMPlexComputeResidual_Patch_Internal(DM dm, PetscSection section, IS cellIS, PetscReal t, Vec locX, Vec locX_t, Vec locF, void *user)
4158: {
4159:   DM_Plex        *mesh       = (DM_Plex *)dm->data;
4160:   const char     *name       = "Residual";
4161:   DM              dmAux      = NULL;
4162:   DMLabel         ghostLabel = NULL;
4163:   PetscDS         prob       = NULL;
4164:   PetscDS         probAux    = NULL;
4165:   PetscBool       useFEM     = PETSC_FALSE;
4166:   PetscBool       isImplicit = (locX_t || t == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE;
4167:   DMField         coordField = NULL;
4168:   Vec             locA;
4169:   PetscScalar    *u = NULL, *u_t, *a, *uL = NULL, *uR = NULL;
4170:   IS              chunkIS;
4171:   const PetscInt *cells;
4172:   PetscInt        cStart, cEnd, numCells;
4173:   PetscInt        Nf, f, totDim, totDimAux, numChunks, cellChunkSize, chunk, fStart, fEnd;
4174:   PetscInt        maxDegree = PETSC_MAX_INT;
4175:   PetscFormKey    key;
4176:   PetscQuadrature affineQuad = NULL, *quads = NULL;
4177:   PetscFEGeom    *affineGeom = NULL, **geoms = NULL;

4179:   PetscFunctionBegin;
4180:   PetscCall(PetscLogEventBegin(DMPLEX_ResidualFEM, dm, 0, 0, 0));
4181:   /* FEM+FVM */
4182:   /* 1: Get sizes from dm and dmAux */
4183:   PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
4184:   PetscCall(DMGetDS(dm, &prob));
4185:   PetscCall(PetscDSGetNumFields(prob, &Nf));
4186:   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
4187:   PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &locA));
4188:   if (locA) {
4189:     PetscCall(VecGetDM(locA, &dmAux));
4190:     PetscCall(DMGetDS(dmAux, &probAux));
4191:     PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
4192:   }
4193:   /* 2: Get geometric data */
4194:   for (f = 0; f < Nf; ++f) {
4195:     PetscObject  obj;
4196:     PetscClassId id;
4197:     PetscBool    fimp;

4199:     PetscCall(PetscDSGetImplicit(prob, f, &fimp));
4200:     if (isImplicit != fimp) continue;
4201:     PetscCall(PetscDSGetDiscretization(prob, f, &obj));
4202:     PetscCall(PetscObjectGetClassId(obj, &id));
4203:     if (id == PETSCFE_CLASSID) useFEM = PETSC_TRUE;
4204:     PetscCheck(id != PETSCFV_CLASSID, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Use of FVM with PCPATCH not yet implemented");
4205:   }
4206:   if (useFEM) {
4207:     PetscCall(DMGetCoordinateField(dm, &coordField));
4208:     PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
4209:     if (maxDegree <= 1) {
4210:       PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &affineQuad));
4211:       if (affineQuad) PetscCall(DMSNESGetFEGeom(coordField, cellIS, affineQuad, PETSC_FALSE, &affineGeom));
4212:     } else {
4213:       PetscCall(PetscCalloc2(Nf, &quads, Nf, &geoms));
4214:       for (f = 0; f < Nf; ++f) {
4215:         PetscObject  obj;
4216:         PetscClassId id;
4217:         PetscBool    fimp;

4219:         PetscCall(PetscDSGetImplicit(prob, f, &fimp));
4220:         if (isImplicit != fimp) continue;
4221:         PetscCall(PetscDSGetDiscretization(prob, f, &obj));
4222:         PetscCall(PetscObjectGetClassId(obj, &id));
4223:         if (id == PETSCFE_CLASSID) {
4224:           PetscFE fe = (PetscFE)obj;

4226:           PetscCall(PetscFEGetQuadrature(fe, &quads[f]));
4227:           PetscCall(PetscObjectReference((PetscObject)quads[f]));
4228:           PetscCall(DMSNESGetFEGeom(coordField, cellIS, quads[f], PETSC_FALSE, &geoms[f]));
4229:         }
4230:       }
4231:     }
4232:   }
4233:   /* Loop over chunks */
4234:   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
4235:   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
4236:   if (useFEM) PetscCall(ISCreate(PETSC_COMM_SELF, &chunkIS));
4237:   numCells      = cEnd - cStart;
4238:   numChunks     = 1;
4239:   cellChunkSize = numCells / numChunks;
4240:   numChunks     = PetscMin(1, numCells);
4241:   key.label     = NULL;
4242:   key.value     = 0;
4243:   key.part      = 0;
4244:   for (chunk = 0; chunk < numChunks; ++chunk) {
4245:     PetscScalar     *elemVec, *fluxL = NULL, *fluxR = NULL;
4246:     PetscReal       *vol   = NULL;
4247:     PetscFVFaceGeom *fgeom = NULL;
4248:     PetscInt         cS = cStart + chunk * cellChunkSize, cE = PetscMin(cS + cellChunkSize, cEnd), numCells = cE - cS, c;
4249:     PetscInt         numFaces = 0;

4251:     /* Extract field coefficients */
4252:     if (useFEM) {
4253:       PetscCall(ISGetPointSubrange(chunkIS, cS, cE, cells));
4254:       PetscCall(DMPlexGetCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a));
4255:       PetscCall(DMGetWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVec));
4256:       PetscCall(PetscArrayzero(elemVec, numCells * totDim));
4257:     }
4258:     /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */
4259:     /* Loop over fields */
4260:     for (f = 0; f < Nf; ++f) {
4261:       PetscObject  obj;
4262:       PetscClassId id;
4263:       PetscBool    fimp;
4264:       PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;

4266:       key.field = f;
4267:       PetscCall(PetscDSGetImplicit(prob, f, &fimp));
4268:       if (isImplicit != fimp) continue;
4269:       PetscCall(PetscDSGetDiscretization(prob, f, &obj));
4270:       PetscCall(PetscObjectGetClassId(obj, &id));
4271:       if (id == PETSCFE_CLASSID) {
4272:         PetscFE         fe        = (PetscFE)obj;
4273:         PetscFEGeom    *geom      = affineGeom ? affineGeom : geoms[f];
4274:         PetscFEGeom    *chunkGeom = NULL;
4275:         PetscQuadrature quad      = affineQuad ? affineQuad : quads[f];
4276:         PetscInt        Nq, Nb;

4278:         PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
4279:         PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
4280:         PetscCall(PetscFEGetDimension(fe, &Nb));
4281:         blockSize = Nb;
4282:         batchSize = numBlocks * blockSize;
4283:         PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
4284:         numChunks = numCells / (numBatches * batchSize);
4285:         Ne        = numChunks * numBatches * batchSize;
4286:         Nr        = numCells % (numBatches * batchSize);
4287:         offset    = numCells - Nr;
4288:         /* Integrate FE residual to get elemVec (need fields at quadrature points) */
4289:         /*   For FV, I think we use a P0 basis and the cell coefficients (for subdivided cells, we can tweak the basis tabulation to be the indicator function) */
4290:         PetscCall(PetscFEGeomGetChunk(geom, 0, offset, &chunkGeom));
4291:         PetscCall(PetscFEIntegrateResidual(prob, key, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec));
4292:         PetscCall(PetscFEGeomGetChunk(geom, offset, numCells, &chunkGeom));
4293:         PetscCall(PetscFEIntegrateResidual(prob, key, Nr, chunkGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, &a[offset * totDimAux], t, &elemVec[offset * totDim]));
4294:         PetscCall(PetscFEGeomRestoreChunk(geom, offset, numCells, &chunkGeom));
4295:       } else if (id == PETSCFV_CLASSID) {
4296:         PetscFV fv = (PetscFV)obj;

4298:         Ne = numFaces;
4299:         /* Riemann solve over faces (need fields at face centroids) */
4300:         /*   We need to evaluate FE fields at those coordinates */
4301:         PetscCall(PetscFVIntegrateRHSFunction(fv, prob, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR));
4302:       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
4303:     }
4304:     /* Loop over domain */
4305:     if (useFEM) {
4306:       /* Add elemVec to locX */
4307:       for (c = cS; c < cE; ++c) {
4308:         const PetscInt cell = cells ? cells[c] : c;
4309:         const PetscInt cind = c - cStart;

4311:         if (mesh->printFEM > 1) PetscCall(DMPrintCellVector(cell, name, totDim, &elemVec[cind * totDim]));
4312:         if (ghostLabel) {
4313:           PetscInt ghostVal;

4315:           PetscCall(DMLabelGetValue(ghostLabel, cell, &ghostVal));
4316:           if (ghostVal > 0) continue;
4317:         }
4318:         PetscCall(DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cind * totDim], ADD_ALL_VALUES));
4319:       }
4320:     }
4321:     /* Handle time derivative */
4322:     if (locX_t) {
4323:       PetscScalar *x_t, *fa;

4325:       PetscCall(VecGetArray(locF, &fa));
4326:       PetscCall(VecGetArray(locX_t, &x_t));
4327:       for (f = 0; f < Nf; ++f) {
4328:         PetscFV      fv;
4329:         PetscObject  obj;
4330:         PetscClassId id;
4331:         PetscInt     pdim, d;

4333:         PetscCall(PetscDSGetDiscretization(prob, f, &obj));
4334:         PetscCall(PetscObjectGetClassId(obj, &id));
4335:         if (id != PETSCFV_CLASSID) continue;
4336:         fv = (PetscFV)obj;
4337:         PetscCall(PetscFVGetNumComponents(fv, &pdim));
4338:         for (c = cS; c < cE; ++c) {
4339:           const PetscInt cell = cells ? cells[c] : c;
4340:           PetscScalar   *u_t, *r;

4342:           if (ghostLabel) {
4343:             PetscInt ghostVal;

4345:             PetscCall(DMLabelGetValue(ghostLabel, cell, &ghostVal));
4346:             if (ghostVal > 0) continue;
4347:           }
4348:           PetscCall(DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t));
4349:           PetscCall(DMPlexPointLocalFieldRef(dm, cell, f, fa, &r));
4350:           for (d = 0; d < pdim; ++d) r[d] += u_t[d];
4351:         }
4352:       }
4353:       PetscCall(VecRestoreArray(locX_t, &x_t));
4354:       PetscCall(VecRestoreArray(locF, &fa));
4355:     }
4356:     if (useFEM) {
4357:       PetscCall(DMPlexRestoreCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a));
4358:       PetscCall(DMRestoreWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVec));
4359:     }
4360:   }
4361:   if (useFEM) PetscCall(ISDestroy(&chunkIS));
4362:   PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
4363:   /* TODO Could include boundary residual here (see DMPlexComputeResidual_Internal) */
4364:   if (useFEM) {
4365:     if (maxDegree <= 1) {
4366:       PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, affineQuad, PETSC_FALSE, &affineGeom));
4367:       PetscCall(PetscQuadratureDestroy(&affineQuad));
4368:     } else {
4369:       for (f = 0; f < Nf; ++f) {
4370:         PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, quads[f], PETSC_FALSE, &geoms[f]));
4371:         PetscCall(PetscQuadratureDestroy(&quads[f]));
4372:       }
4373:       PetscCall(PetscFree2(quads, geoms));
4374:     }
4375:   }
4376:   PetscCall(PetscLogEventEnd(DMPLEX_ResidualFEM, dm, 0, 0, 0));
4377:   PetscFunctionReturn(PETSC_SUCCESS);
4378: }

4380: /*
4381:   We always assemble JacP, and if the matrix is different from Jac and two different sets of point functions are provided, we also assemble Jac

4383:   X   - The local solution vector
4384:   X_t - The local solution time derivative vector, or NULL
4385: */
4386: PetscErrorCode DMPlexComputeJacobian_Patch_Internal(DM dm, PetscSection section, PetscSection globalSection, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Mat Jac, Mat JacP, void *ctx)
4387: {
4388:   DM_Plex        *mesh = (DM_Plex *)dm->data;
4389:   const char     *name = "Jacobian", *nameP = "JacobianPre";
4390:   DM              dmAux = NULL;
4391:   PetscDS         prob, probAux = NULL;
4392:   PetscSection    sectionAux = NULL;
4393:   Vec             A;
4394:   DMField         coordField;
4395:   PetscFEGeom    *cgeomFEM;
4396:   PetscQuadrature qGeom = NULL;
4397:   Mat             J = Jac, JP = JacP;
4398:   PetscScalar    *work, *u = NULL, *u_t = NULL, *a = NULL, *elemMat = NULL, *elemMatP = NULL, *elemMatD = NULL;
4399:   PetscBool       hasJac, hasPrec, hasDyn, assembleJac, *isFE, hasFV = PETSC_FALSE;
4400:   const PetscInt *cells;
4401:   PetscFormKey    key;
4402:   PetscInt        Nf, fieldI, fieldJ, maxDegree, numCells, cStart, cEnd, numChunks, chunkSize, chunk, totDim, totDimAux = 0, sz, wsz, off = 0, offCell = 0;

4404:   PetscFunctionBegin;
4405:   PetscCall(ISGetLocalSize(cellIS, &numCells));
4406:   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
4407:   PetscCall(PetscLogEventBegin(DMPLEX_JacobianFEM, dm, 0, 0, 0));
4408:   PetscCall(DMGetDS(dm, &prob));
4409:   PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &A));
4410:   if (A) {
4411:     PetscCall(VecGetDM(A, &dmAux));
4412:     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
4413:     PetscCall(DMGetDS(dmAux, &probAux));
4414:   }
4415:   /* Get flags */
4416:   PetscCall(PetscDSGetNumFields(prob, &Nf));
4417:   PetscCall(DMGetWorkArray(dm, Nf, MPIU_BOOL, &isFE));
4418:   for (fieldI = 0; fieldI < Nf; ++fieldI) {
4419:     PetscObject  disc;
4420:     PetscClassId id;
4421:     PetscCall(PetscDSGetDiscretization(prob, fieldI, &disc));
4422:     PetscCall(PetscObjectGetClassId(disc, &id));
4423:     if (id == PETSCFE_CLASSID) {
4424:       isFE[fieldI] = PETSC_TRUE;
4425:     } else if (id == PETSCFV_CLASSID) {
4426:       hasFV        = PETSC_TRUE;
4427:       isFE[fieldI] = PETSC_FALSE;
4428:     }
4429:   }
4430:   PetscCall(PetscDSHasJacobian(prob, &hasJac));
4431:   PetscCall(PetscDSHasJacobianPreconditioner(prob, &hasPrec));
4432:   PetscCall(PetscDSHasDynamicJacobian(prob, &hasDyn));
4433:   assembleJac = hasJac && hasPrec && (Jac != JacP) ? PETSC_TRUE : PETSC_FALSE;
4434:   hasDyn      = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
4435:   if (hasFV) PetscCall(MatSetOption(JP, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); /* No allocated space for FV stuff, so ignore the zero entries */
4436:   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
4437:   if (probAux) PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
4438:   /* Compute batch sizes */
4439:   if (isFE[0]) {
4440:     PetscFE         fe;
4441:     PetscQuadrature q;
4442:     PetscInt        numQuadPoints, numBatches, batchSize, numBlocks, blockSize, Nb;

4444:     PetscCall(PetscDSGetDiscretization(prob, 0, (PetscObject *)&fe));
4445:     PetscCall(PetscFEGetQuadrature(fe, &q));
4446:     PetscCall(PetscQuadratureGetData(q, NULL, NULL, &numQuadPoints, NULL, NULL));
4447:     PetscCall(PetscFEGetDimension(fe, &Nb));
4448:     PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
4449:     blockSize = Nb * numQuadPoints;
4450:     batchSize = numBlocks * blockSize;
4451:     chunkSize = numBatches * batchSize;
4452:     numChunks = numCells / chunkSize + numCells % chunkSize;
4453:     PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
4454:   } else {
4455:     chunkSize = numCells;
4456:     numChunks = 1;
4457:   }
4458:   /* Get work space */
4459:   wsz = (((X ? 1 : 0) + (X_t ? 1 : 0) + (dmAux ? 1 : 0)) * totDim + ((hasJac ? 1 : 0) + (hasPrec ? 1 : 0) + (hasDyn ? 1 : 0)) * totDim * totDim) * chunkSize;
4460:   PetscCall(DMGetWorkArray(dm, wsz, MPIU_SCALAR, &work));
4461:   PetscCall(PetscArrayzero(work, wsz));
4462:   off      = 0;
4463:   u        = X ? (sz = chunkSize * totDim, off += sz, work + off - sz) : NULL;
4464:   u_t      = X_t ? (sz = chunkSize * totDim, off += sz, work + off - sz) : NULL;
4465:   a        = dmAux ? (sz = chunkSize * totDimAux, off += sz, work + off - sz) : NULL;
4466:   elemMat  = hasJac ? (sz = chunkSize * totDim * totDim, off += sz, work + off - sz) : NULL;
4467:   elemMatP = hasPrec ? (sz = chunkSize * totDim * totDim, off += sz, work + off - sz) : NULL;
4468:   elemMatD = hasDyn ? (sz = chunkSize * totDim * totDim, off += sz, work + off - sz) : NULL;
4469:   PetscCheck(off == wsz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error is workspace size %" PetscInt_FMT " should be %" PetscInt_FMT, off, wsz);
4470:   /* Setup geometry */
4471:   PetscCall(DMGetCoordinateField(dm, &coordField));
4472:   PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
4473:   if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &qGeom));
4474:   if (!qGeom) {
4475:     PetscFE fe;

4477:     PetscCall(PetscDSGetDiscretization(prob, 0, (PetscObject *)&fe));
4478:     PetscCall(PetscFEGetQuadrature(fe, &qGeom));
4479:     PetscCall(PetscObjectReference((PetscObject)qGeom));
4480:   }
4481:   PetscCall(DMSNESGetFEGeom(coordField, cellIS, qGeom, PETSC_FALSE, &cgeomFEM));
4482:   /* Compute volume integrals */
4483:   if (assembleJac) PetscCall(MatZeroEntries(J));
4484:   PetscCall(MatZeroEntries(JP));
4485:   key.label = NULL;
4486:   key.value = 0;
4487:   key.part  = 0;
4488:   for (chunk = 0; chunk < numChunks; ++chunk, offCell += chunkSize) {
4489:     const PetscInt Ncell = PetscMin(chunkSize, numCells - offCell);
4490:     PetscInt       c;

4492:     /* Extract values */
4493:     for (c = 0; c < Ncell; ++c) {
4494:       const PetscInt cell = cells ? cells[c + offCell] : c + offCell;
4495:       PetscScalar   *x = NULL, *x_t = NULL;
4496:       PetscInt       i;

4498:       if (X) {
4499:         PetscCall(DMPlexVecGetClosure(dm, section, X, cell, NULL, &x));
4500:         for (i = 0; i < totDim; ++i) u[c * totDim + i] = x[i];
4501:         PetscCall(DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x));
4502:       }
4503:       if (X_t) {
4504:         PetscCall(DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t));
4505:         for (i = 0; i < totDim; ++i) u_t[c * totDim + i] = x_t[i];
4506:         PetscCall(DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t));
4507:       }
4508:       if (dmAux) {
4509:         PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, A, cell, NULL, &x));
4510:         for (i = 0; i < totDimAux; ++i) a[c * totDimAux + i] = x[i];
4511:         PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, A, cell, NULL, &x));
4512:       }
4513:     }
4514:     for (fieldI = 0; fieldI < Nf; ++fieldI) {
4515:       PetscFE fe;
4516:       PetscCall(PetscDSGetDiscretization(prob, fieldI, (PetscObject *)&fe));
4517:       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
4518:         key.field = fieldI * Nf + fieldJ;
4519:         if (hasJac) PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, key, Ncell, cgeomFEM, u, u_t, probAux, a, t, X_tShift, elemMat));
4520:         if (hasPrec) PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_PRE, key, Ncell, cgeomFEM, u, u_t, probAux, a, t, X_tShift, elemMatP));
4521:         if (hasDyn) PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, key, Ncell, cgeomFEM, u, u_t, probAux, a, t, X_tShift, elemMatD));
4522:       }
4523:       /* For finite volume, add the identity */
4524:       if (!isFE[fieldI]) {
4525:         PetscFV  fv;
4526:         PetscInt eOffset = 0, Nc, fc, foff;

4528:         PetscCall(PetscDSGetFieldOffset(prob, fieldI, &foff));
4529:         PetscCall(PetscDSGetDiscretization(prob, fieldI, (PetscObject *)&fv));
4530:         PetscCall(PetscFVGetNumComponents(fv, &Nc));
4531:         for (c = 0; c < chunkSize; ++c, eOffset += totDim * totDim) {
4532:           for (fc = 0; fc < Nc; ++fc) {
4533:             const PetscInt i = foff + fc;
4534:             if (hasJac) elemMat[eOffset + i * totDim + i] = 1.0;
4535:             if (hasPrec) elemMatP[eOffset + i * totDim + i] = 1.0;
4536:           }
4537:         }
4538:       }
4539:     }
4540:     /*   Add contribution from X_t */
4541:     if (hasDyn) {
4542:       for (c = 0; c < chunkSize * totDim * totDim; ++c) elemMat[c] += X_tShift * elemMatD[c];
4543:     }
4544:     /* Insert values into matrix */
4545:     for (c = 0; c < Ncell; ++c) {
4546:       const PetscInt cell = cells ? cells[c + offCell] : c + offCell;
4547:       if (mesh->printFEM > 1) {
4548:         if (hasJac) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[(c - cStart) * totDim * totDim]));
4549:         if (hasPrec) PetscCall(DMPrintCellMatrix(cell, nameP, totDim, totDim, &elemMatP[(c - cStart) * totDim * totDim]));
4550:       }
4551:       if (assembleJac) PetscCall(DMPlexMatSetClosure_Internal(dm, section, globalSection, mesh->useMatClPerm, Jac, cell, &elemMat[(c - cStart) * totDim * totDim], ADD_VALUES));
4552:       PetscCall(DMPlexMatSetClosure_Internal(dm, section, globalSection, mesh->useMatClPerm, JP, cell, &elemMat[(c - cStart) * totDim * totDim], ADD_VALUES));
4553:     }
4554:   }
4555:   /* Cleanup */
4556:   PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, qGeom, PETSC_FALSE, &cgeomFEM));
4557:   PetscCall(PetscQuadratureDestroy(&qGeom));
4558:   if (hasFV) PetscCall(MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE));
4559:   PetscCall(DMRestoreWorkArray(dm, Nf, MPIU_BOOL, &isFE));
4560:   PetscCall(DMRestoreWorkArray(dm, ((1 + (X_t ? 1 : 0) + (dmAux ? 1 : 0)) * totDim + ((hasJac ? 1 : 0) + (hasPrec ? 1 : 0) + (hasDyn ? 1 : 0)) * totDim * totDim) * chunkSize, MPIU_SCALAR, &work));
4561:   /* Compute boundary integrals */
4562:   /* PetscCall(DMPlexComputeBdJacobian_Internal(dm, X, X_t, t, X_tShift, Jac, JacP, ctx)); */
4563:   /* Assemble matrix */
4564:   if (assembleJac) {
4565:     PetscCall(MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY));
4566:     PetscCall(MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY));
4567:   }
4568:   PetscCall(MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY));
4569:   PetscCall(MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY));
4570:   PetscCall(PetscLogEventEnd(DMPLEX_JacobianFEM, dm, 0, 0, 0));
4571:   PetscFunctionReturn(PETSC_SUCCESS);
4572: }

4574: /******** FEM Assembly Function ********/

4576: static PetscErrorCode DMConvertPlex_Internal(DM dm, DM *plex, PetscBool copy)
4577: {
4578:   PetscBool isPlex;

4580:   PetscFunctionBegin;
4581:   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
4582:   if (isPlex) {
4583:     *plex = dm;
4584:     PetscCall(PetscObjectReference((PetscObject)dm));
4585:   } else {
4586:     PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex));
4587:     if (!*plex) {
4588:       PetscCall(DMConvert(dm, DMPLEX, plex));
4589:       PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex));
4590:     } else {
4591:       PetscCall(PetscObjectReference((PetscObject)*plex));
4592:     }
4593:     if (copy) PetscCall(DMCopyAuxiliaryVec(dm, *plex));
4594:   }
4595:   PetscFunctionReturn(PETSC_SUCCESS);
4596: }

4598: /*@
4599:   DMPlexGetGeometryFVM - Return precomputed geometric data

4601:   Collective

4603:   Input Parameter:
4604: . dm - The `DM`

4606:   Output Parameters:
4607: + facegeom  - The values precomputed from face geometry
4608: . cellgeom  - The values precomputed from cell geometry
4609: - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell

4611:   Level: developer

4613: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMTSSetRHSFunctionLocal()`
4614: @*/
4615: PetscErrorCode DMPlexGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
4616: {
4617:   DM plex;

4619:   PetscFunctionBegin;
4621:   PetscCall(DMConvertPlex_Internal(dm, &plex, PETSC_TRUE));
4622:   PetscCall(DMPlexGetDataFVM(plex, NULL, cellgeom, facegeom, NULL));
4623:   if (minRadius) PetscCall(DMPlexGetMinRadius(plex, minRadius));
4624:   PetscCall(DMDestroy(&plex));
4625:   PetscFunctionReturn(PETSC_SUCCESS);
4626: }

4628: /*@
4629:   DMPlexGetGradientDM - Return gradient data layout

4631:   Collective

4633:   Input Parameters:
4634: + dm - The `DM`
4635: - fv - The `PetscFV`

4637:   Output Parameter:
4638: . dmGrad - The layout for gradient values

4640:   Level: developer

4642: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetGeometryFVM()`
4643: @*/
4644: PetscErrorCode DMPlexGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
4645: {
4646:   DM        plex;
4647:   PetscBool computeGradients;

4649:   PetscFunctionBegin;
4652:   PetscAssertPointer(dmGrad, 3);
4653:   PetscCall(PetscFVGetComputeGradients(fv, &computeGradients));
4654:   if (!computeGradients) {
4655:     *dmGrad = NULL;
4656:     PetscFunctionReturn(PETSC_SUCCESS);
4657:   }
4658:   PetscCall(DMConvertPlex_Internal(dm, &plex, PETSC_TRUE));
4659:   PetscCall(DMPlexGetDataFVM(plex, fv, NULL, NULL, dmGrad));
4660:   PetscCall(DMDestroy(&plex));
4661:   PetscFunctionReturn(PETSC_SUCCESS);
4662: }

4664: static PetscErrorCode DMPlexComputeBdResidual_Single_Internal(DM dm, PetscReal t, PetscWeakForm wf, PetscFormKey key, Vec locX, Vec locX_t, Vec locF, DMField coordField, IS facetIS)
4665: {
4666:   DM_Plex        *mesh = (DM_Plex *)dm->data;
4667:   DM              plex = NULL, plexA = NULL;
4668:   const char     *name = "BdResidual";
4669:   DMEnclosureType encAux;
4670:   PetscDS         prob, probAux       = NULL;
4671:   PetscSection    section, sectionAux = NULL;
4672:   Vec             locA = NULL;
4673:   PetscScalar    *u = NULL, *u_t = NULL, *a = NULL, *elemVec = NULL;
4674:   PetscInt        totDim, totDimAux = 0;

4676:   PetscFunctionBegin;
4677:   PetscCall(DMConvert(dm, DMPLEX, &plex));
4678:   PetscCall(DMGetLocalSection(dm, &section));
4679:   PetscCall(DMGetDS(dm, &prob));
4680:   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
4681:   PetscCall(DMGetAuxiliaryVec(dm, key.label, key.value, key.part, &locA));
4682:   if (locA) {
4683:     DM dmAux;

4685:     PetscCall(VecGetDM(locA, &dmAux));
4686:     PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
4687:     PetscCall(DMConvert(dmAux, DMPLEX, &plexA));
4688:     PetscCall(DMGetDS(plexA, &probAux));
4689:     PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
4690:     PetscCall(DMGetLocalSection(plexA, &sectionAux));
4691:   }
4692:   {
4693:     PetscFEGeom    *fgeom;
4694:     PetscInt        maxDegree;
4695:     PetscQuadrature qGeom = NULL;
4696:     IS              pointIS;
4697:     const PetscInt *points;
4698:     PetscInt        numFaces, face, Nq;

4700:     PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
4701:     if (!pointIS) goto end; /* No points with that id on this process */
4702:     {
4703:       IS isectIS;

4705:       /* TODO: Special cases of ISIntersect where it is quick to check a priori if one is a superset of the other */
4706:       PetscCall(ISIntersect_Caching_Internal(facetIS, pointIS, &isectIS));
4707:       PetscCall(ISDestroy(&pointIS));
4708:       pointIS = isectIS;
4709:     }
4710:     PetscCall(ISGetLocalSize(pointIS, &numFaces));
4711:     PetscCall(ISGetIndices(pointIS, &points));
4712:     PetscCall(PetscMalloc4(numFaces * totDim, &u, locX_t ? numFaces * totDim : 0, &u_t, numFaces * totDim, &elemVec, locA ? numFaces * totDimAux : 0, &a));
4713:     PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
4714:     if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &qGeom));
4715:     if (!qGeom) {
4716:       PetscFE fe;

4718:       PetscCall(PetscDSGetDiscretization(prob, key.field, (PetscObject *)&fe));
4719:       PetscCall(PetscFEGetFaceQuadrature(fe, &qGeom));
4720:       PetscCall(PetscObjectReference((PetscObject)qGeom));
4721:     }
4722:     PetscCall(PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL));
4723:     PetscCall(DMSNESGetFEGeom(coordField, pointIS, qGeom, PETSC_TRUE, &fgeom));
4724:     for (face = 0; face < numFaces; ++face) {
4725:       const PetscInt point = points[face], *support;
4726:       PetscScalar   *x     = NULL;
4727:       PetscInt       i;

4729:       PetscCall(DMPlexGetSupport(dm, point, &support));
4730:       PetscCall(DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x));
4731:       for (i = 0; i < totDim; ++i) u[face * totDim + i] = x[i];
4732:       PetscCall(DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x));
4733:       if (locX_t) {
4734:         PetscCall(DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x));
4735:         for (i = 0; i < totDim; ++i) u_t[face * totDim + i] = x[i];
4736:         PetscCall(DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x));
4737:       }
4738:       if (locA) {
4739:         PetscInt subp;

4741:         PetscCall(DMGetEnclosurePoint(plexA, dm, encAux, support[0], &subp));
4742:         PetscCall(DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x));
4743:         for (i = 0; i < totDimAux; ++i) a[face * totDimAux + i] = x[i];
4744:         PetscCall(DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x));
4745:       }
4746:     }
4747:     PetscCall(PetscArrayzero(elemVec, numFaces * totDim));
4748:     {
4749:       PetscFE      fe;
4750:       PetscInt     Nb;
4751:       PetscFEGeom *chunkGeom = NULL;
4752:       /* Conforming batches */
4753:       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
4754:       /* Remainder */
4755:       PetscInt Nr, offset;

4757:       PetscCall(PetscDSGetDiscretization(prob, key.field, (PetscObject *)&fe));
4758:       PetscCall(PetscFEGetDimension(fe, &Nb));
4759:       PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
4760:       /* TODO: documentation is unclear about what is going on with these numbers: how should Nb / Nq factor in ? */
4761:       blockSize = Nb;
4762:       batchSize = numBlocks * blockSize;
4763:       PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
4764:       numChunks = numFaces / (numBatches * batchSize);
4765:       Ne        = numChunks * numBatches * batchSize;
4766:       Nr        = numFaces % (numBatches * batchSize);
4767:       offset    = numFaces - Nr;
4768:       PetscCall(PetscFEGeomGetChunk(fgeom, 0, offset, &chunkGeom));
4769:       PetscCall(PetscFEIntegrateBdResidual(prob, wf, key, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec));
4770:       PetscCall(PetscFEGeomRestoreChunk(fgeom, 0, offset, &chunkGeom));
4771:       PetscCall(PetscFEGeomGetChunk(fgeom, offset, numFaces, &chunkGeom));
4772:       PetscCall(PetscFEIntegrateBdResidual(prob, wf, key, Nr, chunkGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, &elemVec[offset * totDim]));
4773:       PetscCall(PetscFEGeomRestoreChunk(fgeom, offset, numFaces, &chunkGeom));
4774:     }
4775:     for (face = 0; face < numFaces; ++face) {
4776:       const PetscInt point = points[face], *support;

4778:       if (mesh->printFEM > 1) PetscCall(DMPrintCellVector(point, name, totDim, &elemVec[face * totDim]));
4779:       PetscCall(DMPlexGetSupport(plex, point, &support));
4780:       PetscCall(DMPlexVecSetClosure(plex, NULL, locF, support[0], &elemVec[face * totDim], ADD_ALL_VALUES));
4781:     }
4782:     PetscCall(DMSNESRestoreFEGeom(coordField, pointIS, qGeom, PETSC_TRUE, &fgeom));
4783:     PetscCall(PetscQuadratureDestroy(&qGeom));
4784:     PetscCall(ISRestoreIndices(pointIS, &points));
4785:     PetscCall(ISDestroy(&pointIS));
4786:     PetscCall(PetscFree4(u, u_t, elemVec, a));
4787:   }
4788: end:
4789:   if (mesh->printFEM) {
4790:     PetscSection s;
4791:     Vec          locFbc;
4792:     PetscInt     pStart, pEnd, maxDof;
4793:     PetscScalar *zeroes;

4795:     PetscCall(DMGetLocalSection(dm, &s));
4796:     PetscCall(VecDuplicate(locF, &locFbc));
4797:     PetscCall(VecCopy(locF, locFbc));
4798:     PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
4799:     PetscCall(PetscSectionGetMaxDof(s, &maxDof));
4800:     PetscCall(PetscCalloc1(maxDof, &zeroes));
4801:     for (PetscInt p = pStart; p < pEnd; p++) PetscCall(VecSetValuesSection(locFbc, s, p, zeroes, INSERT_BC_VALUES));
4802:     PetscCall(PetscFree(zeroes));
4803:     PetscCall(DMPrintLocalVec(dm, name, mesh->printTol, locFbc));
4804:     PetscCall(VecDestroy(&locFbc));
4805:   }
4806:   PetscCall(DMDestroy(&plex));
4807:   PetscCall(DMDestroy(&plexA));
4808:   PetscFunctionReturn(PETSC_SUCCESS);
4809: }

4811: PetscErrorCode DMPlexComputeBdResidualSingle(DM dm, PetscReal t, PetscWeakForm wf, PetscFormKey key, Vec locX, Vec locX_t, Vec locF)
4812: {
4813:   DMField  coordField;
4814:   DMLabel  depthLabel;
4815:   IS       facetIS;
4816:   PetscInt dim;

4818:   PetscFunctionBegin;
4819:   PetscCall(DMGetDimension(dm, &dim));
4820:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
4821:   PetscCall(DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS));
4822:   PetscCall(DMGetCoordinateField(dm, &coordField));
4823:   PetscCall(DMPlexComputeBdResidual_Single_Internal(dm, t, wf, key, locX, locX_t, locF, coordField, facetIS));
4824:   PetscCall(ISDestroy(&facetIS));
4825:   PetscFunctionReturn(PETSC_SUCCESS);
4826: }

4828: static PetscErrorCode DMPlexComputeBdResidual_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
4829: {
4830:   PetscDS  prob;
4831:   PetscInt numBd, bd;
4832:   DMField  coordField = NULL;
4833:   IS       facetIS    = NULL;
4834:   DMLabel  depthLabel;
4835:   PetscInt dim;

4837:   PetscFunctionBegin;
4838:   PetscCall(DMGetDS(dm, &prob));
4839:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
4840:   PetscCall(DMGetDimension(dm, &dim));
4841:   PetscCall(DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS));
4842:   PetscCall(PetscDSGetNumBoundary(prob, &numBd));
4843:   for (bd = 0; bd < numBd; ++bd) {
4844:     PetscWeakForm           wf;
4845:     DMBoundaryConditionType type;
4846:     DMLabel                 label;
4847:     const PetscInt         *values;
4848:     PetscInt                field, numValues, v;
4849:     PetscObject             obj;
4850:     PetscClassId            id;
4851:     PetscFormKey            key;

4853:     PetscCall(PetscDSGetBoundary(prob, bd, &wf, &type, NULL, &label, &numValues, &values, &field, NULL, NULL, NULL, NULL, NULL));
4854:     if (type & DM_BC_ESSENTIAL) continue;
4855:     PetscCall(PetscDSGetDiscretization(prob, field, &obj));
4856:     PetscCall(PetscObjectGetClassId(obj, &id));
4857:     if (id != PETSCFE_CLASSID) continue;
4858:     if (!facetIS) {
4859:       DMLabel  depthLabel;
4860:       PetscInt dim;

4862:       PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
4863:       PetscCall(DMGetDimension(dm, &dim));
4864:       PetscCall(DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS));
4865:     }
4866:     PetscCall(DMGetCoordinateField(dm, &coordField));
4867:     for (v = 0; v < numValues; ++v) {
4868:       key.label = label;
4869:       key.value = values[v];
4870:       key.field = field;
4871:       key.part  = 0;
4872:       PetscCall(DMPlexComputeBdResidual_Single_Internal(dm, t, wf, key, locX, locX_t, locF, coordField, facetIS));
4873:     }
4874:   }
4875:   PetscCall(ISDestroy(&facetIS));
4876:   PetscFunctionReturn(PETSC_SUCCESS);
4877: }

4879: PetscErrorCode DMPlexComputeResidual_Internal(DM dm, PetscFormKey key, IS cellIS, PetscReal time, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
4880: {
4881:   DM_Plex        *mesh       = (DM_Plex *)dm->data;
4882:   const char     *name       = "Residual";
4883:   DM              dmAux      = NULL;
4884:   DM              dmGrad     = NULL;
4885:   DMLabel         ghostLabel = NULL;
4886:   PetscDS         ds         = NULL;
4887:   PetscDS         dsAux      = NULL;
4888:   PetscSection    section    = NULL;
4889:   PetscBool       useFEM     = PETSC_FALSE;
4890:   PetscBool       useFVM     = PETSC_FALSE;
4891:   PetscBool       isImplicit = (locX_t || time == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE;
4892:   PetscFV         fvm        = NULL;
4893:   DMField         coordField = NULL;
4894:   Vec             locA, cellGeometryFVM = NULL, faceGeometryFVM = NULL, locGrad = NULL;
4895:   PetscScalar    *u = NULL, *u_t, *a, *uL, *uR;
4896:   IS              chunkIS;
4897:   const PetscInt *cells;
4898:   PetscInt        cStart, cEnd, numCells;
4899:   PetscInt        Nf, f, totDim, totDimAux, numChunks, cellChunkSize, faceChunkSize, chunk, fStart, fEnd;
4900:   PetscInt        maxDegree  = PETSC_MAX_INT;
4901:   PetscQuadrature affineQuad = NULL, *quads = NULL;
4902:   PetscFEGeom    *affineGeom = NULL, **geoms = NULL;

4904:   PetscFunctionBegin;
4905:   PetscCall(PetscLogEventBegin(DMPLEX_ResidualFEM, dm, 0, 0, 0));
4906:   if (!cellIS) goto end;
4907:   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
4908:   if (cStart >= cEnd) goto end;
4909:   /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */
4910:   /* TODO The FVM geometry is over-manipulated. Make the precalc functions return exactly what we need */
4911:   /* FEM+FVM */
4912:   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
4913:   /* 1: Get sizes from dm and dmAux */
4914:   PetscCall(DMGetLocalSection(dm, &section));
4915:   PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
4916:   PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &ds, NULL));
4917:   PetscCall(PetscDSGetNumFields(ds, &Nf));
4918:   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
4919:   PetscCall(DMGetAuxiliaryVec(dm, key.label, key.value, key.part, &locA));
4920:   if (locA) {
4921:     PetscInt subcell;
4922:     PetscCall(VecGetDM(locA, &dmAux));
4923:     PetscCall(DMGetEnclosurePoint(dmAux, dm, DM_ENC_UNKNOWN, cells ? cells[cStart] : cStart, &subcell));
4924:     PetscCall(DMGetCellDS(dmAux, subcell, &dsAux, NULL));
4925:     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
4926:   }
4927:   /* 2: Get geometric data */
4928:   for (f = 0; f < Nf; ++f) {
4929:     PetscObject  obj;
4930:     PetscClassId id;
4931:     PetscBool    fimp;

4933:     PetscCall(PetscDSGetImplicit(ds, f, &fimp));
4934:     if (isImplicit != fimp) continue;
4935:     PetscCall(PetscDSGetDiscretization(ds, f, &obj));
4936:     PetscCall(PetscObjectGetClassId(obj, &id));
4937:     if (id == PETSCFE_CLASSID) useFEM = PETSC_TRUE;
4938:     if (id == PETSCFV_CLASSID) {
4939:       useFVM = PETSC_TRUE;
4940:       fvm    = (PetscFV)obj;
4941:     }
4942:   }
4943:   if (useFEM) {
4944:     PetscCall(DMGetCoordinateField(dm, &coordField));
4945:     PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
4946:     if (maxDegree <= 1) {
4947:       PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &affineQuad));
4948:       if (affineQuad) PetscCall(DMSNESGetFEGeom(coordField, cellIS, affineQuad, PETSC_FALSE, &affineGeom));
4949:     } else {
4950:       PetscCall(PetscCalloc2(Nf, &quads, Nf, &geoms));
4951:       for (f = 0; f < Nf; ++f) {
4952:         PetscObject  obj;
4953:         PetscClassId id;
4954:         PetscBool    fimp;

4956:         PetscCall(PetscDSGetImplicit(ds, f, &fimp));
4957:         if (isImplicit != fimp) continue;
4958:         PetscCall(PetscDSGetDiscretization(ds, f, &obj));
4959:         PetscCall(PetscObjectGetClassId(obj, &id));
4960:         if (id == PETSCFE_CLASSID) {
4961:           PetscFE fe = (PetscFE)obj;

4963:           PetscCall(PetscFEGetQuadrature(fe, &quads[f]));
4964:           PetscCall(PetscObjectReference((PetscObject)quads[f]));
4965:           PetscCall(DMSNESGetFEGeom(coordField, cellIS, quads[f], PETSC_FALSE, &geoms[f]));
4966:         }
4967:       }
4968:     }
4969:   }
4970:   // Handle non-essential (e.g. outflow) boundary values
4971:   if (useFVM) {
4972:     PetscCall(DMPlexInsertBoundaryValuesFVM(dm, fvm, locX, time, &locGrad));
4973:     PetscCall(DMPlexGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL));
4974:     PetscCall(DMPlexGetGradientDM(dm, fvm, &dmGrad));
4975:   }
4976:   /* Loop over chunks */
4977:   if (useFEM) PetscCall(ISCreate(PETSC_COMM_SELF, &chunkIS));
4978:   numCells      = cEnd - cStart;
4979:   numChunks     = 1;
4980:   cellChunkSize = numCells / numChunks;
4981:   faceChunkSize = (fEnd - fStart) / numChunks;
4982:   numChunks     = PetscMin(1, numCells);
4983:   for (chunk = 0; chunk < numChunks; ++chunk) {
4984:     PetscScalar     *elemVec, *fluxL, *fluxR;
4985:     PetscReal       *vol;
4986:     PetscFVFaceGeom *fgeom;
4987:     PetscInt         cS = cStart + chunk * cellChunkSize, cE = PetscMin(cS + cellChunkSize, cEnd), numCells = cE - cS, c;
4988:     PetscInt         fS = fStart + chunk * faceChunkSize, fE = PetscMin(fS + faceChunkSize, fEnd), numFaces = 0, face;

4990:     /* Extract field coefficients */
4991:     if (useFEM) {
4992:       PetscCall(ISGetPointSubrange(chunkIS, cS, cE, cells));
4993:       PetscCall(DMPlexGetCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a));
4994:       PetscCall(DMGetWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVec));
4995:       PetscCall(PetscArrayzero(elemVec, numCells * totDim));
4996:     }
4997:     if (useFVM) {
4998:       PetscCall(DMPlexGetFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR));
4999:       PetscCall(DMPlexGetFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol));
5000:       PetscCall(DMGetWorkArray(dm, numFaces * totDim, MPIU_SCALAR, &fluxL));
5001:       PetscCall(DMGetWorkArray(dm, numFaces * totDim, MPIU_SCALAR, &fluxR));
5002:       PetscCall(PetscArrayzero(fluxL, numFaces * totDim));
5003:       PetscCall(PetscArrayzero(fluxR, numFaces * totDim));
5004:     }
5005:     /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */
5006:     /* Loop over fields */
5007:     for (f = 0; f < Nf; ++f) {
5008:       PetscObject  obj;
5009:       PetscClassId id;
5010:       PetscBool    fimp;
5011:       PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;

5013:       key.field = f;
5014:       PetscCall(PetscDSGetImplicit(ds, f, &fimp));
5015:       if (isImplicit != fimp) continue;
5016:       PetscCall(PetscDSGetDiscretization(ds, f, &obj));
5017:       PetscCall(PetscObjectGetClassId(obj, &id));
5018:       if (id == PETSCFE_CLASSID) {
5019:         PetscFE         fe        = (PetscFE)obj;
5020:         PetscFEGeom    *geom      = affineGeom ? affineGeom : geoms[f];
5021:         PetscFEGeom    *chunkGeom = NULL;
5022:         PetscQuadrature quad      = affineQuad ? affineQuad : quads[f];
5023:         PetscInt        Nq, Nb;

5025:         PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
5026:         PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
5027:         PetscCall(PetscFEGetDimension(fe, &Nb));
5028:         blockSize = Nb;
5029:         batchSize = numBlocks * blockSize;
5030:         PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
5031:         numChunks = numCells / (numBatches * batchSize);
5032:         Ne        = numChunks * numBatches * batchSize;
5033:         Nr        = numCells % (numBatches * batchSize);
5034:         offset    = numCells - Nr;
5035:         /* Integrate FE residual to get elemVec (need fields at quadrature points) */
5036:         /*   For FV, I think we use a P0 basis and the cell coefficients (for subdivided cells, we can tweak the basis tabulation to be the indicator function) */
5037:         PetscCall(PetscFEGeomGetChunk(geom, 0, offset, &chunkGeom));
5038:         PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, u_t, dsAux, a, t, elemVec));
5039:         PetscCall(PetscFEGeomGetChunk(geom, offset, numCells, &chunkGeom));
5040:         PetscCall(PetscFEIntegrateResidual(ds, key, Nr, chunkGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), dsAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, &elemVec[offset * totDim]));
5041:         PetscCall(PetscFEGeomRestoreChunk(geom, offset, numCells, &chunkGeom));
5042:       } else if (id == PETSCFV_CLASSID) {
5043:         PetscFV fv = (PetscFV)obj;

5045:         Ne = numFaces;
5046:         /* Riemann solve over faces (need fields at face centroids) */
5047:         /*   We need to evaluate FE fields at those coordinates */
5048:         PetscCall(PetscFVIntegrateRHSFunction(fv, ds, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR));
5049:       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
5050:     }
5051:     /* Loop over domain */
5052:     if (useFEM) {
5053:       /* Add elemVec to locX */
5054:       for (c = cS; c < cE; ++c) {
5055:         const PetscInt cell = cells ? cells[c] : c;
5056:         const PetscInt cind = c - cStart;

5058:         if (mesh->printFEM > 1) PetscCall(DMPrintCellVector(cell, name, totDim, &elemVec[cind * totDim]));
5059:         if (ghostLabel) {
5060:           PetscInt ghostVal;

5062:           PetscCall(DMLabelGetValue(ghostLabel, cell, &ghostVal));
5063:           if (ghostVal > 0) continue;
5064:         }
5065:         PetscCall(DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cind * totDim], ADD_ALL_VALUES));
5066:       }
5067:     }
5068:     if (useFVM) {
5069:       PetscScalar *fa;
5070:       PetscInt     iface;

5072:       PetscCall(VecGetArray(locF, &fa));
5073:       for (f = 0; f < Nf; ++f) {
5074:         PetscFV      fv;
5075:         PetscObject  obj;
5076:         PetscClassId id;
5077:         PetscInt     cdim, foff, pdim;

5079:         PetscCall(DMGetCoordinateDim(dm, &cdim));
5080:         PetscCall(PetscDSGetDiscretization(ds, f, &obj));
5081:         PetscCall(PetscDSGetFieldOffset(ds, f, &foff));
5082:         PetscCall(PetscObjectGetClassId(obj, &id));
5083:         if (id != PETSCFV_CLASSID) continue;
5084:         fv = (PetscFV)obj;
5085:         PetscCall(PetscFVGetNumComponents(fv, &pdim));
5086:         /* Accumulate fluxes to cells */
5087:         for (face = fS, iface = 0; face < fE; ++face) {
5088:           const PetscInt *scells;
5089:           PetscScalar    *fL = NULL, *fR = NULL;
5090:           PetscInt        ghost, d, nsupp, nchild;

5092:           PetscCall(DMLabelGetValue(ghostLabel, face, &ghost));
5093:           PetscCall(DMPlexGetSupportSize(dm, face, &nsupp));
5094:           PetscCall(DMPlexGetTreeChildren(dm, face, &nchild, NULL));
5095:           if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
5096:           PetscCall(DMPlexGetSupport(dm, face, &scells));
5097:           PetscCall(DMLabelGetValue(ghostLabel, scells[0], &ghost));
5098:           if (ghost <= 0) PetscCall(DMPlexPointLocalFieldRef(dm, scells[0], f, fa, &fL));
5099:           PetscCall(DMLabelGetValue(ghostLabel, scells[1], &ghost));
5100:           if (ghost <= 0) PetscCall(DMPlexPointLocalFieldRef(dm, scells[1], f, fa, &fR));
5101:           if (mesh->printFVM > 1) {
5102:             PetscCall(DMPrintCellVectorReal(face, "Residual: normal", cdim, fgeom[iface].normal));
5103:             PetscCall(DMPrintCellVector(face, "Residual: left state", pdim, &uL[iface * totDim + foff]));
5104:             PetscCall(DMPrintCellVector(face, "Residual: right state", pdim, &uR[iface * totDim + foff]));
5105:             PetscCall(DMPrintCellVector(face, "Residual: left flux", pdim, &fluxL[iface * totDim + foff]));
5106:             PetscCall(DMPrintCellVector(face, "Residual: right flux", pdim, &fluxR[iface * totDim + foff]));
5107:           }
5108:           for (d = 0; d < pdim; ++d) {
5109:             if (fL) fL[d] -= fluxL[iface * totDim + foff + d];
5110:             if (fR) fR[d] += fluxR[iface * totDim + foff + d];
5111:           }
5112:           ++iface;
5113:         }
5114:       }
5115:       PetscCall(VecRestoreArray(locF, &fa));
5116:     }
5117:     /* Handle time derivative */
5118:     if (locX_t) {
5119:       PetscScalar *x_t, *fa;

5121:       PetscCall(VecGetArray(locF, &fa));
5122:       PetscCall(VecGetArray(locX_t, &x_t));
5123:       for (f = 0; f < Nf; ++f) {
5124:         PetscFV      fv;
5125:         PetscObject  obj;
5126:         PetscClassId id;
5127:         PetscInt     pdim, d;

5129:         PetscCall(PetscDSGetDiscretization(ds, f, &obj));
5130:         PetscCall(PetscObjectGetClassId(obj, &id));
5131:         if (id != PETSCFV_CLASSID) continue;
5132:         fv = (PetscFV)obj;
5133:         PetscCall(PetscFVGetNumComponents(fv, &pdim));
5134:         for (c = cS; c < cE; ++c) {
5135:           const PetscInt cell = cells ? cells[c] : c;
5136:           PetscScalar   *u_t, *r;

5138:           if (ghostLabel) {
5139:             PetscInt ghostVal;

5141:             PetscCall(DMLabelGetValue(ghostLabel, cell, &ghostVal));
5142:             if (ghostVal > 0) continue;
5143:           }
5144:           PetscCall(DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t));
5145:           PetscCall(DMPlexPointLocalFieldRef(dm, cell, f, fa, &r));
5146:           for (d = 0; d < pdim; ++d) r[d] += u_t[d];
5147:         }
5148:       }
5149:       PetscCall(VecRestoreArray(locX_t, &x_t));
5150:       PetscCall(VecRestoreArray(locF, &fa));
5151:     }
5152:     if (useFEM) {
5153:       PetscCall(DMPlexRestoreCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a));
5154:       PetscCall(DMRestoreWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVec));
5155:     }
5156:     if (useFVM) {
5157:       PetscCall(DMPlexRestoreFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR));
5158:       PetscCall(DMPlexRestoreFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol));
5159:       PetscCall(DMRestoreWorkArray(dm, numFaces * totDim, MPIU_SCALAR, &fluxL));
5160:       PetscCall(DMRestoreWorkArray(dm, numFaces * totDim, MPIU_SCALAR, &fluxR));
5161:       if (dmGrad) PetscCall(DMRestoreLocalVector(dmGrad, &locGrad));
5162:     }
5163:   }
5164:   if (useFEM) PetscCall(ISDestroy(&chunkIS));
5165:   PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));

5167:   if (useFEM) {
5168:     PetscCall(DMPlexComputeBdResidual_Internal(dm, locX, locX_t, t, locF, user));

5170:     if (maxDegree <= 1) {
5171:       PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, affineQuad, PETSC_FALSE, &affineGeom));
5172:       PetscCall(PetscQuadratureDestroy(&affineQuad));
5173:     } else {
5174:       for (f = 0; f < Nf; ++f) {
5175:         PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, quads[f], PETSC_FALSE, &geoms[f]));
5176:         PetscCall(PetscQuadratureDestroy(&quads[f]));
5177:       }
5178:       PetscCall(PetscFree2(quads, geoms));
5179:     }
5180:   }

5182:   /* FEM */
5183:   /* 1: Get sizes from dm and dmAux */
5184:   /* 2: Get geometric data */
5185:   /* 3: Handle boundary values */
5186:   /* 4: Loop over domain */
5187:   /*   Extract coefficients */
5188:   /* Loop over fields */
5189:   /*   Set tiling for FE*/
5190:   /*   Integrate FE residual to get elemVec */
5191:   /*     Loop over subdomain */
5192:   /*       Loop over quad points */
5193:   /*         Transform coords to real space */
5194:   /*         Evaluate field and aux fields at point */
5195:   /*         Evaluate residual at point */
5196:   /*         Transform residual to real space */
5197:   /*       Add residual to elemVec */
5198:   /* Loop over domain */
5199:   /*   Add elemVec to locX */

5201:   /* FVM */
5202:   /* Get geometric data */
5203:   /* If using gradients */
5204:   /*   Compute gradient data */
5205:   /*   Loop over domain faces */
5206:   /*     Count computational faces */
5207:   /*     Reconstruct cell gradient */
5208:   /*   Loop over domain cells */
5209:   /*     Limit cell gradients */
5210:   /* Handle boundary values */
5211:   /* Loop over domain faces */
5212:   /*   Read out field, centroid, normal, volume for each side of face */
5213:   /* Riemann solve over faces */
5214:   /* Loop over domain faces */
5215:   /*   Accumulate fluxes to cells */
5216:   /* TODO Change printFEM to printDisc here */
5217:   if (mesh->printFEM) {
5218:     Vec          locFbc;
5219:     PetscInt     pStart, pEnd, p, maxDof;
5220:     PetscScalar *zeroes;

5222:     PetscCall(VecDuplicate(locF, &locFbc));
5223:     PetscCall(VecCopy(locF, locFbc));
5224:     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
5225:     PetscCall(PetscSectionGetMaxDof(section, &maxDof));
5226:     PetscCall(PetscCalloc1(maxDof, &zeroes));
5227:     for (p = pStart; p < pEnd; p++) PetscCall(VecSetValuesSection(locFbc, section, p, zeroes, INSERT_BC_VALUES));
5228:     PetscCall(PetscFree(zeroes));
5229:     PetscCall(DMPrintLocalVec(dm, name, mesh->printTol, locFbc));
5230:     PetscCall(VecDestroy(&locFbc));
5231:   }
5232: end:
5233:   PetscCall(PetscLogEventEnd(DMPLEX_ResidualFEM, dm, 0, 0, 0));
5234:   PetscFunctionReturn(PETSC_SUCCESS);
5235: }

5237: /*
5238:   1) Allow multiple kernels for BdResidual for hybrid DS

5240:   DONE 2) Get out dsAux for either side at the same time as cohesive cell dsAux

5242:   DONE 3) Change DMGetCellFields() to get different aux data a[] for each side
5243:      - I think I just need to replace a[] with the closure from each face

5245:   4) Run both kernels for each non-hybrid field with correct dsAux, and then hybrid field as before
5246: */
5247: PetscErrorCode DMPlexComputeResidual_Hybrid_Internal(DM dm, PetscFormKey key[], IS cellIS, PetscReal time, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
5248: {
5249:   DM_Plex        *mesh       = (DM_Plex *)dm->data;
5250:   const char     *name       = "Hybrid Residual";
5251:   DM              dmAux[3]   = {NULL, NULL, NULL};
5252:   DMLabel         ghostLabel = NULL;
5253:   PetscDS         ds         = NULL;
5254:   PetscDS         dsIn       = NULL;
5255:   PetscDS         dsAux[3]   = {NULL, NULL, NULL};
5256:   Vec             locA[3]    = {NULL, NULL, NULL};
5257:   DM              dmScale[3] = {NULL, NULL, NULL};
5258:   PetscDS         dsScale[3] = {NULL, NULL, NULL};
5259:   Vec             locS[3]    = {NULL, NULL, NULL};
5260:   PetscSection    section    = NULL;
5261:   DMField         coordField = NULL;
5262:   PetscScalar    *a[3]       = {NULL, NULL, NULL};
5263:   PetscScalar    *s[3]       = {NULL, NULL, NULL};
5264:   PetscScalar    *u          = NULL, *u_t;
5265:   PetscScalar    *elemVecNeg, *elemVecPos, *elemVecCoh;
5266:   IS              chunkIS;
5267:   const PetscInt *cells;
5268:   PetscInt       *faces;
5269:   PetscInt        cStart, cEnd, numCells;
5270:   PetscInt        Nf, f, totDim, totDimIn, totDimAux[3], totDimScale[3], numChunks, cellChunkSize, chunk;
5271:   PetscInt        maxDegree  = PETSC_MAX_INT;
5272:   PetscQuadrature affineQuad = NULL, *quads = NULL;
5273:   PetscFEGeom    *affineGeom = NULL, **geoms = NULL;

5275:   PetscFunctionBegin;
5276:   PetscCall(PetscLogEventBegin(DMPLEX_ResidualFEM, dm, 0, 0, 0));
5277:   if (!cellIS) goto end;
5278:   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
5279:   PetscCall(ISGetLocalSize(cellIS, &numCells));
5280:   if (cStart >= cEnd) goto end;
5281:   if ((key[0].label == key[1].label) && (key[0].value == key[1].value) && (key[0].part == key[1].part)) {
5282:     const char *name;
5283:     PetscCall(PetscObjectGetName((PetscObject)key[0].label, &name));
5284:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Form keys for each side of a cohesive surface must be different (%s, %" PetscInt_FMT ", %" PetscInt_FMT ")", name, key[0].value, key[0].part);
5285:   }
5286:   /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */
5287:   /* FEM */
5288:   /* 1: Get sizes from dm and dmAux */
5289:   PetscCall(DMGetSection(dm, &section));
5290:   PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
5291:   PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &ds, &dsIn));
5292:   PetscCall(PetscDSGetNumFields(ds, &Nf));
5293:   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
5294:   PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn));
5295:   PetscCall(DMGetAuxiliaryVec(dm, key[2].label, key[2].value, key[2].part, &locA[2]));
5296:   if (locA[2]) {
5297:     const PetscInt cellStart = cells ? cells[cStart] : cStart;

5299:     PetscCall(VecGetDM(locA[2], &dmAux[2]));
5300:     PetscCall(DMGetCellDS(dmAux[2], cellStart, &dsAux[2], NULL));
5301:     PetscCall(PetscDSGetTotalDimension(dsAux[2], &totDimAux[2]));
5302:     {
5303:       const PetscInt *cone;
5304:       PetscInt        c;

5306:       PetscCall(DMPlexGetCone(dm, cellStart, &cone));
5307:       for (c = 0; c < 2; ++c) {
5308:         const PetscInt *support;
5309:         PetscInt        ssize, s;

5311:         PetscCall(DMPlexGetSupport(dm, cone[c], &support));
5312:         PetscCall(DMPlexGetSupportSize(dm, cone[c], &ssize));
5313:         PetscCheck(ssize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " from cell %" PetscInt_FMT " has support size %" PetscInt_FMT " != 2", cone[c], cellStart, ssize);
5314:         if (support[0] == cellStart) s = 1;
5315:         else if (support[1] == cellStart) s = 0;
5316:         else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " does not have cell %" PetscInt_FMT " in its support", cone[c], cellStart);
5317:         PetscCall(DMGetAuxiliaryVec(dm, key[c].label, key[c].value, key[c].part, &locA[c]));
5318:         PetscCheck(locA[c], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must have auxiliary vector for (%p, %" PetscInt_FMT ", %" PetscInt_FMT ")", (void *)key[c].label, key[c].value, key[c].part);
5319:         if (locA[c]) PetscCall(VecGetDM(locA[c], &dmAux[c]));
5320:         else dmAux[c] = dmAux[2];
5321:         PetscCall(DMGetCellDS(dmAux[c], support[s], &dsAux[c], NULL));
5322:         PetscCall(PetscDSGetTotalDimension(dsAux[c], &totDimAux[c]));
5323:       }
5324:     }
5325:   }
5326:   /* Handle mass matrix scaling
5327:        The field in key[2] is the field to be scaled, and the scaling field is the first in the dsScale */
5328:   PetscCall(DMGetAuxiliaryVec(dm, key[2].label, -key[2].value, key[2].part, &locS[2]));
5329:   if (locS[2]) {
5330:     const PetscInt cellStart = cells ? cells[cStart] : cStart;
5331:     PetscInt       Nb, Nbs;

5333:     PetscCall(VecGetDM(locS[2], &dmScale[2]));
5334:     PetscCall(DMGetCellDS(dmScale[2], cellStart, &dsScale[2], NULL));
5335:     PetscCall(PetscDSGetTotalDimension(dsScale[2], &totDimScale[2]));
5336:     // BRAD: This is not set correctly
5337:     key[2].field = 2;
5338:     PetscCall(PetscDSGetFieldSize(ds, key[2].field, &Nb));
5339:     PetscCall(PetscDSGetFieldSize(dsScale[2], 0, &Nbs));
5340:     PetscCheck(Nb == Nbs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field %" PetscInt_FMT " of size %" PetscInt_FMT " cannot be scaled by field of size %" PetscInt_FMT, key[2].field, Nb, Nbs);
5341:     {
5342:       const PetscInt *cone;
5343:       PetscInt        c;

5345:       locS[1] = locS[0] = locS[2];
5346:       dmScale[1] = dmScale[0] = dmScale[2];
5347:       PetscCall(DMPlexGetCone(dm, cellStart, &cone));
5348:       for (c = 0; c < 2; ++c) {
5349:         const PetscInt *support;
5350:         PetscInt        ssize, s;

5352:         PetscCall(DMPlexGetSupport(dm, cone[c], &support));
5353:         PetscCall(DMPlexGetSupportSize(dm, cone[c], &ssize));
5354:         PetscCheck(ssize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " from cell %" PetscInt_FMT " has support size %" PetscInt_FMT " != 2", cone[c], cellStart, ssize);
5355:         if (support[0] == cellStart) s = 1;
5356:         else if (support[1] == cellStart) s = 0;
5357:         else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " does not have cell %" PetscInt_FMT " in its support", cone[c], cellStart);
5358:         PetscCall(DMGetCellDS(dmScale[c], support[s], &dsScale[c], NULL));
5359:         PetscCall(PetscDSGetTotalDimension(dsScale[c], &totDimScale[c]));
5360:       }
5361:     }
5362:   }
5363:   /* 2: Setup geometric data */
5364:   PetscCall(DMGetCoordinateField(dm, &coordField));
5365:   PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
5366:   if (maxDegree > 1) {
5367:     PetscCall(PetscCalloc2(Nf, &quads, Nf, &geoms));
5368:     for (f = 0; f < Nf; ++f) {
5369:       PetscFE fe;

5371:       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
5372:       if (fe) {
5373:         PetscCall(PetscFEGetQuadrature(fe, &quads[f]));
5374:         PetscCall(PetscObjectReference((PetscObject)quads[f]));
5375:       }
5376:     }
5377:   }
5378:   /* Loop over chunks */
5379:   cellChunkSize = numCells;
5380:   numChunks     = !numCells ? 0 : PetscCeilReal(((PetscReal)numCells) / cellChunkSize);
5381:   PetscCall(PetscCalloc1(2 * cellChunkSize, &faces));
5382:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2 * cellChunkSize, faces, PETSC_USE_POINTER, &chunkIS));
5383:   /* Extract field coefficients */
5384:   /* NOTE This needs the end cap faces to have identical orientations */
5385:   PetscCall(DMPlexGetHybridCellFields(dm, cellIS, locX, locX_t, locA[2], &u, &u_t, &a[2]));
5386:   PetscCall(DMPlexGetHybridFields(dm, dmAux, dsAux, cellIS, locA, PETSC_TRUE, a));
5387:   PetscCall(DMPlexGetHybridFields(dm, dmScale, dsScale, cellIS, locS, PETSC_TRUE, s));
5388:   PetscCall(DMGetWorkArray(dm, cellChunkSize * totDim, MPIU_SCALAR, &elemVecNeg));
5389:   PetscCall(DMGetWorkArray(dm, cellChunkSize * totDim, MPIU_SCALAR, &elemVecPos));
5390:   PetscCall(DMGetWorkArray(dm, cellChunkSize * totDim, MPIU_SCALAR, &elemVecCoh));
5391:   for (chunk = 0; chunk < numChunks; ++chunk) {
5392:     PetscInt cS = cStart + chunk * cellChunkSize, cE = PetscMin(cS + cellChunkSize, cEnd), numCells = cE - cS, c;

5394:     PetscCall(PetscArrayzero(elemVecNeg, cellChunkSize * totDim));
5395:     PetscCall(PetscArrayzero(elemVecPos, cellChunkSize * totDim));
5396:     PetscCall(PetscArrayzero(elemVecCoh, cellChunkSize * totDim));
5397:     /* Get faces */
5398:     for (c = cS; c < cE; ++c) {
5399:       const PetscInt  cell = cells ? cells[c] : c;
5400:       const PetscInt *cone;
5401:       PetscCall(DMPlexGetCone(dm, cell, &cone));
5402:       faces[(c - cS) * 2 + 0] = cone[0];
5403:       faces[(c - cS) * 2 + 1] = cone[1];
5404:     }
5405:     PetscCall(ISGeneralSetIndices(chunkIS, 2 * cellChunkSize, faces, PETSC_USE_POINTER));
5406:     /* Get geometric data */
5407:     if (maxDegree <= 1) {
5408:       if (!affineQuad) PetscCall(DMFieldCreateDefaultQuadrature(coordField, chunkIS, &affineQuad));
5409:       if (affineQuad) PetscCall(DMSNESGetFEGeom(coordField, chunkIS, affineQuad, PETSC_TRUE, &affineGeom));
5410:     } else {
5411:       for (f = 0; f < Nf; ++f) {
5412:         if (quads[f]) PetscCall(DMSNESGetFEGeom(coordField, chunkIS, quads[f], PETSC_TRUE, &geoms[f]));
5413:       }
5414:     }
5415:     /* Loop over fields */
5416:     for (f = 0; f < Nf; ++f) {
5417:       PetscFE         fe;
5418:       PetscFEGeom    *geom      = affineGeom ? affineGeom : geoms[f];
5419:       PetscFEGeom    *chunkGeom = NULL, *remGeom = NULL;
5420:       PetscQuadrature quad = affineQuad ? affineQuad : quads[f];
5421:       PetscInt        numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset, Nq, Nb;
5422:       PetscBool       isCohesiveField;

5424:       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
5425:       if (!fe) continue;
5426:       PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
5427:       PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
5428:       PetscCall(PetscFEGetDimension(fe, &Nb));
5429:       blockSize = Nb;
5430:       batchSize = numBlocks * blockSize;
5431:       PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
5432:       numChunks = numCells / (numBatches * batchSize);
5433:       Ne        = numChunks * numBatches * batchSize;
5434:       Nr        = numCells % (numBatches * batchSize);
5435:       offset    = numCells - Nr;
5436:       PetscCall(PetscFEGeomGetChunk(geom, 0, offset * 2, &chunkGeom));
5437:       PetscCall(PetscFEGeomGetChunk(geom, offset * 2, numCells * 2, &remGeom));
5438:       PetscCall(PetscDSGetCohesive(ds, f, &isCohesiveField));
5439:       chunkGeom->isCohesive = remGeom->isCohesive = PETSC_TRUE;
5440:       key[0].field                                = f;
5441:       key[1].field                                = f;
5442:       key[2].field                                = f;
5443:       PetscCall(PetscFEIntegrateHybridResidual(ds, dsIn, key[0], 0, Ne, chunkGeom, u, u_t, dsAux[0], a[0], t, elemVecNeg));
5444:       PetscCall(PetscFEIntegrateHybridResidual(ds, dsIn, key[0], 0, Nr, remGeom, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[0], PetscSafePointerPlusOffset(a[0], offset * totDimAux[0]), t, &elemVecNeg[offset * totDim]));
5445:       PetscCall(PetscFEIntegrateHybridResidual(ds, dsIn, key[1], 1, Ne, chunkGeom, u, u_t, dsAux[1], a[1], t, elemVecPos));
5446:       PetscCall(PetscFEIntegrateHybridResidual(ds, dsIn, key[1], 1, Nr, remGeom, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[1], PetscSafePointerPlusOffset(a[1], offset * totDimAux[1]), t, &elemVecPos[offset * totDim]));
5447:       PetscCall(PetscFEIntegrateHybridResidual(ds, dsIn, key[2], 2, Ne, chunkGeom, u, u_t, dsAux[2], a[2], t, elemVecCoh));
5448:       PetscCall(PetscFEIntegrateHybridResidual(ds, dsIn, key[2], 2, Nr, remGeom, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[2], PetscSafePointerPlusOffset(a[2], offset * totDimAux[2]), t, &elemVecCoh[offset * totDim]));
5449:       PetscCall(PetscFEGeomRestoreChunk(geom, offset, numCells, &remGeom));
5450:       PetscCall(PetscFEGeomRestoreChunk(geom, 0, offset, &chunkGeom));
5451:     }
5452:     /* Add elemVec to locX */
5453:     for (c = cS; c < cE; ++c) {
5454:       const PetscInt cell = cells ? cells[c] : c;
5455:       const PetscInt cind = c - cStart;
5456:       PetscInt       i;

5458:       /* Scale element values */
5459:       if (locS[0]) {
5460:         PetscInt  Nb, off = cind * totDim, soff = cind * totDimScale[0];
5461:         PetscBool cohesive;

5463:         for (f = 0; f < Nf; ++f) {
5464:           PetscCall(PetscDSGetFieldSize(ds, f, &Nb));
5465:           PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
5466:           if (f == key[2].field) {
5467:             PetscCheck(cohesive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Scaling should not happen for face fields");
5468:             // No cohesive scaling field is currently input
5469:             for (i = 0; i < Nb; ++i) elemVecCoh[off + i] += s[0][soff + i] * elemVecNeg[off + i] + s[1][soff + i] * elemVecPos[off + i];
5470:             off += Nb;
5471:           } else {
5472:             const PetscInt N = cohesive ? Nb : Nb * 2;

5474:             for (i = 0; i < N; ++i) elemVecCoh[off + i] += elemVecNeg[off + i] + elemVecPos[off + i];
5475:             off += N;
5476:           }
5477:         }
5478:       } else {
5479:         for (i = cind * totDim; i < (cind + 1) * totDim; ++i) elemVecCoh[i] += elemVecNeg[i] + elemVecPos[i];
5480:       }
5481:       if (mesh->printFEM > 1) PetscCall(DMPrintCellVector(cell, name, totDim, &elemVecCoh[cind * totDim]));
5482:       if (ghostLabel) {
5483:         PetscInt ghostVal;

5485:         PetscCall(DMLabelGetValue(ghostLabel, cell, &ghostVal));
5486:         if (ghostVal > 0) continue;
5487:       }
5488:       PetscCall(DMPlexVecSetClosure(dm, section, locF, cell, &elemVecCoh[cind * totDim], ADD_ALL_VALUES));
5489:     }
5490:   }
5491:   PetscCall(DMPlexRestoreCellFields(dm, cellIS, locX, locX_t, locA[2], &u, &u_t, &a[2]));
5492:   PetscCall(DMPlexRestoreHybridFields(dm, dmAux, dsAux, cellIS, locA, PETSC_TRUE, a));
5493:   PetscCall(DMPlexRestoreHybridFields(dm, dmScale, dsScale, cellIS, locS, PETSC_TRUE, s));
5494:   PetscCall(DMRestoreWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVecNeg));
5495:   PetscCall(DMRestoreWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVecPos));
5496:   PetscCall(DMRestoreWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVecCoh));
5497:   PetscCall(PetscFree(faces));
5498:   PetscCall(ISDestroy(&chunkIS));
5499:   PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
5500:   if (maxDegree <= 1) {
5501:     PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, affineQuad, PETSC_FALSE, &affineGeom));
5502:     PetscCall(PetscQuadratureDestroy(&affineQuad));
5503:   } else {
5504:     for (f = 0; f < Nf; ++f) {
5505:       if (geoms) PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, quads[f], PETSC_FALSE, &geoms[f]));
5506:       if (quads) PetscCall(PetscQuadratureDestroy(&quads[f]));
5507:     }
5508:     PetscCall(PetscFree2(quads, geoms));
5509:   }
5510:   if (mesh->printFEM) {
5511:     Vec          locFbc;
5512:     PetscInt     pStart, pEnd, p, maxDof;
5513:     PetscScalar *zeroes;

5515:     PetscCall(VecDuplicate(locF, &locFbc));
5516:     PetscCall(VecCopy(locF, locFbc));
5517:     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
5518:     PetscCall(PetscSectionGetMaxDof(section, &maxDof));
5519:     PetscCall(PetscCalloc1(maxDof, &zeroes));
5520:     for (p = pStart; p < pEnd; p++) PetscCall(VecSetValuesSection(locFbc, section, p, zeroes, INSERT_BC_VALUES));
5521:     PetscCall(PetscFree(zeroes));
5522:     PetscCall(DMPrintLocalVec(dm, name, mesh->printTol, locFbc));
5523:     PetscCall(VecDestroy(&locFbc));
5524:   }
5525: end:
5526:   PetscCall(PetscLogEventEnd(DMPLEX_ResidualFEM, dm, 0, 0, 0));
5527:   PetscFunctionReturn(PETSC_SUCCESS);
5528: }

5530: static PetscErrorCode DMPlexComputeBdJacobian_Single_Internal(DM dm, PetscReal t, PetscWeakForm wf, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt fieldI, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, DMField coordField, IS facetIS)
5531: {
5532:   DM_Plex        *mesh = (DM_Plex *)dm->data;
5533:   DM              plex = NULL, plexA = NULL, tdm;
5534:   DMEnclosureType encAux;
5535:   PetscDS         ds, dsAux           = NULL;
5536:   PetscSection    section, sectionAux = NULL;
5537:   PetscSection    globalSection;
5538:   Vec             locA = NULL, tv;
5539:   PetscScalar    *u = NULL, *u_t = NULL, *a = NULL, *elemMat = NULL, *elemMatP = NULL;
5540:   PetscInt        v;
5541:   PetscInt        Nf, totDim, totDimAux = 0;
5542:   PetscBool       hasJac = PETSC_FALSE, hasPrec = PETSC_FALSE, transform;

5544:   PetscFunctionBegin;
5545:   PetscCall(DMHasBasisTransform(dm, &transform));
5546:   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
5547:   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
5548:   PetscCall(DMGetLocalSection(dm, &section));
5549:   PetscCall(DMGetDS(dm, &ds));
5550:   PetscCall(PetscDSGetNumFields(ds, &Nf));
5551:   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
5552:   PetscCall(PetscWeakFormHasBdJacobian(wf, &hasJac));
5553:   PetscCall(PetscWeakFormHasBdJacobianPreconditioner(wf, &hasPrec));
5554:   if (!hasJac && !hasPrec) PetscFunctionReturn(PETSC_SUCCESS);
5555:   PetscCall(DMConvert(dm, DMPLEX, &plex));
5556:   PetscCall(DMGetAuxiliaryVec(dm, label, values[0], 0, &locA));
5557:   if (locA) {
5558:     DM dmAux;

5560:     PetscCall(VecGetDM(locA, &dmAux));
5561:     PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
5562:     PetscCall(DMConvert(dmAux, DMPLEX, &plexA));
5563:     PetscCall(DMGetDS(plexA, &dsAux));
5564:     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
5565:     PetscCall(DMGetLocalSection(plexA, &sectionAux));
5566:   }

5568:   PetscCall(DMGetGlobalSection(dm, &globalSection));
5569:   for (v = 0; v < numValues; ++v) {
5570:     PetscFEGeom    *fgeom;
5571:     PetscInt        maxDegree;
5572:     PetscQuadrature qGeom = NULL;
5573:     IS              pointIS;
5574:     const PetscInt *points;
5575:     PetscFormKey    key;
5576:     PetscInt        numFaces, face, Nq;

5578:     key.label = label;
5579:     key.value = values[v];
5580:     key.part  = 0;
5581:     PetscCall(DMLabelGetStratumIS(label, values[v], &pointIS));
5582:     if (!pointIS) continue; /* No points with that id on this process */
5583:     {
5584:       IS isectIS;

5586:       /* TODO: Special cases of ISIntersect where it is quick to check a prior if one is a superset of the other */
5587:       PetscCall(ISIntersect_Caching_Internal(facetIS, pointIS, &isectIS));
5588:       PetscCall(ISDestroy(&pointIS));
5589:       pointIS = isectIS;
5590:     }
5591:     PetscCall(ISGetLocalSize(pointIS, &numFaces));
5592:     PetscCall(ISGetIndices(pointIS, &points));
5593:     PetscCall(PetscMalloc5(numFaces * totDim, &u, locX_t ? numFaces * totDim : 0, &u_t, hasJac ? numFaces * totDim * totDim : 0, &elemMat, hasPrec ? numFaces * totDim * totDim : 0, &elemMatP, locA ? numFaces * totDimAux : 0, &a));
5594:     PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
5595:     if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &qGeom));
5596:     if (!qGeom) {
5597:       PetscFE fe;

5599:       PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&fe));
5600:       PetscCall(PetscFEGetFaceQuadrature(fe, &qGeom));
5601:       PetscCall(PetscObjectReference((PetscObject)qGeom));
5602:     }
5603:     PetscCall(PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL));
5604:     PetscCall(DMSNESGetFEGeom(coordField, pointIS, qGeom, PETSC_TRUE, &fgeom));
5605:     for (face = 0; face < numFaces; ++face) {
5606:       const PetscInt point = points[face], *support;
5607:       PetscScalar   *x     = NULL;
5608:       PetscInt       i;

5610:       PetscCall(DMPlexGetSupport(dm, point, &support));
5611:       PetscCall(DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x));
5612:       for (i = 0; i < totDim; ++i) u[face * totDim + i] = x[i];
5613:       PetscCall(DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x));
5614:       if (locX_t) {
5615:         PetscCall(DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x));
5616:         for (i = 0; i < totDim; ++i) u_t[face * totDim + i] = x[i];
5617:         PetscCall(DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x));
5618:       }
5619:       if (locA) {
5620:         PetscInt subp;
5621:         PetscCall(DMGetEnclosurePoint(plexA, dm, encAux, support[0], &subp));
5622:         PetscCall(DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x));
5623:         for (i = 0; i < totDimAux; ++i) a[face * totDimAux + i] = x[i];
5624:         PetscCall(DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x));
5625:       }
5626:     }
5627:     if (elemMat) PetscCall(PetscArrayzero(elemMat, numFaces * totDim * totDim));
5628:     if (elemMatP) PetscCall(PetscArrayzero(elemMatP, numFaces * totDim * totDim));
5629:     {
5630:       PetscFE  fe;
5631:       PetscInt Nb;
5632:       /* Conforming batches */
5633:       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
5634:       /* Remainder */
5635:       PetscFEGeom *chunkGeom = NULL;
5636:       PetscInt     fieldJ, Nr, offset;

5638:       PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&fe));
5639:       PetscCall(PetscFEGetDimension(fe, &Nb));
5640:       PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
5641:       blockSize = Nb;
5642:       batchSize = numBlocks * blockSize;
5643:       PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
5644:       numChunks = numFaces / (numBatches * batchSize);
5645:       Ne        = numChunks * numBatches * batchSize;
5646:       Nr        = numFaces % (numBatches * batchSize);
5647:       offset    = numFaces - Nr;
5648:       PetscCall(PetscFEGeomGetChunk(fgeom, 0, offset, &chunkGeom));
5649:       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
5650:         key.field = fieldI * Nf + fieldJ;
5651:         if (hasJac) PetscCall(PetscFEIntegrateBdJacobian(ds, wf, PETSCFE_JACOBIAN, key, Ne, chunkGeom, u, u_t, dsAux, a, t, X_tShift, elemMat));
5652:         if (hasPrec) PetscCall(PetscFEIntegrateBdJacobian(ds, wf, PETSCFE_JACOBIAN_PRE, key, Ne, chunkGeom, u, u_t, dsAux, a, t, X_tShift, elemMatP));
5653:       }
5654:       PetscCall(PetscFEGeomGetChunk(fgeom, offset, numFaces, &chunkGeom));
5655:       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
5656:         key.field = fieldI * Nf + fieldJ;
5657:         if (hasJac)
5658:           PetscCall(PetscFEIntegrateBdJacobian(ds, wf, PETSCFE_JACOBIAN, key, Nr, chunkGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), dsAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMat[offset * totDim * totDim]));
5659:         if (hasPrec)
5660:           PetscCall(PetscFEIntegrateBdJacobian(ds, wf, PETSCFE_JACOBIAN_PRE, key, Nr, chunkGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), dsAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMatP[offset * totDim * totDim]));
5661:       }
5662:       PetscCall(PetscFEGeomRestoreChunk(fgeom, offset, numFaces, &chunkGeom));
5663:     }
5664:     for (face = 0; face < numFaces; ++face) {
5665:       const PetscInt point = points[face], *support;

5667:       /* Transform to global basis before insertion in Jacobian */
5668:       PetscCall(DMPlexGetSupport(plex, point, &support));
5669:       if (hasJac && transform) PetscCall(DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, support[0], PETSC_TRUE, totDim, &elemMat[face * totDim * totDim]));
5670:       if (hasPrec && transform) PetscCall(DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, support[0], PETSC_TRUE, totDim, &elemMatP[face * totDim * totDim]));
5671:       if (hasPrec) {
5672:         if (hasJac) {
5673:           if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(point, "BdJacobian", totDim, totDim, &elemMat[face * totDim * totDim]));
5674:           PetscCall(DMPlexMatSetClosure_Internal(plex, section, globalSection, mesh->useMatClPerm, Jac, support[0], &elemMat[face * totDim * totDim], ADD_VALUES));
5675:         }
5676:         if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(point, "BdJacobian", totDim, totDim, &elemMatP[face * totDim * totDim]));
5677:         PetscCall(DMPlexMatSetClosure_Internal(plex, section, globalSection, mesh->useMatClPerm, JacP, support[0], &elemMatP[face * totDim * totDim], ADD_VALUES));
5678:       } else {
5679:         if (hasJac) {
5680:           if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(point, "BdJacobian", totDim, totDim, &elemMat[face * totDim * totDim]));
5681:           PetscCall(DMPlexMatSetClosure_Internal(plex, section, globalSection, mesh->useMatClPerm, Jac, support[0], &elemMat[face * totDim * totDim], ADD_VALUES));
5682:         }
5683:       }
5684:     }
5685:     PetscCall(DMSNESRestoreFEGeom(coordField, pointIS, qGeom, PETSC_TRUE, &fgeom));
5686:     PetscCall(PetscQuadratureDestroy(&qGeom));
5687:     PetscCall(ISRestoreIndices(pointIS, &points));
5688:     PetscCall(ISDestroy(&pointIS));
5689:     PetscCall(PetscFree5(u, u_t, elemMat, elemMatP, a));
5690:   }
5691:   if (plex) PetscCall(DMDestroy(&plex));
5692:   if (plexA) PetscCall(DMDestroy(&plexA));
5693:   PetscFunctionReturn(PETSC_SUCCESS);
5694: }

5696: PetscErrorCode DMPlexComputeBdJacobianSingle(DM dm, PetscReal t, PetscWeakForm wf, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt field, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP)
5697: {
5698:   DMField  coordField;
5699:   DMLabel  depthLabel;
5700:   IS       facetIS;
5701:   PetscInt dim;

5703:   PetscFunctionBegin;
5704:   PetscCall(DMGetDimension(dm, &dim));
5705:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
5706:   PetscCall(DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS));
5707:   PetscCall(DMGetCoordinateField(dm, &coordField));
5708:   PetscCall(DMPlexComputeBdJacobian_Single_Internal(dm, t, wf, label, numValues, values, field, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS));
5709:   PetscCall(ISDestroy(&facetIS));
5710:   PetscFunctionReturn(PETSC_SUCCESS);
5711: }

5713: static PetscErrorCode DMPlexComputeBdJacobian_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
5714: {
5715:   PetscDS  prob;
5716:   PetscInt dim, numBd, bd;
5717:   DMLabel  depthLabel;
5718:   DMField  coordField = NULL;
5719:   IS       facetIS;

5721:   PetscFunctionBegin;
5722:   PetscCall(DMGetDS(dm, &prob));
5723:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
5724:   PetscCall(DMGetDimension(dm, &dim));
5725:   PetscCall(DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS));
5726:   PetscCall(PetscDSGetNumBoundary(prob, &numBd));
5727:   PetscCall(DMGetCoordinateField(dm, &coordField));
5728:   for (bd = 0; bd < numBd; ++bd) {
5729:     PetscWeakForm           wf;
5730:     DMBoundaryConditionType type;
5731:     DMLabel                 label;
5732:     const PetscInt         *values;
5733:     PetscInt                fieldI, numValues;
5734:     PetscObject             obj;
5735:     PetscClassId            id;

5737:     PetscCall(PetscDSGetBoundary(prob, bd, &wf, &type, NULL, &label, &numValues, &values, &fieldI, NULL, NULL, NULL, NULL, NULL));
5738:     if (type & DM_BC_ESSENTIAL) continue;
5739:     PetscCall(PetscDSGetDiscretization(prob, fieldI, &obj));
5740:     PetscCall(PetscObjectGetClassId(obj, &id));
5741:     if (id != PETSCFE_CLASSID) continue;
5742:     PetscCall(DMPlexComputeBdJacobian_Single_Internal(dm, t, wf, label, numValues, values, fieldI, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS));
5743:   }
5744:   PetscCall(ISDestroy(&facetIS));
5745:   PetscFunctionReturn(PETSC_SUCCESS);
5746: }

5748: PetscErrorCode DMPlexComputeJacobian_Internal(DM dm, PetscFormKey key, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Mat Jac, Mat JacP, void *user)
5749: {
5750:   DM_Plex        *mesh  = (DM_Plex *)dm->data;
5751:   const char     *name  = "Jacobian";
5752:   DM              dmAux = NULL, plex, tdm;
5753:   DMEnclosureType encAux;
5754:   Vec             A, tv;
5755:   DMField         coordField;
5756:   PetscDS         prob, probAux = NULL;
5757:   PetscSection    section, globalSection, sectionAux;
5758:   PetscScalar    *elemMat, *elemMatP, *elemMatD, *u, *u_t, *a = NULL;
5759:   const PetscInt *cells;
5760:   PetscInt        Nf, fieldI, fieldJ;
5761:   PetscInt        totDim, totDimAux = 0, cStart, cEnd, numCells, c;
5762:   PetscBool       hasJac = PETSC_FALSE, hasPrec = PETSC_FALSE, hasDyn, hasFV = PETSC_FALSE, transform;

5764:   PetscFunctionBegin;
5765:   PetscCall(PetscLogEventBegin(DMPLEX_JacobianFEM, dm, 0, 0, 0));
5766:   if (!cellIS) goto end;
5767:   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
5768:   PetscCall(ISGetLocalSize(cellIS, &numCells));
5769:   if (cStart >= cEnd) goto end;
5770:   PetscCall(DMHasBasisTransform(dm, &transform));
5771:   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
5772:   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
5773:   PetscCall(DMGetLocalSection(dm, &section));
5774:   PetscCall(DMGetGlobalSection(dm, &globalSection));
5775:   PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob, NULL));
5776:   PetscCall(PetscDSGetNumFields(prob, &Nf));
5777:   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
5778:   PetscCall(PetscDSHasJacobian(prob, &hasJac));
5779:   PetscCall(PetscDSHasJacobianPreconditioner(prob, &hasPrec));
5780:   /* user passed in the same matrix, avoid double contributions and
5781:      only assemble the Jacobian */
5782:   if (hasJac && Jac == JacP) hasPrec = PETSC_FALSE;
5783:   PetscCall(PetscDSHasDynamicJacobian(prob, &hasDyn));
5784:   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
5785:   PetscCall(DMGetAuxiliaryVec(dm, key.label, key.value, key.part, &A));
5786:   if (A) {
5787:     PetscCall(VecGetDM(A, &dmAux));
5788:     PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
5789:     PetscCall(DMConvert(dmAux, DMPLEX, &plex));
5790:     PetscCall(DMGetLocalSection(plex, &sectionAux));
5791:     PetscCall(DMGetDS(dmAux, &probAux));
5792:     PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
5793:   }
5794:   PetscCall(PetscMalloc5(numCells * totDim, &u, X_t ? numCells * totDim : 0, &u_t, hasJac ? numCells * totDim * totDim : 0, &elemMat, hasPrec ? numCells * totDim * totDim : 0, &elemMatP, hasDyn ? numCells * totDim * totDim : 0, &elemMatD));
5795:   if (dmAux) PetscCall(PetscMalloc1(numCells * totDimAux, &a));
5796:   PetscCall(DMGetCoordinateField(dm, &coordField));
5797:   for (c = cStart; c < cEnd; ++c) {
5798:     const PetscInt cell = cells ? cells[c] : c;
5799:     const PetscInt cind = c - cStart;
5800:     PetscScalar   *x = NULL, *x_t = NULL;
5801:     PetscInt       i;

5803:     PetscCall(DMPlexVecGetClosure(dm, section, X, cell, NULL, &x));
5804:     for (i = 0; i < totDim; ++i) u[cind * totDim + i] = x[i];
5805:     PetscCall(DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x));
5806:     if (X_t) {
5807:       PetscCall(DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t));
5808:       for (i = 0; i < totDim; ++i) u_t[cind * totDim + i] = x_t[i];
5809:       PetscCall(DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t));
5810:     }
5811:     if (dmAux) {
5812:       PetscInt subcell;
5813:       PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell));
5814:       PetscCall(DMPlexVecGetClosure(plex, sectionAux, A, subcell, NULL, &x));
5815:       for (i = 0; i < totDimAux; ++i) a[cind * totDimAux + i] = x[i];
5816:       PetscCall(DMPlexVecRestoreClosure(plex, sectionAux, A, subcell, NULL, &x));
5817:     }
5818:   }
5819:   if (hasJac) PetscCall(PetscArrayzero(elemMat, numCells * totDim * totDim));
5820:   if (hasPrec) PetscCall(PetscArrayzero(elemMatP, numCells * totDim * totDim));
5821:   if (hasDyn) PetscCall(PetscArrayzero(elemMatD, numCells * totDim * totDim));
5822:   for (fieldI = 0; fieldI < Nf; ++fieldI) {
5823:     PetscClassId    id;
5824:     PetscFE         fe;
5825:     PetscQuadrature qGeom = NULL;
5826:     PetscInt        Nb;
5827:     /* Conforming batches */
5828:     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
5829:     /* Remainder */
5830:     PetscInt     Nr, offset, Nq;
5831:     PetscInt     maxDegree;
5832:     PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;

5834:     PetscCall(PetscDSGetDiscretization(prob, fieldI, (PetscObject *)&fe));
5835:     PetscCall(PetscObjectGetClassId((PetscObject)fe, &id));
5836:     if (id == PETSCFV_CLASSID) {
5837:       hasFV = PETSC_TRUE;
5838:       continue;
5839:     }
5840:     PetscCall(PetscFEGetDimension(fe, &Nb));
5841:     PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
5842:     PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
5843:     if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &qGeom));
5844:     if (!qGeom) {
5845:       PetscCall(PetscFEGetQuadrature(fe, &qGeom));
5846:       PetscCall(PetscObjectReference((PetscObject)qGeom));
5847:     }
5848:     PetscCall(PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL));
5849:     PetscCall(DMSNESGetFEGeom(coordField, cellIS, qGeom, PETSC_FALSE, &cgeomFEM));
5850:     blockSize = Nb;
5851:     batchSize = numBlocks * blockSize;
5852:     PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
5853:     numChunks = numCells / (numBatches * batchSize);
5854:     Ne        = numChunks * numBatches * batchSize;
5855:     Nr        = numCells % (numBatches * batchSize);
5856:     offset    = numCells - Nr;
5857:     PetscCall(PetscFEGeomGetChunk(cgeomFEM, 0, offset, &chunkGeom));
5858:     PetscCall(PetscFEGeomGetChunk(cgeomFEM, offset, numCells, &remGeom));
5859:     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
5860:       key.field = fieldI * Nf + fieldJ;
5861:       if (hasJac) {
5862:         PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat));
5863:         PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMat[offset * totDim * totDim]));
5864:       }
5865:       if (hasPrec) {
5866:         PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_PRE, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatP));
5867:         PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_PRE, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMatP[offset * totDim * totDim]));
5868:       }
5869:       if (hasDyn) {
5870:         PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD));
5871:         PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMatD[offset * totDim * totDim]));
5872:       }
5873:     }
5874:     PetscCall(PetscFEGeomRestoreChunk(cgeomFEM, offset, numCells, &remGeom));
5875:     PetscCall(PetscFEGeomRestoreChunk(cgeomFEM, 0, offset, &chunkGeom));
5876:     PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, qGeom, PETSC_FALSE, &cgeomFEM));
5877:     PetscCall(PetscQuadratureDestroy(&qGeom));
5878:   }
5879:   /*   Add contribution from X_t */
5880:   if (hasDyn) {
5881:     for (c = 0; c < numCells * totDim * totDim; ++c) elemMat[c] += X_tShift * elemMatD[c];
5882:   }
5883:   if (hasFV) {
5884:     PetscClassId id;
5885:     PetscFV      fv;
5886:     PetscInt     offsetI, NcI, NbI = 1, fc, f;

5888:     for (fieldI = 0; fieldI < Nf; ++fieldI) {
5889:       PetscCall(PetscDSGetDiscretization(prob, fieldI, (PetscObject *)&fv));
5890:       PetscCall(PetscDSGetFieldOffset(prob, fieldI, &offsetI));
5891:       PetscCall(PetscObjectGetClassId((PetscObject)fv, &id));
5892:       if (id != PETSCFV_CLASSID) continue;
5893:       /* Put in the weighted identity */
5894:       PetscCall(PetscFVGetNumComponents(fv, &NcI));
5895:       for (c = cStart; c < cEnd; ++c) {
5896:         const PetscInt cind    = c - cStart;
5897:         const PetscInt eOffset = cind * totDim * totDim;
5898:         PetscReal      vol;

5900:         PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL));
5901:         for (fc = 0; fc < NcI; ++fc) {
5902:           for (f = 0; f < NbI; ++f) {
5903:             const PetscInt i = offsetI + f * NcI + fc;
5904:             if (hasPrec) {
5905:               if (hasJac) elemMat[eOffset + i * totDim + i] = vol;
5906:               elemMatP[eOffset + i * totDim + i] = vol;
5907:             } else {
5908:               elemMat[eOffset + i * totDim + i] = vol;
5909:             }
5910:           }
5911:         }
5912:       }
5913:     }
5914:     /* No allocated space for FV stuff, so ignore the zero entries */
5915:     PetscCall(MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
5916:   }
5917:   /* Insert values into matrix */
5918:   for (c = cStart; c < cEnd; ++c) {
5919:     const PetscInt cell = cells ? cells[c] : c;
5920:     const PetscInt cind = c - cStart;

5922:     /* Transform to global basis before insertion in Jacobian */
5923:     if (transform) PetscCall(DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, cell, PETSC_TRUE, totDim, &elemMat[cind * totDim * totDim]));
5924:     if (hasPrec) {
5925:       if (hasJac) {
5926:         if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind * totDim * totDim]));
5927:         PetscCall(DMPlexMatSetClosure_Internal(dm, section, globalSection, mesh->useMatClPerm, Jac, cell, &elemMat[cind * totDim * totDim], ADD_VALUES));
5928:       }
5929:       if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatP[cind * totDim * totDim]));
5930:       PetscCall(DMPlexMatSetClosure_Internal(dm, section, globalSection, mesh->useMatClPerm, JacP, cell, &elemMatP[cind * totDim * totDim], ADD_VALUES));
5931:     } else {
5932:       if (hasJac) {
5933:         if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind * totDim * totDim]));
5934:         PetscCall(DMPlexMatSetClosure_Internal(dm, section, globalSection, mesh->useMatClPerm, JacP, cell, &elemMat[cind * totDim * totDim], ADD_VALUES));
5935:       }
5936:     }
5937:   }
5938:   PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
5939:   if (hasFV) PetscCall(MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE));
5940:   PetscCall(PetscFree5(u, u_t, elemMat, elemMatP, elemMatD));
5941:   if (dmAux) {
5942:     PetscCall(PetscFree(a));
5943:     PetscCall(DMDestroy(&plex));
5944:   }
5945:   /* Compute boundary integrals */
5946:   PetscCall(DMPlexComputeBdJacobian_Internal(dm, X, X_t, t, X_tShift, Jac, JacP, user));
5947:   /* Assemble matrix */
5948: end: {
5949:   PetscBool assOp = hasJac && hasPrec ? PETSC_TRUE : PETSC_FALSE, gassOp;

5951:   PetscCall(MPIU_Allreduce(&assOp, &gassOp, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm)));
5952:   if (hasJac && hasPrec) {
5953:     PetscCall(MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY));
5954:     PetscCall(MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY));
5955:   }
5956: }
5957:   PetscCall(MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY));
5958:   PetscCall(MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY));
5959:   PetscCall(PetscLogEventEnd(DMPLEX_JacobianFEM, dm, 0, 0, 0));
5960:   PetscFunctionReturn(PETSC_SUCCESS);
5961: }

5963: PetscErrorCode DMPlexComputeJacobian_Hybrid_Internal(DM dm, PetscFormKey key[], IS cellIS, PetscReal t, PetscReal X_tShift, Vec locX, Vec locX_t, Mat Jac, Mat JacP, void *user)
5964: {
5965:   DM_Plex        *mesh          = (DM_Plex *)dm->data;
5966:   const char     *name          = "Hybrid Jacobian";
5967:   DM              dmAux[3]      = {NULL, NULL, NULL};
5968:   DMLabel         ghostLabel    = NULL;
5969:   DM              plex          = NULL;
5970:   DM              plexA         = NULL;
5971:   PetscDS         ds            = NULL;
5972:   PetscDS         dsIn          = NULL;
5973:   PetscDS         dsAux[3]      = {NULL, NULL, NULL};
5974:   Vec             locA[3]       = {NULL, NULL, NULL};
5975:   DM              dmScale[3]    = {NULL, NULL, NULL};
5976:   PetscDS         dsScale[3]    = {NULL, NULL, NULL};
5977:   Vec             locS[3]       = {NULL, NULL, NULL};
5978:   PetscSection    section       = NULL;
5979:   PetscSection    sectionAux[3] = {NULL, NULL, NULL};
5980:   DMField         coordField    = NULL;
5981:   PetscScalar    *a[3]          = {NULL, NULL, NULL};
5982:   PetscScalar    *s[3]          = {NULL, NULL, NULL};
5983:   PetscScalar    *u             = NULL, *u_t;
5984:   PetscScalar    *elemMatNeg, *elemMatPos, *elemMatCoh;
5985:   PetscScalar    *elemMatNegP, *elemMatPosP, *elemMatCohP;
5986:   PetscSection    globalSection;
5987:   IS              chunkIS;
5988:   const PetscInt *cells;
5989:   PetscInt       *faces;
5990:   PetscInt        cStart, cEnd, numCells;
5991:   PetscInt        Nf, fieldI, fieldJ, totDim, totDimIn, totDimAux[3], totDimScale[3], numChunks, cellChunkSize, chunk;
5992:   PetscInt        maxDegree  = PETSC_MAX_INT;
5993:   PetscQuadrature affineQuad = NULL, *quads = NULL;
5994:   PetscFEGeom    *affineGeom = NULL, **geoms = NULL;
5995:   PetscBool       hasBdJac, hasBdPrec;

5997:   PetscFunctionBegin;
5998:   PetscCall(PetscLogEventBegin(DMPLEX_JacobianFEM, dm, 0, 0, 0));
5999:   if (!cellIS) goto end;
6000:   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
6001:   PetscCall(ISGetLocalSize(cellIS, &numCells));
6002:   if (cStart >= cEnd) goto end;
6003:   if ((key[0].label == key[1].label) && (key[0].value == key[1].value) && (key[0].part == key[1].part)) {
6004:     const char *name;
6005:     PetscCall(PetscObjectGetName((PetscObject)key[0].label, &name));
6006:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Form keys for each side of a cohesive surface must be different (%s, %" PetscInt_FMT ", %" PetscInt_FMT ")", name, key[0].value, key[0].part);
6007:   }
6008:   PetscCall(DMConvert(dm, DMPLEX, &plex));
6009:   PetscCall(DMGetSection(dm, &section));
6010:   PetscCall(DMGetGlobalSection(dm, &globalSection));
6011:   PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
6012:   PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &ds, &dsIn));
6013:   PetscCall(PetscDSGetNumFields(ds, &Nf));
6014:   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
6015:   PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn));
6016:   PetscCall(PetscDSHasBdJacobian(ds, &hasBdJac));
6017:   PetscCall(PetscDSHasBdJacobianPreconditioner(ds, &hasBdPrec));
6018:   PetscCall(DMGetAuxiliaryVec(dm, key[2].label, key[2].value, key[2].part, &locA[2]));
6019:   if (locA[2]) {
6020:     const PetscInt cellStart = cells ? cells[cStart] : cStart;

6022:     PetscCall(VecGetDM(locA[2], &dmAux[2]));
6023:     PetscCall(DMConvert(dmAux[2], DMPLEX, &plexA));
6024:     PetscCall(DMGetSection(dmAux[2], &sectionAux[2]));
6025:     PetscCall(DMGetCellDS(dmAux[2], cellStart, &dsAux[2], NULL));
6026:     PetscCall(PetscDSGetTotalDimension(dsAux[2], &totDimAux[2]));
6027:     {
6028:       const PetscInt *cone;
6029:       PetscInt        c;

6031:       PetscCall(DMPlexGetCone(dm, cellStart, &cone));
6032:       for (c = 0; c < 2; ++c) {
6033:         const PetscInt *support;
6034:         PetscInt        ssize, s;

6036:         PetscCall(DMPlexGetSupport(dm, cone[c], &support));
6037:         PetscCall(DMPlexGetSupportSize(dm, cone[c], &ssize));
6038:         PetscCheck(ssize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " from cell %" PetscInt_FMT " has support size %" PetscInt_FMT " != 2", cone[c], cellStart, ssize);
6039:         if (support[0] == cellStart) s = 1;
6040:         else if (support[1] == cellStart) s = 0;
6041:         else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " does not have cell %" PetscInt_FMT " in its support", cone[c], cellStart);
6042:         PetscCall(DMGetAuxiliaryVec(dm, key[c].label, key[c].value, key[c].part, &locA[c]));
6043:         if (locA[c]) PetscCall(VecGetDM(locA[c], &dmAux[c]));
6044:         else dmAux[c] = dmAux[2];
6045:         PetscCall(DMGetCellDS(dmAux[c], support[s], &dsAux[c], NULL));
6046:         PetscCall(PetscDSGetTotalDimension(dsAux[c], &totDimAux[c]));
6047:       }
6048:     }
6049:   }
6050:   /* Handle mass matrix scaling
6051:        The field in key[2] is the field to be scaled, and the scaling field is the first in the dsScale */
6052:   PetscCall(DMGetAuxiliaryVec(dm, key[2].label, -key[2].value, key[2].part, &locS[2]));
6053:   if (locS[2]) {
6054:     const PetscInt cellStart = cells ? cells[cStart] : cStart;
6055:     PetscInt       Nb, Nbs;

6057:     PetscCall(VecGetDM(locS[2], &dmScale[2]));
6058:     PetscCall(DMGetCellDS(dmScale[2], cells ? cells[cStart] : cStart, &dsScale[2], NULL));
6059:     PetscCall(PetscDSGetTotalDimension(dsScale[2], &totDimScale[2]));
6060:     // BRAD: This is not set correctly
6061:     key[2].field = 2;
6062:     PetscCall(PetscDSGetFieldSize(ds, key[2].field, &Nb));
6063:     PetscCall(PetscDSGetFieldSize(dsScale[2], 0, &Nbs));
6064:     PetscCheck(Nb == Nbs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field %" PetscInt_FMT " of size %" PetscInt_FMT " cannot be scaled by field of size %" PetscInt_FMT, key[2].field, Nb, Nbs);
6065:     {
6066:       const PetscInt *cone;
6067:       PetscInt        c;

6069:       locS[1] = locS[0] = locS[2];
6070:       dmScale[1] = dmScale[0] = dmScale[2];
6071:       PetscCall(DMPlexGetCone(dm, cellStart, &cone));
6072:       for (c = 0; c < 2; ++c) {
6073:         const PetscInt *support;
6074:         PetscInt        ssize, s;

6076:         PetscCall(DMPlexGetSupport(dm, cone[c], &support));
6077:         PetscCall(DMPlexGetSupportSize(dm, cone[c], &ssize));
6078:         PetscCheck(ssize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " from cell %" PetscInt_FMT " has support size %" PetscInt_FMT " != 2", cone[c], cellStart, ssize);
6079:         if (support[0] == cellStart) s = 1;
6080:         else if (support[1] == cellStart) s = 0;
6081:         else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " does not have cell %" PetscInt_FMT " in its support", cone[c], cellStart);
6082:         PetscCall(DMGetCellDS(dmScale[c], support[s], &dsScale[c], NULL));
6083:         PetscCall(PetscDSGetTotalDimension(dsScale[c], &totDimScale[c]));
6084:       }
6085:     }
6086:   }
6087:   /* 2: Setup geometric data */
6088:   PetscCall(DMGetCoordinateField(dm, &coordField));
6089:   PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
6090:   if (maxDegree > 1) {
6091:     PetscInt f;
6092:     PetscCall(PetscCalloc2(Nf, &quads, Nf, &geoms));
6093:     for (f = 0; f < Nf; ++f) {
6094:       PetscFE fe;

6096:       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
6097:       if (fe) {
6098:         PetscCall(PetscFEGetQuadrature(fe, &quads[f]));
6099:         PetscCall(PetscObjectReference((PetscObject)quads[f]));
6100:       }
6101:     }
6102:   }
6103:   /* Loop over chunks */
6104:   cellChunkSize = numCells;
6105:   numChunks     = !numCells ? 0 : PetscCeilReal(((PetscReal)numCells) / cellChunkSize);
6106:   PetscCall(PetscCalloc1(2 * cellChunkSize, &faces));
6107:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1 * cellChunkSize, faces, PETSC_USE_POINTER, &chunkIS));
6108:   /* Extract field coefficients */
6109:   /* NOTE This needs the end cap faces to have identical orientations */
6110:   PetscCall(DMPlexGetHybridCellFields(dm, cellIS, locX, locX_t, locA[2], &u, &u_t, &a[2]));
6111:   PetscCall(DMPlexGetHybridFields(dm, dmAux, dsAux, cellIS, locA, PETSC_TRUE, a));
6112:   PetscCall(DMPlexGetHybridFields(dm, dmScale, dsScale, cellIS, locS, PETSC_TRUE, s));
6113:   PetscCall(DMGetWorkArray(dm, hasBdJac ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatNeg));
6114:   PetscCall(DMGetWorkArray(dm, hasBdJac ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatPos));
6115:   PetscCall(DMGetWorkArray(dm, hasBdJac ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatCoh));
6116:   PetscCall(DMGetWorkArray(dm, hasBdPrec ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatNegP));
6117:   PetscCall(DMGetWorkArray(dm, hasBdPrec ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatPosP));
6118:   PetscCall(DMGetWorkArray(dm, hasBdPrec ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatCohP));
6119:   for (chunk = 0; chunk < numChunks; ++chunk) {
6120:     PetscInt cS = cStart + chunk * cellChunkSize, cE = PetscMin(cS + cellChunkSize, cEnd), numCells = cE - cS, c;

6122:     if (hasBdJac) {
6123:       PetscCall(PetscArrayzero(elemMatNeg, cellChunkSize * totDim * totDim));
6124:       PetscCall(PetscArrayzero(elemMatPos, cellChunkSize * totDim * totDim));
6125:       PetscCall(PetscArrayzero(elemMatCoh, cellChunkSize * totDim * totDim));
6126:     }
6127:     if (hasBdPrec) {
6128:       PetscCall(PetscArrayzero(elemMatNegP, cellChunkSize * totDim * totDim));
6129:       PetscCall(PetscArrayzero(elemMatPosP, cellChunkSize * totDim * totDim));
6130:       PetscCall(PetscArrayzero(elemMatCohP, cellChunkSize * totDim * totDim));
6131:     }
6132:     /* Get faces */
6133:     for (c = cS; c < cE; ++c) {
6134:       const PetscInt  cell = cells ? cells[c] : c;
6135:       const PetscInt *cone;
6136:       PetscCall(DMPlexGetCone(plex, cell, &cone));
6137:       faces[(c - cS) * 2 + 0] = cone[0];
6138:       faces[(c - cS) * 2 + 1] = cone[1];
6139:     }
6140:     PetscCall(ISGeneralSetIndices(chunkIS, 2 * cellChunkSize, faces, PETSC_USE_POINTER));
6141:     if (maxDegree <= 1) {
6142:       if (!affineQuad) PetscCall(DMFieldCreateDefaultQuadrature(coordField, chunkIS, &affineQuad));
6143:       if (affineQuad) PetscCall(DMSNESGetFEGeom(coordField, chunkIS, affineQuad, PETSC_TRUE, &affineGeom));
6144:     } else {
6145:       PetscInt f;
6146:       for (f = 0; f < Nf; ++f) {
6147:         if (quads[f]) PetscCall(DMSNESGetFEGeom(coordField, chunkIS, quads[f], PETSC_TRUE, &geoms[f]));
6148:       }
6149:     }

6151:     for (fieldI = 0; fieldI < Nf; ++fieldI) {
6152:       PetscFE         feI;
6153:       PetscFEGeom    *geom      = affineGeom ? affineGeom : geoms[fieldI];
6154:       PetscFEGeom    *chunkGeom = NULL, *remGeom = NULL;
6155:       PetscQuadrature quad = affineQuad ? affineQuad : quads[fieldI];
6156:       PetscInt        numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset, Nq, Nb;
6157:       PetscBool       isCohesiveField;

6159:       PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
6160:       if (!feI) continue;
6161:       PetscCall(PetscFEGetTileSizes(feI, NULL, &numBlocks, NULL, &numBatches));
6162:       PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
6163:       PetscCall(PetscFEGetDimension(feI, &Nb));
6164:       blockSize = Nb;
6165:       batchSize = numBlocks * blockSize;
6166:       PetscCall(PetscFESetTileSizes(feI, blockSize, numBlocks, batchSize, numBatches));
6167:       numChunks = numCells / (numBatches * batchSize);
6168:       Ne        = numChunks * numBatches * batchSize;
6169:       Nr        = numCells % (numBatches * batchSize);
6170:       offset    = numCells - Nr;
6171:       PetscCall(PetscFEGeomGetChunk(geom, 0, offset * 2, &chunkGeom));
6172:       PetscCall(PetscFEGeomGetChunk(geom, offset * 2, numCells * 2, &remGeom));
6173:       PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveField));
6174:       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
6175:         PetscFE feJ;

6177:         PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
6178:         if (!feJ) continue;
6179:         key[0].field = fieldI * Nf + fieldJ;
6180:         key[1].field = fieldI * Nf + fieldJ;
6181:         key[2].field = fieldI * Nf + fieldJ;
6182:         if (hasBdJac) {
6183:           PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN, key[0], 0, Ne, chunkGeom, u, u_t, dsAux[0], a[0], t, X_tShift, elemMatNeg));
6184:           PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN, key[0], 0, Nr, remGeom, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[0], PetscSafePointerPlusOffset(a[0], offset * totDimAux[0]), t, X_tShift, &elemMatNeg[offset * totDim * totDim]));
6185:           PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN, key[1], 1, Ne, chunkGeom, u, u_t, dsAux[1], a[1], t, X_tShift, elemMatPos));
6186:           PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN, key[1], 1, Nr, remGeom, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[1], PetscSafePointerPlusOffset(a[1], offset * totDimAux[1]), t, X_tShift, &elemMatPos[offset * totDim * totDim]));
6187:         }
6188:         if (hasBdPrec) {
6189:           PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN_PRE, key[0], 0, Ne, chunkGeom, u, u_t, dsAux[0], a[0], t, X_tShift, elemMatNegP));
6190:           PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN_PRE, key[0], 0, Nr, remGeom, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[0], &a[0][offset * totDimAux[0]], t, X_tShift, &elemMatNegP[offset * totDim * totDim]));
6191:           PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN_PRE, key[1], 1, Ne, chunkGeom, u, u_t, dsAux[1], a[1], t, X_tShift, elemMatPosP));
6192:           PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN_PRE, key[1], 1, Nr, remGeom, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[1], &a[1][offset * totDimAux[1]], t, X_tShift, &elemMatPosP[offset * totDim * totDim]));
6193:         }
6194:         if (hasBdJac) {
6195:           PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN, key[2], 2, Ne, chunkGeom, u, u_t, dsAux[2], a[2], t, X_tShift, elemMatCoh));
6196:           PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN, key[2], 2, Nr, remGeom, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[2], PetscSafePointerPlusOffset(a[2], offset * totDimAux[2]), t, X_tShift, &elemMatCoh[offset * totDim * totDim]));
6197:         }
6198:         if (hasBdPrec) {
6199:           PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN_PRE, key[2], 2, Ne, chunkGeom, u, u_t, dsAux[2], a[2], t, X_tShift, elemMatCohP));
6200:           PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN_PRE, key[2], 2, Nr, remGeom, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[2], &a[2][offset * totDimAux[2]], t, X_tShift, &elemMatCohP[offset * totDim * totDim]));
6201:         }
6202:       }
6203:       PetscCall(PetscFEGeomRestoreChunk(geom, offset, numCells, &remGeom));
6204:       PetscCall(PetscFEGeomRestoreChunk(geom, 0, offset, &chunkGeom));
6205:     }
6206:     /* Insert values into matrix */
6207:     for (c = cS; c < cE; ++c) {
6208:       const PetscInt cell = cells ? cells[c] : c;
6209:       const PetscInt cind = c - cS, coff = cind * totDim * totDim;
6210:       PetscInt       i, j;

6212:       /* Scale element values */
6213:       if (locS[0]) {
6214:         PetscInt  Nb, soff = cind * totDimScale[0], off = 0;
6215:         PetscBool cohesive;

6217:         for (fieldI = 0; fieldI < Nf; ++fieldI) {
6218:           PetscCall(PetscDSGetFieldSize(ds, fieldI, &Nb));
6219:           PetscCall(PetscDSGetCohesive(ds, fieldI, &cohesive));

6221:           if (fieldI == key[2].field) {
6222:             PetscCheck(cohesive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Scaling should not happen for face fields");
6223:             for (i = 0; i < Nb; ++i) {
6224:               for (j = 0; j < totDim; ++j) elemMatCoh[coff + (off + i) * totDim + j] += s[0][soff + i] * elemMatNeg[coff + (off + i) * totDim + j] + s[1][soff + i] * elemMatPos[coff + (off + i) * totDim + j];
6225:               if (hasBdPrec)
6226:                 for (j = 0; j < totDim; ++j) elemMatCohP[coff + (off + i) * totDim + j] += s[0][soff + i] * elemMatNegP[coff + (off + i) * totDim + j] + s[1][soff + i] * elemMatPosP[coff + (off + i) * totDim + j];
6227:             }
6228:             off += Nb;
6229:           } else {
6230:             const PetscInt N = cohesive ? Nb : Nb * 2;

6232:             for (i = 0; i < N; ++i) {
6233:               for (j = 0; j < totDim; ++j) elemMatCoh[coff + (off + i) * totDim + j] += elemMatNeg[coff + (off + i) * totDim + j] + elemMatPos[coff + (off + i) * totDim + j];
6234:               if (hasBdPrec)
6235:                 for (j = 0; j < totDim; ++j) elemMatCohP[coff + (off + i) * totDim + j] += elemMatNegP[coff + (off + i) * totDim + j] + elemMatPosP[coff + (off + i) * totDim + j];
6236:             }
6237:             off += N;
6238:           }
6239:         }
6240:       } else {
6241:         for (i = 0; i < totDim * totDim; ++i) elemMatCoh[coff + i] += elemMatNeg[coff + i] + elemMatPos[coff + i];
6242:         if (hasBdPrec)
6243:           for (i = 0; i < totDim * totDim; ++i) elemMatCohP[coff + i] += elemMatNegP[coff + i] + elemMatPosP[coff + i];
6244:       }
6245:       if (hasBdPrec) {
6246:         if (hasBdJac) {
6247:           if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatCoh[cind * totDim * totDim]));
6248:           PetscCall(DMPlexMatSetClosure_Internal(plex, section, globalSection, mesh->useMatClPerm, Jac, cell, &elemMatCoh[cind * totDim * totDim], ADD_VALUES));
6249:         }
6250:         if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatCohP[cind * totDim * totDim]));
6251:         PetscCall(DMPlexMatSetClosure(plex, section, globalSection, JacP, cell, &elemMatCohP[cind * totDim * totDim], ADD_VALUES));
6252:       } else if (hasBdJac) {
6253:         if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatCoh[cind * totDim * totDim]));
6254:         PetscCall(DMPlexMatSetClosure_Internal(plex, section, globalSection, mesh->useMatClPerm, JacP, cell, &elemMatCoh[cind * totDim * totDim], ADD_VALUES));
6255:       }
6256:     }
6257:   }
6258:   PetscCall(DMPlexRestoreCellFields(dm, cellIS, locX, locX_t, locA[2], &u, &u_t, &a[2]));
6259:   PetscCall(DMPlexRestoreHybridFields(dm, dmAux, dsAux, cellIS, locA, PETSC_TRUE, a));
6260:   PetscCall(DMRestoreWorkArray(dm, hasBdJac ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatNeg));
6261:   PetscCall(DMRestoreWorkArray(dm, hasBdJac ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatPos));
6262:   PetscCall(DMRestoreWorkArray(dm, hasBdJac ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatCoh));
6263:   PetscCall(DMRestoreWorkArray(dm, hasBdPrec ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatNegP));
6264:   PetscCall(DMRestoreWorkArray(dm, hasBdPrec ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatPosP));
6265:   PetscCall(DMRestoreWorkArray(dm, hasBdPrec ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatCohP));
6266:   PetscCall(PetscFree(faces));
6267:   PetscCall(ISDestroy(&chunkIS));
6268:   PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
6269:   if (maxDegree <= 1) {
6270:     PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, affineQuad, PETSC_FALSE, &affineGeom));
6271:     PetscCall(PetscQuadratureDestroy(&affineQuad));
6272:   } else {
6273:     PetscInt f;
6274:     for (f = 0; f < Nf; ++f) {
6275:       if (geoms) PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, quads[f], PETSC_FALSE, &geoms[f]));
6276:       if (quads) PetscCall(PetscQuadratureDestroy(&quads[f]));
6277:     }
6278:     PetscCall(PetscFree2(quads, geoms));
6279:   }
6280:   if (dmAux[2]) PetscCall(DMDestroy(&plexA));
6281:   PetscCall(DMDestroy(&plex));
6282: end:
6283:   PetscCall(PetscLogEventEnd(DMPLEX_JacobianFEM, dm, 0, 0, 0));
6284:   PetscFunctionReturn(PETSC_SUCCESS);
6285: }

6287: /*
6288:   DMPlexComputeJacobian_Action_Internal - Form the local portion of the Jacobian action Z = J(X) Y at the local solution X using pointwise functions specified by the user.

6290:   Input Parameters:
6291: + dm     - The mesh
6292: . key    - The PetscWeakFormKey indicating where integration should happen
6293: . cellIS - The cells to integrate over
6294: . t      - The time
6295: . X_tShift - The multiplier for the Jacobian with respect to X_t
6296: . X      - Local solution vector
6297: . X_t    - Time-derivative of the local solution vector
6298: . Y      - Local input vector
6299: - user   - the user context

6301:   Output Parameter:
6302: . Z - Local output vector

6304:   Note:
6305:   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
6306:   like a GPU, or vectorize on a multicore machine.
6307: */
6308: PetscErrorCode DMPlexComputeJacobian_Action_Internal(DM dm, PetscFormKey key, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Vec Y, Vec Z, void *user)
6309: {
6310:   DM_Plex        *mesh  = (DM_Plex *)dm->data;
6311:   const char     *name  = "Jacobian";
6312:   DM              dmAux = NULL, plex, plexAux = NULL;
6313:   DMEnclosureType encAux;
6314:   Vec             A;
6315:   DMField         coordField;
6316:   PetscDS         prob, probAux = NULL;
6317:   PetscQuadrature quad;
6318:   PetscSection    section, globalSection, sectionAux;
6319:   PetscScalar    *elemMat, *elemMatD, *u, *u_t, *a = NULL, *y, *z;
6320:   const PetscInt *cells;
6321:   PetscInt        Nf, fieldI, fieldJ;
6322:   PetscInt        totDim, totDimAux = 0, cStart, cEnd, numCells, c;
6323:   PetscBool       hasDyn;

6325:   PetscFunctionBegin;
6326:   if (!cellIS) PetscFunctionReturn(PETSC_SUCCESS);
6327:   PetscCall(PetscLogEventBegin(DMPLEX_JacobianFEM, dm, 0, 0, 0));
6328:   PetscCall(DMConvert(dm, DMPLEX, &plex));
6329:   PetscCall(ISGetLocalSize(cellIS, &numCells));
6330:   PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
6331:   PetscCall(DMGetLocalSection(dm, &section));
6332:   PetscCall(DMGetGlobalSection(dm, &globalSection));
6333:   PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob, NULL));
6334:   PetscCall(PetscDSGetNumFields(prob, &Nf));
6335:   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
6336:   PetscCall(PetscDSHasDynamicJacobian(prob, &hasDyn));
6337:   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
6338:   PetscCall(DMGetAuxiliaryVec(dm, key.label, key.value, key.part, &A));
6339:   if (A) {
6340:     PetscCall(VecGetDM(A, &dmAux));
6341:     PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
6342:     PetscCall(DMConvert(dmAux, DMPLEX, &plexAux));
6343:     PetscCall(DMGetLocalSection(plexAux, &sectionAux));
6344:     PetscCall(DMGetDS(dmAux, &probAux));
6345:     PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
6346:   }
6347:   PetscCall(VecSet(Z, 0.0));
6348:   PetscCall(PetscMalloc6(numCells * totDim, &u, X_t ? numCells * totDim : 0, &u_t, numCells * totDim * totDim, &elemMat, hasDyn ? numCells * totDim * totDim : 0, &elemMatD, numCells * totDim, &y, totDim, &z));
6349:   if (dmAux) PetscCall(PetscMalloc1(numCells * totDimAux, &a));
6350:   PetscCall(DMGetCoordinateField(dm, &coordField));
6351:   for (c = cStart; c < cEnd; ++c) {
6352:     const PetscInt cell = cells ? cells[c] : c;
6353:     const PetscInt cind = c - cStart;
6354:     PetscScalar   *x = NULL, *x_t = NULL;
6355:     PetscInt       i;

6357:     PetscCall(DMPlexVecGetClosure(plex, section, X, cell, NULL, &x));
6358:     for (i = 0; i < totDim; ++i) u[cind * totDim + i] = x[i];
6359:     PetscCall(DMPlexVecRestoreClosure(plex, section, X, cell, NULL, &x));
6360:     if (X_t) {
6361:       PetscCall(DMPlexVecGetClosure(plex, section, X_t, cell, NULL, &x_t));
6362:       for (i = 0; i < totDim; ++i) u_t[cind * totDim + i] = x_t[i];
6363:       PetscCall(DMPlexVecRestoreClosure(plex, section, X_t, cell, NULL, &x_t));
6364:     }
6365:     if (dmAux) {
6366:       PetscInt subcell;
6367:       PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell));
6368:       PetscCall(DMPlexVecGetClosure(plexAux, sectionAux, A, subcell, NULL, &x));
6369:       for (i = 0; i < totDimAux; ++i) a[cind * totDimAux + i] = x[i];
6370:       PetscCall(DMPlexVecRestoreClosure(plexAux, sectionAux, A, subcell, NULL, &x));
6371:     }
6372:     PetscCall(DMPlexVecGetClosure(plex, section, Y, cell, NULL, &x));
6373:     for (i = 0; i < totDim; ++i) y[cind * totDim + i] = x[i];
6374:     PetscCall(DMPlexVecRestoreClosure(plex, section, Y, cell, NULL, &x));
6375:   }
6376:   PetscCall(PetscArrayzero(elemMat, numCells * totDim * totDim));
6377:   if (hasDyn) PetscCall(PetscArrayzero(elemMatD, numCells * totDim * totDim));
6378:   for (fieldI = 0; fieldI < Nf; ++fieldI) {
6379:     PetscFE  fe;
6380:     PetscInt Nb;
6381:     /* Conforming batches */
6382:     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
6383:     /* Remainder */
6384:     PetscInt        Nr, offset, Nq;
6385:     PetscQuadrature qGeom = NULL;
6386:     PetscInt        maxDegree;
6387:     PetscFEGeom    *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;

6389:     PetscCall(PetscDSGetDiscretization(prob, fieldI, (PetscObject *)&fe));
6390:     PetscCall(PetscFEGetQuadrature(fe, &quad));
6391:     PetscCall(PetscFEGetDimension(fe, &Nb));
6392:     PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
6393:     PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
6394:     if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &qGeom));
6395:     if (!qGeom) {
6396:       PetscCall(PetscFEGetQuadrature(fe, &qGeom));
6397:       PetscCall(PetscObjectReference((PetscObject)qGeom));
6398:     }
6399:     PetscCall(PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL));
6400:     PetscCall(DMSNESGetFEGeom(coordField, cellIS, qGeom, PETSC_FALSE, &cgeomFEM));
6401:     blockSize = Nb;
6402:     batchSize = numBlocks * blockSize;
6403:     PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
6404:     numChunks = numCells / (numBatches * batchSize);
6405:     Ne        = numChunks * numBatches * batchSize;
6406:     Nr        = numCells % (numBatches * batchSize);
6407:     offset    = numCells - Nr;
6408:     PetscCall(PetscFEGeomGetChunk(cgeomFEM, 0, offset, &chunkGeom));
6409:     PetscCall(PetscFEGeomGetChunk(cgeomFEM, offset, numCells, &remGeom));
6410:     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
6411:       key.field = fieldI * Nf + fieldJ;
6412:       PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat));
6413:       PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMat[offset * totDim * totDim]));
6414:       if (hasDyn) {
6415:         PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD));
6416:         PetscCall(PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, &a[offset * totDimAux], t, X_tShift, &elemMatD[offset * totDim * totDim]));
6417:       }
6418:     }
6419:     PetscCall(PetscFEGeomRestoreChunk(cgeomFEM, offset, numCells, &remGeom));
6420:     PetscCall(PetscFEGeomRestoreChunk(cgeomFEM, 0, offset, &chunkGeom));
6421:     PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, qGeom, PETSC_FALSE, &cgeomFEM));
6422:     PetscCall(PetscQuadratureDestroy(&qGeom));
6423:   }
6424:   if (hasDyn) {
6425:     for (c = 0; c < numCells * totDim * totDim; ++c) elemMat[c] += X_tShift * elemMatD[c];
6426:   }
6427:   for (c = cStart; c < cEnd; ++c) {
6428:     const PetscInt     cell = cells ? cells[c] : c;
6429:     const PetscInt     cind = c - cStart;
6430:     const PetscBLASInt M = totDim, one = 1;
6431:     const PetscScalar  a = 1.0, b = 0.0;

6433:     PetscCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[cind * totDim * totDim], &M, &y[cind * totDim], &one, &b, z, &one));
6434:     if (mesh->printFEM > 1) {
6435:       PetscCall(DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[cind * totDim * totDim]));
6436:       PetscCall(DMPrintCellVector(c, "Y", totDim, &y[cind * totDim]));
6437:       PetscCall(DMPrintCellVector(c, "Z", totDim, z));
6438:     }
6439:     PetscCall(DMPlexVecSetClosure(dm, section, Z, cell, z, ADD_VALUES));
6440:   }
6441:   PetscCall(PetscFree6(u, u_t, elemMat, elemMatD, y, z));
6442:   if (mesh->printFEM) {
6443:     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)Z), "Z:\n"));
6444:     PetscCall(VecView(Z, NULL));
6445:   }
6446:   PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
6447:   PetscCall(PetscFree(a));
6448:   PetscCall(DMDestroy(&plexAux));
6449:   PetscCall(DMDestroy(&plex));
6450:   PetscCall(PetscLogEventEnd(DMPLEX_JacobianFEM, dm, 0, 0, 0));
6451:   PetscFunctionReturn(PETSC_SUCCESS);
6452: }