Actual source code: fe.c

petsc-master 2020-01-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:    PetscFEViewFromOptions - View from Options

169:    Collective on PetscFE

171:    Input Parameters:
172: +  A - the PetscFE object
173: .  obj - Optional object
174: -  name - command line option

176:    Level: intermediate
177: .seealso:  PetscFE(), PetscFEView(), PetscObjectViewFromOptions(), PetscFECreate()
178: @*/
179: PetscErrorCode  PetscFEViewFromOptions(PetscFE A,PetscObject obj,const char name[])
180: {

185:   PetscObjectViewFromOptions((PetscObject)A,obj,name);
186:   return(0);
187: }

189: /*@C
190:   PetscFEView - Views a PetscFE

192:   Collective on fem

194:   Input Parameter:
195: + fem - the PetscFE object to view
196: - viewer   - the viewer

198:   Level: beginner

200: .seealso PetscFEDestroy()
201: @*/
202: PetscErrorCode PetscFEView(PetscFE fem, PetscViewer viewer)
203: {
204:   PetscBool      iascii;

210:   if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &viewer);}
211:   PetscObjectPrintClassNamePrefixType((PetscObject)fem, viewer);
212:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
213:   if (fem->ops->view) {(*fem->ops->view)(fem, viewer);}
214:   return(0);
215: }

217: /*@
218:   PetscFESetFromOptions - sets parameters in a PetscFE from the options database

220:   Collective on fem

222:   Input Parameter:
223: . fem - the PetscFE object to set options for

225:   Options Database:
226: + -petscfe_num_blocks  - the number of cell blocks to integrate concurrently
227: - -petscfe_num_batches - the number of cell batches to integrate serially

229:   Level: intermediate

231: .seealso PetscFEView()
232: @*/
233: PetscErrorCode PetscFESetFromOptions(PetscFE fem)
234: {
235:   const char    *defaultType;
236:   char           name[256];
237:   PetscBool      flg;

242:   if (!((PetscObject) fem)->type_name) {
243:     defaultType = PETSCFEBASIC;
244:   } else {
245:     defaultType = ((PetscObject) fem)->type_name;
246:   }
247:   if (!PetscFERegisterAllCalled) {PetscFERegisterAll();}

249:   PetscObjectOptionsBegin((PetscObject) fem);
250:   PetscOptionsFList("-petscfe_type", "Finite element space", "PetscFESetType", PetscFEList, defaultType, name, 256, &flg);
251:   if (flg) {
252:     PetscFESetType(fem, name);
253:   } else if (!((PetscObject) fem)->type_name) {
254:     PetscFESetType(fem, defaultType);
255:   }
256:   PetscOptionsBoundedInt("-petscfe_num_blocks", "The number of cell blocks to integrate concurrently", "PetscSpaceSetTileSizes", fem->numBlocks, &fem->numBlocks, NULL,1);
257:   PetscOptionsBoundedInt("-petscfe_num_batches", "The number of cell batches to integrate serially", "PetscSpaceSetTileSizes", fem->numBatches, &fem->numBatches, NULL,1);
258:   if (fem->ops->setfromoptions) {
259:     (*fem->ops->setfromoptions)(PetscOptionsObject,fem);
260:   }
261:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
262:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) fem);
263:   PetscOptionsEnd();
264:   PetscFEViewFromOptions(fem, NULL, "-petscfe_view");
265:   return(0);
266: }

268: /*@C
269:   PetscFESetUp - Construct data structures for the PetscFE

271:   Collective on fem

273:   Input Parameter:
274: . fem - the PetscFE object to setup

276:   Level: intermediate

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

286:   if (fem->setupcalled) return(0);
287:   fem->setupcalled = PETSC_TRUE;
288:   if (fem->ops->setup) {(*fem->ops->setup)(fem);}
289:   return(0);
290: }

292: /*@
293:   PetscFEDestroy - Destroys a PetscFE object

295:   Collective on fem

297:   Input Parameter:
298: . fem - the PetscFE object to destroy

300:   Level: beginner

302: .seealso PetscFEView()
303: @*/
304: PetscErrorCode PetscFEDestroy(PetscFE *fem)
305: {

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

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

315:   if ((*fem)->subspaces) {
316:     PetscInt dim, d;

318:     PetscDualSpaceGetDimension((*fem)->dualSpace, &dim);
319:     for (d = 0; d < dim; ++d) {PetscFEDestroy(&(*fem)->subspaces[d]);}
320:   }
321:   PetscFree((*fem)->subspaces);
322:   PetscFree((*fem)->invV);
323:   PetscTabulationDestroy(&(*fem)->T);
324:   PetscTabulationDestroy(&(*fem)->Tf);
325:   PetscTabulationDestroy(&(*fem)->Tc);
326:   PetscSpaceDestroy(&(*fem)->basisSpace);
327:   PetscDualSpaceDestroy(&(*fem)->dualSpace);
328:   PetscQuadratureDestroy(&(*fem)->quadrature);
329:   PetscQuadratureDestroy(&(*fem)->faceQuadrature);

331:   if ((*fem)->ops->destroy) {(*(*fem)->ops->destroy)(*fem);}
332:   PetscHeaderDestroy(fem);
333:   return(0);
334: }

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

339:   Collective

341:   Input Parameter:
342: . comm - The communicator for the PetscFE object

344:   Output Parameter:
345: . fem - The PetscFE object

347:   Level: beginner

349: .seealso: PetscFESetType(), PETSCFEGALERKIN
350: @*/
351: PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem)
352: {
353:   PetscFE        f;

358:   PetscCitationsRegister(FECitation,&FEcite);
359:   *fem = NULL;
360:   PetscFEInitializePackage();

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

364:   f->basisSpace    = NULL;
365:   f->dualSpace     = NULL;
366:   f->numComponents = 1;
367:   f->subspaces     = NULL;
368:   f->invV          = NULL;
369:   f->T             = NULL;
370:   f->Tf            = NULL;
371:   f->Tc            = NULL;
372:   PetscArrayzero(&f->quadrature, 1);
373:   PetscArrayzero(&f->faceQuadrature, 1);
374:   f->blockSize     = 0;
375:   f->numBlocks     = 1;
376:   f->batchSize     = 0;
377:   f->numBatches    = 1;

379:   *fem = f;
380:   return(0);
381: }

383: /*@
384:   PetscFEGetSpatialDimension - Returns the spatial dimension of the element

386:   Not collective

388:   Input Parameter:
389: . fem - The PetscFE object

391:   Output Parameter:
392: . dim - The spatial dimension

394:   Level: intermediate

396: .seealso: PetscFECreate()
397: @*/
398: PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim)
399: {
400:   DM             dm;

406:   PetscDualSpaceGetDM(fem->dualSpace, &dm);
407:   DMGetDimension(dm, dim);
408:   return(0);
409: }

411: /*@
412:   PetscFESetNumComponents - Sets the number of components in the element

414:   Not collective

416:   Input Parameters:
417: + fem - The PetscFE object
418: - comp - The number of field components

420:   Level: intermediate

422: .seealso: PetscFECreate()
423: @*/
424: PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp)
425: {
428:   fem->numComponents = comp;
429:   return(0);
430: }

432: /*@
433:   PetscFEGetNumComponents - Returns the number of components in the element

435:   Not collective

437:   Input Parameter:
438: . fem - The PetscFE object

440:   Output Parameter:
441: . comp - The number of field components

443:   Level: intermediate

445: .seealso: PetscFECreate()
446: @*/
447: PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp)
448: {
452:   *comp = fem->numComponents;
453:   return(0);
454: }

456: /*@
457:   PetscFESetTileSizes - Sets the tile sizes for evaluation

459:   Not collective

461:   Input Parameters:
462: + fem - The PetscFE object
463: . blockSize - The number of elements in a block
464: . numBlocks - The number of blocks in a batch
465: . batchSize - The number of elements in a batch
466: - numBatches - The number of batches in a chunk

468:   Level: intermediate

470: .seealso: PetscFECreate()
471: @*/
472: PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches)
473: {
476:   fem->blockSize  = blockSize;
477:   fem->numBlocks  = numBlocks;
478:   fem->batchSize  = batchSize;
479:   fem->numBatches = numBatches;
480:   return(0);
481: }

483: /*@
484:   PetscFEGetTileSizes - Returns the tile sizes for evaluation

486:   Not collective

488:   Input Parameter:
489: . fem - The PetscFE object

491:   Output Parameters:
492: + blockSize - The number of elements in a block
493: . numBlocks - The number of blocks in a batch
494: . batchSize - The number of elements in a batch
495: - numBatches - The number of batches in a chunk

497:   Level: intermediate

499: .seealso: PetscFECreate()
500: @*/
501: PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PetscInt *blockSize, PetscInt *numBlocks, PetscInt *batchSize, PetscInt *numBatches)
502: {
509:   if (blockSize)  *blockSize  = fem->blockSize;
510:   if (numBlocks)  *numBlocks  = fem->numBlocks;
511:   if (batchSize)  *batchSize  = fem->batchSize;
512:   if (numBatches) *numBatches = fem->numBatches;
513:   return(0);
514: }

516: /*@
517:   PetscFEGetBasisSpace - Returns the PetscSpace used for approximation of the solution

519:   Not collective

521:   Input Parameter:
522: . fem - The PetscFE object

524:   Output Parameter:
525: . sp - The PetscSpace object

527:   Level: intermediate

529: .seealso: PetscFECreate()
530: @*/
531: PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp)
532: {
536:   *sp = fem->basisSpace;
537:   return(0);
538: }

540: /*@
541:   PetscFESetBasisSpace - Sets the PetscSpace used for approximation of the solution

543:   Not collective

545:   Input Parameters:
546: + fem - The PetscFE object
547: - sp - The PetscSpace object

549:   Level: intermediate

551: .seealso: PetscFECreate()
552: @*/
553: PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp)
554: {

560:   PetscSpaceDestroy(&fem->basisSpace);
561:   fem->basisSpace = sp;
562:   PetscObjectReference((PetscObject) fem->basisSpace);
563:   return(0);
564: }

566: /*@
567:   PetscFEGetDualSpace - Returns the PetscDualSpace used to define the inner product

569:   Not collective

571:   Input Parameter:
572: . fem - The PetscFE object

574:   Output Parameter:
575: . sp - The PetscDualSpace object

577:   Level: intermediate

579: .seealso: PetscFECreate()
580: @*/
581: PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp)
582: {
586:   *sp = fem->dualSpace;
587:   return(0);
588: }

590: /*@
591:   PetscFESetDualSpace - Sets the PetscDualSpace used to define the inner product

593:   Not collective

595:   Input Parameters:
596: + fem - The PetscFE object
597: - sp - The PetscDualSpace object

599:   Level: intermediate

601: .seealso: PetscFECreate()
602: @*/
603: PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp)
604: {

610:   PetscDualSpaceDestroy(&fem->dualSpace);
611:   fem->dualSpace = sp;
612:   PetscObjectReference((PetscObject) fem->dualSpace);
613:   return(0);
614: }

616: /*@
617:   PetscFEGetQuadrature - Returns the PetscQuadrature used to calculate inner products

619:   Not collective

621:   Input Parameter:
622: . fem - The PetscFE object

624:   Output Parameter:
625: . q - The PetscQuadrature object

627:   Level: intermediate

629: .seealso: PetscFECreate()
630: @*/
631: PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q)
632: {
636:   *q = fem->quadrature;
637:   return(0);
638: }

640: /*@
641:   PetscFESetQuadrature - Sets the PetscQuadrature used to calculate inner products

643:   Not collective

645:   Input Parameters:
646: + fem - The PetscFE object
647: - q - The PetscQuadrature object

649:   Level: intermediate

651: .seealso: PetscFECreate()
652: @*/
653: PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q)
654: {
655:   PetscInt       Nc, qNc;

660:   PetscFEGetNumComponents(fem, &Nc);
661:   PetscQuadratureGetNumComponents(q, &qNc);
662:   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);
663:   PetscTabulationDestroy(&fem->T);
664:   PetscTabulationDestroy(&fem->Tc);
665:   PetscQuadratureDestroy(&fem->quadrature);
666:   fem->quadrature = q;
667:   PetscObjectReference((PetscObject) q);
668:   return(0);
669: }

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

674:   Not collective

676:   Input Parameter:
677: . fem - The PetscFE object

679:   Output Parameter:
680: . q - The PetscQuadrature object

682:   Level: intermediate

684: .seealso: PetscFECreate()
685: @*/
686: PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q)
687: {
691:   *q = fem->faceQuadrature;
692:   return(0);
693: }

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

698:   Not collective

700:   Input Parameters:
701: + fem - The PetscFE object
702: - q - The PetscQuadrature object

704:   Level: intermediate

706: .seealso: PetscFECreate()
707: @*/
708: PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q)
709: {
710:   PetscInt       Nc, qNc;

715:   PetscFEGetNumComponents(fem, &Nc);
716:   PetscQuadratureGetNumComponents(q, &qNc);
717:   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);
718:   PetscTabulationDestroy(&fem->Tf);
719:   PetscQuadratureDestroy(&fem->faceQuadrature);
720:   fem->faceQuadrature = q;
721:   PetscObjectReference((PetscObject) q);
722:   return(0);
723: }

725: /*@
726:   PetscFECopyQuadrature - Copy both volumetric and surface quadrature

728:   Not collective

730:   Input Parameters:
731: + sfe - The PetscFE source for the quadratures
732: - tfe - The PetscFE target for the quadratures

734:   Level: intermediate

736: .seealso: PetscFECreate(), PetscFESetQuadrature(), PetscFESetFaceQuadrature()
737: @*/
738: PetscErrorCode PetscFECopyQuadrature(PetscFE sfe, PetscFE tfe)
739: {
740:   PetscQuadrature q;
741:   PetscErrorCode  ierr;

746:   PetscFEGetQuadrature(sfe, &q);
747:   PetscFESetQuadrature(tfe,  q);
748:   PetscFEGetFaceQuadrature(sfe, &q);
749:   PetscFESetFaceQuadrature(tfe,  q);
750:   return(0);
751: }

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

756:   Not collective

758:   Input Parameter:
759: . fem - The PetscFE object

761:   Output Parameter:
762: . numDof - Array with the number of dofs per dimension

764:   Level: intermediate

766: .seealso: PetscFECreate()
767: @*/
768: PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt **numDof)
769: {

775:   PetscDualSpaceGetNumDof(fem->dualSpace, numDof);
776:   return(0);
777: }

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

782:   Not collective

784:   Input Parameter:
785: . fem - The PetscFE object

787:   Output Parameter:
788: . T - The basis function values and derivatives at quadrature points

790:   Note:
791: $ T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
792: $ T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
793: $ T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e

795:   Level: intermediate

797: .seealso: PetscFECreateTabulation(), PetscTabulationDestroy()
798: @*/
799: PetscErrorCode PetscFEGetCellTabulation(PetscFE fem, PetscTabulation *T)
800: {
801:   PetscInt         npoints;
802:   const PetscReal *points;
803:   PetscErrorCode   ierr;

808:   PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL);
809:   if (!fem->T) {PetscFECreateTabulation(fem, 1, npoints, points, 1, &fem->T);}
810:   *T = fem->T;
811:   return(0);
812: }

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

817:   Not collective

819:   Input Parameter:
820: . fem - The PetscFE object

822:   Output Parameters:
823: . Tf - The basis function values and derviatives at face quadrature points

825:   Note:
826: $ T->T[0] = Bf[((f*Nq + q)*pdim + i)*Nc + c] is the value at point f,q for basis function i and component c
827: $ T->T[1] = Df[(((f*Nq + q)*pdim + i)*Nc + c)*dim + d] is the derivative value at point f,q for basis function i, component c, in direction d
828: $ T->T[2] = Hf[((((f*Nq + q)*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point f,q for basis function i, component c, in directions d and e

830:   Level: intermediate

832: .seealso: PetscFEGetCellTabulation(), PetscFECreateTabulation(), PetscTabulationDestroy()
833: @*/
834: PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscTabulation *Tf)
835: {
836:   PetscErrorCode   ierr;

841:   if (!fem->Tf) {
842:     const PetscReal  xi0[3] = {-1., -1., -1.};
843:     PetscReal        v0[3], J[9], detJ;
844:     PetscQuadrature  fq;
845:     PetscDualSpace   sp;
846:     DM               dm;
847:     const PetscInt  *faces;
848:     PetscInt         dim, numFaces, f, npoints, q;
849:     const PetscReal *points;
850:     PetscReal       *facePoints;

852:     PetscFEGetDualSpace(fem, &sp);
853:     PetscDualSpaceGetDM(sp, &dm);
854:     DMGetDimension(dm, &dim);
855:     DMPlexGetConeSize(dm, 0, &numFaces);
856:     DMPlexGetCone(dm, 0, &faces);
857:     PetscFEGetFaceQuadrature(fem, &fq);
858:     if (fq) {
859:       PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL);
860:       PetscMalloc1(numFaces*npoints*dim, &facePoints);
861:       for (f = 0; f < numFaces; ++f) {
862:         DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ);
863:         for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim-1, xi0, v0, J, &points[q*(dim-1)], &facePoints[(f*npoints+q)*dim]);
864:       }
865:       PetscFECreateTabulation(fem, numFaces, npoints, facePoints, 1, &fem->Tf);
866:       PetscFree(facePoints);
867:     }
868:   }
869:   *Tf = fem->Tf;
870:   return(0);
871: }

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

876:   Not collective

878:   Input Parameter:
879: . fem - The PetscFE object

881:   Output Parameters:
882: . Tc - The basis function values at face centroid points

884:   Note:
885: $ T->T[0] = Bf[(f*pdim + i)*Nc + c] is the value at point f for basis function i and component c

887:   Level: intermediate

889: .seealso: PetscFEGetFaceTabulation(), PetscFEGetCellTabulation(), PetscFECreateTabulation(), PetscTabulationDestroy()
890: @*/
891: PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscTabulation *Tc)
892: {
893:   PetscErrorCode   ierr;

898:   if (!fem->Tc) {
899:     PetscDualSpace  sp;
900:     DM              dm;
901:     const PetscInt *cone;
902:     PetscReal      *centroids;
903:     PetscInt        dim, numFaces, f;

905:     PetscFEGetDualSpace(fem, &sp);
906:     PetscDualSpaceGetDM(sp, &dm);
907:     DMGetDimension(dm, &dim);
908:     DMPlexGetConeSize(dm, 0, &numFaces);
909:     DMPlexGetCone(dm, 0, &cone);
910:     PetscMalloc1(numFaces*dim, &centroids);
911:     for (f = 0; f < numFaces; ++f) {DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, &centroids[f*dim], NULL);}
912:     PetscFECreateTabulation(fem, 1, numFaces, centroids, 0, &fem->Tc);
913:     PetscFree(centroids);
914:   }
915:   *Tc = fem->Tc;
916:   return(0);
917: }

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

922:   Not collective

924:   Input Parameters:
925: + fem     - The PetscFE object
926: . nrepl   - The number of replicas
927: . npoints - The number of tabulation points in a replica
928: . points  - The tabulation point coordinates
929: - K       - The number of derivatives calculated

931:   Output Parameter:
932: . T - The basis function values and derivatives at tabulation points

934:   Note:
935: $ T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
936: $ T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
937: $ T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e

939:   Level: intermediate

941: .seealso: PetscFEGetCellTabulation(), PetscTabulationDestroy()
942: @*/
943: PetscErrorCode PetscFECreateTabulation(PetscFE fem, PetscInt nrepl, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation *T)
944: {
945:   DM               dm;
946:   PetscDualSpace   Q;
947:   PetscInt         Nb;   /* Dimension of FE space P */
948:   PetscInt         Nc;   /* Field components */
949:   PetscInt         cdim; /* Reference coordinate dimension */
950:   PetscInt         k;
951:   PetscErrorCode   ierr;

954:   if (!npoints || !fem->dualSpace || K < 0) {
955:     *T = NULL;
956:     return(0);
957:   }
961:   PetscFEGetDualSpace(fem, &Q);
962:   PetscDualSpaceGetDM(Q, &dm);
963:   DMGetDimension(dm, &cdim);
964:   PetscDualSpaceGetDimension(Q, &Nb);
965:   PetscFEGetNumComponents(fem, &Nc);
966:   PetscMalloc1(1, T);
967:   (*T)->K    = !cdim ? 0 : K;
968:   (*T)->Nr   = nrepl;
969:   (*T)->Np   = npoints;
970:   (*T)->Nb   = Nb;
971:   (*T)->Nc   = Nc;
972:   (*T)->cdim = cdim;
973:   PetscMalloc1((*T)->K+1, &(*T)->T);
974:   for (k = 0; k <= (*T)->K; ++k) {
975:     PetscMalloc1(nrepl*npoints*Nb*Nc*PetscPowInt(cdim, k), &(*T)->T[k]);
976:   }
977:   (*fem->ops->createtabulation)(fem, nrepl*npoints, points, K, *T);
978:   return(0);
979: }

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

984:   Not collective

986:   Input Parameters:
987: + fem     - The PetscFE object
988: . npoints - The number of tabulation points
989: . points  - The tabulation point coordinates
990: . K       - The number of derivatives calculated
991: - T       - An existing tabulation object with enough allocated space

993:   Output Parameter:
994: . T - The basis function values and derivatives at tabulation points

996:   Note:
997: $ T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
998: $ T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
999: $ T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e

1001:   Level: intermediate

1003: .seealso: PetscFEGetCellTabulation(), PetscTabulationDestroy()
1004: @*/
1005: PetscErrorCode PetscFEComputeTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
1006: {

1010:   if (!npoints || !fem->dualSpace || K < 0) return(0);
1014: #ifdef PETSC_USE_DEBUG
1015:   {
1016:     DM               dm;
1017:     PetscDualSpace   Q;
1018:     PetscInt         Nb;   /* Dimension of FE space P */
1019:     PetscInt         Nc;   /* Field components */
1020:     PetscInt         cdim; /* Reference coordinate dimension */

1022:     PetscFEGetDualSpace(fem, &Q);
1023:     PetscDualSpaceGetDM(Q, &dm);
1024:     DMGetDimension(dm, &cdim);
1025:     PetscDualSpaceGetDimension(Q, &Nb);
1026:     PetscFEGetNumComponents(fem, &Nc);
1027:     if (T->K    != (!cdim ? 0 : K)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation K %D must match requested K %D", T->K, !cdim ? 0 : K);
1028:     if (T->Nb   != Nb)              SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nb %D must match requested Nb %D", T->Nb, Nb);
1029:     if (T->Nc   != Nc)              SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nc %D must match requested Nc %D", T->Nc, Nc);
1030:     if (T->cdim != cdim)            SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation cdim %D must match requested cdim %D", T->cdim, cdim);
1031:   }
1032: #endif
1033:   T->Nr = 1;
1034:   T->Np = npoints;
1035:   (*fem->ops->createtabulation)(fem, npoints, points, K, T);
1036:   return(0);
1037: }

1039: /*@C
1040:   PetscTabulationDestroy - Frees memory from the associated tabulation.

1042:   Not collective

1044:   Input Parameter:
1045: . T - The tabulation

1047:   Level: intermediate

1049: .seealso: PetscFECreateTabulation(), PetscFEGetCellTabulation()
1050: @*/
1051: PetscErrorCode PetscTabulationDestroy(PetscTabulation *T)
1052: {
1053:   PetscInt       k;

1058:   if (!T || !(*T)) return(0);
1059:   for (k = 0; k <= (*T)->K; ++k) {PetscFree((*T)->T[k]);}
1060:   PetscFree((*T)->T);
1061:   PetscFree(*T);
1062:   *T = NULL;
1063:   return(0);
1064: }

1066: PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
1067: {
1068:   PetscSpace     bsp, bsubsp;
1069:   PetscDualSpace dsp, dsubsp;
1070:   PetscInt       dim, depth, numComp, i, j, coneSize, order;
1071:   PetscFEType    type;
1072:   DM             dm;
1073:   DMLabel        label;
1074:   PetscReal      *xi, *v, *J, detJ;
1075:   const char     *name;
1076:   PetscQuadrature origin, fullQuad, subQuad;

1082:   PetscFEGetBasisSpace(fe,&bsp);
1083:   PetscFEGetDualSpace(fe,&dsp);
1084:   PetscDualSpaceGetDM(dsp,&dm);
1085:   DMGetDimension(dm,&dim);
1086:   DMPlexGetDepthLabel(dm,&label);
1087:   DMLabelGetValue(label,refPoint,&depth);
1088:   PetscCalloc1(depth,&xi);
1089:   PetscMalloc1(dim,&v);
1090:   PetscMalloc1(dim*dim,&J);
1091:   for (i = 0; i < depth; i++) xi[i] = 0.;
1092:   PetscQuadratureCreate(PETSC_COMM_SELF,&origin);
1093:   PetscQuadratureSetData(origin,depth,0,1,xi,NULL);
1094:   DMPlexComputeCellGeometryFEM(dm,refPoint,origin,v,J,NULL,&detJ);
1095:   /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */
1096:   for (i = 1; i < dim; i++) {
1097:     for (j = 0; j < depth; j++) {
1098:       J[i * depth + j] = J[i * dim + j];
1099:     }
1100:   }
1101:   PetscQuadratureDestroy(&origin);
1102:   PetscDualSpaceGetPointSubspace(dsp,refPoint,&dsubsp);
1103:   PetscSpaceCreateSubspace(bsp,dsubsp,v,J,NULL,NULL,PETSC_OWN_POINTER,&bsubsp);
1104:   PetscSpaceSetUp(bsubsp);
1105:   PetscFECreate(PetscObjectComm((PetscObject)fe),trFE);
1106:   PetscFEGetType(fe,&type);
1107:   PetscFESetType(*trFE,type);
1108:   PetscFEGetNumComponents(fe,&numComp);
1109:   PetscFESetNumComponents(*trFE,numComp);
1110:   PetscFESetBasisSpace(*trFE,bsubsp);
1111:   PetscFESetDualSpace(*trFE,dsubsp);
1112:   PetscObjectGetName((PetscObject) fe, &name);
1113:   if (name) {PetscFESetName(*trFE, name);}
1114:   PetscFEGetQuadrature(fe,&fullQuad);
1115:   PetscQuadratureGetOrder(fullQuad,&order);
1116:   DMPlexGetConeSize(dm,refPoint,&coneSize);
1117:   if (coneSize == 2 * depth) {
1118:     PetscDTGaussTensorQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);
1119:   } else {
1120:     PetscDTGaussJacobiQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);
1121:   }
1122:   PetscFESetQuadrature(*trFE,subQuad);
1123:   PetscFESetUp(*trFE);
1124:   PetscQuadratureDestroy(&subQuad);
1125:   PetscSpaceDestroy(&bsubsp);
1126:   return(0);
1127: }

1129: PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE)
1130: {
1131:   PetscInt       hStart, hEnd;
1132:   PetscDualSpace dsp;
1133:   DM             dm;

1139:   *trFE = NULL;
1140:   PetscFEGetDualSpace(fe,&dsp);
1141:   PetscDualSpaceGetDM(dsp,&dm);
1142:   DMPlexGetHeightStratum(dm,height,&hStart,&hEnd);
1143:   if (hEnd <= hStart) return(0);
1144:   PetscFECreatePointTrace(fe,hStart,trFE);
1145:   return(0);
1146: }


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

1152:   Not collective

1154:   Input Parameter:
1155: . fe - The PetscFE

1157:   Output Parameter:
1158: . dim - The dimension

1160:   Level: intermediate

1162: .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
1163: @*/
1164: PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
1165: {

1171:   if (fem->ops->getdimension) {(*fem->ops->getdimension)(fem, dim);}
1172:   return(0);
1173: }

1175: /*@C
1176:   PetscFEPushforward - Map the reference element function to real space

1178:   Input Parameters:
1179: + fe     - The PetscFE
1180: . fegeom - The cell geometry
1181: . Nv     - The number of function values
1182: - vals   - The function values

1184:   Output Parameter:
1185: . vals   - The transformed function values

1187:   Level: advanced

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

1191: .seealso: PetscDualSpacePushforward()
1192: @*/
1193: PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1194: {

1198:   PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);
1199:   return(0);
1200: }

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

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

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

1214:   Level: advanced

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

1218: .seealso: PetscFEPushforward(), PetscDualSpacePushforwardGradient(), PetscDualSpacePushforward()
1219: @*/
1220: PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1221: {

1225:   PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);
1226:   return(0);
1227: }

1229: /*
1230: Purpose: Compute element vector for chunk of elements

1232: Input:
1233:   Sizes:
1234:      Ne:  number of elements
1235:      Nf:  number of fields
1236:      PetscFE
1237:        dim: spatial dimension
1238:        Nb:  number of basis functions
1239:        Nc:  number of field components
1240:        PetscQuadrature
1241:          Nq:  number of quadrature points

1243:   Geometry:
1244:      PetscFEGeom[Ne] possibly *Nq
1245:        PetscReal v0s[dim]
1246:        PetscReal n[dim]
1247:        PetscReal jacobians[dim*dim]
1248:        PetscReal jacobianInverses[dim*dim]
1249:        PetscReal jacobianDeterminants
1250:   FEM:
1251:      PetscFE
1252:        PetscQuadrature
1253:          PetscReal   quadPoints[Nq*dim]
1254:          PetscReal   quadWeights[Nq]
1255:        PetscReal   basis[Nq*Nb*Nc]
1256:        PetscReal   basisDer[Nq*Nb*Nc*dim]
1257:      PetscScalar coefficients[Ne*Nb*Nc]
1258:      PetscScalar elemVec[Ne*Nb*Nc]

1260:   Problem:
1261:      PetscInt f: the active field
1262:      f0, f1

1264:   Work Space:
1265:      PetscFE
1266:        PetscScalar f0[Nq*dim];
1267:        PetscScalar f1[Nq*dim*dim];
1268:        PetscScalar u[Nc];
1269:        PetscScalar gradU[Nc*dim];
1270:        PetscReal   x[dim];
1271:        PetscScalar realSpaceDer[dim];

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

1275: Input:
1276:   Sizes:
1277:      N_cb: Number of serial cell batches

1279:   Geometry:
1280:      PetscReal v0s[Ne*dim]
1281:      PetscReal jacobians[Ne*dim*dim]        possibly *Nq
1282:      PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq
1283:      PetscReal jacobianDeterminants[Ne]     possibly *Nq
1284:   FEM:
1285:      static PetscReal   quadPoints[Nq*dim]
1286:      static PetscReal   quadWeights[Nq]
1287:      static PetscReal   basis[Nq*Nb*Nc]
1288:      static PetscReal   basisDer[Nq*Nb*Nc*dim]
1289:      PetscScalar coefficients[Ne*Nb*Nc]
1290:      PetscScalar elemVec[Ne*Nb*Nc]

1292: ex62.c:
1293:   PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[],
1294:                                                const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[],
1295:                                                void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]),
1296:                                                void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[])

1298: ex52.c:
1299:   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)
1300:   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)

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

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

1309: ex52_integrateElementOpenCL.c:
1310: PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[],
1311:                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1312:                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)

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

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

1320:   Not collective

1322:   Input Parameters:
1323: + fem          - The PetscFE object for the field being integrated
1324: . prob         - The PetscDS specifying the discretizations and continuum functions
1325: . field        - The field being integrated
1326: . Ne           - The number of elements in the chunk
1327: . cgeom        - The cell geometry for each cell in the chunk
1328: . coefficients - The array of FEM basis coefficients for the elements
1329: . probAux      - The PetscDS specifying the auxiliary discretizations
1330: - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements

1332:   Output Parameter
1333: . integral     - the integral for this field

1335:   Level: intermediate

1337: .seealso: PetscFEIntegrateResidual()
1338: @*/
1339: PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1340:                                 const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1341: {
1342:   PetscFE        fe;

1347:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1348:   if (fe->ops->integrate) {(*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral);}
1349:   return(0);
1350: }

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

1355:   Not collective

1357:   Input Parameters:
1358: + fem          - The PetscFE object for the field being integrated
1359: . prob         - The PetscDS specifying the discretizations and continuum functions
1360: . field        - The field being integrated
1361: . obj_func     - The function to be integrated
1362: . Ne           - The number of elements in the chunk
1363: . fgeom        - The face geometry for each face in the chunk
1364: . coefficients - The array of FEM basis coefficients for the elements
1365: . probAux      - The PetscDS specifying the auxiliary discretizations
1366: - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements

1368:   Output Parameter
1369: . integral     - the integral for this field

1371:   Level: intermediate

1373: .seealso: PetscFEIntegrateResidual()
1374: @*/
1375: PetscErrorCode PetscFEIntegrateBd(PetscDS prob, PetscInt field,
1376:                                   void (*obj_func)(PetscInt, PetscInt, PetscInt,
1377:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1378:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1379:                                                    PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1380:                                   PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1381: {
1382:   PetscFE        fe;

1387:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1388:   if (fe->ops->integratebd) {(*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral);}
1389:   return(0);
1390: }

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

1395:   Not collective

1397:   Input Parameters:
1398: + fem          - The PetscFE object for the field being integrated
1399: . prob         - The PetscDS specifying the discretizations and continuum functions
1400: . field        - The field being integrated
1401: . Ne           - The number of elements in the chunk
1402: . cgeom        - The cell geometry for each cell in the chunk
1403: . coefficients - The array of FEM basis coefficients for the elements
1404: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1405: . probAux      - The PetscDS specifying the auxiliary discretizations
1406: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1407: - t            - The time

1409:   Output Parameter
1410: . elemVec      - the element residual vectors from each element

1412:   Note:
1413: $ Loop over batch of elements (e):
1414: $   Loop over quadrature points (q):
1415: $     Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q
1416: $     Call f_0 and f_1
1417: $   Loop over element vector entries (f,fc --> i):
1418: $     elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u)

1420:   Level: intermediate

1422: .seealso: PetscFEIntegrateResidual()
1423: @*/
1424: PetscErrorCode PetscFEIntegrateResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1425:                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1426: {
1427:   PetscFE        fe;

1432:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1433:   if (fe->ops->integrateresidual) {(*fe->ops->integrateresidual)(prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);}
1434:   return(0);
1435: }

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

1440:   Not collective

1442:   Input Parameters:
1443: + fem          - The PetscFE object for the field being integrated
1444: . prob         - The PetscDS specifying the discretizations and continuum functions
1445: . field        - The field being integrated
1446: . Ne           - The number of elements in the chunk
1447: . fgeom        - The face geometry for each cell in the chunk
1448: . coefficients - The array of FEM basis coefficients for the elements
1449: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1450: . probAux      - The PetscDS specifying the auxiliary discretizations
1451: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1452: - t            - The time

1454:   Output Parameter
1455: . elemVec      - the element residual vectors from each element

1457:   Level: intermediate

1459: .seealso: PetscFEIntegrateResidual()
1460: @*/
1461: PetscErrorCode PetscFEIntegrateBdResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom,
1462:                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1463: {
1464:   PetscFE        fe;

1469:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1470:   if (fe->ops->integratebdresidual) {(*fe->ops->integratebdresidual)(prob, field, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);}
1471:   return(0);
1472: }

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

1477:   Not collective

1479:   Input Parameters:
1480: + fem          - The PetscFE object for the field being integrated
1481: . prob         - The PetscDS specifying the discretizations and continuum functions
1482: . jtype        - The type of matrix pointwise functions that should be used
1483: . fieldI       - The test field being integrated
1484: . fieldJ       - The basis field being integrated
1485: . Ne           - The number of elements in the chunk
1486: . cgeom        - The cell geometry for each cell in the chunk
1487: . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1488: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1489: . probAux      - The PetscDS specifying the auxiliary discretizations
1490: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1491: . t            - The time
1492: - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)

1494:   Output Parameter
1495: . elemMat      - the element matrices for the Jacobian from each element

1497:   Note:
1498: $ Loop over batch of elements (e):
1499: $   Loop over element matrix entries (f,fc,g,gc --> i,j):
1500: $     Loop over quadrature points (q):
1501: $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1502: $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1503: $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1504: $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1505: $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1506:   Level: intermediate

1508: .seealso: PetscFEIntegrateResidual()
1509: @*/
1510: PetscErrorCode PetscFEIntegrateJacobian(PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom,
1511:                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1512: {
1513:   PetscFE        fe;

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

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

1526:   Not collective

1528:   Input Parameters:
1529: . prob         - The PetscDS specifying the discretizations and continuum functions
1530: . fieldI       - The test field being integrated
1531: . fieldJ       - The basis field being integrated
1532: . Ne           - The number of elements in the chunk
1533: . fgeom        - The face geometry for each cell in the chunk
1534: . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1535: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1536: . probAux      - The PetscDS specifying the auxiliary discretizations
1537: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1538: . t            - The time
1539: - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)

1541:   Output Parameter
1542: . elemMat              - the element matrices for the Jacobian from each element

1544:   Note:
1545: $ Loop over batch of elements (e):
1546: $   Loop over element matrix entries (f,fc,g,gc --> i,j):
1547: $     Loop over quadrature points (q):
1548: $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1549: $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1550: $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1551: $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1552: $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1553:   Level: intermediate

1555: .seealso: PetscFEIntegrateJacobian(), PetscFEIntegrateResidual()
1556: @*/
1557: PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom,
1558:                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1559: {
1560:   PetscFE        fe;

1565:   PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
1566:   if (fe->ops->integratebdjacobian) {(*fe->ops->integratebdjacobian)(prob, fieldI, fieldJ, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);}
1567:   return(0);
1568: }

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

1573:   Input Parameters:
1574: + fe     - The finite element space
1575: - height - The height of the Plex point

1577:   Output Parameter:
1578: . subfe  - The subspace of this FE space

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

1582:   Level: advanced

1584: .seealso: PetscFECreateDefault()
1585: @*/
1586: PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe)
1587: {
1588:   PetscSpace      P, subP;
1589:   PetscDualSpace  Q, subQ;
1590:   PetscQuadrature subq;
1591:   PetscFEType     fetype;
1592:   PetscInt        dim, Nc;
1593:   PetscErrorCode  ierr;

1598:   if (height == 0) {
1599:     *subfe = fe;
1600:     return(0);
1601:   }
1602:   PetscFEGetBasisSpace(fe, &P);
1603:   PetscFEGetDualSpace(fe, &Q);
1604:   PetscFEGetNumComponents(fe, &Nc);
1605:   PetscFEGetFaceQuadrature(fe, &subq);
1606:   PetscDualSpaceGetDimension(Q, &dim);
1607:   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);}
1608:   if (!fe->subspaces) {PetscCalloc1(dim, &fe->subspaces);}
1609:   if (height <= dim) {
1610:     if (!fe->subspaces[height-1]) {
1611:       PetscFE     sub;
1612:       const char *name;

1614:       PetscSpaceGetHeightSubspace(P, height, &subP);
1615:       PetscDualSpaceGetHeightSubspace(Q, height, &subQ);
1616:       PetscFECreate(PetscObjectComm((PetscObject) fe), &sub);
1617:       PetscObjectGetName((PetscObject) fe,  &name);
1618:       PetscObjectSetName((PetscObject) sub,  name);
1619:       PetscFEGetType(fe, &fetype);
1620:       PetscFESetType(sub, fetype);
1621:       PetscFESetBasisSpace(sub, subP);
1622:       PetscFESetDualSpace(sub, subQ);
1623:       PetscFESetNumComponents(sub, Nc);
1624:       PetscFESetUp(sub);
1625:       PetscFESetQuadrature(sub, subq);
1626:       fe->subspaces[height-1] = sub;
1627:     }
1628:     *subfe = fe->subspaces[height-1];
1629:   } else {
1630:     *subfe = NULL;
1631:   }
1632:   return(0);
1633: }

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

1640:   Collective on fem

1642:   Input Parameter:
1643: . fe - The initial PetscFE

1645:   Output Parameter:
1646: . feRef - The refined PetscFE

1648:   Level: advanced

1650: .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
1651: @*/
1652: PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)
1653: {
1654:   PetscSpace       P, Pref;
1655:   PetscDualSpace   Q, Qref;
1656:   DM               K, Kref;
1657:   PetscQuadrature  q, qref;
1658:   const PetscReal *v0, *jac;
1659:   PetscInt         numComp, numSubelements;
1660:   PetscErrorCode   ierr;

1663:   PetscFEGetBasisSpace(fe, &P);
1664:   PetscFEGetDualSpace(fe, &Q);
1665:   PetscFEGetQuadrature(fe, &q);
1666:   PetscDualSpaceGetDM(Q, &K);
1667:   /* Create space */
1668:   PetscObjectReference((PetscObject) P);
1669:   Pref = P;
1670:   /* Create dual space */
1671:   PetscDualSpaceDuplicate(Q, &Qref);
1672:   DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref);
1673:   PetscDualSpaceSetDM(Qref, Kref);
1674:   DMDestroy(&Kref);
1675:   PetscDualSpaceSetUp(Qref);
1676:   /* Create element */
1677:   PetscFECreate(PetscObjectComm((PetscObject) fe), feRef);
1678:   PetscFESetType(*feRef, PETSCFECOMPOSITE);
1679:   PetscFESetBasisSpace(*feRef, Pref);
1680:   PetscFESetDualSpace(*feRef, Qref);
1681:   PetscFEGetNumComponents(fe,    &numComp);
1682:   PetscFESetNumComponents(*feRef, numComp);
1683:   PetscFESetUp(*feRef);
1684:   PetscSpaceDestroy(&Pref);
1685:   PetscDualSpaceDestroy(&Qref);
1686:   /* Create quadrature */
1687:   PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL);
1688:   PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);
1689:   PetscFESetQuadrature(*feRef, qref);
1690:   PetscQuadratureDestroy(&qref);
1691:   return(0);
1692: }

1694: /*@C
1695:   PetscFECreateDefault - Create a PetscFE for basic FEM computation

1697:   Collective

1699:   Input Parameters:
1700: + comm      - The MPI comm
1701: . dim       - The spatial dimension
1702: . Nc        - The number of components
1703: . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
1704: . prefix    - The options prefix, or NULL
1705: - qorder    - The quadrature order or PETSC_DETERMINE to use PetscSpace polynomial degree

1707:   Output Parameter:
1708: . fem - The PetscFE object

1710:   Level: beginner

1712: .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
1713: @*/
1714: PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem)
1715: {
1716:   PetscQuadrature q, fq;
1717:   DM              K;
1718:   PetscSpace      P;
1719:   PetscDualSpace  Q;
1720:   PetscInt        order, quadPointsPerEdge;
1721:   PetscBool       tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
1722:   PetscErrorCode  ierr;

1725:   /* Create space */
1726:   PetscSpaceCreate(comm, &P);
1727:   PetscObjectSetOptionsPrefix((PetscObject) P, prefix);
1728:   PetscSpacePolynomialSetTensor(P, tensor);
1729:   PetscSpaceSetNumComponents(P, Nc);
1730:   PetscSpaceSetNumVariables(P, dim);
1731:   PetscSpaceSetFromOptions(P);
1732:   PetscSpaceSetUp(P);
1733:   PetscSpaceGetDegree(P, &order, NULL);
1734:   PetscSpacePolynomialGetTensor(P, &tensor);
1735:   /* Create dual space */
1736:   PetscDualSpaceCreate(comm, &Q);
1737:   PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);
1738:   PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);
1739:   PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);
1740:   PetscDualSpaceSetDM(Q, K);
1741:   DMDestroy(&K);
1742:   PetscDualSpaceSetNumComponents(Q, Nc);
1743:   PetscDualSpaceSetOrder(Q, order);
1744:   PetscDualSpaceLagrangeSetTensor(Q, tensor);
1745:   PetscDualSpaceSetFromOptions(Q);
1746:   PetscDualSpaceSetUp(Q);
1747:   /* Create element */
1748:   PetscFECreate(comm, fem);
1749:   PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);
1750:   PetscFESetBasisSpace(*fem, P);
1751:   PetscFESetDualSpace(*fem, Q);
1752:   PetscFESetNumComponents(*fem, Nc);
1753:   PetscFESetFromOptions(*fem);
1754:   PetscFESetUp(*fem);
1755:   PetscSpaceDestroy(&P);
1756:   PetscDualSpaceDestroy(&Q);
1757:   /* Create quadrature (with specified order if given) */
1758:   qorder = qorder >= 0 ? qorder : order;
1759:   PetscObjectOptionsBegin((PetscObject)*fem);
1760:   PetscOptionsBoundedInt("-petscfe_default_quadrature_order","Quadrature order is one less than quadrature points per edge","PetscFECreateDefault",qorder,&qorder,NULL,0);
1761:   PetscOptionsEnd();
1762:   quadPointsPerEdge = PetscMax(qorder + 1,1);
1763:   if (isSimplex) {
1764:     PetscDTGaussJacobiQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);
1765:     PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);
1766:   } else {
1767:     PetscDTGaussTensorQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);
1768:     PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);
1769:   }
1770:   PetscFESetQuadrature(*fem, q);
1771:   PetscFESetFaceQuadrature(*fem, fq);
1772:   PetscQuadratureDestroy(&q);
1773:   PetscQuadratureDestroy(&fq);
1774:   return(0);
1775: }

1777: /*@C
1778:   PetscFESetName - Names the FE and its subobjects

1780:   Not collective

1782:   Input Parameters:
1783: + fe   - The PetscFE
1784: - name - The name

1786:   Level: intermediate

1788: .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
1789: @*/
1790: PetscErrorCode PetscFESetName(PetscFE fe, const char name[])
1791: {
1792:   PetscSpace     P;
1793:   PetscDualSpace Q;

1797:   PetscFEGetBasisSpace(fe, &P);
1798:   PetscFEGetDualSpace(fe, &Q);
1799:   PetscObjectSetName((PetscObject) fe, name);
1800:   PetscObjectSetName((PetscObject) P,  name);
1801:   PetscObjectSetName((PetscObject) Q,  name);
1802:   return(0);
1803: }

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

1810:   for (f = 0; f < Nf; ++f) {
1811:     PetscFE          fe;
1812:     const PetscInt   cdim = T[f]->cdim;
1813:     const PetscInt   Nq   = T[f]->Np;
1814:     const PetscInt   Nbf  = T[f]->Nb;
1815:     const PetscInt   Ncf  = T[f]->Nc;
1816:     const PetscReal *Bq   = &T[f]->T[0][(r*Nq+q)*Nbf*Ncf];
1817:     const PetscReal *Dq   = &T[f]->T[1][(r*Nq+q)*Nbf*Ncf*cdim];
1818:     PetscInt         b, c, d;

1820:     PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);
1821:     for (c = 0; c < Ncf; ++c) u[fOffset+c] = 0.0;
1822:     for (d = 0; d < cdim*Ncf; ++d) u_x[fOffset*cdim+d] = 0.0;
1823:     for (b = 0; b < Nbf; ++b) {
1824:       for (c = 0; c < Ncf; ++c) {
1825:         const PetscInt cidx = b*Ncf+c;

1827:         u[fOffset+c] += Bq[cidx]*coefficients[dOffset+b];
1828:         for (d = 0; d < cdim; ++d) u_x[(fOffset+c)*cdim+d] += Dq[cidx*cdim+d]*coefficients[dOffset+b];
1829:       }
1830:     }
1831:     PetscFEPushforward(fe, fegeom, 1, &u[fOffset]);
1832:     PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset*cdim]);
1833:     if (u_t) {
1834:       for (c = 0; c < Ncf; ++c) u_t[fOffset+c] = 0.0;
1835:       for (b = 0; b < Nbf; ++b) {
1836:         for (c = 0; c < Ncf; ++c) {
1837:           const PetscInt cidx = b*Ncf+c;

1839:           u_t[fOffset+c] += Bq[cidx]*coefficients_t[dOffset+b];
1840:         }
1841:       }
1842:       PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]);
1843:     }
1844:     fOffset += Ncf;
1845:     dOffset += Nbf;
1846:   }
1847:   return 0;
1848: }

1850: PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[])
1851: {
1852:   PetscFE           fe;
1853:   PetscTabulation Tc;
1854:   PetscInt          b, c;
1855:   PetscErrorCode    ierr;

1857:   if (!prob) return 0;
1858:   PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1859:   PetscFEGetFaceCentroidTabulation(fe, &Tc);
1860:   {
1861:     const PetscReal *faceBasis = Tc->T[0];
1862:     const PetscInt   Nb        = Tc->Nb;
1863:     const PetscInt   Nc        = Tc->Nc;

1865:     for (c = 0; c < Nc; ++c) {u[c] = 0.0;}
1866:     for (b = 0; b < Nb; ++b) {
1867:       for (c = 0; c < Nc; ++c) {
1868:         const PetscInt cidx = b*Nc+c;

1870:         u[c] += coefficients[cidx]*faceBasis[faceLoc*Nb*Nc+cidx];
1871:       }
1872:     }
1873:   }
1874:   return 0;
1875: }

1877: PetscErrorCode PetscFEUpdateElementVec_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[])
1878: {
1879:   const PetscInt   dim      = T->cdim;
1880:   const PetscInt   Nq       = T->Np;
1881:   const PetscInt   Nb       = T->Nb;
1882:   const PetscInt   Nc       = T->Nc;
1883:   const PetscReal *basis    = &T->T[0][r*Nq*Nb*Nc];
1884:   const PetscReal *basisDer = &T->T[1][r*Nq*Nb*Nc*dim];
1885:   PetscInt         q, b, c, d;
1886:   PetscErrorCode   ierr;

1888:   for (b = 0; b < Nb; ++b) elemVec[b] = 0.0;
1889:   for (q = 0; q < Nq; ++q) {
1890:     for (b = 0; b < Nb; ++b) {
1891:       for (c = 0; c < Nc; ++c) {
1892:         const PetscInt bcidx = b*Nc+c;

1894:         tmpBasis[bcidx] = basis[q*Nb*Nc+bcidx];
1895:         for (d = 0; d < dim; ++d) tmpBasisDer[bcidx*dim+d] = basisDer[q*Nb*Nc*dim+bcidx*dim+d];
1896:       }
1897:     }
1898:     PetscFEPushforward(fe, fegeom, Nb, tmpBasis);
1899:     PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer);
1900:     for (b = 0; b < Nb; ++b) {
1901:       for (c = 0; c < Nc; ++c) {
1902:         const PetscInt bcidx = b*Nc+c;
1903:         const PetscInt qcidx = q*Nc+c;

1905:         elemVec[b] += tmpBasis[bcidx]*f0[qcidx];
1906:         for (d = 0; d < dim; ++d) elemVec[b] += tmpBasisDer[bcidx*dim+d]*f1[qcidx*dim+d];
1907:       }
1908:     }
1909:   }
1910:   return(0);
1911: }

1913: PetscErrorCode PetscFEUpdateElementMat_Internal(PetscFE feI, PetscFE feJ, PetscInt r, PetscInt q, PetscTabulation TI, PetscScalar tmpBasisI[], PetscScalar tmpBasisDerI[], PetscTabulation TJ, PetscScalar tmpBasisJ[], PetscScalar tmpBasisDerJ[], PetscFEGeom *fegeom, const PetscScalar g0[], const PetscScalar g1[], const PetscScalar g2[], const PetscScalar g3[], PetscInt eOffset, PetscInt totDim, PetscInt offsetI, PetscInt offsetJ, PetscScalar elemMat[])
1914: {
1915:   const PetscInt   dim       = TI->cdim;
1916:   const PetscInt   NqI       = TI->Np;
1917:   const PetscInt   NbI       = TI->Nb;
1918:   const PetscInt   NcI       = TI->Nc;
1919:   const PetscReal *basisI    = &TI->T[0][(r*NqI+q)*NbI*NcI];
1920:   const PetscReal *basisDerI = &TI->T[1][(r*NqI+q)*NbI*NcI*dim];
1921:   const PetscInt   NqJ       = TJ->Np;
1922:   const PetscInt   NbJ       = TJ->Nb;
1923:   const PetscInt   NcJ       = TJ->Nc;
1924:   const PetscReal *basisJ    = &TJ->T[0][(r*NqJ+q)*NbJ*NcJ];
1925:   const PetscReal *basisDerJ = &TJ->T[1][(r*NqJ+q)*NbJ*NcJ*dim];
1926:   PetscInt         f, fc, g, gc, df, dg;
1927:   PetscErrorCode   ierr;

1929:   for (f = 0; f < NbI; ++f) {
1930:     for (fc = 0; fc < NcI; ++fc) {
1931:       const PetscInt fidx = f*NcI+fc; /* Test function basis index */

1933:       tmpBasisI[fidx] = basisI[fidx];
1934:       for (df = 0; df < dim; ++df) tmpBasisDerI[fidx*dim+df] = basisDerI[fidx*dim+df];
1935:     }
1936:   }
1937:   PetscFEPushforward(feI, fegeom, NbI, tmpBasisI);
1938:   PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI);
1939:   for (g = 0; g < NbJ; ++g) {
1940:     for (gc = 0; gc < NcJ; ++gc) {
1941:       const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */

1943:       tmpBasisJ[gidx] = basisJ[gidx];
1944:       for (dg = 0; dg < dim; ++dg) tmpBasisDerJ[gidx*dim+dg] = basisDerJ[gidx*dim+dg];
1945:     }
1946:   }
1947:   PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ);
1948:   PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ);
1949:   for (f = 0; f < NbI; ++f) {
1950:     for (fc = 0; fc < NcI; ++fc) {
1951:       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
1952:       const PetscInt i    = offsetI+f; /* Element matrix row */
1953:       for (g = 0; g < NbJ; ++g) {
1954:         for (gc = 0; gc < NcJ; ++gc) {
1955:           const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
1956:           const PetscInt j    = offsetJ+g; /* Element matrix column */
1957:           const PetscInt fOff = eOffset+i*totDim+j;

1959:           elemMat[fOff] += tmpBasisI[fidx]*g0[fc*NcJ+gc]*tmpBasisJ[gidx];
1960:           for (df = 0; df < dim; ++df) {
1961:             elemMat[fOff] += tmpBasisI[fidx]*g1[(fc*NcJ+gc)*dim+df]*tmpBasisDerJ[gidx*dim+df];
1962:             elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g2[(fc*NcJ+gc)*dim+df]*tmpBasisJ[gidx];
1963:             for (dg = 0; dg < dim; ++dg) {
1964:               elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g3[((fc*NcJ+gc)*dim+df)*dim+dg]*tmpBasisDerJ[gidx*dim+dg];
1965:             }
1966:           }
1967:         }
1968:       }
1969:     }
1970:   }
1971:   return(0);
1972: }

1974: PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom)
1975: {
1976:   PetscDualSpace  dsp;
1977:   DM              dm;
1978:   PetscQuadrature quadDef;
1979:   PetscInt        dim, cdim, Nq;
1980:   PetscErrorCode  ierr;

1983:   PetscFEGetDualSpace(fe, &dsp);
1984:   PetscDualSpaceGetDM(dsp, &dm);
1985:   DMGetDimension(dm, &dim);
1986:   DMGetCoordinateDim(dm, &cdim);
1987:   PetscFEGetQuadrature(fe, &quadDef);
1988:   quad = quad ? quad : quadDef;
1989:   PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);
1990:   PetscMalloc1(Nq*cdim,      &cgeom->v);
1991:   PetscMalloc1(Nq*cdim*cdim, &cgeom->J);
1992:   PetscMalloc1(Nq*cdim*cdim, &cgeom->invJ);
1993:   PetscMalloc1(Nq,           &cgeom->detJ);
1994:   cgeom->dim       = dim;
1995:   cgeom->dimEmbed  = cdim;
1996:   cgeom->numCells  = 1;
1997:   cgeom->numPoints = Nq;
1998:   DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ);
1999:   return(0);
2000: }

2002: PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom)
2003: {

2007:   PetscFree(cgeom->v);
2008:   PetscFree(cgeom->J);
2009:   PetscFree(cgeom->invJ);
2010:   PetscFree(cgeom->detJ);
2011:   return(0);
2012: }