Actual source code: space.c
1: #include <petsc/private/petscfeimpl.h>
2: #include <petscdmshell.h>
4: PetscClassId PETSCSPACE_CLASSID = 0;
6: PetscFunctionList PetscSpaceList = NULL;
7: PetscBool PetscSpaceRegisterAllCalled = PETSC_FALSE;
9: /*@C
10: PetscSpaceRegister - Adds a new `PetscSpace` implementation
12: Not Collective
14: Input Parameters:
15: + sname - The name of a new user-defined creation routine
16: - function - The creation routine for the implementation type
18: Example Usage:
19: .vb
20: PetscSpaceRegister("my_space", MyPetscSpaceCreate);
21: .ve
23: Then, your PetscSpace type can be chosen with the procedural interface via
24: .vb
25: PetscSpaceCreate(MPI_Comm, PetscSpace *);
26: PetscSpaceSetType(PetscSpace, "my_space");
27: .ve
28: or at runtime via the option
29: .vb
30: -petscspace_type my_space
31: .ve
33: Level: advanced
35: Note:
36: `PetscSpaceRegister()` may be called multiple times to add several user-defined types of `PetscSpace`. The creation function is called
37: when the type is set to 'name'.
39: .seealso: `PetscSpace`, `PetscSpaceRegisterAll()`, `PetscSpaceRegisterDestroy()`
40: @*/
41: PetscErrorCode PetscSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscSpace))
42: {
43: PetscFunctionBegin;
44: PetscCall(PetscFunctionListAdd(&PetscSpaceList, sname, function));
45: PetscFunctionReturn(PETSC_SUCCESS);
46: }
48: /*@C
49: PetscSpaceSetType - Builds a particular `PetscSpace`
51: Collective
53: Input Parameters:
54: + sp - The `PetscSpace` object
55: - name - The kind of space
57: Options Database Key:
58: . -petscspace_type <type> - Sets the `PetscSpace` type; use -help for a list of available types
60: Level: intermediate
62: .seealso: `PetscSpace`, `PetscSpaceType`, `PetscSpaceGetType()`, `PetscSpaceCreate()`
63: @*/
64: PetscErrorCode PetscSpaceSetType(PetscSpace sp, PetscSpaceType name)
65: {
66: PetscErrorCode (*r)(PetscSpace);
67: PetscBool match;
69: PetscFunctionBegin;
71: PetscCall(PetscObjectTypeCompare((PetscObject)sp, name, &match));
72: if (match) PetscFunctionReturn(PETSC_SUCCESS);
74: PetscCall(PetscSpaceRegisterAll());
75: PetscCall(PetscFunctionListFind(PetscSpaceList, name, &r));
76: PetscCheck(r, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSpace type: %s", name);
78: PetscTryTypeMethod(sp, destroy);
79: sp->ops->destroy = NULL;
81: sp->dim = PETSC_DETERMINE;
82: PetscCall((*r)(sp));
83: PetscCall(PetscObjectChangeTypeName((PetscObject)sp, name));
84: PetscFunctionReturn(PETSC_SUCCESS);
85: }
87: /*@C
88: PetscSpaceGetType - Gets the `PetscSpaceType` (as a string) from the object.
90: Not Collective
92: Input Parameter:
93: . sp - The `PetscSpace`
95: Output Parameter:
96: . name - The `PetscSpace` type name
98: Level: intermediate
100: .seealso: `PetscSpaceType`, `PetscSpace`, `PetscSpaceSetType()`, `PetscSpaceCreate()`
101: @*/
102: PetscErrorCode PetscSpaceGetType(PetscSpace sp, PetscSpaceType *name)
103: {
104: PetscFunctionBegin;
106: PetscAssertPointer(name, 2);
107: if (!PetscSpaceRegisterAllCalled) PetscCall(PetscSpaceRegisterAll());
108: *name = ((PetscObject)sp)->type_name;
109: PetscFunctionReturn(PETSC_SUCCESS);
110: }
112: /*@C
113: PetscSpaceViewFromOptions - View a `PetscSpace` based on values in the options database
115: Collective
117: Input Parameters:
118: + A - the `PetscSpace` object
119: . obj - Optional object that provides the options name prefix
120: - name - command line option name
122: Level: intermediate
124: .seealso: `PetscSpace`, `PetscSpaceView()`, `PetscObjectViewFromOptions()`, `PetscSpaceCreate()`
125: @*/
126: PetscErrorCode PetscSpaceViewFromOptions(PetscSpace A, PetscObject obj, const char name[])
127: {
128: PetscFunctionBegin;
130: PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
131: PetscFunctionReturn(PETSC_SUCCESS);
132: }
134: /*@C
135: PetscSpaceView - Views a `PetscSpace`
137: Collective
139: Input Parameters:
140: + sp - the `PetscSpace` object to view
141: - v - the viewer
143: Level: beginner
145: .seealso: `PetscSpace`, `PetscViewer`, `PetscSpaceViewFromOptions()`, `PetscSpaceDestroy()`
146: @*/
147: PetscErrorCode PetscSpaceView(PetscSpace sp, PetscViewer v)
148: {
149: PetscInt pdim;
150: PetscBool iascii;
152: PetscFunctionBegin;
155: if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sp), &v));
156: PetscCall(PetscSpaceGetDimension(sp, &pdim));
157: PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)sp, v));
158: PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii));
159: PetscCall(PetscViewerASCIIPushTab(v));
160: if (iascii) PetscCall(PetscViewerASCIIPrintf(v, "Space in %" PetscInt_FMT " variables with %" PetscInt_FMT " components, size %" PetscInt_FMT "\n", sp->Nv, sp->Nc, pdim));
161: PetscTryTypeMethod(sp, view, v);
162: PetscCall(PetscViewerASCIIPopTab(v));
163: PetscFunctionReturn(PETSC_SUCCESS);
164: }
166: /*@
167: PetscSpaceSetFromOptions - sets parameters in a `PetscSpace` from the options database
169: Collective
171: Input Parameter:
172: . sp - the `PetscSpace` object to set options for
174: Options Database Keys:
175: + -petscspace_degree <deg> - the approximation order of the space
176: . -petscspace_variables <n> - the number of different variables, e.g. x and y
177: - -petscspace_components <c> - the number of components, say d for a vector field
179: Level: intermediate
181: .seealso: `PetscSpace`, `PetscSpaceView()`
182: @*/
183: PetscErrorCode PetscSpaceSetFromOptions(PetscSpace sp)
184: {
185: const char *defaultType;
186: char name[256];
187: PetscBool flg;
189: PetscFunctionBegin;
191: if (!((PetscObject)sp)->type_name) {
192: defaultType = PETSCSPACEPOLYNOMIAL;
193: } else {
194: defaultType = ((PetscObject)sp)->type_name;
195: }
196: if (!PetscSpaceRegisterAllCalled) PetscCall(PetscSpaceRegisterAll());
198: PetscObjectOptionsBegin((PetscObject)sp);
199: PetscCall(PetscOptionsFList("-petscspace_type", "Linear space", "PetscSpaceSetType", PetscSpaceList, defaultType, name, 256, &flg));
200: if (flg) {
201: PetscCall(PetscSpaceSetType(sp, name));
202: } else if (!((PetscObject)sp)->type_name) {
203: PetscCall(PetscSpaceSetType(sp, defaultType));
204: }
205: {
206: PetscCall(PetscOptionsDeprecated("-petscspace_order", "-petscspace_degree", "3.11", NULL));
207: }
208: PetscCall(PetscOptionsBoundedInt("-petscspace_degree", "The (maximally included) polynomial degree", "PetscSpaceSetDegree", sp->degree, &sp->degree, NULL, 0));
209: PetscCall(PetscOptionsBoundedInt("-petscspace_variables", "The number of different variables, e.g. x and y", "PetscSpaceSetNumVariables", sp->Nv, &sp->Nv, NULL, 0));
210: PetscCall(PetscOptionsBoundedInt("-petscspace_components", "The number of components", "PetscSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL, 0));
211: PetscTryTypeMethod(sp, setfromoptions, PetscOptionsObject);
212: /* process any options handlers added with PetscObjectAddOptionsHandler() */
213: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)sp, PetscOptionsObject));
214: PetscOptionsEnd();
215: PetscCall(PetscSpaceViewFromOptions(sp, NULL, "-petscspace_view"));
216: PetscFunctionReturn(PETSC_SUCCESS);
217: }
219: /*@C
220: PetscSpaceSetUp - Construct data structures for the `PetscSpace`
222: Collective
224: Input Parameter:
225: . sp - the `PetscSpace` object to setup
227: Level: intermediate
229: .seealso: `PetscSpace`, `PetscSpaceView()`, `PetscSpaceDestroy()`
230: @*/
231: PetscErrorCode PetscSpaceSetUp(PetscSpace sp)
232: {
233: PetscFunctionBegin;
235: PetscTryTypeMethod(sp, setup);
236: PetscFunctionReturn(PETSC_SUCCESS);
237: }
239: /*@
240: PetscSpaceDestroy - Destroys a `PetscSpace` object
242: Collective
244: Input Parameter:
245: . sp - the `PetscSpace` object to destroy
247: Level: beginner
249: .seealso: `PetscSpace`, `PetscSpaceCreate()`
250: @*/
251: PetscErrorCode PetscSpaceDestroy(PetscSpace *sp)
252: {
253: PetscFunctionBegin;
254: if (!*sp) PetscFunctionReturn(PETSC_SUCCESS);
257: if (--((PetscObject)*sp)->refct > 0) {
258: *sp = NULL;
259: PetscFunctionReturn(PETSC_SUCCESS);
260: }
261: ((PetscObject)*sp)->refct = 0;
262: PetscCall(DMDestroy(&(*sp)->dm));
264: PetscUseTypeMethod(*sp, destroy);
265: PetscCall(PetscHeaderDestroy(sp));
266: PetscFunctionReturn(PETSC_SUCCESS);
267: }
269: /*@
270: PetscSpaceCreate - Creates an empty `PetscSpace` object. The type can then be set with `PetscSpaceSetType()`.
272: Collective
274: Input Parameter:
275: . comm - The communicator for the `PetscSpace` object
277: Output Parameter:
278: . sp - The `PetscSpace` object
280: Level: beginner
282: .seealso: `PetscSpace`, `PetscSpaceSetType()`, `PETSCSPACEPOLYNOMIAL`
283: @*/
284: PetscErrorCode PetscSpaceCreate(MPI_Comm comm, PetscSpace *sp)
285: {
286: PetscSpace s;
288: PetscFunctionBegin;
289: PetscAssertPointer(sp, 2);
290: PetscCall(PetscCitationsRegister(FECitation, &FEcite));
291: *sp = NULL;
292: PetscCall(PetscFEInitializePackage());
294: PetscCall(PetscHeaderCreate(s, PETSCSPACE_CLASSID, "PetscSpace", "Linear Space", "PetscSpace", comm, PetscSpaceDestroy, PetscSpaceView));
296: s->degree = 0;
297: s->maxDegree = PETSC_DETERMINE;
298: s->Nc = 1;
299: s->Nv = 0;
300: s->dim = PETSC_DETERMINE;
301: PetscCall(DMShellCreate(comm, &s->dm));
302: PetscCall(PetscSpaceSetType(s, PETSCSPACEPOLYNOMIAL));
304: *sp = s;
305: PetscFunctionReturn(PETSC_SUCCESS);
306: }
308: /*@
309: PetscSpaceGetDimension - Return the dimension of this space, i.e. the number of basis vectors
311: Input Parameter:
312: . sp - The `PetscSpace`
314: Output Parameter:
315: . dim - The dimension
317: Level: intermediate
319: .seealso: `PetscSpace`, `PetscSpaceGetDegree()`, `PetscSpaceCreate()`
320: @*/
321: PetscErrorCode PetscSpaceGetDimension(PetscSpace sp, PetscInt *dim)
322: {
323: PetscFunctionBegin;
325: PetscAssertPointer(dim, 2);
326: if (sp->dim == PETSC_DETERMINE) PetscTryTypeMethod(sp, getdimension, &sp->dim);
327: *dim = sp->dim;
328: PetscFunctionReturn(PETSC_SUCCESS);
329: }
331: /*@
332: PetscSpaceGetDegree - Return the polynomial degrees that characterize this space
334: Input Parameter:
335: . sp - The `PetscSpace`
337: Output Parameters:
338: + minDegree - The degree of the largest polynomial space contained in the space
339: - maxDegree - The degree of the smallest polynomial space containing the space
341: Level: intermediate
343: .seealso: `PetscSpace`, `PetscSpaceSetDegree()`, `PetscSpaceGetDimension()`, `PetscSpaceCreate()`
344: @*/
345: PetscErrorCode PetscSpaceGetDegree(PetscSpace sp, PetscInt *minDegree, PetscInt *maxDegree)
346: {
347: PetscFunctionBegin;
349: if (minDegree) PetscAssertPointer(minDegree, 2);
350: if (maxDegree) PetscAssertPointer(maxDegree, 3);
351: if (minDegree) *minDegree = sp->degree;
352: if (maxDegree) *maxDegree = sp->maxDegree;
353: PetscFunctionReturn(PETSC_SUCCESS);
354: }
356: /*@
357: PetscSpaceSetDegree - Set the degree of approximation for this space.
359: Input Parameters:
360: + sp - The `PetscSpace`
361: . degree - The degree of the largest polynomial space contained in the space
362: - maxDegree - The degree of the largest polynomial space containing the space. One of degree and maxDegree can be `PETSC_DETERMINE`.
364: Level: intermediate
366: .seealso: `PetscSpace`, `PetscSpaceGetDegree()`, `PetscSpaceCreate()`
367: @*/
368: PetscErrorCode PetscSpaceSetDegree(PetscSpace sp, PetscInt degree, PetscInt maxDegree)
369: {
370: PetscFunctionBegin;
372: sp->degree = degree;
373: sp->maxDegree = maxDegree;
374: PetscFunctionReturn(PETSC_SUCCESS);
375: }
377: /*@
378: PetscSpaceGetNumComponents - Return the number of components for this space
380: Input Parameter:
381: . sp - The `PetscSpace`
383: Output Parameter:
384: . Nc - The number of components
386: Level: intermediate
388: Note:
389: A vector space, for example, will have d components, where d is the spatial dimension
391: .seealso: `PetscSpace`, `PetscSpaceSetNumComponents()`, `PetscSpaceGetNumVariables()`, `PetscSpaceGetDimension()`, `PetscSpaceCreate()`
392: @*/
393: PetscErrorCode PetscSpaceGetNumComponents(PetscSpace sp, PetscInt *Nc)
394: {
395: PetscFunctionBegin;
397: PetscAssertPointer(Nc, 2);
398: *Nc = sp->Nc;
399: PetscFunctionReturn(PETSC_SUCCESS);
400: }
402: /*@
403: PetscSpaceSetNumComponents - Set the number of components for this space
405: Input Parameters:
406: + sp - The `PetscSpace`
407: - Nc - The number of components
409: Level: intermediate
411: .seealso: `PetscSpace`, `PetscSpaceGetNumComponents()`, `PetscSpaceSetNumVariables()`, `PetscSpaceCreate()`
412: @*/
413: PetscErrorCode PetscSpaceSetNumComponents(PetscSpace sp, PetscInt Nc)
414: {
415: PetscFunctionBegin;
417: sp->Nc = Nc;
418: PetscFunctionReturn(PETSC_SUCCESS);
419: }
421: /*@
422: PetscSpaceSetNumVariables - Set the number of variables for this space
424: Input Parameters:
425: + sp - The `PetscSpace`
426: - n - The number of variables, e.g. x, y, z...
428: Level: intermediate
430: .seealso: `PetscSpace`, `PetscSpaceGetNumVariables()`, `PetscSpaceSetNumComponents()`, `PetscSpaceCreate()`
431: @*/
432: PetscErrorCode PetscSpaceSetNumVariables(PetscSpace sp, PetscInt n)
433: {
434: PetscFunctionBegin;
436: sp->Nv = n;
437: PetscFunctionReturn(PETSC_SUCCESS);
438: }
440: /*@
441: PetscSpaceGetNumVariables - Return the number of variables for this space
443: Input Parameter:
444: . sp - The `PetscSpace`
446: Output Parameter:
447: . n - The number of variables, e.g. x, y, z...
449: Level: intermediate
451: .seealso: `PetscSpace`, `PetscSpaceSetNumVariables()`, `PetscSpaceGetNumComponents()`, `PetscSpaceGetDimension()`, `PetscSpaceCreate()`
452: @*/
453: PetscErrorCode PetscSpaceGetNumVariables(PetscSpace sp, PetscInt *n)
454: {
455: PetscFunctionBegin;
457: PetscAssertPointer(n, 2);
458: *n = sp->Nv;
459: PetscFunctionReturn(PETSC_SUCCESS);
460: }
462: /*@C
463: PetscSpaceEvaluate - Evaluate the basis functions and their derivatives (jet) at each point
465: Input Parameters:
466: + sp - The `PetscSpace`
467: . npoints - The number of evaluation points, in reference coordinates
468: - points - The point coordinates
470: Output Parameters:
471: + B - The function evaluations in a `npoints` x `nfuncs` array
472: . D - The derivative evaluations in a `npoints` x `nfuncs` x `dim` array
473: - H - The second derivative evaluations in a `npoints` x `nfuncs` x `dim` x `dim` array
475: Level: beginner
477: Note:
478: Above `nfuncs` is the dimension of the space, and `dim` is the spatial dimension. The coordinates are given
479: on the reference cell, not in real space.
481: .seealso: `PetscSpace`, `PetscFECreateTabulation()`, `PetscFEGetCellTabulation()`, `PetscSpaceCreate()`
482: @*/
483: PetscErrorCode PetscSpaceEvaluate(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
484: {
485: PetscFunctionBegin;
486: if (!npoints) PetscFunctionReturn(PETSC_SUCCESS);
488: if (sp->Nv) PetscAssertPointer(points, 3);
489: if (B) PetscAssertPointer(B, 4);
490: if (D) PetscAssertPointer(D, 5);
491: if (H) PetscAssertPointer(H, 6);
492: PetscTryTypeMethod(sp, evaluate, npoints, points, B, D, H);
493: PetscFunctionReturn(PETSC_SUCCESS);
494: }
496: /*@
497: PetscSpaceGetHeightSubspace - Get the subset of the primal space basis that is supported on a mesh point of a given height.
499: Not Collective
501: Input Parameters:
502: + sp - the `PetscSpace` object
503: - height - the height of the mesh point for which the subspace is desired
505: Output Parameter:
506: . subsp - the subspace
508: Level: advanced
510: Notes:
511: If the space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
512: pointwise values are not defined on the element boundaries), or if the implementation of `PetscSpace` does not
513: support extracting subspaces, then NULL is returned.
515: This does not increment the reference count on the returned space, and the user should not destroy it.
517: .seealso: `PetscDualSpaceGetHeightSubspace()`, `PetscSpace`
518: @*/
519: PetscErrorCode PetscSpaceGetHeightSubspace(PetscSpace sp, PetscInt height, PetscSpace *subsp)
520: {
521: PetscFunctionBegin;
523: PetscAssertPointer(subsp, 3);
524: *subsp = NULL;
525: PetscTryTypeMethod(sp, getheightsubspace, height, subsp);
526: PetscFunctionReturn(PETSC_SUCCESS);
527: }