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