Actual source code: stagutils.c
1: /* Additional functions in the DMStag API, which are not part of the general DM API. */
2: #include <petsc/private/dmstagimpl.h>
3: #include <petscdmproduct.h>
5: PetscErrorCode DMRestrictHook_Coordinates(DM, DM, void *);
7: /*@C
8: DMStagGetBoundaryTypes - get boundary types
10: Not Collective
12: Input Parameter:
13: . dm - the `DMSTAG` object
15: Output Parameters:
16: + boundaryTypeX - boundary type for x direction
17: . boundaryTypeY - boundary type for y direction, not set for one dimensional problems
18: - boundaryTypeZ - boundary type for z direction, not set for one and two dimensional problems
20: Level: intermediate
22: .seealso: [](ch_stag), `DMSTAG`, `DMBoundaryType`
23: @*/
24: PetscErrorCode DMStagGetBoundaryTypes(DM dm, DMBoundaryType *boundaryTypeX, DMBoundaryType *boundaryTypeY, DMBoundaryType *boundaryTypeZ)
25: {
26: const DM_Stag *const stag = (DM_Stag *)dm->data;
27: PetscInt dim;
29: PetscFunctionBegin;
31: PetscCall(DMGetDimension(dm, &dim));
32: if (boundaryTypeX) *boundaryTypeX = stag->boundaryType[0];
33: if (boundaryTypeY && dim > 1) *boundaryTypeY = stag->boundaryType[1];
34: if (boundaryTypeZ && dim > 2) *boundaryTypeZ = stag->boundaryType[2];
35: PetscFunctionReturn(PETSC_SUCCESS);
36: }
38: static PetscErrorCode DMStagGetProductCoordinateArrays_Private(DM dm, void *arrX, void *arrY, void *arrZ, PetscBool read)
39: {
40: PetscInt dim, d, dofCheck[DMSTAG_MAX_STRATA], s;
41: DM dmCoord;
42: void *arr[DMSTAG_MAX_DIM];
43: PetscBool checkDof;
45: PetscFunctionBegin;
47: PetscCall(DMGetDimension(dm, &dim));
48: PetscCheck(dim <= DMSTAG_MAX_DIM, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented for %" PetscInt_FMT " dimensions", dim);
49: arr[0] = arrX;
50: arr[1] = arrY;
51: arr[2] = arrZ;
52: PetscCall(DMGetCoordinateDM(dm, &dmCoord));
53: PetscCheck(dmCoord, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM does not have a coordinate DM");
54: {
55: PetscBool isProduct;
56: DMType dmType;
57: PetscCall(DMGetType(dmCoord, &dmType));
58: PetscCall(PetscStrcmp(DMPRODUCT, dmType, &isProduct));
59: PetscCheck(isProduct, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Coordinate DM is not of type DMPRODUCT");
60: }
61: for (s = 0; s < DMSTAG_MAX_STRATA; ++s) dofCheck[s] = 0;
62: checkDof = PETSC_FALSE;
63: for (d = 0; d < dim; ++d) {
64: DM subDM;
65: DMType dmType;
66: PetscBool isStag;
67: PetscInt dof[DMSTAG_MAX_STRATA], subDim;
68: Vec coord1d_local;
70: /* Ignore unrequested arrays */
71: if (!arr[d]) continue;
73: PetscCall(DMProductGetDM(dmCoord, d, &subDM));
74: PetscCheck(subDM, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Coordinate DM is missing sub DM %" PetscInt_FMT, d);
75: PetscCall(DMGetDimension(subDM, &subDim));
76: PetscCheck(subDim == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Coordinate sub-DM is not of dimension 1");
77: PetscCall(DMGetType(subDM, &dmType));
78: PetscCall(PetscStrcmp(DMSTAG, dmType, &isStag));
79: PetscCheck(isStag, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Coordinate sub-DM is not of type DMSTAG");
80: PetscCall(DMStagGetDOF(subDM, &dof[0], &dof[1], &dof[2], &dof[3]));
81: if (!checkDof) {
82: for (s = 0; s < DMSTAG_MAX_STRATA; ++s) dofCheck[s] = dof[s];
83: checkDof = PETSC_TRUE;
84: } else {
85: for (s = 0; s < DMSTAG_MAX_STRATA; ++s) PetscCheck(dofCheck[s] == dof[s], PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Coordinate sub-DMs have different dofs");
86: }
87: PetscCall(DMGetCoordinatesLocal(subDM, &coord1d_local));
88: if (read) {
89: PetscCall(DMStagVecGetArrayRead(subDM, coord1d_local, arr[d]));
90: } else {
91: PetscCall(DMStagVecGetArray(subDM, coord1d_local, arr[d]));
92: }
93: }
94: PetscFunctionReturn(PETSC_SUCCESS);
95: }
97: /*@C
98: DMStagGetProductCoordinateArrays - extract local product coordinate arrays, one per dimension
100: Logically Collective
102: Input Parameter:
103: . dm - the `DMSTAG` object
105: Output Parameters:
106: + arrX - local 1D coordinate arrays for x direction
107: . arrY - local 1D coordinate arrays for y direction, not set for one dimensional problems
108: - arrZ - local 1D coordinate arrays for z direction, not set for one and two dimensional problems
110: Level: intermediate
112: Notes:
113: A high-level helper function to quickly extract local coordinate arrays.
115: Note that 2-dimensional arrays are returned. See
116: `DMStagVecGetArray()`, which is called internally to produce these arrays
117: representing coordinates on elements and vertices (element boundaries)
118: for a 1-dimensional `DMSTAG` in each coordinate direction.
120: One should use `DMStagGetProductCoordinateLocationSlot()` to determine appropriate
121: indices for the second dimension in these returned arrays. This function
122: checks that the coordinate array is a suitable product of 1-dimensional
123: `DMSTAG` objects.
125: .seealso: [](ch_stag), `DMSTAG`, `DMPRODUCT`, `DMStagGetProductCoordinateArraysRead()`, `DMStagSetUniformCoordinates()`, `DMStagSetUniformCoordinatesProduct()`, `DMStagGetProductCoordinateLocationSlot()`
126: @*/
127: PetscErrorCode DMStagGetProductCoordinateArrays(DM dm, void *arrX, void *arrY, void *arrZ)
128: {
129: PetscFunctionBegin;
130: PetscCall(DMStagGetProductCoordinateArrays_Private(dm, arrX, arrY, arrZ, PETSC_FALSE));
131: PetscFunctionReturn(PETSC_SUCCESS);
132: }
134: /*@C
135: DMStagGetProductCoordinateArraysRead - extract product coordinate arrays, read-only
137: Logically Collective
139: Input Parameter:
140: . dm - the `DMSTAG` object
142: Output Parameters:
143: + arrX - local 1D coordinate arrays for x direction
144: . arrY - local 1D coordinate arrays for y direction, not set for one dimensional problems
145: - arrZ - local 1D coordinate arrays for z direction, not set for one and two dimensional problems
147: Level: intermediate
149: Note:
150: See `DMStagGetProductCoordinateArrays()` for more information.
152: .seealso: [](ch_stag), `DMSTAG`, `DMPRODUCT`, `DMStagGetProductCoordinateArrays()`, `DMStagSetUniformCoordinates()`, `DMStagSetUniformCoordinatesProduct()`, `DMStagGetProductCoordinateLocationSlot()`
153: @*/
154: PetscErrorCode DMStagGetProductCoordinateArraysRead(DM dm, void *arrX, void *arrY, void *arrZ)
155: {
156: PetscFunctionBegin;
157: PetscCall(DMStagGetProductCoordinateArrays_Private(dm, arrX, arrY, arrZ, PETSC_TRUE));
158: PetscFunctionReturn(PETSC_SUCCESS);
159: }
161: /*@C
162: DMStagGetProductCoordinateLocationSlot - get slot for use with local product coordinate arrays
164: Not Collective
166: Input Parameters:
167: + dm - the `DMSTAG` object
168: - loc - the grid location
170: Output Parameter:
171: . slot - the index to use in local arrays
173: Level: intermediate
175: Notes:
176: High-level helper function to get slot indices for 1D coordinate `DM`s,
177: for use with `DMStagGetProductCoordinateArrays()` and related functions.
179: For `loc`, one should use `DMSTAG_LEFT`, `DMSTAG_ELEMENT`, or `DMSTAG_RIGHT` for "previous", "center" and "next"
180: locations, respectively, in each dimension.
181: One can equivalently use `DMSTAG_DOWN` or `DMSTAG_BACK` in place of `DMSTAG_LEFT`,
182: and `DMSTAG_UP` or `DMSTACK_FRONT` in place of `DMSTAG_RIGHT`;
184: This function checks that the coordinates are actually set up so that using the
185: slots from any of the 1D coordinate sub-`DM`s are valid for all the 1D coordinate sub-`DM`s.
187: .seealso: [](ch_stag), `DMSTAG`, `DMPRODUCT`, `DMStagGetProductCoordinateArrays()`, `DMStagGetProductCoordinateArraysRead()`, `DMStagSetUniformCoordinates()`
188: @*/
189: PETSC_EXTERN PetscErrorCode DMStagGetProductCoordinateLocationSlot(DM dm, DMStagStencilLocation loc, PetscInt *slot)
190: {
191: DM dmCoord;
192: PetscInt dim, dofCheck[DMSTAG_MAX_STRATA], s, d;
194: PetscFunctionBegin;
196: PetscCall(DMGetDimension(dm, &dim));
197: PetscCall(DMGetCoordinateDM(dm, &dmCoord));
198: PetscCheck(dmCoord, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM does not have a coordinate DM");
199: {
200: PetscBool isProduct;
201: DMType dmType;
202: PetscCall(DMGetType(dmCoord, &dmType));
203: PetscCall(PetscStrcmp(DMPRODUCT, dmType, &isProduct));
204: PetscCheck(isProduct, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Coordinate DM is not of type DMPRODUCT");
205: }
206: for (s = 0; s < DMSTAG_MAX_STRATA; ++s) dofCheck[s] = 0;
207: for (d = 0; d < dim; ++d) {
208: DM subDM;
209: DMType dmType;
210: PetscBool isStag;
211: PetscInt dof[DMSTAG_MAX_STRATA], subDim;
212: PetscCall(DMProductGetDM(dmCoord, d, &subDM));
213: PetscCheck(subDM, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Coordinate DM is missing sub DM %" PetscInt_FMT, d);
214: PetscCall(DMGetDimension(subDM, &subDim));
215: PetscCheck(subDim == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Coordinate sub-DM is not of dimension 1");
216: PetscCall(DMGetType(subDM, &dmType));
217: PetscCall(PetscStrcmp(DMSTAG, dmType, &isStag));
218: PetscCheck(isStag, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Coordinate sub-DM is not of type DMSTAG");
219: PetscCall(DMStagGetDOF(subDM, &dof[0], &dof[1], &dof[2], &dof[3]));
220: if (d == 0) {
221: const PetscInt component = 0;
222: for (s = 0; s < DMSTAG_MAX_STRATA; ++s) dofCheck[s] = dof[s];
223: PetscCall(DMStagGetLocationSlot(subDM, loc, component, slot));
224: } else {
225: for (s = 0; s < DMSTAG_MAX_STRATA; ++s) PetscCheck(dofCheck[s] == dof[s], PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Coordinate sub-DMs have different dofs");
226: }
227: }
228: PetscFunctionReturn(PETSC_SUCCESS);
229: }
231: /*@C
232: DMStagGetCorners - return global element indices of the local region (excluding ghost points)
234: Not Collective
236: Input Parameter:
237: . dm - the `DMSTAG` object
239: Output Parameters:
240: + x - starting element index in first direction
241: . y - starting element index in second direction
242: . z - starting element index in third direction
243: . m - element width in first direction
244: . n - element width in second direction
245: . p - element width in third direction
246: . nExtrax - number of extra partial elements in first direction
247: . nExtray - number of extra partial elements in second direction
248: - nExtraz - number of extra partial elements in third direction
250: Level: beginner
252: Notes:
253: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
255: The number of extra partial elements is either 1 or 0.
256: The value is 1 on right, top, and front non-periodic domain ("physical") boundaries,
257: in the x, y, and z directions respectively, and otherwise 0.
259: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetGhostCorners()`, `DMDAGetCorners()`
260: @*/
261: PetscErrorCode DMStagGetCorners(DM dm, PetscInt *x, PetscInt *y, PetscInt *z, PetscInt *m, PetscInt *n, PetscInt *p, PetscInt *nExtrax, PetscInt *nExtray, PetscInt *nExtraz)
262: {
263: const DM_Stag *const stag = (DM_Stag *)dm->data;
265: PetscFunctionBegin;
267: if (x) *x = stag->start[0];
268: if (y) *y = stag->start[1];
269: if (z) *z = stag->start[2];
270: if (m) *m = stag->n[0];
271: if (n) *n = stag->n[1];
272: if (p) *p = stag->n[2];
273: if (nExtrax) *nExtrax = stag->boundaryType[0] != DM_BOUNDARY_PERIODIC && stag->lastRank[0] ? 1 : 0;
274: if (nExtray) *nExtray = stag->boundaryType[1] != DM_BOUNDARY_PERIODIC && stag->lastRank[1] ? 1 : 0;
275: if (nExtraz) *nExtraz = stag->boundaryType[2] != DM_BOUNDARY_PERIODIC && stag->lastRank[2] ? 1 : 0;
276: PetscFunctionReturn(PETSC_SUCCESS);
277: }
279: /*@C
280: DMStagGetDOF - get number of DOF associated with each stratum of the grid
282: Not Collective
284: Input Parameter:
285: . dm - the `DMSTAG` object
287: Output Parameters:
288: + dof0 - the number of points per 0-cell (vertex/node)
289: . dof1 - the number of points per 1-cell (element in 1D, edge in 2D and 3D)
290: . dof2 - the number of points per 2-cell (element in 2D, face in 3D)
291: - dof3 - the number of points per 3-cell (element in 3D)
293: Level: beginner
295: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetCorners()`, `DMStagGetGhostCorners()`, `DMStagGetGlobalSizes()`, `DMStagGetStencilWidth()`, `DMStagGetBoundaryTypes()`, `DMStagGetLocationDOF()`, `DMDAGetDof()`
296: @*/
297: PetscErrorCode DMStagGetDOF(DM dm, PetscInt *dof0, PetscInt *dof1, PetscInt *dof2, PetscInt *dof3)
298: {
299: const DM_Stag *const stag = (DM_Stag *)dm->data;
301: PetscFunctionBegin;
303: if (dof0) *dof0 = stag->dof[0];
304: if (dof1) *dof1 = stag->dof[1];
305: if (dof2) *dof2 = stag->dof[2];
306: if (dof3) *dof3 = stag->dof[3];
307: PetscFunctionReturn(PETSC_SUCCESS);
308: }
310: /*@C
311: DMStagGetGhostCorners - return global element indices of the local region, including ghost points
313: Not Collective
315: Input Parameter:
316: . dm - the `DMSTAG` object
318: Output Parameters:
319: + x - the starting element index in the first direction
320: . y - the starting element index in the second direction
321: . z - the starting element index in the third direction
322: . m - the element width in the first direction
323: . n - the element width in the second direction
324: - p - the element width in the third direction
326: Level: beginner
328: Note:
329: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to `NULL` in this case.
331: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetCorners()`, `DMDAGetGhostCorners()`
332: @*/
333: PetscErrorCode DMStagGetGhostCorners(DM dm, PetscInt *x, PetscInt *y, PetscInt *z, PetscInt *m, PetscInt *n, PetscInt *p)
334: {
335: const DM_Stag *const stag = (DM_Stag *)dm->data;
337: PetscFunctionBegin;
339: if (x) *x = stag->startGhost[0];
340: if (y) *y = stag->startGhost[1];
341: if (z) *z = stag->startGhost[2];
342: if (m) *m = stag->nGhost[0];
343: if (n) *n = stag->nGhost[1];
344: if (p) *p = stag->nGhost[2];
345: PetscFunctionReturn(PETSC_SUCCESS);
346: }
348: /*@C
349: DMStagGetGlobalSizes - get global element counts
351: Not Collective
353: Input Parameter:
354: . dm - the `DMSTAG` object
356: Output Parameters:
357: + M - global element counts in the x direction
358: . N - global element counts in the y direction
359: - P - global element counts in the z direction
361: Level: beginner
363: Note:
364: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to `NULL` in this case.
366: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetLocalSizes()`, `DMDAGetInfo()`
367: @*/
368: PetscErrorCode DMStagGetGlobalSizes(DM dm, PetscInt *M, PetscInt *N, PetscInt *P)
369: {
370: const DM_Stag *const stag = (DM_Stag *)dm->data;
372: PetscFunctionBegin;
374: if (M) *M = stag->N[0];
375: if (N) *N = stag->N[1];
376: if (P) *P = stag->N[2];
377: PetscFunctionReturn(PETSC_SUCCESS);
378: }
380: /*@C
381: DMStagGetIsFirstRank - get boolean value for whether this rank is first in each direction in the rank grid
383: Not Collective
385: Input Parameter:
386: . dm - the `DMSTAG` object
388: Output Parameters:
389: + isFirstRank0 - whether this rank is first in the x direction
390: . isFirstRank1 - whether this rank is first in the y direction
391: - isFirstRank2 - whether this rank is first in the z direction
393: Level: intermediate
395: Note:
396: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to `NULL` in this case.
398: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetIsLastRank()`
399: @*/
400: PetscErrorCode DMStagGetIsFirstRank(DM dm, PetscBool *isFirstRank0, PetscBool *isFirstRank1, PetscBool *isFirstRank2)
401: {
402: const DM_Stag *const stag = (DM_Stag *)dm->data;
404: PetscFunctionBegin;
406: if (isFirstRank0) *isFirstRank0 = stag->firstRank[0];
407: if (isFirstRank1) *isFirstRank1 = stag->firstRank[1];
408: if (isFirstRank2) *isFirstRank2 = stag->firstRank[2];
409: PetscFunctionReturn(PETSC_SUCCESS);
410: }
412: /*@C
413: DMStagGetIsLastRank - get boolean value for whether this rank is last in each direction in the rank grid
415: Not Collective
417: Input Parameter:
418: . dm - the `DMSTAG` object
420: Output Parameters:
421: + isLastRank0 - whether this rank is last in the x direction
422: . isLastRank1 - whether this rank is last in the y direction
423: - isLastRank2 - whether this rank is last in the z direction
425: Level: intermediate
427: Note:
428: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to `NULL` in this case.
430: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetIsFirstRank()`
431: @*/
432: PetscErrorCode DMStagGetIsLastRank(DM dm, PetscBool *isLastRank0, PetscBool *isLastRank1, PetscBool *isLastRank2)
433: {
434: const DM_Stag *const stag = (DM_Stag *)dm->data;
436: PetscFunctionBegin;
438: if (isLastRank0) *isLastRank0 = stag->lastRank[0];
439: if (isLastRank1) *isLastRank1 = stag->lastRank[1];
440: if (isLastRank2) *isLastRank2 = stag->lastRank[2];
441: PetscFunctionReturn(PETSC_SUCCESS);
442: }
444: /*@C
445: DMStagGetLocalSizes - get local elementwise sizes
447: Not Collective
449: Input Parameter:
450: . dm - the `DMSTAG` object
452: Output Parameters:
453: + m - local element counts (excluding ghosts) in the x direction
454: . n - local element counts (excluding ghosts) in the y direction
455: - p - local element counts (excluding ghosts) in the z direction
457: Level: beginner
459: Note:
460: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to `NULL` in this case.
462: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetGlobalSizes()`, `DMStagGetDOF()`, `DMStagGetNumRanks()`, `DMDAGetLocalInfo()`
463: @*/
464: PetscErrorCode DMStagGetLocalSizes(DM dm, PetscInt *m, PetscInt *n, PetscInt *p)
465: {
466: const DM_Stag *const stag = (DM_Stag *)dm->data;
468: PetscFunctionBegin;
470: if (m) *m = stag->n[0];
471: if (n) *n = stag->n[1];
472: if (p) *p = stag->n[2];
473: PetscFunctionReturn(PETSC_SUCCESS);
474: }
476: /*@C
477: DMStagGetNumRanks - get number of ranks in each direction in the global grid decomposition
479: Not Collective
481: Input Parameter:
482: . dm - the `DMSTAG` object
484: Output Parameters:
485: + nRanks0 - number of ranks in the x direction in the grid decomposition
486: . nRanks1 - number of ranks in the y direction in the grid decomposition
487: - nRanks2 - number of ranks in the z direction in the grid decomposition
489: Level: intermediate
491: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetGlobalSizes()`, `DMStagGetLocalSize()`, `DMStagSetNumRanks()`, `DMDAGetInfo()`
492: @*/
493: PetscErrorCode DMStagGetNumRanks(DM dm, PetscInt *nRanks0, PetscInt *nRanks1, PetscInt *nRanks2)
494: {
495: const DM_Stag *const stag = (DM_Stag *)dm->data;
497: PetscFunctionBegin;
499: if (nRanks0) *nRanks0 = stag->nRanks[0];
500: if (nRanks1) *nRanks1 = stag->nRanks[1];
501: if (nRanks2) *nRanks2 = stag->nRanks[2];
502: PetscFunctionReturn(PETSC_SUCCESS);
503: }
505: /*@C
506: DMStagGetEntries - get number of native entries in the global representation
508: Not Collective
510: Input Parameter:
511: . dm - the `DMSTAG` object
513: Output Parameter:
514: . entries - number of rank-native entries in the global representation
516: Level: developer
518: Note:
519: This is the number of entries on this rank for a global vector associated with `dm`.
520: That is, it is value of `size` returned by `VecGetLocalSize(vec,&size)` when
521: `DMCreateGlobalVector(dm,&vec) is used to create a `Vec`. Users would typically
522: use these functions.
524: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetDOF()`, `DMStagGetEntriesLocal()`, `DMStagGetEntriesPerElement()`, `DMCreateLocalVector()`
525: @*/
526: PetscErrorCode DMStagGetEntries(DM dm, PetscInt *entries)
527: {
528: const DM_Stag *const stag = (DM_Stag *)dm->data;
530: PetscFunctionBegin;
532: if (entries) *entries = stag->entries;
533: PetscFunctionReturn(PETSC_SUCCESS);
534: }
536: /*@C
537: DMStagGetEntriesLocal - get number of entries in the local representation
539: Not Collective
541: Input Parameter:
542: . dm - the `DMSTAG` object
544: Output Parameter:
545: . entries - number of entries in the local representation
547: Level: developer
549: Note:
550: This is the number of entries on this rank in the local representation.
551: That is, it is value of `size` returned by `VecGetSize(vec,&size)` or
552: `VecGetLocalSize(vec,&size)` when `DMCreateLocalVector(dm,&vec)` is used to
553: create a `Vec`. Users would typically use these functions.
555: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetDOF()`, `DMStagGetEntries()`, `DMStagGetEntriesPerElement()`, `DMCreateLocalVector()`
556: @*/
557: PetscErrorCode DMStagGetEntriesLocal(DM dm, PetscInt *entries)
558: {
559: const DM_Stag *const stag = (DM_Stag *)dm->data;
561: PetscFunctionBegin;
563: if (entries) *entries = stag->entriesGhost;
564: PetscFunctionReturn(PETSC_SUCCESS);
565: }
567: /*@C
568: DMStagGetEntriesPerElement - get number of entries per element in the local representation
570: Not Collective
572: Input Parameter:
573: . dm - the `DMSTAG` object
575: Output Parameter:
576: . entriesPerElement - number of entries associated with each element in the local representation
578: Level: developer
580: Notes:
581: This is the natural block size for most local operations. In 1D it is equal to `dof0` $+$ `dof1`,
582: in 2D it is equal to `dof0` $+ 2$`dof1` $+$ `dof2`, and in 3D it is equal to `dof0` $+ 3$`dof1` $+ 3$`dof2` $+$ `dof3`
584: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetDOF()`
585: @*/
586: PetscErrorCode DMStagGetEntriesPerElement(DM dm, PetscInt *entriesPerElement)
587: {
588: const DM_Stag *const stag = (DM_Stag *)dm->data;
590: PetscFunctionBegin;
592: if (entriesPerElement) *entriesPerElement = stag->entriesPerElement;
593: PetscFunctionReturn(PETSC_SUCCESS);
594: }
596: /*@C
597: DMStagGetStencilType - get elementwise ghost/halo stencil type
599: Not Collective
601: Input Parameter:
602: . dm - the `DMSTAG` object
604: Output Parameter:
605: . stencilType - the elementwise ghost stencil type: `DMSTAG_STENCIL_BOX`, `DMSTAG_STENCIL_STAR`, or `DMSTAG_STENCIL_NONE`
607: Level: beginner
609: .seealso: [](ch_stag), `DMSTAG`, `DMStagSetStencilType()`, `DMStagGetStencilWidth`, `DMStagStencilType`
610: @*/
611: PetscErrorCode DMStagGetStencilType(DM dm, DMStagStencilType *stencilType)
612: {
613: DM_Stag *const stag = (DM_Stag *)dm->data;
615: PetscFunctionBegin;
617: *stencilType = stag->stencilType;
618: PetscFunctionReturn(PETSC_SUCCESS);
619: }
621: /*@C
622: DMStagGetStencilWidth - get elementwise stencil width
624: Not Collective
626: Input Parameter:
627: . dm - the `DMSTAG` object
629: Output Parameter:
630: . stencilWidth - stencil/halo/ghost width in elements
632: Level: beginner
634: .seealso: [](ch_stag), `DMSTAG`, `DMStagSetStencilWidth()`, `DMStagGetStencilType()`, `DMDAGetStencilType()`
635: @*/
636: PetscErrorCode DMStagGetStencilWidth(DM dm, PetscInt *stencilWidth)
637: {
638: const DM_Stag *const stag = (DM_Stag *)dm->data;
640: PetscFunctionBegin;
642: if (stencilWidth) *stencilWidth = stag->stencilWidth;
643: PetscFunctionReturn(PETSC_SUCCESS);
644: }
646: /*@C
647: DMStagGetOwnershipRanges - get elements per rank in each direction
649: Not Collective
651: Input Parameter:
652: . dm - the `DMSTAG` object
654: Output Parameters:
655: + lx - ownership along x direction (optional)
656: . ly - ownership along y direction (optional)
657: - lz - ownership along z direction (optional)
659: Level: intermediate
661: Notes:
662: These correspond to the optional final arguments passed to `DMStagCreate1d()`, `DMStagCreate2d()`, and `DMStagCreate3d()`.
664: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to `NULL` in this case.
666: In C you should not free these arrays, nor change the values in them.
667: They will only have valid values while the `DMSTAG` they came from still exists (has not been destroyed).
669: .seealso: [](ch_stag), `DMSTAG`, `DMStagSetGlobalSizes()`, `DMStagSetOwnershipRanges()`, `DMStagCreate1d()`, `DMStagCreate2d()`, `DMStagCreate3d()`, `DMDAGetOwnershipRanges()`
670: @*/
671: PetscErrorCode DMStagGetOwnershipRanges(DM dm, const PetscInt *lx[], const PetscInt *ly[], const PetscInt *lz[])
672: {
673: const DM_Stag *const stag = (DM_Stag *)dm->data;
675: PetscFunctionBegin;
677: if (lx) *lx = stag->l[0];
678: if (ly) *ly = stag->l[1];
679: if (lz) *lz = stag->l[2];
680: PetscFunctionReturn(PETSC_SUCCESS);
681: }
683: /*@C
684: DMStagCreateCompatibleDMStag - create a compatible `DMSTAG` with different dof/stratum
686: Collective
688: Input Parameters:
689: + dm - the `DMSTAG` object
690: . dof0 - number of dof on the first stratum in the new `DMSTAG`
691: . dof1 - number of dof on the second stratum in the new `DMSTAG`
692: . dof2 - number of dof on the third stratum in the new `DMSTAG`
693: - dof3 - number of dof on the fourth stratum in the new `DMSTAG`
695: Output Parameter:
696: . newdm - the new, compatible `DMSTAG`
698: Level: intermediate
700: Notes:
701: DOF supplied for strata too big for the dimension are ignored; these may be set to `0`.
702: For example, for a 2-dimensional `DMSTAG`, `dof2` sets the number of dof per element,
703: and `dof3` is unused. For a 3-dimensional `DMSTAG`, `dof3` sets the number of DOF per element.
705: In contrast to `DMDACreateCompatibleDMDA()`, coordinates are not reused.
707: .seealso: [](ch_stag), `DMSTAG`, `DMDACreateCompatibleDMDA()`, `DMGetCompatibility()`, `DMStagMigrateVec()`
708: @*/
709: PetscErrorCode DMStagCreateCompatibleDMStag(DM dm, PetscInt dof0, PetscInt dof1, PetscInt dof2, PetscInt dof3, DM *newdm)
710: {
711: PetscFunctionBegin;
713: PetscCall(DMStagDuplicateWithoutSetup(dm, PetscObjectComm((PetscObject)dm), newdm));
714: PetscCall(DMStagSetDOF(*newdm, dof0, dof1, dof2, dof3));
715: PetscCall(DMSetUp(*newdm));
716: PetscFunctionReturn(PETSC_SUCCESS);
717: }
719: /*@C
720: DMStagGetLocationSlot - get index to use in accessing raw local arrays
722: Not Collective
724: Input Parameters:
725: + dm - the `DMSTAG` object
726: . loc - location relative to an element
727: - c - component
729: Output Parameter:
730: . slot - index to use
732: Level: beginner
734: Notes:
735: Provides an appropriate index to use with `DMStagVecGetArray()` and friends.
736: This is required so that the user doesn't need to know about the ordering of
737: dof associated with each local element.
739: .seealso: [](ch_stag), `DMSTAG`, `DMStagVecGetArray()`, `DMStagVecGetArrayRead()`, `DMStagGetDOF()`, `DMStagGetEntriesPerElement()`
740: @*/
741: PetscErrorCode DMStagGetLocationSlot(DM dm, DMStagStencilLocation loc, PetscInt c, PetscInt *slot)
742: {
743: DM_Stag *const stag = (DM_Stag *)dm->data;
745: PetscFunctionBegin;
747: if (PetscDefined(USE_DEBUG)) {
748: PetscInt dof;
749: PetscCall(DMStagGetLocationDOF(dm, loc, &dof));
750: PetscCheck(dof >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Location %s has no dof attached", DMStagStencilLocations[loc]);
751: PetscCheck(c <= dof - 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Supplied component number (%" PetscInt_FMT ") for location %s is too big (maximum %" PetscInt_FMT ")", c, DMStagStencilLocations[loc], dof - 1);
752: }
753: *slot = stag->locationOffsets[loc] + c;
754: PetscFunctionReturn(PETSC_SUCCESS);
755: }
757: /*@C
758: DMStagMigrateVec - transfer a vector associated with a `DMSTAG` to a vector associated with a compatible `DMSTAG`
760: Collective
762: Input Parameters:
763: + dm - the source `DMSTAG` object
764: . vec - the source vector, compatible with `dm`
765: . dmTo - the compatible destination `DMSTAG` object
766: - vecTo - the destination vector, compatible with `dmTo`
768: Level: advanced
770: Notes:
771: Extra dof are ignored, and unfilled dof are zeroed.
772: Currently only implemented to migrate global vectors to global vectors.
773: For the definition of compatibility of `DM`s, see `DMGetCompatibility()`.
775: .seealso: [](ch_stag), `DMSTAG`, `DMStagCreateCompatibleDMStag()`, `DMGetCompatibility()`, `DMStagVecSplitToDMDA()`
776: @*/
777: PetscErrorCode DMStagMigrateVec(DM dm, Vec vec, DM dmTo, Vec vecTo)
778: {
779: DM_Stag *const stag = (DM_Stag *)dm->data;
780: DM_Stag *const stagTo = (DM_Stag *)dmTo->data;
781: PetscInt nLocalTo, nLocal, dim, i, j, k;
782: PetscInt start[DMSTAG_MAX_DIM], startGhost[DMSTAG_MAX_DIM], n[DMSTAG_MAX_DIM], nExtra[DMSTAG_MAX_DIM], offset[DMSTAG_MAX_DIM];
783: Vec vecToLocal, vecLocal;
784: PetscBool compatible, compatibleSet;
785: const PetscScalar *arr;
786: PetscScalar *arrTo;
787: const PetscInt epe = stag->entriesPerElement;
788: const PetscInt epeTo = stagTo->entriesPerElement;
790: PetscFunctionBegin;
795: PetscCall(DMGetCompatibility(dm, dmTo, &compatible, &compatibleSet));
796: PetscCheck(compatibleSet && compatible, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "DMStag objects must be shown to be compatible");
797: PetscCall(DMGetDimension(dm, &dim));
798: PetscCall(VecGetLocalSize(vec, &nLocal));
799: PetscCall(VecGetLocalSize(vecTo, &nLocalTo));
800: PetscCheck(nLocal == stag->entries && nLocalTo == stagTo->entries, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Vector migration only implemented for global vector to global vector.");
801: PetscCall(DMStagGetCorners(dm, &start[0], &start[1], &start[2], &n[0], &n[1], &n[2], &nExtra[0], &nExtra[1], &nExtra[2]));
802: PetscCall(DMStagGetGhostCorners(dm, &startGhost[0], &startGhost[1], &startGhost[2], NULL, NULL, NULL));
803: for (i = 0; i < DMSTAG_MAX_DIM; ++i) offset[i] = start[i] - startGhost[i];
805: /* Proceed by transferring to a local vector, copying, and transferring back to a global vector */
806: PetscCall(DMGetLocalVector(dm, &vecLocal));
807: PetscCall(DMGetLocalVector(dmTo, &vecToLocal));
808: PetscCall(DMGlobalToLocalBegin(dm, vec, INSERT_VALUES, vecLocal));
809: PetscCall(DMGlobalToLocalEnd(dm, vec, INSERT_VALUES, vecLocal));
810: PetscCall(VecGetArrayRead(vecLocal, &arr));
811: PetscCall(VecGetArray(vecToLocal, &arrTo));
812: /* Note that some superfluous copying of entries on partial dummy elements is done */
813: if (dim == 1) {
814: for (i = offset[0]; i < offset[0] + n[0] + nExtra[0]; ++i) {
815: PetscInt d = 0, dTo = 0, b = 0, bTo = 0;
816: PetscInt si;
817: for (si = 0; si < 2; ++si) {
818: b += stag->dof[si];
819: bTo += stagTo->dof[si];
820: for (; d < b && dTo < bTo; ++d, ++dTo) arrTo[i * epeTo + dTo] = arr[i * epe + d];
821: for (; dTo < bTo; ++dTo) arrTo[i * epeTo + dTo] = 0.0;
822: d = b;
823: }
824: }
825: } else if (dim == 2) {
826: const PetscInt epr = stag->nGhost[0] * epe;
827: const PetscInt eprTo = stagTo->nGhost[0] * epeTo;
828: for (j = offset[1]; j < offset[1] + n[1] + nExtra[1]; ++j) {
829: for (i = offset[0]; i < offset[0] + n[0] + nExtra[0]; ++i) {
830: const PetscInt base = j * epr + i * epe;
831: const PetscInt baseTo = j * eprTo + i * epeTo;
832: PetscInt d = 0, dTo = 0, b = 0, bTo = 0;
833: const PetscInt s[4] = {0, 1, 1, 2}; /* Dimensions of points, in order */
834: PetscInt si;
835: for (si = 0; si < 4; ++si) {
836: b += stag->dof[s[si]];
837: bTo += stagTo->dof[s[si]];
838: for (; d < b && dTo < bTo; ++d, ++dTo) arrTo[baseTo + dTo] = arr[base + d];
839: for (; dTo < bTo; ++dTo) arrTo[baseTo + dTo] = 0.0;
840: d = b;
841: }
842: }
843: }
844: } else if (dim == 3) {
845: const PetscInt epr = stag->nGhost[0] * epe;
846: const PetscInt eprTo = stagTo->nGhost[0] * epeTo;
847: const PetscInt epl = stag->nGhost[1] * epr;
848: const PetscInt eplTo = stagTo->nGhost[1] * eprTo;
849: for (k = offset[2]; k < offset[2] + n[2] + nExtra[2]; ++k) {
850: for (j = offset[1]; j < offset[1] + n[1] + nExtra[1]; ++j) {
851: for (i = offset[0]; i < offset[0] + n[0] + nExtra[0]; ++i) {
852: PetscInt d = 0, dTo = 0, b = 0, bTo = 0;
853: const PetscInt base = k * epl + j * epr + i * epe;
854: const PetscInt baseTo = k * eplTo + j * eprTo + i * epeTo;
855: const PetscInt s[8] = {0, 1, 1, 2, 1, 2, 2, 3}; /* dimensions of points, in order */
856: PetscInt is;
857: for (is = 0; is < 8; ++is) {
858: b += stag->dof[s[is]];
859: bTo += stagTo->dof[s[is]];
860: for (; d < b && dTo < bTo; ++d, ++dTo) arrTo[baseTo + dTo] = arr[base + d];
861: for (; dTo < bTo; ++dTo) arrTo[baseTo + dTo] = 0.0;
862: d = b;
863: }
864: }
865: }
866: }
867: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %" PetscInt_FMT, dim);
868: PetscCall(VecRestoreArrayRead(vecLocal, &arr));
869: PetscCall(VecRestoreArray(vecToLocal, &arrTo));
870: PetscCall(DMRestoreLocalVector(dm, &vecLocal));
871: PetscCall(DMLocalToGlobalBegin(dmTo, vecToLocal, INSERT_VALUES, vecTo));
872: PetscCall(DMLocalToGlobalEnd(dmTo, vecToLocal, INSERT_VALUES, vecTo));
873: PetscCall(DMRestoreLocalVector(dmTo, &vecToLocal));
874: PetscFunctionReturn(PETSC_SUCCESS);
875: }
877: /*@C
878: DMStagPopulateLocalToGlobalInjective - populate an internal 1-to-1 local-to-global map
880: Collective
882: Creates an internal object which explicitly maps a single local degree of
883: freedom to each global degree of freedom. This is used, if populated,
884: instead of SCATTER_REVERSE_LOCAL with the (1-to-many, in general)
885: global-to-local map, when DMLocalToGlobal() is called with INSERT_VALUES.
886: This allows usage, for example, even in the periodic, 1-rank case, where
887: the inverse of the global-to-local map, even when restricted to on-rank
888: communication, is non-injective. This is at the cost of storing an additional
889: VecScatter object inside each `DMSTAG` object.
891: Input Parameter:
892: . dm - the `DMSTAG` object
894: Level: developer
896: Notes:
897: In normal usage, library users shouldn't be concerned with this function,
898: as it is called during `DMSetUp()`, when required.
900: Returns immediately if the internal map is already populated.
902: Developer Notes:
903: This could, if desired, be moved up to a general `DM` routine. It would allow,
904: for example, `DMDA` to support `DMLocalToGlobal()` with `INSERT_VALUES`,
905: even in the single-rank periodic case.
907: .seealso: [](ch_stag), `DMSTAG`, `DMLocalToGlobal()`, `VecScatter`
908: @*/
909: PetscErrorCode DMStagPopulateLocalToGlobalInjective(DM dm)
910: {
911: PetscInt dim;
912: DM_Stag *const stag = (DM_Stag *)dm->data;
914: PetscFunctionBegin;
916: if (stag->ltog_injective) PetscFunctionReturn(PETSC_SUCCESS); /* Don't re-populate */
917: PetscCall(DMGetDimension(dm, &dim));
918: switch (dim) {
919: case 1:
920: PetscCall(DMStagPopulateLocalToGlobalInjective_1d(dm));
921: break;
922: case 2:
923: PetscCall(DMStagPopulateLocalToGlobalInjective_2d(dm));
924: break;
925: case 3:
926: PetscCall(DMStagPopulateLocalToGlobalInjective_3d(dm));
927: break;
928: default:
929: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %" PetscInt_FMT, dim);
930: }
931: PetscFunctionReturn(PETSC_SUCCESS);
932: }
934: static PetscErrorCode DMStagRestoreProductCoordinateArrays_Private(DM dm, void *arrX, void *arrY, void *arrZ, PetscBool read)
935: {
936: PetscInt dim, d;
937: void *arr[DMSTAG_MAX_DIM];
938: DM dmCoord;
940: PetscFunctionBegin;
942: PetscCall(DMGetDimension(dm, &dim));
943: PetscCheck(dim <= DMSTAG_MAX_DIM, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented for %" PetscInt_FMT " dimensions", dim);
944: arr[0] = arrX;
945: arr[1] = arrY;
946: arr[2] = arrZ;
947: PetscCall(DMGetCoordinateDM(dm, &dmCoord));
948: for (d = 0; d < dim; ++d) {
949: DM subDM;
950: Vec coord1d_local;
952: /* Ignore unrequested arrays */
953: if (!arr[d]) continue;
955: PetscCall(DMProductGetDM(dmCoord, d, &subDM));
956: PetscCall(DMGetCoordinatesLocal(subDM, &coord1d_local));
957: if (read) {
958: PetscCall(DMStagVecRestoreArrayRead(subDM, coord1d_local, arr[d]));
959: } else {
960: PetscCall(DMStagVecRestoreArray(subDM, coord1d_local, arr[d]));
961: }
962: }
963: PetscFunctionReturn(PETSC_SUCCESS);
964: }
966: /*@C
967: DMStagRestoreProductCoordinateArrays - restore local array access
969: Logically Collective
971: Input Parameter:
972: . dm - the `DMSTAG` object
974: Output Parameters:
975: + arrX - local 1D coordinate arrays for x direction
976: . arrY - local 1D coordinate arrays for y direction
977: - arrZ - local 1D coordinate arrays for z direction
979: Level: intermediate
981: Notes:
982: This function does not automatically perform a local->global scatter to populate global coordinates from the local coordinates.
983: Thus, it may be required to explicitly perform these operations in some situations, as in the following partial example\:
984: .vb
985: PetscCall(DMGetCoordinateDM(dm, &cdm));
986: for (PetscInt d = 0; d < 3; ++d) {
987: DM subdm;
988: Vec coor, coor_local;
990: PetscCall(DMProductGetDM(cdm, d, &subdm));
991: PetscCall(DMGetCoordinates(subdm, &coor));
992: PetscCall(DMGetCoordinatesLocal(subdm, &coor_local));
993: PetscCall(DMLocalToGlobal(subdm, coor_local, INSERT_VALUES, coor));
994: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Coordinates dim %" PetscInt_FMT ":\n", d));
995: PetscCall(VecView(coor, PETSC_VIEWER_STDOUT_WORLD));
996: }
997: .ve
999: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetProductCoordinateArrays()`, `DMStagGetProductCoordinateArraysRead()`
1000: @*/
1001: PetscErrorCode DMStagRestoreProductCoordinateArrays(DM dm, void *arrX, void *arrY, void *arrZ)
1002: {
1003: PetscFunctionBegin;
1004: PetscCall(DMStagRestoreProductCoordinateArrays_Private(dm, arrX, arrY, arrZ, PETSC_FALSE));
1005: PetscFunctionReturn(PETSC_SUCCESS);
1006: }
1008: /*@C
1009: DMStagRestoreProductCoordinateArraysRead - restore local product array access, read-only
1011: Logically Collective
1013: Input Parameters:
1014: + dm - the `DMSTAG` object
1015: . arrX - local 1D coordinate arrays for x direction
1016: . arrY - local 1D coordinate arrays for y direction
1017: - arrZ - local 1D coordinate arrays for z direction
1019: Level: intermediate
1021: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetProductCoordinateArrays()`, `DMStagGetProductCoordinateArraysRead()`
1022: @*/
1023: PetscErrorCode DMStagRestoreProductCoordinateArraysRead(DM dm, void *arrX, void *arrY, void *arrZ)
1024: {
1025: PetscFunctionBegin;
1026: PetscCall(DMStagRestoreProductCoordinateArrays_Private(dm, arrX, arrY, arrZ, PETSC_TRUE));
1027: PetscFunctionReturn(PETSC_SUCCESS);
1028: }
1030: /*@C
1031: DMStagSetBoundaryTypes - set `DMSTAG` boundary types
1033: Logically Collective; boundaryType0, boundaryType1, and boundaryType2 must contain common values
1035: Input Parameters:
1036: + dm - the `DMSTAG` object
1037: . boundaryType2 - boundary type for x direction
1038: . boundaryType1 - boundary type for y direction, not set for one dimensional problems
1039: - boundaryType0 - boundary type for z direction, not set for one and two dimensional problems
1041: Level: advanced
1043: Note:
1044: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1046: .seealso: [](ch_stag), `DMSTAG`, `DMBoundaryType`, `DMStagCreate1d()`, `DMStagCreate2d()`, `DMStagCreate3d()`, `DMDASetBoundaryType()`
1047: @*/
1048: PetscErrorCode DMStagSetBoundaryTypes(DM dm, DMBoundaryType boundaryType0, DMBoundaryType boundaryType1, DMBoundaryType boundaryType2)
1049: {
1050: DM_Stag *const stag = (DM_Stag *)dm->data;
1051: PetscInt dim;
1053: PetscFunctionBegin;
1054: PetscCall(DMGetDimension(dm, &dim));
1059: PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
1060: stag->boundaryType[0] = boundaryType0;
1061: if (dim > 1) stag->boundaryType[1] = boundaryType1;
1062: if (dim > 2) stag->boundaryType[2] = boundaryType2;
1063: PetscFunctionReturn(PETSC_SUCCESS);
1064: }
1066: /*@C
1067: DMStagSetCoordinateDMType - set DM type to store coordinates
1069: Logically Collective; `dmtype` must contain common value
1071: Input Parameters:
1072: + dm - the `DMSTAG` object
1073: - dmtype - DMtype for coordinates, either `DMSTAG` or `DMPRODUCT`
1075: Level: advanced
1077: .seealso: [](ch_stag), `DMSTAG`, `DMPRODUCT`, `DMGetCoordinateDM()`, `DMStagSetUniformCoordinates()`, `DMStagSetUniformCoordinatesExplicit()`, `DMStagSetUniformCoordinatesProduct()`, `DMType`
1078: @*/
1079: PetscErrorCode DMStagSetCoordinateDMType(DM dm, DMType dmtype)
1080: {
1081: DM_Stag *const stag = (DM_Stag *)dm->data;
1083: PetscFunctionBegin;
1085: PetscCall(PetscFree(stag->coordinateDMType));
1086: PetscCall(PetscStrallocpy(dmtype, (char **)&stag->coordinateDMType));
1087: PetscFunctionReturn(PETSC_SUCCESS);
1088: }
1090: /*@C
1091: DMStagSetDOF - set dof/stratum
1093: Logically Collective; `dof0`, `dof1`, `dof2`, and `dof3` must contain common values
1095: Input Parameters:
1096: + dm - the `DMSTAG` object
1097: . dof0 - the number of points per 0-cell (vertex/node)
1098: . dof1 - the number of points per 1-cell (element in 1D, edge in 2D and 3D)
1099: . dof2 - the number of points per 2-cell (element in 2D, face in 3D)
1100: - dof3 - the number of points per 3-cell (element in 3D)
1102: Level: advanced
1104: Note:
1105: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1107: .seealso: [](ch_stag), `DMSTAG`, `DMDASetDof()`
1108: @*/
1109: PetscErrorCode DMStagSetDOF(DM dm, PetscInt dof0, PetscInt dof1, PetscInt dof2, PetscInt dof3)
1110: {
1111: DM_Stag *const stag = (DM_Stag *)dm->data;
1112: PetscInt dim;
1114: PetscFunctionBegin;
1120: PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
1121: PetscCall(DMGetDimension(dm, &dim));
1122: PetscCheck(dof0 >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "dof0 cannot be negative");
1123: PetscCheck(dof1 >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "dof1 cannot be negative");
1124: PetscCheck(dim <= 1 || dof2 >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "dof2 cannot be negative");
1125: PetscCheck(dim <= 2 || dof3 >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "dof3 cannot be negative");
1126: stag->dof[0] = dof0;
1127: stag->dof[1] = dof1;
1128: if (dim > 1) stag->dof[2] = dof2;
1129: if (dim > 2) stag->dof[3] = dof3;
1130: PetscFunctionReturn(PETSC_SUCCESS);
1131: }
1133: /*@C
1134: DMStagSetNumRanks - set ranks in each direction in the global rank grid
1136: Logically Collective; `nRanks0`, `nRanks1`, and `nRanks2` must contain common values
1138: Input Parameters:
1139: + dm - the `DMSTAG` object
1140: . nRanks0 - number of ranks in the x direction
1141: . nRanks1 - number of ranks in the y direction
1142: - nRanks2 - number of ranks in the z direction
1144: Level: developer
1146: Note:
1147: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1149: .seealso: [](ch_stag), `DMSTAG`, `DMDASetNumProcs()`
1150: @*/
1151: PetscErrorCode DMStagSetNumRanks(DM dm, PetscInt nRanks0, PetscInt nRanks1, PetscInt nRanks2)
1152: {
1153: DM_Stag *const stag = (DM_Stag *)dm->data;
1154: PetscInt dim;
1156: PetscFunctionBegin;
1161: PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
1162: PetscCall(DMGetDimension(dm, &dim));
1163: PetscCheck(nRanks0 == PETSC_DECIDE || nRanks0 >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "number of ranks in X direction cannot be less than 1");
1164: PetscCheck(dim <= 1 || nRanks1 == PETSC_DECIDE || nRanks1 >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "number of ranks in Y direction cannot be less than 1");
1165: PetscCheck(dim <= 2 || nRanks2 == PETSC_DECIDE || nRanks2 >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "number of ranks in Z direction cannot be less than 1");
1166: if (nRanks0) stag->nRanks[0] = nRanks0;
1167: if (dim > 1 && nRanks1) stag->nRanks[1] = nRanks1;
1168: if (dim > 2 && nRanks2) stag->nRanks[2] = nRanks2;
1169: PetscFunctionReturn(PETSC_SUCCESS);
1170: }
1172: /*@C
1173: DMStagSetStencilType - set elementwise ghost/halo stencil type
1175: Logically Collective; `stencilType` must contain common value
1177: Input Parameters:
1178: + dm - the `DMSTAG` object
1179: - stencilType - the elementwise ghost stencil type: `DMSTAG_STENCIL_BOX`, `DMSTAG_STENCIL_STAR`, or `DMSTAG_STENCIL_NONE`
1181: Level: beginner
1183: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetStencilType()`, `DMStagSetStencilWidth()`, `DMStagStencilType`
1184: @*/
1185: PetscErrorCode DMStagSetStencilType(DM dm, DMStagStencilType stencilType)
1186: {
1187: DM_Stag *const stag = (DM_Stag *)dm->data;
1189: PetscFunctionBegin;
1192: PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
1193: stag->stencilType = stencilType;
1194: PetscFunctionReturn(PETSC_SUCCESS);
1195: }
1197: /*@C
1198: DMStagSetStencilWidth - set elementwise stencil width
1200: Logically Collective; `stencilWidth` must contain common value
1202: Input Parameters:
1203: + dm - the `DMSTAG` object
1204: - stencilWidth - stencil/halo/ghost width in elements
1206: Level: beginner
1208: Note:
1209: The width value is not used when `DMSTAG_STENCIL_NONE` is specified.
1211: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetStencilWidth()`, `DMStagGetStencilType()`, `DMStagStencilType`
1212: @*/
1213: PetscErrorCode DMStagSetStencilWidth(DM dm, PetscInt stencilWidth)
1214: {
1215: DM_Stag *const stag = (DM_Stag *)dm->data;
1217: PetscFunctionBegin;
1220: PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
1221: PetscCheck(stencilWidth >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil width must be non-negative");
1222: stag->stencilWidth = stencilWidth;
1223: PetscFunctionReturn(PETSC_SUCCESS);
1224: }
1226: /*@C
1227: DMStagSetGlobalSizes - set global element counts in each direction
1229: Logically Collective; `N0`, `N1`, and `N2` must contain common values
1231: Input Parameters:
1232: + dm - the `DMSTAG` object
1233: . N0 - global elementwise size in the x direction
1234: . N1 - global elementwise size in the y direction
1235: - N2 - global elementwise size in the z direction
1237: Level: advanced
1239: Note:
1240: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1242: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetGlobalSizes()`, `DMDASetSizes()`
1243: @*/
1244: PetscErrorCode DMStagSetGlobalSizes(DM dm, PetscInt N0, PetscInt N1, PetscInt N2)
1245: {
1246: DM_Stag *const stag = (DM_Stag *)dm->data;
1247: PetscInt dim;
1249: PetscFunctionBegin;
1251: PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
1252: PetscCall(DMGetDimension(dm, &dim));
1253: PetscCheck(N0 >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Number of elements in X direction must be positive");
1254: PetscCheck(dim <= 1 || N1 >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Number of elements in Y direction must be positive");
1255: PetscCheck(dim <= 2 || N2 >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Number of elements in Z direction must be positive");
1256: if (N0) stag->N[0] = N0;
1257: if (N1) stag->N[1] = N1;
1258: if (N2) stag->N[2] = N2;
1259: PetscFunctionReturn(PETSC_SUCCESS);
1260: }
1262: /*@C
1263: DMStagSetOwnershipRanges - set elements per rank in each direction
1265: Logically Collective; `lx`, `ly`, and `lz` must contain common values
1267: Input Parameters:
1268: + dm - the `DMSTAG` object
1269: . lx - element counts for each rank in the x direction
1270: . ly - element counts for each rank in the y direction
1271: - lz - element counts for each rank in the z direction
1273: Level: developer
1275: Note:
1276: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to `NULL` in this case.
1278: .seealso: [](ch_stag), `DMSTAG`, `DMStagSetGlobalSizes()`, `DMStagGetOwnershipRanges()`, `DMDASetOwnershipRanges()`
1279: @*/
1280: PetscErrorCode DMStagSetOwnershipRanges(DM dm, PetscInt const *lx, PetscInt const *ly, PetscInt const *lz)
1281: {
1282: DM_Stag *const stag = (DM_Stag *)dm->data;
1283: const PetscInt *lin[3];
1284: PetscInt d, dim;
1286: PetscFunctionBegin;
1288: PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
1289: lin[0] = lx;
1290: lin[1] = ly;
1291: lin[2] = lz;
1292: PetscCall(DMGetDimension(dm, &dim));
1293: for (d = 0; d < dim; ++d) {
1294: if (lin[d]) {
1295: PetscCheck(stag->nRanks[d] >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of ranks");
1296: if (!stag->l[d]) PetscCall(PetscMalloc1(stag->nRanks[d], &stag->l[d]));
1297: PetscCall(PetscArraycpy(stag->l[d], lin[d], stag->nRanks[d]));
1298: }
1299: }
1300: PetscFunctionReturn(PETSC_SUCCESS);
1301: }
1303: /*@C
1304: DMStagSetUniformCoordinates - set `DMSTAG` coordinates to be a uniform grid
1306: Collective
1308: Input Parameters:
1309: + dm - the `DMSTAG` object
1310: . xmin - minimum global coordinate value in the x direction
1311: . xmax - maximum global coordinate values in the x direction
1312: . ymin - minimum global coordinate value in the y direction
1313: . ymax - maximum global coordinate value in the y direction
1314: . zmin - minimum global coordinate value in the z direction
1315: - zmax - maximum global coordinate value in the z direction
1317: Level: advanced
1319: Notes:
1320: `DMSTAG` supports 2 different types of coordinate `DM`: `DMSTAG` and `DMPRODUCT`.
1321: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1323: Local coordinates are populated (using `DMSetCoordinatesLocal()`), linearly
1324: extrapolated to ghost cells, including those outside the physical domain.
1325: This is also done in case of periodic boundaries, meaning that the same
1326: global point may have different coordinates in different local representations,
1327: which are equivalent assuming a periodicity implied by the arguments to this function,
1328: i.e. two points are equivalent if their difference is a multiple of $($`xmax` $-$ `xmin` $)$
1329: in the x direction, $($ `ymax` $-$ `ymin` $)$ in the y direction, and $($ `zmax` $-$ `zmin` $)$ in the z direction.
1331: .seealso: [](ch_stag), `DMSTAG`, `DMPRODUCT`, `DMStagSetUniformCoordinatesExplicit()`, `DMStagSetUniformCoordinatesProduct()`, `DMStagSetCoordinateDMType()`, `DMGetCoordinateDM()`, `DMGetCoordinates()`, `DMDASetUniformCoordinates()`, `DMBoundaryType`
1332: @*/
1333: PetscErrorCode DMStagSetUniformCoordinates(DM dm, PetscReal xmin, PetscReal xmax, PetscReal ymin, PetscReal ymax, PetscReal zmin, PetscReal zmax)
1334: {
1335: DM_Stag *const stag = (DM_Stag *)dm->data;
1336: PetscBool flg_stag, flg_product;
1338: PetscFunctionBegin;
1340: PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This function must be called after DMSetUp()");
1341: PetscCheck(stag->coordinateDMType, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "You must first call DMStagSetCoordinateDMType()");
1342: PetscCall(PetscStrcmp(stag->coordinateDMType, DMSTAG, &flg_stag));
1343: PetscCall(PetscStrcmp(stag->coordinateDMType, DMPRODUCT, &flg_product));
1344: if (flg_stag) {
1345: PetscCall(DMStagSetUniformCoordinatesExplicit(dm, xmin, xmax, ymin, ymax, zmin, zmax));
1346: } else if (flg_product) {
1347: PetscCall(DMStagSetUniformCoordinatesProduct(dm, xmin, xmax, ymin, ymax, zmin, zmax));
1348: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported DM Type %s", stag->coordinateDMType);
1349: PetscFunctionReturn(PETSC_SUCCESS);
1350: }
1352: /*@C
1353: DMStagSetUniformCoordinatesExplicit - set `DMSTAG` coordinates to be a uniform grid, storing all values
1355: Collective
1357: Input Parameters:
1358: + dm - the `DMSTAG` object
1359: . xmin - minimum global coordinate value in the x direction
1360: . xmax - maximum global coordinate values in the x direction
1361: . ymin - minimum global coordinate value in the y direction
1362: . ymax - maximum global coordinate value in the y direction
1363: . zmin - minimum global coordinate value in the z direction
1364: - zmax - maximum global coordinate value in the z direction
1366: Level: beginner
1368: Notes:
1369: `DMSTAG` supports 2 different types of coordinate `DM`: either another `DMSTAG`, or a `DMPRODUCT`.
1370: If the grid is orthogonal, using `DMPRODUCT` should be more efficient.
1372: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1374: See the manual page for `DMStagSetUniformCoordinates()` for information on how
1375: coordinates for dummy cells outside the physical domain boundary are populated.
1377: .seealso: [](ch_stag), `DMSTAG`, `DMStagSetUniformCoordinates()`, `DMStagSetUniformCoordinatesProduct()`, `DMStagSetCoordinateDMType()`
1378: @*/
1379: PetscErrorCode DMStagSetUniformCoordinatesExplicit(DM dm, PetscReal xmin, PetscReal xmax, PetscReal ymin, PetscReal ymax, PetscReal zmin, PetscReal zmax)
1380: {
1381: DM_Stag *const stag = (DM_Stag *)dm->data;
1382: PetscInt dim;
1383: PetscBool flg;
1385: PetscFunctionBegin;
1387: PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This function must be called after DMSetUp()");
1388: PetscCall(PetscStrcmp(stag->coordinateDMType, DMSTAG, &flg));
1389: PetscCheck(!stag->coordinateDMType || flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Refusing to change an already-set DM coordinate type");
1390: PetscCall(DMStagSetCoordinateDMType(dm, DMSTAG));
1391: PetscCall(DMGetDimension(dm, &dim));
1392: switch (dim) {
1393: case 1:
1394: PetscCall(DMStagSetUniformCoordinatesExplicit_1d(dm, xmin, xmax));
1395: break;
1396: case 2:
1397: PetscCall(DMStagSetUniformCoordinatesExplicit_2d(dm, xmin, xmax, ymin, ymax));
1398: break;
1399: case 3:
1400: PetscCall(DMStagSetUniformCoordinatesExplicit_3d(dm, xmin, xmax, ymin, ymax, zmin, zmax));
1401: break;
1402: default:
1403: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %" PetscInt_FMT, dim);
1404: }
1405: PetscCall(DMCoarsenHookRemove(dm, DMRestrictHook_Coordinates, NULL, NULL));
1406: PetscFunctionReturn(PETSC_SUCCESS);
1407: }
1409: /*@C
1410: DMStagSetUniformCoordinatesProduct - create uniform coordinates, as a product of 1D arrays
1412: Set the coordinate `DM` to be a `DMPRODUCT` of 1D `DMSTAG` objects, each of which have a coordinate `DM` (also a 1d `DMSTAG`) holding uniform coordinates.
1414: Collective
1416: Input Parameters:
1417: + dm - the `DMSTAG` object
1418: . xmin - minimum global coordinate value in the x direction
1419: . xmax - maximum global coordinate values in the x direction
1420: . ymin - minimum global coordinate value in the y direction
1421: . ymax - maximum global coordinate value in the y direction
1422: . zmin - minimum global coordinate value in the z direction
1423: - zmax - maximum global coordinate value in the z direction
1425: Level: intermediate
1427: Notes:
1428: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1430: The per-dimension 1-dimensional `DMSTAG` objects that comprise the product
1431: always have active 0-cells (vertices, element boundaries) and 1-cells
1432: (element centers).
1434: See the manual page for `DMStagSetUniformCoordinates()` for information on how
1435: coordinates for dummy cells outside the physical domain boundary are populated.
1437: .seealso: [](ch_stag), `DMSTAG`, `DMPRODUCT`, `DMStagSetUniformCoordinates()`, `DMStagSetUniformCoordinatesExplicit()`, `DMStagSetCoordinateDMType()`
1438: @*/
1439: PetscErrorCode DMStagSetUniformCoordinatesProduct(DM dm, PetscReal xmin, PetscReal xmax, PetscReal ymin, PetscReal ymax, PetscReal zmin, PetscReal zmax)
1440: {
1441: DM_Stag *const stag = (DM_Stag *)dm->data;
1442: DM dmc;
1443: PetscInt dim, d, dof0, dof1;
1444: PetscBool flg;
1446: PetscFunctionBegin;
1448: PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This function must be called after DMSetUp()");
1449: PetscCall(PetscStrcmp(stag->coordinateDMType, DMPRODUCT, &flg));
1450: PetscCheck(!stag->coordinateDMType || flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Refusing to change an already-set DM coordinate type");
1451: PetscCall(DMStagSetCoordinateDMType(dm, DMPRODUCT));
1452: PetscCall(DMGetCoordinateDM(dm, &dmc));
1453: PetscCall(DMGetDimension(dm, &dim));
1455: /* Create 1D sub-DMs, living on subcommunicators.
1456: Always include both vertex and element dof, regardless of the active strata of the DMStag */
1457: dof0 = 1;
1458: dof1 = 1;
1460: for (d = 0; d < dim; ++d) {
1461: DM subdm;
1462: MPI_Comm subcomm;
1463: PetscMPIInt color;
1464: const PetscMPIInt key = 0; /* let existing rank break ties */
1466: /* Choose colors based on position in the plane orthogonal to this dim, and split */
1467: switch (d) {
1468: case 0:
1469: color = (dim > 1 ? stag->rank[1] : 0) + (dim > 2 ? stag->nRanks[1] * stag->rank[2] : 0);
1470: break;
1471: case 1:
1472: color = stag->rank[0] + (dim > 2 ? stag->nRanks[0] * stag->rank[2] : 0);
1473: break;
1474: case 2:
1475: color = stag->rank[0] + stag->nRanks[0] * stag->rank[1];
1476: break;
1477: default:
1478: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension index %" PetscInt_FMT, d);
1479: }
1480: PetscCallMPI(MPI_Comm_split(PetscObjectComm((PetscObject)dm), color, key, &subcomm));
1482: /* Create sub-DMs living on these new communicators (which are destroyed by DMProduct) */
1483: PetscCall(DMStagCreate1d(subcomm, stag->boundaryType[d], stag->N[d], dof0, dof1, stag->stencilType, stag->stencilWidth, stag->l[d], &subdm));
1484: /* Copy vector and matrix type information from parent DM */
1485: PetscCall(DMSetVecType(subdm, dm->vectype));
1486: PetscCall(DMSetMatType(subdm, dm->mattype));
1487: PetscCall(DMSetUp(subdm));
1488: switch (d) {
1489: case 0:
1490: PetscCall(DMStagSetUniformCoordinatesExplicit(subdm, xmin, xmax, 0, 0, 0, 0));
1491: break;
1492: case 1:
1493: PetscCall(DMStagSetUniformCoordinatesExplicit(subdm, ymin, ymax, 0, 0, 0, 0));
1494: break;
1495: case 2:
1496: PetscCall(DMStagSetUniformCoordinatesExplicit(subdm, zmin, zmax, 0, 0, 0, 0));
1497: break;
1498: default:
1499: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dimension index %" PetscInt_FMT, d);
1500: }
1501: PetscCall(DMProductSetDM(dmc, d, subdm));
1502: PetscCall(DMProductSetDimensionIndex(dmc, d, 0));
1503: PetscCall(DMDestroy(&subdm));
1504: PetscCallMPI(MPI_Comm_free(&subcomm));
1505: }
1506: PetscCall(DMCoarsenHookRemove(dm, DMRestrictHook_Coordinates, NULL, NULL));
1507: PetscFunctionReturn(PETSC_SUCCESS);
1508: }
1510: /*@C
1511: DMStagVecGetArray - get access to local array
1513: Logically Collective
1515: Input Parameters:
1516: + dm - the `DMSTAG` object
1517: - vec - the `Vec` object
1519: Output Parameter:
1520: . array - the array
1522: Level: beginner
1524: Note:
1525: This function returns a (dim+1)-dimensional array for a dim-dimensional
1526: `DMSTAG`.
1528: The first 1-3 dimensions indicate an element in the global
1529: numbering, using the standard C ordering.
1531: The final dimension in this array corresponds to a degree
1532: of freedom with respect to this element, for example corresponding to
1533: the element or one of its neighboring faces, edges, or vertices.
1535: For example, for a 3D `DMSTAG`, indexing is `array[k][j][i][idx]`, where `k` is the
1536: index in the z-direction, `j` is the index in the y-direction, and `i` is the
1537: index in the x-direction.
1539: `idx` is obtained with `DMStagGetLocationSlot()`, since the correct offset
1540: into the $(d+1)$-dimensional C array for a $d$-dimensional `DMSTAG` depends on the grid size and the number
1541: of DOF stored at each location.
1543: `DMStagVecRestoreArray()` must be called, once finished with the array
1545: .seealso: [](ch_stag), `DMSTAG`, `DMStagVecGetArrayRead()`, `DMStagGetLocationSlot()`, `DMGetLocalVector()`, `DMCreateLocalVector()`, `DMGetGlobalVector()`, `DMCreateGlobalVector()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`
1546: @*/
1547: PetscErrorCode DMStagVecGetArray(DM dm, Vec vec, void *array)
1548: {
1549: DM_Stag *const stag = (DM_Stag *)dm->data;
1550: PetscInt dim;
1551: PetscInt nLocal;
1553: PetscFunctionBegin;
1556: PetscCall(DMGetDimension(dm, &dim));
1557: PetscCall(VecGetLocalSize(vec, &nLocal));
1558: PetscCheck(nLocal == stag->entriesGhost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMStag local size %" PetscInt_FMT, nLocal, stag->entriesGhost);
1559: switch (dim) {
1560: case 1:
1561: PetscCall(VecGetArray2d(vec, stag->nGhost[0], stag->entriesPerElement, stag->startGhost[0], 0, (PetscScalar ***)array));
1562: break;
1563: case 2:
1564: PetscCall(VecGetArray3d(vec, stag->nGhost[1], stag->nGhost[0], stag->entriesPerElement, stag->startGhost[1], stag->startGhost[0], 0, (PetscScalar ****)array));
1565: break;
1566: case 3:
1567: PetscCall(VecGetArray4d(vec, stag->nGhost[2], stag->nGhost[1], stag->nGhost[0], stag->entriesPerElement, stag->startGhost[2], stag->startGhost[1], stag->startGhost[0], 0, (PetscScalar *****)array));
1568: break;
1569: default:
1570: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %" PetscInt_FMT, dim);
1571: }
1572: PetscFunctionReturn(PETSC_SUCCESS);
1573: }
1575: /*@C
1576: DMStagVecGetArrayRead - get read-only access to a local array
1578: Logically Collective
1580: See the man page for `DMStagVecGetArray()` for more information.
1582: Input Parameters:
1583: + dm - the `DMSTAG` object
1584: - vec - the `Vec` object
1586: Output Parameter:
1587: . array - the read-only array
1589: Level: beginner
1591: Note:
1592: `DMStagVecRestoreArrayRead()` must be called, once finished with the array
1594: .seealso: [](ch_stag), `DMSTAG`, `DMStagGetLocationSlot()`, `DMGetLocalVector()`, `DMCreateLocalVector()`, `DMGetGlobalVector()`, `DMCreateGlobalVector()`, `DMDAVecGetArrayRead()`, `DMDAVecGetArrayDOFRead()`
1595: @*/
1596: PetscErrorCode DMStagVecGetArrayRead(DM dm, Vec vec, void *array)
1597: {
1598: DM_Stag *const stag = (DM_Stag *)dm->data;
1599: PetscInt dim;
1600: PetscInt nLocal;
1602: PetscFunctionBegin;
1605: PetscCall(DMGetDimension(dm, &dim));
1606: PetscCall(VecGetLocalSize(vec, &nLocal));
1607: PetscCheck(nLocal == stag->entriesGhost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMSTAG local size %" PetscInt_FMT, nLocal, stag->entriesGhost);
1608: switch (dim) {
1609: case 1:
1610: PetscCall(VecGetArray2dRead(vec, stag->nGhost[0], stag->entriesPerElement, stag->startGhost[0], 0, (PetscScalar ***)array));
1611: break;
1612: case 2:
1613: PetscCall(VecGetArray3dRead(vec, stag->nGhost[1], stag->nGhost[0], stag->entriesPerElement, stag->startGhost[1], stag->startGhost[0], 0, (PetscScalar ****)array));
1614: break;
1615: case 3:
1616: PetscCall(VecGetArray4dRead(vec, stag->nGhost[2], stag->nGhost[1], stag->nGhost[0], stag->entriesPerElement, stag->startGhost[2], stag->startGhost[1], stag->startGhost[0], 0, (PetscScalar *****)array));
1617: break;
1618: default:
1619: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %" PetscInt_FMT, dim);
1620: }
1621: PetscFunctionReturn(PETSC_SUCCESS);
1622: }
1624: /*@C
1625: DMStagVecRestoreArray - restore access to a raw array
1627: Logically Collective
1629: Input Parameters:
1630: + dm - the `DMSTAG` object
1631: - vec - the `Vec` object
1633: Output Parameter:
1634: . array - the array
1636: Level: beginner
1638: .seealso: [](ch_stag), `DMSTAG`, `DMStagVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
1639: @*/
1640: PetscErrorCode DMStagVecRestoreArray(DM dm, Vec vec, void *array)
1641: {
1642: DM_Stag *const stag = (DM_Stag *)dm->data;
1643: PetscInt dim;
1644: PetscInt nLocal;
1646: PetscFunctionBegin;
1649: PetscCall(DMGetDimension(dm, &dim));
1650: PetscCall(VecGetLocalSize(vec, &nLocal));
1651: PetscCheck(nLocal == stag->entriesGhost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMSTAG local size %" PetscInt_FMT, nLocal, stag->entriesGhost);
1652: switch (dim) {
1653: case 1:
1654: PetscCall(VecRestoreArray2d(vec, stag->nGhost[0], stag->entriesPerElement, stag->startGhost[0], 0, (PetscScalar ***)array));
1655: break;
1656: case 2:
1657: PetscCall(VecRestoreArray3d(vec, stag->nGhost[1], stag->nGhost[0], stag->entriesPerElement, stag->startGhost[1], stag->startGhost[0], 0, (PetscScalar ****)array));
1658: break;
1659: case 3:
1660: PetscCall(VecRestoreArray4d(vec, stag->nGhost[2], stag->nGhost[1], stag->nGhost[0], stag->entriesPerElement, stag->startGhost[2], stag->startGhost[1], stag->startGhost[0], 0, (PetscScalar *****)array));
1661: break;
1662: default:
1663: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %" PetscInt_FMT, dim);
1664: }
1665: PetscFunctionReturn(PETSC_SUCCESS);
1666: }
1668: /*@C
1669: DMStagVecRestoreArrayRead - restore read-only access to a raw array
1671: Logically Collective
1673: Input Parameters:
1674: + dm - the `DMSTAG` object
1675: - vec - the Vec object
1677: Output Parameter:
1678: . array - the read-only array
1680: Level: beginner
1682: .seealso: [](ch_stag), `DMSTAG`, `DMStagVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecRestoreArrayDOFRead()`
1683: @*/
1684: PetscErrorCode DMStagVecRestoreArrayRead(DM dm, Vec vec, void *array)
1685: {
1686: DM_Stag *const stag = (DM_Stag *)dm->data;
1687: PetscInt dim;
1688: PetscInt nLocal;
1690: PetscFunctionBegin;
1693: PetscCall(DMGetDimension(dm, &dim));
1694: PetscCall(VecGetLocalSize(vec, &nLocal));
1695: PetscCheck(nLocal == stag->entriesGhost, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMSTAG local size %" PetscInt_FMT, nLocal, stag->entriesGhost);
1696: switch (dim) {
1697: case 1:
1698: PetscCall(VecRestoreArray2dRead(vec, stag->nGhost[0], stag->entriesPerElement, stag->startGhost[0], 0, (PetscScalar ***)array));
1699: break;
1700: case 2:
1701: PetscCall(VecRestoreArray3dRead(vec, stag->nGhost[1], stag->nGhost[0], stag->entriesPerElement, stag->startGhost[1], stag->startGhost[0], 0, (PetscScalar ****)array));
1702: break;
1703: case 3:
1704: PetscCall(VecRestoreArray4dRead(vec, stag->nGhost[2], stag->nGhost[1], stag->nGhost[0], stag->entriesPerElement, stag->startGhost[2], stag->startGhost[1], stag->startGhost[0], 0, (PetscScalar *****)array));
1705: break;
1706: default:
1707: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %" PetscInt_FMT, dim);
1708: }
1709: PetscFunctionReturn(PETSC_SUCCESS);
1710: }