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; 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, 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 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: /*@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 - 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: /*@C
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: /*@C
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: /*@C
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: PetscCall(PetscFEInitializePackage());
452: PetscCall(PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView));
453: s->order = 0;
454: s->Nc = 1;
455: s->k = 0;
456: s->spdim = -1;
457: s->spintdim = -1;
458: s->uniform = PETSC_TRUE;
459: s->setupcalled = PETSC_FALSE;
460: *sp = s;
461: PetscFunctionReturn(PETSC_SUCCESS);
462: }
464: /*@C
465: PetscDualSpaceDuplicate - Creates a duplicate `PetscDualSpace` object that is not setup.
467: Collective
469: Input Parameter:
470: . sp - The original `PetscDualSpace`
472: Output Parameter:
473: . spNew - The duplicate `PetscDualSpace`
475: Level: beginner
477: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`
478: @*/
479: PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
480: {
481: DM dm;
482: PetscDualSpaceType type;
483: const char *name;
485: PetscFunctionBegin;
487: PetscAssertPointer(spNew, 2);
488: PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)sp), spNew));
489: name = ((PetscObject)sp)->name;
490: if (name) { PetscCall(PetscObjectSetName((PetscObject)*spNew, name)); }
491: PetscCall(PetscDualSpaceGetType(sp, &type));
492: PetscCall(PetscDualSpaceSetType(*spNew, type));
493: PetscCall(PetscDualSpaceGetDM(sp, &dm));
494: PetscCall(PetscDualSpaceSetDM(*spNew, dm));
496: (*spNew)->order = sp->order;
497: (*spNew)->k = sp->k;
498: (*spNew)->Nc = sp->Nc;
499: (*spNew)->uniform = sp->uniform;
500: PetscTryTypeMethod(sp, duplicate, *spNew);
501: PetscFunctionReturn(PETSC_SUCCESS);
502: }
504: /*@C
505: PetscDualSpaceGetDM - Get the `DM` representing the reference cell of a `PetscDualSpace`
507: Not Collective
509: Input Parameter:
510: . sp - The `PetscDualSpace`
512: Output Parameter:
513: . dm - The reference cell, that is a `DM` that consists of a single cell
515: Level: intermediate
517: .seealso: `PetscDualSpace`, `PetscDualSpaceSetDM()`, `PetscDualSpaceCreate()`
518: @*/
519: PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
520: {
521: PetscFunctionBegin;
523: PetscAssertPointer(dm, 2);
524: *dm = sp->dm;
525: PetscFunctionReturn(PETSC_SUCCESS);
526: }
528: /*@C
529: PetscDualSpaceSetDM - Get the `DM` representing the reference cell
531: Not Collective
533: Input Parameters:
534: + sp - The `PetscDual`Space
535: - dm - The reference cell
537: Level: intermediate
539: .seealso: `PetscDualSpace`, `DM`, `PetscDualSpaceGetDM()`, `PetscDualSpaceCreate()`
540: @*/
541: PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
542: {
543: PetscFunctionBegin;
546: PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change DM after dualspace is set up");
547: PetscCall(PetscObjectReference((PetscObject)dm));
548: if (sp->dm && sp->dm != dm) PetscCall(PetscDualSpaceClearDMData_Internal(sp, sp->dm));
549: PetscCall(DMDestroy(&sp->dm));
550: sp->dm = dm;
551: PetscFunctionReturn(PETSC_SUCCESS);
552: }
554: /*@C
555: PetscDualSpaceGetOrder - Get the order of the dual space
557: Not Collective
559: Input Parameter:
560: . sp - The `PetscDualSpace`
562: Output Parameter:
563: . order - The order
565: Level: intermediate
567: .seealso: `PetscDualSpace`, `PetscDualSpaceSetOrder()`, `PetscDualSpaceCreate()`
568: @*/
569: PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
570: {
571: PetscFunctionBegin;
573: PetscAssertPointer(order, 2);
574: *order = sp->order;
575: PetscFunctionReturn(PETSC_SUCCESS);
576: }
578: /*@C
579: PetscDualSpaceSetOrder - Set the order of the dual space
581: Not Collective
583: Input Parameters:
584: + sp - The `PetscDualSpace`
585: - order - The order
587: Level: intermediate
589: .seealso: `PetscDualSpace`, `PetscDualSpaceGetOrder()`, `PetscDualSpaceCreate()`
590: @*/
591: PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
592: {
593: PetscFunctionBegin;
595: PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change order after dualspace is set up");
596: sp->order = order;
597: PetscFunctionReturn(PETSC_SUCCESS);
598: }
600: /*@C
601: PetscDualSpaceGetNumComponents - Return the number of components for this space
603: Input Parameter:
604: . sp - The `PetscDualSpace`
606: Output Parameter:
607: . Nc - The number of components
609: Level: intermediate
611: Note:
612: A vector space, for example, will have d components, where d is the spatial dimension
614: .seealso: `PetscDualSpaceSetNumComponents()`, `PetscDualSpaceGetDimension()`, `PetscDualSpaceCreate()`, `PetscDualSpace`
615: @*/
616: PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
617: {
618: PetscFunctionBegin;
620: PetscAssertPointer(Nc, 2);
621: *Nc = sp->Nc;
622: PetscFunctionReturn(PETSC_SUCCESS);
623: }
625: /*@C
626: PetscDualSpaceSetNumComponents - Set the number of components for this space
628: Input Parameters:
629: + sp - The `PetscDualSpace`
630: - Nc - The number of components
632: Level: intermediate
634: .seealso: `PetscDualSpaceGetNumComponents()`, `PetscDualSpaceCreate()`, `PetscDualSpace`
635: @*/
636: PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
637: {
638: PetscFunctionBegin;
640: PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up");
641: sp->Nc = Nc;
642: PetscFunctionReturn(PETSC_SUCCESS);
643: }
645: /*@C
646: PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
648: Not Collective
650: Input Parameters:
651: + sp - The `PetscDualSpace`
652: - i - The basis number
654: Output Parameter:
655: . functional - The basis functional
657: Level: intermediate
659: .seealso: `PetscDualSpace`, `PetscQuadrature`, `PetscDualSpaceGetDimension()`, `PetscDualSpaceCreate()`
660: @*/
661: PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
662: {
663: PetscInt dim;
665: PetscFunctionBegin;
667: PetscAssertPointer(functional, 3);
668: PetscCall(PetscDualSpaceGetDimension(sp, &dim));
669: PetscCheck(!(i < 0) && !(i >= dim), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", i, dim);
670: *functional = sp->functional[i];
671: PetscFunctionReturn(PETSC_SUCCESS);
672: }
674: /*@C
675: PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
677: Not Collective
679: Input Parameter:
680: . sp - The `PetscDualSpace`
682: Output Parameter:
683: . dim - The dimension
685: Level: intermediate
687: .seealso: `PetscDualSpace`, `PetscDualSpaceGetFunctional()`, `PetscDualSpaceCreate()`
688: @*/
689: PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
690: {
691: PetscFunctionBegin;
693: PetscAssertPointer(dim, 2);
694: if (sp->spdim < 0) {
695: PetscSection section;
697: PetscCall(PetscDualSpaceGetSection(sp, §ion));
698: if (section) {
699: PetscCall(PetscSectionGetStorageSize(section, &sp->spdim));
700: } else sp->spdim = 0;
701: }
702: *dim = sp->spdim;
703: PetscFunctionReturn(PETSC_SUCCESS);
704: }
706: /*@C
707: PetscDualSpaceGetInteriorDimension - Get the interior dimension of the dual space, i.e. the number of basis functionals assigned to the interior of the reference domain
709: Not Collective
711: Input Parameter:
712: . sp - The `PetscDualSpace`
714: Output Parameter:
715: . intdim - The dimension
717: Level: intermediate
719: .seealso: `PetscDualSpace`, `PetscDualSpaceGetFunctional()`, `PetscDualSpaceCreate()`
720: @*/
721: PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace sp, PetscInt *intdim)
722: {
723: PetscFunctionBegin;
725: PetscAssertPointer(intdim, 2);
726: if (sp->spintdim < 0) {
727: PetscSection section;
729: PetscCall(PetscDualSpaceGetSection(sp, §ion));
730: if (section) {
731: PetscCall(PetscSectionGetConstrainedStorageSize(section, &sp->spintdim));
732: } else sp->spintdim = 0;
733: }
734: *intdim = sp->spintdim;
735: PetscFunctionReturn(PETSC_SUCCESS);
736: }
738: /*@C
739: PetscDualSpaceGetUniform - Whether this dual space is uniform
741: Not Collective
743: Input Parameter:
744: . sp - A dual space
746: Output Parameter:
747: . uniform - `PETSC_TRUE` if (a) the dual space is the same for each point in a stratum of the reference `DMPLEX`, and
748: (b) every symmetry of each point in the reference `DMPLEX` is also a symmetry of the point's dual space.
750: Level: advanced
752: Note:
753: All of the usual spaces on simplex or tensor-product elements will be uniform, only reference cells
754: with non-uniform strata (like trianguar-prisms) or anisotropic hp dual spaces will not be uniform.
756: .seealso: `PetscDualSpace`, `PetscDualSpaceGetPointSubspace()`, `PetscDualSpaceGetSymmetries()`
757: @*/
758: PetscErrorCode PetscDualSpaceGetUniform(PetscDualSpace sp, PetscBool *uniform)
759: {
760: PetscFunctionBegin;
762: PetscAssertPointer(uniform, 2);
763: *uniform = sp->uniform;
764: PetscFunctionReturn(PETSC_SUCCESS);
765: }
767: /*@CC
768: PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
770: Not Collective
772: Input Parameter:
773: . sp - The `PetscDualSpace`
775: Output Parameter:
776: . numDof - An array of length dim+1 which holds the number of dofs for each dimension
778: Level: intermediate
780: Note:
781: Do not free `numDof`
783: .seealso: `PetscDualSpace`, `PetscDualSpaceGetFunctional()`, `PetscDualSpaceCreate()`
784: @*/
785: PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt *numDof[])
786: {
787: PetscFunctionBegin;
789: PetscAssertPointer(numDof, 2);
790: 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");
791: if (!sp->numDof) {
792: DM dm;
793: PetscInt depth, d;
794: PetscSection section;
796: PetscCall(PetscDualSpaceGetDM(sp, &dm));
797: PetscCall(DMPlexGetDepth(dm, &depth));
798: PetscCall(PetscCalloc1(depth + 1, &sp->numDof));
799: PetscCall(PetscDualSpaceGetSection(sp, §ion));
800: for (d = 0; d <= depth; d++) {
801: PetscInt dStart, dEnd;
803: PetscCall(DMPlexGetDepthStratum(dm, d, &dStart, &dEnd));
804: if (dEnd <= dStart) continue;
805: PetscCall(PetscSectionGetDof(section, dStart, &sp->numDof[d]));
806: }
807: }
808: *numDof = sp->numDof;
809: PetscCheck(*numDof, PetscObjectComm((PetscObject)sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
810: PetscFunctionReturn(PETSC_SUCCESS);
811: }
813: /* create the section of the right size and set a permutation for topological ordering */
814: PetscErrorCode PetscDualSpaceSectionCreate_Internal(PetscDualSpace sp, PetscSection *topSection)
815: {
816: DM dm;
817: PetscInt pStart, pEnd, cStart, cEnd, c, depth, count, i;
818: PetscInt *seen, *perm;
819: PetscSection section;
821: PetscFunctionBegin;
822: dm = sp->dm;
823: PetscCall(PetscSectionCreate(PETSC_COMM_SELF, §ion));
824: PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
825: PetscCall(PetscSectionSetChart(section, pStart, pEnd));
826: PetscCall(PetscCalloc1(pEnd - pStart, &seen));
827: PetscCall(PetscMalloc1(pEnd - pStart, &perm));
828: PetscCall(DMPlexGetDepth(dm, &depth));
829: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
830: for (c = cStart, count = 0; c < cEnd; c++) {
831: PetscInt closureSize = -1, e;
832: PetscInt *closure = NULL;
834: perm[count++] = c;
835: seen[c - pStart] = 1;
836: PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
837: for (e = 0; e < closureSize; e++) {
838: PetscInt point = closure[2 * e];
840: if (seen[point - pStart]) continue;
841: perm[count++] = point;
842: seen[point - pStart] = 1;
843: }
844: PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
845: }
846: PetscCheck(count == pEnd - pStart, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad topological ordering");
847: for (i = 0; i < pEnd - pStart; i++)
848: if (perm[i] != i) break;
849: if (i < pEnd - pStart) {
850: IS permIS;
852: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, perm, PETSC_OWN_POINTER, &permIS));
853: PetscCall(ISSetPermutation(permIS));
854: PetscCall(PetscSectionSetPermutation(section, permIS));
855: PetscCall(ISDestroy(&permIS));
856: } else {
857: PetscCall(PetscFree(perm));
858: }
859: PetscCall(PetscFree(seen));
860: *topSection = section;
861: PetscFunctionReturn(PETSC_SUCCESS);
862: }
864: /* mark boundary points and set up */
865: PetscErrorCode PetscDualSpaceSectionSetUp_Internal(PetscDualSpace sp, PetscSection section)
866: {
867: DM dm;
868: DMLabel boundary;
869: PetscInt pStart, pEnd, p;
871: PetscFunctionBegin;
872: dm = sp->dm;
873: PetscCall(DMLabelCreate(PETSC_COMM_SELF, "boundary", &boundary));
874: PetscCall(PetscDualSpaceGetDM(sp, &dm));
875: PetscCall(DMPlexMarkBoundaryFaces(dm, 1, boundary));
876: PetscCall(DMPlexLabelComplete(dm, boundary));
877: PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
878: for (p = pStart; p < pEnd; p++) {
879: PetscInt bval;
881: PetscCall(DMLabelGetValue(boundary, p, &bval));
882: if (bval == 1) {
883: PetscInt dof;
885: PetscCall(PetscSectionGetDof(section, p, &dof));
886: PetscCall(PetscSectionSetConstraintDof(section, p, dof));
887: }
888: }
889: PetscCall(DMLabelDestroy(&boundary));
890: PetscCall(PetscSectionSetUp(section));
891: PetscFunctionReturn(PETSC_SUCCESS);
892: }
894: /*@C
895: PetscDualSpaceGetSection - Create a `PetscSection` over the reference cell with the layout from this space
897: Collective
899: Input Parameter:
900: . sp - The `PetscDualSpace`
902: Output Parameter:
903: . section - The section
905: Level: advanced
907: .seealso: `PetscDualSpace`, `PetscSection`, `PetscDualSpaceCreate()`, `DMPLEX`
908: @*/
909: PetscErrorCode PetscDualSpaceGetSection(PetscDualSpace sp, PetscSection *section)
910: {
911: PetscInt pStart, pEnd, p;
913: PetscFunctionBegin;
914: if (!sp->dm) {
915: *section = NULL;
916: PetscFunctionReturn(PETSC_SUCCESS);
917: }
918: if (!sp->pointSection) {
919: /* mark the boundary */
920: PetscCall(PetscDualSpaceSectionCreate_Internal(sp, &sp->pointSection));
921: PetscCall(DMPlexGetChart(sp->dm, &pStart, &pEnd));
922: for (p = pStart; p < pEnd; p++) {
923: PetscDualSpace psp;
925: PetscCall(PetscDualSpaceGetPointSubspace(sp, p, &psp));
926: if (psp) {
927: PetscInt dof;
929: PetscCall(PetscDualSpaceGetInteriorDimension(psp, &dof));
930: PetscCall(PetscSectionSetDof(sp->pointSection, p, dof));
931: }
932: }
933: PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, sp->pointSection));
934: }
935: *section = sp->pointSection;
936: PetscFunctionReturn(PETSC_SUCCESS);
937: }
939: /*@C
940: PetscDualSpaceGetInteriorSection - Create a `PetscSection` over the reference cell with the layout from this space
941: for interior degrees of freedom
943: Collective
945: Input Parameter:
946: . sp - The `PetscDualSpace`
948: Output Parameter:
949: . section - The interior section
951: Level: advanced
953: Note:
954: Most reference domains have one cell, in which case the only cell will have
955: all of the interior degrees of freedom in the interior section. But
956: for `PETSCDUALSPACEREFINED` there may be other mesh points in the interior,
957: and this section describes their layout.
959: .seealso: `PetscDualSpace`, `PetscSection`, `PetscDualSpaceCreate()`, `DMPLEX`
960: @*/
961: PetscErrorCode PetscDualSpaceGetInteriorSection(PetscDualSpace sp, PetscSection *section)
962: {
963: PetscInt pStart, pEnd, p;
965: PetscFunctionBegin;
966: if (!sp->dm) {
967: *section = NULL;
968: PetscFunctionReturn(PETSC_SUCCESS);
969: }
970: if (!sp->intPointSection) {
971: PetscSection full_section;
973: PetscCall(PetscDualSpaceGetSection(sp, &full_section));
974: PetscCall(PetscDualSpaceSectionCreate_Internal(sp, &sp->intPointSection));
975: PetscCall(PetscSectionGetChart(full_section, &pStart, &pEnd));
976: for (p = pStart; p < pEnd; p++) {
977: PetscInt dof, cdof;
979: PetscCall(PetscSectionGetDof(full_section, p, &dof));
980: PetscCall(PetscSectionGetConstraintDof(full_section, p, &cdof));
981: PetscCall(PetscSectionSetDof(sp->intPointSection, p, dof - cdof));
982: }
983: PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, sp->intPointSection));
984: }
985: *section = sp->intPointSection;
986: PetscFunctionReturn(PETSC_SUCCESS);
987: }
989: /* this assumes that all of the point dual spaces store their interior dofs first, which is true when the point DMs
990: * have one cell */
991: PetscErrorCode PetscDualSpacePushForwardSubspaces_Internal(PetscDualSpace sp, PetscInt sStart, PetscInt sEnd)
992: {
993: PetscReal *sv0, *v0, *J;
994: PetscSection section;
995: PetscInt dim, s, k;
996: DM dm;
998: PetscFunctionBegin;
999: PetscCall(PetscDualSpaceGetDM(sp, &dm));
1000: PetscCall(DMGetDimension(dm, &dim));
1001: PetscCall(PetscDualSpaceGetSection(sp, §ion));
1002: PetscCall(PetscMalloc3(dim, &v0, dim, &sv0, dim * dim, &J));
1003: PetscCall(PetscDualSpaceGetFormDegree(sp, &k));
1004: for (s = sStart; s < sEnd; s++) {
1005: PetscReal detJ, hdetJ;
1006: PetscDualSpace ssp;
1007: PetscInt dof, off, f, sdim;
1008: PetscInt i, j;
1009: DM sdm;
1011: PetscCall(PetscDualSpaceGetPointSubspace(sp, s, &ssp));
1012: if (!ssp) continue;
1013: PetscCall(PetscSectionGetDof(section, s, &dof));
1014: PetscCall(PetscSectionGetOffset(section, s, &off));
1015: /* get the first vertex of the reference cell */
1016: PetscCall(PetscDualSpaceGetDM(ssp, &sdm));
1017: PetscCall(DMGetDimension(sdm, &sdim));
1018: PetscCall(DMPlexComputeCellGeometryAffineFEM(sdm, 0, sv0, NULL, NULL, &hdetJ));
1019: PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, s, v0, J, NULL, &detJ));
1020: /* compactify Jacobian */
1021: for (i = 0; i < dim; i++)
1022: for (j = 0; j < sdim; j++) J[i * sdim + j] = J[i * dim + j];
1023: for (f = 0; f < dof; f++) {
1024: PetscQuadrature fn;
1026: PetscCall(PetscDualSpaceGetFunctional(ssp, f, &fn));
1027: PetscCall(PetscQuadraturePushForward(fn, dim, sv0, v0, J, k, &sp->functional[off + f]));
1028: }
1029: }
1030: PetscCall(PetscFree3(v0, sv0, J));
1031: PetscFunctionReturn(PETSC_SUCCESS);
1032: }
1034: /*@C
1035: PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
1037: Input Parameters:
1038: + sp - The `PetscDualSpace` object
1039: . f - The basis functional index
1040: . time - The time
1041: . 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)
1042: . numComp - The number of components for the function
1043: . func - The input function
1044: - ctx - A context for the function
1046: Output Parameter:
1047: . value - numComp output values
1049: Calling sequence:
1050: .vb
1051: PetscErrorCode func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt numComponents, PetscScalar values[], void *ctx)
1052: .ve
1054: Level: beginner
1056: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1057: @*/
1058: 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)
1059: {
1060: PetscFunctionBegin;
1062: PetscAssertPointer(cgeom, 4);
1063: PetscAssertPointer(value, 8);
1064: PetscUseTypeMethod(sp, apply, f, time, cgeom, numComp, func, ctx, value);
1065: PetscFunctionReturn(PETSC_SUCCESS);
1066: }
1068: /*@C
1069: PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by `PetscDualSpaceGetAllData()`
1071: Input Parameters:
1072: + sp - The `PetscDualSpace` object
1073: - pointEval - Evaluation at the points returned by `PetscDualSpaceGetAllData()`
1075: Output Parameter:
1076: . spValue - The values of all dual space functionals
1078: Level: advanced
1080: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1081: @*/
1082: PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1083: {
1084: PetscFunctionBegin;
1086: PetscUseTypeMethod(sp, applyall, pointEval, spValue);
1087: PetscFunctionReturn(PETSC_SUCCESS);
1088: }
1090: /*@C
1091: PetscDualSpaceApplyInterior - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by `PetscDualSpaceGetInteriorData()`
1093: Input Parameters:
1094: + sp - The `PetscDualSpace` object
1095: - pointEval - Evaluation at the points returned by `PetscDualSpaceGetInteriorData()`
1097: Output Parameter:
1098: . spValue - The values of interior dual space functionals
1100: Level: advanced
1102: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1103: @*/
1104: PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1105: {
1106: PetscFunctionBegin;
1108: PetscUseTypeMethod(sp, applyint, pointEval, spValue);
1109: PetscFunctionReturn(PETSC_SUCCESS);
1110: }
1112: /*@C
1113: PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
1115: Input Parameters:
1116: + sp - The `PetscDualSpace` object
1117: . f - The basis functional index
1118: . time - The time
1119: . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
1120: . Nc - The number of components for the function
1121: . func - The input function
1122: - ctx - A context for the function
1124: Output Parameter:
1125: . value - The output value
1127: Calling sequence:
1128: .vb
1129: PetscErrorCode func(PetscInt dim, PetscReal time, const PetscReal x[],PetscInt numComponents, PetscScalar values[], void *ctx)
1130: .ve
1132: Level: advanced
1134: Note:
1135: 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.
1137: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1138: @*/
1139: 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)
1140: {
1141: DM dm;
1142: PetscQuadrature n;
1143: const PetscReal *points, *weights;
1144: PetscReal x[3];
1145: PetscScalar *val;
1146: PetscInt dim, dE, qNc, c, Nq, q;
1147: PetscBool isAffine;
1149: PetscFunctionBegin;
1151: PetscAssertPointer(value, 8);
1152: PetscCall(PetscDualSpaceGetDM(sp, &dm));
1153: PetscCall(PetscDualSpaceGetFunctional(sp, f, &n));
1154: PetscCall(PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights));
1155: 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);
1156: PetscCheck(qNc == Nc, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_SIZ, "The quadrature components %" PetscInt_FMT " != function components %" PetscInt_FMT, qNc, Nc);
1157: PetscCall(DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val));
1158: *value = 0.0;
1159: isAffine = cgeom->isAffine;
1160: dE = cgeom->dimEmbed;
1161: for (q = 0; q < Nq; ++q) {
1162: if (isAffine) {
1163: CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q * dim], x);
1164: PetscCall((*func)(dE, time, x, Nc, val, ctx));
1165: } else {
1166: PetscCall((*func)(dE, time, &cgeom->v[dE * q], Nc, val, ctx));
1167: }
1168: for (c = 0; c < Nc; ++c) *value += val[c] * weights[q * Nc + c];
1169: }
1170: PetscCall(DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val));
1171: PetscFunctionReturn(PETSC_SUCCESS);
1172: }
1174: /*@C
1175: PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by `PetscDualSpaceGetAllData()`
1177: Input Parameters:
1178: + sp - The `PetscDualSpace` object
1179: - pointEval - Evaluation at the points returned by `PetscDualSpaceGetAllData()`
1181: Output Parameter:
1182: . spValue - The values of all dual space functionals
1184: Level: advanced
1186: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1187: @*/
1188: PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1189: {
1190: Vec pointValues, dofValues;
1191: Mat allMat;
1193: PetscFunctionBegin;
1195: PetscAssertPointer(pointEval, 2);
1196: PetscAssertPointer(spValue, 3);
1197: PetscCall(PetscDualSpaceGetAllData(sp, NULL, &allMat));
1198: if (!sp->allNodeValues) PetscCall(MatCreateVecs(allMat, &sp->allNodeValues, NULL));
1199: pointValues = sp->allNodeValues;
1200: if (!sp->allDofValues) PetscCall(MatCreateVecs(allMat, NULL, &sp->allDofValues));
1201: dofValues = sp->allDofValues;
1202: PetscCall(VecPlaceArray(pointValues, pointEval));
1203: PetscCall(VecPlaceArray(dofValues, spValue));
1204: PetscCall(MatMult(allMat, pointValues, dofValues));
1205: PetscCall(VecResetArray(dofValues));
1206: PetscCall(VecResetArray(pointValues));
1207: PetscFunctionReturn(PETSC_SUCCESS);
1208: }
1210: /*@C
1211: PetscDualSpaceApplyInteriorDefault - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by `PetscDualSpaceGetInteriorData()`
1213: Input Parameters:
1214: + sp - The `PetscDualSpace` object
1215: - pointEval - Evaluation at the points returned by `PetscDualSpaceGetInteriorData()`
1217: Output Parameter:
1218: . spValue - The values of interior dual space functionals
1220: Level: advanced
1222: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1223: @*/
1224: PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1225: {
1226: Vec pointValues, dofValues;
1227: Mat intMat;
1229: PetscFunctionBegin;
1231: PetscAssertPointer(pointEval, 2);
1232: PetscAssertPointer(spValue, 3);
1233: PetscCall(PetscDualSpaceGetInteriorData(sp, NULL, &intMat));
1234: if (!sp->intNodeValues) PetscCall(MatCreateVecs(intMat, &sp->intNodeValues, NULL));
1235: pointValues = sp->intNodeValues;
1236: if (!sp->intDofValues) PetscCall(MatCreateVecs(intMat, NULL, &sp->intDofValues));
1237: dofValues = sp->intDofValues;
1238: PetscCall(VecPlaceArray(pointValues, pointEval));
1239: PetscCall(VecPlaceArray(dofValues, spValue));
1240: PetscCall(MatMult(intMat, pointValues, dofValues));
1241: PetscCall(VecResetArray(dofValues));
1242: PetscCall(VecResetArray(pointValues));
1243: PetscFunctionReturn(PETSC_SUCCESS);
1244: }
1246: /*@C
1247: PetscDualSpaceGetAllData - Get all quadrature nodes from this space, and the matrix that sends quadrature node values to degree-of-freedom values
1249: Input Parameter:
1250: . sp - The dualspace
1252: Output Parameters:
1253: + allNodes - A `PetscQuadrature` object containing all evaluation nodes, pass `NULL` if not needed
1254: - allMat - A `Mat` for the node-to-dof transformation, pass `NULL` if not needed
1256: Level: advanced
1258: .seealso: `PetscQuadrature`, `PetscDualSpace`, `PetscDualSpaceCreate()`, `Mat`
1259: @*/
1260: PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace sp, PeOp PetscQuadrature *allNodes, PeOp Mat *allMat)
1261: {
1262: PetscFunctionBegin;
1264: if (allNodes) PetscAssertPointer(allNodes, 2);
1265: if (allMat) PetscAssertPointer(allMat, 3);
1266: if ((!sp->allNodes || !sp->allMat) && sp->ops->createalldata) {
1267: PetscQuadrature qpoints;
1268: Mat amat;
1270: PetscUseTypeMethod(sp, createalldata, &qpoints, &amat);
1271: PetscCall(PetscQuadratureDestroy(&sp->allNodes));
1272: PetscCall(MatDestroy(&sp->allMat));
1273: sp->allNodes = qpoints;
1274: sp->allMat = amat;
1275: }
1276: if (allNodes) *allNodes = sp->allNodes;
1277: if (allMat) *allMat = sp->allMat;
1278: PetscFunctionReturn(PETSC_SUCCESS);
1279: }
1281: /*@C
1282: PetscDualSpaceCreateAllDataDefault - Create all evaluation nodes and the node-to-dof matrix by examining functionals
1284: Input Parameter:
1285: . sp - The dualspace
1287: Output Parameters:
1288: + allNodes - A `PetscQuadrature` object containing all evaluation nodes
1289: - allMat - A `Mat` for the node-to-dof transformation
1291: Level: advanced
1293: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`, `Mat`, `PetscQuadrature`
1294: @*/
1295: PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat)
1296: {
1297: PetscInt spdim;
1298: PetscInt numPoints, offset;
1299: PetscReal *points;
1300: PetscInt f, dim;
1301: PetscInt Nc, nrows, ncols;
1302: PetscInt maxNumPoints;
1303: PetscQuadrature q;
1304: Mat A;
1306: PetscFunctionBegin;
1307: PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
1308: PetscCall(PetscDualSpaceGetDimension(sp, &spdim));
1309: if (!spdim) {
1310: PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allNodes));
1311: PetscCall(PetscQuadratureSetData(*allNodes, 0, 0, 0, NULL, NULL));
1312: }
1313: nrows = spdim;
1314: PetscCall(PetscDualSpaceGetFunctional(sp, 0, &q));
1315: PetscCall(PetscQuadratureGetData(q, &dim, NULL, &numPoints, NULL, NULL));
1316: maxNumPoints = numPoints;
1317: for (f = 1; f < spdim; f++) {
1318: PetscInt Np;
1320: PetscCall(PetscDualSpaceGetFunctional(sp, f, &q));
1321: PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Np, NULL, NULL));
1322: numPoints += Np;
1323: maxNumPoints = PetscMax(maxNumPoints, Np);
1324: }
1325: ncols = numPoints * Nc;
1326: PetscCall(PetscMalloc1(dim * numPoints, &points));
1327: PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nrows, ncols, maxNumPoints * Nc, NULL, &A));
1328: for (f = 0, offset = 0; f < spdim; f++) {
1329: const PetscReal *p, *w;
1330: PetscInt Np, i;
1331: PetscInt fnc;
1333: PetscCall(PetscDualSpaceGetFunctional(sp, f, &q));
1334: PetscCall(PetscQuadratureGetData(q, NULL, &fnc, &Np, &p, &w));
1335: PetscCheck(fnc == Nc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "functional component mismatch");
1336: for (i = 0; i < Np * dim; i++) points[offset * dim + i] = p[i];
1337: for (i = 0; i < Np * Nc; i++) PetscCall(MatSetValue(A, f, offset * Nc, w[i], INSERT_VALUES));
1338: offset += Np;
1339: }
1340: PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1341: PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1342: PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allNodes));
1343: PetscCall(PetscQuadratureSetData(*allNodes, dim, 0, numPoints, points, NULL));
1344: *allMat = A;
1345: PetscFunctionReturn(PETSC_SUCCESS);
1346: }
1348: /*@C
1349: PetscDualSpaceGetInteriorData - Get all quadrature points necessary to compute the interior degrees of freedom from
1350: this space, as well as the matrix that computes the degrees of freedom from the quadrature
1351: values.
1353: Input Parameter:
1354: . sp - The dualspace
1356: Output Parameters:
1357: + intNodes - A `PetscQuadrature` object containing all evaluation points needed to evaluate interior degrees of freedom,
1358: pass `NULL` if not needed
1359: - intMat - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is
1360: the size of the constrained layout (`PetscSectionGetConstrainStorageSize()`) of the dual space section,
1361: npoints is the number of points in intNodes and nc is `PetscDualSpaceGetNumComponents()`.
1362: Pass `NULL` if not needed
1364: Level: advanced
1366: Notes:
1367: Degrees of freedom are interior degrees of freedom if they belong (by
1368: `PetscDualSpaceGetSection()`) to interior points in the references, complementary boundary
1369: degrees of freedom are marked as constrained in the section returned by
1370: `PetscDualSpaceGetSection()`).
1372: .seealso: `PetscDualSpace`, `PetscQuadrature`, `Mat`, `PetscDualSpaceCreate()`, `PetscDualSpaceGetDimension()`, `PetscDualSpaceGetNumComponents()`, `PetscQuadratureGetData()`
1373: @*/
1374: PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace sp, PeOp PetscQuadrature *intNodes, PeOp Mat *intMat)
1375: {
1376: PetscFunctionBegin;
1378: if (intNodes) PetscAssertPointer(intNodes, 2);
1379: if (intMat) PetscAssertPointer(intMat, 3);
1380: if ((!sp->intNodes || !sp->intMat) && sp->ops->createintdata) {
1381: PetscQuadrature qpoints;
1382: Mat imat;
1384: PetscUseTypeMethod(sp, createintdata, &qpoints, &imat);
1385: PetscCall(PetscQuadratureDestroy(&sp->intNodes));
1386: PetscCall(MatDestroy(&sp->intMat));
1387: sp->intNodes = qpoints;
1388: sp->intMat = imat;
1389: }
1390: if (intNodes) *intNodes = sp->intNodes;
1391: if (intMat) *intMat = sp->intMat;
1392: PetscFunctionReturn(PETSC_SUCCESS);
1393: }
1395: /*@C
1396: PetscDualSpaceCreateInteriorDataDefault - Create quadrature points by examining interior functionals and create the matrix mapping quadrature point values to interior dual space values
1398: Input Parameter:
1399: . sp - The dualspace
1401: Output Parameters:
1402: + intNodes - A `PetscQuadrature` object containing all evaluation points needed to evaluate interior degrees of freedom
1403: - intMat - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is
1404: the size of the constrained layout (`PetscSectionGetConstrainStorageSize()`) of the dual space section,
1405: npoints is the number of points in allNodes and nc is `PetscDualSpaceGetNumComponents()`.
1407: Level: advanced
1409: .seealso: `PetscDualSpace`, `PetscQuadrature`, `Mat`, `PetscDualSpaceCreate()`, `PetscDualSpaceGetInteriorData()`
1410: @*/
1411: PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat)
1412: {
1413: DM dm;
1414: PetscInt spdim0;
1415: PetscInt Nc;
1416: PetscInt pStart, pEnd, p, f;
1417: PetscSection section;
1418: PetscInt numPoints, offset, matoffset;
1419: PetscReal *points;
1420: PetscInt dim;
1421: PetscInt *nnz;
1422: PetscQuadrature q;
1423: Mat imat;
1425: PetscFunctionBegin;
1427: PetscCall(PetscDualSpaceGetSection(sp, §ion));
1428: PetscCall(PetscSectionGetConstrainedStorageSize(section, &spdim0));
1429: if (!spdim0) {
1430: *intNodes = NULL;
1431: *intMat = NULL;
1432: PetscFunctionReturn(PETSC_SUCCESS);
1433: }
1434: PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
1435: PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1436: PetscCall(PetscDualSpaceGetDM(sp, &dm));
1437: PetscCall(DMGetDimension(dm, &dim));
1438: PetscCall(PetscMalloc1(spdim0, &nnz));
1439: for (p = pStart, f = 0, numPoints = 0; p < pEnd; p++) {
1440: PetscInt dof, cdof, off, d;
1442: PetscCall(PetscSectionGetDof(section, p, &dof));
1443: PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
1444: if (!(dof - cdof)) continue;
1445: PetscCall(PetscSectionGetOffset(section, p, &off));
1446: for (d = 0; d < dof; d++, off++, f++) {
1447: PetscInt Np;
1449: PetscCall(PetscDualSpaceGetFunctional(sp, off, &q));
1450: PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Np, NULL, NULL));
1451: nnz[f] = Np * Nc;
1452: numPoints += Np;
1453: }
1454: }
1455: PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, spdim0, numPoints * Nc, 0, nnz, &imat));
1456: PetscCall(PetscFree(nnz));
1457: PetscCall(PetscMalloc1(dim * numPoints, &points));
1458: for (p = pStart, f = 0, offset = 0, matoffset = 0; p < pEnd; p++) {
1459: PetscInt dof, cdof, off, d;
1461: PetscCall(PetscSectionGetDof(section, p, &dof));
1462: PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
1463: if (!(dof - cdof)) continue;
1464: PetscCall(PetscSectionGetOffset(section, p, &off));
1465: for (d = 0; d < dof; d++, off++, f++) {
1466: const PetscReal *p;
1467: const PetscReal *w;
1468: PetscInt Np, i;
1470: PetscCall(PetscDualSpaceGetFunctional(sp, off, &q));
1471: PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Np, &p, &w));
1472: for (i = 0; i < Np * dim; i++) points[offset + i] = p[i];
1473: for (i = 0; i < Np * Nc; i++) PetscCall(MatSetValue(imat, f, matoffset + i, w[i], INSERT_VALUES));
1474: offset += Np * dim;
1475: matoffset += Np * Nc;
1476: }
1477: }
1478: PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, intNodes));
1479: PetscCall(PetscQuadratureSetData(*intNodes, dim, 0, numPoints, points, NULL));
1480: PetscCall(MatAssemblyBegin(imat, MAT_FINAL_ASSEMBLY));
1481: PetscCall(MatAssemblyEnd(imat, MAT_FINAL_ASSEMBLY));
1482: *intMat = imat;
1483: PetscFunctionReturn(PETSC_SUCCESS);
1484: }
1486: /*@C
1487: PetscDualSpaceEqual - Determine if two dual spaces are equivalent
1489: Input Parameters:
1490: + A - A `PetscDualSpace` object
1491: - B - Another `PetscDualSpace` object
1493: Output Parameter:
1494: . equal - `PETSC_TRUE` if the dual spaces are equivalent
1496: Level: advanced
1498: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1499: @*/
1500: PetscErrorCode PetscDualSpaceEqual(PetscDualSpace A, PetscDualSpace B, PetscBool *equal)
1501: {
1502: PetscInt sizeA, sizeB, dimA, dimB;
1503: const PetscInt *dofA, *dofB;
1504: PetscQuadrature quadA, quadB;
1505: Mat matA, matB;
1507: PetscFunctionBegin;
1510: PetscAssertPointer(equal, 3);
1511: *equal = PETSC_FALSE;
1512: PetscCall(PetscDualSpaceGetDimension(A, &sizeA));
1513: PetscCall(PetscDualSpaceGetDimension(B, &sizeB));
1514: if (sizeB != sizeA) PetscFunctionReturn(PETSC_SUCCESS);
1515: PetscCall(DMGetDimension(A->dm, &dimA));
1516: PetscCall(DMGetDimension(B->dm, &dimB));
1517: if (dimA != dimB) PetscFunctionReturn(PETSC_SUCCESS);
1519: PetscCall(PetscDualSpaceGetNumDof(A, &dofA));
1520: PetscCall(PetscDualSpaceGetNumDof(B, &dofB));
1521: for (PetscInt d = 0; d < dimA; d++) {
1522: if (dofA[d] != dofB[d]) PetscFunctionReturn(PETSC_SUCCESS);
1523: }
1525: PetscCall(PetscDualSpaceGetInteriorData(A, &quadA, &matA));
1526: PetscCall(PetscDualSpaceGetInteriorData(B, &quadB, &matB));
1527: if (!quadA && !quadB) {
1528: *equal = PETSC_TRUE;
1529: } else if (quadA && quadB) {
1530: PetscCall(PetscQuadratureEqual(quadA, quadB, equal));
1531: if (*equal == PETSC_FALSE) PetscFunctionReturn(PETSC_SUCCESS);
1532: if (!matA && !matB) PetscFunctionReturn(PETSC_SUCCESS);
1533: if (matA && matB) PetscCall(MatEqual(matA, matB, equal));
1534: else *equal = PETSC_FALSE;
1535: }
1536: PetscFunctionReturn(PETSC_SUCCESS);
1537: }
1539: /*@C
1540: PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
1542: Input Parameters:
1543: + sp - The `PetscDualSpace` object
1544: . f - The basis functional index
1545: . time - The time
1546: . cgeom - A context with geometric information for this cell, we currently just use the centroid
1547: . Nc - The number of components for the function
1548: . func - The input function
1549: - ctx - A context for the function
1551: Output Parameter:
1552: . value - The output value (scalar)
1554: Calling sequence:
1555: .vb
1556: PetscErrorCode func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt numComponents, PetscScalar values[], void *ctx)
1557: .ve
1559: Level: advanced
1561: Note:
1562: 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.
1564: .seealso: `PetscDualSpace`, `PetscDualSpaceCreate()`
1565: @*/
1566: 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)
1567: {
1568: DM dm;
1569: PetscQuadrature n;
1570: const PetscReal *points, *weights;
1571: PetscScalar *val;
1572: PetscInt dimEmbed, qNc, c, Nq, q;
1574: PetscFunctionBegin;
1576: PetscAssertPointer(value, 8);
1577: PetscCall(PetscDualSpaceGetDM(sp, &dm));
1578: PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
1579: PetscCall(PetscDualSpaceGetFunctional(sp, f, &n));
1580: PetscCall(PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights));
1581: PetscCheck(qNc == Nc, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_SIZ, "The quadrature components %" PetscInt_FMT " != function components %" PetscInt_FMT, qNc, Nc);
1582: PetscCall(DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val));
1583: *value = 0.;
1584: for (q = 0; q < Nq; ++q) {
1585: PetscCall((*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx));
1586: for (c = 0; c < Nc; ++c) *value += val[c] * weights[q * Nc + c];
1587: }
1588: PetscCall(DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val));
1589: PetscFunctionReturn(PETSC_SUCCESS);
1590: }
1592: /*@C
1593: PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
1594: given height. This assumes that the reference cell is symmetric over points of this height.
1596: Not Collective
1598: Input Parameters:
1599: + sp - the `PetscDualSpace` object
1600: - height - the height of the mesh point for which the subspace is desired
1602: Output Parameter:
1603: . subsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the
1604: point, which will be of lesser dimension if height > 0.
1606: Level: advanced
1608: Notes:
1609: If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
1610: pointwise values are not defined on the element boundaries), or if the implementation of `PetscDualSpace` does not
1611: support extracting subspaces, then `NULL` is returned.
1613: This does not increment the reference count on the returned dual space, and the user should not destroy it.
1615: .seealso: `PetscDualSpace`, `PetscSpaceGetHeightSubspace()`, `PetscDualSpaceGetPointSubspace()`
1616: @*/
1617: PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
1618: {
1619: PetscInt depth = -1, cStart, cEnd;
1620: DM dm;
1622: PetscFunctionBegin;
1624: PetscAssertPointer(subsp, 3);
1625: 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");
1626: *subsp = NULL;
1627: dm = sp->dm;
1628: PetscCall(DMPlexGetDepth(dm, &depth));
1629: PetscCheck(height >= 0 && height <= depth, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid height");
1630: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1631: if (height == 0 && cEnd == cStart + 1) {
1632: *subsp = sp;
1633: PetscFunctionReturn(PETSC_SUCCESS);
1634: }
1635: if (!sp->heightSpaces) {
1636: PetscInt h;
1637: PetscCall(PetscCalloc1(depth + 1, &sp->heightSpaces));
1639: for (h = 0; h <= depth; h++) {
1640: if (h == 0 && cEnd == cStart + 1) continue;
1641: if (sp->ops->createheightsubspace) PetscUseTypeMethod(sp, createheightsubspace, height, &sp->heightSpaces[h]);
1642: else if (sp->pointSpaces) {
1643: PetscInt hStart, hEnd;
1645: PetscCall(DMPlexGetHeightStratum(dm, h, &hStart, &hEnd));
1646: if (hEnd > hStart) {
1647: const char *name;
1649: PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[hStart]));
1650: if (sp->pointSpaces[hStart]) {
1651: PetscCall(PetscObjectGetName((PetscObject)sp, &name));
1652: PetscCall(PetscObjectSetName((PetscObject)sp->pointSpaces[hStart], name));
1653: }
1654: sp->heightSpaces[h] = sp->pointSpaces[hStart];
1655: }
1656: }
1657: }
1658: }
1659: *subsp = sp->heightSpaces[height];
1660: PetscFunctionReturn(PETSC_SUCCESS);
1661: }
1663: /*@C
1664: PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.
1666: Not Collective
1668: Input Parameters:
1669: + sp - the `PetscDualSpace` object
1670: - point - the point (in the dual space's DM) for which the subspace is desired
1672: Output Parameter:
1673: . bdsp - the subspace.
1675: Level: advanced
1677: Notes:
1678: The functionals in the subspace are with respect to the intrinsic geometry of the point,
1679: which will be of lesser dimension if height > 0.
1681: If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
1682: defined on the element boundaries), or if the implementation of `PetscDualSpace` does not support extracting
1683: subspaces, then `NULL` is returned.
1685: This does not increment the reference count on the returned dual space, and the user should not destroy it.
1687: .seealso: `PetscDualSpace`, `PetscDualSpaceGetHeightSubspace()`
1688: @*/
1689: PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
1690: {
1691: PetscInt pStart = 0, pEnd = 0, cStart, cEnd;
1692: DM dm;
1694: PetscFunctionBegin;
1696: PetscAssertPointer(bdsp, 3);
1697: *bdsp = NULL;
1698: dm = sp->dm;
1699: PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1700: PetscCheck(point >= pStart && point <= pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point");
1701: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1702: 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 */
1703: *bdsp = sp;
1704: PetscFunctionReturn(PETSC_SUCCESS);
1705: }
1706: if (!sp->pointSpaces) {
1707: PetscInt p;
1708: PetscCall(PetscCalloc1(pEnd - pStart, &sp->pointSpaces));
1710: for (p = 0; p < pEnd - pStart; p++) {
1711: if (p + pStart == cStart && cEnd == cStart + 1) continue;
1712: if (sp->ops->createpointsubspace) PetscUseTypeMethod(sp, createpointsubspace, p + pStart, &sp->pointSpaces[p]);
1713: else if (sp->heightSpaces || sp->ops->createheightsubspace) {
1714: PetscInt dim, depth, height;
1715: DMLabel label;
1717: PetscCall(DMPlexGetDepth(dm, &dim));
1718: PetscCall(DMPlexGetDepthLabel(dm, &label));
1719: PetscCall(DMLabelGetValue(label, p + pStart, &depth));
1720: height = dim - depth;
1721: PetscCall(PetscDualSpaceGetHeightSubspace(sp, height, &sp->pointSpaces[p]));
1722: PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[p]));
1723: }
1724: }
1725: }
1726: *bdsp = sp->pointSpaces[point - pStart];
1727: PetscFunctionReturn(PETSC_SUCCESS);
1728: }
1730: /*@C
1731: PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis
1733: Not Collective
1735: Input Parameter:
1736: . sp - the `PetscDualSpace` object
1738: Output Parameters:
1739: + perms - Permutations of the interior degrees of freedom, parameterized by the point orientation
1740: - flips - Sign reversal of the interior degrees of freedom, parameterized by the point orientation
1742: Level: developer
1744: Note:
1745: The permutation and flip arrays are organized in the following way
1746: .vb
1747: perms[p][ornt][dof # on point] = new local dof #
1748: flips[p][ornt][dof # on point] = reversal or not
1749: .ve
1751: .seealso: `PetscDualSpace`
1752: @*/
1753: PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
1754: {
1755: PetscFunctionBegin;
1757: if (perms) {
1758: PetscAssertPointer(perms, 2);
1759: *perms = NULL;
1760: }
1761: if (flips) {
1762: PetscAssertPointer(flips, 3);
1763: *flips = NULL;
1764: }
1765: PetscTryTypeMethod(sp, getsymmetries, perms, flips);
1766: PetscFunctionReturn(PETSC_SUCCESS);
1767: }
1769: /*@C
1770: PetscDualSpaceGetFormDegree - Get the form degree k for the k-form the describes the pushforwards/pullbacks of this
1771: dual space's functionals.
1773: Input Parameter:
1774: . dsp - The `PetscDualSpace`
1776: Output Parameter:
1777: . k - The *signed* degree k of the k. If k >= 0, this means that the degrees of freedom are k-forms, and are stored
1778: in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms. So for example,
1779: 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).
1780: If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the
1781: Hodge star map. So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms
1782: but are stored as 1-forms.
1784: Level: developer
1786: .seealso: `PetscDualSpace`, `PetscDTAltV`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceTransformType`
1787: @*/
1788: PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace dsp, PetscInt *k)
1789: {
1790: PetscFunctionBeginHot;
1792: PetscAssertPointer(k, 2);
1793: *k = dsp->k;
1794: PetscFunctionReturn(PETSC_SUCCESS);
1795: }
1797: /*@C
1798: PetscDualSpaceSetFormDegree - Set the form degree k for the k-form the describes the pushforwards/pullbacks of this
1799: dual space's functionals.
1801: Input Parameters:
1802: + dsp - The `PetscDualSpace`
1803: - k - The *signed* degree k of the k. If k >= 0, this means that the degrees of freedom are k-forms, and are stored
1804: in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms. So for example,
1805: 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).
1806: If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the
1807: Hodge star map. So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms
1808: but are stored as 1-forms.
1810: Level: developer
1812: .seealso: `PetscDualSpace`, `PetscDTAltV`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceTransformType`
1813: @*/
1814: PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace dsp, PetscInt k)
1815: {
1816: PetscInt dim;
1818: PetscFunctionBeginHot;
1820: PetscCheck(!dsp->setupcalled, PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up");
1821: dim = dsp->dm->dim;
1822: 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);
1823: dsp->k = k;
1824: PetscFunctionReturn(PETSC_SUCCESS);
1825: }
1827: /*@C
1828: PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space
1830: Input Parameter:
1831: . dsp - The `PetscDualSpace`
1833: Output Parameter:
1834: . k - The simplex dimension
1836: Level: developer
1838: Note:
1839: Currently supported values are
1840: .vb
1841: 0: These are H_1 methods that only transform coordinates
1842: 1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM)
1843: 2: These are the same as 1
1844: 3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM)
1845: .ve
1847: .seealso: `PetscDualSpace`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceTransformType`
1848: @*/
1849: PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k)
1850: {
1851: PetscInt dim;
1853: PetscFunctionBeginHot;
1855: PetscAssertPointer(k, 2);
1856: dim = dsp->dm->dim;
1857: if (!dsp->k) *k = IDENTITY_TRANSFORM;
1858: else if (dsp->k == 1) *k = COVARIANT_PIOLA_TRANSFORM;
1859: else if (dsp->k == -(dim - 1)) *k = CONTRAVARIANT_PIOLA_TRANSFORM;
1860: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported transformation");
1861: PetscFunctionReturn(PETSC_SUCCESS);
1862: }
1864: /*@C
1865: PetscDualSpaceTransform - Transform the function values
1867: Input Parameters:
1868: + dsp - The `PetscDualSpace`
1869: . trans - The type of transform
1870: . isInverse - Flag to invert the transform
1871: . fegeom - The cell geometry
1872: . Nv - The number of function samples
1873: . Nc - The number of function components
1874: - vals - The function values
1876: Output Parameter:
1877: . vals - The transformed function values
1879: Level: intermediate
1881: Note:
1882: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
1884: .seealso: `PetscDualSpace`, `PetscDualSpaceTransformGradient()`, `PetscDualSpaceTransformHessian()`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransformType`
1885: @*/
1886: PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1887: {
1888: PetscReal Jstar[9] = {0};
1889: PetscInt dim, v, c, Nk;
1891: PetscFunctionBeginHot;
1893: PetscAssertPointer(fegeom, 4);
1894: PetscAssertPointer(vals, 7);
1895: /* TODO: not handling dimEmbed != dim right now */
1896: dim = dsp->dm->dim;
1897: /* No change needed for 0-forms */
1898: if (!dsp->k) PetscFunctionReturn(PETSC_SUCCESS);
1899: PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(dsp->k), &Nk));
1900: /* TODO: use fegeom->isAffine */
1901: PetscCall(PetscDTAltVPullbackMatrix(dim, dim, isInverse ? fegeom->J : fegeom->invJ, dsp->k, Jstar));
1902: for (v = 0; v < Nv; ++v) {
1903: switch (Nk) {
1904: case 1:
1905: for (c = 0; c < Nc; c++) vals[v * Nc + c] *= Jstar[0];
1906: break;
1907: case 2:
1908: for (c = 0; c < Nc; c += 2) DMPlex_Mult2DReal_Internal(Jstar, 1, &vals[v * Nc + c], &vals[v * Nc + c]);
1909: break;
1910: case 3:
1911: for (c = 0; c < Nc; c += 3) DMPlex_Mult3DReal_Internal(Jstar, 1, &vals[v * Nc + c], &vals[v * Nc + c]);
1912: break;
1913: default:
1914: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported form size %" PetscInt_FMT " for transformation", Nk);
1915: }
1916: }
1917: PetscFunctionReturn(PETSC_SUCCESS);
1918: }
1920: /*@C
1921: PetscDualSpaceTransformGradient - Transform the function gradient values
1923: Input Parameters:
1924: + dsp - The `PetscDualSpace`
1925: . trans - The type of transform
1926: . isInverse - Flag to invert the transform
1927: . fegeom - The cell geometry
1928: . Nv - The number of function gradient samples
1929: . Nc - The number of function components
1930: - vals - The function gradient values
1932: Output Parameter:
1933: . vals - The transformed function gradient values
1935: Level: intermediate
1937: Note:
1938: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
1940: .seealso: `PetscDualSpace`, `PetscDualSpaceTransform()`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransformType`
1941: @*/
1942: PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1943: {
1944: const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed;
1945: PetscInt v, c, d;
1947: PetscFunctionBeginHot;
1949: PetscAssertPointer(fegeom, 4);
1950: PetscAssertPointer(vals, 7);
1951: PetscAssert(dE > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %" PetscInt_FMT, dE);
1952: /* Transform gradient */
1953: if (dim == dE) {
1954: for (v = 0; v < Nv; ++v) {
1955: for (c = 0; c < Nc; ++c) {
1956: switch (dim) {
1957: case 1:
1958: vals[(v * Nc + c) * dim] *= fegeom->invJ[0];
1959: break;
1960: case 2:
1961: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[(v * Nc + c) * dim], &vals[(v * Nc + c) * dim]);
1962: break;
1963: case 3:
1964: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[(v * Nc + c) * dim], &vals[(v * Nc + c) * dim]);
1965: break;
1966: default:
1967: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
1968: }
1969: }
1970: }
1971: } else {
1972: for (v = 0; v < Nv; ++v) {
1973: for (c = 0; c < Nc; ++c) DMPlex_MultTransposeReal_Internal(fegeom->invJ, dim, dE, 1, &vals[(v * Nc + c) * dE], &vals[(v * Nc + c) * dE]);
1974: }
1975: }
1976: /* Assume its a vector, otherwise assume its a bunch of scalars */
1977: if (Nc == 1 || Nc != dim) PetscFunctionReturn(PETSC_SUCCESS);
1978: switch (trans) {
1979: case IDENTITY_TRANSFORM:
1980: break;
1981: case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
1982: if (isInverse) {
1983: for (v = 0; v < Nv; ++v) {
1984: for (d = 0; d < dim; ++d) {
1985: switch (dim) {
1986: case 2:
1987: DMPlex_MultTranspose2DReal_Internal(fegeom->J, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
1988: break;
1989: case 3:
1990: DMPlex_MultTranspose3DReal_Internal(fegeom->J, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
1991: break;
1992: default:
1993: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
1994: }
1995: }
1996: }
1997: } else {
1998: for (v = 0; v < Nv; ++v) {
1999: for (d = 0; d < dim; ++d) {
2000: switch (dim) {
2001: case 2:
2002: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2003: break;
2004: case 3:
2005: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2006: break;
2007: default:
2008: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2009: }
2010: }
2011: }
2012: }
2013: break;
2014: case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
2015: if (isInverse) {
2016: for (v = 0; v < Nv; ++v) {
2017: for (d = 0; d < dim; ++d) {
2018: switch (dim) {
2019: case 2:
2020: DMPlex_Mult2DReal_Internal(fegeom->invJ, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2021: break;
2022: case 3:
2023: DMPlex_Mult3DReal_Internal(fegeom->invJ, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2024: break;
2025: default:
2026: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2027: }
2028: for (c = 0; c < Nc; ++c) vals[(v * Nc + c) * dim + d] *= fegeom->detJ[0];
2029: }
2030: }
2031: } else {
2032: for (v = 0; v < Nv; ++v) {
2033: for (d = 0; d < dim; ++d) {
2034: switch (dim) {
2035: case 2:
2036: DMPlex_Mult2DReal_Internal(fegeom->J, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2037: break;
2038: case 3:
2039: DMPlex_Mult3DReal_Internal(fegeom->J, dim, &vals[v * Nc * dim + d], &vals[v * Nc * dim + d]);
2040: break;
2041: default:
2042: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2043: }
2044: for (c = 0; c < Nc; ++c) vals[(v * Nc + c) * dim + d] /= fegeom->detJ[0];
2045: }
2046: }
2047: }
2048: break;
2049: }
2050: PetscFunctionReturn(PETSC_SUCCESS);
2051: }
2053: /*@C
2054: PetscDualSpaceTransformHessian - Transform the function Hessian values
2056: Input Parameters:
2057: + dsp - The `PetscDualSpace`
2058: . trans - The type of transform
2059: . isInverse - Flag to invert the transform
2060: . fegeom - The cell geometry
2061: . Nv - The number of function Hessian samples
2062: . Nc - The number of function components
2063: - vals - The function gradient values
2065: Output Parameter:
2066: . vals - The transformed function Hessian values
2068: Level: intermediate
2070: Note:
2071: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2073: .seealso: `PetscDualSpace`, `PetscDualSpaceTransform()`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransformType`
2074: @*/
2075: PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
2076: {
2077: const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed;
2078: PetscInt v, c;
2080: PetscFunctionBeginHot;
2082: PetscAssertPointer(fegeom, 4);
2083: PetscAssertPointer(vals, 7);
2084: PetscAssert(dE > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %" PetscInt_FMT, dE);
2085: /* Transform Hessian: J^{-T}_{ik} J^{-T}_{jl} H(f)_{kl} = J^{-T}_{ik} H(f)_{kl} J^{-1}_{lj} */
2086: if (dim == dE) {
2087: for (v = 0; v < Nv; ++v) {
2088: for (c = 0; c < Nc; ++c) {
2089: switch (dim) {
2090: case 1:
2091: vals[(v * Nc + c) * dim * dim] *= PetscSqr(fegeom->invJ[0]);
2092: break;
2093: case 2:
2094: DMPlex_PTAP2DReal_Internal(fegeom->invJ, &vals[(v * Nc + c) * dim * dim], &vals[(v * Nc + c) * dim * dim]);
2095: break;
2096: case 3:
2097: DMPlex_PTAP3DReal_Internal(fegeom->invJ, &vals[(v * Nc + c) * dim * dim], &vals[(v * Nc + c) * dim * dim]);
2098: break;
2099: default:
2100: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2101: }
2102: }
2103: }
2104: } else {
2105: for (v = 0; v < Nv; ++v) {
2106: 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]);
2107: }
2108: }
2109: /* Assume its a vector, otherwise assume its a bunch of scalars */
2110: if (Nc == 1 || Nc != dim) PetscFunctionReturn(PETSC_SUCCESS);
2111: switch (trans) {
2112: case IDENTITY_TRANSFORM:
2113: break;
2114: case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
2115: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Piola mapping for Hessians not yet supported");
2116: case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
2117: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Piola mapping for Hessians not yet supported");
2118: }
2119: PetscFunctionReturn(PETSC_SUCCESS);
2120: }
2122: /*@C
2123: 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.
2125: Input Parameters:
2126: + dsp - The `PetscDualSpace`
2127: . fegeom - The geometry for this cell
2128: . Nq - The number of function samples
2129: . Nc - The number of function components
2130: - pointEval - The function values
2132: Output Parameter:
2133: . pointEval - The transformed function values
2135: Level: advanced
2137: Notes:
2138: 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.
2140: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2142: .seealso: `PetscDualSpace`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2143: @*/
2144: PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2145: {
2146: PetscDualSpaceTransformType trans;
2147: PetscInt k;
2149: PetscFunctionBeginHot;
2151: PetscAssertPointer(fegeom, 2);
2152: PetscAssertPointer(pointEval, 5);
2153: /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2154: This determines their transformation properties. */
2155: PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2156: switch (k) {
2157: case 0: /* H^1 point evaluations */
2158: trans = IDENTITY_TRANSFORM;
2159: break;
2160: case 1: /* Hcurl preserves tangential edge traces */
2161: trans = COVARIANT_PIOLA_TRANSFORM;
2162: break;
2163: case 2:
2164: case 3: /* Hdiv preserve normal traces */
2165: trans = CONTRAVARIANT_PIOLA_TRANSFORM;
2166: break;
2167: default:
2168: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2169: }
2170: PetscCall(PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, pointEval));
2171: PetscFunctionReturn(PETSC_SUCCESS);
2172: }
2174: /*@C
2175: 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.
2177: Input Parameters:
2178: + dsp - The `PetscDualSpace`
2179: . fegeom - The geometry for this cell
2180: . Nq - The number of function samples
2181: . Nc - The number of function components
2182: - pointEval - The function values
2184: Output Parameter:
2185: . pointEval - The transformed function values
2187: Level: advanced
2189: Notes:
2190: 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.
2192: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2194: .seealso: `PetscDualSpace`, `PetscDualSpacePullback()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2195: @*/
2196: PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2197: {
2198: PetscDualSpaceTransformType trans;
2199: PetscInt k;
2201: PetscFunctionBeginHot;
2203: PetscAssertPointer(fegeom, 2);
2204: PetscAssertPointer(pointEval, 5);
2205: /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2206: This determines their transformation properties. */
2207: PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2208: switch (k) {
2209: case 0: /* H^1 point evaluations */
2210: trans = IDENTITY_TRANSFORM;
2211: break;
2212: case 1: /* Hcurl preserves tangential edge traces */
2213: trans = COVARIANT_PIOLA_TRANSFORM;
2214: break;
2215: case 2:
2216: case 3: /* Hdiv preserve normal traces */
2217: trans = CONTRAVARIANT_PIOLA_TRANSFORM;
2218: break;
2219: default:
2220: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2221: }
2222: PetscCall(PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval));
2223: PetscFunctionReturn(PETSC_SUCCESS);
2224: }
2226: /*@C
2227: 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.
2229: Input Parameters:
2230: + dsp - The `PetscDualSpace`
2231: . fegeom - The geometry for this cell
2232: . Nq - The number of function gradient samples
2233: . Nc - The number of function components
2234: - pointEval - The function gradient values
2236: Output Parameter:
2237: . pointEval - The transformed function gradient values
2239: Level: advanced
2241: Notes:
2242: 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.
2244: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2246: .seealso: `PetscDualSpace`, `PetscDualSpacePushforward()`, `PetscDualSpacePullback()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2247: @*/
2248: PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2249: {
2250: PetscDualSpaceTransformType trans;
2251: PetscInt k;
2253: PetscFunctionBeginHot;
2255: PetscAssertPointer(fegeom, 2);
2256: PetscAssertPointer(pointEval, 5);
2257: /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2258: This determines their transformation properties. */
2259: PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2260: switch (k) {
2261: case 0: /* H^1 point evaluations */
2262: trans = IDENTITY_TRANSFORM;
2263: break;
2264: case 1: /* Hcurl preserves tangential edge traces */
2265: trans = COVARIANT_PIOLA_TRANSFORM;
2266: break;
2267: case 2:
2268: case 3: /* Hdiv preserve normal traces */
2269: trans = CONTRAVARIANT_PIOLA_TRANSFORM;
2270: break;
2271: default:
2272: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2273: }
2274: PetscCall(PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval));
2275: PetscFunctionReturn(PETSC_SUCCESS);
2276: }
2278: /*@C
2279: 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.
2281: Input Parameters:
2282: + dsp - The `PetscDualSpace`
2283: . fegeom - The geometry for this cell
2284: . Nq - The number of function Hessian samples
2285: . Nc - The number of function components
2286: - pointEval - The function gradient values
2288: Output Parameter:
2289: . pointEval - The transformed function Hessian values
2291: Level: advanced
2293: Notes:
2294: 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.
2296: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2298: .seealso: `PetscDualSpace`, `PetscDualSpacePushforward()`, `PetscDualSpacePullback()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2299: @*/
2300: PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2301: {
2302: PetscDualSpaceTransformType trans;
2303: PetscInt k;
2305: PetscFunctionBeginHot;
2307: PetscAssertPointer(fegeom, 2);
2308: PetscAssertPointer(pointEval, 5);
2309: /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2310: This determines their transformation properties. */
2311: PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2312: switch (k) {
2313: case 0: /* H^1 point evaluations */
2314: trans = IDENTITY_TRANSFORM;
2315: break;
2316: case 1: /* Hcurl preserves tangential edge traces */
2317: trans = COVARIANT_PIOLA_TRANSFORM;
2318: break;
2319: case 2:
2320: case 3: /* Hdiv preserve normal traces */
2321: trans = CONTRAVARIANT_PIOLA_TRANSFORM;
2322: break;
2323: default:
2324: SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2325: }
2326: PetscCall(PetscDualSpaceTransformHessian(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval));
2327: PetscFunctionReturn(PETSC_SUCCESS);
2328: }