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