Actual source code: fe.c

petsc-master 2019-08-18
Report Typos and Errors
  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: PetscFunctionList PetscFEList              = NULL;
 54: PetscBool         PetscFERegisterAllCalled = PETSC_FALSE;

 56: /*@C
 57:   PetscFERegister - Adds a new PetscFE implementation

 59:   Not Collective

 61:   Input Parameters:
 62: + name        - The name of a new user-defined creation routine
 63: - create_func - The creation routine itself

 65:   Notes:
 66:   PetscFERegister() may be called multiple times to add several user-defined PetscFEs

 68:   Sample usage:
 69: .vb
 70:     PetscFERegister("my_fe", MyPetscFECreate);
 71: .ve

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

 83:   Level: advanced

 85: .seealso: PetscFERegisterAll(), PetscFERegisterDestroy()

 87: @*/
 88: PetscErrorCode PetscFERegister(const char sname[], PetscErrorCode (*function)(PetscFE))
 89: {

 93:   PetscFunctionListAdd(&PetscFEList, sname, function);
 94:   return(0);
 95: }

 97: /*@C
 98:   PetscFESetType - Builds a particular PetscFE

100:   Collective on fem

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

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

109:   Level: intermediate

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

121:   PetscObjectTypeCompare((PetscObject) fem, name, &match);
122:   if (match) return(0);

124:   if (!PetscFERegisterAllCalled) {PetscFERegisterAll();}
125:   PetscFunctionListFind(PetscFEList, name, &r);
126:   if (!r) SETERRQ1(PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFE type: %s", name);

128:   if (fem->ops->destroy) {
129:     (*fem->ops->destroy)(fem);
130:     fem->ops->destroy = NULL;
131:   }
132:   (*r)(fem);
133:   PetscObjectChangeTypeName((PetscObject) fem, name);
134:   return(0);
135: }

137: /*@C
138:   PetscFEGetType - Gets the PetscFE type name (as a string) from the object.

140:   Not Collective

142:   Input Parameter:
143: . fem  - The PetscFE

145:   Output Parameter:
146: . name - The PetscFE type name

148:   Level: intermediate

150: .seealso: PetscFESetType(), PetscFECreate()
151: @*/
152: PetscErrorCode PetscFEGetType(PetscFE fem, PetscFEType *name)
153: {

159:   if (!PetscFERegisterAllCalled) {
160:     PetscFERegisterAll();
161:   }
162:   *name = ((PetscObject) fem)->type_name;
163:   return(0);
164: }

166: /*@C
167:   PetscFEView - Views a PetscFE

169:   Collective on fem

171:   Input Parameter:
172: + fem - the PetscFE object to view
173: - viewer   - the viewer

175:   Level: developer

177: .seealso PetscFEDestroy()
178: @*/
179: PetscErrorCode PetscFEView(PetscFE fem, PetscViewer viewer)
180: {
181:   PetscBool      iascii;

187:   if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &viewer);}
188:   PetscObjectPrintClassNamePrefixType((PetscObject)fem, viewer);
189:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
190:   if (fem->ops->view) {(*fem->ops->view)(fem, viewer);}
191:   return(0);
192: }

194: /*@
195:   PetscFESetFromOptions - sets parameters in a PetscFE from the options database

197:   Collective on fem

199:   Input Parameter:
200: . fem - the PetscFE object to set options for

202:   Options Database:
203: + -petscfe_num_blocks  - the number of cell blocks to integrate concurrently
204: - -petscfe_num_batches - the number of cell batches to integrate serially

206:   Level: developer

208: .seealso PetscFEView()
209: @*/
210: PetscErrorCode PetscFESetFromOptions(PetscFE fem)
211: {
212:   const char    *defaultType;
213:   char           name[256];
214:   PetscBool      flg;

219:   if (!((PetscObject) fem)->type_name) {
220:     defaultType = PETSCFEBASIC;
221:   } else {
222:     defaultType = ((PetscObject) fem)->type_name;
223:   }
224:   if (!PetscFERegisterAllCalled) {PetscFERegisterAll();}

226:   PetscObjectOptionsBegin((PetscObject) fem);
227:   PetscOptionsFList("-petscfe_type", "Finite element space", "PetscFESetType", PetscFEList, defaultType, name, 256, &flg);
228:   if (flg) {
229:     PetscFESetType(fem, name);
230:   } else if (!((PetscObject) fem)->type_name) {
231:     PetscFESetType(fem, defaultType);
232:   }
233:   PetscOptionsBoundedInt("-petscfe_num_blocks", "The number of cell blocks to integrate concurrently", "PetscSpaceSetTileSizes", fem->numBlocks, &fem->numBlocks, NULL,1);
234:   PetscOptionsBoundedInt("-petscfe_num_batches", "The number of cell batches to integrate serially", "PetscSpaceSetTileSizes", fem->numBatches, &fem->numBatches, NULL,1);
235:   if (fem->ops->setfromoptions) {
236:     (*fem->ops->setfromoptions)(PetscOptionsObject,fem);
237:   }
238:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
239:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) fem);
240:   PetscOptionsEnd();
241:   PetscFEViewFromOptions(fem, NULL, "-petscfe_view");
242:   return(0);
243: }

245: /*@C
246:   PetscFESetUp - Construct data structures for the PetscFE

248:   Collective on fem

250:   Input Parameter:
251: . fem - the PetscFE object to setup

253:   Level: developer

255: .seealso PetscFEView(), PetscFEDestroy()
256: @*/
257: PetscErrorCode PetscFESetUp(PetscFE fem)
258: {

263:   if (fem->setupcalled) return(0);
264:   fem->setupcalled = PETSC_TRUE;
265:   if (fem->ops->setup) {(*fem->ops->setup)(fem);}
266:   return(0);
267: }

269: /*@
270:   PetscFEDestroy - Destroys a PetscFE object

272:   Collective on fem

274:   Input Parameter:
275: . fem - the PetscFE object to destroy

277:   Level: developer

279: .seealso PetscFEView()
280: @*/
281: PetscErrorCode PetscFEDestroy(PetscFE *fem)
282: {

286:   if (!*fem) return(0);

289:   if (--((PetscObject)(*fem))->refct > 0) {*fem = 0; return(0);}
290:   ((PetscObject) (*fem))->refct = 0;

292:   if ((*fem)->subspaces) {
293:     PetscInt dim, d;

295:     PetscDualSpaceGetDimension((*fem)->dualSpace, &dim);
296:     for (d = 0; d < dim; ++d) {PetscFEDestroy(&(*fem)->subspaces[d]);}
297:   }
298:   PetscFree((*fem)->subspaces);
299:   PetscFree((*fem)->invV);
300:   PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->B, &(*fem)->D, NULL /*&(*fem)->H*/);
301:   PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->Bf, &(*fem)->Df, NULL /*&(*fem)->Hf*/);
302:   PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->F, NULL, NULL);
303:   PetscSpaceDestroy(&(*fem)->basisSpace);
304:   PetscDualSpaceDestroy(&(*fem)->dualSpace);
305:   PetscQuadratureDestroy(&(*fem)->quadrature);
306:   PetscQuadratureDestroy(&(*fem)->faceQuadrature);

308:   if ((*fem)->ops->destroy) {(*(*fem)->ops->destroy)(*fem);}
309:   PetscHeaderDestroy(fem);
310:   return(0);
311: }

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

316:   Collective

318:   Input Parameter:
319: . comm - The communicator for the PetscFE object

321:   Output Parameter:
322: . fem - The PetscFE object

324:   Level: beginner

326: .seealso: PetscFESetType(), PETSCFEGALERKIN
327: @*/
328: PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem)
329: {
330:   PetscFE        f;

335:   PetscCitationsRegister(FECitation,&FEcite);
336:   *fem = NULL;
337:   PetscFEInitializePackage();

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

341:   f->basisSpace    = NULL;
342:   f->dualSpace     = NULL;
343:   f->numComponents = 1;
344:   f->subspaces     = NULL;
345:   f->invV          = NULL;
346:   f->B             = NULL;
347:   f->D             = NULL;
348:   f->H             = NULL;
349:   f->Bf            = NULL;
350:   f->Df            = NULL;
351:   f->Hf            = NULL;
352:   PetscArrayzero(&f->quadrature, 1);
353:   PetscArrayzero(&f->faceQuadrature, 1);
354:   f->blockSize     = 0;
355:   f->numBlocks     = 1;
356:   f->batchSize     = 0;
357:   f->numBatches    = 1;

359:   *fem = f;
360:   return(0);
361: }

363: /*@
364:   PetscFEGetSpatialDimension - Returns the spatial dimension of the element

366:   Not collective

368:   Input Parameter:
369: . fem - The PetscFE object

371:   Output Parameter:
372: . dim - The spatial dimension

374:   Level: intermediate

376: .seealso: PetscFECreate()
377: @*/
378: PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim)
379: {
380:   DM             dm;

386:   PetscDualSpaceGetDM(fem->dualSpace, &dm);
387:   DMGetDimension(dm, dim);
388:   return(0);
389: }

391: /*@
392:   PetscFESetNumComponents - Sets the number of components in the element

394:   Not collective

396:   Input Parameters:
397: + fem - The PetscFE object
398: - comp - The number of field components

400:   Level: intermediate

402: .seealso: PetscFECreate()
403: @*/
404: PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp)
405: {
408:   fem->numComponents = comp;
409:   return(0);
410: }

412: /*@
413:   PetscFEGetNumComponents - Returns the number of components in the element

415:   Not collective

417:   Input Parameter:
418: . fem - The PetscFE object

420:   Output Parameter:
421: . comp - The number of field components

423:   Level: intermediate

425: .seealso: PetscFECreate()
426: @*/
427: PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp)
428: {
432:   *comp = fem->numComponents;
433:   return(0);
434: }

436: /*@
437:   PetscFESetTileSizes - Sets the tile sizes for evaluation

439:   Not collective

441:   Input Parameters:
442: + fem - The PetscFE object
443: . blockSize - The number of elements in a block
444: . numBlocks - The number of blocks in a batch
445: . batchSize - The number of elements in a batch
446: - numBatches - The number of batches in a chunk

448:   Level: intermediate

450: .seealso: PetscFECreate()
451: @*/
452: PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches)
453: {
456:   fem->blockSize  = blockSize;
457:   fem->numBlocks  = numBlocks;
458:   fem->batchSize  = batchSize;
459:   fem->numBatches = numBatches;
460:   return(0);
461: }

463: /*@
464:   PetscFEGetTileSizes - Returns the tile sizes for evaluation

466:   Not collective

468:   Input Parameter:
469: . fem - The PetscFE object

471:   Output Parameters:
472: + blockSize - The number of elements in a block
473: . numBlocks - The number of blocks in a batch
474: . batchSize - The number of elements in a batch
475: - numBatches - The number of batches in a chunk

477:   Level: intermediate

479: .seealso: PetscFECreate()
480: @*/
481: PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PetscInt *blockSize, PetscInt *numBlocks, PetscInt *batchSize, PetscInt *numBatches)
482: {
489:   if (blockSize)  *blockSize  = fem->blockSize;
490:   if (numBlocks)  *numBlocks  = fem->numBlocks;
491:   if (batchSize)  *batchSize  = fem->batchSize;
492:   if (numBatches) *numBatches = fem->numBatches;
493:   return(0);
494: }

496: /*@
497:   PetscFEGetBasisSpace - Returns the PetscSpace used for approximation of the solution

499:   Not collective

501:   Input Parameter:
502: . fem - The PetscFE object

504:   Output Parameter:
505: . sp - The PetscSpace object

507:   Level: intermediate

509: .seealso: PetscFECreate()
510: @*/
511: PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp)
512: {
516:   *sp = fem->basisSpace;
517:   return(0);
518: }

520: /*@
521:   PetscFESetBasisSpace - Sets the PetscSpace used for approximation of the solution

523:   Not collective

525:   Input Parameters:
526: + fem - The PetscFE object
527: - sp - The PetscSpace object

529:   Level: intermediate

531: .seealso: PetscFECreate()
532: @*/
533: PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp)
534: {

540:   PetscSpaceDestroy(&fem->basisSpace);
541:   fem->basisSpace = sp;
542:   PetscObjectReference((PetscObject) fem->basisSpace);
543:   return(0);
544: }

546: /*@
547:   PetscFEGetDualSpace - Returns the PetscDualSpace used to define the inner product

549:   Not collective

551:   Input Parameter:
552: . fem - The PetscFE object

554:   Output Parameter:
555: . sp - The PetscDualSpace object

557:   Level: intermediate

559: .seealso: PetscFECreate()
560: @*/
561: PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp)
562: {
566:   *sp = fem->dualSpace;
567:   return(0);
568: }

570: /*@
571:   PetscFESetDualSpace - Sets the PetscDualSpace used to define the inner product

573:   Not collective

575:   Input Parameters:
576: + fem - The PetscFE object
577: - sp - The PetscDualSpace object

579:   Level: intermediate

581: .seealso: PetscFECreate()
582: @*/
583: PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp)
584: {

590:   PetscDualSpaceDestroy(&fem->dualSpace);
591:   fem->dualSpace = sp;
592:   PetscObjectReference((PetscObject) fem->dualSpace);
593:   return(0);
594: }

596: /*@
597:   PetscFEGetQuadrature - Returns the PetscQuadrature used to calculate inner products

599:   Not collective

601:   Input Parameter:
602: . fem - The PetscFE object

604:   Output Parameter:
605: . q - The PetscQuadrature object

607:   Level: intermediate

609: .seealso: PetscFECreate()
610: @*/
611: PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q)
612: {
616:   *q = fem->quadrature;
617:   return(0);
618: }

620: /*@
621:   PetscFESetQuadrature - Sets the PetscQuadrature used to calculate inner products

623:   Not collective

625:   Input Parameters:
626: + fem - The PetscFE object
627: - q - The PetscQuadrature object

629:   Level: intermediate

631: .seealso: PetscFECreate()
632: @*/
633: PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q)
634: {
635:   PetscInt       Nc, qNc;

640:   PetscFEGetNumComponents(fem, &Nc);
641:   PetscQuadratureGetNumComponents(q, &qNc);
642:   if ((qNc != 1) && (Nc != qNc)) SETERRQ2(PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_SIZ, "FE components %D != Quadrature components %D and non-scalar quadrature", Nc, qNc);
643:   PetscFERestoreTabulation(fem, 0, NULL, &fem->B, &fem->D, NULL /*&(*fem)->H*/);
644:   PetscQuadratureDestroy(&fem->quadrature);
645:   fem->quadrature = q;
646:   PetscObjectReference((PetscObject) q);
647:   return(0);
648: }

650: /*@
651:   PetscFEGetFaceQuadrature - Returns the PetscQuadrature used to calculate inner products on faces

653:   Not collective

655:   Input Parameter:
656: . fem - The PetscFE object

658:   Output Parameter:
659: . q - The PetscQuadrature object

661:   Level: intermediate

663: .seealso: PetscFECreate()
664: @*/
665: PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q)
666: {
670:   *q = fem->faceQuadrature;
671:   return(0);
672: }

674: /*@
675:   PetscFESetFaceQuadrature - Sets the PetscQuadrature used to calculate inner products on faces

677:   Not collective

679:   Input Parameters:
680: + fem - The PetscFE object
681: - q - The PetscQuadrature object

683:   Level: intermediate

685: .seealso: PetscFECreate()
686: @*/
687: PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q)
688: {

693:   PetscFERestoreTabulation(fem, 0, NULL, &fem->Bf, &fem->Df, NULL /*&(*fem)->Hf*/);
694:   PetscQuadratureDestroy(&fem->faceQuadrature);
695:   fem->faceQuadrature = q;
696:   PetscObjectReference((PetscObject) q);
697:   return(0);
698: }

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

703:   Not collective

705:   Input Parameter:
706: . fem - The PetscFE object

708:   Output Parameter:
709: . numDof - Array with the number of dofs per dimension

711:   Level: intermediate

713: .seealso: PetscFECreate()
714: @*/
715: PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt **numDof)
716: {

722:   PetscDualSpaceGetNumDof(fem->dualSpace, numDof);
723:   return(0);
724: }

726: /*@C
727:   PetscFEGetDefaultTabulation - Returns the tabulation of the basis functions at the quadrature points

729:   Not collective

731:   Input Parameter:
732: . fem - The PetscFE object

734:   Output Parameters:
735: + B - The basis function values at quadrature points
736: . D - The basis function derivatives at quadrature points
737: - H - The basis function second derivatives at quadrature points

739:   Note:
740: $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
741: $ D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
742: $ 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

744:   Level: intermediate

746: .seealso: PetscFEGetTabulation(), PetscFERestoreTabulation()
747: @*/
748: PetscErrorCode PetscFEGetDefaultTabulation(PetscFE fem, PetscReal **B, PetscReal **D, PetscReal **H)
749: {
750:   PetscInt         npoints;
751:   const PetscReal *points;
752:   PetscErrorCode   ierr;

759:   PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL);
760:   if (!fem->B) {PetscFEGetTabulation(fem, npoints, points, &fem->B, &fem->D, NULL/*&fem->H*/);}
761:   if (B) *B = fem->B;
762:   if (D) *D = fem->D;
763:   if (H) *H = fem->H;
764:   return(0);
765: }

767: PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscReal **Bf, PetscReal **Df, PetscReal **Hf)
768: {
769:   PetscErrorCode   ierr;

776:   if (!fem->Bf) {
777:     const PetscReal  xi0[3] = {-1., -1., -1.};
778:     PetscReal        v0[3], J[9], detJ;
779:     PetscQuadrature  fq;
780:     PetscDualSpace   sp;
781:     DM               dm;
782:     const PetscInt  *faces;
783:     PetscInt         dim, numFaces, f, npoints, q;
784:     const PetscReal *points;
785:     PetscReal       *facePoints;

787:     PetscFEGetDualSpace(fem, &sp);
788:     PetscDualSpaceGetDM(sp, &dm);
789:     DMGetDimension(dm, &dim);
790:     DMPlexGetConeSize(dm, 0, &numFaces);
791:     DMPlexGetCone(dm, 0, &faces);
792:     PetscFEGetFaceQuadrature(fem, &fq);
793:     if (fq) {
794:       PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL);
795:       PetscMalloc1(numFaces*npoints*dim, &facePoints);
796:       for (f = 0; f < numFaces; ++f) {
797:         DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ);
798:         for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim-1, xi0, v0, J, &points[q*(dim-1)], &facePoints[(f*npoints+q)*dim]);
799:       }
800:       PetscFEGetTabulation(fem, numFaces*npoints, facePoints, &fem->Bf, &fem->Df, NULL/*&fem->Hf*/);
801:       PetscFree(facePoints);
802:     }
803:   }
804:   if (Bf) *Bf = fem->Bf;
805:   if (Df) *Df = fem->Df;
806:   if (Hf) *Hf = fem->Hf;
807:   return(0);
808: }

810: PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscReal **F)
811: {
812:   PetscErrorCode   ierr;

817:   if (!fem->F) {
818:     PetscDualSpace  sp;
819:     DM              dm;
820:     const PetscInt *cone;
821:     PetscReal      *centroids;
822:     PetscInt        dim, numFaces, f;

824:     PetscFEGetDualSpace(fem, &sp);
825:     PetscDualSpaceGetDM(sp, &dm);
826:     DMGetDimension(dm, &dim);
827:     DMPlexGetConeSize(dm, 0, &numFaces);
828:     DMPlexGetCone(dm, 0, &cone);
829:     PetscMalloc1(numFaces*dim, &centroids);
830:     for (f = 0; f < numFaces; ++f) {DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, &centroids[f*dim], NULL);}
831:     PetscFEGetTabulation(fem, numFaces, centroids, &fem->F, NULL, NULL);
832:     PetscFree(centroids);
833:   }
834:   *F = fem->F;
835:   return(0);
836: }

838: /*@C
839:   PetscFEGetTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided.

841:   Not collective

843:   Input Parameters:
844: + fem     - The PetscFE object
845: . npoints - The number of tabulation points
846: - points  - The tabulation point coordinates

848:   Output Parameters:
849: + B - The basis function values at tabulation points
850: . D - The basis function derivatives at tabulation points
851: - H - The basis function second derivatives at tabulation points

853:   Note:
854: $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
855: $ D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
856: $ 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

858:   Level: intermediate

860: .seealso: PetscFERestoreTabulation(), PetscFEGetDefaultTabulation()
861: @*/
862: PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
863: {
864:   DM               dm;
865:   PetscInt         pdim; /* Dimension of FE space P */
866:   PetscInt         dim;  /* Spatial dimension */
867:   PetscInt         comp; /* Field components */
868:   PetscErrorCode   ierr;

871:   if (!npoints) {
872:     if (B) *B = NULL;
873:     if (D) *D = NULL;
874:     if (H) *H = NULL;
875:     return(0);
876:   }
882:   PetscDualSpaceGetDM(fem->dualSpace, &dm);
883:   DMGetDimension(dm, &dim);
884:   PetscDualSpaceGetDimension(fem->dualSpace, &pdim);
885:   PetscFEGetNumComponents(fem, &comp);
886:   if (B) {DMGetWorkArray(dm, npoints*pdim*comp, MPIU_REAL, B);}
887:   if (!dim) {
888:     if (D) *D = NULL;
889:     if (H) *H = NULL;
890:   } else {
891:     if (D) {DMGetWorkArray(dm, npoints*pdim*comp*dim, MPIU_REAL, D);}
892:     if (H) {DMGetWorkArray(dm, npoints*pdim*comp*dim*dim, MPIU_REAL, H);}
893:   }
894:   (*fem->ops->gettabulation)(fem, npoints, points, B ? *B : NULL, D ? *D : NULL, H ? *H : NULL);
895:   return(0);
896: }

898: PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
899: {
900:   DM             dm;

905:   PetscDualSpaceGetDM(fem->dualSpace, &dm);
906:   if (B && *B) {DMRestoreWorkArray(dm, 0, MPIU_REAL, B);}
907:   if (D && *D) {DMRestoreWorkArray(dm, 0, MPIU_REAL, D);}
908:   if (H && *H) {DMRestoreWorkArray(dm, 0, MPIU_REAL, H);}
909:   return(0);
910: }

912: PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
913: {
914:   PetscSpace     bsp, bsubsp;
915:   PetscDualSpace dsp, dsubsp;
916:   PetscInt       dim, depth, numComp, i, j, coneSize, order;
917:   PetscFEType    type;
918:   DM             dm;
919:   DMLabel        label;
920:   PetscReal      *xi, *v, *J, detJ;
921:   const char     *name;
922:   PetscQuadrature origin, fullQuad, subQuad;

928:   PetscFEGetBasisSpace(fe,&bsp);
929:   PetscFEGetDualSpace(fe,&dsp);
930:   PetscDualSpaceGetDM(dsp,&dm);
931:   DMGetDimension(dm,&dim);
932:   DMPlexGetDepthLabel(dm,&label);
933:   DMLabelGetValue(label,refPoint,&depth);
934:   PetscCalloc1(depth,&xi);
935:   PetscMalloc1(dim,&v);
936:   PetscMalloc1(dim*dim,&J);
937:   for (i = 0; i < depth; i++) xi[i] = 0.;
938:   PetscQuadratureCreate(PETSC_COMM_SELF,&origin);
939:   PetscQuadratureSetData(origin,depth,0,1,xi,NULL);
940:   DMPlexComputeCellGeometryFEM(dm,refPoint,origin,v,J,NULL,&detJ);
941:   /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */
942:   for (i = 1; i < dim; i++) {
943:     for (j = 0; j < depth; j++) {
944:       J[i * depth + j] = J[i * dim + j];
945:     }
946:   }
947:   PetscQuadratureDestroy(&origin);
948:   PetscDualSpaceGetPointSubspace(dsp,refPoint,&dsubsp);
949:   PetscSpaceCreateSubspace(bsp,dsubsp,v,J,NULL,NULL,PETSC_OWN_POINTER,&bsubsp);
950:   PetscSpaceSetUp(bsubsp);
951:   PetscFECreate(PetscObjectComm((PetscObject)fe),trFE);
952:   PetscFEGetType(fe,&type);
953:   PetscFESetType(*trFE,type);
954:   PetscFEGetNumComponents(fe,&numComp);
955:   PetscFESetNumComponents(*trFE,numComp);
956:   PetscFESetBasisSpace(*trFE,bsubsp);
957:   PetscFESetDualSpace(*trFE,dsubsp);
958:   PetscObjectGetName((PetscObject) fe, &name);
959:   if (name) {PetscFESetName(*trFE, name);}
960:   PetscFEGetQuadrature(fe,&fullQuad);
961:   PetscQuadratureGetOrder(fullQuad,&order);
962:   DMPlexGetConeSize(dm,refPoint,&coneSize);
963:   if (coneSize == 2 * depth) {
964:     PetscDTGaussTensorQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);
965:   } else {
966:     PetscDTGaussJacobiQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);
967:   }
968:   PetscFESetQuadrature(*trFE,subQuad);
969:   PetscFESetUp(*trFE);
970:   PetscQuadratureDestroy(&subQuad);
971:   PetscSpaceDestroy(&bsubsp);
972:   return(0);
973: }

975: PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE)
976: {
977:   PetscInt       hStart, hEnd;
978:   PetscDualSpace dsp;
979:   DM             dm;

985:   *trFE = NULL;
986:   PetscFEGetDualSpace(fe,&dsp);
987:   PetscDualSpaceGetDM(dsp,&dm);
988:   DMPlexGetHeightStratum(dm,height,&hStart,&hEnd);
989:   if (hEnd <= hStart) return(0);
990:   PetscFECreatePointTrace(fe,hStart,trFE);
991:   return(0);
992: }


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

998:   Not collective

1000:   Input Parameter:
1001: . fe - The PetscFE

1003:   Output Parameter:
1004: . dim - The dimension

1006:   Level: intermediate

1008: .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
1009: @*/
1010: PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
1011: {

1017:   if (fem->ops->getdimension) {(*fem->ops->getdimension)(fem, dim);}
1018:   return(0);
1019: }

1021: /*@C
1022:   PetscFEPushforward - Map the reference element function to real space

1024:   Input Parameters:
1025: + fe     - The PetscFE
1026: . fegeom - The cell geometry
1027: . Nv     - The number of function values
1028: - vals   - The function values

1030:   Output Parameter:
1031: . vals   - The transformed function values

1033:   Level: advanced

1035:   Note: This just forwards the call onto PetscDualSpacePushforward().

1037: .seealso: PetscDualSpacePushforward()
1038: @*/
1039: PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1040: {

1044:   PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);
1045:   return(0);
1046: }

1048: /*@C
1049:   PetscFEPushforwardGradient - Map the reference element function gradient to real space

1051:   Input Parameters:
1052: + fe     - The PetscFE
1053: . fegeom - The cell geometry
1054: . Nv     - The number of function gradient values
1055: - vals   - The function gradient values

1057:   Output Parameter:
1058: . vals   - The transformed function gradient values

1060:   Level: advanced

1062:   Note: This just forwards the call onto PetscDualSpacePushforwardGradient().

1064: .seealso: PetscFEPushforward(), PetscDualSpacePushforwardGradient(), PetscDualSpacePushforward()
1065: @*/
1066: PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1067: {

1071:   PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);
1072:   return(0);
1073: }

1075: /*
1076: Purpose: Compute element vector for chunk of elements

1078: Input:
1079:   Sizes:
1080:      Ne:  number of elements
1081:      Nf:  number of fields
1082:      PetscFE
1083:        dim: spatial dimension
1084:        Nb:  number of basis functions
1085:        Nc:  number of field components
1086:        PetscQuadrature
1087:          Nq:  number of quadrature points

1089:   Geometry:
1090:      PetscFEGeom[Ne] possibly *Nq
1091:        PetscReal v0s[dim]
1092:        PetscReal n[dim]
1093:        PetscReal jacobians[dim*dim]
1094:        PetscReal jacobianInverses[dim*dim]
1095:        PetscReal jacobianDeterminants
1096:   FEM:
1097:      PetscFE
1098:        PetscQuadrature
1099:          PetscReal   quadPoints[Nq*dim]
1100:          PetscReal   quadWeights[Nq]
1101:        PetscReal   basis[Nq*Nb*Nc]
1102:        PetscReal   basisDer[Nq*Nb*Nc*dim]
1103:      PetscScalar coefficients[Ne*Nb*Nc]
1104:      PetscScalar elemVec[Ne*Nb*Nc]

1106:   Problem:
1107:      PetscInt f: the active field
1108:      f0, f1

1110:   Work Space:
1111:      PetscFE
1112:        PetscScalar f0[Nq*dim];
1113:        PetscScalar f1[Nq*dim*dim];
1114:        PetscScalar u[Nc];
1115:        PetscScalar gradU[Nc*dim];
1116:        PetscReal   x[dim];
1117:        PetscScalar realSpaceDer[dim];

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

1121: Input:
1122:   Sizes:
1123:      N_cb: Number of serial cell batches

1125:   Geometry:
1126:      PetscReal v0s[Ne*dim]
1127:      PetscReal jacobians[Ne*dim*dim]        possibly *Nq
1128:      PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq
1129:      PetscReal jacobianDeterminants[Ne]     possibly *Nq
1130:   FEM:
1131:      static PetscReal   quadPoints[Nq*dim]
1132:      static PetscReal   quadWeights[Nq]
1133:      static PetscReal   basis[Nq*Nb*Nc]
1134:      static PetscReal   basisDer[Nq*Nb*Nc*dim]
1135:      PetscScalar coefficients[Ne*Nb*Nc]
1136:      PetscScalar elemVec[Ne*Nb*Nc]

1138: ex62.c:
1139:   PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[],
1140:                                                const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[],
1141:                                                void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]),
1142:                                                void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[])

1144: ex52.c:
1145:   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)
1146:   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)

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

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

1155: ex52_integrateElementOpenCL.c:
1156: PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[],
1157:                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1158:                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)

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

1163: /*@C
1164:   PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration

1166:   Not collective

1168:   Input Parameters:
1169: + fem          - The PetscFE object for the field being integrated
1170: . prob         - The PetscDS specifying the discretizations and continuum functions
1171: . field        - The field being integrated
1172: . Ne           - The number of elements in the chunk
1173: . cgeom        - The cell geometry for each cell in the chunk
1174: . coefficients - The array of FEM basis coefficients for the elements
1175: . probAux      - The PetscDS specifying the auxiliary discretizations
1176: - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements

1178:   Output Parameter
1179: . integral     - the integral for this field

1181:   Level: developer

1183: .seealso: PetscFEIntegrateResidual()
1184: @*/
1185: PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1186:                                 const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1187: {
1188:   PetscFE        fe;

1193:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1194:   if (fe->ops->integrate) {(*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral);}
1195:   return(0);
1196: }

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

1201:   Not collective

1203:   Input Parameters:
1204: + fem          - The PetscFE object for the field being integrated
1205: . prob         - The PetscDS specifying the discretizations and continuum functions
1206: . field        - The field being integrated
1207: . obj_func     - The function to be integrated
1208: . Ne           - The number of elements in the chunk
1209: . fgeom        - The face geometry for each face in the chunk
1210: . coefficients - The array of FEM basis coefficients for the elements
1211: . probAux      - The PetscDS specifying the auxiliary discretizations
1212: - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements

1214:   Output Parameter
1215: . integral     - the integral for this field

1217:   Level: developer

1219: .seealso: PetscFEIntegrateResidual()
1220: @*/
1221: PetscErrorCode PetscFEIntegrateBd(PetscDS prob, PetscInt field,
1222:                                   void (*obj_func)(PetscInt, PetscInt, PetscInt,
1223:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1224:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1225:                                                    PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1226:                                   PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1227: {
1228:   PetscFE        fe;

1233:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1234:   if (fe->ops->integratebd) {(*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral);}
1235:   return(0);
1236: }

1238: /*@C
1239:   PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration

1241:   Not collective

1243:   Input Parameters:
1244: + fem          - The PetscFE object for the field being integrated
1245: . prob         - The PetscDS specifying the discretizations and continuum functions
1246: . field        - The field being integrated
1247: . Ne           - The number of elements in the chunk
1248: . cgeom        - The cell geometry for each cell in the chunk
1249: . coefficients - The array of FEM basis coefficients for the elements
1250: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1251: . probAux      - The PetscDS specifying the auxiliary discretizations
1252: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1253: - t            - The time

1255:   Output Parameter
1256: . elemVec      - the element residual vectors from each element

1258:   Note:
1259: $ Loop over batch of elements (e):
1260: $   Loop over quadrature points (q):
1261: $     Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q
1262: $     Call f_0 and f_1
1263: $   Loop over element vector entries (f,fc --> i):
1264: $     elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u)

1266:   Level: developer

1268: .seealso: PetscFEIntegrateResidual()
1269: @*/
1270: PetscErrorCode PetscFEIntegrateResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1271:                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1272: {
1273:   PetscFE        fe;

1278:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1279:   if (fe->ops->integrateresidual) {(*fe->ops->integrateresidual)(prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);}
1280:   return(0);
1281: }

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

1286:   Not collective

1288:   Input Parameters:
1289: + fem          - The PetscFE object for the field being integrated
1290: . prob         - The PetscDS specifying the discretizations and continuum functions
1291: . field        - The field being integrated
1292: . Ne           - The number of elements in the chunk
1293: . fgeom        - The face geometry for each cell in the chunk
1294: . coefficients - The array of FEM basis coefficients for the elements
1295: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1296: . probAux      - The PetscDS specifying the auxiliary discretizations
1297: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1298: - t            - The time

1300:   Output Parameter
1301: . elemVec      - the element residual vectors from each element

1303:   Level: developer

1305: .seealso: PetscFEIntegrateResidual()
1306: @*/
1307: PetscErrorCode PetscFEIntegrateBdResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom,
1308:                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1309: {
1310:   PetscFE        fe;

1315:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1316:   if (fe->ops->integratebdresidual) {(*fe->ops->integratebdresidual)(prob, field, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);}
1317:   return(0);
1318: }

1320: /*@C
1321:   PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration

1323:   Not collective

1325:   Input Parameters:
1326: + fem          - The PetscFE object for the field being integrated
1327: . prob         - The PetscDS specifying the discretizations and continuum functions
1328: . jtype        - The type of matrix pointwise functions that should be used
1329: . fieldI       - The test field being integrated
1330: . fieldJ       - The basis field being integrated
1331: . Ne           - The number of elements in the chunk
1332: . cgeom        - The cell geometry for each cell in the chunk
1333: . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1334: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1335: . probAux      - The PetscDS specifying the auxiliary discretizations
1336: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1337: . t            - The time
1338: - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)

1340:   Output Parameter
1341: . elemMat      - the element matrices for the Jacobian from each element

1343:   Note:
1344: $ Loop over batch of elements (e):
1345: $   Loop over element matrix entries (f,fc,g,gc --> i,j):
1346: $     Loop over quadrature points (q):
1347: $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1348: $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1349: $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1350: $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1351: $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1352:   Level: developer

1354: .seealso: PetscFEIntegrateResidual()
1355: @*/
1356: PetscErrorCode PetscFEIntegrateJacobian(PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom,
1357:                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1358: {
1359:   PetscFE        fe;

1364:   PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
1365:   if (fe->ops->integratejacobian) {(*fe->ops->integratejacobian)(prob, jtype, fieldI, fieldJ, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);}
1366:   return(0);
1367: }

1369: /*@C
1370:   PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration

1372:   Not collective

1374:   Input Parameters:
1375: . prob         - The PetscDS specifying the discretizations and continuum functions
1376: . fieldI       - The test field being integrated
1377: . fieldJ       - The basis field being integrated
1378: . Ne           - The number of elements in the chunk
1379: . fgeom        - The face geometry for each cell in the chunk
1380: . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1381: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1382: . probAux      - The PetscDS specifying the auxiliary discretizations
1383: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1384: . t            - The time
1385: - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)

1387:   Output Parameter
1388: . elemMat              - the element matrices for the Jacobian from each element

1390:   Note:
1391: $ Loop over batch of elements (e):
1392: $   Loop over element matrix entries (f,fc,g,gc --> i,j):
1393: $     Loop over quadrature points (q):
1394: $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1395: $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1396: $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1397: $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1398: $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1399:   Level: developer

1401: .seealso: PetscFEIntegrateJacobian(), PetscFEIntegrateResidual()
1402: @*/
1403: PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom,
1404:                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1405: {
1406:   PetscFE        fe;

1411:   PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
1412:   if (fe->ops->integratebdjacobian) {(*fe->ops->integratebdjacobian)(prob, fieldI, fieldJ, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);}
1413:   return(0);
1414: }

1416: PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe)
1417: {
1418:   PetscSpace      P, subP;
1419:   PetscDualSpace  Q, subQ;
1420:   PetscQuadrature subq;
1421:   PetscFEType     fetype;
1422:   PetscInt        dim, Nc;
1423:   PetscErrorCode  ierr;

1428:   if (height == 0) {
1429:     *subfe = fe;
1430:     return(0);
1431:   }
1432:   PetscFEGetBasisSpace(fe, &P);
1433:   PetscFEGetDualSpace(fe, &Q);
1434:   PetscFEGetNumComponents(fe, &Nc);
1435:   PetscFEGetFaceQuadrature(fe, &subq);
1436:   PetscDualSpaceGetDimension(Q, &dim);
1437:   if (height > dim || height < 0) {SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %D for dimension %D space", height, dim);}
1438:   if (!fe->subspaces) {PetscCalloc1(dim, &fe->subspaces);}
1439:   if (height <= dim) {
1440:     if (!fe->subspaces[height-1]) {
1441:       PetscFE     sub;
1442:       const char *name;

1444:       PetscSpaceGetHeightSubspace(P, height, &subP);
1445:       PetscDualSpaceGetHeightSubspace(Q, height, &subQ);
1446:       PetscFECreate(PetscObjectComm((PetscObject) fe), &sub);
1447:       PetscObjectGetName((PetscObject) fe,  &name);
1448:       PetscObjectSetName((PetscObject) sub,  name);
1449:       PetscFEGetType(fe, &fetype);
1450:       PetscFESetType(sub, fetype);
1451:       PetscFESetBasisSpace(sub, subP);
1452:       PetscFESetDualSpace(sub, subQ);
1453:       PetscFESetNumComponents(sub, Nc);
1454:       PetscFESetUp(sub);
1455:       PetscFESetQuadrature(sub, subq);
1456:       fe->subspaces[height-1] = sub;
1457:     }
1458:     *subfe = fe->subspaces[height-1];
1459:   } else {
1460:     *subfe = NULL;
1461:   }
1462:   return(0);
1463: }

1465: /*@
1466:   PetscFERefine - Create a "refined" PetscFE object that refines the reference cell into smaller copies. This is typically used
1467:   to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more
1468:   sparsity). It is also used to create an interpolation between regularly refined meshes.

1470:   Collective on fem

1472:   Input Parameter:
1473: . fe - The initial PetscFE

1475:   Output Parameter:
1476: . feRef - The refined PetscFE

1478:   Level: developer

1480: .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
1481: @*/
1482: PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)
1483: {
1484:   PetscSpace       P, Pref;
1485:   PetscDualSpace   Q, Qref;
1486:   DM               K, Kref;
1487:   PetscQuadrature  q, qref;
1488:   const PetscReal *v0, *jac;
1489:   PetscInt         numComp, numSubelements;
1490:   PetscErrorCode   ierr;

1493:   PetscFEGetBasisSpace(fe, &P);
1494:   PetscFEGetDualSpace(fe, &Q);
1495:   PetscFEGetQuadrature(fe, &q);
1496:   PetscDualSpaceGetDM(Q, &K);
1497:   /* Create space */
1498:   PetscObjectReference((PetscObject) P);
1499:   Pref = P;
1500:   /* Create dual space */
1501:   PetscDualSpaceDuplicate(Q, &Qref);
1502:   DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref);
1503:   PetscDualSpaceSetDM(Qref, Kref);
1504:   DMDestroy(&Kref);
1505:   PetscDualSpaceSetUp(Qref);
1506:   /* Create element */
1507:   PetscFECreate(PetscObjectComm((PetscObject) fe), feRef);
1508:   PetscFESetType(*feRef, PETSCFECOMPOSITE);
1509:   PetscFESetBasisSpace(*feRef, Pref);
1510:   PetscFESetDualSpace(*feRef, Qref);
1511:   PetscFEGetNumComponents(fe,    &numComp);
1512:   PetscFESetNumComponents(*feRef, numComp);
1513:   PetscFESetUp(*feRef);
1514:   PetscSpaceDestroy(&Pref);
1515:   PetscDualSpaceDestroy(&Qref);
1516:   /* Create quadrature */
1517:   PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL);
1518:   PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);
1519:   PetscFESetQuadrature(*feRef, qref);
1520:   PetscQuadratureDestroy(&qref);
1521:   return(0);
1522: }

1524: /*@C
1525:   PetscFECreateDefault - Create a PetscFE for basic FEM computation

1527:   Collective

1529:   Input Parameters:
1530: + comm      - The MPI comm
1531: . dim       - The spatial dimension
1532: . Nc        - The number of components
1533: . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
1534: . prefix    - The options prefix, or NULL
1535: - qorder    - The quadrature order

1537:   Output Parameter:
1538: . fem - The PetscFE object

1540:   Level: beginner

1542: .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
1543: @*/
1544: PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem)
1545: {
1546:   PetscQuadrature q, fq;
1547:   DM              K;
1548:   PetscSpace      P;
1549:   PetscDualSpace  Q;
1550:   PetscInt        order, quadPointsPerEdge;
1551:   PetscBool       tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
1552:   PetscErrorCode  ierr;

1555:   /* Create space */
1556:   PetscSpaceCreate(comm, &P);
1557:   PetscObjectSetOptionsPrefix((PetscObject) P, prefix);
1558:   PetscSpacePolynomialSetTensor(P, tensor);
1559:   PetscSpaceSetNumComponents(P, Nc);
1560:   PetscSpaceSetNumVariables(P, dim);
1561:   PetscSpaceSetFromOptions(P);
1562:   PetscSpaceSetUp(P);
1563:   PetscSpaceGetDegree(P, &order, NULL);
1564:   PetscSpacePolynomialGetTensor(P, &tensor);
1565:   /* Create dual space */
1566:   PetscDualSpaceCreate(comm, &Q);
1567:   PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);
1568:   PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);
1569:   PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);
1570:   PetscDualSpaceSetDM(Q, K);
1571:   DMDestroy(&K);
1572:   PetscDualSpaceSetNumComponents(Q, Nc);
1573:   PetscDualSpaceSetOrder(Q, order);
1574:   PetscDualSpaceLagrangeSetTensor(Q, tensor);
1575:   PetscDualSpaceSetFromOptions(Q);
1576:   PetscDualSpaceSetUp(Q);
1577:   /* Create element */
1578:   PetscFECreate(comm, fem);
1579:   PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);
1580:   PetscFESetBasisSpace(*fem, P);
1581:   PetscFESetDualSpace(*fem, Q);
1582:   PetscFESetNumComponents(*fem, Nc);
1583:   PetscFESetFromOptions(*fem);
1584:   PetscFESetUp(*fem);
1585:   PetscSpaceDestroy(&P);
1586:   PetscDualSpaceDestroy(&Q);
1587:   /* Create quadrature (with specified order if given) */
1588:   qorder = qorder >= 0 ? qorder : order;
1589:   PetscObjectOptionsBegin((PetscObject)*fem);
1590:   PetscOptionsBoundedInt("-petscfe_default_quadrature_order","Quadrature order is one less than quadrature points per edge","PetscFECreateDefault",qorder,&qorder,NULL,0);
1591:   PetscOptionsEnd();
1592:   quadPointsPerEdge = PetscMax(qorder + 1,1);
1593:   if (isSimplex) {
1594:     PetscDTGaussJacobiQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);
1595:     PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);
1596:   } else {
1597:     PetscDTGaussTensorQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);
1598:     PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);
1599:   }
1600:   PetscFESetQuadrature(*fem, q);
1601:   PetscFESetFaceQuadrature(*fem, fq);
1602:   PetscQuadratureDestroy(&q);
1603:   PetscQuadratureDestroy(&fq);
1604:   return(0);
1605: }

1607: /*@C
1608:   PetscFESetName - Names the FE and its subobjects

1610:   Not collective

1612:   Input Parameters:
1613: + fe   - The PetscFE
1614: - name - The name

1616:   Level: beginner

1618: .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
1619: @*/
1620: PetscErrorCode PetscFESetName(PetscFE fe, const char name[])
1621: {
1622:   PetscSpace     P;
1623:   PetscDualSpace Q;

1627:   PetscFEGetBasisSpace(fe, &P);
1628:   PetscFEGetDualSpace(fe, &Q);
1629:   PetscObjectSetName((PetscObject) fe, name);
1630:   PetscObjectSetName((PetscObject) P,  name);
1631:   PetscObjectSetName((PetscObject) Q,  name);
1632:   return(0);
1633: }

1635: PetscErrorCode PetscFEEvaluateFieldJets_Internal(PetscDS ds, PetscInt dim, PetscInt Nf, const PetscInt Nb[], const PetscInt Nc[], PetscInt q, PetscReal *basisField[], PetscReal *basisFieldDer[], PetscFEGeom *fegeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscScalar u[], PetscScalar u_x[], PetscScalar u_t[])
1636: {
1637:   PetscInt       dOffset = 0, fOffset = 0, f;

1640:   for (f = 0; f < Nf; ++f) {
1641:     PetscFE          fe;
1642:     const PetscInt   Nbf = Nb[f], Ncf = Nc[f];
1643:     const PetscReal *Bq = &basisField[f][q*Nbf*Ncf];
1644:     const PetscReal *Dq = &basisFieldDer[f][q*Nbf*Ncf*dim];
1645:     PetscInt         b, c, d;

1647:     PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);
1648:     for (c = 0; c < Ncf; ++c)     u[fOffset+c] = 0.0;
1649:     for (d = 0; d < dim*Ncf; ++d) u_x[fOffset*dim+d] = 0.0;
1650:     for (b = 0; b < Nbf; ++b) {
1651:       for (c = 0; c < Ncf; ++c) {
1652:         const PetscInt cidx = b*Ncf+c;

1654:         u[fOffset+c] += Bq[cidx]*coefficients[dOffset+b];
1655:         for (d = 0; d < dim; ++d) u_x[(fOffset+c)*dim+d] += Dq[cidx*dim+d]*coefficients[dOffset+b];
1656:       }
1657:     }
1658:     PetscFEPushforward(fe, fegeom, 1, &u[fOffset]);
1659:     PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset*dim]);
1660:     if (u_t) {
1661:       for (c = 0; c < Ncf; ++c) u_t[fOffset+c] = 0.0;
1662:       for (b = 0; b < Nbf; ++b) {
1663:         for (c = 0; c < Ncf; ++c) {
1664:           const PetscInt cidx = b*Ncf+c;

1666:           u_t[fOffset+c] += Bq[cidx]*coefficients_t[dOffset+b];
1667:         }
1668:       }
1669:       PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]);
1670:     }
1671:     fOffset += Ncf;
1672:     dOffset += Nbf;
1673:   }
1674:   return 0;
1675: }

1677: PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[])
1678: {
1679:   PetscFE        fe;
1680:   PetscReal     *faceBasis;
1681:   PetscInt       Nb, Nc, b, c;

1684:   if (!prob) return 0;
1685:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1686:   PetscFEGetDimension(fe, &Nb);
1687:   PetscFEGetNumComponents(fe, &Nc);
1688:   PetscFEGetFaceCentroidTabulation(fe, &faceBasis);
1689:   for (c = 0; c < Nc; ++c) {u[c] = 0.0;}
1690:   for (b = 0; b < Nb; ++b) {
1691:     for (c = 0; c < Nc; ++c) {
1692:       const PetscInt cidx = b*Nc+c;

1694:       u[c] += coefficients[cidx]*faceBasis[faceLoc*Nb*Nc+cidx];
1695:     }
1696:   }
1697:   return 0;
1698: }

1700: PetscErrorCode PetscFEUpdateElementVec_Internal(PetscFE fe, PetscInt dim, PetscInt Nq, PetscInt Nb, PetscInt Nc, PetscReal basis[], PetscReal basisDer[], PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[])
1701: {
1702:   PetscInt       q, b, c, d;

1705:   for (b = 0; b < Nb; ++b) elemVec[b] = 0.0;
1706:   for (q = 0; q < Nq; ++q) {
1707:     for (b = 0; b < Nb; ++b) {
1708:       for (c = 0; c < Nc; ++c) {
1709:         const PetscInt bcidx = b*Nc+c;

1711:         tmpBasis[bcidx] = basis[q*Nb*Nc+bcidx];
1712:         for (d = 0; d < dim; ++d) tmpBasisDer[bcidx*dim+d] = basisDer[q*Nb*Nc*dim+bcidx*dim+d];
1713:       }
1714:     }
1715:     PetscFEPushforward(fe, fegeom, Nb, tmpBasis);
1716:     PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer);
1717:     for (b = 0; b < Nb; ++b) {
1718:       for (c = 0; c < Nc; ++c) {
1719:         const PetscInt bcidx = b*Nc+c;
1720:         const PetscInt qcidx = q*Nc+c;

1722:         elemVec[b] += tmpBasis[bcidx]*f0[qcidx];
1723:         for (d = 0; d < dim; ++d) elemVec[b] += tmpBasisDer[bcidx*dim+d]*f1[qcidx*dim+d];
1724:       }
1725:     }
1726:   }
1727:   return(0);
1728: }

1730: PetscErrorCode PetscFEUpdateElementMat_Internal(PetscFE feI, PetscFE feJ, PetscInt dim, PetscInt NbI, PetscInt NcI, const PetscReal basisI[], const PetscReal basisDerI[], PetscScalar tmpBasisI[], PetscScalar tmpBasisDerI[], PetscInt NbJ, PetscInt NcJ, const PetscReal basisJ[], const PetscReal basisDerJ[], 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[])
1731: {
1732:   PetscInt       f, fc, g, gc, df, dg;

1735:   for (f = 0; f < NbI; ++f) {
1736:     for (fc = 0; fc < NcI; ++fc) {
1737:       const PetscInt fidx = f*NcI+fc; /* Test function basis index */

1739:       tmpBasisI[fidx] = basisI[fidx];
1740:       for (df = 0; df < dim; ++df) tmpBasisDerI[fidx*dim+df] = basisDerI[fidx*dim+df];
1741:     }
1742:   }
1743:   PetscFEPushforward(feI, fegeom, NbI, tmpBasisI);
1744:   PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI);
1745:   for (g = 0; g < NbJ; ++g) {
1746:     for (gc = 0; gc < NcJ; ++gc) {
1747:       const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */

1749:       tmpBasisJ[gidx] = basisJ[gidx];
1750:       for (dg = 0; dg < dim; ++dg) tmpBasisDerJ[gidx*dim+dg] = basisDerJ[gidx*dim+dg];
1751:     }
1752:   }
1753:   PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ);
1754:   PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ);
1755:   for (f = 0; f < NbI; ++f) {
1756:     for (fc = 0; fc < NcI; ++fc) {
1757:       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
1758:       const PetscInt i    = offsetI+f; /* Element matrix row */
1759:       for (g = 0; g < NbJ; ++g) {
1760:         for (gc = 0; gc < NcJ; ++gc) {
1761:           const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
1762:           const PetscInt j    = offsetJ+g; /* Element matrix column */
1763:           const PetscInt fOff = eOffset+i*totDim+j;

1765:           elemMat[fOff] += tmpBasisI[fidx]*g0[fc*NcJ+gc]*tmpBasisJ[gidx];
1766:           for (df = 0; df < dim; ++df) {
1767:             elemMat[fOff] += tmpBasisI[fidx]*g1[(fc*NcJ+gc)*dim+df]*tmpBasisDerJ[gidx*dim+df];
1768:             elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g2[(fc*NcJ+gc)*dim+df]*tmpBasisJ[gidx];
1769:             for (dg = 0; dg < dim; ++dg) {
1770:               elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g3[((fc*NcJ+gc)*dim+df)*dim+dg]*tmpBasisDerJ[gidx*dim+dg];
1771:             }
1772:           }
1773:         }
1774:       }
1775:     }
1776:   }
1777:   return(0);
1778: }

1780: PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom)
1781: {
1782:   PetscDualSpace  dsp;
1783:   DM              dm;
1784:   PetscQuadrature quadDef;
1785:   PetscInt        dim, cdim, Nq;
1786:   PetscErrorCode  ierr;

1789:   PetscFEGetDualSpace(fe, &dsp);
1790:   PetscDualSpaceGetDM(dsp, &dm);
1791:   DMGetDimension(dm, &dim);
1792:   DMGetCoordinateDim(dm, &cdim);
1793:   PetscFEGetQuadrature(fe, &quadDef);
1794:   quad = quad ? quad : quadDef;
1795:   PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);
1796:   PetscMalloc1(Nq*cdim,      &cgeom->v);
1797:   PetscMalloc1(Nq*cdim*cdim, &cgeom->J);
1798:   PetscMalloc1(Nq*cdim*cdim, &cgeom->invJ);
1799:   PetscMalloc1(Nq,           &cgeom->detJ);
1800:   cgeom->dim       = dim;
1801:   cgeom->dimEmbed  = cdim;
1802:   cgeom->numCells  = 1;
1803:   cgeom->numPoints = Nq;
1804:   DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ);
1805:   return(0);
1806: }

1808: PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom)
1809: {

1813:   PetscFree(cgeom->v);
1814:   PetscFree(cgeom->J);
1815:   PetscFree(cgeom->invJ);
1816:   PetscFree(cgeom->detJ);
1817:   return(0);
1818: }