Actual source code: dualspace.c

  1: #include <petsc/private/petscfeimpl.h>
  2: #include <petscdmplex.h>

  4: PetscClassId PETSCDUALSPACE_CLASSID = 0;

  6: PetscLogEvent PETSCDUALSPACE_SetUp;

  8: PetscFunctionList PetscDualSpaceList              = NULL;
  9: PetscBool         PetscDualSpaceRegisterAllCalled = PETSC_FALSE;

 11: /*
 12:   PetscDualSpaceLatticePointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to at most 'max'.
 13:                                                      Ordering is lexicographic with lowest index as least significant in ordering.
 14:                                                      e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {0,2}.

 16:   Input Parameters:
 17: + len - The length of the tuple
 18: . max - The maximum sum
 19: - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition

 21:   Output Parameter:
 22: . tup - A tuple of `len` integers whose sum is at most `max`

 24:   Level: developer

 26: .seealso: `PetscDualSpaceType`, `PetscDualSpaceTensorPointLexicographic_Internal()`
 27: */
 28: PetscErrorCode PetscDualSpaceLatticePointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
 29: {
 30:   PetscFunctionBegin;
 31:   while (len--) {
 32:     max -= tup[len];
 33:     if (!max) {
 34:       tup[len] = 0;
 35:       break;
 36:     }
 37:   }
 38:   tup[++len]++;
 39:   PetscFunctionReturn(PETSC_SUCCESS);
 40: }

 42: /*
 43:   PetscDualSpaceTensorPointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max'.
 44:                                                     Ordering is lexicographic with lowest index as least significant in ordering.
 45:                                                     e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,1}, {0,2}, {1,2}, {2,2}.

 47:   Input Parameters:
 48: + len - The length of the tuple
 49: . max - The maximum value
 50: - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition

 52:   Output Parameter:
 53: . tup - A tuple of `len` integers whose entries are at most `max`

 55:   Level: developer

 57: .seealso: `PetscDualSpaceType`, `PetscDualSpaceLatticePointLexicographic_Internal()`
 58: */
 59: PetscErrorCode PetscDualSpaceTensorPointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
 60: {
 61:   PetscInt i;

 63:   PetscFunctionBegin;
 64:   for (i = 0; i < len; i++) {
 65:     if (tup[i] < max) {
 66:       break;
 67:     } else {
 68:       tup[i] = 0;
 69:     }
 70:   }
 71:   tup[i]++;
 72:   PetscFunctionReturn(PETSC_SUCCESS);
 73: }

 75: /*@C
 76:   PetscDualSpaceRegister - Adds a new `PetscDualSpaceType`

 78:   Not Collective, No Fortran Support

 80:   Input Parameters:
 81: + sname    - The name of a new user-defined creation routine
 82: - function - The creation routine

 84:   Example Usage:
 85: .vb
 86:     PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate);
 87: .ve

 89:   Then, your PetscDualSpace type can be chosen with the procedural interface via
 90: .vb
 91:     PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
 92:     PetscDualSpaceSetType(PetscDualSpace, "my_dual_space");
 93: .ve
 94:   or at runtime via the option
 95: .vb
 96:     -petscdualspace_type my_dual_space
 97: .ve

 99:   Level: advanced

101:   Note:
102:   `PetscDualSpaceRegister()` may be called multiple times to add several user-defined `PetscDualSpace`

104: .seealso: `PetscDualSpace`, `PetscDualSpaceType`, `PetscDualSpaceRegisterAll()`, `PetscDualSpaceRegisterDestroy()`
105: @*/
106: PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace))
107: {
108:   PetscFunctionBegin;
109:   PetscCall(PetscFunctionListAdd(&PetscDualSpaceList, sname, function));
110:   PetscFunctionReturn(PETSC_SUCCESS);
111: }

113: /*@C
114:   PetscDualSpaceSetType - Builds a particular `PetscDualSpace` based on its `PetscDualSpaceType`

116:   Collective

118:   Input Parameters:
119: + sp   - The `PetscDualSpace` object
120: - name - The kind of space

122:   Options Database Key:
123: . -petscdualspace_type type - Sets the `PetscDualSpace` type; see `PetscDualSpaceType` for the choices

125:   Level: intermediate

127: .seealso: `PetscDualSpace`, `PetscDualSpaceType`, `PetscDualSpaceGetType()`, `PetscDualSpaceCreate()`
128: @*/
129: PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name)
130: {
131:   PetscErrorCode (*r)(PetscDualSpace);
132:   PetscBool match;

134:   PetscFunctionBegin;
136:   PetscCall(PetscObjectTypeCompare((PetscObject)sp, name, &match));
137:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

139:   if (!PetscDualSpaceRegisterAllCalled) PetscCall(PetscDualSpaceRegisterAll());
140:   PetscCall(PetscFunctionListFind(PetscDualSpaceList, name, &r));
141:   PetscCheck(r, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name);

143:   PetscTryTypeMethod(sp, destroy);
144:   sp->ops->destroy = NULL;

146:   PetscCall((*r)(sp));
147:   PetscCall(PetscObjectChangeTypeName((PetscObject)sp, name));
148:   PetscFunctionReturn(PETSC_SUCCESS);
149: }

151: /*@C
152:   PetscDualSpaceGetType - Gets the `PetscDualSpaceType` name (as a string) from the object.

154:   Not Collective

156:   Input Parameter:
157: . sp - The `PetscDualSpace`

159:   Output Parameter:
160: . name - The `PetscDualSpaceType` name

162:   Level: intermediate

164: .seealso: `PetscDualSpace`, `PetscDualSpaceType`, `PetscDualSpaceSetType()`, `PetscDualSpaceCreate()`
165: @*/
166: PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name)
167: {
168:   PetscFunctionBegin;
170:   PetscAssertPointer(name, 2);
171:   if (!PetscDualSpaceRegisterAllCalled) PetscCall(PetscDualSpaceRegisterAll());
172:   *name = ((PetscObject)sp)->type_name;
173:   PetscFunctionReturn(PETSC_SUCCESS);
174: }

176: static PetscErrorCode PetscDualSpaceView_ASCII(PetscDualSpace sp, PetscViewer v)
177: {
178:   PetscViewerFormat format;
179:   PetscInt          pdim, f;

181:   PetscFunctionBegin;
182:   PetscCall(PetscDualSpaceGetDimension(sp, &pdim));
183:   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)sp, v));
184:   PetscCall(PetscViewerASCIIPushTab(v));
185:   if (sp->k != 0 && sp->k != PETSC_FORM_DEGREE_UNDEFINED) {
186:     PetscCall(PetscViewerASCIIPrintf(v, "Dual space for %" PetscInt_FMT "-forms %swith %" PetscInt_FMT " components, size %" PetscInt_FMT "\n", PetscAbsInt(sp->k), sp->k < 0 ? "(stored in dual form) " : "", sp->Nc, pdim));
187:   } else {
188:     PetscCall(PetscViewerASCIIPrintf(v, "Dual space with %" PetscInt_FMT " components, size %" PetscInt_FMT "\n", sp->Nc, pdim));
189:   }
190:   PetscTryTypeMethod(sp, view, v);
191:   PetscCall(PetscViewerGetFormat(v, &format));
192:   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
193:     PetscCall(PetscViewerASCIIPushTab(v));
194:     for (f = 0; f < pdim; ++f) {
195:       PetscCall(PetscViewerASCIIPrintf(v, "Dual basis vector %" PetscInt_FMT "\n", f));
196:       PetscCall(PetscViewerASCIIPushTab(v));
197:       PetscCall(PetscQuadratureView(sp->functional[f], v));
198:       PetscCall(PetscViewerASCIIPopTab(v));
199:     }
200:     PetscCall(PetscViewerASCIIPopTab(v));
201:   }
202:   PetscCall(PetscViewerASCIIPopTab(v));
203:   PetscFunctionReturn(PETSC_SUCCESS);
204: }

206: /*@C
207:   PetscDualSpaceViewFromOptions - View a `PetscDualSpace` based on values in the options database

209:   Collective

211:   Input Parameters:
212: + A    - the `PetscDualSpace` object
213: . obj  - Optional object, provides the options prefix
214: - name - command line option name

216:   Level: intermediate

218:   Note:
219:   See `PetscObjectViewFromOptions()` for possible command line values

221: .seealso: `PetscDualSpace`, `PetscDualSpaceView()`, `PetscObjectViewFromOptions()`, `PetscDualSpaceCreate()`
222: @*/
223: PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace A, PeOp PetscObject obj, const char name[])
224: {
225:   PetscFunctionBegin;
227:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
228:   PetscFunctionReturn(PETSC_SUCCESS);
229: }

231: /*@C
232:   PetscDualSpaceView - Views a `PetscDualSpace`

234:   Collective

236:   Input Parameters:
237: + sp - the `PetscDualSpace` object to view
238: - v  - the viewer

240:   Level: beginner

242: .seealso: `PetscViewer`, `PetscDualSpaceDestroy()`, `PetscDualSpace`
243: @*/
244: PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
245: {
246:   PetscBool isascii;

248:   PetscFunctionBegin;
251:   if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sp), &v));
252:   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &isascii));
253:   if (isascii) PetscCall(PetscDualSpaceView_ASCII(sp, v));
254:   PetscFunctionReturn(PETSC_SUCCESS);
255: }

257: /*@C
258:   PetscDualSpaceSetFromOptions - sets parameters in a `PetscDualSpace` from the options database

260:   Collective

262:   Input Parameter:
263: . sp - the `PetscDualSpace` object to set options for

265:   Options Database Keys:
266: + -petscdualspace_order order                          - the approximation order of the space
267: . -petscdualspace_form_degree deg                      - the form degree, say 0 for point evaluations, or 2 for area integrals
268: . -petscdualspace_components c                         - the number of components, say d for a vector field
269: . -petscdualspace_refcell celltype                     - Reference cell type name
270: . -petscdualspace_lagrange_continuity (true|false)     - Flag for continuous element
271: . -petscdualspace_lagrange_tensor (true|false)         - Flag for tensor dual space
272: . -petscdualspace_lagrange_trimmed (true|false)        - Flag for trimmed dual space
273: . -petscdualspace_lagrange_node_type nodetype          - Lagrange node location type
274: . -petscdualspace_lagrange_node_endpoints (true|false) - Flag for nodes that include endpoints
275: . -petscdualspace_lagrange_node_exponent exponent      - Gauss-Jacobi weight function exponent
276: . -petscdualspace_lagrange_use_moments (true|false)    - Use moments (where appropriate) for functionals
277: - -petscdualspace_lagrange_moment_order order          - Quadrature order for moment functionals

279:   Level: intermediate

281: .seealso: `PetscDualSpaceView()`, `PetscDualSpace`, `PetscObjectSetFromOptions()`
282: @*/
283: PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
284: {
285:   DMPolytopeType refCell = DM_POLYTOPE_TRIANGLE;
286:   const char    *defaultType;
287:   char           name[256];
288:   PetscBool      flg;

290:   PetscFunctionBegin;
292:   if (!((PetscObject)sp)->type_name) defaultType = PETSCDUALSPACELAGRANGE;
293:   else defaultType = ((PetscObject)sp)->type_name;
294:   if (!PetscSpaceRegisterAllCalled) PetscCall(PetscSpaceRegisterAll());

296:   PetscObjectOptionsBegin((PetscObject)sp);
297:   PetscCall(PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg));
298:   if (flg) PetscCall(PetscDualSpaceSetType(sp, name));
299:   else if (!((PetscObject)sp)->type_name) PetscCall(PetscDualSpaceSetType(sp, defaultType));
300:   PetscCall(PetscOptionsBoundedInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL, 0));
301:   PetscCall(PetscOptionsInt("-petscdualspace_form_degree", "The form degree of the dofs", "PetscDualSpaceSetFormDegree", sp->k, &sp->k, NULL));
302:   PetscCall(PetscOptionsBoundedInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL, 1));
303:   PetscTryTypeMethod(sp, setfromoptions, PetscOptionsObject);
304:   PetscCall(PetscOptionsEnum("-petscdualspace_refcell", "Reference cell shape", "PetscDualSpaceSetReferenceCell", DMPolytopeTypes, (PetscEnum)refCell, (PetscEnum *)&refCell, &flg));
305:   if (flg) {
306:     DM K;

308:     PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, refCell, &K));
309:     PetscCall(PetscDualSpaceSetDM(sp, K));
310:     PetscCall(DMDestroy(&K));
311:   }

313:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
314:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)sp, PetscOptionsObject));
315:   PetscOptionsEnd();
316:   sp->setfromoptionscalled = PETSC_TRUE;
317:   PetscFunctionReturn(PETSC_SUCCESS);
318: }

320: /*@C
321:   PetscDualSpaceSetUp - Construct a basis for a `PetscDualSpace`

323:   Collective

325:   Input Parameter:
326: . sp - the `PetscDualSpace` object to setup

328:   Level: intermediate

330: .seealso: `PetscDualSpaceView()`, `PetscDualSpaceDestroy()`, `PetscDualSpace`
331: @*/
332: PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
333: {
334:   PetscFunctionBegin;
336:   if (sp->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
337:   PetscCall(PetscLogEventBegin(PETSCDUALSPACE_SetUp, sp, 0, 0, 0));
338:   sp->setupcalled = PETSC_TRUE;
339:   PetscTryTypeMethod(sp, setup);
340:   PetscCall(PetscLogEventEnd(PETSCDUALSPACE_SetUp, sp, 0, 0, 0));
341:   if (sp->setfromoptionscalled) PetscCall(PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view"));
342:   PetscFunctionReturn(PETSC_SUCCESS);
343: }

345: static PetscErrorCode PetscDualSpaceClearDMData_Internal(PetscDualSpace sp, DM dm)
346: {
347:   PetscInt pStart = -1, pEnd = -1, depth = -1;

349:   PetscFunctionBegin;
350:   if (!dm) PetscFunctionReturn(PETSC_SUCCESS);
351:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
352:   PetscCall(DMPlexGetDepth(dm, &depth));

354:   if (sp->pointSpaces) {
355:     PetscInt i;

357:     for (i = 0; i < pEnd - pStart; i++) PetscCall(PetscDualSpaceDestroy(&sp->pointSpaces[i]));
358:   }
359:   PetscCall(PetscFree(sp->pointSpaces));

361:   if (sp->heightSpaces) {
362:     PetscInt i;

364:     for (i = 0; i <= depth; i++) PetscCall(PetscDualSpaceDestroy(&sp->heightSpaces[i]));
365:   }
366:   PetscCall(PetscFree(sp->heightSpaces));

368:   PetscCall(PetscSectionDestroy(&sp->pointSection));
369:   PetscCall(PetscSectionDestroy(&sp->intPointSection));
370:   PetscCall(PetscQuadratureDestroy(&sp->intNodes));
371:   PetscCall(VecDestroy(&sp->intDofValues));
372:   PetscCall(VecDestroy(&sp->intNodeValues));
373:   PetscCall(MatDestroy(&sp->intMat));
374:   PetscCall(PetscQuadratureDestroy(&sp->allNodes));
375:   PetscCall(VecDestroy(&sp->allDofValues));
376:   PetscCall(VecDestroy(&sp->allNodeValues));
377:   PetscCall(MatDestroy(&sp->allMat));
378:   PetscCall(PetscFree(sp->numDof));
379:   PetscFunctionReturn(PETSC_SUCCESS);
380: }

382: /*@C
383:   PetscDualSpaceDestroy - Destroys a `PetscDualSpace` object

385:   Collective

387:   Input Parameter:
388: . sp - the `PetscDualSpace` object to destroy

390:   Level: beginner

392: .seealso: `PetscDualSpace`, `PetscDualSpaceView()`, `PetscDualSpace()`, `PetscDualSpaceCreate()`
393: @*/
394: PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
395: {
396:   PetscInt dim, f;
397:   DM       dm;

399:   PetscFunctionBegin;
400:   if (!*sp) PetscFunctionReturn(PETSC_SUCCESS);

403:   if (--((PetscObject)*sp)->refct > 0) {
404:     *sp = NULL;
405:     PetscFunctionReturn(PETSC_SUCCESS);
406:   }
407:   ((PetscObject)*sp)->refct = 0;

409:   PetscCall(PetscDualSpaceGetDimension(*sp, &dim));
410:   dm = (*sp)->dm;

412:   PetscTryTypeMethod(*sp, destroy);
413:   PetscCall(PetscDualSpaceClearDMData_Internal(*sp, dm));

415:   for (f = 0; f < dim; ++f) PetscCall(PetscQuadratureDestroy(&(*sp)->functional[f]));
416:   PetscCall(PetscFree((*sp)->functional));
417:   PetscCall(DMDestroy(&(*sp)->dm));
418:   PetscCall(PetscHeaderDestroy(sp));
419:   PetscFunctionReturn(PETSC_SUCCESS);
420: }

422: /*@C
423:   PetscDualSpaceCreate - Creates an empty `PetscDualSpace` object. The type can then be set with `PetscDualSpaceSetType()`.

425:   Collective

427:   Input Parameter:
428: . comm - The communicator for the `PetscDualSpace` object

430:   Output Parameter:
431: . sp - The `PetscDualSpace` object

433:   Level: beginner

435: .seealso: `PetscDualSpace`, `PetscDualSpaceSetType()`, `PETSCDUALSPACELAGRANGE`
436: @*/
437: PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
438: {
439:   PetscDualSpace s;

441:   PetscFunctionBegin;
442:   PetscAssertPointer(sp, 2);
443:   PetscCall(PetscCitationsRegister(FECitation, &FEcite));
444:   PetscCall(PetscFEInitializePackage());

446:   PetscCall(PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView));
447:   s->order       = 0;
448:   s->Nc          = 1;
449:   s->k           = 0;
450:   s->spdim       = -1;
451:   s->spintdim    = -1;
452:   s->uniform     = PETSC_TRUE;
453:   s->setupcalled = PETSC_FALSE;
454:   *sp            = s;
455:   PetscFunctionReturn(PETSC_SUCCESS);
456: }

458: /*@C
459:   PetscDualSpaceDuplicate - Creates a duplicate `PetscDualSpace` object that is not setup.

461:   Collective

463:   Input Parameter:
464: . sp - The original `PetscDualSpace`

466:   Output Parameter:
467: . spNew - The duplicate `PetscDualSpace`

469:   Level: beginner

471: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`
472: @*/
473: PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
474: {
475:   DM                 dm;
476:   PetscDualSpaceType type;
477:   const char        *name;

479:   PetscFunctionBegin;
481:   PetscAssertPointer(spNew, 2);
482:   PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)sp), spNew));
483:   name = ((PetscObject)sp)->name;
484:   if (name) PetscCall(PetscObjectSetName((PetscObject)*spNew, name));
485:   PetscCall(PetscDualSpaceGetType(sp, &type));
486:   PetscCall(PetscDualSpaceSetType(*spNew, type));
487:   PetscCall(PetscDualSpaceGetDM(sp, &dm));
488:   PetscCall(PetscDualSpaceSetDM(*spNew, dm));

490:   (*spNew)->order   = sp->order;
491:   (*spNew)->k       = sp->k;
492:   (*spNew)->Nc      = sp->Nc;
493:   (*spNew)->uniform = sp->uniform;
494:   PetscTryTypeMethod(sp, duplicate, *spNew);
495:   PetscFunctionReturn(PETSC_SUCCESS);
496: }

498: /*@C
499:   PetscDualSpaceGetDM - Get the `DM` representing the reference cell of a `PetscDualSpace`

501:   Not Collective

503:   Input Parameter:
504: . sp - The `PetscDualSpace`

506:   Output Parameter:
507: . dm - The reference cell, that is a `DM` that consists of a single cell

509:   Level: intermediate

511: .seealso: `PetscDualSpace`, `PetscDualSpaceSetDM()`, `PetscDualSpaceCreate()`
512: @*/
513: PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
514: {
515:   PetscFunctionBegin;
517:   PetscAssertPointer(dm, 2);
518:   *dm = sp->dm;
519:   PetscFunctionReturn(PETSC_SUCCESS);
520: }

522: /*@C
523:   PetscDualSpaceSetDM - Get the `DM` representing the reference cell

525:   Not Collective

527:   Input Parameters:
528: + sp - The `PetscDual`Space
529: - dm - The reference cell

531:   Level: intermediate

533: .seealso: `PetscDualSpace`, `DM`, `PetscDualSpaceGetDM()`, `PetscDualSpaceCreate()`
534: @*/
535: PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
536: {
537:   PetscFunctionBegin;
540:   PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change DM after dualspace is set up");
541:   PetscCall(PetscObjectReference((PetscObject)dm));
542:   if (sp->dm && sp->dm != dm) PetscCall(PetscDualSpaceClearDMData_Internal(sp, sp->dm));
543:   PetscCall(DMDestroy(&sp->dm));
544:   sp->dm = dm;
545:   PetscFunctionReturn(PETSC_SUCCESS);
546: }

548: /*@C
549:   PetscDualSpaceGetOrder - Get the order of the dual space

551:   Not Collective

553:   Input Parameter:
554: . sp - The `PetscDualSpace`

556:   Output Parameter:
557: . order - The order

559:   Level: intermediate

561: .seealso: `PetscDualSpace`, `PetscDualSpaceSetOrder()`, `PetscDualSpaceCreate()`
562: @*/
563: PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
564: {
565:   PetscFunctionBegin;
567:   PetscAssertPointer(order, 2);
568:   *order = sp->order;
569:   PetscFunctionReturn(PETSC_SUCCESS);
570: }

572: /*@C
573:   PetscDualSpaceSetOrder - Set the order of the dual space

575:   Not Collective

577:   Input Parameters:
578: + sp    - The `PetscDualSpace`
579: - order - The order

581:   Level: intermediate

583: .seealso: `PetscDualSpace`, `PetscDualSpaceGetOrder()`, `PetscDualSpaceCreate()`
584: @*/
585: PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
586: {
587:   PetscFunctionBegin;
589:   PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change order after dualspace is set up");
590:   sp->order = order;
591:   PetscFunctionReturn(PETSC_SUCCESS);
592: }

594: /*@C
595:   PetscDualSpaceGetNumComponents - Return the number of components for this space

597:   Input Parameter:
598: . sp - The `PetscDualSpace`

600:   Output Parameter:
601: . Nc - The number of components

603:   Level: intermediate

605:   Note:
606:   A vector space, for example, will have d components, where d is the spatial dimension

608: .seealso: `PetscDualSpaceSetNumComponents()`, `PetscDualSpaceGetDimension()`, `PetscDualSpaceCreate()`, `PetscDualSpace`
609: @*/
610: PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
611: {
612:   PetscFunctionBegin;
614:   PetscAssertPointer(Nc, 2);
615:   *Nc = sp->Nc;
616:   PetscFunctionReturn(PETSC_SUCCESS);
617: }

619: /*@C
620:   PetscDualSpaceSetNumComponents - Set the number of components for this space

622:   Input Parameters:
623: + sp - The `PetscDualSpace`
624: - Nc - The number of components

626:   Level: intermediate

628: .seealso: `PetscDualSpaceGetNumComponents()`, `PetscDualSpaceCreate()`, `PetscDualSpace`
629: @*/
630: PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
631: {
632:   PetscFunctionBegin;
634:   PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up");
635:   sp->Nc = Nc;
636:   PetscFunctionReturn(PETSC_SUCCESS);
637: }

639: /*@C
640:   PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space

642:   Not Collective

644:   Input Parameters:
645: + sp - The `PetscDualSpace`
646: - i  - The basis number

648:   Output Parameter:
649: . functional - The basis functional

651:   Level: intermediate

653: .seealso: `PetscDualSpace`, `PetscQuadrature`, `PetscDualSpaceGetDimension()`, `PetscDualSpaceCreate()`
654: @*/
655: PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
656: {
657:   PetscInt dim;

659:   PetscFunctionBegin;
661:   PetscAssertPointer(functional, 3);
662:   PetscCall(PetscDualSpaceGetDimension(sp, &dim));
663:   PetscCheck(!(i < 0) && !(i >= dim), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", i, dim);
664:   *functional = sp->functional[i];
665:   PetscFunctionReturn(PETSC_SUCCESS);
666: }

668: /*@C
669:   PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals

671:   Not Collective

673:   Input Parameter:
674: . sp - The `PetscDualSpace`

676:   Output Parameter:
677: . dim - The dimension

679:   Level: intermediate

681: .seealso: `PetscDualSpace`, `PetscDualSpaceGetFunctional()`, `PetscDualSpaceCreate()`
682: @*/
683: PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
684: {
685:   PetscFunctionBegin;
687:   PetscAssertPointer(dim, 2);
688:   if (sp->spdim < 0) {
689:     PetscSection section;

691:     PetscCall(PetscDualSpaceGetSection(sp, &section));
692:     if (section) PetscCall(PetscSectionGetStorageSize(section, &sp->spdim));
693:     else sp->spdim = 0;
694:   }
695:   *dim = sp->spdim;
696:   PetscFunctionReturn(PETSC_SUCCESS);
697: }

699: /*@C
700:   PetscDualSpaceGetInteriorDimension - Get the interior dimension of the dual space, i.e. the number of basis functionals assigned to the interior of the reference domain

702:   Not Collective

704:   Input Parameter:
705: . sp - The `PetscDualSpace`

707:   Output Parameter:
708: . intdim - The dimension

710:   Level: intermediate

712: .seealso: `PetscDualSpace`, `PetscDualSpaceGetFunctional()`, `PetscDualSpaceCreate()`
713: @*/
714: PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace sp, PetscInt *intdim)
715: {
716:   PetscFunctionBegin;
718:   PetscAssertPointer(intdim, 2);
719:   if (sp->spintdim < 0) {
720:     PetscSection section;

722:     PetscCall(PetscDualSpaceGetSection(sp, &section));
723:     if (section) PetscCall(PetscSectionGetConstrainedStorageSize(section, &sp->spintdim));
724:     else sp->spintdim = 0;
725:   }
726:   *intdim = sp->spintdim;
727:   PetscFunctionReturn(PETSC_SUCCESS);
728: }

730: /*@C
731:   PetscDualSpaceGetUniform - Whether this dual space is uniform

733:   Not Collective

735:   Input Parameter:
736: . sp - A dual space

738:   Output Parameter:
739: . uniform - `PETSC_TRUE` if (a) the dual space is the same for each point in a stratum of the reference `DMPLEX`, and
740:              (b) every symmetry of each point in the reference `DMPLEX` is also a symmetry of the point's dual space.

742:   Level: advanced

744:   Note:
745:   All of the usual spaces on simplex or tensor-product elements will be uniform, only reference cells
746:   with non-uniform strata (like trianguar-prisms) or anisotropic hp dual spaces will not be uniform.

748: .seealso: `PetscDualSpace`, `PetscDualSpaceGetPointSubspace()`, `PetscDualSpaceGetSymmetries()`
749: @*/
750: PetscErrorCode PetscDualSpaceGetUniform(PetscDualSpace sp, PetscBool *uniform)
751: {
752:   PetscFunctionBegin;
754:   PetscAssertPointer(uniform, 2);
755:   *uniform = sp->uniform;
756:   PetscFunctionReturn(PETSC_SUCCESS);
757: }

759: /*@CC
760:   PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension

762:   Not Collective

764:   Input Parameter:
765: . sp - The `PetscDualSpace`

767:   Output Parameter:
768: . numDof - An array of length dim+1 which holds the number of dofs for each dimension

770:   Level: intermediate

772:   Note:
773:   Do not free `numDof`

775: .seealso: `PetscDualSpace`, `PetscDualSpaceGetFunctional()`, `PetscDualSpaceCreate()`
776: @*/
777: PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt *numDof[])
778: {
779:   PetscFunctionBegin;
781:   PetscAssertPointer(numDof, 2);
782:   PetscCheck(sp->uniform, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A non-uniform space does not have a fixed number of dofs for each height");
783:   if (!sp->numDof) {
784:     DM           dm;
785:     PetscInt     depth, d;
786:     PetscSection section;

788:     PetscCall(PetscDualSpaceGetDM(sp, &dm));
789:     PetscCall(DMPlexGetDepth(dm, &depth));
790:     PetscCall(PetscCalloc1(depth + 1, &sp->numDof));
791:     PetscCall(PetscDualSpaceGetSection(sp, &section));
792:     for (d = 0; d <= depth; d++) {
793:       PetscInt dStart, dEnd;

795:       PetscCall(DMPlexGetDepthStratum(dm, d, &dStart, &dEnd));
796:       if (dEnd <= dStart) continue;
797:       PetscCall(PetscSectionGetDof(section, dStart, &sp->numDof[d]));
798:     }
799:   }
800:   *numDof = sp->numDof;
801:   PetscCheck(*numDof, PetscObjectComm((PetscObject)sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
802:   PetscFunctionReturn(PETSC_SUCCESS);
803: }

805: /* create the section of the right size and set a permutation for topological ordering */
806: PetscErrorCode PetscDualSpaceSectionCreate_Internal(PetscDualSpace sp, PetscSection *topSection)
807: {
808:   DM           dm;
809:   PetscInt     pStart, pEnd, cStart, cEnd, c, depth, count, i;
810:   PetscInt    *seen, *perm;
811:   PetscSection section;

813:   PetscFunctionBegin;
814:   dm = sp->dm;
815:   PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &section));
816:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
817:   PetscCall(PetscSectionSetChart(section, pStart, pEnd));
818:   PetscCall(PetscCalloc1(pEnd - pStart, &seen));
819:   PetscCall(PetscMalloc1(pEnd - pStart, &perm));
820:   PetscCall(DMPlexGetDepth(dm, &depth));
821:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
822:   for (c = cStart, count = 0; c < cEnd; c++) {
823:     PetscInt  closureSize = -1, e;
824:     PetscInt *closure     = NULL;

826:     perm[count++]    = c;
827:     seen[c - pStart] = 1;
828:     PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
829:     for (e = 0; e < closureSize; e++) {
830:       PetscInt point = closure[2 * e];

832:       if (seen[point - pStart]) continue;
833:       perm[count++]        = point;
834:       seen[point - pStart] = 1;
835:     }
836:     PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
837:   }
838:   PetscCheck(count == pEnd - pStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad topological ordering");
839:   for (i = 0; i < pEnd - pStart; i++)
840:     if (perm[i] != i) break;
841:   if (i < pEnd - pStart) {
842:     IS permIS;

844:     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, perm, PETSC_OWN_POINTER, &permIS));
845:     PetscCall(ISSetPermutation(permIS));
846:     PetscCall(PetscSectionSetPermutation(section, permIS));
847:     PetscCall(ISDestroy(&permIS));
848:   } else {
849:     PetscCall(PetscFree(perm));
850:   }
851:   PetscCall(PetscFree(seen));
852:   *topSection = section;
853:   PetscFunctionReturn(PETSC_SUCCESS);
854: }

856: /* mark boundary points and set up */
857: PetscErrorCode PetscDualSpaceSectionSetUp_Internal(PetscDualSpace sp, PetscSection section)
858: {
859:   DM       dm;
860:   DMLabel  boundary;
861:   PetscInt pStart, pEnd, p;

863:   PetscFunctionBegin;
864:   dm = sp->dm;
865:   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "boundary", &boundary));
866:   PetscCall(PetscDualSpaceGetDM(sp, &dm));
867:   PetscCall(DMPlexMarkBoundaryFaces(dm, 1, boundary));
868:   PetscCall(DMPlexLabelComplete(dm, boundary));
869:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
870:   for (p = pStart; p < pEnd; p++) {
871:     PetscInt bval;

873:     PetscCall(DMLabelGetValue(boundary, p, &bval));
874:     if (bval == 1) {
875:       PetscInt dof;

877:       PetscCall(PetscSectionGetDof(section, p, &dof));
878:       PetscCall(PetscSectionSetConstraintDof(section, p, dof));
879:     }
880:   }
881:   PetscCall(DMLabelDestroy(&boundary));
882:   PetscCall(PetscSectionSetUp(section));
883:   PetscFunctionReturn(PETSC_SUCCESS);
884: }

886: /*@C
887:   PetscDualSpaceGetSection - Create a `PetscSection` over the reference cell with the layout from this space

889:   Collective

891:   Input Parameter:
892: . sp - The `PetscDualSpace`

894:   Output Parameter:
895: . section - The section

897:   Level: advanced

899: .seealso: `PetscDualSpace`, `PetscSection`, `PetscDualSpaceCreate()`, `DMPLEX`
900: @*/
901: PetscErrorCode PetscDualSpaceGetSection(PetscDualSpace sp, PetscSection *section)
902: {
903:   PetscInt pStart, pEnd, p;

905:   PetscFunctionBegin;
906:   if (!sp->dm) {
907:     *section = NULL;
908:     PetscFunctionReturn(PETSC_SUCCESS);
909:   }
910:   if (!sp->pointSection) {
911:     /* mark the boundary */
912:     PetscCall(PetscDualSpaceSectionCreate_Internal(sp, &sp->pointSection));
913:     PetscCall(DMPlexGetChart(sp->dm, &pStart, &pEnd));
914:     for (p = pStart; p < pEnd; p++) {
915:       PetscDualSpace psp;

917:       PetscCall(PetscDualSpaceGetPointSubspace(sp, p, &psp));
918:       if (psp) {
919:         PetscInt dof;

921:         PetscCall(PetscDualSpaceGetInteriorDimension(psp, &dof));
922:         PetscCall(PetscSectionSetDof(sp->pointSection, p, dof));
923:       }
924:     }
925:     PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, sp->pointSection));
926:   }
927:   *section = sp->pointSection;
928:   PetscFunctionReturn(PETSC_SUCCESS);
929: }

931: /*@C
932:   PetscDualSpaceGetInteriorSection - Create a `PetscSection` over the reference cell with the layout from this space
933:   for interior degrees of freedom

935:   Collective

937:   Input Parameter:
938: . sp - The `PetscDualSpace`

940:   Output Parameter:
941: . section - The interior section

943:   Level: advanced

945:   Note:
946:   Most reference domains have one cell, in which case the only cell will have
947:   all of the interior degrees of freedom in the interior section.  But
948:   for `PETSCDUALSPACEREFINED` there may be other mesh points in the interior,
949:   and this section describes their layout.

951: .seealso: `PetscDualSpace`, `PetscSection`, `PetscDualSpaceCreate()`, `DMPLEX`
952: @*/
953: PetscErrorCode PetscDualSpaceGetInteriorSection(PetscDualSpace sp, PetscSection *section)
954: {
955:   PetscInt pStart, pEnd, p;

957:   PetscFunctionBegin;
958:   if (!sp->dm) {
959:     *section = NULL;
960:     PetscFunctionReturn(PETSC_SUCCESS);
961:   }
962:   if (!sp->intPointSection) {
963:     PetscSection full_section;

965:     PetscCall(PetscDualSpaceGetSection(sp, &full_section));
966:     PetscCall(PetscDualSpaceSectionCreate_Internal(sp, &sp->intPointSection));
967:     PetscCall(PetscSectionGetChart(full_section, &pStart, &pEnd));
968:     for (p = pStart; p < pEnd; p++) {
969:       PetscInt dof, cdof;

971:       PetscCall(PetscSectionGetDof(full_section, p, &dof));
972:       PetscCall(PetscSectionGetConstraintDof(full_section, p, &cdof));
973:       PetscCall(PetscSectionSetDof(sp->intPointSection, p, dof - cdof));
974:     }
975:     PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, sp->intPointSection));
976:   }
977:   *section = sp->intPointSection;
978:   PetscFunctionReturn(PETSC_SUCCESS);
979: }

981: /* this assumes that all of the point dual spaces store their interior dofs first, which is true when the point DMs
982:  * have one cell */
983: PetscErrorCode PetscDualSpacePushForwardSubspaces_Internal(PetscDualSpace sp, PetscInt sStart, PetscInt sEnd)
984: {
985:   PetscReal   *sv0, *v0, *J;
986:   PetscSection section;
987:   PetscInt     dim, s, k;
988:   DM           dm;

990:   PetscFunctionBegin;
991:   PetscCall(PetscDualSpaceGetDM(sp, &dm));
992:   PetscCall(DMGetDimension(dm, &dim));
993:   PetscCall(PetscDualSpaceGetSection(sp, &section));
994:   PetscCall(PetscMalloc3(dim, &v0, dim, &sv0, dim * dim, &J));
995:   PetscCall(PetscDualSpaceGetFormDegree(sp, &k));
996:   for (s = sStart; s < sEnd; s++) {
997:     PetscReal      detJ, hdetJ;
998:     PetscDualSpace ssp;
999:     PetscInt       dof, off, f, sdim;
1000:     PetscInt       i, j;
1001:     DM             sdm;

1003:     PetscCall(PetscDualSpaceGetPointSubspace(sp, s, &ssp));
1004:     if (!ssp) continue;
1005:     PetscCall(PetscSectionGetDof(section, s, &dof));
1006:     PetscCall(PetscSectionGetOffset(section, s, &off));
1007:     /* get the first vertex of the reference cell */
1008:     PetscCall(PetscDualSpaceGetDM(ssp, &sdm));
1009:     PetscCall(DMGetDimension(sdm, &sdim));
1010:     PetscCall(DMPlexComputeCellGeometryAffineFEM(sdm, 0, sv0, NULL, NULL, &hdetJ));
1011:     PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, s, v0, J, NULL, &detJ));
1012:     /* compactify Jacobian */
1013:     for (i = 0; i < dim; i++)
1014:       for (j = 0; j < sdim; j++) J[i * sdim + j] = J[i * dim + j];
1015:     for (f = 0; f < dof; f++) {
1016:       PetscQuadrature fn;

1018:       PetscCall(PetscDualSpaceGetFunctional(ssp, f, &fn));
1019:       PetscCall(PetscQuadraturePushForward(fn, dim, sv0, v0, J, k, &sp->functional[off + f]));
1020:     }
1021:   }
1022:   PetscCall(PetscFree3(v0, sv0, J));
1023:   PetscFunctionReturn(PETSC_SUCCESS);
1024: }

1026: /*@C
1027:   PetscDualSpaceApply - Apply a functional from the dual space basis to an input function

1029:   Input Parameters:
1030: + sp      - The `PetscDualSpace` object
1031: . f       - The basis functional index
1032: . time    - The time
1033: . cgeom   - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian) (or evaluated at the coordinates of the functional)
1034: . numComp - The number of components for the function
1035: . func    - The input function
1036: - ctx     - A context for the function

1038:   Output Parameter:
1039: . value - numComp output values

1041:   Calling sequence:
1042: .vb
1043:   PetscErrorCode func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt numComponents, PetscScalar values[], PetscCtx ctx)
1044: .ve

1046:   Level: beginner

1048: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1049: @*/
1050: PetscErrorCode PetscDualSpaceApply(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFEGeom *cgeom, PetscInt numComp, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), PetscCtx ctx, PetscScalar *value)
1051: {
1052:   PetscFunctionBegin;
1054:   PetscAssertPointer(cgeom, 4);
1055:   PetscAssertPointer(value, 8);
1056:   PetscUseTypeMethod(sp, apply, f, time, cgeom, numComp, func, ctx, value);
1057:   PetscFunctionReturn(PETSC_SUCCESS);
1058: }

1060: /*@C
1061:   PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by `PetscDualSpaceGetAllData()`

1063:   Input Parameters:
1064: + sp        - The `PetscDualSpace` object
1065: - pointEval - Evaluation at the points returned by `PetscDualSpaceGetAllData()`

1067:   Output Parameter:
1068: . spValue - The values of all dual space functionals

1070:   Level: advanced

1072: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1073: @*/
1074: PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1075: {
1076:   PetscFunctionBegin;
1078:   PetscUseTypeMethod(sp, applyall, pointEval, spValue);
1079:   PetscFunctionReturn(PETSC_SUCCESS);
1080: }

1082: /*@C
1083:   PetscDualSpaceApplyInterior - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by `PetscDualSpaceGetInteriorData()`

1085:   Input Parameters:
1086: + sp        - The `PetscDualSpace` object
1087: - pointEval - Evaluation at the points returned by `PetscDualSpaceGetInteriorData()`

1089:   Output Parameter:
1090: . spValue - The values of interior dual space functionals

1092:   Level: advanced

1094: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1095: @*/
1096: PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1097: {
1098:   PetscFunctionBegin;
1100:   PetscUseTypeMethod(sp, applyint, pointEval, spValue);
1101:   PetscFunctionReturn(PETSC_SUCCESS);
1102: }

1104: /*@C
1105:   PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.

1107:   Input Parameters:
1108: + sp    - The `PetscDualSpace` object
1109: . f     - The basis functional index
1110: . time  - The time
1111: . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
1112: . Nc    - The number of components for the function
1113: . func  - The input function
1114: - ctx   - A context for the function

1116:   Output Parameter:
1117: . value - The output value

1119:   Calling sequence:
1120: .vb
1121:    PetscErrorCode func(PetscInt dim, PetscReal time, const PetscReal x[],PetscInt numComponents, PetscScalar values[], PetscCtx ctx)
1122: .ve

1124:   Level: advanced

1126:   Note:
1127:   The idea is to evaluate the functional as an integral $ n(f) = \int dx n(x) . f(x) $ where both n and f have Nc components.

1129: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1130: @*/
1131: PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFEGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), PetscCtx ctx, PetscScalar *value)
1132: {
1133:   DM               dm;
1134:   PetscQuadrature  n;
1135:   const PetscReal *points, *weights;
1136:   PetscReal        x[3];
1137:   PetscScalar     *val;
1138:   PetscInt         dim, dE, qNc, c, Nq, q;
1139:   PetscBool        isAffine;

1141:   PetscFunctionBegin;
1143:   PetscAssertPointer(value, 8);
1144:   PetscCall(PetscDualSpaceGetDM(sp, &dm));
1145:   PetscCall(PetscDualSpaceGetFunctional(sp, f, &n));
1146:   PetscCall(PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights));
1147:   PetscCheck(dim == cgeom->dim, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %" PetscInt_FMT " != cell geometry dimension %" PetscInt_FMT, dim, cgeom->dim);
1148:   PetscCheck(qNc == Nc, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_SIZ, "The quadrature components %" PetscInt_FMT " != function components %" PetscInt_FMT, qNc, Nc);
1149:   PetscCall(DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val));
1150:   *value   = 0.0;
1151:   isAffine = cgeom->isAffine;
1152:   dE       = cgeom->dimEmbed;
1153:   for (q = 0; q < Nq; ++q) {
1154:     if (isAffine) {
1155:       CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q * dim], x);
1156:       PetscCall((*func)(dE, time, x, Nc, val, ctx));
1157:     } else {
1158:       PetscCall((*func)(dE, time, &cgeom->v[dE * q], Nc, val, ctx));
1159:     }
1160:     for (c = 0; c < Nc; ++c) *value += val[c] * weights[q * Nc + c];
1161:   }
1162:   PetscCall(DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val));
1163:   PetscFunctionReturn(PETSC_SUCCESS);
1164: }

1166: /*@C
1167:   PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by `PetscDualSpaceGetAllData()`

1169:   Input Parameters:
1170: + sp        - The `PetscDualSpace` object
1171: - pointEval - Evaluation at the points returned by `PetscDualSpaceGetAllData()`

1173:   Output Parameter:
1174: . spValue - The values of all dual space functionals

1176:   Level: advanced

1178: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1179: @*/
1180: PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1181: {
1182:   Vec pointValues, dofValues;
1183:   Mat allMat;

1185:   PetscFunctionBegin;
1187:   PetscAssertPointer(pointEval, 2);
1188:   PetscAssertPointer(spValue, 3);
1189:   PetscCall(PetscDualSpaceGetAllData(sp, NULL, &allMat));
1190:   if (!sp->allNodeValues) PetscCall(MatCreateVecs(allMat, &sp->allNodeValues, NULL));
1191:   pointValues = sp->allNodeValues;
1192:   if (!sp->allDofValues) PetscCall(MatCreateVecs(allMat, NULL, &sp->allDofValues));
1193:   dofValues = sp->allDofValues;
1194:   PetscCall(VecPlaceArray(pointValues, pointEval));
1195:   PetscCall(VecPlaceArray(dofValues, spValue));
1196:   PetscCall(MatMult(allMat, pointValues, dofValues));
1197:   PetscCall(VecResetArray(dofValues));
1198:   PetscCall(VecResetArray(pointValues));
1199:   PetscFunctionReturn(PETSC_SUCCESS);
1200: }

1202: /*@C
1203:   PetscDualSpaceApplyInteriorDefault - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by `PetscDualSpaceGetInteriorData()`

1205:   Input Parameters:
1206: + sp        - The `PetscDualSpace` object
1207: - pointEval - Evaluation at the points returned by `PetscDualSpaceGetInteriorData()`

1209:   Output Parameter:
1210: . spValue - The values of interior dual space functionals

1212:   Level: advanced

1214: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1215: @*/
1216: PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1217: {
1218:   Vec pointValues, dofValues;
1219:   Mat intMat;

1221:   PetscFunctionBegin;
1223:   PetscAssertPointer(pointEval, 2);
1224:   PetscAssertPointer(spValue, 3);
1225:   PetscCall(PetscDualSpaceGetInteriorData(sp, NULL, &intMat));
1226:   if (!sp->intNodeValues) PetscCall(MatCreateVecs(intMat, &sp->intNodeValues, NULL));
1227:   pointValues = sp->intNodeValues;
1228:   if (!sp->intDofValues) PetscCall(MatCreateVecs(intMat, NULL, &sp->intDofValues));
1229:   dofValues = sp->intDofValues;
1230:   PetscCall(VecPlaceArray(pointValues, pointEval));
1231:   PetscCall(VecPlaceArray(dofValues, spValue));
1232:   PetscCall(MatMult(intMat, pointValues, dofValues));
1233:   PetscCall(VecResetArray(dofValues));
1234:   PetscCall(VecResetArray(pointValues));
1235:   PetscFunctionReturn(PETSC_SUCCESS);
1236: }

1238: /*@C
1239:   PetscDualSpaceGetAllData - Get all quadrature nodes from this space, and the matrix that sends quadrature node values to degree-of-freedom values

1241:   Input Parameter:
1242: . sp - The dualspace

1244:   Output Parameters:
1245: + allNodes - A `PetscQuadrature` object containing all evaluation nodes, pass `NULL` if not needed
1246: - allMat   - A `Mat` for the node-to-dof transformation, pass `NULL` if not needed

1248:   Level: advanced

1250: .seealso: `PetscQuadrature`, `PetscDualSpace`, `PetscDualSpaceCreate()`, `Mat`
1251: @*/
1252: PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace sp, PeOp PetscQuadrature *allNodes, PeOp Mat *allMat)
1253: {
1254:   PetscFunctionBegin;
1256:   if (allNodes) PetscAssertPointer(allNodes, 2);
1257:   if (allMat) PetscAssertPointer(allMat, 3);
1258:   if ((!sp->allNodes || !sp->allMat) && sp->ops->createalldata) {
1259:     PetscQuadrature qpoints;
1260:     Mat             amat;

1262:     PetscUseTypeMethod(sp, createalldata, &qpoints, &amat);
1263:     PetscCall(PetscQuadratureDestroy(&sp->allNodes));
1264:     PetscCall(MatDestroy(&sp->allMat));
1265:     sp->allNodes = qpoints;
1266:     sp->allMat   = amat;
1267:   }
1268:   if (allNodes) *allNodes = sp->allNodes;
1269:   if (allMat) *allMat = sp->allMat;
1270:   PetscFunctionReturn(PETSC_SUCCESS);
1271: }

1273: /*@C
1274:   PetscDualSpaceCreateAllDataDefault - Create all evaluation nodes and the node-to-dof matrix by examining functionals

1276:   Input Parameter:
1277: . sp - The dualspace

1279:   Output Parameters:
1280: + allNodes - A `PetscQuadrature` object containing all evaluation nodes
1281: - allMat   - A `Mat` for the node-to-dof transformation

1283:   Level: advanced

1285: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`, `Mat`, `PetscQuadrature`
1286: @*/
1287: PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat)
1288: {
1289:   PetscInt        spdim;
1290:   PetscInt        numPoints, offset;
1291:   PetscReal      *points;
1292:   PetscInt        f, dim;
1293:   PetscInt        Nc, nrows, ncols;
1294:   PetscInt        maxNumPoints;
1295:   PetscQuadrature q;
1296:   Mat             A;

1298:   PetscFunctionBegin;
1299:   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
1300:   PetscCall(PetscDualSpaceGetDimension(sp, &spdim));
1301:   if (!spdim) {
1302:     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allNodes));
1303:     PetscCall(PetscQuadratureSetData(*allNodes, 0, 0, 0, NULL, NULL));
1304:   }
1305:   nrows = spdim;
1306:   PetscCall(PetscDualSpaceGetFunctional(sp, 0, &q));
1307:   PetscCall(PetscQuadratureGetData(q, &dim, NULL, &numPoints, NULL, NULL));
1308:   maxNumPoints = numPoints;
1309:   for (f = 1; f < spdim; f++) {
1310:     PetscInt Np;

1312:     PetscCall(PetscDualSpaceGetFunctional(sp, f, &q));
1313:     PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Np, NULL, NULL));
1314:     numPoints += Np;
1315:     maxNumPoints = PetscMax(maxNumPoints, Np);
1316:   }
1317:   ncols = numPoints * Nc;
1318:   PetscCall(PetscMalloc1(dim * numPoints, &points));
1319:   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nrows, ncols, maxNumPoints * Nc, NULL, &A));
1320:   for (f = 0, offset = 0; f < spdim; f++) {
1321:     const PetscReal *p, *w;
1322:     PetscInt         Np, i;
1323:     PetscInt         fnc;

1325:     PetscCall(PetscDualSpaceGetFunctional(sp, f, &q));
1326:     PetscCall(PetscQuadratureGetData(q, NULL, &fnc, &Np, &p, &w));
1327:     PetscCheck(fnc == Nc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "functional component mismatch");
1328:     for (i = 0; i < Np * dim; i++) points[offset * dim + i] = p[i];
1329:     for (i = 0; i < Np * Nc; i++) PetscCall(MatSetValue(A, f, offset * Nc, w[i], INSERT_VALUES));
1330:     offset += Np;
1331:   }
1332:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1333:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1334:   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allNodes));
1335:   PetscCall(PetscQuadratureSetData(*allNodes, dim, 0, numPoints, points, NULL));
1336:   *allMat = A;
1337:   PetscFunctionReturn(PETSC_SUCCESS);
1338: }

1340: /*@C
1341:   PetscDualSpaceGetInteriorData - Get all quadrature points necessary to compute the interior degrees of freedom from
1342:   this space, as well as the matrix that computes the degrees of freedom from the quadrature
1343:   values.

1345:   Input Parameter:
1346: . sp - The dualspace

1348:   Output Parameters:
1349: + intNodes - A `PetscQuadrature` object containing all evaluation points needed to evaluate interior degrees of freedom,
1350:              pass `NULL` if not needed
1351: - intMat   - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is
1352:              the size of the constrained layout (`PetscSectionGetConstrainStorageSize()`) of the dual space section,
1353:              npoints is the number of points in intNodes and nc is `PetscDualSpaceGetNumComponents()`.
1354:              Pass `NULL` if not needed

1356:   Level: advanced

1358:   Notes:
1359:   Degrees of freedom are interior degrees of freedom if they belong (by
1360:   `PetscDualSpaceGetSection()`) to interior points in the references, complementary boundary
1361:   degrees of freedom are marked as constrained in the section returned by
1362:   `PetscDualSpaceGetSection()`).

1364: .seealso: `PetscDualSpace`, `PetscQuadrature`, `Mat`, `PetscDualSpaceCreate()`, `PetscDualSpaceGetDimension()`, `PetscDualSpaceGetNumComponents()`, `PetscQuadratureGetData()`
1365: @*/
1366: PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace sp, PeOp PetscQuadrature *intNodes, PeOp Mat *intMat)
1367: {
1368:   PetscFunctionBegin;
1370:   if (intNodes) PetscAssertPointer(intNodes, 2);
1371:   if (intMat) PetscAssertPointer(intMat, 3);
1372:   if ((!sp->intNodes || !sp->intMat) && sp->ops->createintdata) {
1373:     PetscQuadrature qpoints;
1374:     Mat             imat;

1376:     PetscUseTypeMethod(sp, createintdata, &qpoints, &imat);
1377:     PetscCall(PetscQuadratureDestroy(&sp->intNodes));
1378:     PetscCall(MatDestroy(&sp->intMat));
1379:     sp->intNodes = qpoints;
1380:     sp->intMat   = imat;
1381:   }
1382:   if (intNodes) *intNodes = sp->intNodes;
1383:   if (intMat) *intMat = sp->intMat;
1384:   PetscFunctionReturn(PETSC_SUCCESS);
1385: }

1387: /*@C
1388:   PetscDualSpaceCreateInteriorDataDefault - Create quadrature points by examining interior functionals and create the matrix mapping quadrature point values to interior dual space values

1390:   Input Parameter:
1391: . sp - The dualspace

1393:   Output Parameters:
1394: + intNodes - A `PetscQuadrature` object containing all evaluation points needed to evaluate interior degrees of freedom
1395: - intMat   - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is
1396:               the size of the constrained layout (`PetscSectionGetConstrainStorageSize()`) of the dual space section,
1397:               npoints is the number of points in allNodes and nc is `PetscDualSpaceGetNumComponents()`.

1399:   Level: advanced

1401: .seealso: `PetscDualSpace`, `PetscQuadrature`, `Mat`, `PetscDualSpaceCreate()`, `PetscDualSpaceGetInteriorData()`
1402: @*/
1403: PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat)
1404: {
1405:   DM              dm;
1406:   PetscInt        spdim0;
1407:   PetscInt        Nc;
1408:   PetscInt        pStart, pEnd, p, f;
1409:   PetscSection    section;
1410:   PetscInt        numPoints, offset, matoffset;
1411:   PetscReal      *points;
1412:   PetscInt        dim;
1413:   PetscInt       *nnz;
1414:   PetscQuadrature q;
1415:   Mat             imat;

1417:   PetscFunctionBegin;
1419:   PetscCall(PetscDualSpaceGetSection(sp, &section));
1420:   PetscCall(PetscSectionGetConstrainedStorageSize(section, &spdim0));
1421:   if (!spdim0) {
1422:     *intNodes = NULL;
1423:     *intMat   = NULL;
1424:     PetscFunctionReturn(PETSC_SUCCESS);
1425:   }
1426:   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
1427:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1428:   PetscCall(PetscDualSpaceGetDM(sp, &dm));
1429:   PetscCall(DMGetDimension(dm, &dim));
1430:   PetscCall(PetscMalloc1(spdim0, &nnz));
1431:   for (p = pStart, f = 0, numPoints = 0; p < pEnd; p++) {
1432:     PetscInt dof, cdof, off, d;

1434:     PetscCall(PetscSectionGetDof(section, p, &dof));
1435:     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
1436:     if (!(dof - cdof)) continue;
1437:     PetscCall(PetscSectionGetOffset(section, p, &off));
1438:     for (d = 0; d < dof; d++, off++, f++) {
1439:       PetscInt Np;

1441:       PetscCall(PetscDualSpaceGetFunctional(sp, off, &q));
1442:       PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Np, NULL, NULL));
1443:       nnz[f] = Np * Nc;
1444:       numPoints += Np;
1445:     }
1446:   }
1447:   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, spdim0, numPoints * Nc, 0, nnz, &imat));
1448:   PetscCall(PetscFree(nnz));
1449:   PetscCall(PetscMalloc1(dim * numPoints, &points));
1450:   for (p = pStart, f = 0, offset = 0, matoffset = 0; p < pEnd; p++) {
1451:     PetscInt dof, cdof, off, d;

1453:     PetscCall(PetscSectionGetDof(section, p, &dof));
1454:     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
1455:     if (!(dof - cdof)) continue;
1456:     PetscCall(PetscSectionGetOffset(section, p, &off));
1457:     for (d = 0; d < dof; d++, off++, f++) {
1458:       const PetscReal *p;
1459:       const PetscReal *w;
1460:       PetscInt         Np, i;

1462:       PetscCall(PetscDualSpaceGetFunctional(sp, off, &q));
1463:       PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Np, &p, &w));
1464:       for (i = 0; i < Np * dim; i++) points[offset + i] = p[i];
1465:       for (i = 0; i < Np * Nc; i++) PetscCall(MatSetValue(imat, f, matoffset + i, w[i], INSERT_VALUES));
1466:       offset += Np * dim;
1467:       matoffset += Np * Nc;
1468:     }
1469:   }
1470:   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, intNodes));
1471:   PetscCall(PetscQuadratureSetData(*intNodes, dim, 0, numPoints, points, NULL));
1472:   PetscCall(MatAssemblyBegin(imat, MAT_FINAL_ASSEMBLY));
1473:   PetscCall(MatAssemblyEnd(imat, MAT_FINAL_ASSEMBLY));
1474:   *intMat = imat;
1475:   PetscFunctionReturn(PETSC_SUCCESS);
1476: }

1478: /*@C
1479:   PetscDualSpaceEqual - Determine if two dual spaces are equivalent

1481:   Input Parameters:
1482: + A - A `PetscDualSpace` object
1483: - B - Another `PetscDualSpace` object

1485:   Output Parameter:
1486: . equal - `PETSC_TRUE` if the dual spaces are equivalent

1488:   Level: advanced

1490: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1491: @*/
1492: PetscErrorCode PetscDualSpaceEqual(PetscDualSpace A, PetscDualSpace B, PetscBool *equal)
1493: {
1494:   PetscInt        sizeA, sizeB, dimA, dimB;
1495:   const PetscInt *dofA, *dofB;
1496:   PetscQuadrature quadA, quadB;
1497:   Mat             matA, matB;

1499:   PetscFunctionBegin;
1502:   PetscAssertPointer(equal, 3);
1503:   *equal = PETSC_FALSE;
1504:   PetscCall(PetscDualSpaceGetDimension(A, &sizeA));
1505:   PetscCall(PetscDualSpaceGetDimension(B, &sizeB));
1506:   if (sizeB != sizeA) PetscFunctionReturn(PETSC_SUCCESS);
1507:   PetscCall(DMGetDimension(A->dm, &dimA));
1508:   PetscCall(DMGetDimension(B->dm, &dimB));
1509:   if (dimA != dimB) PetscFunctionReturn(PETSC_SUCCESS);

1511:   PetscCall(PetscDualSpaceGetNumDof(A, &dofA));
1512:   PetscCall(PetscDualSpaceGetNumDof(B, &dofB));
1513:   for (PetscInt d = 0; d < dimA; d++) {
1514:     if (dofA[d] != dofB[d]) PetscFunctionReturn(PETSC_SUCCESS);
1515:   }

1517:   PetscCall(PetscDualSpaceGetInteriorData(A, &quadA, &matA));
1518:   PetscCall(PetscDualSpaceGetInteriorData(B, &quadB, &matB));
1519:   if (!quadA && !quadB) {
1520:     *equal = PETSC_TRUE;
1521:   } else if (quadA && quadB) {
1522:     PetscCall(PetscQuadratureEqual(quadA, quadB, equal));
1523:     if (*equal == PETSC_FALSE) PetscFunctionReturn(PETSC_SUCCESS);
1524:     if (!matA && !matB) PetscFunctionReturn(PETSC_SUCCESS);
1525:     if (matA && matB) PetscCall(MatEqual(matA, matB, equal));
1526:     else *equal = PETSC_FALSE;
1527:   }
1528:   PetscFunctionReturn(PETSC_SUCCESS);
1529: }

1531: /*@C
1532:   PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.

1534:   Input Parameters:
1535: + sp    - The `PetscDualSpace` object
1536: . f     - The basis functional index
1537: . time  - The time
1538: . cgeom - A context with geometric information for this cell, we currently just use the centroid
1539: . Nc    - The number of components for the function
1540: . func  - The input function
1541: - ctx   - A context for the function

1543:   Output Parameter:
1544: . value - The output value (scalar)

1546:   Calling sequence:
1547: .vb
1548:   PetscErrorCode func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt numComponents, PetscScalar values[], PetscCtx ctx)
1549: .ve

1551:   Level: advanced

1553:   Note:
1554:   The idea is to evaluate the functional as an integral $ n(f) = \int dx n(x) . f(x)$ where both n and f have Nc components.

1556: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1557: @*/
1558: PetscErrorCode PetscDualSpaceApplyFVM(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFVCellGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), PetscCtx ctx, PetscScalar *value)
1559: {
1560:   DM               dm;
1561:   PetscQuadrature  n;
1562:   const PetscReal *points, *weights;
1563:   PetscScalar     *val;
1564:   PetscInt         dimEmbed, qNc, c, Nq, q;

1566:   PetscFunctionBegin;
1568:   PetscAssertPointer(value, 8);
1569:   PetscCall(PetscDualSpaceGetDM(sp, &dm));
1570:   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
1571:   PetscCall(PetscDualSpaceGetFunctional(sp, f, &n));
1572:   PetscCall(PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights));
1573:   PetscCheck(qNc == Nc, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_SIZ, "The quadrature components %" PetscInt_FMT " != function components %" PetscInt_FMT, qNc, Nc);
1574:   PetscCall(DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val));
1575:   *value = 0.;
1576:   for (q = 0; q < Nq; ++q) {
1577:     PetscCall((*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx));
1578:     for (c = 0; c < Nc; ++c) *value += val[c] * weights[q * Nc + c];
1579:   }
1580:   PetscCall(DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val));
1581:   PetscFunctionReturn(PETSC_SUCCESS);
1582: }

1584: /*@C
1585:   PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
1586:   given height.  This assumes that the reference cell is symmetric over points of this height.

1588:   Not Collective

1590:   Input Parameters:
1591: + sp     - the `PetscDualSpace` object
1592: - height - the height of the mesh point for which the subspace is desired

1594:   Output Parameter:
1595: . subsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
1596:           point, which will be of lesser dimension if height > 0.

1598:   Level: advanced

1600:   Notes:
1601:   If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
1602:   pointwise values are not defined on the element boundaries), or if the implementation of `PetscDualSpace` does not
1603:   support extracting subspaces, then `NULL` is returned.

1605:   This does not increment the reference count on the returned dual space, and the user should not destroy it.

1607: .seealso: `PetscDualSpace`, `PetscSpaceGetHeightSubspace()`, `PetscDualSpaceGetPointSubspace()`
1608: @*/
1609: PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
1610: {
1611:   PetscInt depth = -1, cStart, cEnd;
1612:   DM       dm;

1614:   PetscFunctionBegin;
1616:   PetscAssertPointer(subsp, 3);
1617:   PetscCheck(sp->uniform, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A non-uniform dual space does not have a single dual space at each height");
1618:   *subsp = NULL;
1619:   dm     = sp->dm;
1620:   PetscCall(DMPlexGetDepth(dm, &depth));
1621:   PetscCheck(height >= 0 && height <= depth, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid height");
1622:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1623:   if (height == 0 && cEnd == cStart + 1) {
1624:     *subsp = sp;
1625:     PetscFunctionReturn(PETSC_SUCCESS);
1626:   }
1627:   if (!sp->heightSpaces) {
1628:     PetscInt h;
1629:     PetscCall(PetscCalloc1(depth + 1, &sp->heightSpaces));

1631:     for (h = 0; h <= depth; h++) {
1632:       if (h == 0 && cEnd == cStart + 1) continue;
1633:       if (sp->ops->createheightsubspace) PetscUseTypeMethod(sp, createheightsubspace, height, &sp->heightSpaces[h]);
1634:       else if (sp->pointSpaces) {
1635:         PetscInt hStart, hEnd;

1637:         PetscCall(DMPlexGetHeightStratum(dm, h, &hStart, &hEnd));
1638:         if (hEnd > hStart) {
1639:           const char *name;

1641:           PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[hStart]));
1642:           if (sp->pointSpaces[hStart]) {
1643:             PetscCall(PetscObjectGetName((PetscObject)sp, &name));
1644:             PetscCall(PetscObjectSetName((PetscObject)sp->pointSpaces[hStart], name));
1645:           }
1646:           sp->heightSpaces[h] = sp->pointSpaces[hStart];
1647:         }
1648:       }
1649:     }
1650:   }
1651:   *subsp = sp->heightSpaces[height];
1652:   PetscFunctionReturn(PETSC_SUCCESS);
1653: }

1655: /*@C
1656:   PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.

1658:   Not Collective

1660:   Input Parameters:
1661: + sp    - the `PetscDualSpace` object
1662: - point - the point (in the dual space's DM) for which the subspace is desired

1664:   Output Parameter:
1665: . bdsp - the subspace.

1667:   Level: advanced

1669:   Notes:
1670:   The functionals in the subspace are with respect to the intrinsic geometry of the point,
1671:   which will be of lesser dimension if height > 0.

1673:   If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
1674:   defined on the element boundaries), or if the implementation of `PetscDualSpace` does not support extracting
1675:   subspaces, then `NULL` is returned.

1677:   This does not increment the reference count on the returned dual space, and the user should not destroy it.

1679: .seealso: `PetscDualSpace`, `PetscDualSpaceGetHeightSubspace()`
1680: @*/
1681: PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
1682: {
1683:   PetscInt pStart = 0, pEnd = 0, cStart, cEnd;
1684:   DM       dm;

1686:   PetscFunctionBegin;
1688:   PetscAssertPointer(bdsp, 3);
1689:   *bdsp = NULL;
1690:   dm    = sp->dm;
1691:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1692:   PetscCheck(point >= pStart && point <= pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point");
1693:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1694:   if (point == cStart && cEnd == cStart + 1) { /* the dual space is only equivalent to the dual space on a cell if the reference mesh has just one cell */
1695:     *bdsp = sp;
1696:     PetscFunctionReturn(PETSC_SUCCESS);
1697:   }
1698:   if (!sp->pointSpaces) {
1699:     PetscInt p;
1700:     PetscCall(PetscCalloc1(pEnd - pStart, &sp->pointSpaces));

1702:     for (p = 0; p < pEnd - pStart; p++) {
1703:       if (p + pStart == cStart && cEnd == cStart + 1) continue;
1704:       if (sp->ops->createpointsubspace) PetscUseTypeMethod(sp, createpointsubspace, p + pStart, &sp->pointSpaces[p]);
1705:       else if (sp->heightSpaces || sp->ops->createheightsubspace) {
1706:         PetscInt dim, depth, height;
1707:         DMLabel  label;

1709:         PetscCall(DMPlexGetDepth(dm, &dim));
1710:         PetscCall(DMPlexGetDepthLabel(dm, &label));
1711:         PetscCall(DMLabelGetValue(label, p + pStart, &depth));
1712:         height = dim - depth;
1713:         PetscCall(PetscDualSpaceGetHeightSubspace(sp, height, &sp->pointSpaces[p]));
1714:         PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[p]));
1715:       }
1716:     }
1717:   }
1718:   *bdsp = sp->pointSpaces[point - pStart];
1719:   PetscFunctionReturn(PETSC_SUCCESS);
1720: }

1722: /*@C
1723:   PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis

1725:   Not Collective

1727:   Input Parameter:
1728: . sp - the `PetscDualSpace` object

1730:   Output Parameters:
1731: + perms - Permutations of the interior degrees of freedom, parameterized by the point orientation
1732: - flips - Sign reversal of the interior degrees of freedom, parameterized by the point orientation

1734:   Level: developer

1736:   Note:
1737:   The permutation and flip arrays are organized in the following way
1738: .vb
1739:   perms[p][ornt][dof # on point] = new local dof #
1740:   flips[p][ornt][dof # on point] = reversal or not
1741: .ve

1743: .seealso: `PetscDualSpace`
1744: @*/
1745: PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
1746: {
1747:   PetscFunctionBegin;
1749:   if (perms) {
1750:     PetscAssertPointer(perms, 2);
1751:     *perms = NULL;
1752:   }
1753:   if (flips) {
1754:     PetscAssertPointer(flips, 3);
1755:     *flips = NULL;
1756:   }
1757:   PetscTryTypeMethod(sp, getsymmetries, perms, flips);
1758:   PetscFunctionReturn(PETSC_SUCCESS);
1759: }

1761: /*@C
1762:   PetscDualSpaceGetFormDegree - Get the form degree k for the k-form the describes the pushforwards/pullbacks of this
1763:   dual space's functionals.

1765:   Input Parameter:
1766: . dsp - The `PetscDualSpace`

1768:   Output Parameter:
1769: . k - The *signed* degree k of the k.  If k >= 0, this means that the degrees of freedom are k-forms, and are stored
1770:         in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms.  So for example,
1771:         the 1-form basis in 3-D is (dx, dy, dz), and the 2-form basis in 3-D is (dx wedge dy, dx wedge dz, dy wedge dz).
1772:         If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the
1773:         Hodge star map.  So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms
1774:         but are stored as 1-forms.

1776:   Level: developer

1778: .seealso: `PetscDualSpace`, `PetscDTAltV`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceTransformType`
1779: @*/
1780: PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace dsp, PetscInt *k)
1781: {
1782:   PetscFunctionBeginHot;
1784:   PetscAssertPointer(k, 2);
1785:   *k = dsp->k;
1786:   PetscFunctionReturn(PETSC_SUCCESS);
1787: }

1789: /*@C
1790:   PetscDualSpaceSetFormDegree - Set the form degree k for the k-form the describes the pushforwards/pullbacks of this
1791:   dual space's functionals.

1793:   Input Parameters:
1794: + dsp - The `PetscDualSpace`
1795: - k   - The *signed* degree k of the k.  If k >= 0, this means that the degrees of freedom are k-forms, and are stored
1796:         in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms.  So for example,
1797:         the 1-form basis in 3-D is (dx, dy, dz), and the 2-form basis in 3-D is (dx wedge dy, dx wedge dz, dy wedge dz).
1798:         If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the
1799:         Hodge star map.  So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms
1800:         but are stored as 1-forms.

1802:   Level: developer

1804: .seealso: `PetscDualSpace`, `PetscDTAltV`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceTransformType`
1805: @*/
1806: PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace dsp, PetscInt k)
1807: {
1808:   PetscInt dim;

1810:   PetscFunctionBeginHot;
1812:   PetscCheck(!dsp->setupcalled, PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up");
1813:   dim = dsp->dm->dim;
1814:   PetscCheck((k >= -dim && k <= dim) || k == PETSC_FORM_DEGREE_UNDEFINED, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported %" PetscInt_FMT "-form on %" PetscInt_FMT "-dimensional reference cell", PetscAbsInt(k), dim);
1815:   dsp->k = k;
1816:   PetscFunctionReturn(PETSC_SUCCESS);
1817: }

1819: /*@C
1820:   PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space

1822:   Input Parameter:
1823: . dsp - The `PetscDualSpace`

1825:   Output Parameter:
1826: . k - The simplex dimension

1828:   Level: developer

1830:   Note:
1831:   Currently supported values are
1832: .vb
1833:   0: These are H_1 methods that only transform coordinates
1834:   1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM)
1835:   2: These are the same as 1
1836:   3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM)
1837: .ve

1839: .seealso: `PetscDualSpace`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceTransformType`
1840: @*/
1841: PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k)
1842: {
1843:   PetscInt dim;

1845:   PetscFunctionBeginHot;
1847:   PetscAssertPointer(k, 2);
1848:   dim = dsp->dm->dim;
1849:   if (!dsp->k) *k = IDENTITY_TRANSFORM;
1850:   else if (dsp->k == 1) *k = COVARIANT_PIOLA_TRANSFORM;
1851:   else if (dsp->k == -(dim - 1)) *k = CONTRAVARIANT_PIOLA_TRANSFORM;
1852:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported transformation");
1853:   PetscFunctionReturn(PETSC_SUCCESS);
1854: }

1856: /*@C
1857:   PetscDualSpaceTransform - Transform the function values

1859:   Input Parameters:
1860: + dsp       - The `PetscDualSpace`
1861: . trans     - The type of transform
1862: . isInverse - Flag to invert the transform
1863: . fegeom    - The cell geometry
1864: . Nv        - The number of function samples
1865: . Nc        - The number of function components
1866: - vals      - The function values

1868:   Output Parameter:
1869: . vals - The transformed function values

1871:   Level: intermediate

1873:   Note:
1874:   This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.

1876: .seealso: `PetscDualSpace`, `PetscDualSpaceTransformGradient()`, `PetscDualSpaceTransformHessian()`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransformType`
1877: @*/
1878: PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1879: {
1880:   PetscReal Jstar[9] = {0};
1881:   PetscInt  dim, v, c, Nk;

1883:   PetscFunctionBeginHot;
1885:   PetscAssertPointer(fegeom, 4);
1886:   PetscAssertPointer(vals, 7);
1887:   /* TODO: not handling dimEmbed != dim right now */
1888:   dim = dsp->dm->dim;
1889:   /* No change needed for 0-forms */
1890:   if (!dsp->k) PetscFunctionReturn(PETSC_SUCCESS);
1891:   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(dsp->k), &Nk));
1892:   /* TODO: use fegeom->isAffine */
1893:   PetscCall(PetscDTAltVPullbackMatrix(dim, dim, isInverse ? fegeom->J : fegeom->invJ, dsp->k, Jstar));
1894:   for (v = 0; v < Nv; ++v) {
1895:     switch (Nk) {
1896:     case 1:
1897:       for (c = 0; c < Nc; c++) vals[v * Nc + c] *= Jstar[0];
1898:       break;
1899:     case 2:
1900:       for (c = 0; c < Nc; c += 2) DMPlex_Mult2DReal_Internal(Jstar, 1, &vals[v * Nc + c], &vals[v * Nc + c]);
1901:       break;
1902:     case 3:
1903:       for (c = 0; c < Nc; c += 3) DMPlex_Mult3DReal_Internal(Jstar, 1, &vals[v * Nc + c], &vals[v * Nc + c]);
1904:       break;
1905:     default:
1906:       SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported form size %" PetscInt_FMT " for transformation", Nk);
1907:     }
1908:   }
1909:   PetscFunctionReturn(PETSC_SUCCESS);
1910: }

1912: /*@C
1913:   PetscDualSpaceTransformGradient - Transform the function gradient values

1915:   Input Parameters:
1916: + dsp       - The `PetscDualSpace`
1917: . trans     - The type of transform
1918: . isInverse - Flag to invert the transform
1919: . fegeom    - The cell geometry
1920: . Nv        - The number of function gradient samples
1921: . Nc        - The number of function components
1922: - vals      - The function gradient values

1924:   Output Parameter:
1925: . vals - The transformed function gradient values

1927:   Level: intermediate

1929:   Note:
1930:   This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.

1932: .seealso: `PetscDualSpace`, `PetscDualSpaceTransform()`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransformType`
1933: @*/
1934: PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1935: {
1936:   const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed;
1937:   PetscInt       v, c, d;

1939:   PetscFunctionBeginHot;
1941:   PetscAssertPointer(fegeom, 4);
1942:   PetscAssertPointer(vals, 7);
1943:   PetscAssert(dE > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %" PetscInt_FMT, dE);
1944:   /* Transform gradient */
1945:   if (dim == dE) {
1946:     for (v = 0; v < Nv; ++v) {
1947:       for (c = 0; c < Nc; ++c) {
1948:         switch (dim) {
1949:         case 1:
1950:           vals[(v * Nc + c) * dim] *= fegeom->invJ[0];
1951:           break;
1952:         case 2:
1953:           DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[(v * Nc + c) * dim], &vals[(v * Nc + c) * dim]);
1954:           break;
1955:         case 3:
1956:           DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[(v * Nc + c) * dim], &vals[(v * Nc + c) * dim]);
1957:           break;
1958:         default:
1959:           SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
1960:         }
1961:       }
1962:     }
1963:   } else {
1964:     for (v = 0; v < Nv; ++v) {
1965:       for (c = 0; c < Nc; ++c) DMPlex_MultTransposeReal_Internal(fegeom->invJ, dim, dE, 1, &vals[(v * Nc + c) * dE], &vals[(v * Nc + c) * dE]);
1966:     }
1967:   }
1968:   /* Assume its a vector, otherwise assume its a bunch of scalars */
1969:   if (Nc == 1 || Nc != dim) PetscFunctionReturn(PETSC_SUCCESS);
1970:   switch (trans) {
1971:   case IDENTITY_TRANSFORM:
1972:     break;
1973:   case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
1974:     if (isInverse) {
1975:       for (v = 0; v < Nv; ++v) {
1976:         for (d = 0; d < dim; ++d) {
1977:           switch (dim) {
1978:           case 2:
1979:             DMPlex_MultTranspose2DReal_Internal(fegeom->J, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
1980:             break;
1981:           case 3:
1982:             DMPlex_MultTranspose3DReal_Internal(fegeom->J, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
1983:             break;
1984:           default:
1985:             SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
1986:           }
1987:         }
1988:       }
1989:     } else {
1990:       for (v = 0; v < Nv; ++v) {
1991:         for (d = 0; d < dim; ++d) {
1992:           switch (dim) {
1993:           case 2:
1994:             DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
1995:             break;
1996:           case 3:
1997:             DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
1998:             break;
1999:           default:
2000:             SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2001:           }
2002:         }
2003:       }
2004:     }
2005:     break;
2006:   case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
2007:     if (isInverse) {
2008:       for (v = 0; v < Nv; ++v) {
2009:         for (d = 0; d < dim; ++d) {
2010:           switch (dim) {
2011:           case 2:
2012:             DMPlex_Mult2DReal_Internal(fegeom->invJ, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2013:             break;
2014:           case 3:
2015:             DMPlex_Mult3DReal_Internal(fegeom->invJ, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2016:             break;
2017:           default:
2018:             SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2019:           }
2020:           for (c = 0; c < Nc; ++c) vals[(v * Nc + c) * dim + d] *= fegeom->detJ[0];
2021:         }
2022:       }
2023:     } else {
2024:       for (v = 0; v < Nv; ++v) {
2025:         for (d = 0; d < dim; ++d) {
2026:           switch (dim) {
2027:           case 2:
2028:             DMPlex_Mult2DReal_Internal(fegeom->J, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2029:             break;
2030:           case 3:
2031:             DMPlex_Mult3DReal_Internal(fegeom->J, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2032:             break;
2033:           default:
2034:             SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2035:           }
2036:           for (c = 0; c < Nc; ++c) vals[(v * Nc + c) * dim + d] /= fegeom->detJ[0];
2037:         }
2038:       }
2039:     }
2040:     break;
2041:   }
2042:   PetscFunctionReturn(PETSC_SUCCESS);
2043: }

2045: /*@C
2046:   PetscDualSpaceTransformHessian - Transform the function Hessian values

2048:   Input Parameters:
2049: + dsp       - The `PetscDualSpace`
2050: . trans     - The type of transform
2051: . isInverse - Flag to invert the transform
2052: . fegeom    - The cell geometry
2053: . Nv        - The number of function Hessian samples
2054: . Nc        - The number of function components
2055: - vals      - The function gradient values

2057:   Output Parameter:
2058: . vals - The transformed function Hessian values

2060:   Level: intermediate

2062:   Note:
2063:   This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.

2065: .seealso: `PetscDualSpace`, `PetscDualSpaceTransform()`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransformType`
2066: @*/
2067: PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
2068: {
2069:   const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed;
2070:   PetscInt       v, c;

2072:   PetscFunctionBeginHot;
2074:   PetscAssertPointer(fegeom, 4);
2075:   PetscAssertPointer(vals, 7);
2076:   PetscAssert(dE > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %" PetscInt_FMT, dE);
2077:   /* Transform Hessian: J^{-T}_{ik} J^{-T}_{jl} H(f)_{kl} = J^{-T}_{ik} H(f)_{kl} J^{-1}_{lj} */
2078:   if (dim == dE) {
2079:     for (v = 0; v < Nv; ++v) {
2080:       for (c = 0; c < Nc; ++c) {
2081:         switch (dim) {
2082:         case 1:
2083:           vals[(v * Nc + c) * dim * dim] *= PetscSqr(fegeom->invJ[0]);
2084:           break;
2085:         case 2:
2086:           DMPlex_PTAP2DReal_Internal(fegeom->invJ, &vals[(v * Nc + c) * dim * dim], &vals[(v * Nc + c) * dim * dim]);
2087:           break;
2088:         case 3:
2089:           DMPlex_PTAP3DReal_Internal(fegeom->invJ, &vals[(v * Nc + c) * dim * dim], &vals[(v * Nc + c) * dim * dim]);
2090:           break;
2091:         default:
2092:           SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2093:         }
2094:       }
2095:     }
2096:   } else {
2097:     for (v = 0; v < Nv; ++v) {
2098:       for (c = 0; c < Nc; ++c) DMPlex_PTAPReal_Internal(fegeom->invJ, dim, dE, &vals[(v * Nc + c) * dE * dE], &vals[(v * Nc + c) * dE * dE]);
2099:     }
2100:   }
2101:   /* Assume its a vector, otherwise assume its a bunch of scalars */
2102:   if (Nc == 1 || Nc != dim) PetscFunctionReturn(PETSC_SUCCESS);
2103:   switch (trans) {
2104:   case IDENTITY_TRANSFORM:
2105:     break;
2106:   case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
2107:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Piola mapping for Hessians not yet supported");
2108:   case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
2109:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Piola mapping for Hessians not yet supported");
2110:   }
2111:   PetscFunctionReturn(PETSC_SUCCESS);
2112: }

2114: /*@C
2115:   PetscDualSpacePullback - Transform the given functional so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method.

2117:   Input Parameters:
2118: + dsp       - The `PetscDualSpace`
2119: . fegeom    - The geometry for this cell
2120: . Nq        - The number of function samples
2121: . Nc        - The number of function components
2122: - pointEval - The function values

2124:   Output Parameter:
2125: . pointEval - The transformed function values

2127:   Level: advanced

2129:   Notes:
2130:   Functions transform in a complementary way (pushforward) to functionals, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.

2132:   This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.

2134: .seealso: `PetscDualSpace`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2135: @*/
2136: PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2137: {
2138:   PetscDualSpaceTransformType trans;
2139:   PetscInt                    k;

2141:   PetscFunctionBeginHot;
2143:   PetscAssertPointer(fegeom, 2);
2144:   PetscAssertPointer(pointEval, 5);
2145:   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2146:      This determines their transformation properties. */
2147:   PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2148:   switch (k) {
2149:   case 0: /* H^1 point evaluations */
2150:     trans = IDENTITY_TRANSFORM;
2151:     break;
2152:   case 1: /* Hcurl preserves tangential edge traces  */
2153:     trans = COVARIANT_PIOLA_TRANSFORM;
2154:     break;
2155:   case 2:
2156:   case 3: /* Hdiv preserve normal traces */
2157:     trans = CONTRAVARIANT_PIOLA_TRANSFORM;
2158:     break;
2159:   default:
2160:     SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2161:   }
2162:   PetscCall(PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, pointEval));
2163:   PetscFunctionReturn(PETSC_SUCCESS);
2164: }

2166: /*@C
2167:   PetscDualSpacePushforward - Transform the given function so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method.

2169:   Input Parameters:
2170: + dsp       - The `PetscDualSpace`
2171: . fegeom    - The geometry for this cell
2172: . Nq        - The number of function samples
2173: . Nc        - The number of function components
2174: - pointEval - The function values

2176:   Output Parameter:
2177: . pointEval - The transformed function values

2179:   Level: advanced

2181:   Notes:
2182:   Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.

2184:   This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.

2186: .seealso: `PetscDualSpace`, `PetscDualSpacePullback()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2187: @*/
2188: PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2189: {
2190:   PetscDualSpaceTransformType trans;
2191:   PetscInt                    k;

2193:   PetscFunctionBeginHot;
2195:   PetscAssertPointer(fegeom, 2);
2196:   PetscAssertPointer(pointEval, 5);
2197:   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2198:      This determines their transformation properties. */
2199:   PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2200:   switch (k) {
2201:   case 0: /* H^1 point evaluations */
2202:     trans = IDENTITY_TRANSFORM;
2203:     break;
2204:   case 1: /* Hcurl preserves tangential edge traces  */
2205:     trans = COVARIANT_PIOLA_TRANSFORM;
2206:     break;
2207:   case 2:
2208:   case 3: /* Hdiv preserve normal traces */
2209:     trans = CONTRAVARIANT_PIOLA_TRANSFORM;
2210:     break;
2211:   default:
2212:     SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2213:   }
2214:   PetscCall(PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval));
2215:   PetscFunctionReturn(PETSC_SUCCESS);
2216: }

2218: /*@C
2219:   PetscDualSpacePushforwardGradient - Transform the given function gradient so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method.

2221:   Input Parameters:
2222: + dsp       - The `PetscDualSpace`
2223: . fegeom    - The geometry for this cell
2224: . Nq        - The number of function gradient samples
2225: . Nc        - The number of function components
2226: - pointEval - The function gradient values

2228:   Output Parameter:
2229: . pointEval - The transformed function gradient values

2231:   Level: advanced

2233:   Notes:
2234:   Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.

2236:   This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.

2238: .seealso: `PetscDualSpace`, `PetscDualSpacePushforward()`, `PetscDualSpacePullback()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2239: @*/
2240: PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2241: {
2242:   PetscDualSpaceTransformType trans;
2243:   PetscInt                    k;

2245:   PetscFunctionBeginHot;
2247:   PetscAssertPointer(fegeom, 2);
2248:   PetscAssertPointer(pointEval, 5);
2249:   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2250:      This determines their transformation properties. */
2251:   PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2252:   switch (k) {
2253:   case 0: /* H^1 point evaluations */
2254:     trans = IDENTITY_TRANSFORM;
2255:     break;
2256:   case 1: /* Hcurl preserves tangential edge traces  */
2257:     trans = COVARIANT_PIOLA_TRANSFORM;
2258:     break;
2259:   case 2:
2260:   case 3: /* Hdiv preserve normal traces */
2261:     trans = CONTRAVARIANT_PIOLA_TRANSFORM;
2262:     break;
2263:   default:
2264:     SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2265:   }
2266:   PetscCall(PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval));
2267:   PetscFunctionReturn(PETSC_SUCCESS);
2268: }

2270: /*@C
2271:   PetscDualSpacePushforwardHessian - Transform the given function Hessian so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method.

2273:   Input Parameters:
2274: + dsp       - The `PetscDualSpace`
2275: . fegeom    - The geometry for this cell
2276: . Nq        - The number of function Hessian samples
2277: . Nc        - The number of function components
2278: - pointEval - The function gradient values

2280:   Output Parameter:
2281: . pointEval - The transformed function Hessian values

2283:   Level: advanced

2285:   Notes:
2286:   Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.

2288:   This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.

2290: .seealso: `PetscDualSpace`, `PetscDualSpacePushforward()`, `PetscDualSpacePullback()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2291: @*/
2292: PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2293: {
2294:   PetscDualSpaceTransformType trans;
2295:   PetscInt                    k;

2297:   PetscFunctionBeginHot;
2299:   PetscAssertPointer(fegeom, 2);
2300:   PetscAssertPointer(pointEval, 5);
2301:   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2302:      This determines their transformation properties. */
2303:   PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2304:   switch (k) {
2305:   case 0: /* H^1 point evaluations */
2306:     trans = IDENTITY_TRANSFORM;
2307:     break;
2308:   case 1: /* Hcurl preserves tangential edge traces  */
2309:     trans = COVARIANT_PIOLA_TRANSFORM;
2310:     break;
2311:   case 2:
2312:   case 3: /* Hdiv preserve normal traces */
2313:     trans = CONTRAVARIANT_PIOLA_TRANSFORM;
2314:     break;
2315:   default:
2316:     SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2317:   }
2318:   PetscCall(PetscDualSpaceTransformHessian(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval));
2319:   PetscFunctionReturn(PETSC_SUCCESS);
2320: }