xref: /petsc/src/dm/dt/fe/interface/fe.c (revision 63f3c55c12ae2f190c391f6df6c540efe018a44a)
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 /*@
701   PetscFECopyQuadrature - Copy both volumetric and surface quadrature
702 
703   Not collective
704 
705   Input Parameters:
706 + sfe - The PetscFE source for the quadratures
707 - tfe - The PetscFE target for the quadratures
708 
709   Level: intermediate
710 
711 .seealso: PetscFECreate(), PetscFESetQuadrature(), PetscFESetFaceQuadrature()
712 @*/
713 PetscErrorCode PetscFECopyQuadrature(PetscFE sfe, PetscFE tfe)
714 {
715   PetscQuadrature q;
716   PetscErrorCode  ierr;
717 
718   PetscFunctionBegin;
719   PetscValidHeaderSpecific(sfe, PETSCFE_CLASSID, 1);
720   PetscValidHeaderSpecific(tfe, PETSCFE_CLASSID, 2);
721   ierr = PetscFEGetQuadrature(sfe, &q);CHKERRQ(ierr);
722   ierr = PetscFESetQuadrature(tfe,  q);CHKERRQ(ierr);
723   ierr = PetscFEGetFaceQuadrature(sfe, &q);CHKERRQ(ierr);
724   ierr = PetscFESetFaceQuadrature(tfe,  q);CHKERRQ(ierr);
725   PetscFunctionReturn(0);
726 }
727 
728 /*@C
729   PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension
730 
731   Not collective
732 
733   Input Parameter:
734 . fem - The PetscFE object
735 
736   Output Parameter:
737 . numDof - Array with the number of dofs per dimension
738 
739   Level: intermediate
740 
741 .seealso: PetscFECreate()
742 @*/
743 PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt **numDof)
744 {
745   PetscErrorCode ierr;
746 
747   PetscFunctionBegin;
748   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
749   PetscValidPointer(numDof, 2);
750   ierr = PetscDualSpaceGetNumDof(fem->dualSpace, numDof);CHKERRQ(ierr);
751   PetscFunctionReturn(0);
752 }
753 
754 /*@C
755   PetscFEGetDefaultTabulation - Returns the tabulation of the basis functions at the quadrature points
756 
757   Not collective
758 
759   Input Parameter:
760 . fem - The PetscFE object
761 
762   Output Parameters:
763 + B - The basis function values at quadrature points
764 . D - The basis function derivatives at quadrature points
765 - H - The basis function second derivatives at quadrature points
766 
767   Note:
768 $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
769 $ D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
770 $ 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
771 
772   Level: intermediate
773 
774 .seealso: PetscFEGetTabulation(), PetscFERestoreTabulation()
775 @*/
776 PetscErrorCode PetscFEGetDefaultTabulation(PetscFE fem, PetscReal **B, PetscReal **D, PetscReal **H)
777 {
778   PetscInt         npoints;
779   const PetscReal *points;
780   PetscErrorCode   ierr;
781 
782   PetscFunctionBegin;
783   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
784   if (B) PetscValidPointer(B, 2);
785   if (D) PetscValidPointer(D, 3);
786   if (H) PetscValidPointer(H, 4);
787   ierr = PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr);
788   if (!fem->B) {ierr = PetscFEGetTabulation(fem, npoints, points, &fem->B, &fem->D, NULL/*&fem->H*/);CHKERRQ(ierr);}
789   if (B) *B = fem->B;
790   if (D) *D = fem->D;
791   if (H) *H = fem->H;
792   PetscFunctionReturn(0);
793 }
794 
795 /*@C
796   PetscFEGetFaceTabulation - Returns the tabulation of the basis functions at the face quadrature points
797 
798   Not collective
799 
800   Input Parameter:
801 . fem - The PetscFE object
802 
803   Output Parameters:
804 + B - The basis function values at face quadrature points
805 . D - The basis function derivatives at face quadrature points
806 - H - The basis function second derivatives at face quadrature points
807 
808   Note:
809 $ Bf[((f*Nq + q)*pdim + i)*Nc + c] is the value at point f,q for basis function i and component c
810 $ 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
811 $ 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
812 
813   Level: intermediate
814 
815 .seealso: PetscFEGetDefaultTabulation(), PetscFEGetTabulation(), PetscFERestoreTabulation()
816 @*/
817 PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscReal **Bf, PetscReal **Df, PetscReal **Hf)
818 {
819   PetscErrorCode   ierr;
820 
821   PetscFunctionBegin;
822   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
823   if (Bf) PetscValidPointer(Bf, 2);
824   if (Df) PetscValidPointer(Df, 3);
825   if (Hf) PetscValidPointer(Hf, 4);
826   if (!fem->Bf) {
827     const PetscReal  xi0[3] = {-1., -1., -1.};
828     PetscReal        v0[3], J[9], detJ;
829     PetscQuadrature  fq;
830     PetscDualSpace   sp;
831     DM               dm;
832     const PetscInt  *faces;
833     PetscInt         dim, numFaces, f, npoints, q;
834     const PetscReal *points;
835     PetscReal       *facePoints;
836 
837     ierr = PetscFEGetDualSpace(fem, &sp);CHKERRQ(ierr);
838     ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
839     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
840     ierr = DMPlexGetConeSize(dm, 0, &numFaces);CHKERRQ(ierr);
841     ierr = DMPlexGetCone(dm, 0, &faces);CHKERRQ(ierr);
842     ierr = PetscFEGetFaceQuadrature(fem, &fq);CHKERRQ(ierr);
843     if (fq) {
844       ierr = PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr);
845       ierr = PetscMalloc1(numFaces*npoints*dim, &facePoints);CHKERRQ(ierr);
846       for (f = 0; f < numFaces; ++f) {
847         ierr = DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ);CHKERRQ(ierr);
848         for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim-1, xi0, v0, J, &points[q*(dim-1)], &facePoints[(f*npoints+q)*dim]);
849       }
850       ierr = PetscFEGetTabulation(fem, numFaces*npoints, facePoints, &fem->Bf, &fem->Df, NULL/*&fem->Hf*/);CHKERRQ(ierr);
851       ierr = PetscFree(facePoints);CHKERRQ(ierr);
852     }
853   }
854   if (Bf) *Bf = fem->Bf;
855   if (Df) *Df = fem->Df;
856   if (Hf) *Hf = fem->Hf;
857   PetscFunctionReturn(0);
858 }
859 
860 /*@C
861   PetscFEGetFaceTabulation - Returns the tabulation of the basis functions at the face centroid points
862 
863   Not collective
864 
865   Input Parameter:
866 . fem - The PetscFE object
867 
868   Output Parameters:
869 + B - The basis function values at face centroid points
870 . D - The basis function derivatives at face centroid points
871 - H - The basis function second derivatives at face centroid points
872 
873   Note:
874 $ Bf[(f*pdim + i)*Nc + c] is the value at point f for basis function i and component c
875 $ Df[((f*pdim + i)*Nc + c)*dim + d] is the derivative value at point f for basis function i, component c, in direction d
876 $ 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
877 
878   Level: intermediate
879 
880 .seealso: PetscFEGetFaceTabulation(), PetscFEGetDefaultTabulation(), PetscFEGetTabulation(), PetscFERestoreTabulation()
881 @*/
882 PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscReal **F)
883 {
884   PetscErrorCode   ierr;
885 
886   PetscFunctionBegin;
887   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
888   PetscValidPointer(F, 2);
889   if (!fem->F) {
890     PetscDualSpace  sp;
891     DM              dm;
892     const PetscInt *cone;
893     PetscReal      *centroids;
894     PetscInt        dim, numFaces, f;
895 
896     ierr = PetscFEGetDualSpace(fem, &sp);CHKERRQ(ierr);
897     ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
898     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
899     ierr = DMPlexGetConeSize(dm, 0, &numFaces);CHKERRQ(ierr);
900     ierr = DMPlexGetCone(dm, 0, &cone);CHKERRQ(ierr);
901     ierr = PetscMalloc1(numFaces*dim, &centroids);CHKERRQ(ierr);
902     for (f = 0; f < numFaces; ++f) {ierr = DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, &centroids[f*dim], NULL);CHKERRQ(ierr);}
903     ierr = PetscFEGetTabulation(fem, numFaces, centroids, &fem->F, NULL, NULL);CHKERRQ(ierr);
904     ierr = PetscFree(centroids);CHKERRQ(ierr);
905   }
906   *F = fem->F;
907   PetscFunctionReturn(0);
908 }
909 
910 /*@C
911   PetscFEGetTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided.
912 
913   Not collective
914 
915   Input Parameters:
916 + fem     - The PetscFE object
917 . npoints - The number of tabulation points
918 - points  - The tabulation point coordinates
919 
920   Output Parameters:
921 + B - The basis function values at tabulation points
922 . D - The basis function derivatives at tabulation points
923 - H - The basis function second derivatives at tabulation points
924 
925   Note:
926 $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
927 $ D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
928 $ 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
929 
930   Level: intermediate
931 
932 .seealso: PetscFERestoreTabulation(), PetscFEGetDefaultTabulation()
933 @*/
934 PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
935 {
936   DM               dm;
937   PetscInt         pdim; /* Dimension of FE space P */
938   PetscInt         dim;  /* Spatial dimension */
939   PetscInt         comp; /* Field components */
940   PetscErrorCode   ierr;
941 
942   PetscFunctionBegin;
943   if (!npoints) {
944     if (B) *B = NULL;
945     if (D) *D = NULL;
946     if (H) *H = NULL;
947     PetscFunctionReturn(0);
948   }
949   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
950   PetscValidPointer(points, 3);
951   if (B) PetscValidPointer(B, 4);
952   if (D) PetscValidPointer(D, 5);
953   if (H) PetscValidPointer(H, 6);
954   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
955   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
956   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
957   ierr = PetscFEGetNumComponents(fem, &comp);CHKERRQ(ierr);
958   if (B) {ierr = DMGetWorkArray(dm, npoints*pdim*comp, MPIU_REAL, B);CHKERRQ(ierr);}
959   if (!dim) {
960     if (D) *D = NULL;
961     if (H) *H = NULL;
962   } else {
963     if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*comp*dim, MPIU_REAL, D);CHKERRQ(ierr);}
964     if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*comp*dim*dim, MPIU_REAL, H);CHKERRQ(ierr);}
965   }
966   ierr = (*fem->ops->gettabulation)(fem, npoints, points, B ? *B : NULL, D ? *D : NULL, H ? *H : NULL);CHKERRQ(ierr);
967   PetscFunctionReturn(0);
968 }
969 
970 /*@C
971   PetscFERestoreTabulation - Frees memory from the associated tabulation.
972 
973   Not collective
974 
975   Input Parameters:
976 + fem     - The PetscFE object
977 . npoints - The number of tabulation points
978 . points  - The tabulation point coordinates
979 . B - The basis function values at tabulation points
980 . D - The basis function derivatives at tabulation points
981 - H - The basis function second derivatives at tabulation points
982 
983   Note:
984 $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
985 $ D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
986 $ 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
987 
988   Level: intermediate
989 
990 .seealso: PetscFEGetTabulation(), PetscFEGetDefaultTabulation()
991 @*/
992 PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
993 {
994   DM             dm;
995   PetscErrorCode ierr;
996 
997   PetscFunctionBegin;
998   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
999   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
1000   if (B && *B) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, B);CHKERRQ(ierr);}
1001   if (D && *D) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, D);CHKERRQ(ierr);}
1002   if (H && *H) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, H);CHKERRQ(ierr);}
1003   PetscFunctionReturn(0);
1004 }
1005 
1006 PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
1007 {
1008   PetscSpace     bsp, bsubsp;
1009   PetscDualSpace dsp, dsubsp;
1010   PetscInt       dim, depth, numComp, i, j, coneSize, order;
1011   PetscFEType    type;
1012   DM             dm;
1013   DMLabel        label;
1014   PetscReal      *xi, *v, *J, detJ;
1015   const char     *name;
1016   PetscQuadrature origin, fullQuad, subQuad;
1017   PetscErrorCode ierr;
1018 
1019   PetscFunctionBegin;
1020   PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1);
1021   PetscValidPointer(trFE,3);
1022   ierr = PetscFEGetBasisSpace(fe,&bsp);CHKERRQ(ierr);
1023   ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr);
1024   ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr);
1025   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
1026   ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr);
1027   ierr = DMLabelGetValue(label,refPoint,&depth);CHKERRQ(ierr);
1028   ierr = PetscCalloc1(depth,&xi);CHKERRQ(ierr);
1029   ierr = PetscMalloc1(dim,&v);CHKERRQ(ierr);
1030   ierr = PetscMalloc1(dim*dim,&J);CHKERRQ(ierr);
1031   for (i = 0; i < depth; i++) xi[i] = 0.;
1032   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,&origin);CHKERRQ(ierr);
1033   ierr = PetscQuadratureSetData(origin,depth,0,1,xi,NULL);CHKERRQ(ierr);
1034   ierr = DMPlexComputeCellGeometryFEM(dm,refPoint,origin,v,J,NULL,&detJ);CHKERRQ(ierr);
1035   /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */
1036   for (i = 1; i < dim; i++) {
1037     for (j = 0; j < depth; j++) {
1038       J[i * depth + j] = J[i * dim + j];
1039     }
1040   }
1041   ierr = PetscQuadratureDestroy(&origin);CHKERRQ(ierr);
1042   ierr = PetscDualSpaceGetPointSubspace(dsp,refPoint,&dsubsp);CHKERRQ(ierr);
1043   ierr = PetscSpaceCreateSubspace(bsp,dsubsp,v,J,NULL,NULL,PETSC_OWN_POINTER,&bsubsp);CHKERRQ(ierr);
1044   ierr = PetscSpaceSetUp(bsubsp);CHKERRQ(ierr);
1045   ierr = PetscFECreate(PetscObjectComm((PetscObject)fe),trFE);CHKERRQ(ierr);
1046   ierr = PetscFEGetType(fe,&type);CHKERRQ(ierr);
1047   ierr = PetscFESetType(*trFE,type);CHKERRQ(ierr);
1048   ierr = PetscFEGetNumComponents(fe,&numComp);CHKERRQ(ierr);
1049   ierr = PetscFESetNumComponents(*trFE,numComp);CHKERRQ(ierr);
1050   ierr = PetscFESetBasisSpace(*trFE,bsubsp);CHKERRQ(ierr);
1051   ierr = PetscFESetDualSpace(*trFE,dsubsp);CHKERRQ(ierr);
1052   ierr = PetscObjectGetName((PetscObject) fe, &name);CHKERRQ(ierr);
1053   if (name) {ierr = PetscFESetName(*trFE, name);CHKERRQ(ierr);}
1054   ierr = PetscFEGetQuadrature(fe,&fullQuad);CHKERRQ(ierr);
1055   ierr = PetscQuadratureGetOrder(fullQuad,&order);CHKERRQ(ierr);
1056   ierr = DMPlexGetConeSize(dm,refPoint,&coneSize);CHKERRQ(ierr);
1057   if (coneSize == 2 * depth) {
1058     ierr = PetscDTGaussTensorQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr);
1059   } else {
1060     ierr = PetscDTGaussJacobiQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr);
1061   }
1062   ierr = PetscFESetQuadrature(*trFE,subQuad);CHKERRQ(ierr);
1063   ierr = PetscFESetUp(*trFE);CHKERRQ(ierr);
1064   ierr = PetscQuadratureDestroy(&subQuad);CHKERRQ(ierr);
1065   ierr = PetscSpaceDestroy(&bsubsp);CHKERRQ(ierr);
1066   PetscFunctionReturn(0);
1067 }
1068 
1069 PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE)
1070 {
1071   PetscInt       hStart, hEnd;
1072   PetscDualSpace dsp;
1073   DM             dm;
1074   PetscErrorCode ierr;
1075 
1076   PetscFunctionBegin;
1077   PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1);
1078   PetscValidPointer(trFE,3);
1079   *trFE = NULL;
1080   ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr);
1081   ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr);
1082   ierr = DMPlexGetHeightStratum(dm,height,&hStart,&hEnd);CHKERRQ(ierr);
1083   if (hEnd <= hStart) PetscFunctionReturn(0);
1084   ierr = PetscFECreatePointTrace(fe,hStart,trFE);CHKERRQ(ierr);
1085   PetscFunctionReturn(0);
1086 }
1087 
1088 
1089 /*@
1090   PetscFEGetDimension - Get the dimension of the finite element space on a cell
1091 
1092   Not collective
1093 
1094   Input Parameter:
1095 . fe - The PetscFE
1096 
1097   Output Parameter:
1098 . dim - The dimension
1099 
1100   Level: intermediate
1101 
1102 .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
1103 @*/
1104 PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
1105 {
1106   PetscErrorCode ierr;
1107 
1108   PetscFunctionBegin;
1109   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1110   PetscValidPointer(dim, 2);
1111   if (fem->ops->getdimension) {ierr = (*fem->ops->getdimension)(fem, dim);CHKERRQ(ierr);}
1112   PetscFunctionReturn(0);
1113 }
1114 
1115 /*@C
1116   PetscFEPushforward - Map the reference element function to real space
1117 
1118   Input Parameters:
1119 + fe     - The PetscFE
1120 . fegeom - The cell geometry
1121 . Nv     - The number of function values
1122 - vals   - The function values
1123 
1124   Output Parameter:
1125 . vals   - The transformed function values
1126 
1127   Level: advanced
1128 
1129   Note: This just forwards the call onto PetscDualSpacePushforward().
1130 
1131 .seealso: PetscDualSpacePushforward()
1132 @*/
1133 PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1134 {
1135   PetscErrorCode ierr;
1136 
1137   PetscFunctionBeginHot;
1138   ierr = PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);CHKERRQ(ierr);
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 /*@C
1143   PetscFEPushforwardGradient - Map the reference element function gradient to real space
1144 
1145   Input Parameters:
1146 + fe     - The PetscFE
1147 . fegeom - The cell geometry
1148 . Nv     - The number of function gradient values
1149 - vals   - The function gradient values
1150 
1151   Output Parameter:
1152 . vals   - The transformed function gradient values
1153 
1154   Level: advanced
1155 
1156   Note: This just forwards the call onto PetscDualSpacePushforwardGradient().
1157 
1158 .seealso: PetscFEPushforward(), PetscDualSpacePushforwardGradient(), PetscDualSpacePushforward()
1159 @*/
1160 PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1161 {
1162   PetscErrorCode ierr;
1163 
1164   PetscFunctionBeginHot;
1165   ierr = PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);CHKERRQ(ierr);
1166   PetscFunctionReturn(0);
1167 }
1168 
1169 /*
1170 Purpose: Compute element vector for chunk of elements
1171 
1172 Input:
1173   Sizes:
1174      Ne:  number of elements
1175      Nf:  number of fields
1176      PetscFE
1177        dim: spatial dimension
1178        Nb:  number of basis functions
1179        Nc:  number of field components
1180        PetscQuadrature
1181          Nq:  number of quadrature points
1182 
1183   Geometry:
1184      PetscFEGeom[Ne] possibly *Nq
1185        PetscReal v0s[dim]
1186        PetscReal n[dim]
1187        PetscReal jacobians[dim*dim]
1188        PetscReal jacobianInverses[dim*dim]
1189        PetscReal jacobianDeterminants
1190   FEM:
1191      PetscFE
1192        PetscQuadrature
1193          PetscReal   quadPoints[Nq*dim]
1194          PetscReal   quadWeights[Nq]
1195        PetscReal   basis[Nq*Nb*Nc]
1196        PetscReal   basisDer[Nq*Nb*Nc*dim]
1197      PetscScalar coefficients[Ne*Nb*Nc]
1198      PetscScalar elemVec[Ne*Nb*Nc]
1199 
1200   Problem:
1201      PetscInt f: the active field
1202      f0, f1
1203 
1204   Work Space:
1205      PetscFE
1206        PetscScalar f0[Nq*dim];
1207        PetscScalar f1[Nq*dim*dim];
1208        PetscScalar u[Nc];
1209        PetscScalar gradU[Nc*dim];
1210        PetscReal   x[dim];
1211        PetscScalar realSpaceDer[dim];
1212 
1213 Purpose: Compute element vector for N_cb batches of elements
1214 
1215 Input:
1216   Sizes:
1217      N_cb: Number of serial cell batches
1218 
1219   Geometry:
1220      PetscReal v0s[Ne*dim]
1221      PetscReal jacobians[Ne*dim*dim]        possibly *Nq
1222      PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq
1223      PetscReal jacobianDeterminants[Ne]     possibly *Nq
1224   FEM:
1225      static PetscReal   quadPoints[Nq*dim]
1226      static PetscReal   quadWeights[Nq]
1227      static PetscReal   basis[Nq*Nb*Nc]
1228      static PetscReal   basisDer[Nq*Nb*Nc*dim]
1229      PetscScalar coefficients[Ne*Nb*Nc]
1230      PetscScalar elemVec[Ne*Nb*Nc]
1231 
1232 ex62.c:
1233   PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[],
1234                                                const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[],
1235                                                void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]),
1236                                                void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[])
1237 
1238 ex52.c:
1239   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)
1240   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)
1241 
1242 ex52_integrateElement.cu
1243 __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec)
1244 
1245 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[],
1246                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1247                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)
1248 
1249 ex52_integrateElementOpenCL.c:
1250 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[],
1251                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
1252                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)
1253 
1254 __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec)
1255 */
1256 
1257 /*@C
1258   PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration
1259 
1260   Not collective
1261 
1262   Input Parameters:
1263 + fem          - The PetscFE object for the field being integrated
1264 . prob         - The PetscDS specifying the discretizations and continuum functions
1265 . field        - The field being integrated
1266 . Ne           - The number of elements in the chunk
1267 . cgeom        - The cell geometry for each cell in the chunk
1268 . coefficients - The array of FEM basis coefficients for the elements
1269 . probAux      - The PetscDS specifying the auxiliary discretizations
1270 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1271 
1272   Output Parameter
1273 . integral     - the integral for this field
1274 
1275   Level: intermediate
1276 
1277 .seealso: PetscFEIntegrateResidual()
1278 @*/
1279 PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1280                                 const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1281 {
1282   PetscFE        fe;
1283   PetscErrorCode ierr;
1284 
1285   PetscFunctionBegin;
1286   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1287   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1288   if (fe->ops->integrate) {ierr = (*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);}
1289   PetscFunctionReturn(0);
1290 }
1291 
1292 /*@C
1293   PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration
1294 
1295   Not collective
1296 
1297   Input Parameters:
1298 + fem          - The PetscFE object for the field being integrated
1299 . prob         - The PetscDS specifying the discretizations and continuum functions
1300 . field        - The field being integrated
1301 . obj_func     - The function to be integrated
1302 . Ne           - The number of elements in the chunk
1303 . fgeom        - The face geometry for each face in the chunk
1304 . coefficients - The array of FEM basis coefficients for the elements
1305 . probAux      - The PetscDS specifying the auxiliary discretizations
1306 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1307 
1308   Output Parameter
1309 . integral     - the integral for this field
1310 
1311   Level: intermediate
1312 
1313 .seealso: PetscFEIntegrateResidual()
1314 @*/
1315 PetscErrorCode PetscFEIntegrateBd(PetscDS prob, PetscInt field,
1316                                   void (*obj_func)(PetscInt, PetscInt, PetscInt,
1317                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1318                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1319                                                    PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1320                                   PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1321 {
1322   PetscFE        fe;
1323   PetscErrorCode ierr;
1324 
1325   PetscFunctionBegin;
1326   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1327   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1328   if (fe->ops->integratebd) {ierr = (*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);}
1329   PetscFunctionReturn(0);
1330 }
1331 
1332 /*@C
1333   PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration
1334 
1335   Not collective
1336 
1337   Input Parameters:
1338 + fem          - The PetscFE object for the field being integrated
1339 . prob         - The PetscDS specifying the discretizations and continuum functions
1340 . field        - The field being integrated
1341 . Ne           - The number of elements in the chunk
1342 . cgeom        - The cell geometry for each cell in the chunk
1343 . coefficients - The array of FEM basis coefficients for the elements
1344 . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1345 . probAux      - The PetscDS specifying the auxiliary discretizations
1346 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1347 - t            - The time
1348 
1349   Output Parameter
1350 . elemVec      - the element residual vectors from each element
1351 
1352   Note:
1353 $ Loop over batch of elements (e):
1354 $   Loop over quadrature points (q):
1355 $     Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q
1356 $     Call f_0 and f_1
1357 $   Loop over element vector entries (f,fc --> i):
1358 $     elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u)
1359 
1360   Level: intermediate
1361 
1362 .seealso: PetscFEIntegrateResidual()
1363 @*/
1364 PetscErrorCode PetscFEIntegrateResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1365                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1366 {
1367   PetscFE        fe;
1368   PetscErrorCode ierr;
1369 
1370   PetscFunctionBegin;
1371   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1372   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1373   if (fe->ops->integrateresidual) {ierr = (*fe->ops->integrateresidual)(prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);}
1374   PetscFunctionReturn(0);
1375 }
1376 
1377 /*@C
1378   PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary
1379 
1380   Not collective
1381 
1382   Input Parameters:
1383 + fem          - The PetscFE object for the field being integrated
1384 . prob         - The PetscDS specifying the discretizations and continuum functions
1385 . field        - The field being integrated
1386 . Ne           - The number of elements in the chunk
1387 . fgeom        - The face geometry for each cell in the chunk
1388 . coefficients - The array of FEM basis coefficients for the elements
1389 . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1390 . probAux      - The PetscDS specifying the auxiliary discretizations
1391 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1392 - t            - The time
1393 
1394   Output Parameter
1395 . elemVec      - the element residual vectors from each element
1396 
1397   Level: intermediate
1398 
1399 .seealso: PetscFEIntegrateResidual()
1400 @*/
1401 PetscErrorCode PetscFEIntegrateBdResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom,
1402                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
1403 {
1404   PetscFE        fe;
1405   PetscErrorCode ierr;
1406 
1407   PetscFunctionBegin;
1408   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1409   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1410   if (fe->ops->integratebdresidual) {ierr = (*fe->ops->integratebdresidual)(prob, field, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);}
1411   PetscFunctionReturn(0);
1412 }
1413 
1414 /*@C
1415   PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration
1416 
1417   Not collective
1418 
1419   Input Parameters:
1420 + fem          - The PetscFE object for the field being integrated
1421 . prob         - The PetscDS specifying the discretizations and continuum functions
1422 . jtype        - The type of matrix pointwise functions that should be used
1423 . fieldI       - The test field being integrated
1424 . fieldJ       - The basis field being integrated
1425 . Ne           - The number of elements in the chunk
1426 . cgeom        - The cell geometry for each cell in the chunk
1427 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1428 . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1429 . probAux      - The PetscDS specifying the auxiliary discretizations
1430 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1431 . t            - The time
1432 - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)
1433 
1434   Output Parameter
1435 . elemMat      - the element matrices for the Jacobian from each element
1436 
1437   Note:
1438 $ Loop over batch of elements (e):
1439 $   Loop over element matrix entries (f,fc,g,gc --> i,j):
1440 $     Loop over quadrature points (q):
1441 $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1442 $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1443 $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1444 $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1445 $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1446   Level: intermediate
1447 
1448 .seealso: PetscFEIntegrateResidual()
1449 @*/
1450 PetscErrorCode PetscFEIntegrateJacobian(PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom,
1451                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1452 {
1453   PetscFE        fe;
1454   PetscErrorCode ierr;
1455 
1456   PetscFunctionBegin;
1457   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1458   ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
1459   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);}
1460   PetscFunctionReturn(0);
1461 }
1462 
1463 /*@C
1464   PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration
1465 
1466   Not collective
1467 
1468   Input Parameters:
1469 . prob         - The PetscDS specifying the discretizations and continuum functions
1470 . fieldI       - The test field being integrated
1471 . fieldJ       - The basis field being integrated
1472 . Ne           - The number of elements in the chunk
1473 . fgeom        - The face geometry for each cell in the chunk
1474 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
1475 . coefficients_t - The array of FEM basis time derivative coefficients for the elements
1476 . probAux      - The PetscDS specifying the auxiliary discretizations
1477 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1478 . t            - The time
1479 - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)
1480 
1481   Output Parameter
1482 . elemMat              - the element matrices for the Jacobian from each element
1483 
1484   Note:
1485 $ Loop over batch of elements (e):
1486 $   Loop over element matrix entries (f,fc,g,gc --> i,j):
1487 $     Loop over quadrature points (q):
1488 $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
1489 $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
1490 $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1491 $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
1492 $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
1493   Level: intermediate
1494 
1495 .seealso: PetscFEIntegrateJacobian(), PetscFEIntegrateResidual()
1496 @*/
1497 PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom,
1498                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1499 {
1500   PetscFE        fe;
1501   PetscErrorCode ierr;
1502 
1503   PetscFunctionBegin;
1504   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1505   ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
1506   if (fe->ops->integratebdjacobian) {ierr = (*fe->ops->integratebdjacobian)(prob, fieldI, fieldJ, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);CHKERRQ(ierr);}
1507   PetscFunctionReturn(0);
1508 }
1509 
1510 /*@
1511   PetscFEGetHeightSubspace - Get the subspace of this space for a mesh point of a given height
1512 
1513   Input Parameters:
1514 + fe     - The finite element space
1515 - height - The height of the Plex point
1516 
1517   Output Parameter:
1518 . subfe  - The subspace of this FE space
1519 
1520   Note: For example, if we want the subspace of this space for a face, we would choose height = 1.
1521 
1522   Level: advanced
1523 
1524 .seealso: PetscFECreateDefault()
1525 @*/
1526 PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe)
1527 {
1528   PetscSpace      P, subP;
1529   PetscDualSpace  Q, subQ;
1530   PetscQuadrature subq;
1531   PetscFEType     fetype;
1532   PetscInt        dim, Nc;
1533   PetscErrorCode  ierr;
1534 
1535   PetscFunctionBegin;
1536   PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1);
1537   PetscValidPointer(subfe, 3);
1538   if (height == 0) {
1539     *subfe = fe;
1540     PetscFunctionReturn(0);
1541   }
1542   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
1543   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1544   ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1545   ierr = PetscFEGetFaceQuadrature(fe, &subq);CHKERRQ(ierr);
1546   ierr = PetscDualSpaceGetDimension(Q, &dim);CHKERRQ(ierr);
1547   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);}
1548   if (!fe->subspaces) {ierr = PetscCalloc1(dim, &fe->subspaces);CHKERRQ(ierr);}
1549   if (height <= dim) {
1550     if (!fe->subspaces[height-1]) {
1551       PetscFE     sub;
1552       const char *name;
1553 
1554       ierr = PetscSpaceGetHeightSubspace(P, height, &subP);CHKERRQ(ierr);
1555       ierr = PetscDualSpaceGetHeightSubspace(Q, height, &subQ);CHKERRQ(ierr);
1556       ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), &sub);CHKERRQ(ierr);
1557       ierr = PetscObjectGetName((PetscObject) fe,  &name);CHKERRQ(ierr);
1558       ierr = PetscObjectSetName((PetscObject) sub,  name);CHKERRQ(ierr);
1559       ierr = PetscFEGetType(fe, &fetype);CHKERRQ(ierr);
1560       ierr = PetscFESetType(sub, fetype);CHKERRQ(ierr);
1561       ierr = PetscFESetBasisSpace(sub, subP);CHKERRQ(ierr);
1562       ierr = PetscFESetDualSpace(sub, subQ);CHKERRQ(ierr);
1563       ierr = PetscFESetNumComponents(sub, Nc);CHKERRQ(ierr);
1564       ierr = PetscFESetUp(sub);CHKERRQ(ierr);
1565       ierr = PetscFESetQuadrature(sub, subq);CHKERRQ(ierr);
1566       fe->subspaces[height-1] = sub;
1567     }
1568     *subfe = fe->subspaces[height-1];
1569   } else {
1570     *subfe = NULL;
1571   }
1572   PetscFunctionReturn(0);
1573 }
1574 
1575 /*@
1576   PetscFERefine - Create a "refined" PetscFE object that refines the reference cell into smaller copies. This is typically used
1577   to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more
1578   sparsity). It is also used to create an interpolation between regularly refined meshes.
1579 
1580   Collective on fem
1581 
1582   Input Parameter:
1583 . fe - The initial PetscFE
1584 
1585   Output Parameter:
1586 . feRef - The refined PetscFE
1587 
1588   Level: advanced
1589 
1590 .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
1591 @*/
1592 PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)
1593 {
1594   PetscSpace       P, Pref;
1595   PetscDualSpace   Q, Qref;
1596   DM               K, Kref;
1597   PetscQuadrature  q, qref;
1598   const PetscReal *v0, *jac;
1599   PetscInt         numComp, numSubelements;
1600   PetscErrorCode   ierr;
1601 
1602   PetscFunctionBegin;
1603   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
1604   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1605   ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1606   ierr = PetscDualSpaceGetDM(Q, &K);CHKERRQ(ierr);
1607   /* Create space */
1608   ierr = PetscObjectReference((PetscObject) P);CHKERRQ(ierr);
1609   Pref = P;
1610   /* Create dual space */
1611   ierr = PetscDualSpaceDuplicate(Q, &Qref);CHKERRQ(ierr);
1612   ierr = DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref);CHKERRQ(ierr);
1613   ierr = PetscDualSpaceSetDM(Qref, Kref);CHKERRQ(ierr);
1614   ierr = DMDestroy(&Kref);CHKERRQ(ierr);
1615   ierr = PetscDualSpaceSetUp(Qref);CHKERRQ(ierr);
1616   /* Create element */
1617   ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), feRef);CHKERRQ(ierr);
1618   ierr = PetscFESetType(*feRef, PETSCFECOMPOSITE);CHKERRQ(ierr);
1619   ierr = PetscFESetBasisSpace(*feRef, Pref);CHKERRQ(ierr);
1620   ierr = PetscFESetDualSpace(*feRef, Qref);CHKERRQ(ierr);
1621   ierr = PetscFEGetNumComponents(fe,    &numComp);CHKERRQ(ierr);
1622   ierr = PetscFESetNumComponents(*feRef, numComp);CHKERRQ(ierr);
1623   ierr = PetscFESetUp(*feRef);CHKERRQ(ierr);
1624   ierr = PetscSpaceDestroy(&Pref);CHKERRQ(ierr);
1625   ierr = PetscDualSpaceDestroy(&Qref);CHKERRQ(ierr);
1626   /* Create quadrature */
1627   ierr = PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL);CHKERRQ(ierr);
1628   ierr = PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);CHKERRQ(ierr);
1629   ierr = PetscFESetQuadrature(*feRef, qref);CHKERRQ(ierr);
1630   ierr = PetscQuadratureDestroy(&qref);CHKERRQ(ierr);
1631   PetscFunctionReturn(0);
1632 }
1633 
1634 /*@C
1635   PetscFECreateDefault - Create a PetscFE for basic FEM computation
1636 
1637   Collective
1638 
1639   Input Parameters:
1640 + comm      - The MPI comm
1641 . dim       - The spatial dimension
1642 . Nc        - The number of components
1643 . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
1644 . prefix    - The options prefix, or NULL
1645 - qorder    - The quadrature order or PETSC_DETERMINE to use PetscSpace polynomial degree
1646 
1647   Output Parameter:
1648 . fem - The PetscFE object
1649 
1650   Level: beginner
1651 
1652 .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
1653 @*/
1654 PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem)
1655 {
1656   PetscQuadrature q, fq;
1657   DM              K;
1658   PetscSpace      P;
1659   PetscDualSpace  Q;
1660   PetscInt        order, quadPointsPerEdge;
1661   PetscBool       tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
1662   PetscErrorCode  ierr;
1663 
1664   PetscFunctionBegin;
1665   /* Create space */
1666   ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr);
1667   ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr);
1668   ierr = PetscSpacePolynomialSetTensor(P, tensor);CHKERRQ(ierr);
1669   ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr);
1670   ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr);
1671   ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr);
1672   ierr = PetscSpaceSetUp(P);CHKERRQ(ierr);
1673   ierr = PetscSpaceGetDegree(P, &order, NULL);CHKERRQ(ierr);
1674   ierr = PetscSpacePolynomialGetTensor(P, &tensor);CHKERRQ(ierr);
1675   /* Create dual space */
1676   ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr);
1677   ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr);
1678   ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr);
1679   ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr);
1680   ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr);
1681   ierr = DMDestroy(&K);CHKERRQ(ierr);
1682   ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr);
1683   ierr = PetscDualSpaceSetOrder(Q, order);CHKERRQ(ierr);
1684   ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor);CHKERRQ(ierr);
1685   ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr);
1686   ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr);
1687   /* Create element */
1688   ierr = PetscFECreate(comm, fem);CHKERRQ(ierr);
1689   ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr);
1690   ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr);
1691   ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr);
1692   ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr);
1693   ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr);
1694   ierr = PetscFESetUp(*fem);CHKERRQ(ierr);
1695   ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr);
1696   ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr);
1697   /* Create quadrature (with specified order if given) */
1698   qorder = qorder >= 0 ? qorder : order;
1699   ierr = PetscObjectOptionsBegin((PetscObject)*fem);CHKERRQ(ierr);
1700   ierr = PetscOptionsBoundedInt("-petscfe_default_quadrature_order","Quadrature order is one less than quadrature points per edge","PetscFECreateDefault",qorder,&qorder,NULL,0);CHKERRQ(ierr);
1701   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1702   quadPointsPerEdge = PetscMax(qorder + 1,1);
1703   if (isSimplex) {
1704     ierr = PetscDTGaussJacobiQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr);
1705     ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr);
1706   } else {
1707     ierr = PetscDTGaussTensorQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr);
1708     ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr);
1709   }
1710   ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr);
1711   ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr);
1712   ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr);
1713   ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr);
1714   PetscFunctionReturn(0);
1715 }
1716 
1717 /*@C
1718   PetscFESetName - Names the FE and its subobjects
1719 
1720   Not collective
1721 
1722   Input Parameters:
1723 + fe   - The PetscFE
1724 - name - The name
1725 
1726   Level: intermediate
1727 
1728 .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
1729 @*/
1730 PetscErrorCode PetscFESetName(PetscFE fe, const char name[])
1731 {
1732   PetscSpace     P;
1733   PetscDualSpace Q;
1734   PetscErrorCode ierr;
1735 
1736   PetscFunctionBegin;
1737   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
1738   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1739   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
1740   ierr = PetscObjectSetName((PetscObject) P,  name);CHKERRQ(ierr);
1741   ierr = PetscObjectSetName((PetscObject) Q,  name);CHKERRQ(ierr);
1742   PetscFunctionReturn(0);
1743 }
1744 
1745 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[])
1746 {
1747   PetscInt       dOffset = 0, fOffset = 0, f;
1748   PetscErrorCode ierr;
1749 
1750   for (f = 0; f < Nf; ++f) {
1751     PetscFE          fe;
1752     const PetscInt   Nbf = Nb[f], Ncf = Nc[f];
1753     const PetscReal *Bq = &basisField[f][q*Nbf*Ncf];
1754     const PetscReal *Dq = &basisFieldDer[f][q*Nbf*Ncf*dim];
1755     PetscInt         b, c, d;
1756 
1757     ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);CHKERRQ(ierr);
1758     for (c = 0; c < Ncf; ++c)     u[fOffset+c] = 0.0;
1759     for (d = 0; d < dim*Ncf; ++d) u_x[fOffset*dim+d] = 0.0;
1760     for (b = 0; b < Nbf; ++b) {
1761       for (c = 0; c < Ncf; ++c) {
1762         const PetscInt cidx = b*Ncf+c;
1763 
1764         u[fOffset+c] += Bq[cidx]*coefficients[dOffset+b];
1765         for (d = 0; d < dim; ++d) u_x[(fOffset+c)*dim+d] += Dq[cidx*dim+d]*coefficients[dOffset+b];
1766       }
1767     }
1768     ierr = PetscFEPushforward(fe, fegeom, 1, &u[fOffset]);CHKERRQ(ierr);
1769     ierr = PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset*dim]);CHKERRQ(ierr);
1770     if (u_t) {
1771       for (c = 0; c < Ncf; ++c) u_t[fOffset+c] = 0.0;
1772       for (b = 0; b < Nbf; ++b) {
1773         for (c = 0; c < Ncf; ++c) {
1774           const PetscInt cidx = b*Ncf+c;
1775 
1776           u_t[fOffset+c] += Bq[cidx]*coefficients_t[dOffset+b];
1777         }
1778       }
1779       ierr = PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]);CHKERRQ(ierr);
1780     }
1781     fOffset += Ncf;
1782     dOffset += Nbf;
1783   }
1784   return 0;
1785 }
1786 
1787 PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[])
1788 {
1789   PetscFE        fe;
1790   PetscReal     *faceBasis;
1791   PetscInt       Nb, Nc, b, c;
1792   PetscErrorCode ierr;
1793 
1794   if (!prob) return 0;
1795   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1796   ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1797   ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1798   ierr = PetscFEGetFaceCentroidTabulation(fe, &faceBasis);CHKERRQ(ierr);
1799   for (c = 0; c < Nc; ++c) {u[c] = 0.0;}
1800   for (b = 0; b < Nb; ++b) {
1801     for (c = 0; c < Nc; ++c) {
1802       const PetscInt cidx = b*Nc+c;
1803 
1804       u[c] += coefficients[cidx]*faceBasis[faceLoc*Nb*Nc+cidx];
1805     }
1806   }
1807   return 0;
1808 }
1809 
1810 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[])
1811 {
1812   PetscInt       q, b, c, d;
1813   PetscErrorCode ierr;
1814 
1815   for (b = 0; b < Nb; ++b) elemVec[b] = 0.0;
1816   for (q = 0; q < Nq; ++q) {
1817     for (b = 0; b < Nb; ++b) {
1818       for (c = 0; c < Nc; ++c) {
1819         const PetscInt bcidx = b*Nc+c;
1820 
1821         tmpBasis[bcidx] = basis[q*Nb*Nc+bcidx];
1822         for (d = 0; d < dim; ++d) tmpBasisDer[bcidx*dim+d] = basisDer[q*Nb*Nc*dim+bcidx*dim+d];
1823       }
1824     }
1825     ierr = PetscFEPushforward(fe, fegeom, Nb, tmpBasis);CHKERRQ(ierr);
1826     ierr = PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer);CHKERRQ(ierr);
1827     for (b = 0; b < Nb; ++b) {
1828       for (c = 0; c < Nc; ++c) {
1829         const PetscInt bcidx = b*Nc+c;
1830         const PetscInt qcidx = q*Nc+c;
1831 
1832         elemVec[b] += tmpBasis[bcidx]*f0[qcidx];
1833         for (d = 0; d < dim; ++d) elemVec[b] += tmpBasisDer[bcidx*dim+d]*f1[qcidx*dim+d];
1834       }
1835     }
1836   }
1837   return(0);
1838 }
1839 
1840 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[])
1841 {
1842   PetscInt       f, fc, g, gc, df, dg;
1843   PetscErrorCode ierr;
1844 
1845   for (f = 0; f < NbI; ++f) {
1846     for (fc = 0; fc < NcI; ++fc) {
1847       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
1848 
1849       tmpBasisI[fidx] = basisI[fidx];
1850       for (df = 0; df < dim; ++df) tmpBasisDerI[fidx*dim+df] = basisDerI[fidx*dim+df];
1851     }
1852   }
1853   ierr = PetscFEPushforward(feI, fegeom, NbI, tmpBasisI);CHKERRQ(ierr);
1854   ierr = PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI);CHKERRQ(ierr);
1855   for (g = 0; g < NbJ; ++g) {
1856     for (gc = 0; gc < NcJ; ++gc) {
1857       const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
1858 
1859       tmpBasisJ[gidx] = basisJ[gidx];
1860       for (dg = 0; dg < dim; ++dg) tmpBasisDerJ[gidx*dim+dg] = basisDerJ[gidx*dim+dg];
1861     }
1862   }
1863   ierr = PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ);CHKERRQ(ierr);
1864   ierr = PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ);CHKERRQ(ierr);
1865   for (f = 0; f < NbI; ++f) {
1866     for (fc = 0; fc < NcI; ++fc) {
1867       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
1868       const PetscInt i    = offsetI+f; /* Element matrix row */
1869       for (g = 0; g < NbJ; ++g) {
1870         for (gc = 0; gc < NcJ; ++gc) {
1871           const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
1872           const PetscInt j    = offsetJ+g; /* Element matrix column */
1873           const PetscInt fOff = eOffset+i*totDim+j;
1874 
1875           elemMat[fOff] += tmpBasisI[fidx]*g0[fc*NcJ+gc]*tmpBasisJ[gidx];
1876           for (df = 0; df < dim; ++df) {
1877             elemMat[fOff] += tmpBasisI[fidx]*g1[(fc*NcJ+gc)*dim+df]*tmpBasisDerJ[gidx*dim+df];
1878             elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g2[(fc*NcJ+gc)*dim+df]*tmpBasisJ[gidx];
1879             for (dg = 0; dg < dim; ++dg) {
1880               elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g3[((fc*NcJ+gc)*dim+df)*dim+dg]*tmpBasisDerJ[gidx*dim+dg];
1881             }
1882           }
1883         }
1884       }
1885     }
1886   }
1887   return(0);
1888 }
1889 
1890 PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom)
1891 {
1892   PetscDualSpace  dsp;
1893   DM              dm;
1894   PetscQuadrature quadDef;
1895   PetscInt        dim, cdim, Nq;
1896   PetscErrorCode  ierr;
1897 
1898   PetscFunctionBegin;
1899   ierr = PetscFEGetDualSpace(fe, &dsp);CHKERRQ(ierr);
1900   ierr = PetscDualSpaceGetDM(dsp, &dm);CHKERRQ(ierr);
1901   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1902   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1903   ierr = PetscFEGetQuadrature(fe, &quadDef);CHKERRQ(ierr);
1904   quad = quad ? quad : quadDef;
1905   ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1906   ierr = PetscMalloc1(Nq*cdim,      &cgeom->v);CHKERRQ(ierr);
1907   ierr = PetscMalloc1(Nq*cdim*cdim, &cgeom->J);CHKERRQ(ierr);
1908   ierr = PetscMalloc1(Nq*cdim*cdim, &cgeom->invJ);CHKERRQ(ierr);
1909   ierr = PetscMalloc1(Nq,           &cgeom->detJ);CHKERRQ(ierr);
1910   cgeom->dim       = dim;
1911   cgeom->dimEmbed  = cdim;
1912   cgeom->numCells  = 1;
1913   cgeom->numPoints = Nq;
1914   ierr = DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ);CHKERRQ(ierr);
1915   PetscFunctionReturn(0);
1916 }
1917 
1918 PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom)
1919 {
1920   PetscErrorCode ierr;
1921 
1922   PetscFunctionBegin;
1923   ierr = PetscFree(cgeom->v);CHKERRQ(ierr);
1924   ierr = PetscFree(cgeom->J);CHKERRQ(ierr);
1925   ierr = PetscFree(cgeom->invJ);CHKERRQ(ierr);
1926   ierr = PetscFree(cgeom->detJ);CHKERRQ(ierr);
1927   PetscFunctionReturn(0);
1928 }
1929