xref: /petsc/src/dm/dt/fe/interface/fe.c (revision 60dd46ca9875962950ab5c8a8f4d94fdb4c0f129)
1 /* Basis Jet Tabulation
2 
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.
7 
8   \psi_i = \sum_k \alpha_{ki} \phi_k
9 
10 Our nodal basis is defined in terms of the dual basis $n_j$
11 
12   n_j \cdot \psi_i = \delta_{ji}
13 
14 and we may act on the first equation to obtain
15 
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
19 
20 so the coefficients of the nodal basis in the prime basis are
21 
22    \alpha = V^{-1}
23 
24 We will define the dual basis vectors $n_j$ using a quadrature rule.
25 
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$.
29 
30 I will have to count the degrees correctly for the Legendre product when we are on simplices.
31 
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> /*I "petscfe.h" I*/
38 #include <petscdmplex.h>
39 
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";
50 
51 PetscClassId PETSCFE_CLASSID = 0;
52 
53 PetscFunctionList PetscFEList              = NULL;
54 PetscBool         PetscFERegisterAllCalled = PETSC_FALSE;
55 
56 /*@C
57   PetscFERegister - Adds a new PetscFE implementation
58 
59   Not Collective
60 
61   Input Parameters:
62 + name        - The name of a new user-defined creation routine
63 - create_func - The creation routine itself
64 
65   Notes:
66   PetscFERegister() may be called multiple times to add several user-defined PetscFEs
67 
68   Sample usage:
69 .vb
70     PetscFERegister("my_fe", MyPetscFECreate);
71 .ve
72 
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
82 
83   Level: advanced
84 
85 .seealso: PetscFERegisterAll(), PetscFERegisterDestroy()
86 
87 @*/
88 PetscErrorCode PetscFERegister(const char sname[], PetscErrorCode (*function)(PetscFE))
89 {
90   PetscErrorCode ierr;
91 
92   PetscFunctionBegin;
93   ierr = PetscFunctionListAdd(&PetscFEList, sname, function);CHKERRQ(ierr);
94   PetscFunctionReturn(0);
95 }
96 
97 /*@C
98   PetscFESetType - Builds a particular PetscFE
99 
100   Collective on fem
101 
102   Input Parameters:
103 + fem  - The PetscFE object
104 - name - The kind of FEM space
105 
106   Options Database Key:
107 . -petscfe_type <type> - Sets the PetscFE type; use -help for a list of available types
108 
109   Level: intermediate
110 
111 .seealso: PetscFEGetType(), PetscFECreate()
112 @*/
113 PetscErrorCode PetscFESetType(PetscFE fem, PetscFEType name)
114 {
115   PetscErrorCode (*r)(PetscFE);
116   PetscBool      match;
117   PetscErrorCode ierr;
118 
119   PetscFunctionBegin;
120   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
121   ierr = PetscObjectTypeCompare((PetscObject) fem, name, &match);CHKERRQ(ierr);
122   if (match) PetscFunctionReturn(0);
123 
124   if (!PetscFERegisterAllCalled) {ierr = PetscFERegisterAll();CHKERRQ(ierr);}
125   ierr = PetscFunctionListFind(PetscFEList, name, &r);CHKERRQ(ierr);
126   if (!r) SETERRQ1(PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFE type: %s", name);
127 
128   if (fem->ops->destroy) {
129     ierr              = (*fem->ops->destroy)(fem);CHKERRQ(ierr);
130     fem->ops->destroy = NULL;
131   }
132   ierr = (*r)(fem);CHKERRQ(ierr);
133   ierr = PetscObjectChangeTypeName((PetscObject) fem, name);CHKERRQ(ierr);
134   PetscFunctionReturn(0);
135 }
136 
137 /*@C
138   PetscFEGetType - Gets the PetscFE type name (as a string) from the object.
139 
140   Not Collective
141 
142   Input Parameter:
143 . fem  - The PetscFE
144 
145   Output Parameter:
146 . name - The PetscFE type name
147 
148   Level: intermediate
149 
150 .seealso: PetscFESetType(), PetscFECreate()
151 @*/
152 PetscErrorCode PetscFEGetType(PetscFE fem, PetscFEType *name)
153 {
154   PetscErrorCode ierr;
155 
156   PetscFunctionBegin;
157   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
158   PetscValidPointer(name, 2);
159   if (!PetscFERegisterAllCalled) {
160     ierr = PetscFERegisterAll();CHKERRQ(ierr);
161   }
162   *name = ((PetscObject) fem)->type_name;
163   PetscFunctionReturn(0);
164 }
165 
166 /*@C
167    PetscFEViewFromOptions - View from Options
168 
169    Collective on PetscFE
170 
171    Input Parameters:
172 +  A - the PetscFE object
173 .  obj - Optional object
174 -  name - command line option
175 
176    Level: intermediate
177 .seealso:  PetscFE(), PetscFEView(), PetscObjectViewFromOptions(), PetscFECreate()
178 @*/
179 PetscErrorCode  PetscFEViewFromOptions(PetscFE A,PetscObject obj,const char name[])
180 {
181   PetscErrorCode ierr;
182 
183   PetscFunctionBegin;
184   PetscValidHeaderSpecific(A,PETSCFE_CLASSID,1);
185   ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr);
186   PetscFunctionReturn(0);
187 }
188 
189 /*@C
190   PetscFEView - Views a PetscFE
191 
192   Collective on fem
193 
194   Input Parameter:
195 + fem - the PetscFE object to view
196 - viewer   - the viewer
197 
198   Level: beginner
199 
200 .seealso PetscFEDestroy()
201 @*/
202 PetscErrorCode PetscFEView(PetscFE fem, PetscViewer viewer)
203 {
204   PetscBool      iascii;
205   PetscErrorCode ierr;
206 
207   PetscFunctionBegin;
208   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
209   if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
210   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &viewer);CHKERRQ(ierr);}
211   ierr = PetscObjectPrintClassNamePrefixType((PetscObject)fem, viewer);CHKERRQ(ierr);
212   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
213   if (fem->ops->view) {ierr = (*fem->ops->view)(fem, viewer);CHKERRQ(ierr);}
214   PetscFunctionReturn(0);
215 }
216 
217 /*@
218   PetscFESetFromOptions - sets parameters in a PetscFE from the options database
219 
220   Collective on fem
221 
222   Input Parameter:
223 . fem - the PetscFE object to set options for
224 
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
228 
229   Level: intermediate
230 
231 .seealso PetscFEView()
232 @*/
233 PetscErrorCode PetscFESetFromOptions(PetscFE fem)
234 {
235   const char    *defaultType;
236   char           name[256];
237   PetscBool      flg;
238   PetscErrorCode ierr;
239 
240   PetscFunctionBegin;
241   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
242   if (!((PetscObject) fem)->type_name) {
243     defaultType = PETSCFEBASIC;
244   } else {
245     defaultType = ((PetscObject) fem)->type_name;
246   }
247   if (!PetscFERegisterAllCalled) {ierr = PetscFERegisterAll();CHKERRQ(ierr);}
248 
249   ierr = PetscObjectOptionsBegin((PetscObject) fem);CHKERRQ(ierr);
250   ierr = PetscOptionsFList("-petscfe_type", "Finite element space", "PetscFESetType", PetscFEList, defaultType, name, 256, &flg);CHKERRQ(ierr);
251   if (flg) {
252     ierr = PetscFESetType(fem, name);CHKERRQ(ierr);
253   } else if (!((PetscObject) fem)->type_name) {
254     ierr = PetscFESetType(fem, defaultType);CHKERRQ(ierr);
255   }
256   ierr = PetscOptionsBoundedInt("-petscfe_num_blocks", "The number of cell blocks to integrate concurrently", "PetscSpaceSetTileSizes", fem->numBlocks, &fem->numBlocks, NULL,1);CHKERRQ(ierr);
257   ierr = PetscOptionsBoundedInt("-petscfe_num_batches", "The number of cell batches to integrate serially", "PetscSpaceSetTileSizes", fem->numBatches, &fem->numBatches, NULL,1);CHKERRQ(ierr);
258   if (fem->ops->setfromoptions) {
259     ierr = (*fem->ops->setfromoptions)(PetscOptionsObject,fem);CHKERRQ(ierr);
260   }
261   /* process any options handlers added with PetscObjectAddOptionsHandler() */
262   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) fem);CHKERRQ(ierr);
263   ierr = PetscOptionsEnd();CHKERRQ(ierr);
264   ierr = PetscFEViewFromOptions(fem, NULL, "-petscfe_view");CHKERRQ(ierr);
265   PetscFunctionReturn(0);
266 }
267 
268 /*@C
269   PetscFESetUp - Construct data structures for the PetscFE
270 
271   Collective on fem
272 
273   Input Parameter:
274 . fem - the PetscFE object to setup
275 
276   Level: intermediate
277 
278 .seealso PetscFEView(), PetscFEDestroy()
279 @*/
280 PetscErrorCode PetscFESetUp(PetscFE fem)
281 {
282   PetscErrorCode ierr;
283 
284   PetscFunctionBegin;
285   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
286   if (fem->setupcalled) PetscFunctionReturn(0);
287   fem->setupcalled = PETSC_TRUE;
288   if (fem->ops->setup) {ierr = (*fem->ops->setup)(fem);CHKERRQ(ierr);}
289   PetscFunctionReturn(0);
290 }
291 
292 /*@
293   PetscFEDestroy - Destroys a PetscFE object
294 
295   Collective on fem
296 
297   Input Parameter:
298 . fem - the PetscFE object to destroy
299 
300   Level: beginner
301 
302 .seealso PetscFEView()
303 @*/
304 PetscErrorCode PetscFEDestroy(PetscFE *fem)
305 {
306   PetscErrorCode ierr;
307 
308   PetscFunctionBegin;
309   if (!*fem) PetscFunctionReturn(0);
310   PetscValidHeaderSpecific((*fem), PETSCFE_CLASSID, 1);
311 
312   if (--((PetscObject)(*fem))->refct > 0) {*fem = 0; PetscFunctionReturn(0);}
313   ((PetscObject) (*fem))->refct = 0;
314 
315   if ((*fem)->subspaces) {
316     PetscInt dim, d;
317 
318     ierr = PetscDualSpaceGetDimension((*fem)->dualSpace, &dim);CHKERRQ(ierr);
319     for (d = 0; d < dim; ++d) {ierr = PetscFEDestroy(&(*fem)->subspaces[d]);CHKERRQ(ierr);}
320   }
321   ierr = PetscFree((*fem)->subspaces);CHKERRQ(ierr);
322   ierr = PetscFree((*fem)->invV);CHKERRQ(ierr);
323   ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->B, &(*fem)->D, NULL /*&(*fem)->H*/);CHKERRQ(ierr);
324   ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->Bf, &(*fem)->Df, NULL /*&(*fem)->Hf*/);CHKERRQ(ierr);
325   ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->F, NULL, NULL);CHKERRQ(ierr);
326   ierr = PetscSpaceDestroy(&(*fem)->basisSpace);CHKERRQ(ierr);
327   ierr = PetscDualSpaceDestroy(&(*fem)->dualSpace);CHKERRQ(ierr);
328   ierr = PetscQuadratureDestroy(&(*fem)->quadrature);CHKERRQ(ierr);
329   ierr = PetscQuadratureDestroy(&(*fem)->faceQuadrature);CHKERRQ(ierr);
330 
331   if ((*fem)->ops->destroy) {ierr = (*(*fem)->ops->destroy)(*fem);CHKERRQ(ierr);}
332   ierr = PetscHeaderDestroy(fem);CHKERRQ(ierr);
333   PetscFunctionReturn(0);
334 }
335 
336 /*@
337   PetscFECreate - Creates an empty PetscFE object. The type can then be set with PetscFESetType().
338 
339   Collective
340 
341   Input Parameter:
342 . comm - The communicator for the PetscFE object
343 
344   Output Parameter:
345 . fem - The PetscFE object
346 
347   Level: beginner
348 
349 .seealso: PetscFESetType(), PETSCFEGALERKIN
350 @*/
351 PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem)
352 {
353   PetscFE        f;
354   PetscErrorCode ierr;
355 
356   PetscFunctionBegin;
357   PetscValidPointer(fem, 2);
358   ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr);
359   *fem = NULL;
360   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
361 
362   ierr = PetscHeaderCreate(f, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView);CHKERRQ(ierr);
363 
364   f->basisSpace    = NULL;
365   f->dualSpace     = NULL;
366   f->numComponents = 1;
367   f->subspaces     = NULL;
368   f->invV          = NULL;
369   f->B             = NULL;
370   f->D             = NULL;
371   f->H             = NULL;
372   f->Bf            = NULL;
373   f->Df            = NULL;
374   f->Hf            = NULL;
375   ierr = PetscArrayzero(&f->quadrature, 1);CHKERRQ(ierr);
376   ierr = PetscArrayzero(&f->faceQuadrature, 1);CHKERRQ(ierr);
377   f->blockSize     = 0;
378   f->numBlocks     = 1;
379   f->batchSize     = 0;
380   f->numBatches    = 1;
381 
382   *fem = f;
383   PetscFunctionReturn(0);
384 }
385 
386 /*@
387   PetscFEGetSpatialDimension - Returns the spatial dimension of the element
388 
389   Not collective
390 
391   Input Parameter:
392 . fem - The PetscFE object
393 
394   Output Parameter:
395 . dim - The spatial dimension
396 
397   Level: intermediate
398 
399 .seealso: PetscFECreate()
400 @*/
401 PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim)
402 {
403   DM             dm;
404   PetscErrorCode ierr;
405 
406   PetscFunctionBegin;
407   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
408   PetscValidPointer(dim, 2);
409   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
410   ierr = DMGetDimension(dm, dim);CHKERRQ(ierr);
411   PetscFunctionReturn(0);
412 }
413 
414 /*@
415   PetscFESetNumComponents - Sets the number of components in the element
416 
417   Not collective
418 
419   Input Parameters:
420 + fem - The PetscFE object
421 - comp - The number of field components
422 
423   Level: intermediate
424 
425 .seealso: PetscFECreate()
426 @*/
427 PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp)
428 {
429   PetscFunctionBegin;
430   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
431   fem->numComponents = comp;
432   PetscFunctionReturn(0);
433 }
434 
435 /*@
436   PetscFEGetNumComponents - Returns the number of components in the element
437 
438   Not collective
439 
440   Input Parameter:
441 . fem - The PetscFE object
442 
443   Output Parameter:
444 . comp - The number of field components
445 
446   Level: intermediate
447 
448 .seealso: PetscFECreate()
449 @*/
450 PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp)
451 {
452   PetscFunctionBegin;
453   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
454   PetscValidPointer(comp, 2);
455   *comp = fem->numComponents;
456   PetscFunctionReturn(0);
457 }
458 
459 /*@
460   PetscFESetTileSizes - Sets the tile sizes for evaluation
461 
462   Not collective
463 
464   Input Parameters:
465 + fem - The PetscFE object
466 . blockSize - The number of elements in a block
467 . numBlocks - The number of blocks in a batch
468 . batchSize - The number of elements in a batch
469 - numBatches - The number of batches in a chunk
470 
471   Level: intermediate
472 
473 .seealso: PetscFECreate()
474 @*/
475 PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches)
476 {
477   PetscFunctionBegin;
478   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
479   fem->blockSize  = blockSize;
480   fem->numBlocks  = numBlocks;
481   fem->batchSize  = batchSize;
482   fem->numBatches = numBatches;
483   PetscFunctionReturn(0);
484 }
485 
486 /*@
487   PetscFEGetTileSizes - Returns the tile sizes for evaluation
488 
489   Not collective
490 
491   Input Parameter:
492 . fem - The PetscFE object
493 
494   Output Parameters:
495 + blockSize - The number of elements in a block
496 . numBlocks - The number of blocks in a batch
497 . batchSize - The number of elements in a batch
498 - numBatches - The number of batches in a chunk
499 
500   Level: intermediate
501 
502 .seealso: PetscFECreate()
503 @*/
504 PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PetscInt *blockSize, PetscInt *numBlocks, PetscInt *batchSize, PetscInt *numBatches)
505 {
506   PetscFunctionBegin;
507   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
508   if (blockSize)  PetscValidPointer(blockSize,  2);
509   if (numBlocks)  PetscValidPointer(numBlocks,  3);
510   if (batchSize)  PetscValidPointer(batchSize,  4);
511   if (numBatches) PetscValidPointer(numBatches, 5);
512   if (blockSize)  *blockSize  = fem->blockSize;
513   if (numBlocks)  *numBlocks  = fem->numBlocks;
514   if (batchSize)  *batchSize  = fem->batchSize;
515   if (numBatches) *numBatches = fem->numBatches;
516   PetscFunctionReturn(0);
517 }
518 
519 /*@
520   PetscFEGetBasisSpace - Returns the PetscSpace used for approximation of the solution
521 
522   Not collective
523 
524   Input Parameter:
525 . fem - The PetscFE object
526 
527   Output Parameter:
528 . sp - The PetscSpace object
529 
530   Level: intermediate
531 
532 .seealso: PetscFECreate()
533 @*/
534 PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp)
535 {
536   PetscFunctionBegin;
537   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
538   PetscValidPointer(sp, 2);
539   *sp = fem->basisSpace;
540   PetscFunctionReturn(0);
541 }
542 
543 /*@
544   PetscFESetBasisSpace - Sets the PetscSpace used for approximation of the solution
545 
546   Not collective
547 
548   Input Parameters:
549 + fem - The PetscFE object
550 - sp - The PetscSpace object
551 
552   Level: intermediate
553 
554 .seealso: PetscFECreate()
555 @*/
556 PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp)
557 {
558   PetscErrorCode ierr;
559 
560   PetscFunctionBegin;
561   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
562   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2);
563   ierr = PetscSpaceDestroy(&fem->basisSpace);CHKERRQ(ierr);
564   fem->basisSpace = sp;
565   ierr = PetscObjectReference((PetscObject) fem->basisSpace);CHKERRQ(ierr);
566   PetscFunctionReturn(0);
567 }
568 
569 /*@
570   PetscFEGetDualSpace - Returns the PetscDualSpace used to define the inner product
571 
572   Not collective
573 
574   Input Parameter:
575 . fem - The PetscFE object
576 
577   Output Parameter:
578 . sp - The PetscDualSpace object
579 
580   Level: intermediate
581 
582 .seealso: PetscFECreate()
583 @*/
584 PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp)
585 {
586   PetscFunctionBegin;
587   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
588   PetscValidPointer(sp, 2);
589   *sp = fem->dualSpace;
590   PetscFunctionReturn(0);
591 }
592 
593 /*@
594   PetscFESetDualSpace - Sets the PetscDualSpace used to define the inner product
595 
596   Not collective
597 
598   Input Parameters:
599 + fem - The PetscFE object
600 - sp - The PetscDualSpace object
601 
602   Level: intermediate
603 
604 .seealso: PetscFECreate()
605 @*/
606 PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp)
607 {
608   PetscErrorCode ierr;
609 
610   PetscFunctionBegin;
611   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
612   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2);
613   ierr = PetscDualSpaceDestroy(&fem->dualSpace);CHKERRQ(ierr);
614   fem->dualSpace = sp;
615   ierr = PetscObjectReference((PetscObject) fem->dualSpace);CHKERRQ(ierr);
616   PetscFunctionReturn(0);
617 }
618 
619 /*@
620   PetscFEGetQuadrature - Returns the PetscQuadrature used to calculate inner products
621 
622   Not collective
623 
624   Input Parameter:
625 . fem - The PetscFE object
626 
627   Output Parameter:
628 . q - The PetscQuadrature object
629 
630   Level: intermediate
631 
632 .seealso: PetscFECreate()
633 @*/
634 PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q)
635 {
636   PetscFunctionBegin;
637   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
638   PetscValidPointer(q, 2);
639   *q = fem->quadrature;
640   PetscFunctionReturn(0);
641 }
642 
643 /*@
644   PetscFESetQuadrature - Sets the PetscQuadrature used to calculate inner products
645 
646   Not collective
647 
648   Input Parameters:
649 + fem - The PetscFE object
650 - q - The PetscQuadrature object
651 
652   Level: intermediate
653 
654 .seealso: PetscFECreate()
655 @*/
656 PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q)
657 {
658   PetscInt       Nc, qNc;
659   PetscErrorCode ierr;
660 
661   PetscFunctionBegin;
662   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
663   ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr);
664   ierr = PetscQuadratureGetNumComponents(q, &qNc);CHKERRQ(ierr);
665   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);
666   ierr = PetscFERestoreTabulation(fem, 0, NULL, &fem->B, &fem->D, NULL /*&(*fem)->H*/);CHKERRQ(ierr);
667   ierr = PetscQuadratureDestroy(&fem->quadrature);CHKERRQ(ierr);
668   fem->quadrature = q;
669   ierr = PetscObjectReference((PetscObject) q);CHKERRQ(ierr);
670   PetscFunctionReturn(0);
671 }
672 
673 /*@
674   PetscFEGetFaceQuadrature - Returns the PetscQuadrature used to calculate inner products on faces
675 
676   Not collective
677 
678   Input Parameter:
679 . fem - The PetscFE object
680 
681   Output Parameter:
682 . q - The PetscQuadrature object
683 
684   Level: intermediate
685 
686 .seealso: PetscFECreate()
687 @*/
688 PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q)
689 {
690   PetscFunctionBegin;
691   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
692   PetscValidPointer(q, 2);
693   *q = fem->faceQuadrature;
694   PetscFunctionReturn(0);
695 }
696 
697 /*@
698   PetscFESetFaceQuadrature - Sets the PetscQuadrature used to calculate inner products on faces
699 
700   Not collective
701 
702   Input Parameters:
703 + fem - The PetscFE object
704 - q - The PetscQuadrature object
705 
706   Level: intermediate
707 
708 .seealso: PetscFECreate()
709 @*/
710 PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q)
711 {
712   PetscErrorCode ierr;
713 
714   PetscFunctionBegin;
715   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
716   ierr = PetscFERestoreTabulation(fem, 0, NULL, &fem->Bf, &fem->Df, NULL /*&(*fem)->Hf*/);CHKERRQ(ierr);
717   ierr = PetscQuadratureDestroy(&fem->faceQuadrature);CHKERRQ(ierr);
718   fem->faceQuadrature = q;
719   ierr = PetscObjectReference((PetscObject) q);CHKERRQ(ierr);
720   PetscFunctionReturn(0);
721 }
722 
723 /*@
724   PetscFECopyQuadrature - Copy both volumetric and surface quadrature
725 
726   Not collective
727 
728   Input Parameters:
729 + sfe - The PetscFE source for the quadratures
730 - tfe - The PetscFE target for the quadratures
731 
732   Level: intermediate
733 
734 .seealso: PetscFECreate(), PetscFESetQuadrature(), PetscFESetFaceQuadrature()
735 @*/
736 PetscErrorCode PetscFECopyQuadrature(PetscFE sfe, PetscFE tfe)
737 {
738   PetscQuadrature q;
739   PetscErrorCode  ierr;
740 
741   PetscFunctionBegin;
742   PetscValidHeaderSpecific(sfe, PETSCFE_CLASSID, 1);
743   PetscValidHeaderSpecific(tfe, PETSCFE_CLASSID, 2);
744   ierr = PetscFEGetQuadrature(sfe, &q);CHKERRQ(ierr);
745   ierr = PetscFESetQuadrature(tfe,  q);CHKERRQ(ierr);
746   ierr = PetscFEGetFaceQuadrature(sfe, &q);CHKERRQ(ierr);
747   ierr = PetscFESetFaceQuadrature(tfe,  q);CHKERRQ(ierr);
748   PetscFunctionReturn(0);
749 }
750 
751 /*@C
752   PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension
753 
754   Not collective
755 
756   Input Parameter:
757 . fem - The PetscFE object
758 
759   Output Parameter:
760 . numDof - Array with the number of dofs per dimension
761 
762   Level: intermediate
763 
764 .seealso: PetscFECreate()
765 @*/
766 PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt **numDof)
767 {
768   PetscErrorCode ierr;
769 
770   PetscFunctionBegin;
771   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
772   PetscValidPointer(numDof, 2);
773   ierr = PetscDualSpaceGetNumDof(fem->dualSpace, numDof);CHKERRQ(ierr);
774   PetscFunctionReturn(0);
775 }
776 
777 /*@C
778   PetscFEGetDefaultTabulation - Returns the tabulation of the basis functions at the quadrature points
779 
780   Not collective
781 
782   Input Parameter:
783 . fem - The PetscFE object
784 
785   Output Parameters:
786 + B - The basis function values at quadrature points
787 . D - The basis function derivatives at quadrature points
788 - H - The basis function second derivatives at quadrature points
789 
790   Note:
791 $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
792 $ 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 $ 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
794 
795   Level: intermediate
796 
797 .seealso: PetscFEGetTabulation(), PetscFERestoreTabulation()
798 @*/
799 PetscErrorCode PetscFEGetDefaultTabulation(PetscFE fem, PetscReal **B, PetscReal **D, PetscReal **H)
800 {
801   PetscInt         npoints;
802   const PetscReal *points;
803   PetscErrorCode   ierr;
804 
805   PetscFunctionBegin;
806   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
807   if (B) PetscValidPointer(B, 2);
808   if (D) PetscValidPointer(D, 3);
809   if (H) PetscValidPointer(H, 4);
810   ierr = PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr);
811   if (!fem->B) {ierr = PetscFEGetTabulation(fem, npoints, points, &fem->B, &fem->D, NULL/*&fem->H*/);CHKERRQ(ierr);}
812   if (B) *B = fem->B;
813   if (D) *D = fem->D;
814   if (H) *H = fem->H;
815   PetscFunctionReturn(0);
816 }
817 
818 /*@C
819   PetscFEGetFaceTabulation - Returns the tabulation of the basis functions at the face quadrature points
820 
821   Not collective
822 
823   Input Parameter:
824 . fem - The PetscFE object
825 
826   Output Parameters:
827 + B - The basis function values at face quadrature points
828 . D - The basis function derivatives at face quadrature points
829 - H - The basis function second derivatives at face quadrature points
830 
831   Note:
832 $ Bf[((f*Nq + q)*pdim + i)*Nc + c] is the value at point f,q for basis function i and component c
833 $ 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
834 $ 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
835 
836   Level: intermediate
837 
838 .seealso: PetscFEGetDefaultTabulation(), PetscFEGetTabulation(), PetscFERestoreTabulation()
839 @*/
840 PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscReal **Bf, PetscReal **Df, PetscReal **Hf)
841 {
842   PetscErrorCode   ierr;
843 
844   PetscFunctionBegin;
845   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
846   if (Bf) PetscValidPointer(Bf, 2);
847   if (Df) PetscValidPointer(Df, 3);
848   if (Hf) PetscValidPointer(Hf, 4);
849   if (!fem->Bf) {
850     const PetscReal  xi0[3] = {-1., -1., -1.};
851     PetscReal        v0[3], J[9], detJ;
852     PetscQuadrature  fq;
853     PetscDualSpace   sp;
854     DM               dm;
855     const PetscInt  *faces;
856     PetscInt         dim, numFaces, f, npoints, q;
857     const PetscReal *points;
858     PetscReal       *facePoints;
859 
860     ierr = PetscFEGetDualSpace(fem, &sp);CHKERRQ(ierr);
861     ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
862     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
863     ierr = DMPlexGetConeSize(dm, 0, &numFaces);CHKERRQ(ierr);
864     ierr = DMPlexGetCone(dm, 0, &faces);CHKERRQ(ierr);
865     ierr = PetscFEGetFaceQuadrature(fem, &fq);CHKERRQ(ierr);
866     if (fq) {
867       ierr = PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr);
868       ierr = PetscMalloc1(numFaces*npoints*dim, &facePoints);CHKERRQ(ierr);
869       for (f = 0; f < numFaces; ++f) {
870         ierr = DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ);CHKERRQ(ierr);
871         for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim-1, xi0, v0, J, &points[q*(dim-1)], &facePoints[(f*npoints+q)*dim]);
872       }
873       ierr = PetscFEGetTabulation(fem, numFaces*npoints, facePoints, &fem->Bf, &fem->Df, NULL/*&fem->Hf*/);CHKERRQ(ierr);
874       ierr = PetscFree(facePoints);CHKERRQ(ierr);
875     }
876   }
877   if (Bf) *Bf = fem->Bf;
878   if (Df) *Df = fem->Df;
879   if (Hf) *Hf = fem->Hf;
880   PetscFunctionReturn(0);
881 }
882 
883 /*@C
884   PetscFEGetFaceTabulation - Returns the tabulation of the basis functions at the face centroid points
885 
886   Not collective
887 
888   Input Parameter:
889 . fem - The PetscFE object
890 
891   Output Parameters:
892 + B - The basis function values at face centroid points
893 . D - The basis function derivatives at face centroid points
894 - H - The basis function second derivatives at face centroid points
895 
896   Note:
897 $ Bf[(f*pdim + i)*Nc + c] is the value at point f for basis function i and component c
898 $ Df[((f*pdim + i)*Nc + c)*dim + d] is the derivative value at point f for basis function i, component c, in direction d
899 $ 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
900 
901   Level: intermediate
902 
903 .seealso: PetscFEGetFaceTabulation(), PetscFEGetDefaultTabulation(), PetscFEGetTabulation(), PetscFERestoreTabulation()
904 @*/
905 PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscReal **F)
906 {
907   PetscErrorCode   ierr;
908 
909   PetscFunctionBegin;
910   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
911   PetscValidPointer(F, 2);
912   if (!fem->F) {
913     PetscDualSpace  sp;
914     DM              dm;
915     const PetscInt *cone;
916     PetscReal      *centroids;
917     PetscInt        dim, numFaces, f;
918 
919     ierr = PetscFEGetDualSpace(fem, &sp);CHKERRQ(ierr);
920     ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
921     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
922     ierr = DMPlexGetConeSize(dm, 0, &numFaces);CHKERRQ(ierr);
923     ierr = DMPlexGetCone(dm, 0, &cone);CHKERRQ(ierr);
924     ierr = PetscMalloc1(numFaces*dim, &centroids);CHKERRQ(ierr);
925     for (f = 0; f < numFaces; ++f) {ierr = DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, &centroids[f*dim], NULL);CHKERRQ(ierr);}
926     ierr = PetscFEGetTabulation(fem, numFaces, centroids, &fem->F, NULL, NULL);CHKERRQ(ierr);
927     ierr = PetscFree(centroids);CHKERRQ(ierr);
928   }
929   *F = fem->F;
930   PetscFunctionReturn(0);
931 }
932 
933 /*@C
934   PetscFEGetTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided.
935 
936   Not collective
937 
938   Input Parameters:
939 + fem     - The PetscFE object
940 . npoints - The number of tabulation points
941 - points  - The tabulation point coordinates
942 
943   Output Parameters:
944 + B - The basis function values at tabulation points
945 . D - The basis function derivatives at tabulation points
946 - H - The basis function second derivatives at tabulation points
947 
948   Note:
949 $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
950 $ D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
951 $ 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
952 
953   Level: intermediate
954 
955 .seealso: PetscFERestoreTabulation(), PetscFEGetDefaultTabulation()
956 @*/
957 PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
958 {
959   DM               dm;
960   PetscInt         pdim; /* Dimension of FE space P */
961   PetscInt         dim;  /* Spatial dimension */
962   PetscInt         comp; /* Field components */
963   PetscErrorCode   ierr;
964 
965   PetscFunctionBegin;
966   if (!npoints) {
967     if (B) *B = NULL;
968     if (D) *D = NULL;
969     if (H) *H = NULL;
970     PetscFunctionReturn(0);
971   }
972   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
973   PetscValidPointer(points, 3);
974   if (B) PetscValidPointer(B, 4);
975   if (D) PetscValidPointer(D, 5);
976   if (H) PetscValidPointer(H, 6);
977   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
978   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
979   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
980   ierr = PetscFEGetNumComponents(fem, &comp);CHKERRQ(ierr);
981   if (B) {ierr = DMGetWorkArray(dm, npoints*pdim*comp, MPIU_REAL, B);CHKERRQ(ierr);}
982   if (!dim) {
983     if (D) *D = NULL;
984     if (H) *H = NULL;
985   } else {
986     if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*comp*dim, MPIU_REAL, D);CHKERRQ(ierr);}
987     if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*comp*dim*dim, MPIU_REAL, H);CHKERRQ(ierr);}
988   }
989   ierr = (*fem->ops->gettabulation)(fem, npoints, points, B ? *B : NULL, D ? *D : NULL, H ? *H : NULL);CHKERRQ(ierr);
990   PetscFunctionReturn(0);
991 }
992 
993 /*@C
994   PetscFERestoreTabulation - Frees memory from the associated tabulation.
995 
996   Not collective
997 
998   Input Parameters:
999 + fem     - The PetscFE object
1000 . npoints - The number of tabulation points
1001 . points  - The tabulation point coordinates
1002 . B - The basis function values at tabulation points
1003 . D - The basis function derivatives at tabulation points
1004 - H - The basis function second derivatives at tabulation points
1005 
1006   Note:
1007 $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
1008 $ D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
1009 $ 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
1010 
1011   Level: intermediate
1012 
1013 .seealso: PetscFEGetTabulation(), PetscFEGetDefaultTabulation()
1014 @*/
1015 PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
1016 {
1017   DM             dm;
1018   PetscErrorCode ierr;
1019 
1020   PetscFunctionBegin;
1021   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1022   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
1023   if (B && *B) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, B);CHKERRQ(ierr);}
1024   if (D && *D) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, D);CHKERRQ(ierr);}
1025   if (H && *H) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, H);CHKERRQ(ierr);}
1026   PetscFunctionReturn(0);
1027 }
1028 
1029 PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
1030 {
1031   PetscSpace     bsp, bsubsp;
1032   PetscDualSpace dsp, dsubsp;
1033   PetscInt       dim, depth, numComp, i, j, coneSize, order;
1034   PetscFEType    type;
1035   DM             dm;
1036   DMLabel        label;
1037   PetscReal      *xi, *v, *J, detJ;
1038   const char     *name;
1039   PetscQuadrature origin, fullQuad, subQuad;
1040   PetscErrorCode ierr;
1041 
1042   PetscFunctionBegin;
1043   PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1);
1044   PetscValidPointer(trFE,3);
1045   ierr = PetscFEGetBasisSpace(fe,&bsp);CHKERRQ(ierr);
1046   ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr);
1047   ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr);
1048   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
1049   ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr);
1050   ierr = DMLabelGetValue(label,refPoint,&depth);CHKERRQ(ierr);
1051   ierr = PetscCalloc1(depth,&xi);CHKERRQ(ierr);
1052   ierr = PetscMalloc1(dim,&v);CHKERRQ(ierr);
1053   ierr = PetscMalloc1(dim*dim,&J);CHKERRQ(ierr);
1054   for (i = 0; i < depth; i++) xi[i] = 0.;
1055   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,&origin);CHKERRQ(ierr);
1056   ierr = PetscQuadratureSetData(origin,depth,0,1,xi,NULL);CHKERRQ(ierr);
1057   ierr = DMPlexComputeCellGeometryFEM(dm,refPoint,origin,v,J,NULL,&detJ);CHKERRQ(ierr);
1058   /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */
1059   for (i = 1; i < dim; i++) {
1060     for (j = 0; j < depth; j++) {
1061       J[i * depth + j] = J[i * dim + j];
1062     }
1063   }
1064   ierr = PetscQuadratureDestroy(&origin);CHKERRQ(ierr);
1065   ierr = PetscDualSpaceGetPointSubspace(dsp,refPoint,&dsubsp);CHKERRQ(ierr);
1066   ierr = PetscSpaceCreateSubspace(bsp,dsubsp,v,J,NULL,NULL,PETSC_OWN_POINTER,&bsubsp);CHKERRQ(ierr);
1067   ierr = PetscSpaceSetUp(bsubsp);CHKERRQ(ierr);
1068   ierr = PetscFECreate(PetscObjectComm((PetscObject)fe),trFE);CHKERRQ(ierr);
1069   ierr = PetscFEGetType(fe,&type);CHKERRQ(ierr);
1070   ierr = PetscFESetType(*trFE,type);CHKERRQ(ierr);
1071   ierr = PetscFEGetNumComponents(fe,&numComp);CHKERRQ(ierr);
1072   ierr = PetscFESetNumComponents(*trFE,numComp);CHKERRQ(ierr);
1073   ierr = PetscFESetBasisSpace(*trFE,bsubsp);CHKERRQ(ierr);
1074   ierr = PetscFESetDualSpace(*trFE,dsubsp);CHKERRQ(ierr);
1075   ierr = PetscObjectGetName((PetscObject) fe, &name);CHKERRQ(ierr);
1076   if (name) {ierr = PetscFESetName(*trFE, name);CHKERRQ(ierr);}
1077   ierr = PetscFEGetQuadrature(fe,&fullQuad);CHKERRQ(ierr);
1078   ierr = PetscQuadratureGetOrder(fullQuad,&order);CHKERRQ(ierr);
1079   ierr = DMPlexGetConeSize(dm,refPoint,&coneSize);CHKERRQ(ierr);
1080   if (coneSize == 2 * depth) {
1081     ierr = PetscDTGaussTensorQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr);
1082   } else {
1083     ierr = PetscDTGaussJacobiQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr);
1084   }
1085   ierr = PetscFESetQuadrature(*trFE,subQuad);CHKERRQ(ierr);
1086   ierr = PetscFESetUp(*trFE);CHKERRQ(ierr);
1087   ierr = PetscQuadratureDestroy(&subQuad);CHKERRQ(ierr);
1088   ierr = PetscSpaceDestroy(&bsubsp);CHKERRQ(ierr);
1089   PetscFunctionReturn(0);
1090 }
1091 
1092 PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE)
1093 {
1094   PetscInt       hStart, hEnd;
1095   PetscDualSpace dsp;
1096   DM             dm;
1097   PetscErrorCode ierr;
1098 
1099   PetscFunctionBegin;
1100   PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1);
1101   PetscValidPointer(trFE,3);
1102   *trFE = NULL;
1103   ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr);
1104   ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr);
1105   ierr = DMPlexGetHeightStratum(dm,height,&hStart,&hEnd);CHKERRQ(ierr);
1106   if (hEnd <= hStart) PetscFunctionReturn(0);
1107   ierr = PetscFECreatePointTrace(fe,hStart,trFE);CHKERRQ(ierr);
1108   PetscFunctionReturn(0);
1109 }
1110 
1111 
1112 /*@
1113   PetscFEGetDimension - Get the dimension of the finite element space on a cell
1114 
1115   Not collective
1116 
1117   Input Parameter:
1118 . fe - The PetscFE
1119 
1120   Output Parameter:
1121 . dim - The dimension
1122 
1123   Level: intermediate
1124 
1125 .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
1126 @*/
1127 PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
1128 {
1129   PetscErrorCode ierr;
1130 
1131   PetscFunctionBegin;
1132   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1133   PetscValidPointer(dim, 2);
1134   if (fem->ops->getdimension) {ierr = (*fem->ops->getdimension)(fem, dim);CHKERRQ(ierr);}
1135   PetscFunctionReturn(0);
1136 }
1137 
1138 /*@C
1139   PetscFEPushforward - Map the reference element function to real space
1140 
1141   Input Parameters:
1142 + fe     - The PetscFE
1143 . fegeom - The cell geometry
1144 . Nv     - The number of function values
1145 - vals   - The function values
1146 
1147   Output Parameter:
1148 . vals   - The transformed function values
1149 
1150   Level: advanced
1151 
1152   Note: This just forwards the call onto PetscDualSpacePushforward().
1153 
1154 .seealso: PetscDualSpacePushforward()
1155 @*/
1156 PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1157 {
1158   PetscErrorCode ierr;
1159 
1160   PetscFunctionBeginHot;
1161   ierr = PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);CHKERRQ(ierr);
1162   PetscFunctionReturn(0);
1163 }
1164 
1165 /*@C
1166   PetscFEPushforwardGradient - Map the reference element function gradient to real space
1167 
1168   Input Parameters:
1169 + fe     - The PetscFE
1170 . fegeom - The cell geometry
1171 . Nv     - The number of function gradient values
1172 - vals   - The function gradient values
1173 
1174   Output Parameter:
1175 . vals   - The transformed function gradient values
1176 
1177   Level: advanced
1178 
1179   Note: This just forwards the call onto PetscDualSpacePushforwardGradient().
1180 
1181 .seealso: PetscFEPushforward(), PetscDualSpacePushforwardGradient(), PetscDualSpacePushforward()
1182 @*/
1183 PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1184 {
1185   PetscErrorCode ierr;
1186 
1187   PetscFunctionBeginHot;
1188   ierr = PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);CHKERRQ(ierr);
1189   PetscFunctionReturn(0);
1190 }
1191 
1192 /*
1193 Purpose: Compute element vector for chunk of elements
1194 
1195 Input:
1196   Sizes:
1197      Ne:  number of elements
1198      Nf:  number of fields
1199      PetscFE
1200        dim: spatial dimension
1201        Nb:  number of basis functions
1202        Nc:  number of field components
1203        PetscQuadrature
1204          Nq:  number of quadrature points
1205 
1206   Geometry:
1207      PetscFEGeom[Ne] possibly *Nq
1208        PetscReal v0s[dim]
1209        PetscReal n[dim]
1210        PetscReal jacobians[dim*dim]
1211        PetscReal jacobianInverses[dim*dim]
1212        PetscReal jacobianDeterminants
1213   FEM:
1214      PetscFE
1215        PetscQuadrature
1216          PetscReal   quadPoints[Nq*dim]
1217          PetscReal   quadWeights[Nq]
1218        PetscReal   basis[Nq*Nb*Nc]
1219        PetscReal   basisDer[Nq*Nb*Nc*dim]
1220      PetscScalar coefficients[Ne*Nb*Nc]
1221      PetscScalar elemVec[Ne*Nb*Nc]
1222 
1223   Problem:
1224      PetscInt f: the active field
1225      f0, f1
1226 
1227   Work Space:
1228      PetscFE
1229        PetscScalar f0[Nq*dim];
1230        PetscScalar f1[Nq*dim*dim];
1231        PetscScalar u[Nc];
1232        PetscScalar gradU[Nc*dim];
1233        PetscReal   x[dim];
1234        PetscScalar realSpaceDer[dim];
1235 
1236 Purpose: Compute element vector for N_cb batches of elements
1237 
1238 Input:
1239   Sizes:
1240      N_cb: Number of serial cell batches
1241 
1242   Geometry:
1243      PetscReal v0s[Ne*dim]
1244      PetscReal jacobians[Ne*dim*dim]        possibly *Nq
1245      PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq
1246      PetscReal jacobianDeterminants[Ne]     possibly *Nq
1247   FEM:
1248      static PetscReal   quadPoints[Nq*dim]
1249      static PetscReal   quadWeights[Nq]
1250      static PetscReal   basis[Nq*Nb*Nc]
1251      static PetscReal   basisDer[Nq*Nb*Nc*dim]
1252      PetscScalar coefficients[Ne*Nb*Nc]
1253      PetscScalar elemVec[Ne*Nb*Nc]
1254 
1255 ex62.c:
1256   PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[],
1257                                                const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[],
1258                                                void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]),
1259                                                void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[])
1260 
1261 ex52.c:
1262   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)
1263   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)
1264 
1265 ex52_integrateElement.cu
1266 __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec)
1267 
1268 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[],
1269                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1270                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)
1271 
1272 ex52_integrateElementOpenCL.c:
1273 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[],
1274                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1275                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)
1276 
1277 __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec)
1278 */
1279 
1280 /*@C
1281   PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration
1282 
1283   Not collective
1284 
1285   Input Parameters:
1286 + fem          - The PetscFE object for the field being integrated
1287 . prob         - The PetscDS specifying the discretizations and continuum functions
1288 . field        - The field being integrated
1289 . Ne           - The number of elements in the chunk
1290 . cgeom        - The cell geometry for each cell in the chunk
1291 . coefficients - The array of FEM basis coefficients for the elements
1292 . probAux      - The PetscDS specifying the auxiliary discretizations
1293 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1294 
1295   Output Parameter
1296 . integral     - the integral for this field
1297 
1298   Level: intermediate
1299 
1300 .seealso: PetscFEIntegrateResidual()
1301 @*/
1302 PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1303                                 const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1304 {
1305   PetscFE        fe;
1306   PetscErrorCode ierr;
1307 
1308   PetscFunctionBegin;
1309   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1310   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1311   if (fe->ops->integrate) {ierr = (*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);}
1312   PetscFunctionReturn(0);
1313 }
1314 
1315 /*@C
1316   PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration
1317 
1318   Not collective
1319 
1320   Input Parameters:
1321 + fem          - The PetscFE object for the field being integrated
1322 . prob         - The PetscDS specifying the discretizations and continuum functions
1323 . field        - The field being integrated
1324 . obj_func     - The function to be integrated
1325 . Ne           - The number of elements in the chunk
1326 . fgeom        - The face geometry for each face in the chunk
1327 . coefficients - The array of FEM basis coefficients for the elements
1328 . probAux      - The PetscDS specifying the auxiliary discretizations
1329 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1330 
1331   Output Parameter
1332 . integral     - the integral for this field
1333 
1334   Level: intermediate
1335 
1336 .seealso: PetscFEIntegrateResidual()
1337 @*/
1338 PetscErrorCode PetscFEIntegrateBd(PetscDS prob, PetscInt field,
1339                                   void (*obj_func)(PetscInt, PetscInt, PetscInt,
1340                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1341                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1342                                                    PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1343                                   PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1344 {
1345   PetscFE        fe;
1346   PetscErrorCode ierr;
1347 
1348   PetscFunctionBegin;
1349   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1350   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1351   if (fe->ops->integratebd) {ierr = (*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);}
1352   PetscFunctionReturn(0);
1353 }
1354 
1355 /*@C
1356   PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration
1357 
1358   Not collective
1359 
1360   Input Parameters:
1361 + fem          - The PetscFE object for the field being integrated
1362 . prob         - The PetscDS specifying the discretizations and continuum functions
1363 . field        - The field being integrated
1364 . Ne           - The number of elements in the chunk
1365 . cgeom        - The cell geometry for each cell in the chunk
1366 . coefficients - The array of FEM basis coefficients for the elements
1367 . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1368 . probAux      - The PetscDS specifying the auxiliary discretizations
1369 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1370 - t            - The time
1371 
1372   Output Parameter
1373 . elemVec      - the element residual vectors from each element
1374 
1375   Note:
1376 $ Loop over batch of elements (e):
1377 $   Loop over quadrature points (q):
1378 $     Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q
1379 $     Call f_0 and f_1
1380 $   Loop over element vector entries (f,fc --> i):
1381 $     elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u)
1382 
1383   Level: intermediate
1384 
1385 .seealso: PetscFEIntegrateResidual()
1386 @*/
1387 PetscErrorCode PetscFEIntegrateResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1388                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1389 {
1390   PetscFE        fe;
1391   PetscErrorCode ierr;
1392 
1393   PetscFunctionBegin;
1394   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1395   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1396   if (fe->ops->integrateresidual) {ierr = (*fe->ops->integrateresidual)(prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);}
1397   PetscFunctionReturn(0);
1398 }
1399 
1400 /*@C
1401   PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary
1402 
1403   Not collective
1404 
1405   Input Parameters:
1406 + fem          - The PetscFE object for the field being integrated
1407 . prob         - The PetscDS specifying the discretizations and continuum functions
1408 . field        - The field being integrated
1409 . Ne           - The number of elements in the chunk
1410 . fgeom        - The face geometry for each cell in the chunk
1411 . coefficients - The array of FEM basis coefficients for the elements
1412 . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1413 . probAux      - The PetscDS specifying the auxiliary discretizations
1414 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1415 - t            - The time
1416 
1417   Output Parameter
1418 . elemVec      - the element residual vectors from each element
1419 
1420   Level: intermediate
1421 
1422 .seealso: PetscFEIntegrateResidual()
1423 @*/
1424 PetscErrorCode PetscFEIntegrateBdResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom,
1425                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1426 {
1427   PetscFE        fe;
1428   PetscErrorCode ierr;
1429 
1430   PetscFunctionBegin;
1431   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1432   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1433   if (fe->ops->integratebdresidual) {ierr = (*fe->ops->integratebdresidual)(prob, field, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);}
1434   PetscFunctionReturn(0);
1435 }
1436 
1437 /*@C
1438   PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration
1439 
1440   Not collective
1441 
1442   Input Parameters:
1443 + fem          - The PetscFE object for the field being integrated
1444 . prob         - The PetscDS specifying the discretizations and continuum functions
1445 . jtype        - The type of matrix pointwise functions that should be used
1446 . fieldI       - The test field being integrated
1447 . fieldJ       - The basis field being integrated
1448 . Ne           - The number of elements in the chunk
1449 . cgeom        - The cell geometry for each cell in the chunk
1450 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1451 . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1452 . probAux      - The PetscDS specifying the auxiliary discretizations
1453 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1454 . t            - The time
1455 - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)
1456 
1457   Output Parameter
1458 . elemMat      - the element matrices for the Jacobian from each element
1459 
1460   Note:
1461 $ Loop over batch of elements (e):
1462 $   Loop over element matrix entries (f,fc,g,gc --> i,j):
1463 $     Loop over quadrature points (q):
1464 $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1465 $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1466 $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1467 $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1468 $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1469   Level: intermediate
1470 
1471 .seealso: PetscFEIntegrateResidual()
1472 @*/
1473 PetscErrorCode PetscFEIntegrateJacobian(PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom,
1474                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1475 {
1476   PetscFE        fe;
1477   PetscErrorCode ierr;
1478 
1479   PetscFunctionBegin;
1480   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1481   ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
1482   if (fe->ops->integratejacobian) {ierr = (*fe->ops->integratejacobian)(prob, jtype, fieldI, fieldJ, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);CHKERRQ(ierr);}
1483   PetscFunctionReturn(0);
1484 }
1485 
1486 /*@C
1487   PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration
1488 
1489   Not collective
1490 
1491   Input Parameters:
1492 . prob         - The PetscDS specifying the discretizations and continuum functions
1493 . fieldI       - The test field being integrated
1494 . fieldJ       - The basis field being integrated
1495 . Ne           - The number of elements in the chunk
1496 . fgeom        - The face geometry for each cell in the chunk
1497 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1498 . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1499 . probAux      - The PetscDS specifying the auxiliary discretizations
1500 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1501 . t            - The time
1502 - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)
1503 
1504   Output Parameter
1505 . elemMat              - the element matrices for the Jacobian from each element
1506 
1507   Note:
1508 $ Loop over batch of elements (e):
1509 $   Loop over element matrix entries (f,fc,g,gc --> i,j):
1510 $     Loop over quadrature points (q):
1511 $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1512 $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1513 $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1514 $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1515 $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1516   Level: intermediate
1517 
1518 .seealso: PetscFEIntegrateJacobian(), PetscFEIntegrateResidual()
1519 @*/
1520 PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom,
1521                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1522 {
1523   PetscFE        fe;
1524   PetscErrorCode ierr;
1525 
1526   PetscFunctionBegin;
1527   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1528   ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
1529   if (fe->ops->integratebdjacobian) {ierr = (*fe->ops->integratebdjacobian)(prob, fieldI, fieldJ, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);CHKERRQ(ierr);}
1530   PetscFunctionReturn(0);
1531 }
1532 
1533 /*@
1534   PetscFEGetHeightSubspace - Get the subspace of this space for a mesh point of a given height
1535 
1536   Input Parameters:
1537 + fe     - The finite element space
1538 - height - The height of the Plex point
1539 
1540   Output Parameter:
1541 . subfe  - The subspace of this FE space
1542 
1543   Note: For example, if we want the subspace of this space for a face, we would choose height = 1.
1544 
1545   Level: advanced
1546 
1547 .seealso: PetscFECreateDefault()
1548 @*/
1549 PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe)
1550 {
1551   PetscSpace      P, subP;
1552   PetscDualSpace  Q, subQ;
1553   PetscQuadrature subq;
1554   PetscFEType     fetype;
1555   PetscInt        dim, Nc;
1556   PetscErrorCode  ierr;
1557 
1558   PetscFunctionBegin;
1559   PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1);
1560   PetscValidPointer(subfe, 3);
1561   if (height == 0) {
1562     *subfe = fe;
1563     PetscFunctionReturn(0);
1564   }
1565   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
1566   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1567   ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1568   ierr = PetscFEGetFaceQuadrature(fe, &subq);CHKERRQ(ierr);
1569   ierr = PetscDualSpaceGetDimension(Q, &dim);CHKERRQ(ierr);
1570   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);}
1571   if (!fe->subspaces) {ierr = PetscCalloc1(dim, &fe->subspaces);CHKERRQ(ierr);}
1572   if (height <= dim) {
1573     if (!fe->subspaces[height-1]) {
1574       PetscFE     sub;
1575       const char *name;
1576 
1577       ierr = PetscSpaceGetHeightSubspace(P, height, &subP);CHKERRQ(ierr);
1578       ierr = PetscDualSpaceGetHeightSubspace(Q, height, &subQ);CHKERRQ(ierr);
1579       ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), &sub);CHKERRQ(ierr);
1580       ierr = PetscObjectGetName((PetscObject) fe,  &name);CHKERRQ(ierr);
1581       ierr = PetscObjectSetName((PetscObject) sub,  name);CHKERRQ(ierr);
1582       ierr = PetscFEGetType(fe, &fetype);CHKERRQ(ierr);
1583       ierr = PetscFESetType(sub, fetype);CHKERRQ(ierr);
1584       ierr = PetscFESetBasisSpace(sub, subP);CHKERRQ(ierr);
1585       ierr = PetscFESetDualSpace(sub, subQ);CHKERRQ(ierr);
1586       ierr = PetscFESetNumComponents(sub, Nc);CHKERRQ(ierr);
1587       ierr = PetscFESetUp(sub);CHKERRQ(ierr);
1588       ierr = PetscFESetQuadrature(sub, subq);CHKERRQ(ierr);
1589       fe->subspaces[height-1] = sub;
1590     }
1591     *subfe = fe->subspaces[height-1];
1592   } else {
1593     *subfe = NULL;
1594   }
1595   PetscFunctionReturn(0);
1596 }
1597 
1598 /*@
1599   PetscFERefine - Create a "refined" PetscFE object that refines the reference cell into smaller copies. This is typically used
1600   to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more
1601   sparsity). It is also used to create an interpolation between regularly refined meshes.
1602 
1603   Collective on fem
1604 
1605   Input Parameter:
1606 . fe - The initial PetscFE
1607 
1608   Output Parameter:
1609 . feRef - The refined PetscFE
1610 
1611   Level: advanced
1612 
1613 .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
1614 @*/
1615 PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)
1616 {
1617   PetscSpace       P, Pref;
1618   PetscDualSpace   Q, Qref;
1619   DM               K, Kref;
1620   PetscQuadrature  q, qref;
1621   const PetscReal *v0, *jac;
1622   PetscInt         numComp, numSubelements;
1623   PetscErrorCode   ierr;
1624 
1625   PetscFunctionBegin;
1626   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
1627   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1628   ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1629   ierr = PetscDualSpaceGetDM(Q, &K);CHKERRQ(ierr);
1630   /* Create space */
1631   ierr = PetscObjectReference((PetscObject) P);CHKERRQ(ierr);
1632   Pref = P;
1633   /* Create dual space */
1634   ierr = PetscDualSpaceDuplicate(Q, &Qref);CHKERRQ(ierr);
1635   ierr = DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref);CHKERRQ(ierr);
1636   ierr = PetscDualSpaceSetDM(Qref, Kref);CHKERRQ(ierr);
1637   ierr = DMDestroy(&Kref);CHKERRQ(ierr);
1638   ierr = PetscDualSpaceSetUp(Qref);CHKERRQ(ierr);
1639   /* Create element */
1640   ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), feRef);CHKERRQ(ierr);
1641   ierr = PetscFESetType(*feRef, PETSCFECOMPOSITE);CHKERRQ(ierr);
1642   ierr = PetscFESetBasisSpace(*feRef, Pref);CHKERRQ(ierr);
1643   ierr = PetscFESetDualSpace(*feRef, Qref);CHKERRQ(ierr);
1644   ierr = PetscFEGetNumComponents(fe,    &numComp);CHKERRQ(ierr);
1645   ierr = PetscFESetNumComponents(*feRef, numComp);CHKERRQ(ierr);
1646   ierr = PetscFESetUp(*feRef);CHKERRQ(ierr);
1647   ierr = PetscSpaceDestroy(&Pref);CHKERRQ(ierr);
1648   ierr = PetscDualSpaceDestroy(&Qref);CHKERRQ(ierr);
1649   /* Create quadrature */
1650   ierr = PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL);CHKERRQ(ierr);
1651   ierr = PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);CHKERRQ(ierr);
1652   ierr = PetscFESetQuadrature(*feRef, qref);CHKERRQ(ierr);
1653   ierr = PetscQuadratureDestroy(&qref);CHKERRQ(ierr);
1654   PetscFunctionReturn(0);
1655 }
1656 
1657 /*@C
1658   PetscFECreateDefault - Create a PetscFE for basic FEM computation
1659 
1660   Collective
1661 
1662   Input Parameters:
1663 + comm      - The MPI comm
1664 . dim       - The spatial dimension
1665 . Nc        - The number of components
1666 . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
1667 . prefix    - The options prefix, or NULL
1668 - qorder    - The quadrature order or PETSC_DETERMINE to use PetscSpace polynomial degree
1669 
1670   Output Parameter:
1671 . fem - The PetscFE object
1672 
1673   Level: beginner
1674 
1675 .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
1676 @*/
1677 PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem)
1678 {
1679   PetscQuadrature q, fq;
1680   DM              K;
1681   PetscSpace      P;
1682   PetscDualSpace  Q;
1683   PetscInt        order, quadPointsPerEdge;
1684   PetscBool       tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
1685   PetscErrorCode  ierr;
1686 
1687   PetscFunctionBegin;
1688   /* Create space */
1689   ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr);
1690   ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr);
1691   ierr = PetscSpacePolynomialSetTensor(P, tensor);CHKERRQ(ierr);
1692   ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr);
1693   ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr);
1694   ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr);
1695   ierr = PetscSpaceSetUp(P);CHKERRQ(ierr);
1696   ierr = PetscSpaceGetDegree(P, &order, NULL);CHKERRQ(ierr);
1697   ierr = PetscSpacePolynomialGetTensor(P, &tensor);CHKERRQ(ierr);
1698   /* Create dual space */
1699   ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr);
1700   ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr);
1701   ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr);
1702   ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr);
1703   ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr);
1704   ierr = DMDestroy(&K);CHKERRQ(ierr);
1705   ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr);
1706   ierr = PetscDualSpaceSetOrder(Q, order);CHKERRQ(ierr);
1707   ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor);CHKERRQ(ierr);
1708   ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr);
1709   ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr);
1710   /* Create element */
1711   ierr = PetscFECreate(comm, fem);CHKERRQ(ierr);
1712   ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr);
1713   ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr);
1714   ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr);
1715   ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr);
1716   ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr);
1717   ierr = PetscFESetUp(*fem);CHKERRQ(ierr);
1718   ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr);
1719   ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr);
1720   /* Create quadrature (with specified order if given) */
1721   qorder = qorder >= 0 ? qorder : order;
1722   ierr = PetscObjectOptionsBegin((PetscObject)*fem);CHKERRQ(ierr);
1723   ierr = PetscOptionsBoundedInt("-petscfe_default_quadrature_order","Quadrature order is one less than quadrature points per edge","PetscFECreateDefault",qorder,&qorder,NULL,0);CHKERRQ(ierr);
1724   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1725   quadPointsPerEdge = PetscMax(qorder + 1,1);
1726   if (isSimplex) {
1727     ierr = PetscDTGaussJacobiQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr);
1728     ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr);
1729   } else {
1730     ierr = PetscDTGaussTensorQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr);
1731     ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr);
1732   }
1733   ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr);
1734   ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr);
1735   ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr);
1736   ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr);
1737   PetscFunctionReturn(0);
1738 }
1739 
1740 /*@C
1741   PetscFESetName - Names the FE and its subobjects
1742 
1743   Not collective
1744 
1745   Input Parameters:
1746 + fe   - The PetscFE
1747 - name - The name
1748 
1749   Level: intermediate
1750 
1751 .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
1752 @*/
1753 PetscErrorCode PetscFESetName(PetscFE fe, const char name[])
1754 {
1755   PetscSpace     P;
1756   PetscDualSpace Q;
1757   PetscErrorCode ierr;
1758 
1759   PetscFunctionBegin;
1760   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
1761   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1762   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
1763   ierr = PetscObjectSetName((PetscObject) P,  name);CHKERRQ(ierr);
1764   ierr = PetscObjectSetName((PetscObject) Q,  name);CHKERRQ(ierr);
1765   PetscFunctionReturn(0);
1766 }
1767 
1768 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[])
1769 {
1770   PetscInt       dOffset = 0, fOffset = 0, f;
1771   PetscErrorCode ierr;
1772 
1773   for (f = 0; f < Nf; ++f) {
1774     PetscFE          fe;
1775     const PetscInt   Nbf = Nb[f], Ncf = Nc[f];
1776     const PetscReal *Bq = &basisField[f][q*Nbf*Ncf];
1777     const PetscReal *Dq = &basisFieldDer[f][q*Nbf*Ncf*dim];
1778     PetscInt         b, c, d;
1779 
1780     ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);CHKERRQ(ierr);
1781     for (c = 0; c < Ncf; ++c)     u[fOffset+c] = 0.0;
1782     for (d = 0; d < dim*Ncf; ++d) u_x[fOffset*dim+d] = 0.0;
1783     for (b = 0; b < Nbf; ++b) {
1784       for (c = 0; c < Ncf; ++c) {
1785         const PetscInt cidx = b*Ncf+c;
1786 
1787         u[fOffset+c] += Bq[cidx]*coefficients[dOffset+b];
1788         for (d = 0; d < dim; ++d) u_x[(fOffset+c)*dim+d] += Dq[cidx*dim+d]*coefficients[dOffset+b];
1789       }
1790     }
1791     ierr = PetscFEPushforward(fe, fegeom, 1, &u[fOffset]);CHKERRQ(ierr);
1792     ierr = PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset*dim]);CHKERRQ(ierr);
1793     if (u_t) {
1794       for (c = 0; c < Ncf; ++c) u_t[fOffset+c] = 0.0;
1795       for (b = 0; b < Nbf; ++b) {
1796         for (c = 0; c < Ncf; ++c) {
1797           const PetscInt cidx = b*Ncf+c;
1798 
1799           u_t[fOffset+c] += Bq[cidx]*coefficients_t[dOffset+b];
1800         }
1801       }
1802       ierr = PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]);CHKERRQ(ierr);
1803     }
1804     fOffset += Ncf;
1805     dOffset += Nbf;
1806   }
1807   return 0;
1808 }
1809 
1810 PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[])
1811 {
1812   PetscFE        fe;
1813   PetscReal     *faceBasis;
1814   PetscInt       Nb, Nc, b, c;
1815   PetscErrorCode ierr;
1816 
1817   if (!prob) return 0;
1818   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1819   ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1820   ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1821   ierr = PetscFEGetFaceCentroidTabulation(fe, &faceBasis);CHKERRQ(ierr);
1822   for (c = 0; c < Nc; ++c) {u[c] = 0.0;}
1823   for (b = 0; b < Nb; ++b) {
1824     for (c = 0; c < Nc; ++c) {
1825       const PetscInt cidx = b*Nc+c;
1826 
1827       u[c] += coefficients[cidx]*faceBasis[faceLoc*Nb*Nc+cidx];
1828     }
1829   }
1830   return 0;
1831 }
1832 
1833 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[])
1834 {
1835   PetscInt       q, b, c, d;
1836   PetscErrorCode ierr;
1837 
1838   for (b = 0; b < Nb; ++b) elemVec[b] = 0.0;
1839   for (q = 0; q < Nq; ++q) {
1840     for (b = 0; b < Nb; ++b) {
1841       for (c = 0; c < Nc; ++c) {
1842         const PetscInt bcidx = b*Nc+c;
1843 
1844         tmpBasis[bcidx] = basis[q*Nb*Nc+bcidx];
1845         for (d = 0; d < dim; ++d) tmpBasisDer[bcidx*dim+d] = basisDer[q*Nb*Nc*dim+bcidx*dim+d];
1846       }
1847     }
1848     ierr = PetscFEPushforward(fe, fegeom, Nb, tmpBasis);CHKERRQ(ierr);
1849     ierr = PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer);CHKERRQ(ierr);
1850     for (b = 0; b < Nb; ++b) {
1851       for (c = 0; c < Nc; ++c) {
1852         const PetscInt bcidx = b*Nc+c;
1853         const PetscInt qcidx = q*Nc+c;
1854 
1855         elemVec[b] += tmpBasis[bcidx]*f0[qcidx];
1856         for (d = 0; d < dim; ++d) elemVec[b] += tmpBasisDer[bcidx*dim+d]*f1[qcidx*dim+d];
1857       }
1858     }
1859   }
1860   return(0);
1861 }
1862 
1863 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[])
1864 {
1865   PetscInt       f, fc, g, gc, df, dg;
1866   PetscErrorCode ierr;
1867 
1868   for (f = 0; f < NbI; ++f) {
1869     for (fc = 0; fc < NcI; ++fc) {
1870       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
1871 
1872       tmpBasisI[fidx] = basisI[fidx];
1873       for (df = 0; df < dim; ++df) tmpBasisDerI[fidx*dim+df] = basisDerI[fidx*dim+df];
1874     }
1875   }
1876   ierr = PetscFEPushforward(feI, fegeom, NbI, tmpBasisI);CHKERRQ(ierr);
1877   ierr = PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI);CHKERRQ(ierr);
1878   for (g = 0; g < NbJ; ++g) {
1879     for (gc = 0; gc < NcJ; ++gc) {
1880       const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
1881 
1882       tmpBasisJ[gidx] = basisJ[gidx];
1883       for (dg = 0; dg < dim; ++dg) tmpBasisDerJ[gidx*dim+dg] = basisDerJ[gidx*dim+dg];
1884     }
1885   }
1886   ierr = PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ);CHKERRQ(ierr);
1887   ierr = PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ);CHKERRQ(ierr);
1888   for (f = 0; f < NbI; ++f) {
1889     for (fc = 0; fc < NcI; ++fc) {
1890       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
1891       const PetscInt i    = offsetI+f; /* Element matrix row */
1892       for (g = 0; g < NbJ; ++g) {
1893         for (gc = 0; gc < NcJ; ++gc) {
1894           const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
1895           const PetscInt j    = offsetJ+g; /* Element matrix column */
1896           const PetscInt fOff = eOffset+i*totDim+j;
1897 
1898           elemMat[fOff] += tmpBasisI[fidx]*g0[fc*NcJ+gc]*tmpBasisJ[gidx];
1899           for (df = 0; df < dim; ++df) {
1900             elemMat[fOff] += tmpBasisI[fidx]*g1[(fc*NcJ+gc)*dim+df]*tmpBasisDerJ[gidx*dim+df];
1901             elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g2[(fc*NcJ+gc)*dim+df]*tmpBasisJ[gidx];
1902             for (dg = 0; dg < dim; ++dg) {
1903               elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g3[((fc*NcJ+gc)*dim+df)*dim+dg]*tmpBasisDerJ[gidx*dim+dg];
1904             }
1905           }
1906         }
1907       }
1908     }
1909   }
1910   return(0);
1911 }
1912 
1913 PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom)
1914 {
1915   PetscDualSpace  dsp;
1916   DM              dm;
1917   PetscQuadrature quadDef;
1918   PetscInt        dim, cdim, Nq;
1919   PetscErrorCode  ierr;
1920 
1921   PetscFunctionBegin;
1922   ierr = PetscFEGetDualSpace(fe, &dsp);CHKERRQ(ierr);
1923   ierr = PetscDualSpaceGetDM(dsp, &dm);CHKERRQ(ierr);
1924   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1925   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1926   ierr = PetscFEGetQuadrature(fe, &quadDef);CHKERRQ(ierr);
1927   quad = quad ? quad : quadDef;
1928   ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1929   ierr = PetscMalloc1(Nq*cdim,      &cgeom->v);CHKERRQ(ierr);
1930   ierr = PetscMalloc1(Nq*cdim*cdim, &cgeom->J);CHKERRQ(ierr);
1931   ierr = PetscMalloc1(Nq*cdim*cdim, &cgeom->invJ);CHKERRQ(ierr);
1932   ierr = PetscMalloc1(Nq,           &cgeom->detJ);CHKERRQ(ierr);
1933   cgeom->dim       = dim;
1934   cgeom->dimEmbed  = cdim;
1935   cgeom->numCells  = 1;
1936   cgeom->numPoints = Nq;
1937   ierr = DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ);CHKERRQ(ierr);
1938   PetscFunctionReturn(0);
1939 }
1940 
1941 PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom)
1942 {
1943   PetscErrorCode ierr;
1944 
1945   PetscFunctionBegin;
1946   ierr = PetscFree(cgeom->v);CHKERRQ(ierr);
1947   ierr = PetscFree(cgeom->J);CHKERRQ(ierr);
1948   ierr = PetscFree(cgeom->invJ);CHKERRQ(ierr);
1949   ierr = PetscFree(cgeom->detJ);CHKERRQ(ierr);
1950   PetscFunctionReturn(0);
1951 }
1952