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