Actual source code: fe.c

petsc-master 2019-09-15
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: beginner

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: intermediate

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: intermediate

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: beginner

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: /*@C
768:   PetscFEGetFaceTabulation - Returns the tabulation of the basis functions at the face quadrature points

770:   Not collective

772:   Input Parameter:
773: . fem - The PetscFE object

775:   Output Parameters:
776: + B - The basis function values at face quadrature points
777: . D - The basis function derivatives at face quadrature points
778: - H - The basis function second derivatives at face quadrature points

780:   Note:
781: $ Bf[((f*Nq + q)*pdim + i)*Nc + c] is the value at point f,q for basis function i and component c
782: $ 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
783: $ 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

785:   Level: intermediate

787: .seealso: PetscFEGetDefaultTabulation(), PetscFEGetTabulation(), PetscFERestoreTabulation()
788: @*/
789: PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscReal **Bf, PetscReal **Df, PetscReal **Hf)
790: {
791:   PetscErrorCode   ierr;

798:   if (!fem->Bf) {
799:     const PetscReal  xi0[3] = {-1., -1., -1.};
800:     PetscReal        v0[3], J[9], detJ;
801:     PetscQuadrature  fq;
802:     PetscDualSpace   sp;
803:     DM               dm;
804:     const PetscInt  *faces;
805:     PetscInt         dim, numFaces, f, npoints, q;
806:     const PetscReal *points;
807:     PetscReal       *facePoints;

809:     PetscFEGetDualSpace(fem, &sp);
810:     PetscDualSpaceGetDM(sp, &dm);
811:     DMGetDimension(dm, &dim);
812:     DMPlexGetConeSize(dm, 0, &numFaces);
813:     DMPlexGetCone(dm, 0, &faces);
814:     PetscFEGetFaceQuadrature(fem, &fq);
815:     if (fq) {
816:       PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL);
817:       PetscMalloc1(numFaces*npoints*dim, &facePoints);
818:       for (f = 0; f < numFaces; ++f) {
819:         DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ);
820:         for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim-1, xi0, v0, J, &points[q*(dim-1)], &facePoints[(f*npoints+q)*dim]);
821:       }
822:       PetscFEGetTabulation(fem, numFaces*npoints, facePoints, &fem->Bf, &fem->Df, NULL/*&fem->Hf*/);
823:       PetscFree(facePoints);
824:     }
825:   }
826:   if (Bf) *Bf = fem->Bf;
827:   if (Df) *Df = fem->Df;
828:   if (Hf) *Hf = fem->Hf;
829:   return(0);
830: }

832: /*@C
833:   PetscFEGetFaceTabulation - Returns the tabulation of the basis functions at the face centroid points

835:   Not collective

837:   Input Parameter:
838: . fem - The PetscFE object

840:   Output Parameters:
841: + B - The basis function values at face centroid points
842: . D - The basis function derivatives at face centroid points
843: - H - The basis function second derivatives at face centroid points

845:   Note:
846: $ Bf[(f*pdim + i)*Nc + c] is the value at point f for basis function i and component c
847: $ Df[((f*pdim + i)*Nc + c)*dim + d] is the derivative value at point f for basis function i, component c, in direction d
848: $ Hf[(((f*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point f for basis function i, component c, in directions d and e

850:   Level: intermediate

852: .seealso: PetscFEGetFaceTabulation(), PetscFEGetDefaultTabulation(), PetscFEGetTabulation(), PetscFERestoreTabulation()
853: @*/
854: PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscReal **F)
855: {
856:   PetscErrorCode   ierr;

861:   if (!fem->F) {
862:     PetscDualSpace  sp;
863:     DM              dm;
864:     const PetscInt *cone;
865:     PetscReal      *centroids;
866:     PetscInt        dim, numFaces, f;

868:     PetscFEGetDualSpace(fem, &sp);
869:     PetscDualSpaceGetDM(sp, &dm);
870:     DMGetDimension(dm, &dim);
871:     DMPlexGetConeSize(dm, 0, &numFaces);
872:     DMPlexGetCone(dm, 0, &cone);
873:     PetscMalloc1(numFaces*dim, &centroids);
874:     for (f = 0; f < numFaces; ++f) {DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, &centroids[f*dim], NULL);}
875:     PetscFEGetTabulation(fem, numFaces, centroids, &fem->F, NULL, NULL);
876:     PetscFree(centroids);
877:   }
878:   *F = fem->F;
879:   return(0);
880: }

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

885:   Not collective

887:   Input Parameters:
888: + fem     - The PetscFE object
889: . npoints - The number of tabulation points
890: - points  - The tabulation point coordinates

892:   Output Parameters:
893: + B - The basis function values at tabulation points
894: . D - The basis function derivatives at tabulation points
895: - H - The basis function second derivatives at tabulation points

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

902:   Level: intermediate

904: .seealso: PetscFERestoreTabulation(), PetscFEGetDefaultTabulation()
905: @*/
906: PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
907: {
908:   DM               dm;
909:   PetscInt         pdim; /* Dimension of FE space P */
910:   PetscInt         dim;  /* Spatial dimension */
911:   PetscInt         comp; /* Field components */
912:   PetscErrorCode   ierr;

915:   if (!npoints) {
916:     if (B) *B = NULL;
917:     if (D) *D = NULL;
918:     if (H) *H = NULL;
919:     return(0);
920:   }
926:   PetscDualSpaceGetDM(fem->dualSpace, &dm);
927:   DMGetDimension(dm, &dim);
928:   PetscDualSpaceGetDimension(fem->dualSpace, &pdim);
929:   PetscFEGetNumComponents(fem, &comp);
930:   if (B) {DMGetWorkArray(dm, npoints*pdim*comp, MPIU_REAL, B);}
931:   if (!dim) {
932:     if (D) *D = NULL;
933:     if (H) *H = NULL;
934:   } else {
935:     if (D) {DMGetWorkArray(dm, npoints*pdim*comp*dim, MPIU_REAL, D);}
936:     if (H) {DMGetWorkArray(dm, npoints*pdim*comp*dim*dim, MPIU_REAL, H);}
937:   }
938:   (*fem->ops->gettabulation)(fem, npoints, points, B ? *B : NULL, D ? *D : NULL, H ? *H : NULL);
939:   return(0);
940: }

942: /*@C
943:   PetscFERestoreTabulation - Frees memory from the associated tabulation.

945:   Not collective

947:   Input Parameters:
948: + fem     - The PetscFE object
949: . npoints - The number of tabulation points
950: . points  - The tabulation point coordinates
951: . B - The basis function values at tabulation points
952: . D - The basis function derivatives at tabulation points
953: - H - The basis function second derivatives at tabulation points

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

960:   Level: intermediate

962: .seealso: PetscFEGetTabulation(), PetscFEGetDefaultTabulation()
963: @*/
964: PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
965: {
966:   DM             dm;

971:   PetscDualSpaceGetDM(fem->dualSpace, &dm);
972:   if (B && *B) {DMRestoreWorkArray(dm, 0, MPIU_REAL, B);}
973:   if (D && *D) {DMRestoreWorkArray(dm, 0, MPIU_REAL, D);}
974:   if (H && *H) {DMRestoreWorkArray(dm, 0, MPIU_REAL, H);}
975:   return(0);
976: }

978: PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
979: {
980:   PetscSpace     bsp, bsubsp;
981:   PetscDualSpace dsp, dsubsp;
982:   PetscInt       dim, depth, numComp, i, j, coneSize, order;
983:   PetscFEType    type;
984:   DM             dm;
985:   DMLabel        label;
986:   PetscReal      *xi, *v, *J, detJ;
987:   const char     *name;
988:   PetscQuadrature origin, fullQuad, subQuad;

994:   PetscFEGetBasisSpace(fe,&bsp);
995:   PetscFEGetDualSpace(fe,&dsp);
996:   PetscDualSpaceGetDM(dsp,&dm);
997:   DMGetDimension(dm,&dim);
998:   DMPlexGetDepthLabel(dm,&label);
999:   DMLabelGetValue(label,refPoint,&depth);
1000:   PetscCalloc1(depth,&xi);
1001:   PetscMalloc1(dim,&v);
1002:   PetscMalloc1(dim*dim,&J);
1003:   for (i = 0; i < depth; i++) xi[i] = 0.;
1004:   PetscQuadratureCreate(PETSC_COMM_SELF,&origin);
1005:   PetscQuadratureSetData(origin,depth,0,1,xi,NULL);
1006:   DMPlexComputeCellGeometryFEM(dm,refPoint,origin,v,J,NULL,&detJ);
1007:   /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */
1008:   for (i = 1; i < dim; i++) {
1009:     for (j = 0; j < depth; j++) {
1010:       J[i * depth + j] = J[i * dim + j];
1011:     }
1012:   }
1013:   PetscQuadratureDestroy(&origin);
1014:   PetscDualSpaceGetPointSubspace(dsp,refPoint,&dsubsp);
1015:   PetscSpaceCreateSubspace(bsp,dsubsp,v,J,NULL,NULL,PETSC_OWN_POINTER,&bsubsp);
1016:   PetscSpaceSetUp(bsubsp);
1017:   PetscFECreate(PetscObjectComm((PetscObject)fe),trFE);
1018:   PetscFEGetType(fe,&type);
1019:   PetscFESetType(*trFE,type);
1020:   PetscFEGetNumComponents(fe,&numComp);
1021:   PetscFESetNumComponents(*trFE,numComp);
1022:   PetscFESetBasisSpace(*trFE,bsubsp);
1023:   PetscFESetDualSpace(*trFE,dsubsp);
1024:   PetscObjectGetName((PetscObject) fe, &name);
1025:   if (name) {PetscFESetName(*trFE, name);}
1026:   PetscFEGetQuadrature(fe,&fullQuad);
1027:   PetscQuadratureGetOrder(fullQuad,&order);
1028:   DMPlexGetConeSize(dm,refPoint,&coneSize);
1029:   if (coneSize == 2 * depth) {
1030:     PetscDTGaussTensorQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);
1031:   } else {
1032:     PetscDTGaussJacobiQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);
1033:   }
1034:   PetscFESetQuadrature(*trFE,subQuad);
1035:   PetscFESetUp(*trFE);
1036:   PetscQuadratureDestroy(&subQuad);
1037:   PetscSpaceDestroy(&bsubsp);
1038:   return(0);
1039: }

1041: PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE)
1042: {
1043:   PetscInt       hStart, hEnd;
1044:   PetscDualSpace dsp;
1045:   DM             dm;

1051:   *trFE = NULL;
1052:   PetscFEGetDualSpace(fe,&dsp);
1053:   PetscDualSpaceGetDM(dsp,&dm);
1054:   DMPlexGetHeightStratum(dm,height,&hStart,&hEnd);
1055:   if (hEnd <= hStart) return(0);
1056:   PetscFECreatePointTrace(fe,hStart,trFE);
1057:   return(0);
1058: }


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

1064:   Not collective

1066:   Input Parameter:
1067: . fe - The PetscFE

1069:   Output Parameter:
1070: . dim - The dimension

1072:   Level: intermediate

1074: .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
1075: @*/
1076: PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
1077: {

1083:   if (fem->ops->getdimension) {(*fem->ops->getdimension)(fem, dim);}
1084:   return(0);
1085: }

1087: /*@C
1088:   PetscFEPushforward - Map the reference element function to real space

1090:   Input Parameters:
1091: + fe     - The PetscFE
1092: . fegeom - The cell geometry
1093: . Nv     - The number of function values
1094: - vals   - The function values

1096:   Output Parameter:
1097: . vals   - The transformed function values

1099:   Level: advanced

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

1103: .seealso: PetscDualSpacePushforward()
1104: @*/
1105: PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1106: {

1110:   PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);
1111:   return(0);
1112: }

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

1117:   Input Parameters:
1118: + fe     - The PetscFE
1119: . fegeom - The cell geometry
1120: . Nv     - The number of function gradient values
1121: - vals   - The function gradient values

1123:   Output Parameter:
1124: . vals   - The transformed function gradient values

1126:   Level: advanced

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

1130: .seealso: PetscFEPushforward(), PetscDualSpacePushforwardGradient(), PetscDualSpacePushforward()
1131: @*/
1132: PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1133: {

1137:   PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);
1138:   return(0);
1139: }

1141: /*
1142: Purpose: Compute element vector for chunk of elements

1144: Input:
1145:   Sizes:
1146:      Ne:  number of elements
1147:      Nf:  number of fields
1148:      PetscFE
1149:        dim: spatial dimension
1150:        Nb:  number of basis functions
1151:        Nc:  number of field components
1152:        PetscQuadrature
1153:          Nq:  number of quadrature points

1155:   Geometry:
1156:      PetscFEGeom[Ne] possibly *Nq
1157:        PetscReal v0s[dim]
1158:        PetscReal n[dim]
1159:        PetscReal jacobians[dim*dim]
1160:        PetscReal jacobianInverses[dim*dim]
1161:        PetscReal jacobianDeterminants
1162:   FEM:
1163:      PetscFE
1164:        PetscQuadrature
1165:          PetscReal   quadPoints[Nq*dim]
1166:          PetscReal   quadWeights[Nq]
1167:        PetscReal   basis[Nq*Nb*Nc]
1168:        PetscReal   basisDer[Nq*Nb*Nc*dim]
1169:      PetscScalar coefficients[Ne*Nb*Nc]
1170:      PetscScalar elemVec[Ne*Nb*Nc]

1172:   Problem:
1173:      PetscInt f: the active field
1174:      f0, f1

1176:   Work Space:
1177:      PetscFE
1178:        PetscScalar f0[Nq*dim];
1179:        PetscScalar f1[Nq*dim*dim];
1180:        PetscScalar u[Nc];
1181:        PetscScalar gradU[Nc*dim];
1182:        PetscReal   x[dim];
1183:        PetscScalar realSpaceDer[dim];

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

1187: Input:
1188:   Sizes:
1189:      N_cb: Number of serial cell batches

1191:   Geometry:
1192:      PetscReal v0s[Ne*dim]
1193:      PetscReal jacobians[Ne*dim*dim]        possibly *Nq
1194:      PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq
1195:      PetscReal jacobianDeterminants[Ne]     possibly *Nq
1196:   FEM:
1197:      static PetscReal   quadPoints[Nq*dim]
1198:      static PetscReal   quadWeights[Nq]
1199:      static PetscReal   basis[Nq*Nb*Nc]
1200:      static PetscReal   basisDer[Nq*Nb*Nc*dim]
1201:      PetscScalar coefficients[Ne*Nb*Nc]
1202:      PetscScalar elemVec[Ne*Nb*Nc]

1204: ex62.c:
1205:   PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[],
1206:                                                const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[],
1207:                                                void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]),
1208:                                                void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[])

1210: ex52.c:
1211:   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)
1212:   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)

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

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

1221: ex52_integrateElementOpenCL.c:
1222: PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[],
1223:                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1224:                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)

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

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

1232:   Not collective

1234:   Input Parameters:
1235: + fem          - The PetscFE object for the field being integrated
1236: . prob         - The PetscDS specifying the discretizations and continuum functions
1237: . field        - The field being integrated
1238: . Ne           - The number of elements in the chunk
1239: . cgeom        - The cell geometry for each cell in the chunk
1240: . coefficients - The array of FEM basis coefficients for the elements
1241: . probAux      - The PetscDS specifying the auxiliary discretizations
1242: - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements

1244:   Output Parameter
1245: . integral     - the integral for this field

1247:   Level: intermediate

1249: .seealso: PetscFEIntegrateResidual()
1250: @*/
1251: PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1252:                                 const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1253: {
1254:   PetscFE        fe;

1259:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1260:   if (fe->ops->integrate) {(*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral);}
1261:   return(0);
1262: }

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

1267:   Not collective

1269:   Input Parameters:
1270: + fem          - The PetscFE object for the field being integrated
1271: . prob         - The PetscDS specifying the discretizations and continuum functions
1272: . field        - The field being integrated
1273: . obj_func     - The function to be integrated
1274: . Ne           - The number of elements in the chunk
1275: . fgeom        - The face geometry for each face in the chunk
1276: . coefficients - The array of FEM basis coefficients for the elements
1277: . probAux      - The PetscDS specifying the auxiliary discretizations
1278: - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements

1280:   Output Parameter
1281: . integral     - the integral for this field

1283:   Level: intermediate

1285: .seealso: PetscFEIntegrateResidual()
1286: @*/
1287: PetscErrorCode PetscFEIntegrateBd(PetscDS prob, PetscInt field,
1288:                                   void (*obj_func)(PetscInt, PetscInt, PetscInt,
1289:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1290:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1291:                                                    PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1292:                                   PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1293: {
1294:   PetscFE        fe;

1299:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1300:   if (fe->ops->integratebd) {(*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral);}
1301:   return(0);
1302: }

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

1307:   Not collective

1309:   Input Parameters:
1310: + fem          - The PetscFE object for the field being integrated
1311: . prob         - The PetscDS specifying the discretizations and continuum functions
1312: . field        - The field being integrated
1313: . Ne           - The number of elements in the chunk
1314: . cgeom        - The cell geometry for each cell in the chunk
1315: . coefficients - The array of FEM basis coefficients for the elements
1316: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1317: . probAux      - The PetscDS specifying the auxiliary discretizations
1318: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1319: - t            - The time

1321:   Output Parameter
1322: . elemVec      - the element residual vectors from each element

1324:   Note:
1325: $ Loop over batch of elements (e):
1326: $   Loop over quadrature points (q):
1327: $     Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q
1328: $     Call f_0 and f_1
1329: $   Loop over element vector entries (f,fc --> i):
1330: $     elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u)

1332:   Level: intermediate

1334: .seealso: PetscFEIntegrateResidual()
1335: @*/
1336: PetscErrorCode PetscFEIntegrateResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1337:                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1338: {
1339:   PetscFE        fe;

1344:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1345:   if (fe->ops->integrateresidual) {(*fe->ops->integrateresidual)(prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);}
1346:   return(0);
1347: }

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

1352:   Not collective

1354:   Input Parameters:
1355: + fem          - The PetscFE object for the field being integrated
1356: . prob         - The PetscDS specifying the discretizations and continuum functions
1357: . field        - The field being integrated
1358: . Ne           - The number of elements in the chunk
1359: . fgeom        - The face geometry for each cell in the chunk
1360: . coefficients - The array of FEM basis coefficients for the elements
1361: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1362: . probAux      - The PetscDS specifying the auxiliary discretizations
1363: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1364: - t            - The time

1366:   Output Parameter
1367: . elemVec      - the element residual vectors from each element

1369:   Level: intermediate

1371: .seealso: PetscFEIntegrateResidual()
1372: @*/
1373: PetscErrorCode PetscFEIntegrateBdResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom,
1374:                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1375: {
1376:   PetscFE        fe;

1381:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1382:   if (fe->ops->integratebdresidual) {(*fe->ops->integratebdresidual)(prob, field, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);}
1383:   return(0);
1384: }

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

1389:   Not collective

1391:   Input Parameters:
1392: + fem          - The PetscFE object for the field being integrated
1393: . prob         - The PetscDS specifying the discretizations and continuum functions
1394: . jtype        - The type of matrix pointwise functions that should be used
1395: . fieldI       - The test field being integrated
1396: . fieldJ       - The basis field being integrated
1397: . Ne           - The number of elements in the chunk
1398: . cgeom        - The cell geometry for each cell in the chunk
1399: . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1400: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1401: . probAux      - The PetscDS specifying the auxiliary discretizations
1402: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1403: . t            - The time
1404: - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)

1406:   Output Parameter
1407: . elemMat      - the element matrices for the Jacobian from each element

1409:   Note:
1410: $ Loop over batch of elements (e):
1411: $   Loop over element matrix entries (f,fc,g,gc --> i,j):
1412: $     Loop over quadrature points (q):
1413: $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1414: $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1415: $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1416: $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1417: $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1418:   Level: intermediate

1420: .seealso: PetscFEIntegrateResidual()
1421: @*/
1422: PetscErrorCode PetscFEIntegrateJacobian(PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom,
1423:                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1424: {
1425:   PetscFE        fe;

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

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

1438:   Not collective

1440:   Input Parameters:
1441: . prob         - The PetscDS specifying the discretizations and continuum functions
1442: . fieldI       - The test field being integrated
1443: . fieldJ       - The basis field being integrated
1444: . Ne           - The number of elements in the chunk
1445: . fgeom        - The face geometry for each cell in the chunk
1446: . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1447: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1448: . probAux      - The PetscDS specifying the auxiliary discretizations
1449: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1450: . t            - The time
1451: - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)

1453:   Output Parameter
1454: . elemMat              - the element matrices for the Jacobian from each element

1456:   Note:
1457: $ Loop over batch of elements (e):
1458: $   Loop over element matrix entries (f,fc,g,gc --> i,j):
1459: $     Loop over quadrature points (q):
1460: $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1461: $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1462: $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1463: $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1464: $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1465:   Level: intermediate

1467: .seealso: PetscFEIntegrateJacobian(), PetscFEIntegrateResidual()
1468: @*/
1469: PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom,
1470:                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1471: {
1472:   PetscFE        fe;

1477:   PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
1478:   if (fe->ops->integratebdjacobian) {(*fe->ops->integratebdjacobian)(prob, fieldI, fieldJ, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);}
1479:   return(0);
1480: }

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

1485:   Input Parameters:
1486: + fe     - The finite element space
1487: - height - The height of the Plex point

1489:   Output Parameter:
1490: . subfe  - The subspace of this FE space

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

1494:   Level: advanced

1496: .seealso: PetscFECreateDefault()
1497: @*/
1498: PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe)
1499: {
1500:   PetscSpace      P, subP;
1501:   PetscDualSpace  Q, subQ;
1502:   PetscQuadrature subq;
1503:   PetscFEType     fetype;
1504:   PetscInt        dim, Nc;
1505:   PetscErrorCode  ierr;

1510:   if (height == 0) {
1511:     *subfe = fe;
1512:     return(0);
1513:   }
1514:   PetscFEGetBasisSpace(fe, &P);
1515:   PetscFEGetDualSpace(fe, &Q);
1516:   PetscFEGetNumComponents(fe, &Nc);
1517:   PetscFEGetFaceQuadrature(fe, &subq);
1518:   PetscDualSpaceGetDimension(Q, &dim);
1519:   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);}
1520:   if (!fe->subspaces) {PetscCalloc1(dim, &fe->subspaces);}
1521:   if (height <= dim) {
1522:     if (!fe->subspaces[height-1]) {
1523:       PetscFE     sub;
1524:       const char *name;

1526:       PetscSpaceGetHeightSubspace(P, height, &subP);
1527:       PetscDualSpaceGetHeightSubspace(Q, height, &subQ);
1528:       PetscFECreate(PetscObjectComm((PetscObject) fe), &sub);
1529:       PetscObjectGetName((PetscObject) fe,  &name);
1530:       PetscObjectSetName((PetscObject) sub,  name);
1531:       PetscFEGetType(fe, &fetype);
1532:       PetscFESetType(sub, fetype);
1533:       PetscFESetBasisSpace(sub, subP);
1534:       PetscFESetDualSpace(sub, subQ);
1535:       PetscFESetNumComponents(sub, Nc);
1536:       PetscFESetUp(sub);
1537:       PetscFESetQuadrature(sub, subq);
1538:       fe->subspaces[height-1] = sub;
1539:     }
1540:     *subfe = fe->subspaces[height-1];
1541:   } else {
1542:     *subfe = NULL;
1543:   }
1544:   return(0);
1545: }

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

1552:   Collective on fem

1554:   Input Parameter:
1555: . fe - The initial PetscFE

1557:   Output Parameter:
1558: . feRef - The refined PetscFE

1560:   Level: advanced

1562: .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
1563: @*/
1564: PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)
1565: {
1566:   PetscSpace       P, Pref;
1567:   PetscDualSpace   Q, Qref;
1568:   DM               K, Kref;
1569:   PetscQuadrature  q, qref;
1570:   const PetscReal *v0, *jac;
1571:   PetscInt         numComp, numSubelements;
1572:   PetscErrorCode   ierr;

1575:   PetscFEGetBasisSpace(fe, &P);
1576:   PetscFEGetDualSpace(fe, &Q);
1577:   PetscFEGetQuadrature(fe, &q);
1578:   PetscDualSpaceGetDM(Q, &K);
1579:   /* Create space */
1580:   PetscObjectReference((PetscObject) P);
1581:   Pref = P;
1582:   /* Create dual space */
1583:   PetscDualSpaceDuplicate(Q, &Qref);
1584:   DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref);
1585:   PetscDualSpaceSetDM(Qref, Kref);
1586:   DMDestroy(&Kref);
1587:   PetscDualSpaceSetUp(Qref);
1588:   /* Create element */
1589:   PetscFECreate(PetscObjectComm((PetscObject) fe), feRef);
1590:   PetscFESetType(*feRef, PETSCFECOMPOSITE);
1591:   PetscFESetBasisSpace(*feRef, Pref);
1592:   PetscFESetDualSpace(*feRef, Qref);
1593:   PetscFEGetNumComponents(fe,    &numComp);
1594:   PetscFESetNumComponents(*feRef, numComp);
1595:   PetscFESetUp(*feRef);
1596:   PetscSpaceDestroy(&Pref);
1597:   PetscDualSpaceDestroy(&Qref);
1598:   /* Create quadrature */
1599:   PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL);
1600:   PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);
1601:   PetscFESetQuadrature(*feRef, qref);
1602:   PetscQuadratureDestroy(&qref);
1603:   return(0);
1604: }

1606: /*@C
1607:   PetscFECreateDefault - Create a PetscFE for basic FEM computation

1609:   Collective

1611:   Input Parameters:
1612: + comm      - The MPI comm
1613: . dim       - The spatial dimension
1614: . Nc        - The number of components
1615: . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
1616: . prefix    - The options prefix, or NULL
1617: - qorder    - The quadrature order

1619:   Output Parameter:
1620: . fem - The PetscFE object

1622:   Level: beginner

1624: .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
1625: @*/
1626: PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem)
1627: {
1628:   PetscQuadrature q, fq;
1629:   DM              K;
1630:   PetscSpace      P;
1631:   PetscDualSpace  Q;
1632:   PetscInt        order, quadPointsPerEdge;
1633:   PetscBool       tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
1634:   PetscErrorCode  ierr;

1637:   /* Create space */
1638:   PetscSpaceCreate(comm, &P);
1639:   PetscObjectSetOptionsPrefix((PetscObject) P, prefix);
1640:   PetscSpacePolynomialSetTensor(P, tensor);
1641:   PetscSpaceSetNumComponents(P, Nc);
1642:   PetscSpaceSetNumVariables(P, dim);
1643:   PetscSpaceSetFromOptions(P);
1644:   PetscSpaceSetUp(P);
1645:   PetscSpaceGetDegree(P, &order, NULL);
1646:   PetscSpacePolynomialGetTensor(P, &tensor);
1647:   /* Create dual space */
1648:   PetscDualSpaceCreate(comm, &Q);
1649:   PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);
1650:   PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);
1651:   PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);
1652:   PetscDualSpaceSetDM(Q, K);
1653:   DMDestroy(&K);
1654:   PetscDualSpaceSetNumComponents(Q, Nc);
1655:   PetscDualSpaceSetOrder(Q, order);
1656:   PetscDualSpaceLagrangeSetTensor(Q, tensor);
1657:   PetscDualSpaceSetFromOptions(Q);
1658:   PetscDualSpaceSetUp(Q);
1659:   /* Create element */
1660:   PetscFECreate(comm, fem);
1661:   PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);
1662:   PetscFESetBasisSpace(*fem, P);
1663:   PetscFESetDualSpace(*fem, Q);
1664:   PetscFESetNumComponents(*fem, Nc);
1665:   PetscFESetFromOptions(*fem);
1666:   PetscFESetUp(*fem);
1667:   PetscSpaceDestroy(&P);
1668:   PetscDualSpaceDestroy(&Q);
1669:   /* Create quadrature (with specified order if given) */
1670:   qorder = qorder >= 0 ? qorder : order;
1671:   PetscObjectOptionsBegin((PetscObject)*fem);
1672:   PetscOptionsBoundedInt("-petscfe_default_quadrature_order","Quadrature order is one less than quadrature points per edge","PetscFECreateDefault",qorder,&qorder,NULL,0);
1673:   PetscOptionsEnd();
1674:   quadPointsPerEdge = PetscMax(qorder + 1,1);
1675:   if (isSimplex) {
1676:     PetscDTGaussJacobiQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);
1677:     PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);
1678:   } else {
1679:     PetscDTGaussTensorQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);
1680:     PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);
1681:   }
1682:   PetscFESetQuadrature(*fem, q);
1683:   PetscFESetFaceQuadrature(*fem, fq);
1684:   PetscQuadratureDestroy(&q);
1685:   PetscQuadratureDestroy(&fq);
1686:   return(0);
1687: }

1689: /*@C
1690:   PetscFESetName - Names the FE and its subobjects

1692:   Not collective

1694:   Input Parameters:
1695: + fe   - The PetscFE
1696: - name - The name

1698:   Level: intermediate

1700: .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
1701: @*/
1702: PetscErrorCode PetscFESetName(PetscFE fe, const char name[])
1703: {
1704:   PetscSpace     P;
1705:   PetscDualSpace Q;

1709:   PetscFEGetBasisSpace(fe, &P);
1710:   PetscFEGetDualSpace(fe, &Q);
1711:   PetscObjectSetName((PetscObject) fe, name);
1712:   PetscObjectSetName((PetscObject) P,  name);
1713:   PetscObjectSetName((PetscObject) Q,  name);
1714:   return(0);
1715: }

1717: 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[])
1718: {
1719:   PetscInt       dOffset = 0, fOffset = 0, f;

1722:   for (f = 0; f < Nf; ++f) {
1723:     PetscFE          fe;
1724:     const PetscInt   Nbf = Nb[f], Ncf = Nc[f];
1725:     const PetscReal *Bq = &basisField[f][q*Nbf*Ncf];
1726:     const PetscReal *Dq = &basisFieldDer[f][q*Nbf*Ncf*dim];
1727:     PetscInt         b, c, d;

1729:     PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);
1730:     for (c = 0; c < Ncf; ++c)     u[fOffset+c] = 0.0;
1731:     for (d = 0; d < dim*Ncf; ++d) u_x[fOffset*dim+d] = 0.0;
1732:     for (b = 0; b < Nbf; ++b) {
1733:       for (c = 0; c < Ncf; ++c) {
1734:         const PetscInt cidx = b*Ncf+c;

1736:         u[fOffset+c] += Bq[cidx]*coefficients[dOffset+b];
1737:         for (d = 0; d < dim; ++d) u_x[(fOffset+c)*dim+d] += Dq[cidx*dim+d]*coefficients[dOffset+b];
1738:       }
1739:     }
1740:     PetscFEPushforward(fe, fegeom, 1, &u[fOffset]);
1741:     PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset*dim]);
1742:     if (u_t) {
1743:       for (c = 0; c < Ncf; ++c) u_t[fOffset+c] = 0.0;
1744:       for (b = 0; b < Nbf; ++b) {
1745:         for (c = 0; c < Ncf; ++c) {
1746:           const PetscInt cidx = b*Ncf+c;

1748:           u_t[fOffset+c] += Bq[cidx]*coefficients_t[dOffset+b];
1749:         }
1750:       }
1751:       PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]);
1752:     }
1753:     fOffset += Ncf;
1754:     dOffset += Nbf;
1755:   }
1756:   return 0;
1757: }

1759: PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[])
1760: {
1761:   PetscFE        fe;
1762:   PetscReal     *faceBasis;
1763:   PetscInt       Nb, Nc, b, c;

1766:   if (!prob) return 0;
1767:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1768:   PetscFEGetDimension(fe, &Nb);
1769:   PetscFEGetNumComponents(fe, &Nc);
1770:   PetscFEGetFaceCentroidTabulation(fe, &faceBasis);
1771:   for (c = 0; c < Nc; ++c) {u[c] = 0.0;}
1772:   for (b = 0; b < Nb; ++b) {
1773:     for (c = 0; c < Nc; ++c) {
1774:       const PetscInt cidx = b*Nc+c;

1776:       u[c] += coefficients[cidx]*faceBasis[faceLoc*Nb*Nc+cidx];
1777:     }
1778:   }
1779:   return 0;
1780: }

1782: 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[])
1783: {
1784:   PetscInt       q, b, c, d;

1787:   for (b = 0; b < Nb; ++b) elemVec[b] = 0.0;
1788:   for (q = 0; q < Nq; ++q) {
1789:     for (b = 0; b < Nb; ++b) {
1790:       for (c = 0; c < Nc; ++c) {
1791:         const PetscInt bcidx = b*Nc+c;

1793:         tmpBasis[bcidx] = basis[q*Nb*Nc+bcidx];
1794:         for (d = 0; d < dim; ++d) tmpBasisDer[bcidx*dim+d] = basisDer[q*Nb*Nc*dim+bcidx*dim+d];
1795:       }
1796:     }
1797:     PetscFEPushforward(fe, fegeom, Nb, tmpBasis);
1798:     PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer);
1799:     for (b = 0; b < Nb; ++b) {
1800:       for (c = 0; c < Nc; ++c) {
1801:         const PetscInt bcidx = b*Nc+c;
1802:         const PetscInt qcidx = q*Nc+c;

1804:         elemVec[b] += tmpBasis[bcidx]*f0[qcidx];
1805:         for (d = 0; d < dim; ++d) elemVec[b] += tmpBasisDer[bcidx*dim+d]*f1[qcidx*dim+d];
1806:       }
1807:     }
1808:   }
1809:   return(0);
1810: }

1812: 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[])
1813: {
1814:   PetscInt       f, fc, g, gc, df, dg;

1817:   for (f = 0; f < NbI; ++f) {
1818:     for (fc = 0; fc < NcI; ++fc) {
1819:       const PetscInt fidx = f*NcI+fc; /* Test function basis index */

1821:       tmpBasisI[fidx] = basisI[fidx];
1822:       for (df = 0; df < dim; ++df) tmpBasisDerI[fidx*dim+df] = basisDerI[fidx*dim+df];
1823:     }
1824:   }
1825:   PetscFEPushforward(feI, fegeom, NbI, tmpBasisI);
1826:   PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI);
1827:   for (g = 0; g < NbJ; ++g) {
1828:     for (gc = 0; gc < NcJ; ++gc) {
1829:       const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */

1831:       tmpBasisJ[gidx] = basisJ[gidx];
1832:       for (dg = 0; dg < dim; ++dg) tmpBasisDerJ[gidx*dim+dg] = basisDerJ[gidx*dim+dg];
1833:     }
1834:   }
1835:   PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ);
1836:   PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ);
1837:   for (f = 0; f < NbI; ++f) {
1838:     for (fc = 0; fc < NcI; ++fc) {
1839:       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
1840:       const PetscInt i    = offsetI+f; /* Element matrix row */
1841:       for (g = 0; g < NbJ; ++g) {
1842:         for (gc = 0; gc < NcJ; ++gc) {
1843:           const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
1844:           const PetscInt j    = offsetJ+g; /* Element matrix column */
1845:           const PetscInt fOff = eOffset+i*totDim+j;

1847:           elemMat[fOff] += tmpBasisI[fidx]*g0[fc*NcJ+gc]*tmpBasisJ[gidx];
1848:           for (df = 0; df < dim; ++df) {
1849:             elemMat[fOff] += tmpBasisI[fidx]*g1[(fc*NcJ+gc)*dim+df]*tmpBasisDerJ[gidx*dim+df];
1850:             elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g2[(fc*NcJ+gc)*dim+df]*tmpBasisJ[gidx];
1851:             for (dg = 0; dg < dim; ++dg) {
1852:               elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g3[((fc*NcJ+gc)*dim+df)*dim+dg]*tmpBasisDerJ[gidx*dim+dg];
1853:             }
1854:           }
1855:         }
1856:       }
1857:     }
1858:   }
1859:   return(0);
1860: }

1862: PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom)
1863: {
1864:   PetscDualSpace  dsp;
1865:   DM              dm;
1866:   PetscQuadrature quadDef;
1867:   PetscInt        dim, cdim, Nq;
1868:   PetscErrorCode  ierr;

1871:   PetscFEGetDualSpace(fe, &dsp);
1872:   PetscDualSpaceGetDM(dsp, &dm);
1873:   DMGetDimension(dm, &dim);
1874:   DMGetCoordinateDim(dm, &cdim);
1875:   PetscFEGetQuadrature(fe, &quadDef);
1876:   quad = quad ? quad : quadDef;
1877:   PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);
1878:   PetscMalloc1(Nq*cdim,      &cgeom->v);
1879:   PetscMalloc1(Nq*cdim*cdim, &cgeom->J);
1880:   PetscMalloc1(Nq*cdim*cdim, &cgeom->invJ);
1881:   PetscMalloc1(Nq,           &cgeom->detJ);
1882:   cgeom->dim       = dim;
1883:   cgeom->dimEmbed  = cdim;
1884:   cgeom->numCells  = 1;
1885:   cgeom->numPoints = Nq;
1886:   DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ);
1887:   return(0);
1888: }

1890: PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom)
1891: {

1895:   PetscFree(cgeom->v);
1896:   PetscFree(cgeom->J);
1897:   PetscFree(cgeom->invJ);
1898:   PetscFree(cgeom->detJ);
1899:   return(0);
1900: }