Actual source code: dualspace.c
petsc-3.14.6 2021-03-30
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: const char *const PetscDualSpaceReferenceCells[] = {"SIMPLEX", "TENSOR", "PetscDualSpaceReferenceCell", "PETSCDUALSPACE_REFCELL_", NULL};
13: /*
14: PetscDualSpaceLatticePointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to at most 'max'.
15: Ordering is lexicographic with lowest index as least significant in ordering.
16: 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}.
18: Input Parameters:
19: + len - The length of the tuple
20: . max - The maximum sum
21: - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
23: Output Parameter:
24: . tup - A tuple of len integers whos sum is at most 'max'
26: Level: developer
28: .seealso: PetscDualSpaceTensorPointLexicographic_Internal()
29: */
30: PetscErrorCode PetscDualSpaceLatticePointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
31: {
33: while (len--) {
34: max -= tup[len];
35: if (!max) {
36: tup[len] = 0;
37: break;
38: }
39: }
40: tup[++len]++;
41: return(0);
42: }
44: /*
45: PetscDualSpaceTensorPointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max'.
46: Ordering is lexicographic with lowest index as least significant in ordering.
47: 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}.
49: Input Parameters:
50: + len - The length of the tuple
51: . max - The maximum value
52: - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
54: Output Parameter:
55: . tup - A tuple of len integers whos sum is at most 'max'
57: Level: developer
59: .seealso: PetscDualSpaceLatticePointLexicographic_Internal()
60: */
61: PetscErrorCode PetscDualSpaceTensorPointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
62: {
63: PetscInt i;
66: for (i = 0; i < len; i++) {
67: if (tup[i] < max) {
68: break;
69: } else {
70: tup[i] = 0;
71: }
72: }
73: tup[i]++;
74: return(0);
75: }
77: /*@C
78: PetscDualSpaceRegister - Adds a new PetscDualSpace implementation
80: Not Collective
82: Input Parameters:
83: + name - The name of a new user-defined creation routine
84: - create_func - The creation routine itself
86: Notes:
87: PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces
89: Sample usage:
90: .vb
91: PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate);
92: .ve
94: Then, your PetscDualSpace type can be chosen with the procedural interface via
95: .vb
96: PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
97: PetscDualSpaceSetType(PetscDualSpace, "my_dual_space");
98: .ve
99: or at runtime via the option
100: .vb
101: -petscdualspace_type my_dual_space
102: .ve
104: Level: advanced
106: .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy()
108: @*/
109: PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace))
110: {
114: PetscFunctionListAdd(&PetscDualSpaceList, sname, function);
115: return(0);
116: }
118: /*@C
119: PetscDualSpaceSetType - Builds a particular PetscDualSpace
121: Collective on sp
123: Input Parameters:
124: + sp - The PetscDualSpace object
125: - name - The kind of space
127: Options Database Key:
128: . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types
130: Level: intermediate
132: .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate()
133: @*/
134: PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name)
135: {
136: PetscErrorCode (*r)(PetscDualSpace);
137: PetscBool match;
142: PetscObjectTypeCompare((PetscObject) sp, name, &match);
143: if (match) return(0);
145: if (!PetscDualSpaceRegisterAllCalled) {PetscDualSpaceRegisterAll();}
146: PetscFunctionListFind(PetscDualSpaceList, name, &r);
147: if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name);
149: if (sp->ops->destroy) {
150: (*sp->ops->destroy)(sp);
151: sp->ops->destroy = NULL;
152: }
153: (*r)(sp);
154: PetscObjectChangeTypeName((PetscObject) sp, name);
155: return(0);
156: }
158: /*@C
159: PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object.
161: Not Collective
163: Input Parameter:
164: . sp - The PetscDualSpace
166: Output Parameter:
167: . name - The PetscDualSpace type name
169: Level: intermediate
171: .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate()
172: @*/
173: PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name)
174: {
180: if (!PetscDualSpaceRegisterAllCalled) {
181: PetscDualSpaceRegisterAll();
182: }
183: *name = ((PetscObject) sp)->type_name;
184: return(0);
185: }
187: static PetscErrorCode PetscDualSpaceView_ASCII(PetscDualSpace sp, PetscViewer v)
188: {
189: PetscViewerFormat format;
190: PetscInt pdim, f;
191: PetscErrorCode ierr;
194: PetscDualSpaceGetDimension(sp, &pdim);
195: PetscObjectPrintClassNamePrefixType((PetscObject) sp, v);
196: PetscViewerASCIIPushTab(v);
197: if (sp->k) {
198: PetscViewerASCIIPrintf(v, "Dual space for %D-forms %swith %D components, size %D\n", PetscAbsInt(sp->k), sp->k < 0 ? "(stored in dual form) ": "", sp->Nc, pdim);
199: } else {
200: PetscViewerASCIIPrintf(v, "Dual space with %D components, size %D\n", sp->Nc, pdim);
201: }
202: if (sp->ops->view) {(*sp->ops->view)(sp, v);}
203: PetscViewerGetFormat(v, &format);
204: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
205: PetscViewerASCIIPushTab(v);
206: for (f = 0; f < pdim; ++f) {
207: PetscViewerASCIIPrintf(v, "Dual basis vector %D\n", f);
208: PetscViewerASCIIPushTab(v);
209: PetscQuadratureView(sp->functional[f], v);
210: PetscViewerASCIIPopTab(v);
211: }
212: PetscViewerASCIIPopTab(v);
213: }
214: PetscViewerASCIIPopTab(v);
215: return(0);
216: }
218: /*@C
219: PetscDualSpaceViewFromOptions - View from Options
221: Collective on PetscDualSpace
223: Input Parameters:
224: + A - the PetscDualSpace object
225: . obj - Optional object, proivides prefix
226: - name - command line option
228: Level: intermediate
229: .seealso: PetscDualSpace, PetscDualSpaceView(), PetscObjectViewFromOptions(), PetscDualSpaceCreate()
230: @*/
231: PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace A,PetscObject obj,const char name[])
232: {
237: PetscObjectViewFromOptions((PetscObject)A,obj,name);
238: return(0);
239: }
241: /*@
242: PetscDualSpaceView - Views a PetscDualSpace
244: Collective on sp
246: Input Parameter:
247: + sp - the PetscDualSpace object to view
248: - v - the viewer
250: Level: beginner
252: .seealso PetscDualSpaceDestroy(), PetscDualSpace
253: @*/
254: PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
255: {
256: PetscBool iascii;
262: if (!v) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);}
263: PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);
264: if (iascii) {PetscDualSpaceView_ASCII(sp, v);}
265: return(0);
266: }
268: /*@
269: PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database
271: Collective on sp
273: Input Parameter:
274: . sp - the PetscDualSpace object to set options for
276: Options Database:
277: . -petscspace_degree the approximation order of the space
279: Level: intermediate
281: .seealso PetscDualSpaceView(), PetscDualSpace, PetscObjectSetFromOptions()
282: @*/
283: PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
284: {
285: PetscDualSpaceReferenceCell refCell = PETSCDUALSPACE_REFCELL_SIMPLEX;
286: PetscInt refDim = 0;
287: PetscBool flg;
288: const char *defaultType;
289: char name[256];
290: PetscErrorCode ierr;
294: if (!((PetscObject) sp)->type_name) {
295: defaultType = PETSCDUALSPACELAGRANGE;
296: } else {
297: defaultType = ((PetscObject) sp)->type_name;
298: }
299: if (!PetscSpaceRegisterAllCalled) {PetscSpaceRegisterAll();}
301: PetscObjectOptionsBegin((PetscObject) sp);
302: PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);
303: if (flg) {
304: PetscDualSpaceSetType(sp, name);
305: } else if (!((PetscObject) sp)->type_name) {
306: PetscDualSpaceSetType(sp, defaultType);
307: }
308: PetscOptionsBoundedInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL,0);
309: PetscOptionsInt("-petscdualspace_form_degree", "The form degree of the dofs", "PetscDualSpaceSetFormDegree", sp->k, &sp->k, NULL);
310: PetscOptionsBoundedInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL,1);
311: if (sp->ops->setfromoptions) {
312: (*sp->ops->setfromoptions)(PetscOptionsObject,sp);
313: }
314: PetscOptionsBoundedInt("-petscdualspace_refdim", "The spatial dimension of the reference cell", "PetscDualSpaceSetReferenceCell", refDim, &refDim, NULL,0);
315: PetscOptionsEnum("-petscdualspace_refcell", "Reference cell", "PetscDualSpaceSetReferenceCell", PetscDualSpaceReferenceCells, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg);
316: if (flg) {
317: DM K;
319: if (!refDim) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_INCOMP, "Reference cell specified without a dimension. Use -petscdualspace_refdim.");
320: PetscDualSpaceCreateReferenceCell(sp, refDim, refCell == PETSCDUALSPACE_REFCELL_SIMPLEX ? PETSC_TRUE : PETSC_FALSE, &K);
321: PetscDualSpaceSetDM(sp, K);
322: DMDestroy(&K);
323: }
325: /* process any options handlers added with PetscObjectAddOptionsHandler() */
326: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);
327: PetscOptionsEnd();
328: sp->setfromoptionscalled = PETSC_TRUE;
329: return(0);
330: }
332: /*@
333: PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace
335: Collective on sp
337: Input Parameter:
338: . sp - the PetscDualSpace object to setup
340: Level: intermediate
342: .seealso PetscDualSpaceView(), PetscDualSpaceDestroy(), PetscDualSpace
343: @*/
344: PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
345: {
350: if (sp->setupcalled) return(0);
351: PetscLogEventBegin(PETSCDUALSPACE_SetUp, sp, 0, 0, 0);
352: sp->setupcalled = PETSC_TRUE;
353: if (sp->ops->setup) {(*sp->ops->setup)(sp);}
354: PetscLogEventEnd(PETSCDUALSPACE_SetUp, sp, 0, 0, 0);
355: if (sp->setfromoptionscalled) {PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");}
356: return(0);
357: }
359: static PetscErrorCode PetscDualSpaceClearDMData_Internal(PetscDualSpace sp, DM dm)
360: {
361: PetscInt pStart = -1, pEnd = -1, depth = -1;
365: if (!dm) return(0);
366: DMPlexGetChart(dm, &pStart, &pEnd);
367: DMPlexGetDepth(dm, &depth);
369: if (sp->pointSpaces) {
370: PetscInt i;
372: for (i = 0; i < pEnd - pStart; i++) {
373: PetscDualSpaceDestroy(&(sp->pointSpaces[i]));
374: }
375: }
376: PetscFree(sp->pointSpaces);
378: if (sp->heightSpaces) {
379: PetscInt i;
381: for (i = 0; i <= depth; i++) {
382: PetscDualSpaceDestroy(&(sp->heightSpaces[i]));
383: }
384: }
385: PetscFree(sp->heightSpaces);
387: PetscSectionDestroy(&(sp->pointSection));
388: PetscQuadratureDestroy(&(sp->intNodes));
389: VecDestroy(&(sp->intDofValues));
390: VecDestroy(&(sp->intNodeValues));
391: MatDestroy(&(sp->intMat));
392: PetscQuadratureDestroy(&(sp->allNodes));
393: VecDestroy(&(sp->allDofValues));
394: VecDestroy(&(sp->allNodeValues));
395: MatDestroy(&(sp->allMat));
396: PetscFree(sp->numDof);
397: return(0);
398: }
401: /*@
402: PetscDualSpaceDestroy - Destroys a PetscDualSpace object
404: Collective on sp
406: Input Parameter:
407: . sp - the PetscDualSpace object to destroy
409: Level: beginner
411: .seealso PetscDualSpaceView(), PetscDualSpace(), PetscDualSpaceCreate()
412: @*/
413: PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
414: {
415: PetscInt dim, f;
416: DM dm;
420: if (!*sp) return(0);
423: if (--((PetscObject)(*sp))->refct > 0) {*sp = NULL; return(0);}
424: ((PetscObject) (*sp))->refct = 0;
426: PetscDualSpaceGetDimension(*sp, &dim);
427: dm = (*sp)->dm;
429: if ((*sp)->ops->destroy) {(*(*sp)->ops->destroy)(*sp);}
430: PetscDualSpaceClearDMData_Internal(*sp, dm);
432: for (f = 0; f < dim; ++f) {
433: PetscQuadratureDestroy(&(*sp)->functional[f]);
434: }
435: PetscFree((*sp)->functional);
436: DMDestroy(&(*sp)->dm);
437: PetscHeaderDestroy(sp);
438: return(0);
439: }
441: /*@
442: PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType().
444: Collective
446: Input Parameter:
447: . comm - The communicator for the PetscDualSpace object
449: Output Parameter:
450: . sp - The PetscDualSpace object
452: Level: beginner
454: .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE
455: @*/
456: PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
457: {
458: PetscDualSpace s;
463: PetscCitationsRegister(FECitation,&FEcite);
464: *sp = NULL;
465: PetscFEInitializePackage();
467: PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);
469: s->order = 0;
470: s->Nc = 1;
471: s->k = 0;
472: s->spdim = -1;
473: s->spintdim = -1;
474: s->uniform = PETSC_TRUE;
475: s->setupcalled = PETSC_FALSE;
477: *sp = s;
478: return(0);
479: }
481: /*@
482: PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup.
484: Collective on sp
486: Input Parameter:
487: . sp - The original PetscDualSpace
489: Output Parameter:
490: . spNew - The duplicate PetscDualSpace
492: Level: beginner
494: .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType()
495: @*/
496: PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
497: {
498: DM dm;
499: PetscDualSpaceType type;
500: const char *name;
506: PetscDualSpaceCreate(PetscObjectComm((PetscObject)sp), spNew);
507: PetscObjectGetName((PetscObject) sp, &name);
508: PetscObjectSetName((PetscObject) *spNew, name);
509: PetscDualSpaceGetType(sp, &type);
510: PetscDualSpaceSetType(*spNew, type);
511: PetscDualSpaceGetDM(sp, &dm);
512: PetscDualSpaceSetDM(*spNew, dm);
514: (*spNew)->order = sp->order;
515: (*spNew)->k = sp->k;
516: (*spNew)->Nc = sp->Nc;
517: (*spNew)->uniform = sp->uniform;
518: if (sp->ops->duplicate) {(*sp->ops->duplicate)(sp, *spNew);}
519: return(0);
520: }
522: /*@
523: PetscDualSpaceGetDM - Get the DM representing the reference cell
525: Not collective
527: Input Parameter:
528: . sp - The PetscDualSpace
530: Output Parameter:
531: . dm - The reference cell
533: Level: intermediate
535: .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate()
536: @*/
537: PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
538: {
542: *dm = sp->dm;
543: return(0);
544: }
546: /*@
547: PetscDualSpaceSetDM - Get the DM representing the reference cell
549: Not collective
551: Input Parameters:
552: + sp - The PetscDualSpace
553: - dm - The reference cell
555: Level: intermediate
557: .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate()
558: @*/
559: PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
560: {
566: if (sp->setupcalled) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change DM after dualspace is set up");
567: PetscObjectReference((PetscObject) dm);
568: if (sp->dm && sp->dm != dm) {
569: PetscDualSpaceClearDMData_Internal(sp, sp->dm);
570: }
571: DMDestroy(&sp->dm);
572: sp->dm = dm;
573: return(0);
574: }
576: /*@
577: PetscDualSpaceGetOrder - Get the order of the dual space
579: Not collective
581: Input Parameter:
582: . sp - The PetscDualSpace
584: Output Parameter:
585: . order - The order
587: Level: intermediate
589: .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate()
590: @*/
591: PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
592: {
596: *order = sp->order;
597: return(0);
598: }
600: /*@
601: PetscDualSpaceSetOrder - Set the order of the dual space
603: Not collective
605: Input Parameters:
606: + sp - The PetscDualSpace
607: - order - The order
609: Level: intermediate
611: .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate()
612: @*/
613: PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
614: {
617: if (sp->setupcalled) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change order after dualspace is set up");
618: sp->order = order;
619: return(0);
620: }
622: /*@
623: PetscDualSpaceGetNumComponents - Return the number of components for this space
625: Input Parameter:
626: . sp - The PetscDualSpace
628: Output Parameter:
629: . Nc - The number of components
631: Note: A vector space, for example, will have d components, where d is the spatial dimension
633: Level: intermediate
635: .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace
636: @*/
637: PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
638: {
642: *Nc = sp->Nc;
643: return(0);
644: }
646: /*@
647: PetscDualSpaceSetNumComponents - Set the number of components for this space
649: Input Parameters:
650: + sp - The PetscDualSpace
651: - order - The number of components
653: Level: intermediate
655: .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace
656: @*/
657: PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
658: {
661: if (sp->setupcalled) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up");
662: sp->Nc = Nc;
663: return(0);
664: }
666: /*@
667: PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
669: Not collective
671: Input Parameters:
672: + sp - The PetscDualSpace
673: - i - The basis number
675: Output Parameter:
676: . functional - The basis functional
678: Level: intermediate
680: .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate()
681: @*/
682: PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
683: {
684: PetscInt dim;
690: PetscDualSpaceGetDimension(sp, &dim);
691: if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim);
692: *functional = sp->functional[i];
693: return(0);
694: }
696: /*@
697: PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
699: Not collective
701: Input Parameter:
702: . sp - The PetscDualSpace
704: Output Parameter:
705: . dim - The dimension
707: Level: intermediate
709: .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
710: @*/
711: PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
712: {
718: if (sp->spdim < 0) {
719: PetscSection section;
721: PetscDualSpaceGetSection(sp, §ion);
722: if (section) {
723: PetscSectionGetStorageSize(section, &(sp->spdim));
724: } else sp->spdim = 0;
725: }
726: *dim = sp->spdim;
727: return(0);
728: }
730: /*@
731: PetscDualSpaceGetInteriorDimension - Get the interior dimension of the dual space, i.e. the number of basis functionals assigned to the interior of the reference domain
733: Not collective
735: Input Parameter:
736: . sp - The PetscDualSpace
738: Output Parameter:
739: . dim - The dimension
741: Level: intermediate
743: .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
744: @*/
745: PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace sp, PetscInt *intdim)
746: {
752: if (sp->spintdim < 0) {
753: PetscSection section;
755: PetscDualSpaceGetSection(sp, §ion);
756: if (section) {
757: PetscSectionGetConstrainedStorageSize(section, &(sp->spintdim));
758: } else sp->spintdim = 0;
759: }
760: *intdim = sp->spintdim;
761: return(0);
762: }
764: /*@
765: PetscDualSpaceGetUniform - Whether this dual space is uniform
767: Not collective
769: Input Parameters:
770: . sp - A dual space
772: Output Parameters:
773: . uniform - PETSC_TRUE if (a) the dual space is the same for each point in a stratum of the reference DMPlex, and
774: (b) every symmetry of each point in the reference DMPlex is also a symmetry of the point's dual space.
777: Level: advanced
779: Note: all of the usual spaces on simplex or tensor-product elements will be uniform, only reference cells
780: with non-uniform strata (like trianguar-prisms) or anisotropic hp dual spaces will not be uniform.
782: .seealso: PetscDualSpaceGetPointSubspace(), PetscDualSpaceGetSymmetries()
783: @*/
784: PetscErrorCode PetscDualSpaceGetUniform(PetscDualSpace sp, PetscBool *uniform)
785: {
789: *uniform = sp->uniform;
790: return(0);
791: }
794: /*@C
795: PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
797: Not collective
799: Input Parameter:
800: . sp - The PetscDualSpace
802: Output Parameter:
803: . numDof - An array of length dim+1 which holds the number of dofs for each dimension
805: Level: intermediate
807: .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
808: @*/
809: PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof)
810: {
816: if (!sp->uniform) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A non-uniform space does not have a fixed number of dofs for each height");
817: if (!sp->numDof) {
818: DM dm;
819: PetscInt depth, d;
820: PetscSection section;
822: PetscDualSpaceGetDM(sp, &dm);
823: DMPlexGetDepth(dm, &depth);
824: PetscCalloc1(depth+1,&(sp->numDof));
825: PetscDualSpaceGetSection(sp, §ion);
826: for (d = 0; d <= depth; d++) {
827: PetscInt dStart, dEnd;
829: DMPlexGetDepthStratum(dm, d, &dStart, &dEnd);
830: if (dEnd <= dStart) continue;
831: PetscSectionGetDof(section, dStart, &(sp->numDof[d]));
833: }
834: }
835: *numDof = sp->numDof;
836: if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
837: return(0);
838: }
840: /* create the section of the right size and set a permutation for topological ordering */
841: PetscErrorCode PetscDualSpaceSectionCreate_Internal(PetscDualSpace sp, PetscSection *topSection)
842: {
843: DM dm;
844: PetscInt pStart, pEnd, cStart, cEnd, c, depth, count, i;
845: PetscInt *seen, *perm;
846: PetscSection section;
850: dm = sp->dm;
851: PetscSectionCreate(PETSC_COMM_SELF, §ion);
852: DMPlexGetChart(dm, &pStart, &pEnd);
853: PetscSectionSetChart(section, pStart, pEnd);
854: PetscCalloc1(pEnd - pStart, &seen);
855: PetscMalloc1(pEnd - pStart, &perm);
856: DMPlexGetDepth(dm, &depth);
857: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
858: for (c = cStart, count = 0; c < cEnd; c++) {
859: PetscInt closureSize = -1, e;
860: PetscInt *closure = NULL;
862: perm[count++] = c;
863: seen[c-pStart] = 1;
864: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
865: for (e = 0; e < closureSize; e++) {
866: PetscInt point = closure[2*e];
868: if (seen[point-pStart]) continue;
869: perm[count++] = point;
870: seen[point-pStart] = 1;
871: }
872: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
873: }
874: if (count != pEnd - pStart) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad topological ordering");
875: for (i = 0; i < pEnd - pStart; i++) if (perm[i] != i) break;
876: if (i < pEnd - pStart) {
877: IS permIS;
879: ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, perm, PETSC_OWN_POINTER, &permIS);
880: ISSetPermutation(permIS);
881: PetscSectionSetPermutation(section, permIS);
882: ISDestroy(&permIS);
883: } else {
884: PetscFree(perm);
885: }
886: PetscFree(seen);
887: *topSection = section;
888: return(0);
889: }
891: /* mark boundary points and set up */
892: PetscErrorCode PetscDualSpaceSectionSetUp_Internal(PetscDualSpace sp, PetscSection section)
893: {
894: DM dm;
895: DMLabel boundary;
896: PetscInt pStart, pEnd, p;
900: dm = sp->dm;
901: DMLabelCreate(PETSC_COMM_SELF,"boundary",&boundary);
902: PetscDualSpaceGetDM(sp,&dm);
903: DMPlexMarkBoundaryFaces(dm,1,boundary);
904: DMPlexLabelComplete(dm,boundary);
905: DMPlexGetChart(dm, &pStart, &pEnd);
906: for (p = pStart; p < pEnd; p++) {
907: PetscInt bval;
909: DMLabelGetValue(boundary, p, &bval);
910: if (bval == 1) {
911: PetscInt dof;
913: PetscSectionGetDof(section, p, &dof);
914: PetscSectionSetConstraintDof(section, p, dof);
915: }
916: }
917: DMLabelDestroy(&boundary);
918: PetscSectionSetUp(section);
919: return(0);
920: }
922: /*@
923: PetscDualSpaceGetSection - Create a PetscSection over the reference cell with the layout from this space
925: Collective on sp
927: Input Parameters:
928: . sp - The PetscDualSpace
930: Output Parameter:
931: . section - The section
933: Level: advanced
935: .seealso: PetscDualSpaceCreate(), DMPLEX
936: @*/
937: PetscErrorCode PetscDualSpaceGetSection(PetscDualSpace sp, PetscSection *section)
938: {
939: PetscInt pStart, pEnd, p;
943: if (!sp->pointSection) {
944: /* mark the boundary */
945: PetscDualSpaceSectionCreate_Internal(sp, &(sp->pointSection));
946: DMPlexGetChart(sp->dm,&pStart,&pEnd);
947: for (p = pStart; p < pEnd; p++) {
948: PetscDualSpace psp;
950: PetscDualSpaceGetPointSubspace(sp, p, &psp);
951: if (psp) {
952: PetscInt dof;
954: PetscDualSpaceGetInteriorDimension(psp, &dof);
955: PetscSectionSetDof(sp->pointSection,p,dof);
956: }
957: }
958: PetscDualSpaceSectionSetUp_Internal(sp,sp->pointSection);
959: }
960: *section = sp->pointSection;
961: return(0);
962: }
964: /* this assumes that all of the point dual spaces store their interior dofs first, which is true when the point DMs
965: * have one cell */
966: PetscErrorCode PetscDualSpacePushForwardSubspaces_Internal(PetscDualSpace sp, PetscInt sStart, PetscInt sEnd)
967: {
968: PetscReal *sv0, *v0, *J;
969: PetscSection section;
970: PetscInt dim, s, k;
971: DM dm;
975: PetscDualSpaceGetDM(sp, &dm);
976: DMGetDimension(dm, &dim);
977: PetscDualSpaceGetSection(sp, §ion);
978: PetscMalloc3(dim, &v0, dim, &sv0, dim*dim, &J);
979: PetscDualSpaceGetFormDegree(sp, &k);
980: for (s = sStart; s < sEnd; s++) {
981: PetscReal detJ, hdetJ;
982: PetscDualSpace ssp;
983: PetscInt dof, off, f, sdim;
984: PetscInt i, j;
985: DM sdm;
987: PetscDualSpaceGetPointSubspace(sp, s, &ssp);
988: if (!ssp) continue;
989: PetscSectionGetDof(section, s, &dof);
990: PetscSectionGetOffset(section, s, &off);
991: /* get the first vertex of the reference cell */
992: PetscDualSpaceGetDM(ssp, &sdm);
993: DMGetDimension(sdm, &sdim);
994: DMPlexComputeCellGeometryAffineFEM(sdm, 0, sv0, NULL, NULL, &hdetJ);
995: DMPlexComputeCellGeometryAffineFEM(dm, s, v0, J, NULL, &detJ);
996: /* compactify Jacobian */
997: for (i = 0; i < dim; i++) for (j = 0; j < sdim; j++) J[i* sdim + j] = J[i * dim + j];
998: for (f = 0; f < dof; f++) {
999: PetscQuadrature fn;
1001: PetscDualSpaceGetFunctional(ssp, f, &fn);
1002: PetscQuadraturePushForward(fn, dim, sv0, v0, J, k, &(sp->functional[off+f]));
1003: }
1004: }
1005: PetscFree3(v0, sv0, J);
1006: return(0);
1007: }
1009: /*@
1010: PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
1012: Collective on sp
1014: Input Parameters:
1015: + sp - The PetscDualSpace
1016: . dim - The spatial dimension
1017: - simplex - Flag for simplex, otherwise use a tensor-product cell
1019: Output Parameter:
1020: . refdm - The reference cell
1022: Level: intermediate
1024: .seealso: PetscDualSpaceCreate(), DMPLEX
1025: @*/
1026: PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm)
1027: {
1031: DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);
1032: return(0);
1033: }
1035: /*@C
1036: PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
1038: Input Parameters:
1039: + sp - The PetscDualSpace object
1040: . f - The basis functional index
1041: . time - The time
1042: . 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)
1043: . numComp - The number of components for the function
1044: . func - The input function
1045: - ctx - A context for the function
1047: Output Parameter:
1048: . value - numComp output values
1050: Note: The calling sequence for the callback func is given by:
1052: $ func(PetscInt dim, PetscReal time, const PetscReal x[],
1053: $ PetscInt numComponents, PetscScalar values[], void *ctx)
1055: Level: beginner
1057: .seealso: PetscDualSpaceCreate()
1058: @*/
1059: 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)
1060: {
1067: (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);
1068: return(0);
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: beginner
1083: .seealso: PetscDualSpaceCreate()
1084: @*/
1085: PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1086: {
1091: (*sp->ops->applyall)(sp, pointEval, spValue);
1092: return(0);
1093: }
1095: /*@C
1096: PetscDualSpaceApplyInterior - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetInteriorData()
1098: Input Parameters:
1099: + sp - The PetscDualSpace object
1100: - pointEval - Evaluation at the points returned by PetscDualSpaceGetInteriorData()
1102: Output Parameter:
1103: . spValue - The values of interior dual space functionals
1105: Level: beginner
1107: .seealso: PetscDualSpaceCreate()
1108: @*/
1109: PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1110: {
1115: (*sp->ops->applyint)(sp, pointEval, spValue);
1116: return(0);
1117: }
1119: /*@C
1120: PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
1122: Input Parameters:
1123: + sp - The PetscDualSpace object
1124: . f - The basis functional index
1125: . time - The time
1126: . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
1127: . Nc - The number of components for the function
1128: . func - The input function
1129: - ctx - A context for the function
1131: Output Parameter:
1132: . value - The output value
1134: Note: The calling sequence for the callback func is given by:
1136: $ func(PetscInt dim, PetscReal time, const PetscReal x[],
1137: $ PetscInt numComponents, PetscScalar values[], void *ctx)
1139: and the idea is to evaluate the functional as an integral
1141: $ n(f) = int dx n(x) . f(x)
1143: where both n and f have Nc components.
1145: Level: beginner
1147: .seealso: PetscDualSpaceCreate()
1148: @*/
1149: 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)
1150: {
1151: DM dm;
1152: PetscQuadrature n;
1153: const PetscReal *points, *weights;
1154: PetscReal x[3];
1155: PetscScalar *val;
1156: PetscInt dim, dE, qNc, c, Nq, q;
1157: PetscBool isAffine;
1158: PetscErrorCode ierr;
1163: PetscDualSpaceGetDM(sp, &dm);
1164: PetscDualSpaceGetFunctional(sp, f, &n);
1165: PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);
1166: if (dim != cgeom->dim) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %D != cell geometry dimension %D", dim, cgeom->dim);
1167: if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
1168: DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);
1169: *value = 0.0;
1170: isAffine = cgeom->isAffine;
1171: dE = cgeom->dimEmbed;
1172: for (q = 0; q < Nq; ++q) {
1173: if (isAffine) {
1174: CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x);
1175: (*func)(dE, time, x, Nc, val, ctx);
1176: } else {
1177: (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);
1178: }
1179: for (c = 0; c < Nc; ++c) {
1180: *value += val[c]*weights[q*Nc+c];
1181: }
1182: }
1183: DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);
1184: return(0);
1185: }
1187: /*@C
1188: PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllData()
1190: Input Parameters:
1191: + sp - The PetscDualSpace object
1192: - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllData()
1194: Output Parameter:
1195: . spValue - The values of all dual space functionals
1197: Level: beginner
1199: .seealso: PetscDualSpaceCreate()
1200: @*/
1201: PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1202: {
1203: Vec pointValues, dofValues;
1204: Mat allMat;
1205: PetscErrorCode ierr;
1211: PetscDualSpaceGetAllData(sp, NULL, &allMat);
1212: if (!(sp->allNodeValues)) {
1213: MatCreateVecs(allMat, &(sp->allNodeValues), NULL);
1214: }
1215: pointValues = sp->allNodeValues;
1216: if (!(sp->allDofValues)) {
1217: MatCreateVecs(allMat, NULL, &(sp->allDofValues));
1218: }
1219: dofValues = sp->allDofValues;
1220: VecPlaceArray(pointValues, pointEval);
1221: VecPlaceArray(dofValues, spValue);
1222: MatMult(allMat, pointValues, dofValues);
1223: VecResetArray(dofValues);
1224: VecResetArray(pointValues);
1225: return(0);
1226: }
1228: /*@C
1229: PetscDualSpaceApplyInteriorDefault - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetInteriorData()
1231: Input Parameters:
1232: + sp - The PetscDualSpace object
1233: - pointEval - Evaluation at the points returned by PetscDualSpaceGetInteriorData()
1235: Output Parameter:
1236: . spValue - The values of interior dual space functionals
1238: Level: beginner
1240: .seealso: PetscDualSpaceCreate()
1241: @*/
1242: PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1243: {
1244: Vec pointValues, dofValues;
1245: Mat intMat;
1246: PetscErrorCode ierr;
1252: PetscDualSpaceGetInteriorData(sp, NULL, &intMat);
1253: if (!(sp->intNodeValues)) {
1254: MatCreateVecs(intMat, &(sp->intNodeValues), NULL);
1255: }
1256: pointValues = sp->intNodeValues;
1257: if (!(sp->intDofValues)) {
1258: MatCreateVecs(intMat, NULL, &(sp->intDofValues));
1259: }
1260: dofValues = sp->intDofValues;
1261: VecPlaceArray(pointValues, pointEval);
1262: VecPlaceArray(dofValues, spValue);
1263: MatMult(intMat, pointValues, dofValues);
1264: VecResetArray(dofValues);
1265: VecResetArray(pointValues);
1266: return(0);
1267: }
1269: /*@
1270: PetscDualSpaceGetAllData - Get all quadrature nodes from this space, and the matrix that sends quadrature node values to degree-of-freedom values
1272: Input Parameter:
1273: . sp - The dualspace
1275: Output Parameter:
1276: + allNodes - A PetscQuadrature object containing all evaluation nodes
1277: - allMat - A Mat for the node-to-dof transformation
1279: Level: advanced
1281: .seealso: PetscDualSpaceCreate()
1282: @*/
1283: PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat)
1284: {
1291: if ((!sp->allNodes || !sp->allMat) && sp->ops->createalldata) {
1292: PetscQuadrature qpoints;
1293: Mat amat;
1295: (*sp->ops->createalldata)(sp,&qpoints,&amat);
1296: PetscQuadratureDestroy(&(sp->allNodes));
1297: MatDestroy(&(sp->allMat));
1298: sp->allNodes = qpoints;
1299: sp->allMat = amat;
1300: }
1301: if (allNodes) *allNodes = sp->allNodes;
1302: if (allMat) *allMat = sp->allMat;
1303: return(0);
1304: }
1306: /*@
1307: PetscDualSpaceCreateAllDataDefault - Create all evaluation nodes and the node-to-dof matrix by examining functionals
1309: Input Parameter:
1310: . sp - The dualspace
1312: Output Parameter:
1313: + allNodes - A PetscQuadrature object containing all evaluation nodes
1314: - allMat - A Mat for the node-to-dof transformation
1316: Level: advanced
1318: .seealso: PetscDualSpaceCreate()
1319: @*/
1320: PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat)
1321: {
1322: PetscInt spdim;
1323: PetscInt numPoints, offset;
1324: PetscReal *points;
1325: PetscInt f, dim;
1326: PetscInt Nc, nrows, ncols;
1327: PetscInt maxNumPoints;
1328: PetscQuadrature q;
1329: Mat A;
1330: PetscErrorCode ierr;
1333: PetscDualSpaceGetNumComponents(sp, &Nc);
1334: PetscDualSpaceGetDimension(sp,&spdim);
1335: if (!spdim) {
1336: PetscQuadratureCreate(PETSC_COMM_SELF,allNodes);
1337: PetscQuadratureSetData(*allNodes,0,0,0,NULL,NULL);
1338: }
1339: nrows = spdim;
1340: PetscDualSpaceGetFunctional(sp,0,&q);
1341: PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);
1342: maxNumPoints = numPoints;
1343: for (f = 1; f < spdim; f++) {
1344: PetscInt Np;
1346: PetscDualSpaceGetFunctional(sp,f,&q);
1347: PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);
1348: numPoints += Np;
1349: maxNumPoints = PetscMax(maxNumPoints,Np);
1350: }
1351: ncols = numPoints * Nc;
1352: PetscMalloc1(dim*numPoints,&points);
1353: MatCreateSeqAIJ(PETSC_COMM_SELF, nrows, ncols, maxNumPoints * Nc, NULL, &A);
1354: for (f = 0, offset = 0; f < spdim; f++) {
1355: const PetscReal *p, *w;
1356: PetscInt Np, i;
1357: PetscInt fnc;
1359: PetscDualSpaceGetFunctional(sp,f,&q);
1360: PetscQuadratureGetData(q,NULL,&fnc,&Np,&p,&w);
1361: if (fnc != Nc) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "functional component mismatch");
1362: for (i = 0; i < Np * dim; i++) {
1363: points[offset* dim + i] = p[i];
1364: }
1365: for (i = 0; i < Np * Nc; i++) {
1366: MatSetValue(A, f, offset * Nc, w[i], INSERT_VALUES);
1367: }
1368: offset += Np;
1369: }
1370: MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
1371: MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
1372: PetscQuadratureCreate(PETSC_COMM_SELF,allNodes);
1373: PetscQuadratureSetData(*allNodes,dim,0,numPoints,points,NULL);
1374: *allMat = A;
1375: return(0);
1376: }
1378: /*@
1379: PetscDualSpaceGetInteriorData - Get all quadrature points necessary to compute the interior degrees of freedom from
1380: this space, as well as the matrix that computes the degrees of freedom from the quadrature values. Degrees of
1381: freedom are interior degrees of freedom if they belong (by PetscDualSpaceGetSection()) to interior points in the
1382: reference DMPlex: complementary boundary degrees of freedom are marked as constrained in the section returned by
1383: PetscDualSpaceGetSection()).
1385: Input Parameter:
1386: . sp - The dualspace
1388: Output Parameter:
1389: + intNodes - A PetscQuadrature object containing all evaluation points needed to evaluate interior degrees of freedom
1390: - intMat - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is
1391: the size of the constrained layout (PetscSectionGetConstrainStorageSize()) of the dual space section,
1392: npoints is the number of points in intNodes and nc is PetscDualSpaceGetNumComponents().
1394: Level: advanced
1396: .seealso: PetscDualSpaceCreate(), PetscDualSpaceGetDimension(), PetscDualSpaceGetNumComponents(), PetscQuadratureGetData()
1397: @*/
1398: PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat)
1399: {
1406: if ((!sp->intNodes || !sp->intMat) && sp->ops->createintdata) {
1407: PetscQuadrature qpoints;
1408: Mat imat;
1410: (*sp->ops->createintdata)(sp,&qpoints,&imat);
1411: PetscQuadratureDestroy(&(sp->intNodes));
1412: MatDestroy(&(sp->intMat));
1413: sp->intNodes = qpoints;
1414: sp->intMat = imat;
1415: }
1416: if (intNodes) *intNodes = sp->intNodes;
1417: if (intMat) *intMat = sp->intMat;
1418: return(0);
1419: }
1421: /*@
1422: PetscDualSpaceCreateInteriorDataDefault - Create quadrature points by examining interior functionals and create the matrix mapping quadrature point values to interior dual space values
1424: Input Parameter:
1425: . sp - The dualspace
1427: Output Parameter:
1428: + intNodes - A PetscQuadrature object containing all evaluation points needed to evaluate interior degrees of freedom
1429: - intMat - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is
1430: the size of the constrained layout (PetscSectionGetConstrainStorageSize()) of the dual space section,
1431: npoints is the number of points in allNodes and nc is PetscDualSpaceGetNumComponents().
1433: Level: advanced
1435: .seealso: PetscDualSpaceCreate(), PetscDualSpaceGetInteriorData()
1436: @*/
1437: PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat)
1438: {
1439: DM dm;
1440: PetscInt spdim0;
1441: PetscInt Nc;
1442: PetscInt pStart, pEnd, p, f;
1443: PetscSection section;
1444: PetscInt numPoints, offset, matoffset;
1445: PetscReal *points;
1446: PetscInt dim;
1447: PetscInt *nnz;
1448: PetscQuadrature q;
1449: Mat imat;
1450: PetscErrorCode ierr;
1454: PetscDualSpaceGetSection(sp, §ion);
1455: PetscSectionGetConstrainedStorageSize(section, &spdim0);
1456: if (!spdim0) {
1457: *intNodes = NULL;
1458: *intMat = NULL;
1459: return(0);
1460: }
1461: PetscDualSpaceGetNumComponents(sp, &Nc);
1462: PetscSectionGetChart(section, &pStart, &pEnd);
1463: PetscDualSpaceGetDM(sp, &dm);
1464: DMGetDimension(dm, &dim);
1465: PetscMalloc1(spdim0, &nnz);
1466: for (p = pStart, f = 0, numPoints = 0; p < pEnd; p++) {
1467: PetscInt dof, cdof, off, d;
1469: PetscSectionGetDof(section, p, &dof);
1470: PetscSectionGetConstraintDof(section, p, &cdof);
1471: if (!(dof - cdof)) continue;
1472: PetscSectionGetOffset(section, p, &off);
1473: for (d = 0; d < dof; d++, off++, f++) {
1474: PetscInt Np;
1476: PetscDualSpaceGetFunctional(sp,off,&q);
1477: PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);
1478: nnz[f] = Np * Nc;
1479: numPoints += Np;
1480: }
1481: }
1482: MatCreateSeqAIJ(PETSC_COMM_SELF, spdim0, numPoints * Nc, 0, nnz, &imat);
1483: PetscFree(nnz);
1484: PetscMalloc1(dim*numPoints,&points);
1485: for (p = pStart, f = 0, offset = 0, matoffset = 0; p < pEnd; p++) {
1486: PetscInt dof, cdof, off, d;
1488: PetscSectionGetDof(section, p, &dof);
1489: PetscSectionGetConstraintDof(section, p, &cdof);
1490: if (!(dof - cdof)) continue;
1491: PetscSectionGetOffset(section, p, &off);
1492: for (d = 0; d < dof; d++, off++, f++) {
1493: const PetscReal *p;
1494: const PetscReal *w;
1495: PetscInt Np, i;
1497: PetscDualSpaceGetFunctional(sp,off,&q);
1498: PetscQuadratureGetData(q,NULL,NULL,&Np,&p,&w);
1499: for (i = 0; i < Np * dim; i++) {
1500: points[offset + i] = p[i];
1501: }
1502: for (i = 0; i < Np * Nc; i++) {
1503: MatSetValue(imat, f, matoffset + i, w[i],INSERT_VALUES);
1504: }
1505: offset += Np * dim;
1506: matoffset += Np * Nc;
1507: }
1508: }
1509: PetscQuadratureCreate(PETSC_COMM_SELF,intNodes);
1510: PetscQuadratureSetData(*intNodes,dim,0,numPoints,points,NULL);
1511: MatAssemblyBegin(imat, MAT_FINAL_ASSEMBLY);
1512: MatAssemblyEnd(imat, MAT_FINAL_ASSEMBLY);
1513: *intMat = imat;
1514: return(0);
1515: }
1517: /*@C
1518: PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
1520: Input Parameters:
1521: + sp - The PetscDualSpace object
1522: . f - The basis functional index
1523: . time - The time
1524: . cgeom - A context with geometric information for this cell, we currently just use the centroid
1525: . Nc - The number of components for the function
1526: . func - The input function
1527: - ctx - A context for the function
1529: Output Parameter:
1530: . value - The output value (scalar)
1532: Note: The calling sequence for the callback func is given by:
1534: $ func(PetscInt dim, PetscReal time, const PetscReal x[],
1535: $ PetscInt numComponents, PetscScalar values[], void *ctx)
1537: and the idea is to evaluate the functional as an integral
1539: $ n(f) = int dx n(x) . f(x)
1541: where both n and f have Nc components.
1543: Level: beginner
1545: .seealso: PetscDualSpaceCreate()
1546: @*/
1547: 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)
1548: {
1549: DM dm;
1550: PetscQuadrature n;
1551: const PetscReal *points, *weights;
1552: PetscScalar *val;
1553: PetscInt dimEmbed, qNc, c, Nq, q;
1554: PetscErrorCode ierr;
1559: PetscDualSpaceGetDM(sp, &dm);
1560: DMGetCoordinateDim(dm, &dimEmbed);
1561: PetscDualSpaceGetFunctional(sp, f, &n);
1562: PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);
1563: if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
1564: DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);
1565: *value = 0.;
1566: for (q = 0; q < Nq; ++q) {
1567: (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);
1568: for (c = 0; c < Nc; ++c) {
1569: *value += val[c]*weights[q*Nc+c];
1570: }
1571: }
1572: DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);
1573: return(0);
1574: }
1576: /*@
1577: PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
1578: given height. This assumes that the reference cell is symmetric over points of this height.
1580: If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
1581: pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not
1582: support extracting subspaces, then NULL is returned.
1584: This does not increment the reference count on the returned dual space, and the user should not destroy it.
1586: Not collective
1588: Input Parameters:
1589: + sp - the PetscDualSpace object
1590: - height - the height of the mesh point for which the subspace is desired
1592: Output Parameter:
1593: . subsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the
1594: point, which will be of lesser dimension if height > 0.
1596: Level: advanced
1598: .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace
1599: @*/
1600: PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
1601: {
1602: PetscInt depth = -1, cStart, cEnd;
1603: DM dm;
1609: if (!(sp->uniform)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A non-uniform dual space does not have a single dual space at each height");
1610: *subsp = NULL;
1611: dm = sp->dm;
1612: DMPlexGetDepth(dm, &depth);
1613: if (height < 0 || height > depth) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid height");
1614: DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);
1615: if (height == 0 && cEnd == cStart + 1) {
1616: *subsp = sp;
1617: return(0);
1618: }
1619: if (!sp->heightSpaces) {
1620: PetscInt h;
1621: PetscCalloc1(depth+1, &(sp->heightSpaces));
1623: for (h = 0; h <= depth; h++) {
1624: if (h == 0 && cEnd == cStart + 1) continue;
1625: if (sp->ops->createheightsubspace) {(*sp->ops->createheightsubspace)(sp,height,&(sp->heightSpaces[h]));}
1626: else if (sp->pointSpaces) {
1627: PetscInt hStart, hEnd;
1629: DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);
1630: if (hEnd > hStart) {
1631: const char *name;
1633: PetscObjectReference((PetscObject)(sp->pointSpaces[hStart]));
1634: if (sp->pointSpaces[hStart]) {
1635: PetscObjectGetName((PetscObject) sp, &name);
1636: PetscObjectSetName((PetscObject) sp->pointSpaces[hStart], name);
1637: }
1638: sp->heightSpaces[h] = sp->pointSpaces[hStart];
1639: }
1640: }
1641: }
1642: }
1643: *subsp = sp->heightSpaces[height];
1644: return(0);
1645: }
1647: /*@
1648: PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.
1650: If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
1651: defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting
1652: subspaces, then NULL is returned.
1654: This does not increment the reference count on the returned dual space, and the user should not destroy it.
1656: Not collective
1658: Input Parameters:
1659: + sp - the PetscDualSpace object
1660: - point - the point (in the dual space's DM) for which the subspace is desired
1662: Output Parameters:
1663: bdsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the
1664: point, which will be of lesser dimension if height > 0.
1666: Level: advanced
1668: .seealso: PetscDualSpace
1669: @*/
1670: PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
1671: {
1672: PetscInt pStart = 0, pEnd = 0, cStart, cEnd;
1673: DM dm;
1679: *bdsp = NULL;
1680: dm = sp->dm;
1681: DMPlexGetChart(dm, &pStart, &pEnd);
1682: if (point < pStart || point > pEnd) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point");
1683: DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);
1684: 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 */
1685: *bdsp = sp;
1686: return(0);
1687: }
1688: if (!sp->pointSpaces) {
1689: PetscInt p;
1690: PetscCalloc1(pEnd - pStart, &(sp->pointSpaces));
1692: for (p = 0; p < pEnd - pStart; p++) {
1693: if (p + pStart == cStart && cEnd == cStart + 1) continue;
1694: if (sp->ops->createpointsubspace) {(*sp->ops->createpointsubspace)(sp,p+pStart,&(sp->pointSpaces[p]));}
1695: else if (sp->heightSpaces || sp->ops->createheightsubspace) {
1696: PetscInt dim, depth, height;
1697: DMLabel label;
1699: DMPlexGetDepth(dm,&dim);
1700: DMPlexGetDepthLabel(dm,&label);
1701: DMLabelGetValue(label,p+pStart,&depth);
1702: height = dim - depth;
1703: PetscDualSpaceGetHeightSubspace(sp, height, &(sp->pointSpaces[p]));
1704: PetscObjectReference((PetscObject)sp->pointSpaces[p]);
1705: }
1706: }
1707: }
1708: *bdsp = sp->pointSpaces[point - pStart];
1709: return(0);
1710: }
1712: /*@C
1713: PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis
1715: Not collective
1717: Input Parameter:
1718: . sp - the PetscDualSpace object
1720: Output Parameters:
1721: + perms - Permutations of the interior degrees of freedom, parameterized by the point orientation
1722: - flips - Sign reversal of the interior degrees of freedom, parameterized by the point orientation
1724: Note: The permutation and flip arrays are organized in the following way
1725: $ perms[p][ornt][dof # on point] = new local dof #
1726: $ flips[p][ornt][dof # on point] = reversal or not
1728: Level: developer
1730: @*/
1731: PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
1732: {
1739: if (sp->ops->getsymmetries) {(sp->ops->getsymmetries)(sp,perms,flips);}
1740: return(0);
1741: }
1743: /*@
1744: PetscDualSpaceGetFormDegree - Get the form degree k for the k-form the describes the pushforwards/pullbacks of this
1745: dual space's functionals.
1747: Input Parameter:
1748: . dsp - The PetscDualSpace
1750: Output Parameter:
1751: . k - The *signed* degree k of the k. If k >= 0, this means that the degrees of freedom are k-forms, and are stored
1752: in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms. So for example,
1753: 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).
1754: If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the
1755: Hodge star map. So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms
1756: but are stored as 1-forms.
1758: Level: developer
1760: .seealso: PetscDTAltV, PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType
1761: @*/
1762: PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace dsp, PetscInt *k)
1763: {
1767: *k = dsp->k;
1768: return(0);
1769: }
1771: /*@
1772: PetscDualSpaceSetFormDegree - Set the form degree k for the k-form the describes the pushforwards/pullbacks of this
1773: dual space's functionals.
1775: Input Parameter:
1776: + dsp - The PetscDualSpace
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: PetscDTAltV, PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType
1787: @*/
1788: PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace dsp, PetscInt k)
1789: {
1790: PetscInt dim;
1794: if (dsp->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up");
1795: dim = dsp->dm->dim;
1796: if (k < -dim || k > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported %D-form on %D-dimensional reference cell", PetscAbsInt(k), dim);
1797: dsp->k = k;
1798: return(0);
1799: }
1801: /*@
1802: PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space
1804: Input Parameter:
1805: . dsp - The PetscDualSpace
1807: Output Parameter:
1808: . k - The simplex dimension
1810: Level: developer
1812: Note: Currently supported values are
1813: $ 0: These are H_1 methods that only transform coordinates
1814: $ 1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM)
1815: $ 2: These are the same as 1
1816: $ 3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM)
1818: .seealso: PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType
1819: @*/
1820: PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k)
1821: {
1822: PetscInt dim;
1827: dim = dsp->dm->dim;
1828: if (!dsp->k) *k = IDENTITY_TRANSFORM;
1829: else if (dsp->k == 1) *k = COVARIANT_PIOLA_TRANSFORM;
1830: else if (dsp->k == -(dim - 1)) *k = CONTRAVARIANT_PIOLA_TRANSFORM;
1831: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported transformation");
1832: return(0);
1833: }
1835: /*@C
1836: PetscDualSpaceTransform - Transform the function values
1838: Input Parameters:
1839: + dsp - The PetscDualSpace
1840: . trans - The type of transform
1841: . isInverse - Flag to invert the transform
1842: . fegeom - The cell geometry
1843: . Nv - The number of function samples
1844: . Nc - The number of function components
1845: - vals - The function values
1847: Output Parameter:
1848: . vals - The transformed function values
1850: Level: intermediate
1852: Note: This only handles tranformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
1854: .seealso: PetscDualSpaceTransformGradient(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType
1855: @*/
1856: PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1857: {
1858: PetscReal Jstar[9] = {0};
1859: PetscInt dim, v, c, Nk;
1866: /* TODO: not handling dimEmbed != dim right now */
1867: dim = dsp->dm->dim;
1868: /* No change needed for 0-forms */
1869: if (!dsp->k) return(0);
1870: PetscDTBinomialInt(dim, PetscAbsInt(dsp->k), &Nk);
1871: /* TODO: use fegeom->isAffine */
1872: PetscDTAltVPullbackMatrix(dim, dim, isInverse ? fegeom->J : fegeom->invJ, dsp->k, Jstar);
1873: for (v = 0; v < Nv; ++v) {
1874: switch (Nk) {
1875: case 1:
1876: for (c = 0; c < Nc; c++) vals[v*Nc + c] *= Jstar[0];
1877: break;
1878: case 2:
1879: for (c = 0; c < Nc; c += 2) DMPlex_Mult2DReal_Internal(Jstar, 1, &vals[v*Nc + c], &vals[v*Nc + c]);
1880: break;
1881: case 3:
1882: for (c = 0; c < Nc; c += 3) DMPlex_Mult3DReal_Internal(Jstar, 1, &vals[v*Nc + c], &vals[v*Nc + c]);
1883: break;
1884: default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported form size %D for transformation", Nk);
1885: }
1886: }
1887: return(0);
1888: }
1890: /*@C
1891: PetscDualSpaceTransformGradient - Transform the function gradient values
1893: Input Parameters:
1894: + dsp - The PetscDualSpace
1895: . trans - The type of transform
1896: . isInverse - Flag to invert the transform
1897: . fegeom - The cell geometry
1898: . Nv - The number of function gradient samples
1899: . Nc - The number of function components
1900: - vals - The function gradient values
1902: Output Parameter:
1903: . vals - The transformed function values
1905: Level: intermediate
1907: Note: This only handles tranformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
1909: .seealso: PetscDualSpaceTransform(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType
1910: @*/
1911: PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1912: {
1913: const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed;
1914: PetscInt v, c, d;
1920: #ifdef PETSC_USE_DEBUG
1921: if (dE <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %D", dE);
1922: #endif
1923: /* Transform gradient */
1924: if (dim == dE) {
1925: for (v = 0; v < Nv; ++v) {
1926: for (c = 0; c < Nc; ++c) {
1927: switch (dim)
1928: {
1929: case 1: vals[(v*Nc+c)*dim] *= fegeom->invJ[0];break;
1930: case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break;
1931: case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break;
1932: default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1933: }
1934: }
1935: }
1936: } else {
1937: for (v = 0; v < Nv; ++v) {
1938: for (c = 0; c < Nc; ++c) {
1939: DMPlex_MultTransposeReal_Internal(fegeom->invJ, dim, dE, 1, &vals[(v*Nc+c)*dE], &vals[(v*Nc+c)*dE]);
1940: }
1941: }
1942: }
1943: /* Assume its a vector, otherwise assume its a bunch of scalars */
1944: if (Nc == 1 || Nc != dim) return(0);
1945: switch (trans) {
1946: case IDENTITY_TRANSFORM: break;
1947: case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
1948: if (isInverse) {
1949: for (v = 0; v < Nv; ++v) {
1950: for (d = 0; d < dim; ++d) {
1951: switch (dim)
1952: {
1953: case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1954: case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1955: default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1956: }
1957: }
1958: }
1959: } else {
1960: for (v = 0; v < Nv; ++v) {
1961: for (d = 0; d < dim; ++d) {
1962: switch (dim)
1963: {
1964: case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1965: case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1966: default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1967: }
1968: }
1969: }
1970: }
1971: break;
1972: case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
1973: if (isInverse) {
1974: for (v = 0; v < Nv; ++v) {
1975: for (d = 0; d < dim; ++d) {
1976: switch (dim)
1977: {
1978: case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1979: case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1980: default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1981: }
1982: for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] *= fegeom->detJ[0];
1983: }
1984: }
1985: } else {
1986: for (v = 0; v < Nv; ++v) {
1987: for (d = 0; d < dim; ++d) {
1988: switch (dim)
1989: {
1990: case 2: DMPlex_Mult2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1991: case 3: DMPlex_Mult3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1992: default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1993: }
1994: for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] /= fegeom->detJ[0];
1995: }
1996: }
1997: }
1998: break;
1999: }
2000: return(0);
2001: }
2003: /*@C
2004: 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.
2006: Input Parameters:
2007: + dsp - The PetscDualSpace
2008: . fegeom - The geometry for this cell
2009: . Nq - The number of function samples
2010: . Nc - The number of function components
2011: - pointEval - The function values
2013: Output Parameter:
2014: . pointEval - The transformed function values
2016: Level: advanced
2018: Note: 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.
2020: Note: This only handles tranformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2022: .seealso: PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
2023: @*/
2024: PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2025: {
2026: PetscDualSpaceTransformType trans;
2027: PetscInt k;
2028: PetscErrorCode ierr;
2034: /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2035: This determines their transformation properties. */
2036: PetscDualSpaceGetDeRahm(dsp, &k);
2037: switch (k)
2038: {
2039: case 0: /* H^1 point evaluations */
2040: trans = IDENTITY_TRANSFORM;break;
2041: case 1: /* Hcurl preserves tangential edge traces */
2042: trans = COVARIANT_PIOLA_TRANSFORM;break;
2043: case 2:
2044: case 3: /* Hdiv preserve normal traces */
2045: trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
2046: default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k);
2047: }
2048: PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, pointEval);
2049: return(0);
2050: }
2052: /*@C
2053: 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.
2055: Input Parameters:
2056: + dsp - The PetscDualSpace
2057: . fegeom - The geometry for this cell
2058: . Nq - The number of function samples
2059: . Nc - The number of function components
2060: - pointEval - The function values
2062: Output Parameter:
2063: . pointEval - The transformed function values
2065: Level: advanced
2067: Note: 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.
2069: Note: This only handles tranformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2071: .seealso: PetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
2072: @*/
2073: PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2074: {
2075: PetscDualSpaceTransformType trans;
2076: PetscInt k;
2077: PetscErrorCode ierr;
2083: /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2084: This determines their transformation properties. */
2085: PetscDualSpaceGetDeRahm(dsp, &k);
2086: switch (k)
2087: {
2088: case 0: /* H^1 point evaluations */
2089: trans = IDENTITY_TRANSFORM;break;
2090: case 1: /* Hcurl preserves tangential edge traces */
2091: trans = COVARIANT_PIOLA_TRANSFORM;break;
2092: case 2:
2093: case 3: /* Hdiv preserve normal traces */
2094: trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
2095: default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k);
2096: }
2097: PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);
2098: return(0);
2099: }
2101: /*@C
2102: 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.
2104: Input Parameters:
2105: + dsp - The PetscDualSpace
2106: . fegeom - The geometry for this cell
2107: . Nq - The number of function gradient samples
2108: . Nc - The number of function components
2109: - pointEval - The function gradient values
2111: Output Parameter:
2112: . pointEval - The transformed function gradient values
2114: Level: advanced
2116: Note: 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.
2118: Note: This only handles tranformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2120: .seealso: PetscDualSpacePushforward(), PPetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
2121: @*/
2122: PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2123: {
2124: PetscDualSpaceTransformType trans;
2125: PetscInt k;
2126: PetscErrorCode ierr;
2132: /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2133: This determines their transformation properties. */
2134: PetscDualSpaceGetDeRahm(dsp, &k);
2135: switch (k)
2136: {
2137: case 0: /* H^1 point evaluations */
2138: trans = IDENTITY_TRANSFORM;break;
2139: case 1: /* Hcurl preserves tangential edge traces */
2140: trans = COVARIANT_PIOLA_TRANSFORM;break;
2141: case 2:
2142: case 3: /* Hdiv preserve normal traces */
2143: trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
2144: default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k);
2145: }
2146: PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);
2147: return(0);
2148: }