xref: /petsc/src/dm/dt/fe/interface/fe.c (revision b2a1bb3a1fb00847ede9cc2c2c2974a94b6b9cb3)
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 || !fem->dualSpace) {
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   if (fem->dualSpace) {
1023     ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
1024     if (B && *B) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, B);CHKERRQ(ierr);}
1025     if (D && *D) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, D);CHKERRQ(ierr);}
1026     if (H && *H) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, H);CHKERRQ(ierr);}
1027   }
1028   PetscFunctionReturn(0);
1029 }
1030 
1031 PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
1032 {
1033   PetscSpace     bsp, bsubsp;
1034   PetscDualSpace dsp, dsubsp;
1035   PetscInt       dim, depth, numComp, i, j, coneSize, order;
1036   PetscFEType    type;
1037   DM             dm;
1038   DMLabel        label;
1039   PetscReal      *xi, *v, *J, detJ;
1040   const char     *name;
1041   PetscQuadrature origin, fullQuad, subQuad;
1042   PetscErrorCode ierr;
1043 
1044   PetscFunctionBegin;
1045   PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1);
1046   PetscValidPointer(trFE,3);
1047   ierr = PetscFEGetBasisSpace(fe,&bsp);CHKERRQ(ierr);
1048   ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr);
1049   ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr);
1050   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
1051   ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr);
1052   ierr = DMLabelGetValue(label,refPoint,&depth);CHKERRQ(ierr);
1053   ierr = PetscCalloc1(depth,&xi);CHKERRQ(ierr);
1054   ierr = PetscMalloc1(dim,&v);CHKERRQ(ierr);
1055   ierr = PetscMalloc1(dim*dim,&J);CHKERRQ(ierr);
1056   for (i = 0; i < depth; i++) xi[i] = 0.;
1057   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,&origin);CHKERRQ(ierr);
1058   ierr = PetscQuadratureSetData(origin,depth,0,1,xi,NULL);CHKERRQ(ierr);
1059   ierr = DMPlexComputeCellGeometryFEM(dm,refPoint,origin,v,J,NULL,&detJ);CHKERRQ(ierr);
1060   /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */
1061   for (i = 1; i < dim; i++) {
1062     for (j = 0; j < depth; j++) {
1063       J[i * depth + j] = J[i * dim + j];
1064     }
1065   }
1066   ierr = PetscQuadratureDestroy(&origin);CHKERRQ(ierr);
1067   ierr = PetscDualSpaceGetPointSubspace(dsp,refPoint,&dsubsp);CHKERRQ(ierr);
1068   ierr = PetscSpaceCreateSubspace(bsp,dsubsp,v,J,NULL,NULL,PETSC_OWN_POINTER,&bsubsp);CHKERRQ(ierr);
1069   ierr = PetscSpaceSetUp(bsubsp);CHKERRQ(ierr);
1070   ierr = PetscFECreate(PetscObjectComm((PetscObject)fe),trFE);CHKERRQ(ierr);
1071   ierr = PetscFEGetType(fe,&type);CHKERRQ(ierr);
1072   ierr = PetscFESetType(*trFE,type);CHKERRQ(ierr);
1073   ierr = PetscFEGetNumComponents(fe,&numComp);CHKERRQ(ierr);
1074   ierr = PetscFESetNumComponents(*trFE,numComp);CHKERRQ(ierr);
1075   ierr = PetscFESetBasisSpace(*trFE,bsubsp);CHKERRQ(ierr);
1076   ierr = PetscFESetDualSpace(*trFE,dsubsp);CHKERRQ(ierr);
1077   ierr = PetscObjectGetName((PetscObject) fe, &name);CHKERRQ(ierr);
1078   if (name) {ierr = PetscFESetName(*trFE, name);CHKERRQ(ierr);}
1079   ierr = PetscFEGetQuadrature(fe,&fullQuad);CHKERRQ(ierr);
1080   ierr = PetscQuadratureGetOrder(fullQuad,&order);CHKERRQ(ierr);
1081   ierr = DMPlexGetConeSize(dm,refPoint,&coneSize);CHKERRQ(ierr);
1082   if (coneSize == 2 * depth) {
1083     ierr = PetscDTGaussTensorQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr);
1084   } else {
1085     ierr = PetscDTGaussJacobiQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr);
1086   }
1087   ierr = PetscFESetQuadrature(*trFE,subQuad);CHKERRQ(ierr);
1088   ierr = PetscFESetUp(*trFE);CHKERRQ(ierr);
1089   ierr = PetscQuadratureDestroy(&subQuad);CHKERRQ(ierr);
1090   ierr = PetscSpaceDestroy(&bsubsp);CHKERRQ(ierr);
1091   PetscFunctionReturn(0);
1092 }
1093 
1094 PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE)
1095 {
1096   PetscInt       hStart, hEnd;
1097   PetscDualSpace dsp;
1098   DM             dm;
1099   PetscErrorCode ierr;
1100 
1101   PetscFunctionBegin;
1102   PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1);
1103   PetscValidPointer(trFE,3);
1104   *trFE = NULL;
1105   ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr);
1106   ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr);
1107   ierr = DMPlexGetHeightStratum(dm,height,&hStart,&hEnd);CHKERRQ(ierr);
1108   if (hEnd <= hStart) PetscFunctionReturn(0);
1109   ierr = PetscFECreatePointTrace(fe,hStart,trFE);CHKERRQ(ierr);
1110   PetscFunctionReturn(0);
1111 }
1112 
1113 
1114 /*@
1115   PetscFEGetDimension - Get the dimension of the finite element space on a cell
1116 
1117   Not collective
1118 
1119   Input Parameter:
1120 . fe - The PetscFE
1121 
1122   Output Parameter:
1123 . dim - The dimension
1124 
1125   Level: intermediate
1126 
1127 .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
1128 @*/
1129 PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
1130 {
1131   PetscErrorCode ierr;
1132 
1133   PetscFunctionBegin;
1134   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1135   PetscValidPointer(dim, 2);
1136   if (fem->ops->getdimension) {ierr = (*fem->ops->getdimension)(fem, dim);CHKERRQ(ierr);}
1137   PetscFunctionReturn(0);
1138 }
1139 
1140 /*@C
1141   PetscFEPushforward - Map the reference element function to real space
1142 
1143   Input Parameters:
1144 + fe     - The PetscFE
1145 . fegeom - The cell geometry
1146 . Nv     - The number of function values
1147 - vals   - The function values
1148 
1149   Output Parameter:
1150 . vals   - The transformed function values
1151 
1152   Level: advanced
1153 
1154   Note: This just forwards the call onto PetscDualSpacePushforward().
1155 
1156 .seealso: PetscDualSpacePushforward()
1157 @*/
1158 PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1159 {
1160   PetscErrorCode ierr;
1161 
1162   PetscFunctionBeginHot;
1163   ierr = PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);CHKERRQ(ierr);
1164   PetscFunctionReturn(0);
1165 }
1166 
1167 /*@C
1168   PetscFEPushforwardGradient - Map the reference element function gradient to real space
1169 
1170   Input Parameters:
1171 + fe     - The PetscFE
1172 . fegeom - The cell geometry
1173 . Nv     - The number of function gradient values
1174 - vals   - The function gradient values
1175 
1176   Output Parameter:
1177 . vals   - The transformed function gradient values
1178 
1179   Level: advanced
1180 
1181   Note: This just forwards the call onto PetscDualSpacePushforwardGradient().
1182 
1183 .seealso: PetscFEPushforward(), PetscDualSpacePushforwardGradient(), PetscDualSpacePushforward()
1184 @*/
1185 PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1186 {
1187   PetscErrorCode ierr;
1188 
1189   PetscFunctionBeginHot;
1190   ierr = PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);CHKERRQ(ierr);
1191   PetscFunctionReturn(0);
1192 }
1193 
1194 /*
1195 Purpose: Compute element vector for chunk of elements
1196 
1197 Input:
1198   Sizes:
1199      Ne:  number of elements
1200      Nf:  number of fields
1201      PetscFE
1202        dim: spatial dimension
1203        Nb:  number of basis functions
1204        Nc:  number of field components
1205        PetscQuadrature
1206          Nq:  number of quadrature points
1207 
1208   Geometry:
1209      PetscFEGeom[Ne] possibly *Nq
1210        PetscReal v0s[dim]
1211        PetscReal n[dim]
1212        PetscReal jacobians[dim*dim]
1213        PetscReal jacobianInverses[dim*dim]
1214        PetscReal jacobianDeterminants
1215   FEM:
1216      PetscFE
1217        PetscQuadrature
1218          PetscReal   quadPoints[Nq*dim]
1219          PetscReal   quadWeights[Nq]
1220        PetscReal   basis[Nq*Nb*Nc]
1221        PetscReal   basisDer[Nq*Nb*Nc*dim]
1222      PetscScalar coefficients[Ne*Nb*Nc]
1223      PetscScalar elemVec[Ne*Nb*Nc]
1224 
1225   Problem:
1226      PetscInt f: the active field
1227      f0, f1
1228 
1229   Work Space:
1230      PetscFE
1231        PetscScalar f0[Nq*dim];
1232        PetscScalar f1[Nq*dim*dim];
1233        PetscScalar u[Nc];
1234        PetscScalar gradU[Nc*dim];
1235        PetscReal   x[dim];
1236        PetscScalar realSpaceDer[dim];
1237 
1238 Purpose: Compute element vector for N_cb batches of elements
1239 
1240 Input:
1241   Sizes:
1242      N_cb: Number of serial cell batches
1243 
1244   Geometry:
1245      PetscReal v0s[Ne*dim]
1246      PetscReal jacobians[Ne*dim*dim]        possibly *Nq
1247      PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq
1248      PetscReal jacobianDeterminants[Ne]     possibly *Nq
1249   FEM:
1250      static PetscReal   quadPoints[Nq*dim]
1251      static PetscReal   quadWeights[Nq]
1252      static PetscReal   basis[Nq*Nb*Nc]
1253      static PetscReal   basisDer[Nq*Nb*Nc*dim]
1254      PetscScalar coefficients[Ne*Nb*Nc]
1255      PetscScalar elemVec[Ne*Nb*Nc]
1256 
1257 ex62.c:
1258   PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[],
1259                                                const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[],
1260                                                void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]),
1261                                                void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[])
1262 
1263 ex52.c:
1264   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)
1265   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)
1266 
1267 ex52_integrateElement.cu
1268 __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec)
1269 
1270 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[],
1271                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1272                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)
1273 
1274 ex52_integrateElementOpenCL.c:
1275 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[],
1276                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1277                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)
1278 
1279 __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec)
1280 */
1281 
1282 /*@C
1283   PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration
1284 
1285   Not collective
1286 
1287   Input Parameters:
1288 + fem          - The PetscFE object for the field being integrated
1289 . prob         - The PetscDS specifying the discretizations and continuum functions
1290 . field        - The field being integrated
1291 . Ne           - The number of elements in the chunk
1292 . cgeom        - The cell geometry for each cell in the chunk
1293 . coefficients - The array of FEM basis coefficients for the elements
1294 . probAux      - The PetscDS specifying the auxiliary discretizations
1295 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1296 
1297   Output Parameter
1298 . integral     - the integral for this field
1299 
1300   Level: intermediate
1301 
1302 .seealso: PetscFEIntegrateResidual()
1303 @*/
1304 PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1305                                 const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1306 {
1307   PetscFE        fe;
1308   PetscErrorCode ierr;
1309 
1310   PetscFunctionBegin;
1311   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1312   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1313   if (fe->ops->integrate) {ierr = (*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);}
1314   PetscFunctionReturn(0);
1315 }
1316 
1317 /*@C
1318   PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration
1319 
1320   Not collective
1321 
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 . obj_func     - The function to be integrated
1327 . Ne           - The number of elements in the chunk
1328 . fgeom        - The face geometry for each face in the chunk
1329 . coefficients - The array of FEM basis coefficients for the elements
1330 . probAux      - The PetscDS specifying the auxiliary discretizations
1331 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1332 
1333   Output Parameter
1334 . integral     - the integral for this field
1335 
1336   Level: intermediate
1337 
1338 .seealso: PetscFEIntegrateResidual()
1339 @*/
1340 PetscErrorCode PetscFEIntegrateBd(PetscDS prob, PetscInt field,
1341                                   void (*obj_func)(PetscInt, PetscInt, PetscInt,
1342                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1343                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1344                                                    PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1345                                   PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1346 {
1347   PetscFE        fe;
1348   PetscErrorCode ierr;
1349 
1350   PetscFunctionBegin;
1351   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1352   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1353   if (fe->ops->integratebd) {ierr = (*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);}
1354   PetscFunctionReturn(0);
1355 }
1356 
1357 /*@C
1358   PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration
1359 
1360   Not collective
1361 
1362   Input Parameters:
1363 + fem          - The PetscFE object for the field being integrated
1364 . prob         - The PetscDS specifying the discretizations and continuum functions
1365 . field        - The field being integrated
1366 . Ne           - The number of elements in the chunk
1367 . cgeom        - The cell geometry for each cell in the chunk
1368 . coefficients - The array of FEM basis coefficients for the elements
1369 . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1370 . probAux      - The PetscDS specifying the auxiliary discretizations
1371 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1372 - t            - The time
1373 
1374   Output Parameter
1375 . elemVec      - the element residual vectors from each element
1376 
1377   Note:
1378 $ Loop over batch of elements (e):
1379 $   Loop over quadrature points (q):
1380 $     Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q
1381 $     Call f_0 and f_1
1382 $   Loop over element vector entries (f,fc --> i):
1383 $     elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u)
1384 
1385   Level: intermediate
1386 
1387 .seealso: PetscFEIntegrateResidual()
1388 @*/
1389 PetscErrorCode PetscFEIntegrateResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1390                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1391 {
1392   PetscFE        fe;
1393   PetscErrorCode ierr;
1394 
1395   PetscFunctionBegin;
1396   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1397   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1398   if (fe->ops->integrateresidual) {ierr = (*fe->ops->integrateresidual)(prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);}
1399   PetscFunctionReturn(0);
1400 }
1401 
1402 /*@C
1403   PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary
1404 
1405   Not collective
1406 
1407   Input Parameters:
1408 + fem          - The PetscFE object for the field being integrated
1409 . prob         - The PetscDS specifying the discretizations and continuum functions
1410 . field        - The field being integrated
1411 . Ne           - The number of elements in the chunk
1412 . fgeom        - The face geometry for each cell in the chunk
1413 . coefficients - The array of FEM basis coefficients for the elements
1414 . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1415 . probAux      - The PetscDS specifying the auxiliary discretizations
1416 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1417 - t            - The time
1418 
1419   Output Parameter
1420 . elemVec      - the element residual vectors from each element
1421 
1422   Level: intermediate
1423 
1424 .seealso: PetscFEIntegrateResidual()
1425 @*/
1426 PetscErrorCode PetscFEIntegrateBdResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom,
1427                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1428 {
1429   PetscFE        fe;
1430   PetscErrorCode ierr;
1431 
1432   PetscFunctionBegin;
1433   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1434   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1435   if (fe->ops->integratebdresidual) {ierr = (*fe->ops->integratebdresidual)(prob, field, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);}
1436   PetscFunctionReturn(0);
1437 }
1438 
1439 /*@C
1440   PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration
1441 
1442   Not collective
1443 
1444   Input Parameters:
1445 + fem          - The PetscFE object for the field being integrated
1446 . prob         - The PetscDS specifying the discretizations and continuum functions
1447 . jtype        - The type of matrix pointwise functions that should be used
1448 . fieldI       - The test field being integrated
1449 . fieldJ       - The basis field being integrated
1450 . Ne           - The number of elements in the chunk
1451 . cgeom        - The cell geometry for each cell in the chunk
1452 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1453 . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1454 . probAux      - The PetscDS specifying the auxiliary discretizations
1455 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1456 . t            - The time
1457 - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)
1458 
1459   Output Parameter
1460 . elemMat      - the element matrices for the Jacobian from each element
1461 
1462   Note:
1463 $ Loop over batch of elements (e):
1464 $   Loop over element matrix entries (f,fc,g,gc --> i,j):
1465 $     Loop over quadrature points (q):
1466 $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1467 $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1468 $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1469 $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1470 $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1471   Level: intermediate
1472 
1473 .seealso: PetscFEIntegrateResidual()
1474 @*/
1475 PetscErrorCode PetscFEIntegrateJacobian(PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom,
1476                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1477 {
1478   PetscFE        fe;
1479   PetscErrorCode ierr;
1480 
1481   PetscFunctionBegin;
1482   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1483   ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
1484   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);}
1485   PetscFunctionReturn(0);
1486 }
1487 
1488 /*@C
1489   PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration
1490 
1491   Not collective
1492 
1493   Input Parameters:
1494 . prob         - The PetscDS specifying the discretizations and continuum functions
1495 . fieldI       - The test field being integrated
1496 . fieldJ       - The basis field being integrated
1497 . Ne           - The number of elements in the chunk
1498 . fgeom        - The face geometry for each cell in the chunk
1499 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1500 . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1501 . probAux      - The PetscDS specifying the auxiliary discretizations
1502 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1503 . t            - The time
1504 - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)
1505 
1506   Output Parameter
1507 . elemMat              - the element matrices for the Jacobian from each element
1508 
1509   Note:
1510 $ Loop over batch of elements (e):
1511 $   Loop over element matrix entries (f,fc,g,gc --> i,j):
1512 $     Loop over quadrature points (q):
1513 $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1514 $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1515 $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1516 $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1517 $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1518   Level: intermediate
1519 
1520 .seealso: PetscFEIntegrateJacobian(), PetscFEIntegrateResidual()
1521 @*/
1522 PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom,
1523                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1524 {
1525   PetscFE        fe;
1526   PetscErrorCode ierr;
1527 
1528   PetscFunctionBegin;
1529   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1530   ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
1531   if (fe->ops->integratebdjacobian) {ierr = (*fe->ops->integratebdjacobian)(prob, fieldI, fieldJ, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);CHKERRQ(ierr);}
1532   PetscFunctionReturn(0);
1533 }
1534 
1535 /*@
1536   PetscFEGetHeightSubspace - Get the subspace of this space for a mesh point of a given height
1537 
1538   Input Parameters:
1539 + fe     - The finite element space
1540 - height - The height of the Plex point
1541 
1542   Output Parameter:
1543 . subfe  - The subspace of this FE space
1544 
1545   Note: For example, if we want the subspace of this space for a face, we would choose height = 1.
1546 
1547   Level: advanced
1548 
1549 .seealso: PetscFECreateDefault()
1550 @*/
1551 PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe)
1552 {
1553   PetscSpace      P, subP;
1554   PetscDualSpace  Q, subQ;
1555   PetscQuadrature subq;
1556   PetscFEType     fetype;
1557   PetscInt        dim, Nc;
1558   PetscErrorCode  ierr;
1559 
1560   PetscFunctionBegin;
1561   PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1);
1562   PetscValidPointer(subfe, 3);
1563   if (height == 0) {
1564     *subfe = fe;
1565     PetscFunctionReturn(0);
1566   }
1567   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
1568   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1569   ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1570   ierr = PetscFEGetFaceQuadrature(fe, &subq);CHKERRQ(ierr);
1571   ierr = PetscDualSpaceGetDimension(Q, &dim);CHKERRQ(ierr);
1572   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);}
1573   if (!fe->subspaces) {ierr = PetscCalloc1(dim, &fe->subspaces);CHKERRQ(ierr);}
1574   if (height <= dim) {
1575     if (!fe->subspaces[height-1]) {
1576       PetscFE     sub;
1577       const char *name;
1578 
1579       ierr = PetscSpaceGetHeightSubspace(P, height, &subP);CHKERRQ(ierr);
1580       ierr = PetscDualSpaceGetHeightSubspace(Q, height, &subQ);CHKERRQ(ierr);
1581       ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), &sub);CHKERRQ(ierr);
1582       ierr = PetscObjectGetName((PetscObject) fe,  &name);CHKERRQ(ierr);
1583       ierr = PetscObjectSetName((PetscObject) sub,  name);CHKERRQ(ierr);
1584       ierr = PetscFEGetType(fe, &fetype);CHKERRQ(ierr);
1585       ierr = PetscFESetType(sub, fetype);CHKERRQ(ierr);
1586       ierr = PetscFESetBasisSpace(sub, subP);CHKERRQ(ierr);
1587       ierr = PetscFESetDualSpace(sub, subQ);CHKERRQ(ierr);
1588       ierr = PetscFESetNumComponents(sub, Nc);CHKERRQ(ierr);
1589       ierr = PetscFESetUp(sub);CHKERRQ(ierr);
1590       ierr = PetscFESetQuadrature(sub, subq);CHKERRQ(ierr);
1591       fe->subspaces[height-1] = sub;
1592     }
1593     *subfe = fe->subspaces[height-1];
1594   } else {
1595     *subfe = NULL;
1596   }
1597   PetscFunctionReturn(0);
1598 }
1599 
1600 /*@
1601   PetscFERefine - Create a "refined" PetscFE object that refines the reference cell into smaller copies. This is typically used
1602   to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more
1603   sparsity). It is also used to create an interpolation between regularly refined meshes.
1604 
1605   Collective on fem
1606 
1607   Input Parameter:
1608 . fe - The initial PetscFE
1609 
1610   Output Parameter:
1611 . feRef - The refined PetscFE
1612 
1613   Level: advanced
1614 
1615 .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
1616 @*/
1617 PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)
1618 {
1619   PetscSpace       P, Pref;
1620   PetscDualSpace   Q, Qref;
1621   DM               K, Kref;
1622   PetscQuadrature  q, qref;
1623   const PetscReal *v0, *jac;
1624   PetscInt         numComp, numSubelements;
1625   PetscErrorCode   ierr;
1626 
1627   PetscFunctionBegin;
1628   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
1629   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1630   ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1631   ierr = PetscDualSpaceGetDM(Q, &K);CHKERRQ(ierr);
1632   /* Create space */
1633   ierr = PetscObjectReference((PetscObject) P);CHKERRQ(ierr);
1634   Pref = P;
1635   /* Create dual space */
1636   ierr = PetscDualSpaceDuplicate(Q, &Qref);CHKERRQ(ierr);
1637   ierr = DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref);CHKERRQ(ierr);
1638   ierr = PetscDualSpaceSetDM(Qref, Kref);CHKERRQ(ierr);
1639   ierr = DMDestroy(&Kref);CHKERRQ(ierr);
1640   ierr = PetscDualSpaceSetUp(Qref);CHKERRQ(ierr);
1641   /* Create element */
1642   ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), feRef);CHKERRQ(ierr);
1643   ierr = PetscFESetType(*feRef, PETSCFECOMPOSITE);CHKERRQ(ierr);
1644   ierr = PetscFESetBasisSpace(*feRef, Pref);CHKERRQ(ierr);
1645   ierr = PetscFESetDualSpace(*feRef, Qref);CHKERRQ(ierr);
1646   ierr = PetscFEGetNumComponents(fe,    &numComp);CHKERRQ(ierr);
1647   ierr = PetscFESetNumComponents(*feRef, numComp);CHKERRQ(ierr);
1648   ierr = PetscFESetUp(*feRef);CHKERRQ(ierr);
1649   ierr = PetscSpaceDestroy(&Pref);CHKERRQ(ierr);
1650   ierr = PetscDualSpaceDestroy(&Qref);CHKERRQ(ierr);
1651   /* Create quadrature */
1652   ierr = PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL);CHKERRQ(ierr);
1653   ierr = PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);CHKERRQ(ierr);
1654   ierr = PetscFESetQuadrature(*feRef, qref);CHKERRQ(ierr);
1655   ierr = PetscQuadratureDestroy(&qref);CHKERRQ(ierr);
1656   PetscFunctionReturn(0);
1657 }
1658 
1659 /*@C
1660   PetscFECreateDefault - Create a PetscFE for basic FEM computation
1661 
1662   Collective
1663 
1664   Input Parameters:
1665 + comm      - The MPI comm
1666 . dim       - The spatial dimension
1667 . Nc        - The number of components
1668 . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
1669 . prefix    - The options prefix, or NULL
1670 - qorder    - The quadrature order or PETSC_DETERMINE to use PetscSpace polynomial degree
1671 
1672   Output Parameter:
1673 . fem - The PetscFE object
1674 
1675   Level: beginner
1676 
1677 .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
1678 @*/
1679 PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem)
1680 {
1681   PetscQuadrature q, fq;
1682   DM              K;
1683   PetscSpace      P;
1684   PetscDualSpace  Q;
1685   PetscInt        order, quadPointsPerEdge;
1686   PetscBool       tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
1687   PetscErrorCode  ierr;
1688 
1689   PetscFunctionBegin;
1690   /* Create space */
1691   ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr);
1692   ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr);
1693   ierr = PetscSpacePolynomialSetTensor(P, tensor);CHKERRQ(ierr);
1694   ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr);
1695   ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr);
1696   ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr);
1697   ierr = PetscSpaceSetUp(P);CHKERRQ(ierr);
1698   ierr = PetscSpaceGetDegree(P, &order, NULL);CHKERRQ(ierr);
1699   ierr = PetscSpacePolynomialGetTensor(P, &tensor);CHKERRQ(ierr);
1700   /* Create dual space */
1701   ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr);
1702   ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr);
1703   ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr);
1704   ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr);
1705   ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr);
1706   ierr = DMDestroy(&K);CHKERRQ(ierr);
1707   ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr);
1708   ierr = PetscDualSpaceSetOrder(Q, order);CHKERRQ(ierr);
1709   ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor);CHKERRQ(ierr);
1710   ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr);
1711   ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr);
1712   /* Create element */
1713   ierr = PetscFECreate(comm, fem);CHKERRQ(ierr);
1714   ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr);
1715   ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr);
1716   ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr);
1717   ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr);
1718   ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr);
1719   ierr = PetscFESetUp(*fem);CHKERRQ(ierr);
1720   ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr);
1721   ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr);
1722   /* Create quadrature (with specified order if given) */
1723   qorder = qorder >= 0 ? qorder : order;
1724   ierr = PetscObjectOptionsBegin((PetscObject)*fem);CHKERRQ(ierr);
1725   ierr = PetscOptionsBoundedInt("-petscfe_default_quadrature_order","Quadrature order is one less than quadrature points per edge","PetscFECreateDefault",qorder,&qorder,NULL,0);CHKERRQ(ierr);
1726   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1727   quadPointsPerEdge = PetscMax(qorder + 1,1);
1728   if (isSimplex) {
1729     ierr = PetscDTGaussJacobiQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr);
1730     ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr);
1731   } else {
1732     ierr = PetscDTGaussTensorQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr);
1733     ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr);
1734   }
1735   ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr);
1736   ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr);
1737   ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr);
1738   ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr);
1739   PetscFunctionReturn(0);
1740 }
1741 
1742 /*@C
1743   PetscFESetName - Names the FE and its subobjects
1744 
1745   Not collective
1746 
1747   Input Parameters:
1748 + fe   - The PetscFE
1749 - name - The name
1750 
1751   Level: intermediate
1752 
1753 .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
1754 @*/
1755 PetscErrorCode PetscFESetName(PetscFE fe, const char name[])
1756 {
1757   PetscSpace     P;
1758   PetscDualSpace Q;
1759   PetscErrorCode ierr;
1760 
1761   PetscFunctionBegin;
1762   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
1763   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1764   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
1765   ierr = PetscObjectSetName((PetscObject) P,  name);CHKERRQ(ierr);
1766   ierr = PetscObjectSetName((PetscObject) Q,  name);CHKERRQ(ierr);
1767   PetscFunctionReturn(0);
1768 }
1769 
1770 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[])
1771 {
1772   PetscInt       dOffset = 0, fOffset = 0, f;
1773   PetscErrorCode ierr;
1774 
1775   for (f = 0; f < Nf; ++f) {
1776     PetscFE          fe;
1777     const PetscInt   Nbf = Nb[f], Ncf = Nc[f];
1778     const PetscReal *Bq = &basisField[f][q*Nbf*Ncf];
1779     const PetscReal *Dq = &basisFieldDer[f][q*Nbf*Ncf*dim];
1780     PetscInt         b, c, d;
1781 
1782     ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);CHKERRQ(ierr);
1783     for (c = 0; c < Ncf; ++c)     u[fOffset+c] = 0.0;
1784     for (d = 0; d < dim*Ncf; ++d) u_x[fOffset*dim+d] = 0.0;
1785     for (b = 0; b < Nbf; ++b) {
1786       for (c = 0; c < Ncf; ++c) {
1787         const PetscInt cidx = b*Ncf+c;
1788 
1789         u[fOffset+c] += Bq[cidx]*coefficients[dOffset+b];
1790         for (d = 0; d < dim; ++d) u_x[(fOffset+c)*dim+d] += Dq[cidx*dim+d]*coefficients[dOffset+b];
1791       }
1792     }
1793     ierr = PetscFEPushforward(fe, fegeom, 1, &u[fOffset]);CHKERRQ(ierr);
1794     ierr = PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset*dim]);CHKERRQ(ierr);
1795     if (u_t) {
1796       for (c = 0; c < Ncf; ++c) u_t[fOffset+c] = 0.0;
1797       for (b = 0; b < Nbf; ++b) {
1798         for (c = 0; c < Ncf; ++c) {
1799           const PetscInt cidx = b*Ncf+c;
1800 
1801           u_t[fOffset+c] += Bq[cidx]*coefficients_t[dOffset+b];
1802         }
1803       }
1804       ierr = PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]);CHKERRQ(ierr);
1805     }
1806     fOffset += Ncf;
1807     dOffset += Nbf;
1808   }
1809   return 0;
1810 }
1811 
1812 PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[])
1813 {
1814   PetscFE        fe;
1815   PetscReal     *faceBasis;
1816   PetscInt       Nb, Nc, b, c;
1817   PetscErrorCode ierr;
1818 
1819   if (!prob) return 0;
1820   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1821   ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1822   ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1823   ierr = PetscFEGetFaceCentroidTabulation(fe, &faceBasis);CHKERRQ(ierr);
1824   for (c = 0; c < Nc; ++c) {u[c] = 0.0;}
1825   for (b = 0; b < Nb; ++b) {
1826     for (c = 0; c < Nc; ++c) {
1827       const PetscInt cidx = b*Nc+c;
1828 
1829       u[c] += coefficients[cidx]*faceBasis[faceLoc*Nb*Nc+cidx];
1830     }
1831   }
1832   return 0;
1833 }
1834 
1835 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[])
1836 {
1837   PetscInt       q, b, c, d;
1838   PetscErrorCode ierr;
1839 
1840   for (b = 0; b < Nb; ++b) elemVec[b] = 0.0;
1841   for (q = 0; q < Nq; ++q) {
1842     for (b = 0; b < Nb; ++b) {
1843       for (c = 0; c < Nc; ++c) {
1844         const PetscInt bcidx = b*Nc+c;
1845 
1846         tmpBasis[bcidx] = basis[q*Nb*Nc+bcidx];
1847         for (d = 0; d < dim; ++d) tmpBasisDer[bcidx*dim+d] = basisDer[q*Nb*Nc*dim+bcidx*dim+d];
1848       }
1849     }
1850     ierr = PetscFEPushforward(fe, fegeom, Nb, tmpBasis);CHKERRQ(ierr);
1851     ierr = PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer);CHKERRQ(ierr);
1852     for (b = 0; b < Nb; ++b) {
1853       for (c = 0; c < Nc; ++c) {
1854         const PetscInt bcidx = b*Nc+c;
1855         const PetscInt qcidx = q*Nc+c;
1856 
1857         elemVec[b] += tmpBasis[bcidx]*f0[qcidx];
1858         for (d = 0; d < dim; ++d) elemVec[b] += tmpBasisDer[bcidx*dim+d]*f1[qcidx*dim+d];
1859       }
1860     }
1861   }
1862   return(0);
1863 }
1864 
1865 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[])
1866 {
1867   PetscInt       f, fc, g, gc, df, dg;
1868   PetscErrorCode ierr;
1869 
1870   for (f = 0; f < NbI; ++f) {
1871     for (fc = 0; fc < NcI; ++fc) {
1872       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
1873 
1874       tmpBasisI[fidx] = basisI[fidx];
1875       for (df = 0; df < dim; ++df) tmpBasisDerI[fidx*dim+df] = basisDerI[fidx*dim+df];
1876     }
1877   }
1878   ierr = PetscFEPushforward(feI, fegeom, NbI, tmpBasisI);CHKERRQ(ierr);
1879   ierr = PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI);CHKERRQ(ierr);
1880   for (g = 0; g < NbJ; ++g) {
1881     for (gc = 0; gc < NcJ; ++gc) {
1882       const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
1883 
1884       tmpBasisJ[gidx] = basisJ[gidx];
1885       for (dg = 0; dg < dim; ++dg) tmpBasisDerJ[gidx*dim+dg] = basisDerJ[gidx*dim+dg];
1886     }
1887   }
1888   ierr = PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ);CHKERRQ(ierr);
1889   ierr = PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ);CHKERRQ(ierr);
1890   for (f = 0; f < NbI; ++f) {
1891     for (fc = 0; fc < NcI; ++fc) {
1892       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
1893       const PetscInt i    = offsetI+f; /* Element matrix row */
1894       for (g = 0; g < NbJ; ++g) {
1895         for (gc = 0; gc < NcJ; ++gc) {
1896           const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
1897           const PetscInt j    = offsetJ+g; /* Element matrix column */
1898           const PetscInt fOff = eOffset+i*totDim+j;
1899 
1900           elemMat[fOff] += tmpBasisI[fidx]*g0[fc*NcJ+gc]*tmpBasisJ[gidx];
1901           for (df = 0; df < dim; ++df) {
1902             elemMat[fOff] += tmpBasisI[fidx]*g1[(fc*NcJ+gc)*dim+df]*tmpBasisDerJ[gidx*dim+df];
1903             elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g2[(fc*NcJ+gc)*dim+df]*tmpBasisJ[gidx];
1904             for (dg = 0; dg < dim; ++dg) {
1905               elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g3[((fc*NcJ+gc)*dim+df)*dim+dg]*tmpBasisDerJ[gidx*dim+dg];
1906             }
1907           }
1908         }
1909       }
1910     }
1911   }
1912   return(0);
1913 }
1914 
1915 PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom)
1916 {
1917   PetscDualSpace  dsp;
1918   DM              dm;
1919   PetscQuadrature quadDef;
1920   PetscInt        dim, cdim, Nq;
1921   PetscErrorCode  ierr;
1922 
1923   PetscFunctionBegin;
1924   ierr = PetscFEGetDualSpace(fe, &dsp);CHKERRQ(ierr);
1925   ierr = PetscDualSpaceGetDM(dsp, &dm);CHKERRQ(ierr);
1926   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1927   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1928   ierr = PetscFEGetQuadrature(fe, &quadDef);CHKERRQ(ierr);
1929   quad = quad ? quad : quadDef;
1930   ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1931   ierr = PetscMalloc1(Nq*cdim,      &cgeom->v);CHKERRQ(ierr);
1932   ierr = PetscMalloc1(Nq*cdim*cdim, &cgeom->J);CHKERRQ(ierr);
1933   ierr = PetscMalloc1(Nq*cdim*cdim, &cgeom->invJ);CHKERRQ(ierr);
1934   ierr = PetscMalloc1(Nq,           &cgeom->detJ);CHKERRQ(ierr);
1935   cgeom->dim       = dim;
1936   cgeom->dimEmbed  = cdim;
1937   cgeom->numCells  = 1;
1938   cgeom->numPoints = Nq;
1939   ierr = DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ);CHKERRQ(ierr);
1940   PetscFunctionReturn(0);
1941 }
1942 
1943 PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom)
1944 {
1945   PetscErrorCode ierr;
1946 
1947   PetscFunctionBegin;
1948   ierr = PetscFree(cgeom->v);CHKERRQ(ierr);
1949   ierr = PetscFree(cgeom->J);CHKERRQ(ierr);
1950   ierr = PetscFree(cgeom->invJ);CHKERRQ(ierr);
1951   ierr = PetscFree(cgeom->detJ);CHKERRQ(ierr);
1952   PetscFunctionReturn(0);
1953 }
1954