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