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
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; use -help for a list of available types
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, PetscObject obj, const char name[])
224: {
225: PetscFunctionBegin;
227: PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
228: PetscFunctionReturn(PETSC_SUCCESS);
229: }
231: /*@
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 iascii;
248: PetscFunctionBegin;
251: if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sp), &v));
252: PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii));
253: if (iascii) PetscCall(PetscDualSpaceView_ASCII(sp, v));
254: PetscFunctionReturn(PETSC_SUCCESS);
255: }
257: /*@
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 - Flag for continuous element
271: . -petscdualspace_lagrange_tensor - Flag for tensor dual space
272: . -petscdualspace_lagrange_trimmed - Flag for trimmed dual space
273: . -petscdualspace_lagrange_node_type <nodetype> - Lagrange node location type
274: . -petscdualspace_lagrange_node_endpoints - Flag for nodes that include endpoints
275: . -petscdualspace_lagrange_node_exponent - Gauss-Jacobi weight function exponent
276: . -petscdualspace_lagrange_use_moments - 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) {
293: defaultType = PETSCDUALSPACELAGRANGE;
294: } else {
295: defaultType = ((PetscObject)sp)->type_name;
296: }
297: if (!PetscSpaceRegisterAllCalled) PetscCall(PetscSpaceRegisterAll());
299: PetscObjectOptionsBegin((PetscObject)sp);
300: PetscCall(PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg));
301: if (flg) {
302: PetscCall(PetscDualSpaceSetType(sp, name));
303: } else if (!((PetscObject)sp)->type_name) {
304: PetscCall(PetscDualSpaceSetType(sp, defaultType));
305: }
306: PetscCall(PetscOptionsBoundedInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL, 0));
307: PetscCall(PetscOptionsInt("-petscdualspace_form_degree", "The form degree of the dofs", "PetscDualSpaceSetFormDegree", sp->k, &sp->k, NULL));
308: PetscCall(PetscOptionsBoundedInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL, 1));
309: PetscTryTypeMethod(sp, setfromoptions, PetscOptionsObject);
310: PetscCall(PetscOptionsEnum("-petscdualspace_refcell", "Reference cell shape", "PetscDualSpaceSetReferenceCell", DMPolytopeTypes, (PetscEnum)refCell, (PetscEnum *)&refCell, &flg));
311: if (flg) {
312: DM K;
314: PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, refCell, &K));
315: PetscCall(PetscDualSpaceSetDM(sp, K));
316: PetscCall(DMDestroy(&K));
317: }
319: /* process any options handlers added with PetscObjectAddOptionsHandler() */
320: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)sp, PetscOptionsObject));
321: PetscOptionsEnd();
322: sp->setfromoptionscalled = PETSC_TRUE;
323: PetscFunctionReturn(PETSC_SUCCESS);
324: }
326: /*@
327: PetscDualSpaceSetUp - Construct a basis for a `PetscDualSpace`
329: Collective
331: Input Parameter:
332: . sp - the `PetscDualSpace` object to setup
334: Level: intermediate
336: .seealso: `PetscDualSpaceView()`, `PetscDualSpaceDestroy()`, `PetscDualSpace`
337: @*/
338: PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
339: {
340: PetscFunctionBegin;
342: if (sp->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
343: PetscCall(PetscLogEventBegin(PETSCDUALSPACE_SetUp, sp, 0, 0, 0));
344: sp->setupcalled = PETSC_TRUE;
345: PetscTryTypeMethod(sp, setup);
346: PetscCall(PetscLogEventEnd(PETSCDUALSPACE_SetUp, sp, 0, 0, 0));
347: if (sp->setfromoptionscalled) PetscCall(PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view"));
348: PetscFunctionReturn(PETSC_SUCCESS);
349: }
351: static PetscErrorCode PetscDualSpaceClearDMData_Internal(PetscDualSpace sp, DM dm)
352: {
353: PetscInt pStart = -1, pEnd = -1, depth = -1;
355: PetscFunctionBegin;
356: if (!dm) PetscFunctionReturn(PETSC_SUCCESS);
357: PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
358: PetscCall(DMPlexGetDepth(dm, &depth));
360: if (sp->pointSpaces) {
361: PetscInt i;
363: for (i = 0; i < pEnd - pStart; i++) PetscCall(PetscDualSpaceDestroy(&sp->pointSpaces[i]));
364: }
365: PetscCall(PetscFree(sp->pointSpaces));
367: if (sp->heightSpaces) {
368: PetscInt i;
370: for (i = 0; i <= depth; i++) PetscCall(PetscDualSpaceDestroy(&sp->heightSpaces[i]));
371: }
372: PetscCall(PetscFree(sp->heightSpaces));
374: PetscCall(PetscSectionDestroy(&sp->pointSection));
375: PetscCall(PetscSectionDestroy(&sp->intPointSection));
376: PetscCall(PetscQuadratureDestroy(&sp->intNodes));
377: PetscCall(VecDestroy(&sp->intDofValues));
378: PetscCall(VecDestroy(&sp->intNodeValues));
379: PetscCall(MatDestroy(&sp->intMat));
380: PetscCall(PetscQuadratureDestroy(&sp->allNodes));
381: PetscCall(VecDestroy(&sp->allDofValues));
382: PetscCall(VecDestroy(&sp->allNodeValues));
383: PetscCall(MatDestroy(&sp->allMat));
384: PetscCall(PetscFree(sp->numDof));
385: PetscFunctionReturn(PETSC_SUCCESS);
386: }
388: /*@
389: PetscDualSpaceDestroy - Destroys a `PetscDualSpace` object
391: Collective
393: Input Parameter:
394: . sp - the `PetscDualSpace` object to destroy
396: Level: beginner
398: .seealso: `PetscDualSpace`, `PetscDualSpaceView()`, `PetscDualSpace()`, `PetscDualSpaceCreate()`
399: @*/
400: PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
401: {
402: PetscInt dim, f;
403: DM dm;
405: PetscFunctionBegin;
406: if (!*sp) PetscFunctionReturn(PETSC_SUCCESS);
409: if (--((PetscObject)*sp)->refct > 0) {
410: *sp = NULL;
411: PetscFunctionReturn(PETSC_SUCCESS);
412: }
413: ((PetscObject)*sp)->refct = 0;
415: PetscCall(PetscDualSpaceGetDimension(*sp, &dim));
416: dm = (*sp)->dm;
418: PetscTryTypeMethod(*sp, destroy);
419: PetscCall(PetscDualSpaceClearDMData_Internal(*sp, dm));
421: for (f = 0; f < dim; ++f) PetscCall(PetscQuadratureDestroy(&(*sp)->functional[f]));
422: PetscCall(PetscFree((*sp)->functional));
423: PetscCall(DMDestroy(&(*sp)->dm));
424: PetscCall(PetscHeaderDestroy(sp));
425: PetscFunctionReturn(PETSC_SUCCESS);
426: }
428: /*@
429: PetscDualSpaceCreate - Creates an empty `PetscDualSpace` object. The type can then be set with `PetscDualSpaceSetType()`.
431: Collective
433: Input Parameter:
434: . comm - The communicator for the `PetscDualSpace` object
436: Output Parameter:
437: . sp - The `PetscDualSpace` object
439: Level: beginner
441: .seealso: `PetscDualSpace`, `PetscDualSpaceSetType()`, `PETSCDUALSPACELAGRANGE`
442: @*/
443: PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
444: {
445: PetscDualSpace s;
447: PetscFunctionBegin;
448: PetscAssertPointer(sp, 2);
449: PetscCall(PetscCitationsRegister(FECitation, &FEcite));
450: *sp = NULL;
451: PetscCall(PetscFEInitializePackage());
453: PetscCall(PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView));
455: s->order = 0;
456: s->Nc = 1;
457: s->k = 0;
458: s->spdim = -1;
459: s->spintdim = -1;
460: s->uniform = PETSC_TRUE;
461: s->setupcalled = PETSC_FALSE;
463: *sp = s;
464: PetscFunctionReturn(PETSC_SUCCESS);
465: }
467: /*@
468: PetscDualSpaceDuplicate - Creates a duplicate `PetscDualSpace` object that is not setup.
470: Collective
472: Input Parameter:
473: . sp - The original `PetscDualSpace`
475: Output Parameter:
476: . spNew - The duplicate `PetscDualSpace`
478: Level: beginner
480: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`
481: @*/
482: PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
483: {
484: DM dm;
485: PetscDualSpaceType type;
486: const char *name;
488: PetscFunctionBegin;
490: PetscAssertPointer(spNew, 2);
491: PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)sp), spNew));
492: name = ((PetscObject)sp)->name;
493: if (name) { PetscCall(PetscObjectSetName((PetscObject)*spNew, name)); }
494: PetscCall(PetscDualSpaceGetType(sp, &type));
495: PetscCall(PetscDualSpaceSetType(*spNew, type));
496: PetscCall(PetscDualSpaceGetDM(sp, &dm));
497: PetscCall(PetscDualSpaceSetDM(*spNew, dm));
499: (*spNew)->order = sp->order;
500: (*spNew)->k = sp->k;
501: (*spNew)->Nc = sp->Nc;
502: (*spNew)->uniform = sp->uniform;
503: PetscTryTypeMethod(sp, duplicate, *spNew);
504: PetscFunctionReturn(PETSC_SUCCESS);
505: }
507: /*@
508: PetscDualSpaceGetDM - Get the `DM` representing the reference cell of a `PetscDualSpace`
510: Not Collective
512: Input Parameter:
513: . sp - The `PetscDualSpace`
515: Output Parameter:
516: . dm - The reference cell, that is a `DM` that consists of a single cell
518: Level: intermediate
520: .seealso: `PetscDualSpace`, `PetscDualSpaceSetDM()`, `PetscDualSpaceCreate()`
521: @*/
522: PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
523: {
524: PetscFunctionBegin;
526: PetscAssertPointer(dm, 2);
527: *dm = sp->dm;
528: PetscFunctionReturn(PETSC_SUCCESS);
529: }
531: /*@
532: PetscDualSpaceSetDM - Get the `DM` representing the reference cell
534: Not Collective
536: Input Parameters:
537: + sp - The `PetscDual`Space
538: - dm - The reference cell
540: Level: intermediate
542: .seealso: `PetscDualSpace`, `DM`, `PetscDualSpaceGetDM()`, `PetscDualSpaceCreate()`
543: @*/
544: PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
545: {
546: PetscFunctionBegin;
549: PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change DM after dualspace is set up");
550: PetscCall(PetscObjectReference((PetscObject)dm));
551: if (sp->dm && sp->dm != dm) PetscCall(PetscDualSpaceClearDMData_Internal(sp, sp->dm));
552: PetscCall(DMDestroy(&sp->dm));
553: sp->dm = dm;
554: PetscFunctionReturn(PETSC_SUCCESS);
555: }
557: /*@
558: PetscDualSpaceGetOrder - Get the order of the dual space
560: Not Collective
562: Input Parameter:
563: . sp - The `PetscDualSpace`
565: Output Parameter:
566: . order - The order
568: Level: intermediate
570: .seealso: `PetscDualSpace`, `PetscDualSpaceSetOrder()`, `PetscDualSpaceCreate()`
571: @*/
572: PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
573: {
574: PetscFunctionBegin;
576: PetscAssertPointer(order, 2);
577: *order = sp->order;
578: PetscFunctionReturn(PETSC_SUCCESS);
579: }
581: /*@
582: PetscDualSpaceSetOrder - Set the order of the dual space
584: Not Collective
586: Input Parameters:
587: + sp - The `PetscDualSpace`
588: - order - The order
590: Level: intermediate
592: .seealso: `PetscDualSpace`, `PetscDualSpaceGetOrder()`, `PetscDualSpaceCreate()`
593: @*/
594: PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
595: {
596: PetscFunctionBegin;
598: PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change order after dualspace is set up");
599: sp->order = order;
600: PetscFunctionReturn(PETSC_SUCCESS);
601: }
603: /*@
604: PetscDualSpaceGetNumComponents - Return the number of components for this space
606: Input Parameter:
607: . sp - The `PetscDualSpace`
609: Output Parameter:
610: . Nc - The number of components
612: Level: intermediate
614: Note:
615: A vector space, for example, will have d components, where d is the spatial dimension
617: .seealso: `PetscDualSpaceSetNumComponents()`, `PetscDualSpaceGetDimension()`, `PetscDualSpaceCreate()`, `PetscDualSpace`
618: @*/
619: PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
620: {
621: PetscFunctionBegin;
623: PetscAssertPointer(Nc, 2);
624: *Nc = sp->Nc;
625: PetscFunctionReturn(PETSC_SUCCESS);
626: }
628: /*@
629: PetscDualSpaceSetNumComponents - Set the number of components for this space
631: Input Parameters:
632: + sp - The `PetscDualSpace`
633: - Nc - The number of components
635: Level: intermediate
637: .seealso: `PetscDualSpaceGetNumComponents()`, `PetscDualSpaceCreate()`, `PetscDualSpace`
638: @*/
639: PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
640: {
641: PetscFunctionBegin;
643: PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up");
644: sp->Nc = Nc;
645: PetscFunctionReturn(PETSC_SUCCESS);
646: }
648: /*@
649: PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
651: Not Collective
653: Input Parameters:
654: + sp - The `PetscDualSpace`
655: - i - The basis number
657: Output Parameter:
658: . functional - The basis functional
660: Level: intermediate
662: .seealso: `PetscDualSpace`, `PetscQuadrature`, `PetscDualSpaceGetDimension()`, `PetscDualSpaceCreate()`
663: @*/
664: PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
665: {
666: PetscInt dim;
668: PetscFunctionBegin;
670: PetscAssertPointer(functional, 3);
671: PetscCall(PetscDualSpaceGetDimension(sp, &dim));
672: PetscCheck(!(i < 0) && !(i >= dim), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", i, dim);
673: *functional = sp->functional[i];
674: PetscFunctionReturn(PETSC_SUCCESS);
675: }
677: /*@
678: PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
680: Not Collective
682: Input Parameter:
683: . sp - The `PetscDualSpace`
685: Output Parameter:
686: . dim - The dimension
688: Level: intermediate
690: .seealso: `PetscDualSpace`, `PetscDualSpaceGetFunctional()`, `PetscDualSpaceCreate()`
691: @*/
692: PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
693: {
694: PetscFunctionBegin;
696: PetscAssertPointer(dim, 2);
697: if (sp->spdim < 0) {
698: PetscSection section;
700: PetscCall(PetscDualSpaceGetSection(sp, §ion));
701: if (section) {
702: PetscCall(PetscSectionGetStorageSize(section, &sp->spdim));
703: } else sp->spdim = 0;
704: }
705: *dim = sp->spdim;
706: PetscFunctionReturn(PETSC_SUCCESS);
707: }
709: /*@
710: PetscDualSpaceGetInteriorDimension - Get the interior dimension of the dual space, i.e. the number of basis functionals assigned to the interior of the reference domain
712: Not Collective
714: Input Parameter:
715: . sp - The `PetscDualSpace`
717: Output Parameter:
718: . intdim - The dimension
720: Level: intermediate
722: .seealso: `PetscDualSpace`, `PetscDualSpaceGetFunctional()`, `PetscDualSpaceCreate()`
723: @*/
724: PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace sp, PetscInt *intdim)
725: {
726: PetscFunctionBegin;
728: PetscAssertPointer(intdim, 2);
729: if (sp->spintdim < 0) {
730: PetscSection section;
732: PetscCall(PetscDualSpaceGetSection(sp, §ion));
733: if (section) {
734: PetscCall(PetscSectionGetConstrainedStorageSize(section, &sp->spintdim));
735: } else sp->spintdim = 0;
736: }
737: *intdim = sp->spintdim;
738: PetscFunctionReturn(PETSC_SUCCESS);
739: }
741: /*@
742: PetscDualSpaceGetUniform - Whether this dual space is uniform
744: Not Collective
746: Input Parameter:
747: . sp - A dual space
749: Output Parameter:
750: . uniform - `PETSC_TRUE` if (a) the dual space is the same for each point in a stratum of the reference `DMPLEX`, and
751: (b) every symmetry of each point in the reference `DMPLEX` is also a symmetry of the point's dual space.
753: Level: advanced
755: Note:
756: All of the usual spaces on simplex or tensor-product elements will be uniform, only reference cells
757: with non-uniform strata (like trianguar-prisms) or anisotropic hp dual spaces will not be uniform.
759: .seealso: `PetscDualSpace`, `PetscDualSpaceGetPointSubspace()`, `PetscDualSpaceGetSymmetries()`
760: @*/
761: PetscErrorCode PetscDualSpaceGetUniform(PetscDualSpace sp, PetscBool *uniform)
762: {
763: PetscFunctionBegin;
765: PetscAssertPointer(uniform, 2);
766: *uniform = sp->uniform;
767: PetscFunctionReturn(PETSC_SUCCESS);
768: }
770: /*@C
771: PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
773: Not Collective
775: Input Parameter:
776: . sp - The `PetscDualSpace`
778: Output Parameter:
779: . numDof - An array of length dim+1 which holds the number of dofs for each dimension
781: Level: intermediate
783: Note:
784: Do not free `numDof`
786: .seealso: `PetscDualSpace`, `PetscDualSpaceGetFunctional()`, `PetscDualSpaceCreate()`
787: @*/
788: PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt *numDof[])
789: {
790: PetscFunctionBegin;
792: PetscAssertPointer(numDof, 2);
793: 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");
794: if (!sp->numDof) {
795: DM dm;
796: PetscInt depth, d;
797: PetscSection section;
799: PetscCall(PetscDualSpaceGetDM(sp, &dm));
800: PetscCall(DMPlexGetDepth(dm, &depth));
801: PetscCall(PetscCalloc1(depth + 1, &sp->numDof));
802: PetscCall(PetscDualSpaceGetSection(sp, §ion));
803: for (d = 0; d <= depth; d++) {
804: PetscInt dStart, dEnd;
806: PetscCall(DMPlexGetDepthStratum(dm, d, &dStart, &dEnd));
807: if (dEnd <= dStart) continue;
808: PetscCall(PetscSectionGetDof(section, dStart, &sp->numDof[d]));
809: }
810: }
811: *numDof = sp->numDof;
812: PetscCheck(*numDof, PetscObjectComm((PetscObject)sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
813: PetscFunctionReturn(PETSC_SUCCESS);
814: }
816: /* create the section of the right size and set a permutation for topological ordering */
817: PetscErrorCode PetscDualSpaceSectionCreate_Internal(PetscDualSpace sp, PetscSection *topSection)
818: {
819: DM dm;
820: PetscInt pStart, pEnd, cStart, cEnd, c, depth, count, i;
821: PetscInt *seen, *perm;
822: PetscSection section;
824: PetscFunctionBegin;
825: dm = sp->dm;
826: PetscCall(PetscSectionCreate(PETSC_COMM_SELF, §ion));
827: PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
828: PetscCall(PetscSectionSetChart(section, pStart, pEnd));
829: PetscCall(PetscCalloc1(pEnd - pStart, &seen));
830: PetscCall(PetscMalloc1(pEnd - pStart, &perm));
831: PetscCall(DMPlexGetDepth(dm, &depth));
832: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
833: for (c = cStart, count = 0; c < cEnd; c++) {
834: PetscInt closureSize = -1, e;
835: PetscInt *closure = NULL;
837: perm[count++] = c;
838: seen[c - pStart] = 1;
839: PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
840: for (e = 0; e < closureSize; e++) {
841: PetscInt point = closure[2 * e];
843: if (seen[point - pStart]) continue;
844: perm[count++] = point;
845: seen[point - pStart] = 1;
846: }
847: PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
848: }
849: PetscCheck(count == pEnd - pStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad topological ordering");
850: for (i = 0; i < pEnd - pStart; i++)
851: if (perm[i] != i) break;
852: if (i < pEnd - pStart) {
853: IS permIS;
855: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, perm, PETSC_OWN_POINTER, &permIS));
856: PetscCall(ISSetPermutation(permIS));
857: PetscCall(PetscSectionSetPermutation(section, permIS));
858: PetscCall(ISDestroy(&permIS));
859: } else {
860: PetscCall(PetscFree(perm));
861: }
862: PetscCall(PetscFree(seen));
863: *topSection = section;
864: PetscFunctionReturn(PETSC_SUCCESS);
865: }
867: /* mark boundary points and set up */
868: PetscErrorCode PetscDualSpaceSectionSetUp_Internal(PetscDualSpace sp, PetscSection section)
869: {
870: DM dm;
871: DMLabel boundary;
872: PetscInt pStart, pEnd, p;
874: PetscFunctionBegin;
875: dm = sp->dm;
876: PetscCall(DMLabelCreate(PETSC_COMM_SELF, "boundary", &boundary));
877: PetscCall(PetscDualSpaceGetDM(sp, &dm));
878: PetscCall(DMPlexMarkBoundaryFaces(dm, 1, boundary));
879: PetscCall(DMPlexLabelComplete(dm, boundary));
880: PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
881: for (p = pStart; p < pEnd; p++) {
882: PetscInt bval;
884: PetscCall(DMLabelGetValue(boundary, p, &bval));
885: if (bval == 1) {
886: PetscInt dof;
888: PetscCall(PetscSectionGetDof(section, p, &dof));
889: PetscCall(PetscSectionSetConstraintDof(section, p, dof));
890: }
891: }
892: PetscCall(DMLabelDestroy(&boundary));
893: PetscCall(PetscSectionSetUp(section));
894: PetscFunctionReturn(PETSC_SUCCESS);
895: }
897: /*@
898: PetscDualSpaceGetSection - Create a `PetscSection` over the reference cell with the layout from this space
900: Collective
902: Input Parameter:
903: . sp - The `PetscDualSpace`
905: Output Parameter:
906: . section - The section
908: Level: advanced
910: .seealso: `PetscDualSpace`, `PetscSection`, `PetscDualSpaceCreate()`, `DMPLEX`
911: @*/
912: PetscErrorCode PetscDualSpaceGetSection(PetscDualSpace sp, PetscSection *section)
913: {
914: PetscInt pStart, pEnd, p;
916: PetscFunctionBegin;
917: if (!sp->dm) {
918: *section = NULL;
919: PetscFunctionReturn(PETSC_SUCCESS);
920: }
921: if (!sp->pointSection) {
922: /* mark the boundary */
923: PetscCall(PetscDualSpaceSectionCreate_Internal(sp, &sp->pointSection));
924: PetscCall(DMPlexGetChart(sp->dm, &pStart, &pEnd));
925: for (p = pStart; p < pEnd; p++) {
926: PetscDualSpace psp;
928: PetscCall(PetscDualSpaceGetPointSubspace(sp, p, &psp));
929: if (psp) {
930: PetscInt dof;
932: PetscCall(PetscDualSpaceGetInteriorDimension(psp, &dof));
933: PetscCall(PetscSectionSetDof(sp->pointSection, p, dof));
934: }
935: }
936: PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, sp->pointSection));
937: }
938: *section = sp->pointSection;
939: PetscFunctionReturn(PETSC_SUCCESS);
940: }
942: /*@
943: PetscDualSpaceGetInteriorSection - Create a `PetscSection` over the reference cell with the layout from this space
944: for interior degrees of freedom
946: Collective
948: Input Parameter:
949: . sp - The `PetscDualSpace`
951: Output Parameter:
952: . section - The interior section
954: Level: advanced
956: Note:
957: Most reference domains have one cell, in which case the only cell will have
958: all of the interior degrees of freedom in the interior section. But
959: for `PETSCDUALSPACEREFINED` there may be other mesh points in the interior,
960: and this section describes their layout.
962: .seealso: `PetscDualSpace`, `PetscSection`, `PetscDualSpaceCreate()`, `DMPLEX`
963: @*/
964: PetscErrorCode PetscDualSpaceGetInteriorSection(PetscDualSpace sp, PetscSection *section)
965: {
966: PetscInt pStart, pEnd, p;
968: PetscFunctionBegin;
969: if (!sp->dm) {
970: *section = NULL;
971: PetscFunctionReturn(PETSC_SUCCESS);
972: }
973: if (!sp->intPointSection) {
974: PetscSection full_section;
976: PetscCall(PetscDualSpaceGetSection(sp, &full_section));
977: PetscCall(PetscDualSpaceSectionCreate_Internal(sp, &sp->intPointSection));
978: PetscCall(PetscSectionGetChart(full_section, &pStart, &pEnd));
979: for (p = pStart; p < pEnd; p++) {
980: PetscInt dof, cdof;
982: PetscCall(PetscSectionGetDof(full_section, p, &dof));
983: PetscCall(PetscSectionGetConstraintDof(full_section, p, &cdof));
984: PetscCall(PetscSectionSetDof(sp->intPointSection, p, dof - cdof));
985: }
986: PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, sp->intPointSection));
987: }
988: *section = sp->intPointSection;
989: PetscFunctionReturn(PETSC_SUCCESS);
990: }
992: /* this assumes that all of the point dual spaces store their interior dofs first, which is true when the point DMs
993: * have one cell */
994: PetscErrorCode PetscDualSpacePushForwardSubspaces_Internal(PetscDualSpace sp, PetscInt sStart, PetscInt sEnd)
995: {
996: PetscReal *sv0, *v0, *J;
997: PetscSection section;
998: PetscInt dim, s, k;
999: DM dm;
1001: PetscFunctionBegin;
1002: PetscCall(PetscDualSpaceGetDM(sp, &dm));
1003: PetscCall(DMGetDimension(dm, &dim));
1004: PetscCall(PetscDualSpaceGetSection(sp, §ion));
1005: PetscCall(PetscMalloc3(dim, &v0, dim, &sv0, dim * dim, &J));
1006: PetscCall(PetscDualSpaceGetFormDegree(sp, &k));
1007: for (s = sStart; s < sEnd; s++) {
1008: PetscReal detJ, hdetJ;
1009: PetscDualSpace ssp;
1010: PetscInt dof, off, f, sdim;
1011: PetscInt i, j;
1012: DM sdm;
1014: PetscCall(PetscDualSpaceGetPointSubspace(sp, s, &ssp));
1015: if (!ssp) continue;
1016: PetscCall(PetscSectionGetDof(section, s, &dof));
1017: PetscCall(PetscSectionGetOffset(section, s, &off));
1018: /* get the first vertex of the reference cell */
1019: PetscCall(PetscDualSpaceGetDM(ssp, &sdm));
1020: PetscCall(DMGetDimension(sdm, &sdim));
1021: PetscCall(DMPlexComputeCellGeometryAffineFEM(sdm, 0, sv0, NULL, NULL, &hdetJ));
1022: PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, s, v0, J, NULL, &detJ));
1023: /* compactify Jacobian */
1024: for (i = 0; i < dim; i++)
1025: for (j = 0; j < sdim; j++) J[i * sdim + j] = J[i * dim + j];
1026: for (f = 0; f < dof; f++) {
1027: PetscQuadrature fn;
1029: PetscCall(PetscDualSpaceGetFunctional(ssp, f, &fn));
1030: PetscCall(PetscQuadraturePushForward(fn, dim, sv0, v0, J, k, &sp->functional[off + f]));
1031: }
1032: }
1033: PetscCall(PetscFree3(v0, sv0, J));
1034: PetscFunctionReturn(PETSC_SUCCESS);
1035: }
1037: /*@C
1038: PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
1040: Input Parameters:
1041: + sp - The `PetscDualSpace` object
1042: . f - The basis functional index
1043: . time - The time
1044: . 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)
1045: . numComp - The number of components for the function
1046: . func - The input function
1047: - ctx - A context for the function
1049: Output Parameter:
1050: . value - numComp output values
1052: Calling sequence:
1053: .vb
1054: PetscErrorCode func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt numComponents, PetscScalar values[], void *ctx)
1055: .ve
1057: Level: beginner
1059: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1060: @*/
1061: PetscErrorCode PetscDualSpaceApply(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFEGeom *cgeom, PetscInt numComp, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value)
1062: {
1063: PetscFunctionBegin;
1065: PetscAssertPointer(cgeom, 4);
1066: PetscAssertPointer(value, 8);
1067: PetscUseTypeMethod(sp, apply, f, time, cgeom, numComp, func, ctx, value);
1068: PetscFunctionReturn(PETSC_SUCCESS);
1069: }
1071: /*@C
1072: PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by `PetscDualSpaceGetAllData()`
1074: Input Parameters:
1075: + sp - The `PetscDualSpace` object
1076: - pointEval - Evaluation at the points returned by `PetscDualSpaceGetAllData()`
1078: Output Parameter:
1079: . spValue - The values of all dual space functionals
1081: Level: advanced
1083: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1084: @*/
1085: PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1086: {
1087: PetscFunctionBegin;
1089: PetscUseTypeMethod(sp, applyall, pointEval, spValue);
1090: PetscFunctionReturn(PETSC_SUCCESS);
1091: }
1093: /*@C
1094: PetscDualSpaceApplyInterior - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by `PetscDualSpaceGetInteriorData()`
1096: Input Parameters:
1097: + sp - The `PetscDualSpace` object
1098: - pointEval - Evaluation at the points returned by `PetscDualSpaceGetInteriorData()`
1100: Output Parameter:
1101: . spValue - The values of interior dual space functionals
1103: Level: advanced
1105: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1106: @*/
1107: PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1108: {
1109: PetscFunctionBegin;
1111: PetscUseTypeMethod(sp, applyint, pointEval, spValue);
1112: PetscFunctionReturn(PETSC_SUCCESS);
1113: }
1115: /*@C
1116: PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
1118: Input Parameters:
1119: + sp - The `PetscDualSpace` object
1120: . f - The basis functional index
1121: . time - The time
1122: . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
1123: . Nc - The number of components for the function
1124: . func - The input function
1125: - ctx - A context for the function
1127: Output Parameter:
1128: . value - The output value
1130: Calling sequence:
1131: .vb
1132: PetscErrorCode func(PetscInt dim, PetscReal time, const PetscReal x[],PetscInt numComponents, PetscScalar values[], void *ctx)
1133: .ve
1135: Level: advanced
1137: Note:
1138: 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.
1140: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1141: @*/
1142: PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFEGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value)
1143: {
1144: DM dm;
1145: PetscQuadrature n;
1146: const PetscReal *points, *weights;
1147: PetscReal x[3];
1148: PetscScalar *val;
1149: PetscInt dim, dE, qNc, c, Nq, q;
1150: PetscBool isAffine;
1152: PetscFunctionBegin;
1154: PetscAssertPointer(value, 8);
1155: PetscCall(PetscDualSpaceGetDM(sp, &dm));
1156: PetscCall(PetscDualSpaceGetFunctional(sp, f, &n));
1157: PetscCall(PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights));
1158: 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);
1159: PetscCheck(qNc == Nc, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_SIZ, "The quadrature components %" PetscInt_FMT " != function components %" PetscInt_FMT, qNc, Nc);
1160: PetscCall(DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val));
1161: *value = 0.0;
1162: isAffine = cgeom->isAffine;
1163: dE = cgeom->dimEmbed;
1164: for (q = 0; q < Nq; ++q) {
1165: if (isAffine) {
1166: CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q * dim], x);
1167: PetscCall((*func)(dE, time, x, Nc, val, ctx));
1168: } else {
1169: PetscCall((*func)(dE, time, &cgeom->v[dE * q], Nc, val, ctx));
1170: }
1171: for (c = 0; c < Nc; ++c) *value += val[c] * weights[q * Nc + c];
1172: }
1173: PetscCall(DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val));
1174: PetscFunctionReturn(PETSC_SUCCESS);
1175: }
1177: /*@C
1178: PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by `PetscDualSpaceGetAllData()`
1180: Input Parameters:
1181: + sp - The `PetscDualSpace` object
1182: - pointEval - Evaluation at the points returned by `PetscDualSpaceGetAllData()`
1184: Output Parameter:
1185: . spValue - The values of all dual space functionals
1187: Level: advanced
1189: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1190: @*/
1191: PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1192: {
1193: Vec pointValues, dofValues;
1194: Mat allMat;
1196: PetscFunctionBegin;
1198: PetscAssertPointer(pointEval, 2);
1199: PetscAssertPointer(spValue, 3);
1200: PetscCall(PetscDualSpaceGetAllData(sp, NULL, &allMat));
1201: if (!sp->allNodeValues) PetscCall(MatCreateVecs(allMat, &sp->allNodeValues, NULL));
1202: pointValues = sp->allNodeValues;
1203: if (!sp->allDofValues) PetscCall(MatCreateVecs(allMat, NULL, &sp->allDofValues));
1204: dofValues = sp->allDofValues;
1205: PetscCall(VecPlaceArray(pointValues, pointEval));
1206: PetscCall(VecPlaceArray(dofValues, spValue));
1207: PetscCall(MatMult(allMat, pointValues, dofValues));
1208: PetscCall(VecResetArray(dofValues));
1209: PetscCall(VecResetArray(pointValues));
1210: PetscFunctionReturn(PETSC_SUCCESS);
1211: }
1213: /*@C
1214: PetscDualSpaceApplyInteriorDefault - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by `PetscDualSpaceGetInteriorData()`
1216: Input Parameters:
1217: + sp - The `PetscDualSpace` object
1218: - pointEval - Evaluation at the points returned by `PetscDualSpaceGetInteriorData()`
1220: Output Parameter:
1221: . spValue - The values of interior dual space functionals
1223: Level: advanced
1225: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1226: @*/
1227: PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1228: {
1229: Vec pointValues, dofValues;
1230: Mat intMat;
1232: PetscFunctionBegin;
1234: PetscAssertPointer(pointEval, 2);
1235: PetscAssertPointer(spValue, 3);
1236: PetscCall(PetscDualSpaceGetInteriorData(sp, NULL, &intMat));
1237: if (!sp->intNodeValues) PetscCall(MatCreateVecs(intMat, &sp->intNodeValues, NULL));
1238: pointValues = sp->intNodeValues;
1239: if (!sp->intDofValues) PetscCall(MatCreateVecs(intMat, NULL, &sp->intDofValues));
1240: dofValues = sp->intDofValues;
1241: PetscCall(VecPlaceArray(pointValues, pointEval));
1242: PetscCall(VecPlaceArray(dofValues, spValue));
1243: PetscCall(MatMult(intMat, pointValues, dofValues));
1244: PetscCall(VecResetArray(dofValues));
1245: PetscCall(VecResetArray(pointValues));
1246: PetscFunctionReturn(PETSC_SUCCESS);
1247: }
1249: /*@
1250: PetscDualSpaceGetAllData - Get all quadrature nodes from this space, and the matrix that sends quadrature node values to degree-of-freedom values
1252: Input Parameter:
1253: . sp - The dualspace
1255: Output Parameters:
1256: + allNodes - A `PetscQuadrature` object containing all evaluation nodes
1257: - allMat - A `Mat` for the node-to-dof transformation
1259: Level: advanced
1261: .seealso: `PetscQuadrature`, `PetscDualSpace`, `PetscDualSpaceCreate()`, `Mat`
1262: @*/
1263: PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat)
1264: {
1265: PetscFunctionBegin;
1267: if (allNodes) PetscAssertPointer(allNodes, 2);
1268: if (allMat) PetscAssertPointer(allMat, 3);
1269: if ((!sp->allNodes || !sp->allMat) && sp->ops->createalldata) {
1270: PetscQuadrature qpoints;
1271: Mat amat;
1273: PetscUseTypeMethod(sp, createalldata, &qpoints, &amat);
1274: PetscCall(PetscQuadratureDestroy(&sp->allNodes));
1275: PetscCall(MatDestroy(&sp->allMat));
1276: sp->allNodes = qpoints;
1277: sp->allMat = amat;
1278: }
1279: if (allNodes) *allNodes = sp->allNodes;
1280: if (allMat) *allMat = sp->allMat;
1281: PetscFunctionReturn(PETSC_SUCCESS);
1282: }
1284: /*@
1285: PetscDualSpaceCreateAllDataDefault - Create all evaluation nodes and the node-to-dof matrix by examining functionals
1287: Input Parameter:
1288: . sp - The dualspace
1290: Output Parameters:
1291: + allNodes - A `PetscQuadrature` object containing all evaluation nodes
1292: - allMat - A `Mat` for the node-to-dof transformation
1294: Level: advanced
1296: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`, `Mat`, `PetscQuadrature`
1297: @*/
1298: PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat)
1299: {
1300: PetscInt spdim;
1301: PetscInt numPoints, offset;
1302: PetscReal *points;
1303: PetscInt f, dim;
1304: PetscInt Nc, nrows, ncols;
1305: PetscInt maxNumPoints;
1306: PetscQuadrature q;
1307: Mat A;
1309: PetscFunctionBegin;
1310: PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
1311: PetscCall(PetscDualSpaceGetDimension(sp, &spdim));
1312: if (!spdim) {
1313: PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allNodes));
1314: PetscCall(PetscQuadratureSetData(*allNodes, 0, 0, 0, NULL, NULL));
1315: }
1316: nrows = spdim;
1317: PetscCall(PetscDualSpaceGetFunctional(sp, 0, &q));
1318: PetscCall(PetscQuadratureGetData(q, &dim, NULL, &numPoints, NULL, NULL));
1319: maxNumPoints = numPoints;
1320: for (f = 1; f < spdim; f++) {
1321: PetscInt Np;
1323: PetscCall(PetscDualSpaceGetFunctional(sp, f, &q));
1324: PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Np, NULL, NULL));
1325: numPoints += Np;
1326: maxNumPoints = PetscMax(maxNumPoints, Np);
1327: }
1328: ncols = numPoints * Nc;
1329: PetscCall(PetscMalloc1(dim * numPoints, &points));
1330: PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nrows, ncols, maxNumPoints * Nc, NULL, &A));
1331: for (f = 0, offset = 0; f < spdim; f++) {
1332: const PetscReal *p, *w;
1333: PetscInt Np, i;
1334: PetscInt fnc;
1336: PetscCall(PetscDualSpaceGetFunctional(sp, f, &q));
1337: PetscCall(PetscQuadratureGetData(q, NULL, &fnc, &Np, &p, &w));
1338: PetscCheck(fnc == Nc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "functional component mismatch");
1339: for (i = 0; i < Np * dim; i++) points[offset * dim + i] = p[i];
1340: for (i = 0; i < Np * Nc; i++) PetscCall(MatSetValue(A, f, offset * Nc, w[i], INSERT_VALUES));
1341: offset += Np;
1342: }
1343: PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1344: PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1345: PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allNodes));
1346: PetscCall(PetscQuadratureSetData(*allNodes, dim, 0, numPoints, points, NULL));
1347: *allMat = A;
1348: PetscFunctionReturn(PETSC_SUCCESS);
1349: }
1351: /*@
1352: PetscDualSpaceGetInteriorData - Get all quadrature points necessary to compute the interior degrees of freedom from
1353: this space, as well as the matrix that computes the degrees of freedom from the quadrature
1354: values.
1356: Input Parameter:
1357: . sp - The dualspace
1359: Output Parameters:
1360: + intNodes - A `PetscQuadrature` object containing all evaluation points needed to evaluate interior degrees of freedom
1361: - intMat - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is
1362: the size of the constrained layout (`PetscSectionGetConstrainStorageSize()`) of the dual space section,
1363: npoints is the number of points in intNodes and nc is `PetscDualSpaceGetNumComponents()`.
1365: Level: advanced
1367: Notes:
1368: Degrees of freedom are interior degrees of freedom if they belong (by
1369: `PetscDualSpaceGetSection()`) to interior points in the references, complementary boundary
1370: degrees of freedom are marked as constrained in the section returned by
1371: `PetscDualSpaceGetSection()`).
1373: .seealso: `PetscDualSpace`, `PetscQuadrature`, `Mat`, `PetscDualSpaceCreate()`, `PetscDualSpaceGetDimension()`, `PetscDualSpaceGetNumComponents()`, `PetscQuadratureGetData()`
1374: @*/
1375: PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat)
1376: {
1377: PetscFunctionBegin;
1379: if (intNodes) PetscAssertPointer(intNodes, 2);
1380: if (intMat) PetscAssertPointer(intMat, 3);
1381: if ((!sp->intNodes || !sp->intMat) && sp->ops->createintdata) {
1382: PetscQuadrature qpoints;
1383: Mat imat;
1385: PetscUseTypeMethod(sp, createintdata, &qpoints, &imat);
1386: PetscCall(PetscQuadratureDestroy(&sp->intNodes));
1387: PetscCall(MatDestroy(&sp->intMat));
1388: sp->intNodes = qpoints;
1389: sp->intMat = imat;
1390: }
1391: if (intNodes) *intNodes = sp->intNodes;
1392: if (intMat) *intMat = sp->intMat;
1393: PetscFunctionReturn(PETSC_SUCCESS);
1394: }
1396: /*@
1397: PetscDualSpaceCreateInteriorDataDefault - Create quadrature points by examining interior functionals and create the matrix mapping quadrature point values to interior dual space values
1399: Input Parameter:
1400: . sp - The dualspace
1402: Output Parameters:
1403: + intNodes - A `PetscQuadrature` object containing all evaluation points needed to evaluate interior degrees of freedom
1404: - intMat - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is
1405: the size of the constrained layout (`PetscSectionGetConstrainStorageSize()`) of the dual space section,
1406: npoints is the number of points in allNodes and nc is `PetscDualSpaceGetNumComponents()`.
1408: Level: advanced
1410: .seealso: `PetscDualSpace`, `PetscQuadrature`, `Mat`, `PetscDualSpaceCreate()`, `PetscDualSpaceGetInteriorData()`
1411: @*/
1412: PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat)
1413: {
1414: DM dm;
1415: PetscInt spdim0;
1416: PetscInt Nc;
1417: PetscInt pStart, pEnd, p, f;
1418: PetscSection section;
1419: PetscInt numPoints, offset, matoffset;
1420: PetscReal *points;
1421: PetscInt dim;
1422: PetscInt *nnz;
1423: PetscQuadrature q;
1424: Mat imat;
1426: PetscFunctionBegin;
1428: PetscCall(PetscDualSpaceGetSection(sp, §ion));
1429: PetscCall(PetscSectionGetConstrainedStorageSize(section, &spdim0));
1430: if (!spdim0) {
1431: *intNodes = NULL;
1432: *intMat = NULL;
1433: PetscFunctionReturn(PETSC_SUCCESS);
1434: }
1435: PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
1436: PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1437: PetscCall(PetscDualSpaceGetDM(sp, &dm));
1438: PetscCall(DMGetDimension(dm, &dim));
1439: PetscCall(PetscMalloc1(spdim0, &nnz));
1440: for (p = pStart, f = 0, numPoints = 0; p < pEnd; p++) {
1441: PetscInt dof, cdof, off, d;
1443: PetscCall(PetscSectionGetDof(section, p, &dof));
1444: PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
1445: if (!(dof - cdof)) continue;
1446: PetscCall(PetscSectionGetOffset(section, p, &off));
1447: for (d = 0; d < dof; d++, off++, f++) {
1448: PetscInt Np;
1450: PetscCall(PetscDualSpaceGetFunctional(sp, off, &q));
1451: PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Np, NULL, NULL));
1452: nnz[f] = Np * Nc;
1453: numPoints += Np;
1454: }
1455: }
1456: PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, spdim0, numPoints * Nc, 0, nnz, &imat));
1457: PetscCall(PetscFree(nnz));
1458: PetscCall(PetscMalloc1(dim * numPoints, &points));
1459: for (p = pStart, f = 0, offset = 0, matoffset = 0; p < pEnd; p++) {
1460: PetscInt dof, cdof, off, d;
1462: PetscCall(PetscSectionGetDof(section, p, &dof));
1463: PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
1464: if (!(dof - cdof)) continue;
1465: PetscCall(PetscSectionGetOffset(section, p, &off));
1466: for (d = 0; d < dof; d++, off++, f++) {
1467: const PetscReal *p;
1468: const PetscReal *w;
1469: PetscInt Np, i;
1471: PetscCall(PetscDualSpaceGetFunctional(sp, off, &q));
1472: PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Np, &p, &w));
1473: for (i = 0; i < Np * dim; i++) points[offset + i] = p[i];
1474: for (i = 0; i < Np * Nc; i++) PetscCall(MatSetValue(imat, f, matoffset + i, w[i], INSERT_VALUES));
1475: offset += Np * dim;
1476: matoffset += Np * Nc;
1477: }
1478: }
1479: PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, intNodes));
1480: PetscCall(PetscQuadratureSetData(*intNodes, dim, 0, numPoints, points, NULL));
1481: PetscCall(MatAssemblyBegin(imat, MAT_FINAL_ASSEMBLY));
1482: PetscCall(MatAssemblyEnd(imat, MAT_FINAL_ASSEMBLY));
1483: *intMat = imat;
1484: PetscFunctionReturn(PETSC_SUCCESS);
1485: }
1487: /*@
1488: PetscDualSpaceEqual - Determine if two dual spaces are equivalent
1490: Input Parameters:
1491: + A - A `PetscDualSpace` object
1492: - B - Another `PetscDualSpace` object
1494: Output Parameter:
1495: . equal - `PETSC_TRUE` if the dual spaces are equivalent
1497: Level: advanced
1499: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1500: @*/
1501: PetscErrorCode PetscDualSpaceEqual(PetscDualSpace A, PetscDualSpace B, PetscBool *equal)
1502: {
1503: PetscInt sizeA, sizeB, dimA, dimB;
1504: const PetscInt *dofA, *dofB;
1505: PetscQuadrature quadA, quadB;
1506: Mat matA, matB;
1508: PetscFunctionBegin;
1511: PetscAssertPointer(equal, 3);
1512: *equal = PETSC_FALSE;
1513: PetscCall(PetscDualSpaceGetDimension(A, &sizeA));
1514: PetscCall(PetscDualSpaceGetDimension(B, &sizeB));
1515: if (sizeB != sizeA) PetscFunctionReturn(PETSC_SUCCESS);
1516: PetscCall(DMGetDimension(A->dm, &dimA));
1517: PetscCall(DMGetDimension(B->dm, &dimB));
1518: if (dimA != dimB) PetscFunctionReturn(PETSC_SUCCESS);
1520: PetscCall(PetscDualSpaceGetNumDof(A, &dofA));
1521: PetscCall(PetscDualSpaceGetNumDof(B, &dofB));
1522: for (PetscInt d = 0; d < dimA; d++) {
1523: if (dofA[d] != dofB[d]) PetscFunctionReturn(PETSC_SUCCESS);
1524: }
1526: PetscCall(PetscDualSpaceGetInteriorData(A, &quadA, &matA));
1527: PetscCall(PetscDualSpaceGetInteriorData(B, &quadB, &matB));
1528: if (!quadA && !quadB) {
1529: *equal = PETSC_TRUE;
1530: } else if (quadA && quadB) {
1531: PetscCall(PetscQuadratureEqual(quadA, quadB, equal));
1532: if (*equal == PETSC_FALSE) PetscFunctionReturn(PETSC_SUCCESS);
1533: if (!matA && !matB) PetscFunctionReturn(PETSC_SUCCESS);
1534: if (matA && matB) PetscCall(MatEqual(matA, matB, equal));
1535: else *equal = PETSC_FALSE;
1536: }
1537: PetscFunctionReturn(PETSC_SUCCESS);
1538: }
1540: /*@C
1541: PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
1543: Input Parameters:
1544: + sp - The `PetscDualSpace` object
1545: . f - The basis functional index
1546: . time - The time
1547: . cgeom - A context with geometric information for this cell, we currently just use the centroid
1548: . Nc - The number of components for the function
1549: . func - The input function
1550: - ctx - A context for the function
1552: Output Parameter:
1553: . value - The output value (scalar)
1555: Calling sequence:
1556: .vb
1557: PetscErrorCode func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt numComponents, PetscScalar values[], void *ctx)
1558: .ve
1560: Level: advanced
1562: Note:
1563: 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.
1565: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1566: @*/
1567: PetscErrorCode PetscDualSpaceApplyFVM(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFVCellGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value)
1568: {
1569: DM dm;
1570: PetscQuadrature n;
1571: const PetscReal *points, *weights;
1572: PetscScalar *val;
1573: PetscInt dimEmbed, qNc, c, Nq, q;
1575: PetscFunctionBegin;
1577: PetscAssertPointer(value, 8);
1578: PetscCall(PetscDualSpaceGetDM(sp, &dm));
1579: PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
1580: PetscCall(PetscDualSpaceGetFunctional(sp, f, &n));
1581: PetscCall(PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights));
1582: PetscCheck(qNc == Nc, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_SIZ, "The quadrature components %" PetscInt_FMT " != function components %" PetscInt_FMT, qNc, Nc);
1583: PetscCall(DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val));
1584: *value = 0.;
1585: for (q = 0; q < Nq; ++q) {
1586: PetscCall((*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx));
1587: for (c = 0; c < Nc; ++c) *value += val[c] * weights[q * Nc + c];
1588: }
1589: PetscCall(DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val));
1590: PetscFunctionReturn(PETSC_SUCCESS);
1591: }
1593: /*@
1594: PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
1595: given height. This assumes that the reference cell is symmetric over points of this height.
1597: Not Collective
1599: Input Parameters:
1600: + sp - the `PetscDualSpace` object
1601: - height - the height of the mesh point for which the subspace is desired
1603: Output Parameter:
1604: . subsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the
1605: point, which will be of lesser dimension if height > 0.
1607: Level: advanced
1609: Notes:
1610: If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
1611: pointwise values are not defined on the element boundaries), or if the implementation of `PetscDualSpace` does not
1612: support extracting subspaces, then NULL is returned.
1614: This does not increment the reference count on the returned dual space, and the user should not destroy it.
1616: .seealso: `PetscDualSpace`, `PetscSpaceGetHeightSubspace()`, `PetscDualSpaceGetPointSubspace()`
1617: @*/
1618: PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
1619: {
1620: PetscInt depth = -1, cStart, cEnd;
1621: DM dm;
1623: PetscFunctionBegin;
1625: PetscAssertPointer(subsp, 3);
1626: 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");
1627: *subsp = NULL;
1628: dm = sp->dm;
1629: PetscCall(DMPlexGetDepth(dm, &depth));
1630: PetscCheck(height >= 0 && height <= depth, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid height");
1631: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1632: if (height == 0 && cEnd == cStart + 1) {
1633: *subsp = sp;
1634: PetscFunctionReturn(PETSC_SUCCESS);
1635: }
1636: if (!sp->heightSpaces) {
1637: PetscInt h;
1638: PetscCall(PetscCalloc1(depth + 1, &sp->heightSpaces));
1640: for (h = 0; h <= depth; h++) {
1641: if (h == 0 && cEnd == cStart + 1) continue;
1642: if (sp->ops->createheightsubspace) PetscUseTypeMethod(sp, createheightsubspace, height, &sp->heightSpaces[h]);
1643: else if (sp->pointSpaces) {
1644: PetscInt hStart, hEnd;
1646: PetscCall(DMPlexGetHeightStratum(dm, h, &hStart, &hEnd));
1647: if (hEnd > hStart) {
1648: const char *name;
1650: PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[hStart]));
1651: if (sp->pointSpaces[hStart]) {
1652: PetscCall(PetscObjectGetName((PetscObject)sp, &name));
1653: PetscCall(PetscObjectSetName((PetscObject)sp->pointSpaces[hStart], name));
1654: }
1655: sp->heightSpaces[h] = sp->pointSpaces[hStart];
1656: }
1657: }
1658: }
1659: }
1660: *subsp = sp->heightSpaces[height];
1661: PetscFunctionReturn(PETSC_SUCCESS);
1662: }
1664: /*@
1665: PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.
1667: Not Collective
1669: Input Parameters:
1670: + sp - the `PetscDualSpace` object
1671: - point - the point (in the dual space's DM) for which the subspace is desired
1673: Output Parameters:
1674: . bdsp - the subspace.
1676: Level: advanced
1678: Notes:
1679: The functionals in the subspace are with respect to the intrinsic geometry of the point,
1680: which will be of lesser dimension if height > 0.
1682: If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
1683: defined on the element boundaries), or if the implementation of `PetscDualSpace` does not support extracting
1684: subspaces, then `NULL` is returned.
1686: This does not increment the reference count on the returned dual space, and the user should not destroy it.
1688: .seealso: `PetscDualSpace`, `PetscDualSpaceGetHeightSubspace()`
1689: @*/
1690: PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
1691: {
1692: PetscInt pStart = 0, pEnd = 0, cStart, cEnd;
1693: DM dm;
1695: PetscFunctionBegin;
1697: PetscAssertPointer(bdsp, 3);
1698: *bdsp = NULL;
1699: dm = sp->dm;
1700: PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1701: PetscCheck(point >= pStart && point <= pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point");
1702: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1703: 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 */
1704: *bdsp = sp;
1705: PetscFunctionReturn(PETSC_SUCCESS);
1706: }
1707: if (!sp->pointSpaces) {
1708: PetscInt p;
1709: PetscCall(PetscCalloc1(pEnd - pStart, &sp->pointSpaces));
1711: for (p = 0; p < pEnd - pStart; p++) {
1712: if (p + pStart == cStart && cEnd == cStart + 1) continue;
1713: if (sp->ops->createpointsubspace) PetscUseTypeMethod(sp, createpointsubspace, p + pStart, &sp->pointSpaces[p]);
1714: else if (sp->heightSpaces || sp->ops->createheightsubspace) {
1715: PetscInt dim, depth, height;
1716: DMLabel label;
1718: PetscCall(DMPlexGetDepth(dm, &dim));
1719: PetscCall(DMPlexGetDepthLabel(dm, &label));
1720: PetscCall(DMLabelGetValue(label, p + pStart, &depth));
1721: height = dim - depth;
1722: PetscCall(PetscDualSpaceGetHeightSubspace(sp, height, &sp->pointSpaces[p]));
1723: PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[p]));
1724: }
1725: }
1726: }
1727: *bdsp = sp->pointSpaces[point - pStart];
1728: PetscFunctionReturn(PETSC_SUCCESS);
1729: }
1731: /*@C
1732: PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis
1734: Not Collective
1736: Input Parameter:
1737: . sp - the `PetscDualSpace` object
1739: Output Parameters:
1740: + perms - Permutations of the interior degrees of freedom, parameterized by the point orientation
1741: - flips - Sign reversal of the interior degrees of freedom, parameterized by the point orientation
1743: Level: developer
1745: Note:
1746: The permutation and flip arrays are organized in the following way
1747: .vb
1748: perms[p][ornt][dof # on point] = new local dof #
1749: flips[p][ornt][dof # on point] = reversal or not
1750: .ve
1752: .seealso: `PetscDualSpace`
1753: @*/
1754: PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
1755: {
1756: PetscFunctionBegin;
1758: if (perms) {
1759: PetscAssertPointer(perms, 2);
1760: *perms = NULL;
1761: }
1762: if (flips) {
1763: PetscAssertPointer(flips, 3);
1764: *flips = NULL;
1765: }
1766: PetscTryTypeMethod(sp, getsymmetries, perms, flips);
1767: PetscFunctionReturn(PETSC_SUCCESS);
1768: }
1770: /*@
1771: PetscDualSpaceGetFormDegree - Get the form degree k for the k-form the describes the pushforwards/pullbacks of this
1772: dual space's functionals.
1774: Input Parameter:
1775: . dsp - The `PetscDualSpace`
1777: Output Parameter:
1778: . k - The *signed* degree k of the k. If k >= 0, this means that the degrees of freedom are k-forms, and are stored
1779: in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms. So for example,
1780: 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).
1781: If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the
1782: Hodge star map. So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms
1783: but are stored as 1-forms.
1785: Level: developer
1787: .seealso: `PetscDualSpace`, `PetscDTAltV`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceTransformType`
1788: @*/
1789: PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace dsp, PetscInt *k)
1790: {
1791: PetscFunctionBeginHot;
1793: PetscAssertPointer(k, 2);
1794: *k = dsp->k;
1795: PetscFunctionReturn(PETSC_SUCCESS);
1796: }
1798: /*@
1799: PetscDualSpaceSetFormDegree - Set the form degree k for the k-form the describes the pushforwards/pullbacks of this
1800: dual space's functionals.
1802: Input Parameters:
1803: + dsp - The `PetscDualSpace`
1804: - k - The *signed* degree k of the k. If k >= 0, this means that the degrees of freedom are k-forms, and are stored
1805: in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms. So for example,
1806: 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).
1807: If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the
1808: Hodge star map. So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms
1809: but are stored as 1-forms.
1811: Level: developer
1813: .seealso: `PetscDualSpace`, `PetscDTAltV`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceTransformType`
1814: @*/
1815: PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace dsp, PetscInt k)
1816: {
1817: PetscInt dim;
1819: PetscFunctionBeginHot;
1821: PetscCheck(!dsp->setupcalled, PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up");
1822: dim = dsp->dm->dim;
1823: 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);
1824: dsp->k = k;
1825: PetscFunctionReturn(PETSC_SUCCESS);
1826: }
1828: /*@
1829: PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space
1831: Input Parameter:
1832: . dsp - The `PetscDualSpace`
1834: Output Parameter:
1835: . k - The simplex dimension
1837: Level: developer
1839: Note:
1840: Currently supported values are
1841: .vb
1842: 0: These are H_1 methods that only transform coordinates
1843: 1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM)
1844: 2: These are the same as 1
1845: 3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM)
1846: .ve
1848: .seealso: `PetscDualSpace`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceTransformType`
1849: @*/
1850: PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k)
1851: {
1852: PetscInt dim;
1854: PetscFunctionBeginHot;
1856: PetscAssertPointer(k, 2);
1857: dim = dsp->dm->dim;
1858: if (!dsp->k) *k = IDENTITY_TRANSFORM;
1859: else if (dsp->k == 1) *k = COVARIANT_PIOLA_TRANSFORM;
1860: else if (dsp->k == -(dim - 1)) *k = CONTRAVARIANT_PIOLA_TRANSFORM;
1861: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported transformation");
1862: PetscFunctionReturn(PETSC_SUCCESS);
1863: }
1865: /*@C
1866: PetscDualSpaceTransform - Transform the function values
1868: Input Parameters:
1869: + dsp - The `PetscDualSpace`
1870: . trans - The type of transform
1871: . isInverse - Flag to invert the transform
1872: . fegeom - The cell geometry
1873: . Nv - The number of function samples
1874: . Nc - The number of function components
1875: - vals - The function values
1877: Output Parameter:
1878: . vals - The transformed function values
1880: Level: intermediate
1882: Note:
1883: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
1885: .seealso: `PetscDualSpace`, `PetscDualSpaceTransformGradient()`, `PetscDualSpaceTransformHessian()`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransformType`
1886: @*/
1887: PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1888: {
1889: PetscReal Jstar[9] = {0};
1890: PetscInt dim, v, c, Nk;
1892: PetscFunctionBeginHot;
1894: PetscAssertPointer(fegeom, 4);
1895: PetscAssertPointer(vals, 7);
1896: /* TODO: not handling dimEmbed != dim right now */
1897: dim = dsp->dm->dim;
1898: /* No change needed for 0-forms */
1899: if (!dsp->k) PetscFunctionReturn(PETSC_SUCCESS);
1900: PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(dsp->k), &Nk));
1901: /* TODO: use fegeom->isAffine */
1902: PetscCall(PetscDTAltVPullbackMatrix(dim, dim, isInverse ? fegeom->J : fegeom->invJ, dsp->k, Jstar));
1903: for (v = 0; v < Nv; ++v) {
1904: switch (Nk) {
1905: case 1:
1906: for (c = 0; c < Nc; c++) vals[v * Nc + c] *= Jstar[0];
1907: break;
1908: case 2:
1909: for (c = 0; c < Nc; c += 2) DMPlex_Mult2DReal_Internal(Jstar, 1, &vals[v * Nc + c], &vals[v * Nc + c]);
1910: break;
1911: case 3:
1912: for (c = 0; c < Nc; c += 3) DMPlex_Mult3DReal_Internal(Jstar, 1, &vals[v * Nc + c], &vals[v * Nc + c]);
1913: break;
1914: default:
1915: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported form size %" PetscInt_FMT " for transformation", Nk);
1916: }
1917: }
1918: PetscFunctionReturn(PETSC_SUCCESS);
1919: }
1921: /*@C
1922: PetscDualSpaceTransformGradient - Transform the function gradient values
1924: Input Parameters:
1925: + dsp - The `PetscDualSpace`
1926: . trans - The type of transform
1927: . isInverse - Flag to invert the transform
1928: . fegeom - The cell geometry
1929: . Nv - The number of function gradient samples
1930: . Nc - The number of function components
1931: - vals - The function gradient values
1933: Output Parameter:
1934: . vals - The transformed function gradient values
1936: Level: intermediate
1938: Note:
1939: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
1941: .seealso: `PetscDualSpace`, `PetscDualSpaceTransform()`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransformType`
1942: @*/
1943: PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1944: {
1945: const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed;
1946: PetscInt v, c, d;
1948: PetscFunctionBeginHot;
1950: PetscAssertPointer(fegeom, 4);
1951: PetscAssertPointer(vals, 7);
1952: PetscAssert(dE > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %" PetscInt_FMT, dE);
1953: /* Transform gradient */
1954: if (dim == dE) {
1955: for (v = 0; v < Nv; ++v) {
1956: for (c = 0; c < Nc; ++c) {
1957: switch (dim) {
1958: case 1:
1959: vals[(v * Nc + c) * dim] *= fegeom->invJ[0];
1960: break;
1961: case 2:
1962: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[(v * Nc + c) * dim], &vals[(v * Nc + c) * dim]);
1963: break;
1964: case 3:
1965: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[(v * Nc + c) * dim], &vals[(v * Nc + c) * dim]);
1966: break;
1967: default:
1968: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
1969: }
1970: }
1971: }
1972: } else {
1973: for (v = 0; v < Nv; ++v) {
1974: for (c = 0; c < Nc; ++c) DMPlex_MultTransposeReal_Internal(fegeom->invJ, dim, dE, 1, &vals[(v * Nc + c) * dE], &vals[(v * Nc + c) * dE]);
1975: }
1976: }
1977: /* Assume its a vector, otherwise assume its a bunch of scalars */
1978: if (Nc == 1 || Nc != dim) PetscFunctionReturn(PETSC_SUCCESS);
1979: switch (trans) {
1980: case IDENTITY_TRANSFORM:
1981: break;
1982: case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
1983: if (isInverse) {
1984: for (v = 0; v < Nv; ++v) {
1985: for (d = 0; d < dim; ++d) {
1986: switch (dim) {
1987: case 2:
1988: DMPlex_MultTranspose2DReal_Internal(fegeom->J, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
1989: break;
1990: case 3:
1991: DMPlex_MultTranspose3DReal_Internal(fegeom->J, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
1992: break;
1993: default:
1994: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
1995: }
1996: }
1997: }
1998: } else {
1999: for (v = 0; v < Nv; ++v) {
2000: for (d = 0; d < dim; ++d) {
2001: switch (dim) {
2002: case 2:
2003: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2004: break;
2005: case 3:
2006: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2007: break;
2008: default:
2009: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2010: }
2011: }
2012: }
2013: }
2014: break;
2015: case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
2016: if (isInverse) {
2017: for (v = 0; v < Nv; ++v) {
2018: for (d = 0; d < dim; ++d) {
2019: switch (dim) {
2020: case 2:
2021: DMPlex_Mult2DReal_Internal(fegeom->invJ, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2022: break;
2023: case 3:
2024: DMPlex_Mult3DReal_Internal(fegeom->invJ, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2025: break;
2026: default:
2027: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2028: }
2029: for (c = 0; c < Nc; ++c) vals[(v * Nc + c) * dim + d] *= fegeom->detJ[0];
2030: }
2031: }
2032: } else {
2033: for (v = 0; v < Nv; ++v) {
2034: for (d = 0; d < dim; ++d) {
2035: switch (dim) {
2036: case 2:
2037: DMPlex_Mult2DReal_Internal(fegeom->J, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2038: break;
2039: case 3:
2040: DMPlex_Mult3DReal_Internal(fegeom->J, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2041: break;
2042: default:
2043: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2044: }
2045: for (c = 0; c < Nc; ++c) vals[(v * Nc + c) * dim + d] /= fegeom->detJ[0];
2046: }
2047: }
2048: }
2049: break;
2050: }
2051: PetscFunctionReturn(PETSC_SUCCESS);
2052: }
2054: /*@C
2055: PetscDualSpaceTransformHessian - Transform the function Hessian values
2057: Input Parameters:
2058: + dsp - The `PetscDualSpace`
2059: . trans - The type of transform
2060: . isInverse - Flag to invert the transform
2061: . fegeom - The cell geometry
2062: . Nv - The number of function Hessian samples
2063: . Nc - The number of function components
2064: - vals - The function gradient values
2066: Output Parameter:
2067: . vals - The transformed function Hessian values
2069: Level: intermediate
2071: Note:
2072: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2074: .seealso: `PetscDualSpace`, `PetscDualSpaceTransform()`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransformType`
2075: @*/
2076: PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
2077: {
2078: const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed;
2079: PetscInt v, c;
2081: PetscFunctionBeginHot;
2083: PetscAssertPointer(fegeom, 4);
2084: PetscAssertPointer(vals, 7);
2085: PetscAssert(dE > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %" PetscInt_FMT, dE);
2086: /* Transform Hessian: J^{-T}_{ik} J^{-T}_{jl} H(f)_{kl} = J^{-T}_{ik} H(f)_{kl} J^{-1}_{lj} */
2087: if (dim == dE) {
2088: for (v = 0; v < Nv; ++v) {
2089: for (c = 0; c < Nc; ++c) {
2090: switch (dim) {
2091: case 1:
2092: vals[(v * Nc + c) * dim * dim] *= PetscSqr(fegeom->invJ[0]);
2093: break;
2094: case 2:
2095: DMPlex_PTAP2DReal_Internal(fegeom->invJ, &vals[(v * Nc + c) * dim * dim], &vals[(v * Nc + c) * dim * dim]);
2096: break;
2097: case 3:
2098: DMPlex_PTAP3DReal_Internal(fegeom->invJ, &vals[(v * Nc + c) * dim * dim], &vals[(v * Nc + c) * dim * dim]);
2099: break;
2100: default:
2101: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2102: }
2103: }
2104: }
2105: } else {
2106: for (v = 0; v < Nv; ++v) {
2107: 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]);
2108: }
2109: }
2110: /* Assume its a vector, otherwise assume its a bunch of scalars */
2111: if (Nc == 1 || Nc != dim) PetscFunctionReturn(PETSC_SUCCESS);
2112: switch (trans) {
2113: case IDENTITY_TRANSFORM:
2114: break;
2115: case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
2116: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Piola mapping for Hessians not yet supported");
2117: case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
2118: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Piola mapping for Hessians not yet supported");
2119: }
2120: PetscFunctionReturn(PETSC_SUCCESS);
2121: }
2123: /*@C
2124: 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.
2126: Input Parameters:
2127: + dsp - The `PetscDualSpace`
2128: . fegeom - The geometry for this cell
2129: . Nq - The number of function samples
2130: . Nc - The number of function components
2131: - pointEval - The function values
2133: Output Parameter:
2134: . pointEval - The transformed function values
2136: Level: advanced
2138: Notes:
2139: 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.
2141: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2143: .seealso: `PetscDualSpace`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2144: @*/
2145: PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2146: {
2147: PetscDualSpaceTransformType trans;
2148: PetscInt k;
2150: PetscFunctionBeginHot;
2152: PetscAssertPointer(fegeom, 2);
2153: PetscAssertPointer(pointEval, 5);
2154: /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2155: This determines their transformation properties. */
2156: PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2157: switch (k) {
2158: case 0: /* H^1 point evaluations */
2159: trans = IDENTITY_TRANSFORM;
2160: break;
2161: case 1: /* Hcurl preserves tangential edge traces */
2162: trans = COVARIANT_PIOLA_TRANSFORM;
2163: break;
2164: case 2:
2165: case 3: /* Hdiv preserve normal traces */
2166: trans = CONTRAVARIANT_PIOLA_TRANSFORM;
2167: break;
2168: default:
2169: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2170: }
2171: PetscCall(PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, pointEval));
2172: PetscFunctionReturn(PETSC_SUCCESS);
2173: }
2175: /*@C
2176: 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.
2178: Input Parameters:
2179: + dsp - The `PetscDualSpace`
2180: . fegeom - The geometry for this cell
2181: . Nq - The number of function samples
2182: . Nc - The number of function components
2183: - pointEval - The function values
2185: Output Parameter:
2186: . pointEval - The transformed function values
2188: Level: advanced
2190: Notes:
2191: 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.
2193: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2195: .seealso: `PetscDualSpace`, `PetscDualSpacePullback()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2196: @*/
2197: PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2198: {
2199: PetscDualSpaceTransformType trans;
2200: PetscInt k;
2202: PetscFunctionBeginHot;
2204: PetscAssertPointer(fegeom, 2);
2205: PetscAssertPointer(pointEval, 5);
2206: /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2207: This determines their transformation properties. */
2208: PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2209: switch (k) {
2210: case 0: /* H^1 point evaluations */
2211: trans = IDENTITY_TRANSFORM;
2212: break;
2213: case 1: /* Hcurl preserves tangential edge traces */
2214: trans = COVARIANT_PIOLA_TRANSFORM;
2215: break;
2216: case 2:
2217: case 3: /* Hdiv preserve normal traces */
2218: trans = CONTRAVARIANT_PIOLA_TRANSFORM;
2219: break;
2220: default:
2221: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2222: }
2223: PetscCall(PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval));
2224: PetscFunctionReturn(PETSC_SUCCESS);
2225: }
2227: /*@C
2228: 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.
2230: Input Parameters:
2231: + dsp - The `PetscDualSpace`
2232: . fegeom - The geometry for this cell
2233: . Nq - The number of function gradient samples
2234: . Nc - The number of function components
2235: - pointEval - The function gradient values
2237: Output Parameter:
2238: . pointEval - The transformed function gradient values
2240: Level: advanced
2242: Notes:
2243: 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.
2245: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2247: .seealso: `PetscDualSpace`, `PetscDualSpacePushforward()`, `PPetscDualSpacePullback()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2248: @*/
2249: PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2250: {
2251: PetscDualSpaceTransformType trans;
2252: PetscInt k;
2254: PetscFunctionBeginHot;
2256: PetscAssertPointer(fegeom, 2);
2257: PetscAssertPointer(pointEval, 5);
2258: /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2259: This determines their transformation properties. */
2260: PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2261: switch (k) {
2262: case 0: /* H^1 point evaluations */
2263: trans = IDENTITY_TRANSFORM;
2264: break;
2265: case 1: /* Hcurl preserves tangential edge traces */
2266: trans = COVARIANT_PIOLA_TRANSFORM;
2267: break;
2268: case 2:
2269: case 3: /* Hdiv preserve normal traces */
2270: trans = CONTRAVARIANT_PIOLA_TRANSFORM;
2271: break;
2272: default:
2273: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2274: }
2275: PetscCall(PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval));
2276: PetscFunctionReturn(PETSC_SUCCESS);
2277: }
2279: /*@C
2280: 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.
2282: Input Parameters:
2283: + dsp - The `PetscDualSpace`
2284: . fegeom - The geometry for this cell
2285: . Nq - The number of function Hessian samples
2286: . Nc - The number of function components
2287: - pointEval - The function gradient values
2289: Output Parameter:
2290: . pointEval - The transformed function Hessian values
2292: Level: advanced
2294: Notes:
2295: 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.
2297: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2299: .seealso: `PetscDualSpace`, `PetscDualSpacePushforward()`, `PPetscDualSpacePullback()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2300: @*/
2301: PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2302: {
2303: PetscDualSpaceTransformType trans;
2304: PetscInt k;
2306: PetscFunctionBeginHot;
2308: PetscAssertPointer(fegeom, 2);
2309: PetscAssertPointer(pointEval, 5);
2310: /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2311: This determines their transformation properties. */
2312: PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2313: switch (k) {
2314: case 0: /* H^1 point evaluations */
2315: trans = IDENTITY_TRANSFORM;
2316: break;
2317: case 1: /* Hcurl preserves tangential edge traces */
2318: trans = COVARIANT_PIOLA_TRANSFORM;
2319: break;
2320: case 2:
2321: case 3: /* Hdiv preserve normal traces */
2322: trans = CONTRAVARIANT_PIOLA_TRANSFORM;
2323: break;
2324: default:
2325: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2326: }
2327: PetscCall(PetscDualSpaceTransformHessian(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval));
2328: PetscFunctionReturn(PETSC_SUCCESS);
2329: }