Actual source code: fe.c

  1: /* Basis Jet Tabulation

  3: We would like to tabulate the nodal basis functions and derivatives at a set of points, usually quadrature points. We
  4: follow here the derviation in http://www.math.ttu.edu/~kirby/papers/fiat-toms-2004.pdf. The nodal basis $\psi_i$ can
  5: be expressed in terms of a prime basis $\phi_i$ which can be stably evaluated. In PETSc, we will use the Legendre basis
  6: as a prime basis.

  8:   \psi_i = \sum_k \alpha_{ki} \phi_k

 10: Our nodal basis is defined in terms of the dual basis $n_j$

 12:   n_j \cdot \psi_i = \delta_{ji}

 14: and we may act on the first equation to obtain

 16:   n_j \cdot \psi_i = \sum_k \alpha_{ki} n_j \cdot \phi_k
 17:        \delta_{ji} = \sum_k \alpha_{ki} V_{jk}
 18:                  I = V \alpha

 20: so the coefficients of the nodal basis in the prime basis are

 22:    \alpha = V^{-1}

 24: We will define the dual basis vectors $n_j$ using a quadrature rule.

 26: Right now, we will just use the polynomial spaces P^k. I know some elements use the space of symmetric polynomials
 27: (I think Nedelec), but we will neglect this for now. Constraints in the space, e.g. Arnold-Winther elements, can
 28: be implemented exactly as in FIAT using functionals $L_j$.

 30: I will have to count the degrees correctly for the Legendre product when we are on simplices.

 32: We will have three objects:
 33:  - Space, P: this just need point evaluation I think
 34:  - Dual Space, P'+K: This looks like a set of functionals that can act on members of P, each n is defined by a Q
 35:  - FEM: This keeps {P, P', Q}
 36: */
 37: #include <petsc/private/petscfeimpl.h>
 38: #include <petscdmplex.h>

 40: PetscBool  FEcite       = PETSC_FALSE;
 41: const char FECitation[] = "@article{kirby2004,\n"
 42:                           "  title   = {Algorithm 839: FIAT, a New Paradigm for Computing Finite Element Basis Functions},\n"
 43:                           "  journal = {ACM Transactions on Mathematical Software},\n"
 44:                           "  author  = {Robert C. Kirby},\n"
 45:                           "  volume  = {30},\n"
 46:                           "  number  = {4},\n"
 47:                           "  pages   = {502--516},\n"
 48:                           "  doi     = {10.1145/1039813.1039820},\n"
 49:                           "  year    = {2004}\n}\n";

 51: PetscClassId PETSCFE_CLASSID = 0;

 53: PetscLogEvent PETSCFE_SetUp;

 55: PetscFunctionList PetscFEList              = NULL;
 56: PetscBool         PetscFERegisterAllCalled = PETSC_FALSE;

 58: /*@C
 59:   PetscFERegister - Adds a new `PetscFEType`

 61:   Not Collective, No Fortran Support

 63:   Input Parameters:
 64: + sname    - The name of a new user-defined creation routine
 65: - function - The creation routine

 67:   Example Usage:
 68: .vb
 69:     PetscFERegister("my_fe", MyPetscFECreate);
 70: .ve

 72:   Then, your PetscFE type can be chosen with the procedural interface via
 73: .vb
 74:     PetscFECreate(MPI_Comm, PetscFE *);
 75:     PetscFESetType(PetscFE, "my_fe");
 76: .ve
 77:   or at runtime via the option
 78: .vb
 79:     -petscfe_type my_fe
 80: .ve

 82:   Level: advanced

 84:   Note:
 85:   `PetscFERegister()` may be called multiple times to add several user-defined `PetscFE`s

 87: .seealso: `PetscFE`, `PetscFEType`, `PetscFERegisterAll()`, `PetscFERegisterDestroy()`
 88: @*/
 89: PetscErrorCode PetscFERegister(const char sname[], PetscErrorCode (*function)(PetscFE))
 90: {
 91:   PetscFunctionBegin;
 92:   PetscCall(PetscFunctionListAdd(&PetscFEList, sname, function));
 93:   PetscFunctionReturn(PETSC_SUCCESS);
 94: }

 96: /*@
 97:   PetscFESetType - Builds a particular `PetscFE`

 99:   Collective

101:   Input Parameters:
102: + fem  - The `PetscFE` object
103: - name - The kind of FEM space

105:   Options Database Key:
106: . -petscfe_type type - Sets the `PetscFE` type; use -help for a list of available types

108:   Level: intermediate

110: .seealso: `PetscFEType`, `PetscFE`, `PetscFEGetType()`, `PetscFECreate()`
111: @*/
112: PetscErrorCode PetscFESetType(PetscFE fem, PetscFEType name)
113: {
114:   PetscErrorCode (*r)(PetscFE);
115:   PetscBool match;

117:   PetscFunctionBegin;
119:   PetscCall(PetscObjectTypeCompare((PetscObject)fem, name, &match));
120:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

122:   if (!PetscFERegisterAllCalled) PetscCall(PetscFERegisterAll());
123:   PetscCall(PetscFunctionListFind(PetscFEList, name, &r));
124:   PetscCheck(r, PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFE type: %s", name);

126:   PetscTryTypeMethod(fem, destroy);
127:   fem->ops->destroy = NULL;

129:   PetscCall((*r)(fem));
130:   PetscCall(PetscObjectChangeTypeName((PetscObject)fem, name));
131:   PetscFunctionReturn(PETSC_SUCCESS);
132: }

134: /*@
135:   PetscFEGetType - Gets the `PetscFEType` (as a string) from the `PetscFE` object.

137:   Not Collective

139:   Input Parameter:
140: . fem - The `PetscFE`

142:   Output Parameter:
143: . name - The `PetscFEType` name

145:   Level: intermediate

147: .seealso: `PetscFEType`, `PetscFE`, `PetscFESetType()`, `PetscFECreate()`
148: @*/
149: PetscErrorCode PetscFEGetType(PetscFE fem, PetscFEType *name)
150: {
151:   PetscFunctionBegin;
153:   PetscAssertPointer(name, 2);
154:   if (!PetscFERegisterAllCalled) PetscCall(PetscFERegisterAll());
155:   *name = ((PetscObject)fem)->type_name;
156:   PetscFunctionReturn(PETSC_SUCCESS);
157: }

159: /*@
160:   PetscFEViewFromOptions - View from a `PetscFE` based on values in the options database

162:   Collective

164:   Input Parameters:
165: + A    - the `PetscFE` object
166: . obj  - Optional object that provides the options prefix, pass `NULL` to use the options prefix of `A`
167: - name - command line option name

169:   Options Database Key:
170: . -name [viewertype][:...] - option name and values. See `PetscObjectViewFromOptions()` for the possible arguments

172:   Level: intermediate

174: .seealso: `PetscFE`, `PetscFEView()`, `PetscObjectViewFromOptions()`, `PetscFECreate()`
175: @*/
176: PetscErrorCode PetscFEViewFromOptions(PetscFE A, PeOp PetscObject obj, const char name[])
177: {
178:   PetscFunctionBegin;
180:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
181:   PetscFunctionReturn(PETSC_SUCCESS);
182: }

184: /*@
185:   PetscFEView - Views a `PetscFE`

187:   Collective

189:   Input Parameters:
190: + fem    - the `PetscFE` object to view
191: - viewer - the viewer

193:   Level: beginner

195: .seealso: `PetscFE`, `PetscViewer`, `PetscFEDestroy()`, `PetscFEViewFromOptions()`
196: @*/
197: PetscErrorCode PetscFEView(PetscFE fem, PetscViewer viewer)
198: {
199:   PetscBool isascii;

201:   PetscFunctionBegin;
204:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)fem), &viewer));
205:   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)fem, viewer));
206:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
207:   PetscTryTypeMethod(fem, view, viewer);
208:   PetscFunctionReturn(PETSC_SUCCESS);
209: }

211: /*@
212:   PetscFESetFromOptions - sets parameters in a `PetscFE` from the options database

214:   Collective

216:   Input Parameter:
217: . fem - the `PetscFE` object to set options for

219:   Options Database Keys:
220: + -petscfe_num_blocks  - the number of cell blocks to integrate concurrently
221: - -petscfe_num_batches - the number of cell batches to integrate serially

223:   Level: intermediate

225: .seealso: `PetscFE`, `PetscFEView()`
226: @*/
227: PetscErrorCode PetscFESetFromOptions(PetscFE fem)
228: {
229:   const char *defaultType;
230:   char        name[256];
231:   PetscBool   flg;

233:   PetscFunctionBegin;
235:   if (!((PetscObject)fem)->type_name) defaultType = PETSCFEBASIC;
236:   else defaultType = ((PetscObject)fem)->type_name;
237:   if (!PetscFERegisterAllCalled) PetscCall(PetscFERegisterAll());

239:   PetscObjectOptionsBegin((PetscObject)fem);
240:   PetscCall(PetscOptionsFList("-petscfe_type", "Finite element space", "PetscFESetType", PetscFEList, defaultType, name, 256, &flg));
241:   if (flg) PetscCall(PetscFESetType(fem, name));
242:   else if (!((PetscObject)fem)->type_name) PetscCall(PetscFESetType(fem, defaultType));
243:   PetscCall(PetscOptionsBoundedInt("-petscfe_num_blocks", "The number of cell blocks to integrate concurrently", "PetscSpaceSetTileSizes", fem->numBlocks, &fem->numBlocks, NULL, 1));
244:   PetscCall(PetscOptionsBoundedInt("-petscfe_num_batches", "The number of cell batches to integrate serially", "PetscSpaceSetTileSizes", fem->numBatches, &fem->numBatches, NULL, 1));
245:   PetscTryTypeMethod(fem, setfromoptions, PetscOptionsObject);
246:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
247:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)fem, PetscOptionsObject));
248:   PetscOptionsEnd();
249:   PetscCall(PetscFEViewFromOptions(fem, NULL, "-petscfe_view"));
250:   PetscFunctionReturn(PETSC_SUCCESS);
251: }

253: /*@
254:   PetscFESetUp - Construct data structures for the `PetscFE` after the `PetscFEType` has been set

256:   Collective

258:   Input Parameter:
259: . fem - the `PetscFE` object to setup

261:   Level: intermediate

263: .seealso: `PetscFE`, `PetscFEView()`, `PetscFEDestroy()`
264: @*/
265: PetscErrorCode PetscFESetUp(PetscFE fem)
266: {
267:   PetscFunctionBegin;
269:   if (fem->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
270:   PetscCall(PetscLogEventBegin(PETSCFE_SetUp, fem, 0, 0, 0));
271:   fem->setupcalled = PETSC_TRUE;
272:   PetscTryTypeMethod(fem, setup);
273:   PetscCall(PetscLogEventEnd(PETSCFE_SetUp, fem, 0, 0, 0));
274:   PetscFunctionReturn(PETSC_SUCCESS);
275: }

277: /*@
278:   PetscFEDestroy - Destroys a `PetscFE` object

280:   Collective

282:   Input Parameter:
283: . fem - the `PetscFE` object to destroy

285:   Level: beginner

287: .seealso: `PetscFE`, `PetscFEView()`
288: @*/
289: PetscErrorCode PetscFEDestroy(PetscFE *fem)
290: {
291:   PetscFunctionBegin;
292:   if (!*fem) PetscFunctionReturn(PETSC_SUCCESS);

295:   if (--((PetscObject)*fem)->refct > 0) {
296:     *fem = NULL;
297:     PetscFunctionReturn(PETSC_SUCCESS);
298:   }
299:   ((PetscObject)*fem)->refct = 0;

301:   if ((*fem)->subspaces) {
302:     PetscInt dim;

304:     PetscCall(PetscDualSpaceGetDimension((*fem)->dualSpace, &dim));
305:     for (PetscInt d = 0; d < dim; ++d) PetscCall(PetscFEDestroy(&(*fem)->subspaces[d]));
306:   }
307:   PetscCall(PetscFree((*fem)->subspaces));
308:   PetscCall(PetscFree((*fem)->invV));
309:   PetscCall(PetscTabulationDestroy(&(*fem)->T));
310:   PetscCall(PetscTabulationDestroy(&(*fem)->Tf));
311:   PetscCall(PetscTabulationDestroy(&(*fem)->Tc));
312:   PetscCall(PetscSpaceDestroy(&(*fem)->basisSpace));
313:   PetscCall(PetscDualSpaceDestroy(&(*fem)->dualSpace));
314:   PetscCall(PetscQuadratureDestroy(&(*fem)->quadrature));
315:   PetscCall(PetscQuadratureDestroy(&(*fem)->faceQuadrature));
316: #ifdef PETSC_HAVE_LIBCEED
317:   PetscCallCEED(CeedBasisDestroy(&(*fem)->ceedBasis));
318:   PetscCallCEED(CeedDestroy(&(*fem)->ceed));
319: #endif

321:   PetscTryTypeMethod(*fem, destroy);
322:   PetscCall(PetscHeaderDestroy(fem));
323:   PetscFunctionReturn(PETSC_SUCCESS);
324: }

326: /*@
327:   PetscFECreate - Creates an empty `PetscFE` object. The type can then be set with `PetscFESetType()`.

329:   Collective

331:   Input Parameter:
332: . comm - The communicator for the `PetscFE` object

334:   Output Parameter:
335: . fem - The `PetscFE` object

337:   Level: beginner

339: .seealso: `PetscFE`, `PetscFEType`, `PetscFESetType()`, `PetscFECreateDefault()`, `PETSCFEGALERKIN`
340: @*/
341: PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem)
342: {
343:   PetscFE f;

345:   PetscFunctionBegin;
346:   PetscAssertPointer(fem, 2);
347:   PetscCall(PetscCitationsRegister(FECitation, &FEcite));
348:   PetscCall(PetscFEInitializePackage());

350:   PetscCall(PetscHeaderCreate(f, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView));

352:   f->basisSpace    = NULL;
353:   f->dualSpace     = NULL;
354:   f->numComponents = 1;
355:   f->subspaces     = NULL;
356:   f->invV          = NULL;
357:   f->T             = NULL;
358:   f->Tf            = NULL;
359:   f->Tc            = NULL;
360:   PetscCall(PetscArrayzero(&f->quadrature, 1));
361:   PetscCall(PetscArrayzero(&f->faceQuadrature, 1));
362:   f->blockSize  = 0;
363:   f->numBlocks  = 1;
364:   f->batchSize  = 0;
365:   f->numBatches = 1;

367:   *fem = f;
368:   PetscFunctionReturn(PETSC_SUCCESS);
369: }

371: /*@
372:   PetscFEGetSpatialDimension - Returns the spatial dimension of the element

374:   Not Collective

376:   Input Parameter:
377: . fem - The `PetscFE` object

379:   Output Parameter:
380: . dim - The spatial dimension

382:   Level: intermediate

384: .seealso: `PetscFE`, `PetscFECreate()`
385: @*/
386: PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim)
387: {
388:   DM dm;

390:   PetscFunctionBegin;
392:   PetscAssertPointer(dim, 2);
393:   PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm));
394:   PetscCall(DMGetDimension(dm, dim));
395:   PetscFunctionReturn(PETSC_SUCCESS);
396: }

398: /*@
399:   PetscFESetNumComponents - Sets the number of field components in the element

401:   Not Collective

403:   Input Parameters:
404: + fem  - The `PetscFE` object
405: - comp - The number of field components

407:   Level: intermediate

409: .seealso: `PetscFE`, `PetscFECreate()`, `PetscFEGetSpatialDimension()`, `PetscFEGetNumComponents()`
410: @*/
411: PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp)
412: {
413:   PetscFunctionBegin;
415:   fem->numComponents = comp;
416:   PetscFunctionReturn(PETSC_SUCCESS);
417: }

419: /*@
420:   PetscFEGetNumComponents - Returns the number of components in the element

422:   Not Collective

424:   Input Parameter:
425: . fem - The `PetscFE` object

427:   Output Parameter:
428: . comp - The number of field components

430:   Level: intermediate

432: .seealso: `PetscFE`, `PetscFECreate()`, `PetscFEGetSpatialDimension()`
433: @*/
434: PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp)
435: {
436:   PetscFunctionBegin;
438:   PetscAssertPointer(comp, 2);
439:   *comp = fem->numComponents;
440:   PetscFunctionReturn(PETSC_SUCCESS);
441: }

443: /*@
444:   PetscFESetTileSizes - Sets the tile sizes for evaluation

446:   Not Collective

448:   Input Parameters:
449: + fem        - The `PetscFE` object
450: . blockSize  - The number of elements in a block
451: . numBlocks  - The number of blocks in a batch
452: . batchSize  - The number of elements in a batch
453: - numBatches - The number of batches in a chunk

455:   Level: intermediate

457: .seealso: `PetscFE`, `PetscFECreate()`, `PetscFEGetTileSizes()`
458: @*/
459: PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches)
460: {
461:   PetscFunctionBegin;
463:   fem->blockSize  = blockSize;
464:   fem->numBlocks  = numBlocks;
465:   fem->batchSize  = batchSize;
466:   fem->numBatches = numBatches;
467:   PetscFunctionReturn(PETSC_SUCCESS);
468: }

470: /*@
471:   PetscFEGetTileSizes - Returns the tile sizes for evaluation

473:   Not Collective

475:   Input Parameter:
476: . fem - The `PetscFE` object

478:   Output Parameters:
479: + blockSize  - The number of elements in a block, pass `NULL` if not needed
480: . numBlocks  - The number of blocks in a batch, pass `NULL` if not needed
481: . batchSize  - The number of elements in a batch, pass `NULL` if not needed
482: - numBatches - The number of batches in a chunk, pass `NULL` if not needed

484:   Level: intermediate

486: .seealso: `PetscFE`, `PetscFECreate()`, `PetscFESetTileSizes()`
487: @*/
488: PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PeOp PetscInt *blockSize, PeOp PetscInt *numBlocks, PeOp PetscInt *batchSize, PeOp PetscInt *numBatches)
489: {
490:   PetscFunctionBegin;
492:   if (blockSize) PetscAssertPointer(blockSize, 2);
493:   if (numBlocks) PetscAssertPointer(numBlocks, 3);
494:   if (batchSize) PetscAssertPointer(batchSize, 4);
495:   if (numBatches) PetscAssertPointer(numBatches, 5);
496:   if (blockSize) *blockSize = fem->blockSize;
497:   if (numBlocks) *numBlocks = fem->numBlocks;
498:   if (batchSize) *batchSize = fem->batchSize;
499:   if (numBatches) *numBatches = fem->numBatches;
500:   PetscFunctionReturn(PETSC_SUCCESS);
501: }

503: /*@
504:   PetscFEGetBasisSpace - Returns the `PetscSpace` used for the approximation of the solution for the `PetscFE`

506:   Not Collective

508:   Input Parameter:
509: . fem - The `PetscFE` object

511:   Output Parameter:
512: . sp - The `PetscSpace` object

514:   Level: intermediate

516: .seealso: `PetscFE`, `PetscSpace`, `PetscFECreate()`
517: @*/
518: PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp)
519: {
520:   PetscFunctionBegin;
522:   PetscAssertPointer(sp, 2);
523:   *sp = fem->basisSpace;
524:   PetscFunctionReturn(PETSC_SUCCESS);
525: }

527: /*@
528:   PetscFESetBasisSpace - Sets the `PetscSpace` used for the approximation of the solution

530:   Not Collective

532:   Input Parameters:
533: + fem - The `PetscFE` object
534: - sp  - The `PetscSpace` object

536:   Level: intermediate

538:   Developer Notes:
539:   There is `PetscFESetBasisSpace()` but the `PetscFESetDualSpace()`, likely the Basis is unneeded in the function name

541: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()`, `PetscFESetDualSpace()`
542: @*/
543: PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp)
544: {
545:   PetscFunctionBegin;
548:   PetscCall(PetscSpaceDestroy(&fem->basisSpace));
549:   fem->basisSpace = sp;
550:   PetscCall(PetscObjectReference((PetscObject)fem->basisSpace));
551:   PetscFunctionReturn(PETSC_SUCCESS);
552: }

554: /*@
555:   PetscFEGetDualSpace - Returns the `PetscDualSpace` used to define the inner product for a `PetscFE`

557:   Not Collective

559:   Input Parameter:
560: . fem - The `PetscFE` object

562:   Output Parameter:
563: . sp - The `PetscDualSpace` object

565:   Level: intermediate

567: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()`
568: @*/
569: PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp)
570: {
571:   PetscFunctionBegin;
573:   PetscAssertPointer(sp, 2);
574:   *sp = fem->dualSpace;
575:   PetscFunctionReturn(PETSC_SUCCESS);
576: }

578: /*@
579:   PetscFESetDualSpace - Sets the `PetscDualSpace` used to define the inner product

581:   Not Collective

583:   Input Parameters:
584: + fem - The `PetscFE` object
585: - sp  - The `PetscDualSpace` object

587:   Level: intermediate

589: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()`, `PetscFESetBasisSpace()`
590: @*/
591: PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp)
592: {
593:   PetscFunctionBegin;
596:   PetscCall(PetscDualSpaceDestroy(&fem->dualSpace));
597:   fem->dualSpace = sp;
598:   PetscCall(PetscObjectReference((PetscObject)fem->dualSpace));
599:   PetscFunctionReturn(PETSC_SUCCESS);
600: }

602: /*@
603:   PetscFEGetQuadrature - Returns the `PetscQuadrature` used to calculate inner products

605:   Not Collective

607:   Input Parameter:
608: . fem - The `PetscFE` object

610:   Output Parameter:
611: . q - The `PetscQuadrature` object

613:   Level: intermediate

615: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`
616: @*/
617: PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q)
618: {
619:   PetscFunctionBegin;
621:   PetscAssertPointer(q, 2);
622:   *q = fem->quadrature;
623:   PetscFunctionReturn(PETSC_SUCCESS);
624: }

626: /*@
627:   PetscFESetQuadrature - Sets the `PetscQuadrature` used to calculate inner products

629:   Not Collective

631:   Input Parameters:
632: + fem - The `PetscFE` object
633: - q   - The `PetscQuadrature` object

635:   Level: intermediate

637: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFEGetFaceQuadrature()`
638: @*/
639: PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q)
640: {
641:   PetscInt Nc, qNc;

643:   PetscFunctionBegin;
645:   if (q == fem->quadrature) PetscFunctionReturn(PETSC_SUCCESS);
646:   PetscCall(PetscFEGetNumComponents(fem, &Nc));
647:   PetscCall(PetscQuadratureGetNumComponents(q, &qNc));
648:   PetscCheck(!(qNc != 1) || !(Nc != qNc), PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_SIZ, "FE components %" PetscInt_FMT " != Quadrature components %" PetscInt_FMT " and non-scalar quadrature", Nc, qNc);
649:   PetscCall(PetscTabulationDestroy(&fem->T));
650:   PetscCall(PetscTabulationDestroy(&fem->Tc));
651:   PetscCall(PetscObjectReference((PetscObject)q));
652:   PetscCall(PetscQuadratureDestroy(&fem->quadrature));
653:   fem->quadrature = q;
654:   PetscFunctionReturn(PETSC_SUCCESS);
655: }

657: /*@
658:   PetscFEGetFaceQuadrature - Returns the `PetscQuadrature` used to calculate inner products on faces

660:   Not Collective

662:   Input Parameter:
663: . fem - The `PetscFE` object

665:   Output Parameter:
666: . q - The `PetscQuadrature` object

668:   Level: intermediate

670:   Developer Notes:
671:   There is a special face quadrature but not edge, likely this API would benefit from a refactorization

673: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()`, `PetscFESetFaceQuadrature()`
674: @*/
675: PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q)
676: {
677:   PetscFunctionBegin;
679:   PetscAssertPointer(q, 2);
680:   *q = fem->faceQuadrature;
681:   PetscFunctionReturn(PETSC_SUCCESS);
682: }

684: /*@
685:   PetscFESetFaceQuadrature - Sets the `PetscQuadrature` used to calculate inner products on faces

687:   Not Collective

689:   Input Parameters:
690: + fem - The `PetscFE` object
691: - q   - The `PetscQuadrature` object

693:   Level: intermediate

695: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()`
696: @*/
697: PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q)
698: {
699:   PetscInt Nc, qNc;

701:   PetscFunctionBegin;
703:   if (q == fem->faceQuadrature) PetscFunctionReturn(PETSC_SUCCESS);
704:   PetscCall(PetscFEGetNumComponents(fem, &Nc));
705:   PetscCall(PetscQuadratureGetNumComponents(q, &qNc));
706:   PetscCheck(!(qNc != 1) || !(Nc != qNc), PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_SIZ, "FE components %" PetscInt_FMT " != Quadrature components %" PetscInt_FMT " and non-scalar quadrature", Nc, qNc);
707:   PetscCall(PetscTabulationDestroy(&fem->Tf));
708:   PetscCall(PetscObjectReference((PetscObject)q));
709:   PetscCall(PetscQuadratureDestroy(&fem->faceQuadrature));
710:   fem->faceQuadrature = q;
711:   PetscFunctionReturn(PETSC_SUCCESS);
712: }

714: /*@
715:   PetscFECopyQuadrature - Copy both volumetric and surface quadrature to a new `PetscFE`

717:   Not Collective

719:   Input Parameters:
720: + sfe - The `PetscFE` source for the quadratures
721: - tfe - The `PetscFE` target for the quadratures

723:   Level: intermediate

725: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()`, `PetscFESetFaceQuadrature()`
726: @*/
727: PetscErrorCode PetscFECopyQuadrature(PetscFE sfe, PetscFE tfe)
728: {
729:   PetscQuadrature q;

731:   PetscFunctionBegin;
734:   PetscCall(PetscFEGetQuadrature(sfe, &q));
735:   PetscCall(PetscFESetQuadrature(tfe, q));
736:   PetscCall(PetscFEGetFaceQuadrature(sfe, &q));
737:   PetscCall(PetscFESetFaceQuadrature(tfe, q));
738:   PetscFunctionReturn(PETSC_SUCCESS);
739: }

741: /*@C
742:   PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension

744:   Not Collective

746:   Input Parameter:
747: . fem - The `PetscFE` object

749:   Output Parameter:
750: . numDof - Array of length `dim` with the number of dofs in each dimension

752:   Level: intermediate

754: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()`
755: @*/
756: PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt *numDof[])
757: {
758:   PetscFunctionBegin;
760:   PetscAssertPointer(numDof, 2);
761:   PetscCall(PetscDualSpaceGetNumDof(fem->dualSpace, numDof));
762:   PetscFunctionReturn(PETSC_SUCCESS);
763: }

765: /*@C
766:   PetscFEGetCellTabulation - Returns the tabulation of the basis functions at the quadrature points on the reference cell

768:   Not Collective

770:   Input Parameters:
771: + fem - The `PetscFE` object
772: - k   - The highest derivative we need to tabulate, very often 1

774:   Output Parameter:
775: . T - The basis function values and derivatives at quadrature points

777:   Level: intermediate

779:   Note:
780: .vb
781:   T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
782:   T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
783:   T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e
784: .ve

786: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()`
787: @*/
788: PetscErrorCode PetscFEGetCellTabulation(PetscFE fem, PetscInt k, PetscTabulation *T)
789: {
790:   PetscInt         npoints;
791:   const PetscReal *points;

793:   PetscFunctionBegin;
795:   PetscAssertPointer(T, 3);
796:   PetscCall(PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL));
797:   if (!fem->T) PetscCall(PetscFECreateTabulation(fem, 1, npoints, points, k, &fem->T));
798:   PetscCheck(!fem->T || k <= fem->T->K || (!fem->T->cdim && !fem->T->K), PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_OUTOFRANGE, "Requested %" PetscInt_FMT " derivatives, but only tabulated %" PetscInt_FMT, k, fem->T->K);
799:   *T = fem->T;
800:   PetscFunctionReturn(PETSC_SUCCESS);
801: }

803: PetscErrorCode PetscFEExpandFaceQuadrature(PetscFE fe, PetscQuadrature fq, PetscQuadrature *efq)
804: {
805:   DM               dm;
806:   PetscDualSpace   sp;
807:   const PetscInt  *faces;
808:   const PetscReal *points, *weights;
809:   DMPolytopeType   ct;
810:   PetscReal       *facePoints, *faceWeights;
811:   PetscInt         dim, cStart, Nf, Nc, Np, order;

813:   PetscFunctionBegin;
814:   PetscCall(PetscFEGetDualSpace(fe, &sp));
815:   PetscCall(PetscDualSpaceGetDM(sp, &dm));
816:   PetscCall(DMGetDimension(dm, &dim));
817:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
818:   PetscCall(DMPlexGetConeSize(dm, cStart, &Nf));
819:   PetscCall(DMPlexGetCone(dm, cStart, &faces));
820:   PetscCall(PetscQuadratureGetData(fq, NULL, &Nc, &Np, &points, &weights));
821:   PetscCall(PetscMalloc1(Nf * Np * dim, &facePoints));
822:   PetscCall(PetscMalloc1(Nf * Np * Nc, &faceWeights));
823:   for (PetscInt f = 0; f < Nf; ++f) {
824:     const PetscReal xi0[3] = {-1., -1., -1.};
825:     PetscReal       v0[3], J[9], detJ;

827:     PetscCall(DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ));
828:     for (PetscInt q = 0; q < Np; ++q) {
829:       CoordinatesRefToReal(dim, dim - 1, xi0, v0, J, &points[q * (dim - 1)], &facePoints[(f * Np + q) * dim]);
830:       for (PetscInt c = 0; c < Nc; ++c) faceWeights[(f * Np + q) * Nc + c] = weights[q * Nc + c];
831:     }
832:   }
833:   PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)fq), efq));
834:   PetscCall(PetscQuadratureGetCellType(fq, &ct));
835:   PetscCall(PetscQuadratureSetCellType(*efq, ct));
836:   PetscCall(PetscQuadratureGetOrder(fq, &order));
837:   PetscCall(PetscQuadratureSetOrder(*efq, order));
838:   PetscCall(PetscQuadratureSetData(*efq, dim, Nc, Nf * Np, facePoints, faceWeights));
839:   PetscFunctionReturn(PETSC_SUCCESS);
840: }

842: /*@C
843:   PetscFEGetFaceTabulation - Returns the tabulation of the basis functions at the face quadrature points for each face of the reference cell

845:   Not Collective

847:   Input Parameters:
848: + fem - The `PetscFE` object
849: - k   - The highest derivative we need to tabulate, very often 1

851:   Output Parameter:
852: . Tf - The basis function values and derivatives at face quadrature points

854:   Level: intermediate

856:   Note:
857: .vb
858:   T->T[0] = Bf[((f*Nq + q)*pdim + i)*Nc + c] is the value at point f,q for basis function i and component c
859:   T->T[1] = Df[(((f*Nq + q)*pdim + i)*Nc + c)*dim + d] is the derivative value at point f,q for basis function i, component c, in direction d
860:   T->T[2] = Hf[((((f*Nq + q)*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point f,q for basis function i, component c, in directions d and e
861: .ve

863: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()`
864: @*/
865: PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscInt k, PetscTabulation *Tf)
866: {
867:   PetscFunctionBegin;
869:   PetscAssertPointer(Tf, 3);
870:   if (!fem->Tf) {
871:     PetscQuadrature fq;

873:     PetscCall(PetscFEGetFaceQuadrature(fem, &fq));
874:     if (fq) {
875:       PetscQuadrature  efq;
876:       const PetscReal *facePoints;
877:       PetscInt         Np, eNp;

879:       PetscCall(PetscFEExpandFaceQuadrature(fem, fq, &efq));
880:       PetscCall(PetscQuadratureGetData(fq, NULL, NULL, &Np, NULL, NULL));
881:       PetscCall(PetscQuadratureGetData(efq, NULL, NULL, &eNp, &facePoints, NULL));
882:       if (PetscDefined(USE_DEBUG)) {
883:         PetscDualSpace sp;
884:         DM             dm;
885:         PetscInt       cStart, Nf;

887:         PetscCall(PetscFEGetDualSpace(fem, &sp));
888:         PetscCall(PetscDualSpaceGetDM(sp, &dm));
889:         PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
890:         PetscCall(DMPlexGetConeSize(dm, cStart, &Nf));
891:         PetscCheck(Nf == eNp / Np, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of faces %" PetscInt_FMT " != %" PetscInt_FMT " number of quadrature replicas", Nf, eNp / Np);
892:       }
893:       PetscCall(PetscFECreateTabulation(fem, eNp / Np, Np, facePoints, k, &fem->Tf));
894:       PetscCall(PetscQuadratureDestroy(&efq));
895:     }
896:   }
897:   PetscCheck(!fem->Tf || k <= fem->Tf->K, PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_OUTOFRANGE, "Requested %" PetscInt_FMT " derivatives, but only tabulated %" PetscInt_FMT, k, fem->Tf->K);
898:   *Tf = fem->Tf;
899:   PetscFunctionReturn(PETSC_SUCCESS);
900: }

902: /*@C
903:   PetscFEGetFaceCentroidTabulation - Returns the tabulation of the basis functions at the face centroid points

905:   Not Collective

907:   Input Parameter:
908: . fem - The `PetscFE` object

910:   Output Parameter:
911: . Tc - The basis function values at face centroid points

913:   Level: intermediate

915:   Note:
916: .vb
917:   T->T[0] = Bf[(f*pdim + i)*Nc + c] is the value at point f for basis function i and component c
918: .ve

920: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFEGetFaceTabulation()`, `PetscFEGetCellTabulation()`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()`
921: @*/
922: PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscTabulation *Tc)
923: {
924:   PetscFunctionBegin;
926:   PetscAssertPointer(Tc, 2);
927:   if (!fem->Tc) {
928:     PetscDualSpace  sp;
929:     DM              dm;
930:     const PetscInt *cone;
931:     PetscReal      *centroids;
932:     PetscInt        dim, numFaces, f;

934:     PetscCall(PetscFEGetDualSpace(fem, &sp));
935:     PetscCall(PetscDualSpaceGetDM(sp, &dm));
936:     PetscCall(DMGetDimension(dm, &dim));
937:     PetscCall(DMPlexGetConeSize(dm, 0, &numFaces));
938:     PetscCall(DMPlexGetCone(dm, 0, &cone));
939:     PetscCall(PetscMalloc1(numFaces * dim, &centroids));
940:     for (f = 0; f < numFaces; ++f) PetscCall(DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, &centroids[f * dim], NULL));
941:     PetscCall(PetscFECreateTabulation(fem, 1, numFaces, centroids, 0, &fem->Tc));
942:     PetscCall(PetscFree(centroids));
943:   }
944:   *Tc = fem->Tc;
945:   PetscFunctionReturn(PETSC_SUCCESS);
946: }

948: /*@C
949:   PetscFECreateTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided.

951:   Not Collective

953:   Input Parameters:
954: + fem     - The `PetscFE` object
955: . nrepl   - The number of replicas
956: . npoints - The number of tabulation points in a replica
957: . points  - The tabulation point coordinates
958: - K       - The number of derivatives calculated

960:   Output Parameter:
961: . T - The basis function values and derivatives at tabulation points

963:   Level: intermediate

965:   Note:
966: .vb
967:   T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
968:   T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
969:   T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis
970:   T->function i, component c, in directions d and e
971: .ve

973: .seealso: `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscTabulationDestroy()`
974: @*/
975: PetscErrorCode PetscFECreateTabulation(PetscFE fem, PetscInt nrepl, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation *T)
976: {
977:   DM             dm;
978:   PetscDualSpace Q;
979:   PetscInt       Nb;   /* Dimension of FE space P */
980:   PetscInt       Nc;   /* Field components */
981:   PetscInt       cdim; /* Reference coordinate dimension */

983:   PetscFunctionBegin;
984:   if (!npoints || !fem->dualSpace || K < 0) {
985:     *T = NULL;
986:     PetscFunctionReturn(PETSC_SUCCESS);
987:   }
989:   PetscAssertPointer(points, 4);
990:   PetscAssertPointer(T, 6);
991:   PetscCall(PetscFEGetDualSpace(fem, &Q));
992:   PetscCall(PetscDualSpaceGetDM(Q, &dm));
993:   PetscCall(DMGetDimension(dm, &cdim));
994:   PetscCall(PetscDualSpaceGetDimension(Q, &Nb));
995:   PetscCall(PetscFEGetNumComponents(fem, &Nc));
996:   {
997:     PetscSpace sp;
998:     PetscInt   Nv;

1000:     PetscCall(PetscFEGetBasisSpace(fem, &sp));
1001:     PetscCall(PetscSpaceGetNumVariables(sp, &Nv));
1002:     PetscCheck(cdim == Nv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Dual space mesh dim %" PetscInt_FMT " != %" PetscInt_FMT " number of space variables", cdim, Nv);
1003:   }
1004:   PetscCall(PetscMalloc1(1, T));
1005:   (*T)->K    = !cdim ? 0 : K;
1006:   (*T)->Nr   = nrepl;
1007:   (*T)->Np   = npoints;
1008:   (*T)->Nb   = Nb;
1009:   (*T)->Nc   = Nc;
1010:   (*T)->cdim = cdim;
1011:   PetscCall(PetscMalloc1((*T)->K + 1, &(*T)->T));
1012:   for (PetscInt k = 0; k <= (*T)->K; ++k) PetscCall(PetscCalloc1(nrepl * npoints * Nb * Nc * PetscPowInt(cdim, k), &(*T)->T[k]));
1013:   PetscUseTypeMethod(fem, computetabulation, nrepl * npoints, points, K, *T);
1014:   PetscFunctionReturn(PETSC_SUCCESS);
1015: }

1017: /*@C
1018:   PetscFEComputeTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided.

1020:   Not Collective

1022:   Input Parameters:
1023: + fem     - The `PetscFE` object
1024: . npoints - The number of tabulation points
1025: . points  - The tabulation point coordinates
1026: . K       - The number of derivatives calculated
1027: - T       - An existing tabulation object with enough allocated space

1029:   Output Parameter:
1030: . T - The basis function values and derivatives at tabulation points

1032:   Level: intermediate

1034:   Note:
1035: .vb
1036:   T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
1037:   T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
1038:   T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e
1039: .ve

1041: .seealso: `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscTabulationDestroy()`
1042: @*/
1043: PetscErrorCode PetscFEComputeTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
1044: {
1045:   PetscFunctionBeginHot;
1046:   if (!npoints || !fem->dualSpace || K < 0) PetscFunctionReturn(PETSC_SUCCESS);
1048:   PetscAssertPointer(points, 3);
1049:   PetscAssertPointer(T, 5);
1050:   if (PetscDefined(USE_DEBUG)) {
1051:     DM             dm;
1052:     PetscDualSpace Q;
1053:     PetscInt       Nb;   /* Dimension of FE space P */
1054:     PetscInt       Nc;   /* Field components */
1055:     PetscInt       cdim; /* Reference coordinate dimension */

1057:     PetscCall(PetscFEGetDualSpace(fem, &Q));
1058:     PetscCall(PetscDualSpaceGetDM(Q, &dm));
1059:     PetscCall(DMGetDimension(dm, &cdim));
1060:     PetscCall(PetscDualSpaceGetDimension(Q, &Nb));
1061:     PetscCall(PetscFEGetNumComponents(fem, &Nc));
1062:     PetscCheck(T->K == (!cdim ? 0 : K), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation K %" PetscInt_FMT " must match requested K %" PetscInt_FMT, T->K, !cdim ? 0 : K);
1063:     PetscCheck(T->Nb == Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nb %" PetscInt_FMT " must match requested Nb %" PetscInt_FMT, T->Nb, Nb);
1064:     PetscCheck(T->Nc == Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nc %" PetscInt_FMT " must match requested Nc %" PetscInt_FMT, T->Nc, Nc);
1065:     PetscCheck(T->cdim == cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation cdim %" PetscInt_FMT " must match requested cdim %" PetscInt_FMT, T->cdim, cdim);
1066:   }
1067:   T->Nr = 1;
1068:   T->Np = npoints;
1069:   PetscUseTypeMethod(fem, computetabulation, npoints, points, K, T);
1070:   PetscFunctionReturn(PETSC_SUCCESS);
1071: }

1073: /*@
1074:   PetscTabulationDestroy - Frees memory from the associated tabulation.

1076:   Not Collective

1078:   Input Parameter:
1079: . T - The tabulation

1081:   Level: intermediate

1083: .seealso: `PetscTabulation`, `PetscFECreateTabulation()`, `PetscFEGetCellTabulation()`
1084: @*/
1085: PetscErrorCode PetscTabulationDestroy(PetscTabulation *T)
1086: {
1087:   PetscFunctionBegin;
1088:   PetscAssertPointer(T, 1);
1089:   if (!T || !*T) PetscFunctionReturn(PETSC_SUCCESS);
1090:   for (PetscInt k = 0; k <= (*T)->K; ++k) PetscCall(PetscFree((*T)->T[k]));
1091:   PetscCall(PetscFree((*T)->T));
1092:   PetscCall(PetscFree(*T));
1093:   *T = NULL;
1094:   PetscFunctionReturn(PETSC_SUCCESS);
1095: }

1097: static PetscErrorCode PetscFECreatePointTraceDefault_Internal(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
1098: {
1099:   PetscSpace      bsp, bsubsp;
1100:   PetscDualSpace  dsp, dsubsp;
1101:   PetscInt        dim, depth, numComp, i, j, coneSize, order;
1102:   DM              dm;
1103:   DMLabel         label;
1104:   PetscReal      *xi, *v, *J, detJ;
1105:   const char     *name;
1106:   PetscQuadrature origin, fullQuad, subQuad;

1108:   PetscFunctionBegin;
1109:   PetscCall(PetscFEGetBasisSpace(fe, &bsp));
1110:   PetscCall(PetscFEGetDualSpace(fe, &dsp));
1111:   PetscCall(PetscDualSpaceGetDM(dsp, &dm));
1112:   PetscCall(DMGetDimension(dm, &dim));
1113:   PetscCall(DMPlexGetDepthLabel(dm, &label));
1114:   PetscCall(DMLabelGetValue(label, refPoint, &depth));
1115:   PetscCall(PetscCalloc1(depth, &xi));
1116:   PetscCall(PetscMalloc1(dim, &v));
1117:   PetscCall(PetscMalloc1(dim * dim, &J));
1118:   for (i = 0; i < depth; i++) xi[i] = 0.;
1119:   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &origin));
1120:   PetscCall(PetscQuadratureSetData(origin, depth, 0, 1, xi, NULL));
1121:   PetscCall(DMPlexComputeCellGeometryFEM(dm, refPoint, origin, v, J, NULL, &detJ));
1122:   /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */
1123:   for (i = 1; i < dim; i++) {
1124:     for (j = 0; j < depth; j++) J[i * depth + j] = J[i * dim + j];
1125:   }
1126:   PetscCall(PetscQuadratureDestroy(&origin));
1127:   PetscCall(PetscDualSpaceGetPointSubspace(dsp, refPoint, &dsubsp));
1128:   PetscCall(PetscSpaceCreateSubspace(bsp, dsubsp, v, J, NULL, NULL, PETSC_OWN_POINTER, &bsubsp));
1129:   PetscCall(PetscSpaceSetUp(bsubsp));
1130:   PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe), trFE));
1131:   PetscCall(PetscFESetType(*trFE, PETSCFEBASIC));
1132:   PetscCall(PetscFEGetNumComponents(fe, &numComp));
1133:   PetscCall(PetscFESetNumComponents(*trFE, numComp));
1134:   PetscCall(PetscFESetBasisSpace(*trFE, bsubsp));
1135:   PetscCall(PetscFESetDualSpace(*trFE, dsubsp));
1136:   PetscCall(PetscObjectGetName((PetscObject)fe, &name));
1137:   if (name) PetscCall(PetscFESetName(*trFE, name));
1138:   PetscCall(PetscFEGetQuadrature(fe, &fullQuad));
1139:   PetscCall(PetscQuadratureGetOrder(fullQuad, &order));
1140:   PetscCall(DMPlexGetConeSize(dm, refPoint, &coneSize));
1141:   if (coneSize == 2 * depth) PetscCall(PetscDTGaussTensorQuadrature(depth, 1, (order + 2) / 2, -1., 1., &subQuad));
1142:   else PetscCall(PetscDTSimplexQuadrature(depth, order, PETSCDTSIMPLEXQUAD_DEFAULT, &subQuad));
1143:   PetscCall(PetscFESetQuadrature(*trFE, subQuad));
1144:   PetscCall(PetscFESetUp(*trFE));
1145:   PetscCall(PetscQuadratureDestroy(&subQuad));
1146:   PetscCall(PetscSpaceDestroy(&bsubsp));
1147:   PetscFunctionReturn(PETSC_SUCCESS);
1148: }

1150: PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
1151: {
1152:   PetscFunctionBegin;
1154:   PetscAssertPointer(trFE, 3);
1155:   if (fe->ops->createpointtrace) PetscUseTypeMethod(fe, createpointtrace, refPoint, trFE);
1156:   else PetscCall(PetscFECreatePointTraceDefault_Internal(fe, refPoint, trFE));
1157:   PetscFunctionReturn(PETSC_SUCCESS);
1158: }

1160: PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE)
1161: {
1162:   PetscInt       hStart, hEnd;
1163:   PetscDualSpace dsp;
1164:   DM             dm;

1166:   PetscFunctionBegin;
1168:   PetscAssertPointer(trFE, 3);
1169:   *trFE = NULL;
1170:   PetscCall(PetscFEGetDualSpace(fe, &dsp));
1171:   PetscCall(PetscDualSpaceGetDM(dsp, &dm));
1172:   PetscCall(DMPlexGetHeightStratum(dm, height, &hStart, &hEnd));
1173:   if (hEnd <= hStart) PetscFunctionReturn(PETSC_SUCCESS);
1174:   PetscCall(PetscFECreatePointTrace(fe, hStart, trFE));
1175:   PetscFunctionReturn(PETSC_SUCCESS);
1176: }

1178: /*@
1179:   PetscFEGetDimension - Get the dimension of the finite element space on a cell

1181:   Not Collective

1183:   Input Parameter:
1184: . fem - The `PetscFE`

1186:   Output Parameter:
1187: . dim - The dimension

1189:   Level: intermediate

1191: .seealso: `PetscFE`, `PetscFECreate()`, `PetscSpaceGetDimension()`, `PetscDualSpaceGetDimension()`
1192: @*/
1193: PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
1194: {
1195:   PetscFunctionBegin;
1197:   PetscAssertPointer(dim, 2);
1198:   PetscTryTypeMethod(fem, getdimension, dim);
1199:   PetscFunctionReturn(PETSC_SUCCESS);
1200: }

1202: /*@
1203:   PetscFEPushforward - Map the reference element function to real space

1205:   Input Parameters:
1206: + fe     - The `PetscFE`
1207: . fegeom - The cell geometry
1208: . Nv     - The number of function values
1209: - vals   - The function values

1211:   Output Parameter:
1212: . vals - The transformed function values

1214:   Level: advanced

1216:   Notes:
1217:   This just forwards the call onto `PetscDualSpacePushforward()`.

1219:   It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.

1221: .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscDualSpacePushforward()`
1222: @*/
1223: PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1224: {
1225:   PetscFunctionBeginHot;
1226:   PetscCall(PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals));
1227:   PetscFunctionReturn(PETSC_SUCCESS);
1228: }

1230: /*@
1231:   PetscFEPushforwardGradient - Map the reference element function gradient to real space

1233:   Input Parameters:
1234: + fe     - The `PetscFE`
1235: . fegeom - The cell geometry
1236: . Nv     - The number of function gradient values
1237: - vals   - The function gradient values

1239:   Output Parameter:
1240: . vals - The transformed function gradient values

1242:   Level: advanced

1244:   Notes:
1245:   This just forwards the call onto `PetscDualSpacePushforwardGradient()`.

1247:   It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.

1249: .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscFEPushforward()`, `PetscDualSpacePushforwardGradient()`, `PetscDualSpacePushforward()`
1250: @*/
1251: PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1252: {
1253:   PetscFunctionBeginHot;
1254:   PetscCall(PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals));
1255:   PetscFunctionReturn(PETSC_SUCCESS);
1256: }

1258: /*@
1259:   PetscFEPushforwardHessian - Map the reference element function Hessian to real space

1261:   Input Parameters:
1262: + fe     - The `PetscFE`
1263: . fegeom - The cell geometry
1264: . Nv     - The number of function Hessian values
1265: - vals   - The function Hessian values

1267:   Output Parameter:
1268: . vals - The transformed function Hessian values

1270:   Level: advanced

1272:   Notes:
1273:   This just forwards the call onto `PetscDualSpacePushforwardHessian()`.

1275:   It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.

1277:   Developer Notes:
1278:   It is unclear why all these one line convenience routines are desirable

1280: .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscFEPushforward()`, `PetscDualSpacePushforwardHessian()`, `PetscDualSpacePushforward()`
1281: @*/
1282: PetscErrorCode PetscFEPushforwardHessian(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1283: {
1284:   PetscFunctionBeginHot;
1285:   PetscCall(PetscDualSpacePushforwardHessian(fe->dualSpace, fegeom, Nv, fe->numComponents, vals));
1286:   PetscFunctionReturn(PETSC_SUCCESS);
1287: }

1289: /*
1290: Purpose: Compute element vector for chunk of elements

1292: Input:
1293:   Sizes:
1294:      Ne:  number of elements
1295:      Nf:  number of fields
1296:      PetscFE
1297:        dim: spatial dimension
1298:        Nb:  number of basis functions
1299:        Nc:  number of field components
1300:        PetscQuadrature
1301:          Nq:  number of quadrature points

1303:   Geometry:
1304:      PetscFEGeom[Ne] possibly *Nq
1305:        PetscReal v0s[dim]
1306:        PetscReal n[dim]
1307:        PetscReal jacobians[dim*dim]
1308:        PetscReal jacobianInverses[dim*dim]
1309:        PetscReal jacobianDeterminants
1310:   FEM:
1311:      PetscFE
1312:        PetscQuadrature
1313:          PetscReal   quadPoints[Nq*dim]
1314:          PetscReal   quadWeights[Nq]
1315:        PetscReal   basis[Nq*Nb*Nc]
1316:        PetscReal   basisDer[Nq*Nb*Nc*dim]
1317:      PetscScalar coefficients[Ne*Nb*Nc]
1318:      PetscScalar elemVec[Ne*Nb*Nc]

1320:   Problem:
1321:      PetscInt f: the active field
1322:      f0, f1

1324:   Work Space:
1325:      PetscFE
1326:        PetscScalar f0[Nq*dim];
1327:        PetscScalar f1[Nq*dim*dim];
1328:        PetscScalar u[Nc];
1329:        PetscScalar gradU[Nc*dim];
1330:        PetscReal   x[dim];
1331:        PetscScalar realSpaceDer[dim];

1333: Purpose: Compute element vector for N_cb batches of elements

1335: Input:
1336:   Sizes:
1337:      N_cb: Number of serial cell batches

1339:   Geometry:
1340:      PetscReal v0s[Ne*dim]
1341:      PetscReal jacobians[Ne*dim*dim]        possibly *Nq
1342:      PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq
1343:      PetscReal jacobianDeterminants[Ne]     possibly *Nq
1344:   FEM:
1345:      static PetscReal   quadPoints[Nq*dim]
1346:      static PetscReal   quadWeights[Nq]
1347:      static PetscReal   basis[Nq*Nb*Nc]
1348:      static PetscReal   basisDer[Nq*Nb*Nc*dim]
1349:      PetscScalar coefficients[Ne*Nb*Nc]
1350:      PetscScalar elemVec[Ne*Nb*Nc]

1352: ex62.c:
1353:   PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[],
1354:                                                const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[],
1355:                                                void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]),
1356:                                                void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[])

1358: ex52.c:
1359:   PetscErrorCode IntegrateLaplacianBatchCPU(PetscInt Ne, PetscInt Nb, const PetscScalar coefficients[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscInt Nq, const PetscReal quadPoints[], const PetscReal quadWeights[], const PetscReal basisTabulation[], const PetscReal basisDerTabulation[], PetscScalar elemVec[], AppCtx *user)
1360:   PetscErrorCode IntegrateElasticityBatchCPU(PetscInt Ne, PetscInt Nb, PetscInt Ncomp, const PetscScalar coefficients[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscInt Nq, const PetscReal quadPoints[], const PetscReal quadWeights[], const PetscReal basisTabulation[], const PetscReal basisDerTabulation[], PetscScalar elemVec[], AppCtx *user)

1362: ex52_integrateElement.cu
1363: __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec)

1365: PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[],
1366:                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1367:                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)

1369: ex52_integrateElementOpenCL.c:
1370: PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[],
1371:                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1372:                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)

1374: __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec)
1375: */

1377: /*@
1378:   PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration

1380:   Not Collective

1382:   Input Parameters:
1383: + prob            - The `PetscDS` specifying the discretizations and continuum functions
1384: . field           - The field being integrated
1385: . Ne              - The number of elements in the chunk
1386: . cgeom           - The cell geometry for each cell in the chunk
1387: . coefficients    - The array of FEM basis coefficients for the elements
1388: . probAux         - The `PetscDS` specifying the auxiliary discretizations
1389: - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements

1391:   Output Parameter:
1392: . integral - the integral for this field

1394:   Level: intermediate

1396:   Developer Notes:
1397:   The function name begins with `PetscFE` and yet the first argument is `PetscDS` and it has no `PetscFE` arguments.

1399: .seealso: `PetscFE`, `PetscDS`, `PetscFEIntegrateResidual()`, `PetscFEIntegrateBd()`
1400: @*/
1401: PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1402: {
1403:   PetscFE fe;

1405:   PetscFunctionBegin;
1407:   PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe));
1408:   if (fe->ops->integrate) PetscCall((*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral));
1409:   PetscFunctionReturn(PETSC_SUCCESS);
1410: }

1412: /*@C
1413:   PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration

1415:   Not Collective

1417:   Input Parameters:
1418: + prob            - The `PetscDS` specifying the discretizations and continuum functions
1419: . field           - The field being integrated
1420: . obj_func        - The function to be integrated
1421: . Ne              - The number of elements in the chunk
1422: . geom            - The face geometry for each face in the chunk
1423: . coefficients    - The array of FEM basis coefficients for the elements
1424: . probAux         - The `PetscDS` specifying the auxiliary discretizations
1425: - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements

1427:   Output Parameter:
1428: . integral - the integral for this field

1430:   Level: intermediate

1432:   Developer Notes:
1433:   The function name begins with `PetscFE` and yet the first argument is `PetscDS` and it has no `PetscFE` arguments.

1435: .seealso: `PetscFE`, `PetscDS`, `PetscFEIntegrateResidual()`, `PetscFEIntegrate()`
1436: @*/
1437: PetscErrorCode PetscFEIntegrateBd(PetscDS prob, PetscInt field, void (*obj_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[]), PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1438: {
1439:   PetscFE fe;

1441:   PetscFunctionBegin;
1443:   PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe));
1444:   if (fe->ops->integratebd) PetscCall((*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral));
1445:   PetscFunctionReturn(PETSC_SUCCESS);
1446: }

1448: /*@
1449:   PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration

1451:   Not Collective

1453:   Input Parameters:
1454: + ds              - The `PetscDS` specifying the discretizations and continuum functions
1455: . key             - The (label+value, field) being integrated
1456: . Ne              - The number of elements in the chunk
1457: . cgeom           - The cell geometry for each cell in the chunk
1458: . coefficients    - The array of FEM basis coefficients for the elements
1459: . coefficients_t  - The array of FEM basis time derivative coefficients for the elements
1460: . probAux         - The `PetscDS` specifying the auxiliary discretizations
1461: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1462: - t               - The time

1464:   Output Parameter:
1465: . elemVec - the element residual vectors from each element

1467:   Level: intermediate

1469:   Note:
1470: .vb
1471:   Loop over batch of elements (e):
1472:     Loop over quadrature points (q):
1473:       Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q
1474:       Call f_0 and f_1
1475:     Loop over element vector entries (f,fc --> i):
1476:       elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u)
1477: .ve

1479: .seealso: `PetscFEIntegrateBdResidual()`
1480: @*/
1481: PetscErrorCode PetscFEIntegrateResidual(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1482: {
1483:   PetscFE fe;

1485:   PetscFunctionBeginHot;
1487:   PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe));
1488:   if (fe->ops->integrateresidual) PetscCall((*fe->ops->integrateresidual)(ds, key, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec));
1489:   PetscFunctionReturn(PETSC_SUCCESS);
1490: }

1492: /*@
1493:   PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary

1495:   Not Collective

1497:   Input Parameters:
1498: + ds              - The `PetscDS` specifying the discretizations and continuum functions
1499: . wf              - The PetscWeakForm object holding the pointwise functions
1500: . key             - The (label+value, field) being integrated
1501: . Ne              - The number of elements in the chunk
1502: . fgeom           - The face geometry for each cell in the chunk
1503: . coefficients    - The array of FEM basis coefficients for the elements
1504: . coefficients_t  - The array of FEM basis time derivative coefficients for the elements
1505: . probAux         - The `PetscDS` specifying the auxiliary discretizations
1506: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1507: - t               - The time

1509:   Output Parameter:
1510: . elemVec - the element residual vectors from each element

1512:   Level: intermediate

1514: .seealso: `PetscFEIntegrateResidual()`
1515: @*/
1516: PetscErrorCode PetscFEIntegrateBdResidual(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1517: {
1518:   PetscFE fe;

1520:   PetscFunctionBegin;
1522:   PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe));
1523:   if (fe->ops->integratebdresidual) PetscCall((*fe->ops->integratebdresidual)(ds, wf, key, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec));
1524:   PetscFunctionReturn(PETSC_SUCCESS);
1525: }

1527: /*@
1528:   PetscFEIntegrateHybridResidual - Produce the element residual vector for a chunk of hybrid element faces by quadrature integration

1530:   Not Collective

1532:   Input Parameters:
1533: + ds              - The `PetscDS` specifying the discretizations and continuum functions
1534: . dsIn            - The `PetscDS` specifying the discretizations and continuum functions for input
1535: . key             - The (label+value, field) being integrated
1536: . s               - The side of the cell being integrated, 0 for negative and 1 for positive
1537: . Ne              - The number of elements in the chunk
1538: . fgeom           - The face geometry for each cell in the chunk
1539: . cgeom           - The cell geometry for each neighbor cell in the chunk
1540: . coefficients    - The array of FEM basis coefficients for the elements
1541: . coefficients_t  - The array of FEM basis time derivative coefficients for the elements
1542: . probAux         - The `PetscDS` specifying the auxiliary discretizations
1543: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1544: - t               - The time

1546:   Output Parameter:
1547: . elemVec - the element residual vectors from each element

1549:   Level: developer

1551: .seealso: `PetscFEIntegrateResidual()`
1552: @*/
1553: PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS ds, PetscDS dsIn, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1554: {
1555:   PetscFE fe;

1557:   PetscFunctionBegin;
1560:   PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe));
1561:   if (fe->ops->integratehybridresidual) PetscCall((*fe->ops->integratehybridresidual)(ds, dsIn, key, s, Ne, fgeom, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec));
1562:   PetscFunctionReturn(PETSC_SUCCESS);
1563: }

1565: /*@
1566:   PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration

1568:   Not Collective

1570:   Input Parameters:
1571: + rds             - The `PetscDS` specifying the row discretizations and continuum functions
1572: . cds             - The `PetscDS` specifying the column discretizations
1573: . jtype           - The type of matrix pointwise functions that should be used
1574: . key             - The (label+value, fieldI*Nf + fieldJ) being integrated
1575: . Ne              - The number of elements in the chunk
1576: . cgeom           - The cell geometry for each cell in the chunk
1577: . coefficients    - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1578: . coefficients_t  - The array of FEM basis time derivative coefficients for the elements
1579: . dsAux           - The `PetscDS` specifying the auxiliary discretizations
1580: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1581: . t               - The time
1582: - u_tshift        - A multiplier for the dF/du_t term (as opposed to the dF/du term)

1584:   Output Parameter:
1585: . elemMat - the element matrices for the Jacobian from each element

1587:   Level: intermediate

1589:   Note:
1590: .vb
1591:   Loop over batch of elements (e):
1592:     Loop over element matrix entries (f,fc,g,gc --> i,j):
1593:       Loop over quadrature points (q):
1594:         Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1595:           elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1596:                        + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1597:                        + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1598:                        + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1599: .ve

1601: .seealso: `PetscFEIntegrateResidual()`
1602: @*/
1603: PetscErrorCode PetscFEIntegrateJacobian(PetscDS rds, PetscDS cds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1604: {
1605:   PetscFE  fe;
1606:   PetscInt Nf;

1608:   PetscFunctionBegin;
1611:   PetscCall(PetscDSGetNumFields(rds, &Nf));
1612:   PetscCall(PetscDSGetDiscretization(rds, key.field / Nf, (PetscObject *)&fe));
1613:   if (fe->ops->integratejacobian) PetscCall((*fe->ops->integratejacobian)(rds, cds, jtype, key, Ne, cgeom, coefficients, coefficients_t, dsAux, coefficientsAux, t, u_tshift, elemMat));
1614:   PetscFunctionReturn(PETSC_SUCCESS);
1615: }

1617: /*@
1618:   PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration

1620:   Not Collective

1622:   Input Parameters:
1623: + ds              - The `PetscDS` specifying the discretizations and continuum functions
1624: . wf              - The PetscWeakForm holding the pointwise functions
1625: . jtype           - The type of matrix pointwise functions that should be used
1626: . key             - The (label+value, fieldI*Nf + fieldJ) being integrated
1627: . Ne              - The number of elements in the chunk
1628: . fgeom           - The face geometry for each cell in the chunk
1629: . coefficients    - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1630: . coefficients_t  - The array of FEM basis time derivative coefficients for the elements
1631: . probAux         - The `PetscDS` specifying the auxiliary discretizations
1632: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1633: . t               - The time
1634: - u_tshift        - A multiplier for the dF/du_t term (as opposed to the dF/du term)

1636:   Output Parameter:
1637: . elemMat - the element matrices for the Jacobian from each element

1639:   Level: intermediate

1641:   Note:
1642: .vb
1643:   Loop over batch of elements (e):
1644:     Loop over element matrix entries (f,fc,g,gc --> i,j):
1645:       Loop over quadrature points (q):
1646:         Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1647:           elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1648:                        + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1649:                        + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1650:                        + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1651: .ve

1653: .seealso: `PetscFEIntegrateJacobian()`, `PetscFEIntegrateResidual()`
1654: @*/
1655: PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS ds, PetscWeakForm wf, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1656: {
1657:   PetscFE  fe;
1658:   PetscInt Nf;

1660:   PetscFunctionBegin;
1662:   PetscCall(PetscDSGetNumFields(ds, &Nf));
1663:   PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *)&fe));
1664:   if (fe->ops->integratebdjacobian) PetscCall((*fe->ops->integratebdjacobian)(ds, wf, jtype, key, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat));
1665:   PetscFunctionReturn(PETSC_SUCCESS);
1666: }

1668: /*@
1669:   PetscFEIntegrateHybridJacobian - Produce the boundary element Jacobian for a chunk of hybrid elements by quadrature integration

1671:   Not Collective

1673:   Input Parameters:
1674: + ds              - The `PetscDS` specifying the discretizations and continuum functions for the output
1675: . dsIn            - The `PetscDS` specifying the discretizations and continuum functions for the input
1676: . jtype           - The type of matrix pointwise functions that should be used
1677: . key             - The (label+value, fieldI*Nf + fieldJ) being integrated
1678: . s               - The side of the cell being integrated, 0 for negative and 1 for positive
1679: . Ne              - The number of elements in the chunk
1680: . fgeom           - The face geometry for each cell in the chunk
1681: . cgeom           - The cell geometry for each neighbor cell in the chunk
1682: . coefficients    - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1683: . coefficients_t  - The array of FEM basis time derivative coefficients for the elements
1684: . probAux         - The `PetscDS` specifying the auxiliary discretizations
1685: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1686: . t               - The time
1687: - u_tshift        - A multiplier for the dF/du_t term (as opposed to the dF/du term)

1689:   Output Parameter:
1690: . elemMat - the element matrices for the Jacobian from each element

1692:   Level: developer

1694:   Note:
1695: .vb
1696:   Loop over batch of elements (e):
1697:     Loop over element matrix entries (f,fc,g,gc --> i,j):
1698:       Loop over quadrature points (q):
1699:         Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1700:           elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1701:                        + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1702:                        + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1703:                        + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1704: .ve

1706: .seealso: `PetscFEIntegrateJacobian()`, `PetscFEIntegrateResidual()`
1707: @*/
1708: PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS ds, PetscDS dsIn, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1709: {
1710:   PetscFE  fe;
1711:   PetscInt Nf;

1713:   PetscFunctionBegin;
1715:   PetscCall(PetscDSGetNumFields(ds, &Nf));
1716:   PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *)&fe));
1717:   if (fe->ops->integratehybridjacobian) PetscCall((*fe->ops->integratehybridjacobian)(ds, dsIn, jtype, key, s, Ne, fgeom, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat));
1718:   PetscFunctionReturn(PETSC_SUCCESS);
1719: }

1721: /*@
1722:   PetscFEGetHeightSubspace - Get the subspace of this space for a mesh point of a given height

1724:   Input Parameters:
1725: + fe     - The finite element space
1726: - height - The height of the `DMPLEX` point

1728:   Output Parameter:
1729: . subfe - The subspace of this `PetscFE` space

1731:   Level: advanced

1733:   Note:
1734:   For example, if we want the subspace of this space for a face, we would choose height = 1.

1736: .seealso: `PetscFECreateDefault()`
1737: @*/
1738: PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe)
1739: {
1740:   PetscSpace      P, subP;
1741:   PetscDualSpace  Q, subQ;
1742:   PetscQuadrature subq;
1743:   PetscInt        dim, Nc;

1745:   PetscFunctionBegin;
1747:   PetscAssertPointer(subfe, 3);
1748:   if (height == 0) {
1749:     *subfe = fe;
1750:     PetscFunctionReturn(PETSC_SUCCESS);
1751:   }
1752:   PetscCall(PetscFEGetBasisSpace(fe, &P));
1753:   PetscCall(PetscFEGetDualSpace(fe, &Q));
1754:   PetscCall(PetscFEGetNumComponents(fe, &Nc));
1755:   PetscCall(PetscFEGetFaceQuadrature(fe, &subq));
1756:   PetscCall(PetscDualSpaceGetDimension(Q, &dim));
1757:   PetscCheck(height <= dim && height >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %" PetscInt_FMT " for dimension %" PetscInt_FMT " space", height, dim);
1758:   if (!fe->subspaces) PetscCall(PetscCalloc1(dim, &fe->subspaces));
1759:   if (height <= dim) {
1760:     if (!fe->subspaces[height - 1]) {
1761:       PetscFE     sub = NULL;
1762:       const char *name;

1764:       PetscCall(PetscSpaceGetHeightSubspace(P, height, &subP));
1765:       PetscCall(PetscDualSpaceGetHeightSubspace(Q, height, &subQ));
1766:       if (subQ) {
1767:         PetscCall(PetscObjectReference((PetscObject)subP));
1768:         PetscCall(PetscObjectReference((PetscObject)subQ));
1769:         PetscCall(PetscObjectReference((PetscObject)subq));
1770:         PetscCall(PetscFECreateFromSpaces(subP, subQ, subq, NULL, &sub));
1771:       }
1772:       if (sub) {
1773:         PetscCall(PetscObjectGetName((PetscObject)fe, &name));
1774:         if (name) PetscCall(PetscFESetName(sub, name));
1775:       }
1776:       fe->subspaces[height - 1] = sub;
1777:     }
1778:     *subfe = fe->subspaces[height - 1];
1779:   } else {
1780:     *subfe = NULL;
1781:   }
1782:   PetscFunctionReturn(PETSC_SUCCESS);
1783: }

1785: /*@
1786:   PetscFERefine - Create a "refined" `PetscFE` object that refines the reference cell into
1787:   smaller copies.

1789:   Collective

1791:   Input Parameter:
1792: . fe - The initial `PetscFE`

1794:   Output Parameter:
1795: . feRef - The refined `PetscFE`

1797:   Level: advanced

1799:   Notes:
1800:   This is typically used to generate a preconditioner for a higher order method from a lower order method on a
1801:   refined mesh having the same number of dofs (but more sparsity). It is also used to create an
1802:   interpolation between regularly refined meshes.

1804: .seealso: `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`
1805: @*/
1806: PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)
1807: {
1808:   PetscSpace       P, Pref;
1809:   PetscDualSpace   Q, Qref;
1810:   DM               K, Kref;
1811:   PetscQuadrature  q, qref;
1812:   const PetscReal *v0, *jac;
1813:   PetscInt         numComp, numSubelements;
1814:   PetscInt         cStart, cEnd, c;
1815:   PetscDualSpace  *cellSpaces;

1817:   PetscFunctionBegin;
1818:   PetscCall(PetscFEGetBasisSpace(fe, &P));
1819:   PetscCall(PetscFEGetDualSpace(fe, &Q));
1820:   PetscCall(PetscFEGetQuadrature(fe, &q));
1821:   PetscCall(PetscDualSpaceGetDM(Q, &K));
1822:   /* Create space */
1823:   PetscCall(PetscObjectReference((PetscObject)P));
1824:   Pref = P;
1825:   /* Create dual space */
1826:   PetscCall(PetscDualSpaceDuplicate(Q, &Qref));
1827:   PetscCall(PetscDualSpaceSetType(Qref, PETSCDUALSPACEREFINED));
1828:   PetscCall(DMRefine(K, PetscObjectComm((PetscObject)fe), &Kref));
1829:   PetscCall(DMGetCoordinatesLocalSetUp(Kref));
1830:   PetscCall(PetscDualSpaceSetDM(Qref, Kref));
1831:   PetscCall(DMPlexGetHeightStratum(Kref, 0, &cStart, &cEnd));
1832:   PetscCall(PetscMalloc1(cEnd - cStart, &cellSpaces));
1833:   /* TODO: fix for non-uniform refinement */
1834:   for (c = 0; c < cEnd - cStart; c++) cellSpaces[c] = Q;
1835:   PetscCall(PetscDualSpaceRefinedSetCellSpaces(Qref, cellSpaces));
1836:   PetscCall(PetscFree(cellSpaces));
1837:   PetscCall(DMDestroy(&Kref));
1838:   PetscCall(PetscDualSpaceSetUp(Qref));
1839:   /* Create element */
1840:   PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe), feRef));
1841:   PetscCall(PetscFESetType(*feRef, PETSCFECOMPOSITE));
1842:   PetscCall(PetscFESetBasisSpace(*feRef, Pref));
1843:   PetscCall(PetscFESetDualSpace(*feRef, Qref));
1844:   PetscCall(PetscFEGetNumComponents(fe, &numComp));
1845:   PetscCall(PetscFESetNumComponents(*feRef, numComp));
1846:   PetscCall(PetscFESetUp(*feRef));
1847:   PetscCall(PetscSpaceDestroy(&Pref));
1848:   PetscCall(PetscDualSpaceDestroy(&Qref));
1849:   /* Create quadrature */
1850:   PetscCall(PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL));
1851:   PetscCall(PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref));
1852:   PetscCall(PetscFESetQuadrature(*feRef, qref));
1853:   PetscCall(PetscQuadratureDestroy(&qref));
1854:   PetscFunctionReturn(PETSC_SUCCESS);
1855: }

1857: static PetscErrorCode PetscFESetDefaultName_Private(PetscFE fe)
1858: {
1859:   PetscSpace     P;
1860:   PetscDualSpace Q;
1861:   DM             K;
1862:   DMPolytopeType ct;
1863:   PetscInt       degree;
1864:   char           name[64];

1866:   PetscFunctionBegin;
1867:   PetscCall(PetscFEGetBasisSpace(fe, &P));
1868:   PetscCall(PetscSpaceGetDegree(P, &degree, NULL));
1869:   PetscCall(PetscFEGetDualSpace(fe, &Q));
1870:   PetscCall(PetscDualSpaceGetDM(Q, &K));
1871:   PetscCall(DMPlexGetCellType(K, 0, &ct));
1872:   switch (ct) {
1873:   case DM_POLYTOPE_SEGMENT:
1874:   case DM_POLYTOPE_POINT_PRISM_TENSOR:
1875:   case DM_POLYTOPE_QUADRILATERAL:
1876:   case DM_POLYTOPE_SEG_PRISM_TENSOR:
1877:   case DM_POLYTOPE_HEXAHEDRON:
1878:   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
1879:     PetscCall(PetscSNPrintf(name, sizeof(name), "Q%" PetscInt_FMT, degree));
1880:     break;
1881:   case DM_POLYTOPE_TRIANGLE:
1882:   case DM_POLYTOPE_TETRAHEDRON:
1883:     PetscCall(PetscSNPrintf(name, sizeof(name), "P%" PetscInt_FMT, degree));
1884:     break;
1885:   case DM_POLYTOPE_TRI_PRISM:
1886:   case DM_POLYTOPE_TRI_PRISM_TENSOR:
1887:     PetscCall(PetscSNPrintf(name, sizeof(name), "P%" PetscInt_FMT "xQ%" PetscInt_FMT, degree, degree));
1888:     break;
1889:   default:
1890:     PetscCall(PetscSNPrintf(name, sizeof(name), "FE"));
1891:   }
1892:   PetscCall(PetscFESetName(fe, name));
1893:   PetscFunctionReturn(PETSC_SUCCESS);
1894: }

1896: /*@
1897:   PetscFECreateFromSpaces - Create a `PetscFE` from the basis and dual spaces

1899:   Collective

1901:   Input Parameters:
1902: + P  - The basis space
1903: . Q  - The dual space
1904: . q  - The cell quadrature
1905: - fq - The face quadrature

1907:   Output Parameter:
1908: . fem - The `PetscFE` object

1910:   Level: beginner

1912:   Note:
1913:   The `PetscFE` takes ownership of these spaces by calling destroy on each. They should not be used after this call, and for borrowed references from `PetscFEGetSpace()` and the like, the caller must use `PetscObjectReference` before this call.

1915: .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`,
1916:           `PetscFECreateLagrangeByCell()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
1917: @*/
1918: PetscErrorCode PetscFECreateFromSpaces(PetscSpace P, PetscDualSpace Q, PetscQuadrature q, PetscQuadrature fq, PetscFE *fem)
1919: {
1920:   PetscInt    Nc;
1921:   PetscInt    p_Ns = -1, p_Nc = -1, q_Ns = -1, q_Nc = -1;
1922:   PetscBool   p_is_uniform_sum = PETSC_FALSE, p_interleave_basis = PETSC_FALSE, p_interleave_components = PETSC_FALSE;
1923:   PetscBool   q_is_uniform_sum = PETSC_FALSE, q_interleave_basis = PETSC_FALSE, q_interleave_components = PETSC_FALSE;
1924:   const char *prefix;

1926:   PetscFunctionBegin;
1927:   PetscCall(PetscObjectTypeCompare((PetscObject)P, PETSCSPACESUM, &p_is_uniform_sum));
1928:   if (p_is_uniform_sum) {
1929:     PetscSpace subsp_0 = NULL;
1930:     PetscCall(PetscSpaceSumGetNumSubspaces(P, &p_Ns));
1931:     PetscCall(PetscSpaceGetNumComponents(P, &p_Nc));
1932:     PetscCall(PetscSpaceSumGetConcatenate(P, &p_is_uniform_sum));
1933:     PetscCall(PetscSpaceSumGetInterleave(P, &p_interleave_basis, &p_interleave_components));
1934:     for (PetscInt s = 0; s < p_Ns; s++) {
1935:       PetscSpace subsp;

1937:       PetscCall(PetscSpaceSumGetSubspace(P, s, &subsp));
1938:       if (!s) {
1939:         subsp_0 = subsp;
1940:       } else if (subsp != subsp_0) {
1941:         p_is_uniform_sum = PETSC_FALSE;
1942:       }
1943:     }
1944:   }
1945:   PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACESUM, &q_is_uniform_sum));
1946:   if (q_is_uniform_sum) {
1947:     PetscDualSpace subsp_0 = NULL;
1948:     PetscCall(PetscDualSpaceSumGetNumSubspaces(Q, &q_Ns));
1949:     PetscCall(PetscDualSpaceGetNumComponents(Q, &q_Nc));
1950:     PetscCall(PetscDualSpaceSumGetConcatenate(Q, &q_is_uniform_sum));
1951:     PetscCall(PetscDualSpaceSumGetInterleave(Q, &q_interleave_basis, &q_interleave_components));
1952:     for (PetscInt s = 0; s < q_Ns; s++) {
1953:       PetscDualSpace subsp;

1955:       PetscCall(PetscDualSpaceSumGetSubspace(Q, s, &subsp));
1956:       if (!s) {
1957:         subsp_0 = subsp;
1958:       } else if (subsp != subsp_0) {
1959:         q_is_uniform_sum = PETSC_FALSE;
1960:       }
1961:     }
1962:   }
1963:   if (p_is_uniform_sum && q_is_uniform_sum && (p_interleave_basis == q_interleave_basis) && (p_interleave_components == q_interleave_components) && (p_Ns == q_Ns) && (p_Nc == q_Nc)) {
1964:     PetscSpace     scalar_space;
1965:     PetscDualSpace scalar_dspace;
1966:     PetscFE        scalar_fe;

1968:     PetscCall(PetscSpaceSumGetSubspace(P, 0, &scalar_space));
1969:     PetscCall(PetscDualSpaceSumGetSubspace(Q, 0, &scalar_dspace));
1970:     PetscCall(PetscObjectReference((PetscObject)scalar_space));
1971:     PetscCall(PetscObjectReference((PetscObject)scalar_dspace));
1972:     PetscCall(PetscObjectReference((PetscObject)q));
1973:     PetscCall(PetscObjectReference((PetscObject)fq));
1974:     PetscCall(PetscFECreateFromSpaces(scalar_space, scalar_dspace, q, fq, &scalar_fe));
1975:     PetscCall(PetscFECreateVector(scalar_fe, p_Ns, p_interleave_basis, p_interleave_components, fem));
1976:     PetscCall(PetscFEDestroy(&scalar_fe));
1977:   } else {
1978:     PetscCall(PetscFECreate(PetscObjectComm((PetscObject)P), fem));
1979:     PetscCall(PetscFESetType(*fem, PETSCFEBASIC));
1980:   }
1981:   PetscCall(PetscSpaceGetNumComponents(P, &Nc));
1982:   PetscCall(PetscFESetNumComponents(*fem, Nc));
1983:   PetscCall(PetscFESetBasisSpace(*fem, P));
1984:   PetscCall(PetscFESetDualSpace(*fem, Q));
1985:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)P, &prefix));
1986:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*fem, prefix));
1987:   PetscCall(PetscFESetUp(*fem));
1988:   PetscCall(PetscSpaceDestroy(&P));
1989:   PetscCall(PetscDualSpaceDestroy(&Q));
1990:   PetscCall(PetscFESetQuadrature(*fem, q));
1991:   PetscCall(PetscFESetFaceQuadrature(*fem, fq));
1992:   PetscCall(PetscQuadratureDestroy(&q));
1993:   PetscCall(PetscQuadratureDestroy(&fq));
1994:   PetscCall(PetscFESetDefaultName_Private(*fem));
1995:   PetscFunctionReturn(PETSC_SUCCESS);
1996: }

1998: static PetscErrorCode PetscFECreate_Internal(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, const char prefix[], PetscInt degree, PetscInt qorder, PetscBool setFromOptions, PetscFE *fem)
1999: {
2000:   DM              K;
2001:   PetscSpace      P;
2002:   PetscDualSpace  Q;
2003:   PetscQuadrature q, fq;
2004:   PetscBool       tensor;

2006:   PetscFunctionBegin;
2007:   if (prefix) PetscAssertPointer(prefix, 5);
2008:   PetscAssertPointer(fem, 9);
2009:   switch (ct) {
2010:   case DM_POLYTOPE_SEGMENT:
2011:   case DM_POLYTOPE_POINT_PRISM_TENSOR:
2012:   case DM_POLYTOPE_QUADRILATERAL:
2013:   case DM_POLYTOPE_SEG_PRISM_TENSOR:
2014:   case DM_POLYTOPE_HEXAHEDRON:
2015:   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
2016:     tensor = PETSC_TRUE;
2017:     break;
2018:   default:
2019:     tensor = PETSC_FALSE;
2020:   }
2021:   /* Create space */
2022:   PetscCall(PetscSpaceCreate(comm, &P));
2023:   PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL));
2024:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)P, prefix));
2025:   PetscCall(PetscSpacePolynomialSetTensor(P, tensor));
2026:   PetscCall(PetscSpaceSetNumComponents(P, Nc));
2027:   PetscCall(PetscSpaceSetNumVariables(P, dim));
2028:   if (degree >= 0) {
2029:     PetscCall(PetscSpaceSetDegree(P, degree, PETSC_DETERMINE));
2030:     if (ct == DM_POLYTOPE_TRI_PRISM || ct == DM_POLYTOPE_TRI_PRISM_TENSOR) {
2031:       PetscSpace Pend, Pside;

2033:       PetscCall(PetscSpaceSetNumComponents(P, 1));
2034:       PetscCall(PetscSpaceCreate(comm, &Pend));
2035:       PetscCall(PetscSpaceSetType(Pend, PETSCSPACEPOLYNOMIAL));
2036:       PetscCall(PetscSpacePolynomialSetTensor(Pend, PETSC_FALSE));
2037:       PetscCall(PetscSpaceSetNumComponents(Pend, 1));
2038:       PetscCall(PetscSpaceSetNumVariables(Pend, dim - 1));
2039:       PetscCall(PetscSpaceSetDegree(Pend, degree, PETSC_DETERMINE));
2040:       PetscCall(PetscSpaceCreate(comm, &Pside));
2041:       PetscCall(PetscSpaceSetType(Pside, PETSCSPACEPOLYNOMIAL));
2042:       PetscCall(PetscSpacePolynomialSetTensor(Pside, PETSC_FALSE));
2043:       PetscCall(PetscSpaceSetNumComponents(Pside, 1));
2044:       PetscCall(PetscSpaceSetNumVariables(Pside, 1));
2045:       PetscCall(PetscSpaceSetDegree(Pside, degree, PETSC_DETERMINE));
2046:       PetscCall(PetscSpaceSetType(P, PETSCSPACETENSOR));
2047:       PetscCall(PetscSpaceTensorSetNumSubspaces(P, 2));
2048:       PetscCall(PetscSpaceTensorSetSubspace(P, 0, Pend));
2049:       PetscCall(PetscSpaceTensorSetSubspace(P, 1, Pside));
2050:       PetscCall(PetscSpaceDestroy(&Pend));
2051:       PetscCall(PetscSpaceDestroy(&Pside));

2053:       if (Nc > 1) {
2054:         PetscSpace scalar_P = P;

2056:         PetscCall(PetscSpaceCreate(comm, &P));
2057:         PetscCall(PetscSpaceSetNumVariables(P, dim));
2058:         PetscCall(PetscSpaceSetNumComponents(P, Nc));
2059:         PetscCall(PetscSpaceSetType(P, PETSCSPACESUM));
2060:         PetscCall(PetscSpaceSumSetNumSubspaces(P, Nc));
2061:         PetscCall(PetscSpaceSumSetConcatenate(P, PETSC_TRUE));
2062:         PetscCall(PetscSpaceSumSetInterleave(P, PETSC_TRUE, PETSC_FALSE));
2063:         for (PetscInt i = 0; i < Nc; i++) PetscCall(PetscSpaceSumSetSubspace(P, i, scalar_P));
2064:         PetscCall(PetscSpaceDestroy(&scalar_P));
2065:       }
2066:     }
2067:   }
2068:   if (setFromOptions) PetscCall(PetscSpaceSetFromOptions(P));
2069:   PetscCall(PetscSpaceSetUp(P));
2070:   PetscCall(PetscSpaceGetDegree(P, &degree, NULL));
2071:   PetscCall(PetscSpacePolynomialGetTensor(P, &tensor));
2072:   PetscCall(PetscSpaceGetNumComponents(P, &Nc));
2073:   /* Create dual space */
2074:   PetscCall(PetscDualSpaceCreate(comm, &Q));
2075:   PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE));
2076:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)Q, prefix));
2077:   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K));
2078:   PetscCall(PetscDualSpaceSetDM(Q, K));
2079:   PetscCall(DMDestroy(&K));
2080:   PetscCall(PetscDualSpaceSetNumComponents(Q, Nc));
2081:   PetscCall(PetscDualSpaceSetOrder(Q, degree));
2082:   PetscCall(PetscDualSpaceLagrangeSetTensor(Q, (tensor || (ct == DM_POLYTOPE_TRI_PRISM)) ? PETSC_TRUE : PETSC_FALSE));
2083:   if (setFromOptions) PetscCall(PetscDualSpaceSetFromOptions(Q));
2084:   PetscCall(PetscDualSpaceSetUp(Q));
2085:   /* Create quadrature */
2086:   PetscDTSimplexQuadratureType qtype = PETSCDTSIMPLEXQUAD_DEFAULT;

2088:   qorder = qorder >= 0 ? qorder : degree;
2089:   if (setFromOptions) {
2090:     PetscObjectOptionsBegin((PetscObject)P);
2091:     PetscCall(PetscOptionsBoundedInt("-petscfe_default_quadrature_order", "Quadrature order is one less than quadrature points per edge", "PetscFECreateDefault", qorder, &qorder, NULL, 0));
2092:     PetscCall(PetscOptionsEnum("-petscfe_default_quadrature_type", "Simplex quadrature type", "PetscDTSimplexQuadratureType", PetscDTSimplexQuadratureTypes, (PetscEnum)qtype, (PetscEnum *)&qtype, NULL));
2093:     PetscOptionsEnd();
2094:   }
2095:   PetscCall(PetscDTCreateQuadratureByCell(ct, qorder, qtype, &q, &fq));
2096:   /* Create finite element */
2097:   PetscCall(PetscFECreateFromSpaces(P, Q, q, fq, fem));
2098:   if (setFromOptions) PetscCall(PetscFESetFromOptions(*fem));
2099:   PetscFunctionReturn(PETSC_SUCCESS);
2100: }

2102: /*@
2103:   PetscFECreateDefault - Create a `PetscFE` for basic FEM computation

2105:   Collective

2107:   Input Parameters:
2108: + comm      - The MPI comm
2109: . dim       - The spatial dimension
2110: . Nc        - The number of components
2111: . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
2112: . prefix    - The options prefix, or `NULL`
2113: - qorder    - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree

2115:   Output Parameter:
2116: . fem - The `PetscFE` object

2118:   Level: beginner

2120:   Note:
2121:   Each subobject is SetFromOption() during creation, so that the object may be customized from the command line, using the prefix specified above. See the links below for the particular options available.

2123: .seealso: `PetscFECreateLagrange()`, `PetscFECreateByCell()`, `PetscSpaceSetFromOptions()`, `PetscDualSpaceSetFromOptions()`, `PetscFESetFromOptions()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2124: @*/
2125: PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem)
2126: {
2127:   PetscFunctionBegin;
2128:   PetscCall(PetscFECreate_Internal(comm, dim, Nc, DMPolytopeTypeSimpleShape(dim, isSimplex), prefix, PETSC_DECIDE, qorder, PETSC_TRUE, fem));
2129:   PetscFunctionReturn(PETSC_SUCCESS);
2130: }

2132: /*@
2133:   PetscFECreateByCell - Create a `PetscFE` for basic FEM computation

2135:   Collective

2137:   Input Parameters:
2138: + comm   - The MPI comm
2139: . dim    - The spatial dimension
2140: . Nc     - The number of components
2141: . ct     - The celltype of the reference cell
2142: . prefix - The options prefix, or `NULL`
2143: - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree

2145:   Output Parameter:
2146: . fem - The `PetscFE` object

2148:   Level: beginner

2150:   Note:
2151:   Each subobject is SetFromOption() during creation, so that the object may be customized from the command line, using the prefix specified above. See the links below for the particular options available.

2153: .seealso: `PetscFECreateDefault()`, `PetscFECreateLagrange()`, `PetscSpaceSetFromOptions()`, `PetscDualSpaceSetFromOptions()`, `PetscFESetFromOptions()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2154: @*/
2155: PetscErrorCode PetscFECreateByCell(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, const char prefix[], PetscInt qorder, PetscFE *fem)
2156: {
2157:   PetscFunctionBegin;
2158:   PetscCall(PetscFECreate_Internal(comm, dim, Nc, ct, prefix, PETSC_DECIDE, qorder, PETSC_TRUE, fem));
2159:   PetscFunctionReturn(PETSC_SUCCESS);
2160: }

2162: /*@
2163:   PetscFECreateLagrange - Create a `PetscFE` for the basic Lagrange space of degree k

2165:   Collective

2167:   Input Parameters:
2168: + comm      - The MPI comm
2169: . dim       - The spatial dimension
2170: . Nc        - The number of components
2171: . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
2172: . k         - The degree k of the space
2173: - qorder    - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree

2175:   Output Parameter:
2176: . fem - The `PetscFE` object

2178:   Level: beginner

2180:   Note:
2181:   For simplices, this element is the space of maximum polynomial degree k, otherwise it is a tensor product of 1D polynomials, each with maximal degree k.

2183: .seealso: `PetscFECreateLagrangeByCell()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2184: @*/
2185: PetscErrorCode PetscFECreateLagrange(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, PetscInt k, PetscInt qorder, PetscFE *fem)
2186: {
2187:   PetscFunctionBegin;
2188:   PetscCall(PetscFECreate_Internal(comm, dim, Nc, DMPolytopeTypeSimpleShape(dim, isSimplex), NULL, k, qorder, PETSC_FALSE, fem));
2189:   PetscFunctionReturn(PETSC_SUCCESS);
2190: }

2192: /*@
2193:   PetscFECreateLagrangeByCell - Create a `PetscFE` for the basic Lagrange space of degree k

2195:   Collective

2197:   Input Parameters:
2198: + comm   - The MPI comm
2199: . dim    - The spatial dimension
2200: . Nc     - The number of components
2201: . ct     - The celltype of the reference cell
2202: . k      - The degree k of the space
2203: - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree

2205:   Output Parameter:
2206: . fem - The `PetscFE` object

2208:   Level: beginner

2210:   Note:
2211:   For simplices, this element is the space of maximum polynomial degree k, otherwise it is a tensor product of 1D polynomials, each with maximal degree k.

2213: .seealso: `PetscFECreateLagrange()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2214: @*/
2215: PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, PetscInt k, PetscInt qorder, PetscFE *fem)
2216: {
2217:   PetscFunctionBegin;
2218:   PetscCall(PetscFECreate_Internal(comm, dim, Nc, ct, NULL, k, qorder, PETSC_FALSE, fem));
2219:   PetscFunctionReturn(PETSC_SUCCESS);
2220: }

2222: /*@
2223:   PetscFELimitDegree - Copy a `PetscFE` but limit the degree to be in the given range

2225:   Collective

2227:   Input Parameters:
2228: + fe        - The `PetscFE`
2229: . minDegree - The minimum degree, or `PETSC_DETERMINE` for no limit
2230: - maxDegree - The maximum degree, or `PETSC_DETERMINE` for no limit

2232:   Output Parameter:
2233: . newfe - The `PetscFE` object

2235:   Level: advanced

2237:   Note:
2238:   This currently only works for Lagrange elements.

2240: .seealso: `PetscFECreateLagrange()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2241: @*/
2242: PetscErrorCode PetscFELimitDegree(PetscFE fe, PetscInt minDegree, PetscInt maxDegree, PetscFE *newfe)
2243: {
2244:   PetscDualSpace Q;
2245:   PetscBool      islag, issum;
2246:   PetscInt       oldk = 0, k;

2248:   PetscFunctionBegin;
2249:   PetscCall(PetscFEGetDualSpace(fe, &Q));
2250:   PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACELAGRANGE, &islag));
2251:   PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACESUM, &issum));
2252:   if (islag) {
2253:     PetscCall(PetscDualSpaceGetOrder(Q, &oldk));
2254:   } else if (issum) {
2255:     PetscDualSpace subQ;

2257:     PetscCall(PetscDualSpaceSumGetSubspace(Q, 0, &subQ));
2258:     PetscCall(PetscDualSpaceGetOrder(subQ, &oldk));
2259:   } else {
2260:     PetscCall(PetscObjectReference((PetscObject)fe));
2261:     *newfe = fe;
2262:     PetscFunctionReturn(PETSC_SUCCESS);
2263:   }
2264:   k = oldk;
2265:   if (minDegree >= 0) k = PetscMax(k, minDegree);
2266:   if (maxDegree >= 0) k = PetscMin(k, maxDegree);
2267:   if (k != oldk) {
2268:     DM              K;
2269:     PetscSpace      P;
2270:     PetscQuadrature q;
2271:     DMPolytopeType  ct;
2272:     PetscInt        dim, Nc;

2274:     PetscCall(PetscFEGetBasisSpace(fe, &P));
2275:     PetscCall(PetscSpaceGetNumVariables(P, &dim));
2276:     PetscCall(PetscSpaceGetNumComponents(P, &Nc));
2277:     PetscCall(PetscDualSpaceGetDM(Q, &K));
2278:     PetscCall(DMPlexGetCellType(K, 0, &ct));
2279:     PetscCall(PetscFECreateLagrangeByCell(PetscObjectComm((PetscObject)fe), dim, Nc, ct, k, PETSC_DETERMINE, newfe));
2280:     PetscCall(PetscFEGetQuadrature(fe, &q));
2281:     PetscCall(PetscFESetQuadrature(*newfe, q));
2282:   } else {
2283:     PetscCall(PetscObjectReference((PetscObject)fe));
2284:     *newfe = fe;
2285:   }
2286:   PetscFunctionReturn(PETSC_SUCCESS);
2287: }

2289: /*@
2290:   PetscFECreateBrokenElement - Create a discontinuous version of the input `PetscFE`

2292:   Collective

2294:   Input Parameters:
2295: . cgfe - The continuous `PetscFE` object

2297:   Output Parameter:
2298: . dgfe - The discontinuous `PetscFE` object

2300:   Level: advanced

2302:   Note:
2303:   This only works for Lagrange elements.

2305: .seealso: `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`, `PetscFECreateLagrange()`, `PetscFECreateLagrangeByCell()`, `PetscDualSpaceLagrangeSetContinuity()`
2306: @*/
2307: PetscErrorCode PetscFECreateBrokenElement(PetscFE cgfe, PetscFE *dgfe)
2308: {
2309:   PetscSpace      P;
2310:   PetscDualSpace  Q, dgQ;
2311:   PetscQuadrature q, fq;
2312:   PetscBool       is_lagrange, is_sum;

2314:   PetscFunctionBegin;
2315:   PetscCall(PetscFEGetBasisSpace(cgfe, &P));
2316:   PetscCall(PetscObjectReference((PetscObject)P));
2317:   PetscCall(PetscFEGetDualSpace(cgfe, &Q));
2318:   PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACELAGRANGE, &is_lagrange));
2319:   PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACESUM, &is_sum));
2320:   PetscCheck(is_lagrange || is_sum, PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only create broken elements of Lagrange elements");
2321:   PetscCall(PetscDualSpaceDuplicate(Q, &dgQ));
2322:   PetscCall(PetscDualSpaceLagrangeSetContinuity(dgQ, PETSC_FALSE));
2323:   PetscCall(PetscDualSpaceSetUp(dgQ));
2324:   PetscCall(PetscFEGetQuadrature(cgfe, &q));
2325:   PetscCall(PetscObjectReference((PetscObject)q));
2326:   PetscCall(PetscFEGetFaceQuadrature(cgfe, &fq));
2327:   PetscCall(PetscObjectReference((PetscObject)fq));
2328:   PetscCall(PetscFECreateFromSpaces(P, dgQ, q, fq, dgfe));
2329:   PetscFunctionReturn(PETSC_SUCCESS);
2330: }

2332: /*@
2333:   PetscFESetName - Names the `PetscFE` and its subobjects

2335:   Not Collective

2337:   Input Parameters:
2338: + fe   - The `PetscFE`
2339: - name - The name

2341:   Level: intermediate

2343: .seealso: `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2344: @*/
2345: PetscErrorCode PetscFESetName(PetscFE fe, const char name[])
2346: {
2347:   PetscSpace     P;
2348:   PetscDualSpace Q;

2350:   PetscFunctionBegin;
2351:   PetscCall(PetscFEGetBasisSpace(fe, &P));
2352:   PetscCall(PetscFEGetDualSpace(fe, &Q));
2353:   PetscCall(PetscObjectSetName((PetscObject)fe, name));
2354:   PetscCall(PetscObjectSetName((PetscObject)P, name));
2355:   PetscCall(PetscObjectSetName((PetscObject)Q, name));
2356:   PetscFunctionReturn(PETSC_SUCCESS);
2357: }

2359: PetscErrorCode PetscFEEvaluateFieldJets_Internal(PetscDS ds, PetscInt Nf, PetscInt r, PetscInt q, PetscTabulation T[], PetscFEGeom *fegeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscScalar u[], PetscScalar u_x[], PetscScalar u_t[])
2360: {
2361:   PetscInt dOffset = 0, fOffset = 0, f, g;

2363:   for (f = 0; f < Nf; ++f) {
2364:     PetscCheck(r < T[f]->Nr, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Replica number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", r, T[f]->Nr);
2365:     PetscCheck(q < T[f]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", q, T[f]->Np);
2366:     PetscFE          fe;
2367:     const PetscInt   k       = ds->jetDegree[f];
2368:     const PetscInt   cdim    = T[f]->cdim;
2369:     const PetscInt   dE      = fegeom->dimEmbed;
2370:     const PetscInt   Nq      = T[f]->Np;
2371:     const PetscInt   Nbf     = T[f]->Nb;
2372:     const PetscInt   Ncf     = T[f]->Nc;
2373:     const PetscReal *Bq      = &T[f]->T[0][(r * Nq + q) * Nbf * Ncf];
2374:     const PetscReal *Dq      = &T[f]->T[1][(r * Nq + q) * Nbf * Ncf * cdim];
2375:     const PetscReal *Hq      = k > 1 ? &T[f]->T[2][(r * Nq + q) * Nbf * Ncf * cdim * cdim] : NULL;
2376:     PetscInt         hOffset = 0, b, c, d;

2378:     PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
2379:     for (c = 0; c < Ncf; ++c) u[fOffset + c] = 0.0;
2380:     for (d = 0; d < dE * Ncf; ++d) u_x[fOffset * dE + d] = 0.0;
2381:     for (b = 0; b < Nbf; ++b) {
2382:       for (c = 0; c < Ncf; ++c) {
2383:         const PetscInt cidx = b * Ncf + c;

2385:         u[fOffset + c] += Bq[cidx] * coefficients[dOffset + b];
2386:         for (d = 0; d < cdim; ++d) u_x[(fOffset + c) * dE + d] += Dq[cidx * cdim + d] * coefficients[dOffset + b];
2387:       }
2388:     }
2389:     if (k > 1) {
2390:       for (g = 0; g < Nf; ++g) hOffset += T[g]->Nc * dE;
2391:       for (d = 0; d < dE * dE * Ncf; ++d) u_x[hOffset + fOffset * dE * dE + d] = 0.0;
2392:       for (b = 0; b < Nbf; ++b) {
2393:         for (c = 0; c < Ncf; ++c) {
2394:           const PetscInt cidx = b * Ncf + c;

2396:           for (d = 0; d < cdim * cdim; ++d) u_x[hOffset + (fOffset + c) * dE * dE + d] += Hq[cidx * cdim * cdim + d] * coefficients[dOffset + b];
2397:         }
2398:       }
2399:       PetscCall(PetscFEPushforwardHessian(fe, fegeom, 1, &u_x[hOffset + fOffset * dE * dE]));
2400:     }
2401:     PetscCall(PetscFEPushforward(fe, fegeom, 1, &u[fOffset]));
2402:     PetscCall(PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset * dE]));
2403:     if (u_t) {
2404:       for (c = 0; c < Ncf; ++c) u_t[fOffset + c] = 0.0;
2405:       for (b = 0; b < Nbf; ++b) {
2406:         for (c = 0; c < Ncf; ++c) {
2407:           const PetscInt cidx = b * Ncf + c;

2409:           u_t[fOffset + c] += Bq[cidx] * coefficients_t[dOffset + b];
2410:         }
2411:       }
2412:       PetscCall(PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]));
2413:     }
2414:     fOffset += Ncf;
2415:     dOffset += Nbf;
2416:   }
2417:   return PETSC_SUCCESS;
2418: }

2420: PetscErrorCode PetscFEEvaluateFieldJets_Hybrid_Internal(PetscDS ds, PetscInt Nf, PetscInt rc, PetscInt qc, PetscTabulation Tab[], const PetscInt rf[], const PetscInt qf[], PetscTabulation Tabf[], PetscFEGeom *fegeom, PetscFEGeom *fegeomNbr, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscScalar u[], PetscScalar u_x[], PetscScalar u_t[])
2421: {
2422:   PetscInt dOffset = 0, fOffset = 0, f;

2424:   /* f is the field number in the DS */
2425:   for (f = 0; f < Nf; ++f) {
2426:     PetscBool isCohesive;
2427:     PetscInt  Ns;

2429:     if (!Tab[f]) continue;
2430:     PetscCall(PetscDSGetCohesive(ds, f, &isCohesive));
2431:     Ns = isCohesive ? 1 : 2;
2432:     {
2433:       PetscTabulation T   = isCohesive ? Tab[f] : Tabf[f];
2434:       PetscFE         fe  = (PetscFE)ds->disc[f];
2435:       const PetscInt  dEt = T->cdim;
2436:       const PetscInt  dE  = fegeom->dimEmbed;
2437:       const PetscInt  Nq  = T->Np;
2438:       const PetscInt  Nbf = T->Nb;
2439:       const PetscInt  Ncf = T->Nc;

2441:       for (PetscInt s = 0; s < Ns; ++s) {
2442:         const PetscInt   r  = isCohesive ? rc : rf[s];
2443:         const PetscInt   q  = isCohesive ? qc : qf[s];
2444:         const PetscReal *Bq = &T->T[0][(r * Nq + q) * Nbf * Ncf];
2445:         const PetscReal *Dq = &T->T[1][(r * Nq + q) * Nbf * Ncf * dEt];
2446:         PetscInt         b, c, d;

2448:         PetscCheck(r < T->Nr, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " Side %" PetscInt_FMT " Replica number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", f, s, r, T->Nr);
2449:         PetscCheck(q < T->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " Side %" PetscInt_FMT " Point number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", f, s, q, T->Np);
2450:         for (c = 0; c < Ncf; ++c) u[fOffset + c] = 0.0;
2451:         for (d = 0; d < dE * Ncf; ++d) u_x[fOffset * dE + d] = 0.0;
2452:         for (b = 0; b < Nbf; ++b) {
2453:           for (c = 0; c < Ncf; ++c) {
2454:             const PetscInt cidx = b * Ncf + c;

2456:             u[fOffset + c] += Bq[cidx] * coefficients[dOffset + b];
2457:             for (d = 0; d < dEt; ++d) u_x[(fOffset + c) * dE + d] += Dq[cidx * dEt + d] * coefficients[dOffset + b];
2458:           }
2459:         }
2460:         PetscCall(PetscFEPushforward(fe, isCohesive ? fegeom : &fegeomNbr[s], 1, &u[fOffset]));
2461:         PetscCall(PetscFEPushforwardGradient(fe, isCohesive ? fegeom : &fegeomNbr[s], 1, &u_x[fOffset * dE]));
2462:         if (u_t) {
2463:           for (c = 0; c < Ncf; ++c) u_t[fOffset + c] = 0.0;
2464:           for (b = 0; b < Nbf; ++b) {
2465:             for (c = 0; c < Ncf; ++c) {
2466:               const PetscInt cidx = b * Ncf + c;

2468:               u_t[fOffset + c] += Bq[cidx] * coefficients_t[dOffset + b];
2469:             }
2470:           }
2471:           PetscCall(PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]));
2472:         }
2473:         fOffset += Ncf;
2474:         dOffset += Nbf;
2475:       }
2476:     }
2477:   }
2478:   return PETSC_SUCCESS;
2479: }

2481: PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[])
2482: {
2483:   PetscFE         fe;
2484:   PetscTabulation Tc;
2485:   PetscInt        b, c;

2487:   if (!prob) return PETSC_SUCCESS;
2488:   PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe));
2489:   PetscCall(PetscFEGetFaceCentroidTabulation(fe, &Tc));
2490:   {
2491:     const PetscReal *faceBasis = Tc->T[0];
2492:     const PetscInt   Nb        = Tc->Nb;
2493:     const PetscInt   Nc        = Tc->Nc;

2495:     for (c = 0; c < Nc; ++c) u[c] = 0.0;
2496:     for (b = 0; b < Nb; ++b) {
2497:       for (c = 0; c < Nc; ++c) u[c] += coefficients[b] * faceBasis[(faceLoc * Nb + b) * Nc + c];
2498:     }
2499:   }
2500:   return PETSC_SUCCESS;
2501: }

2503: PetscErrorCode PetscFEUpdateElementVec_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscInt e, PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[])
2504: {
2505:   PetscFEGeom      pgeom;
2506:   const PetscInt   dEt      = T->cdim;
2507:   const PetscInt   dE       = fegeom->dimEmbed;
2508:   const PetscInt   Nq       = T->Np;
2509:   const PetscInt   Nb       = T->Nb;
2510:   const PetscInt   Nc       = T->Nc;
2511:   const PetscReal *basis    = &T->T[0][r * Nq * Nb * Nc];
2512:   const PetscReal *basisDer = &T->T[1][r * Nq * Nb * Nc * dEt];
2513:   PetscInt         q, b, c, d;

2515:   for (q = 0; q < Nq; ++q) {
2516:     for (b = 0; b < Nb; ++b) {
2517:       for (c = 0; c < Nc; ++c) {
2518:         const PetscInt bcidx = b * Nc + c;

2520:         tmpBasis[bcidx] = basis[q * Nb * Nc + bcidx];
2521:         for (d = 0; d < dEt; ++d) tmpBasisDer[bcidx * dE + d] = basisDer[q * Nb * Nc * dEt + bcidx * dEt + d];
2522:         for (d = dEt; d < dE; ++d) tmpBasisDer[bcidx * dE + d] = 0.0;
2523:       }
2524:     }
2525:     PetscCall(PetscFEGeomGetCellPoint(fegeom, e, q, &pgeom));
2526:     PetscCall(PetscFEPushforward(fe, &pgeom, Nb, tmpBasis));
2527:     PetscCall(PetscFEPushforwardGradient(fe, &pgeom, Nb, tmpBasisDer));
2528:     for (b = 0; b < Nb; ++b) {
2529:       for (c = 0; c < Nc; ++c) {
2530:         const PetscInt bcidx = b * Nc + c;
2531:         const PetscInt qcidx = q * Nc + c;

2533:         elemVec[b] += tmpBasis[bcidx] * f0[qcidx];
2534:         for (d = 0; d < dE; ++d) elemVec[b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d];
2535:       }
2536:     }
2537:   }
2538:   return PETSC_SUCCESS;
2539: }

2541: PetscErrorCode PetscFEUpdateElementVec_Hybrid_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscInt side, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[])
2542: {
2543:   const PetscInt   dE       = T->cdim;
2544:   const PetscInt   Nq       = T->Np;
2545:   const PetscInt   Nb       = T->Nb;
2546:   const PetscInt   Nc       = T->Nc;
2547:   const PetscReal *basis    = &T->T[0][r * Nq * Nb * Nc];
2548:   const PetscReal *basisDer = &T->T[1][r * Nq * Nb * Nc * dE];

2550:   for (PetscInt q = 0; q < Nq; ++q) {
2551:     for (PetscInt b = 0; b < Nb; ++b) {
2552:       for (PetscInt c = 0; c < Nc; ++c) {
2553:         const PetscInt bcidx = b * Nc + c;

2555:         tmpBasis[bcidx] = basis[q * Nb * Nc + bcidx];
2556:         for (PetscInt d = 0; d < dE; ++d) tmpBasisDer[bcidx * dE + d] = basisDer[q * Nb * Nc * dE + bcidx * dE + d];
2557:       }
2558:     }
2559:     PetscCall(PetscFEPushforward(fe, fegeom, Nb, tmpBasis));
2560:     // TODO This is currently broken since we do not pull the geometry down to the lower dimension
2561:     // PetscCall(PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer));
2562:     if (side == 2) {
2563:       // Integrating over whole cohesive cell, so insert for both sides
2564:       for (PetscInt s = 0; s < 2; ++s) {
2565:         for (PetscInt b = 0; b < Nb; ++b) {
2566:           for (PetscInt c = 0; c < Nc; ++c) {
2567:             const PetscInt bcidx = b * Nc + c;
2568:             const PetscInt qcidx = (q * 2 + s) * Nc + c;

2570:             elemVec[Nb * s + b] += tmpBasis[bcidx] * f0[qcidx];
2571:             for (PetscInt d = 0; d < dE; ++d) elemVec[Nb * s + b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d];
2572:           }
2573:         }
2574:       }
2575:     } else {
2576:       // Integrating over endcaps of cohesive cell, so insert for correct side
2577:       for (PetscInt b = 0; b < Nb; ++b) {
2578:         for (PetscInt c = 0; c < Nc; ++c) {
2579:           const PetscInt bcidx = b * Nc + c;
2580:           const PetscInt qcidx = q * Nc + c;

2582:           elemVec[Nb * side + b] += tmpBasis[bcidx] * f0[qcidx];
2583:           for (PetscInt d = 0; d < dE; ++d) elemVec[Nb * side + b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d];
2584:         }
2585:       }
2586:     }
2587:   }
2588:   return PETSC_SUCCESS;
2589: }

2591: #define petsc_elemmat_kernel_g1(_NbI, _NcI, _NbJ, _NcJ, _dE) \
2592:   do { \
2593:     for (PetscInt fc = 0; fc < (_NcI); ++fc) { \
2594:       for (PetscInt gc = 0; gc < (_NcJ); ++gc) { \
2595:         const PetscScalar *G = g1 + (fc * (_NcJ) + gc) * _dE; \
2596:         for (PetscInt f = 0; f < (_NbI); ++f) { \
2597:           const PetscScalar tBIv = tmpBasisI[f * (_NcI) + fc]; \
2598:           for (PetscInt g = 0; g < (_NbJ); ++g) { \
2599:             const PetscScalar *tBDJ = tmpBasisDerJ + (g * (_NcJ) + gc) * (_dE); \
2600:             PetscScalar        s    = 0.0; \
2601:             for (PetscInt df = 0; df < _dE; ++df) s += G[df] * tBDJ[df]; \
2602:             elemMat[(offsetI + f) * totDim + (offsetJ + g)] += s * tBIv; \
2603:           } \
2604:         } \
2605:       } \
2606:     } \
2607:   } while (0)

2609: #define petsc_elemmat_kernel_g2(_NbI, _NcI, _NbJ, _NcJ, _dE) \
2610:   do { \
2611:     for (PetscInt gc = 0; gc < (_NcJ); ++gc) { \
2612:       for (PetscInt fc = 0; fc < (_NcI); ++fc) { \
2613:         const PetscScalar *G = g2 + (fc * (_NcJ) + gc) * _dE; \
2614:         for (PetscInt g = 0; g < (_NbJ); ++g) { \
2615:           const PetscScalar tBJv = tmpBasisJ[g * (_NcJ) + gc]; \
2616:           for (PetscInt f = 0; f < (_NbI); ++f) { \
2617:             const PetscScalar *tBDI = tmpBasisDerI + (f * (_NcI) + fc) * (_dE); \
2618:             PetscScalar        s    = 0.0; \
2619:             for (PetscInt df = 0; df < _dE; ++df) s += tBDI[df] * G[df]; \
2620:             elemMat[(offsetI + f) * totDim + (offsetJ + g)] += s * tBJv; \
2621:           } \
2622:         } \
2623:       } \
2624:     } \
2625:   } while (0)

2627: #define petsc_elemmat_kernel_g3(_NbI, _NcI, _NbJ, _NcJ, _dE) \
2628:   do { \
2629:     for (PetscInt fc = 0; fc < (_NcI); ++fc) { \
2630:       for (PetscInt gc = 0; gc < (_NcJ); ++gc) { \
2631:         const PetscScalar *G = g3 + (fc * (_NcJ) + gc) * (_dE) * (_dE); \
2632:         for (PetscInt f = 0; f < (_NbI); ++f) { \
2633:           const PetscScalar *tBDI = tmpBasisDerI + (f * (_NcI) + fc) * (_dE); \
2634:           for (PetscInt g = 0; g < (_NbJ); ++g) { \
2635:             PetscScalar        s    = 0.0; \
2636:             const PetscScalar *tBDJ = tmpBasisDerJ + (g * (_NcJ) + gc) * (_dE); \
2637:             for (PetscInt df = 0; df < (_dE); ++df) { \
2638:               for (PetscInt dg = 0; dg < (_dE); ++dg) s += tBDI[df] * G[df * (_dE) + dg] * tBDJ[dg]; \
2639:             } \
2640:             elemMat[(offsetI + f) * totDim + (offsetJ + g)] += s; \
2641:           } \
2642:         } \
2643:       } \
2644:     } \
2645:   } while (0)

2647: PetscErrorCode PetscFEUpdateElementMat_Internal(PetscFE feI, PetscFE feJ, PetscInt r, PetscInt q, PetscTabulation TI, PetscScalar tmpBasisI[], PetscScalar tmpBasisDerI[], PetscTabulation TJ, PetscScalar tmpBasisJ[], PetscScalar tmpBasisDerJ[], PetscFEGeom *fegeom, const PetscScalar g0[], const PetscScalar g1[], const PetscScalar g2[], const PetscScalar g3[], PetscInt totDim, PetscInt offsetI, PetscInt offsetJ, PetscScalar elemMat[])
2648: {
2649:   const PetscInt   cdim      = TI->cdim;
2650:   const PetscInt   dE        = fegeom->dimEmbed;
2651:   const PetscInt   NqI       = TI->Np;
2652:   const PetscInt   NbI       = TI->Nb;
2653:   const PetscInt   NcI       = TI->Nc;
2654:   const PetscReal *basisI    = &TI->T[0][(r * NqI + q) * NbI * NcI];
2655:   const PetscReal *basisDerI = &TI->T[1][(r * NqI + q) * NbI * NcI * cdim];
2656:   const PetscInt   NqJ       = TJ->Np;
2657:   const PetscInt   NbJ       = TJ->Nb;
2658:   const PetscInt   NcJ       = TJ->Nc;
2659:   const PetscReal *basisJ    = &TJ->T[0][(r * NqJ + q) * NbJ * NcJ];
2660:   const PetscReal *basisDerJ = &TJ->T[1][(r * NqJ + q) * NbJ * NcJ * cdim];

2662:   for (PetscInt f = 0; f < NbI; ++f) {
2663:     for (PetscInt fc = 0; fc < NcI; ++fc) {
2664:       const PetscInt fidx = f * NcI + fc; /* Test function basis index */

2666:       tmpBasisI[fidx] = basisI[fidx];
2667:       for (PetscInt df = 0; df < cdim; ++df) tmpBasisDerI[fidx * dE + df] = basisDerI[fidx * cdim + df];
2668:     }
2669:   }
2670:   PetscCall(PetscFEPushforward(feI, fegeom, NbI, tmpBasisI));
2671:   PetscCall(PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI));
2672:   if (feI != feJ) {
2673:     for (PetscInt g = 0; g < NbJ; ++g) {
2674:       for (PetscInt gc = 0; gc < NcJ; ++gc) {
2675:         const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */

2677:         tmpBasisJ[gidx] = basisJ[gidx];
2678:         for (PetscInt dg = 0; dg < cdim; ++dg) tmpBasisDerJ[gidx * dE + dg] = basisDerJ[gidx * cdim + dg];
2679:       }
2680:     }
2681:     PetscCall(PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ));
2682:     PetscCall(PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ));
2683:   } else {
2684:     tmpBasisJ    = tmpBasisI;
2685:     tmpBasisDerJ = tmpBasisDerI;
2686:   }
2687:   if (PetscUnlikely(g0)) {
2688:     for (PetscInt f = 0; f < NbI; ++f) {
2689:       const PetscInt i = offsetI + f; /* Element matrix row */

2691:       for (PetscInt fc = 0; fc < NcI; ++fc) {
2692:         const PetscScalar bI = tmpBasisI[f * NcI + fc]; /* Test function basis value */

2694:         for (PetscInt g = 0; g < NbJ; ++g) {
2695:           const PetscInt j    = offsetJ + g; /* Element matrix column */
2696:           const PetscInt fOff = i * totDim + j;

2698:           for (PetscInt gc = 0; gc < NcJ; ++gc) elemMat[fOff] += bI * g0[fc * NcJ + gc] * tmpBasisJ[g * NcJ + gc];
2699:         }
2700:       }
2701:     }
2702:   }
2703:   if (PetscUnlikely(g1)) {
2704: #if 1
2705:     if (dE == 2) {
2706:       petsc_elemmat_kernel_g1(NbI, NcI, NbJ, NcJ, 2);
2707:     } else if (dE == 3) {
2708:       petsc_elemmat_kernel_g1(NbI, NcI, NbJ, NcJ, 3);
2709:     } else {
2710:       petsc_elemmat_kernel_g1(NbI, NcI, NbJ, NcJ, dE);
2711:     }
2712: #else
2713:     for (PetscInt f = 0; f < NbI; ++f) {
2714:       const PetscInt i = offsetI + f; /* Element matrix row */

2716:       for (PetscInt fc = 0; fc < NcI; ++fc) {
2717:         const PetscScalar bI = tmpBasisI[f * NcI + fc]; /* Test function basis value */

2719:         for (PetscInt g = 0; g < NbJ; ++g) {
2720:           const PetscInt j    = offsetJ + g; /* Element matrix column */
2721:           const PetscInt fOff = i * totDim + j;

2723:           for (PetscInt gc = 0; gc < NcJ; ++gc) {
2724:             const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */

2726:             for (PetscInt df = 0; df < dE; ++df) elemMat[fOff] += bI * g1[(fc * NcJ + gc) * dE + df] * tmpBasisDerJ[gidx * dE + df];
2727:           }
2728:         }
2729:       }
2730:     }
2731: #endif
2732:   }
2733:   if (PetscUnlikely(g2)) {
2734: #if 1
2735:     if (dE == 2) {
2736:       petsc_elemmat_kernel_g2(NbI, NcI, NbJ, NcJ, 2);
2737:     } else if (dE == 3) {
2738:       petsc_elemmat_kernel_g2(NbI, NcI, NbJ, NcJ, 3);
2739:     } else {
2740:       petsc_elemmat_kernel_g2(NbI, NcI, NbJ, NcJ, dE);
2741:     }
2742: #else
2743:     for (PetscInt g = 0; g < NbJ; ++g) {
2744:       const PetscInt j = offsetJ + g; /* Element matrix column */

2746:       for (PetscInt gc = 0; gc < NcJ; ++gc) {
2747:         const PetscScalar bJ = tmpBasisJ[g * NcJ + gc]; /* Trial function basis value */

2749:         for (PetscInt f = 0; f < NbI; ++f) {
2750:           const PetscInt i    = offsetI + f; /* Element matrix row */
2751:           const PetscInt fOff = i * totDim + j;

2753:           for (PetscInt fc = 0; fc < NcI; ++fc) {
2754:             const PetscInt fidx = f * NcI + fc; /* Test function basis index */

2756:             for (PetscInt df = 0; df < dE; ++df) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g2[(fc * NcJ + gc) * dE + df] * bJ;
2757:           }
2758:         }
2759:       }
2760:     }
2761: #endif
2762:   }
2763:   if (PetscUnlikely(g3)) {
2764: #if 1
2765:     if (dE == 2) {
2766:       petsc_elemmat_kernel_g3(NbI, NcI, NbJ, NcJ, 2);
2767:     } else if (dE == 3) {
2768:       petsc_elemmat_kernel_g3(NbI, NcI, NbJ, NcJ, 3);
2769:     } else {
2770:       petsc_elemmat_kernel_g3(NbI, NcI, NbJ, NcJ, dE);
2771:     }
2772: #else
2773:     for (PetscInt f = 0; f < NbI; ++f) {
2774:       const PetscInt i = offsetI + f; /* Element matrix row */

2776:       for (PetscInt fc = 0; fc < NcI; ++fc) {
2777:         const PetscInt fidx = f * NcI + fc; /* Test function basis index */

2779:         for (PetscInt g = 0; g < NbJ; ++g) {
2780:           const PetscInt j    = offsetJ + g; /* Element matrix column */
2781:           const PetscInt fOff = i * totDim + j;

2783:           for (PetscInt gc = 0; gc < NcJ; ++gc) {
2784:             const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */

2786:             for (PetscInt df = 0; df < dE; ++df) {
2787:               for (PetscInt dg = 0; dg < dE; ++dg) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g3[((fc * NcJ + gc) * dE + df) * dE + dg] * tmpBasisDerJ[gidx * dE + dg];
2788:             }
2789:           }
2790:         }
2791:       }
2792:     }
2793: #endif
2794:   }
2795:   return PETSC_SUCCESS;
2796: }

2798: #undef petsc_elemmat_kernel_g1
2799: #undef petsc_elemmat_kernel_g2
2800: #undef petsc_elemmat_kernel_g3

2802: PetscErrorCode PetscFEUpdateElementMat_Hybrid_Internal(PetscFE feI, PetscBool isHybridI, PetscFE feJ, PetscBool isHybridJ, PetscInt r, PetscInt s, PetscInt t, PetscInt q, PetscTabulation TI, PetscScalar tmpBasisI[], PetscScalar tmpBasisDerI[], PetscTabulation TJ, PetscScalar tmpBasisJ[], PetscScalar tmpBasisDerJ[], PetscFEGeom *fegeom, const PetscScalar g0[], const PetscScalar g1[], const PetscScalar g2[], const PetscScalar g3[], PetscInt eOffset, PetscInt totDim, PetscInt offsetI, PetscInt offsetJ, PetscScalar elemMat[])
2803: {
2804:   const PetscInt   dE        = TI->cdim;
2805:   const PetscInt   NqI       = TI->Np;
2806:   const PetscInt   NbI       = TI->Nb;
2807:   const PetscInt   NcI       = TI->Nc;
2808:   const PetscReal *basisI    = &TI->T[0][(r * NqI + q) * NbI * NcI];
2809:   const PetscReal *basisDerI = &TI->T[1][(r * NqI + q) * NbI * NcI * dE];
2810:   const PetscInt   NqJ       = TJ->Np;
2811:   const PetscInt   NbJ       = TJ->Nb;
2812:   const PetscInt   NcJ       = TJ->Nc;
2813:   const PetscReal *basisJ    = &TJ->T[0][(r * NqJ + q) * NbJ * NcJ];
2814:   const PetscReal *basisDerJ = &TJ->T[1][(r * NqJ + q) * NbJ * NcJ * dE];
2815:   const PetscInt   so        = isHybridI ? 0 : s;
2816:   const PetscInt   to        = isHybridJ ? 0 : t;
2817:   PetscInt         f, fc, g, gc, df, dg;

2819:   for (f = 0; f < NbI; ++f) {
2820:     for (fc = 0; fc < NcI; ++fc) {
2821:       const PetscInt fidx = f * NcI + fc; /* Test function basis index */

2823:       tmpBasisI[fidx] = basisI[fidx];
2824:       for (df = 0; df < dE; ++df) tmpBasisDerI[fidx * dE + df] = basisDerI[fidx * dE + df];
2825:     }
2826:   }
2827:   PetscCall(PetscFEPushforward(feI, fegeom, NbI, tmpBasisI));
2828:   PetscCall(PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI));
2829:   for (g = 0; g < NbJ; ++g) {
2830:     for (gc = 0; gc < NcJ; ++gc) {
2831:       const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */

2833:       tmpBasisJ[gidx] = basisJ[gidx];
2834:       for (dg = 0; dg < dE; ++dg) tmpBasisDerJ[gidx * dE + dg] = basisDerJ[gidx * dE + dg];
2835:     }
2836:   }
2837:   PetscCall(PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ));
2838:   // TODO This is currently broken since we do not pull the geometry down to the lower dimension
2839:   // PetscCall(PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ));
2840:   for (f = 0; f < NbI; ++f) {
2841:     for (fc = 0; fc < NcI; ++fc) {
2842:       const PetscInt fidx = f * NcI + fc;           /* Test function basis index */
2843:       const PetscInt i    = offsetI + NbI * so + f; /* Element matrix row */
2844:       for (g = 0; g < NbJ; ++g) {
2845:         for (gc = 0; gc < NcJ; ++gc) {
2846:           const PetscInt gidx = g * NcJ + gc;           /* Trial function basis index */
2847:           const PetscInt j    = offsetJ + NbJ * to + g; /* Element matrix column */
2848:           const PetscInt fOff = eOffset + i * totDim + j;

2850:           elemMat[fOff] += tmpBasisI[fidx] * g0[fc * NcJ + gc] * tmpBasisJ[gidx];
2851:           for (df = 0; df < dE; ++df) {
2852:             elemMat[fOff] += tmpBasisI[fidx] * g1[(fc * NcJ + gc) * dE + df] * tmpBasisDerJ[gidx * dE + df];
2853:             elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g2[(fc * NcJ + gc) * dE + df] * tmpBasisJ[gidx];
2854:             for (dg = 0; dg < dE; ++dg) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g3[((fc * NcJ + gc) * dE + df) * dE + dg] * tmpBasisDerJ[gidx * dE + dg];
2855:           }
2856:         }
2857:       }
2858:     }
2859:   }
2860:   return PETSC_SUCCESS;
2861: }

2863: PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom)
2864: {
2865:   PetscDualSpace  dsp;
2866:   DM              dm;
2867:   PetscQuadrature quadDef;
2868:   PetscInt        dim, cdim, Nq;

2870:   PetscFunctionBegin;
2871:   PetscCall(PetscFEGetDualSpace(fe, &dsp));
2872:   PetscCall(PetscDualSpaceGetDM(dsp, &dm));
2873:   PetscCall(DMGetDimension(dm, &dim));
2874:   PetscCall(DMGetCoordinateDim(dm, &cdim));
2875:   PetscCall(PetscFEGetQuadrature(fe, &quadDef));
2876:   quad = quad ? quad : quadDef;
2877:   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
2878:   PetscCall(PetscMalloc1(Nq * cdim, &cgeom->v));
2879:   PetscCall(PetscMalloc1(Nq * cdim * cdim, &cgeom->J));
2880:   PetscCall(PetscMalloc1(Nq * cdim * cdim, &cgeom->invJ));
2881:   PetscCall(PetscMalloc1(Nq, &cgeom->detJ));
2882:   cgeom->dim       = dim;
2883:   cgeom->dimEmbed  = cdim;
2884:   cgeom->numCells  = 1;
2885:   cgeom->numPoints = Nq;
2886:   PetscCall(DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ));
2887:   PetscFunctionReturn(PETSC_SUCCESS);
2888: }

2890: PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom)
2891: {
2892:   PetscFunctionBegin;
2893:   PetscCall(PetscFree(cgeom->v));
2894:   PetscCall(PetscFree(cgeom->J));
2895:   PetscCall(PetscFree(cgeom->invJ));
2896:   PetscCall(PetscFree(cgeom->detJ));
2897:   PetscFunctionReturn(PETSC_SUCCESS);
2898: }

2900: #if 0
2901: PetscErrorCode PetscFEUpdateElementMat_Internal_SparseIndices(PetscTabulation TI, PetscTabulation TJ, PetscInt dimEmbed, const PetscInt g0[], const PetscInt g1[], const PetscInt g2[], const PetscInt g3[], PetscInt totDim, PetscInt offsetI, PetscInt offsetJ, PetscInt *n_g0, PetscInt **g0_idxs_out, PetscInt *n_g1, PetscInt **g1_idxs_out, PetscInt *n_g2, PetscInt **g2_idxs_out, PetscInt *n_g3, PetscInt **g3_idxs_out)
2902: {
2903:   const PetscInt dE      = dimEmbed;
2904:   const PetscInt NbI     = TI->Nb;
2905:   const PetscInt NcI     = TI->Nc;
2906:   const PetscInt NbJ     = TJ->Nb;
2907:   const PetscInt NcJ     = TJ->Nc;
2908:   PetscBool      has_g0  = g0 ? PETSC_TRUE : PETSC_FALSE;
2909:   PetscBool      has_g1  = g1 ? PETSC_TRUE : PETSC_FALSE;
2910:   PetscBool      has_g2  = g2 ? PETSC_TRUE : PETSC_FALSE;
2911:   PetscBool      has_g3  = g3 ? PETSC_TRUE : PETSC_FALSE;
2912:   PetscInt      *g0_idxs = NULL, *g1_idxs = NULL, *g2_idxs = NULL, *g3_idxs = NULL;
2913:   PetscInt       g0_i, g1_i, g2_i, g3_i;

2915:   PetscFunctionBegin;
2916:   g0_i = g1_i = g2_i = g3_i = 0;
2917:   if (has_g0)
2918:     for (PetscInt i = 0; i < NcI * NcJ; i++)
2919:       if (g0[i]) g0_i += NbI * NbJ;
2920:   if (has_g1)
2921:     for (PetscInt i = 0; i < NcI * NcJ * dE; i++)
2922:       if (g1[i]) g1_i += NbI * NbJ;
2923:   if (has_g2)
2924:     for (PetscInt i = 0; i < NcI * NcJ * dE; i++)
2925:       if (g2[i]) g2_i += NbI * NbJ;
2926:   if (has_g3)
2927:     for (PetscInt i = 0; i < NcI * NcJ * dE * dE; i++)
2928:       if (g3[i]) g3_i += NbI * NbJ;
2929:   if (g0_i == NbI * NbJ * NcI * NcJ) g0_i = 0;
2930:   if (g1_i == NbI * NbJ * NcI * NcJ * dE) g1_i = 0;
2931:   if (g2_i == NbI * NbJ * NcI * NcJ * dE) g2_i = 0;
2932:   if (g3_i == NbI * NbJ * NcI * NcJ * dE * dE) g3_i = 0;
2933:   has_g0 = g0_i ? PETSC_TRUE : PETSC_FALSE;
2934:   has_g1 = g1_i ? PETSC_TRUE : PETSC_FALSE;
2935:   has_g2 = g2_i ? PETSC_TRUE : PETSC_FALSE;
2936:   has_g3 = g3_i ? PETSC_TRUE : PETSC_FALSE;
2937:   if (has_g0) PetscCall(PetscMalloc1(4 * g0_i, &g0_idxs));
2938:   if (has_g1) PetscCall(PetscMalloc1(4 * g1_i, &g1_idxs));
2939:   if (has_g2) PetscCall(PetscMalloc1(4 * g2_i, &g2_idxs));
2940:   if (has_g3) PetscCall(PetscMalloc1(4 * g3_i, &g3_idxs));
2941:   g0_i = g1_i = g2_i = g3_i = 0;

2943:   for (PetscInt f = 0; f < NbI; ++f) {
2944:     const PetscInt i = offsetI + f; /* Element matrix row */
2945:     for (PetscInt fc = 0; fc < NcI; ++fc) {
2946:       const PetscInt fidx = f * NcI + fc; /* Test function basis index */

2948:       for (PetscInt g = 0; g < NbJ; ++g) {
2949:         const PetscInt j    = offsetJ + g; /* Element matrix column */
2950:         const PetscInt fOff = i * totDim + j;
2951:         for (PetscInt gc = 0; gc < NcJ; ++gc) {
2952:           const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */

2954:           if (has_g0) {
2955:             if (g0[fc * NcJ + gc]) {
2956:               g0_idxs[4 * g0_i + 0] = fidx;
2957:               g0_idxs[4 * g0_i + 1] = fc * NcJ + gc;
2958:               g0_idxs[4 * g0_i + 2] = gidx;
2959:               g0_idxs[4 * g0_i + 3] = fOff;
2960:               g0_i++;
2961:             }
2962:           }

2964:           for (PetscInt df = 0; df < dE; ++df) {
2965:             if (has_g1) {
2966:               if (g1[(fc * NcJ + gc) * dE + df]) {
2967:                 g1_idxs[4 * g1_i + 0] = fidx;
2968:                 g1_idxs[4 * g1_i + 1] = (fc * NcJ + gc) * dE + df;
2969:                 g1_idxs[4 * g1_i + 2] = gidx * dE + df;
2970:                 g1_idxs[4 * g1_i + 3] = fOff;
2971:                 g1_i++;
2972:               }
2973:             }
2974:             if (has_g2) {
2975:               if (g2[(fc * NcJ + gc) * dE + df]) {
2976:                 g2_idxs[4 * g2_i + 0] = fidx * dE + df;
2977:                 g2_idxs[4 * g2_i + 1] = (fc * NcJ + gc) * dE + df;
2978:                 g2_idxs[4 * g2_i + 2] = gidx;
2979:                 g2_idxs[4 * g2_i + 3] = fOff;
2980:                 g2_i++;
2981:               }
2982:             }
2983:             if (has_g3) {
2984:               for (PetscInt dg = 0; dg < dE; ++dg) {
2985:                 if (g3[((fc * NcJ + gc) * dE + df) * dE + dg]) {
2986:                   g3_idxs[4 * g3_i + 0] = fidx * dE + df;
2987:                   g3_idxs[4 * g3_i + 1] = ((fc * NcJ + gc) * dE + df) * dE + dg;
2988:                   g3_idxs[4 * g3_i + 2] = gidx * dE + dg;
2989:                   g3_idxs[4 * g3_i + 3] = fOff;
2990:                   g3_i++;
2991:                 }
2992:               }
2993:             }
2994:           }
2995:         }
2996:       }
2997:     }
2998:   }
2999:   *n_g0 = g0_i;
3000:   *n_g1 = g1_i;
3001:   *n_g2 = g2_i;
3002:   *n_g3 = g3_i;

3004:   *g0_idxs_out = g0_idxs;
3005:   *g1_idxs_out = g1_idxs;
3006:   *g2_idxs_out = g2_idxs;
3007:   *g3_idxs_out = g3_idxs;
3008:   PetscFunctionReturn(PETSC_SUCCESS);
3009: }

3011: //example HOW TO USE
3012:       for (PetscInt i = 0; i < g0_sparse_n; i++) {
3013:         PetscInt bM = g0_sparse_idxs[4 * i + 0];
3014:         PetscInt bN = g0_sparse_idxs[4 * i + 1];
3015:         PetscInt bK = g0_sparse_idxs[4 * i + 2];
3016:         PetscInt bO = g0_sparse_idxs[4 * i + 3];
3017:         elemMat[bO] += tmpBasisI[bM] * g0[bN] * tmpBasisJ[bK];
3018:       }
3019: #endif