Actual source code: stagutils.c
petsc-3.13.0 2020-03-29
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>
4: /*@C
5: DMStagGetBoundaryTypes - get boundary types
7: Not Collective
9: Input Parameter:
10: . dm - the DMStag object
12: Output Parameters:
13: . boundaryTypeX,boundaryTypeY,boundaryTypeZ - boundary types
15: Level: intermediate
17: .seealso: DMSTAG, DMDAGetBoundaryTypes()
18: @*/
19: PetscErrorCode DMStagGetBoundaryTypes(DM dm,DMBoundaryType *boundaryTypeX,DMBoundaryType *boundaryTypeY,DMBoundaryType *boundaryTypeZ)
20: {
21: PetscErrorCode ierr;
22: const DM_Stag * const stag = (DM_Stag*)dm->data;
23: PetscInt dim;
27: DMGetDimension(dm,&dim);
28: if (boundaryTypeX ) *boundaryTypeX = stag->boundaryType[0];
29: if (boundaryTypeY && dim > 1) *boundaryTypeY = stag->boundaryType[1];
30: if (boundaryTypeZ && dim > 2) *boundaryTypeZ = stag->boundaryType[2];
31: return(0);
32: }
34: static PetscErrorCode DMStagGetProductCoordinateArrays_Private(DM dm,void* arrX,void* arrY,void* arrZ,PetscBool read)
35: {
37: PetscInt dim,d,dofCheck[DMSTAG_MAX_STRATA],s;
38: DM dmCoord;
39: void* arr[DMSTAG_MAX_DIM];
43: DMGetDimension(dm,&dim);
44: if (dim > DMSTAG_MAX_DIM) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not implemented for %D dimensions",dim);
45: arr[0] = arrX; arr[1] = arrY; arr[2] = arrZ;
46: DMGetCoordinateDM(dm,&dmCoord);
47: if (!dmCoord) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM does not have a coordinate DM");
48: {
49: PetscBool isProduct;
50: DMType dmType;
51: DMGetType(dmCoord,&dmType);
52: PetscStrcmp(DMPRODUCT,dmType,&isProduct);
53: if (!isProduct) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate DM is not of type DMPRODUCT");
54: }
55: for (s=0; s<DMSTAG_MAX_STRATA; ++s) dofCheck[s] = 0;
56: for (d=0; d<dim; ++d) {
57: DM subDM;
58: DMType dmType;
59: PetscBool isStag;
60: PetscInt dof[DMSTAG_MAX_STRATA],subDim;
61: Vec coord1d_local;
63: DMProductGetDM(dmCoord,d,&subDM);
64: if (!subDM) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate DM is missing sub DM %D",d);
65: DMGetDimension(subDM,&subDim);
66: if (subDim != 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate sub-DM is not of dimension 1");
67: DMGetType(subDM,&dmType);
68: PetscStrcmp(DMSTAG,dmType,&isStag);
69: if (!isStag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate sub-DM is not of type DMSTAG");
70: DMStagGetDOF(subDM,&dof[0],&dof[1],&dof[2],&dof[3]);
71: if (d == 0) {
72: for (s=0; s<DMSTAG_MAX_STRATA; ++s) dofCheck[s] = dof[s];
73: } else {
74: for (s=0; s<DMSTAG_MAX_STRATA; ++s) {
75: if (dofCheck[s] != dof[s]) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate sub-DMs have different dofs");
76: }
77: }
78: DMGetCoordinatesLocal(subDM,&coord1d_local);
79: if (read) {
80: DMStagVecGetArrayRead(subDM,coord1d_local,arr[d]);
81: } else {
82: DMStagVecGetArray(subDM,coord1d_local,arr[d]);
83: }
84: }
85: return(0);
86: }
88: /*@C
89: DMStagGetProductCoordinateArrays - extract local product coordinate arrays, one per dimension
91: Logically Collective
93: A high-level helper function to quickly extract local coordinate arrays.
95: Note that 2-dimensional arrays are returned. See
96: DMStagVecGetArray(), which is called internally to produce these arrays
97: representing coordinates on elements and vertices (element boundaries)
98: for a 1-dimensional DMStag in each coordinate direction.
100: One should use DMStagGetProductCoordinateSlot() to determine appropriate
101: indices for the second dimension in these returned arrays. This function
102: checks that the coordinate array is a suitable product of 1-dimensional
103: DMStag objects.
105: Input Parameter:
106: . dm - the DMStag object
108: Output Parameters:
109: . arrX,arrY,arrZ - local 1D coordinate arrays
111: Level: intermediate
113: .seealso: DMSTAG, DMPRODUCT, DMStagGetProductCoordinateArraysRead(), DMStagSetUniformCoordinates(), DMStagSetUniformCoordinatesProduct(), DMStagGetProductCoordinateLocationSlot()
114: @*/
115: PetscErrorCode DMStagGetProductCoordinateArrays(DM dm,void* arrX,void* arrY,void* arrZ)
116: {
120: DMStagGetProductCoordinateArrays_Private(dm,arrX,arrY,arrZ,PETSC_FALSE);
121: return(0);
122: }
124: /*@C
125: DMStagGetProductCoordinateArraysRead - extract product coordinate arrays, read-only
127: Logically Collective
129: See the man page for DMStagGetProductCoordinateArrays() for more information.
131: Input Parameter:
132: . dm - the DMStag object
134: Output Parameters:
135: . arrX,arrY,arrZ - local 1D coordinate arrays
137: Level: intermediate
139: .seealso: DMSTAG, DMPRODUCT, DMStagGetProductCoordinateArrays(), DMStagSetUniformCoordinates(), DMStagSetUniformCoordinatesProduct(), DMStagGetProductCoordinateLocationSlot()
140: @*/
141: PetscErrorCode DMStagGetProductCoordinateArraysRead(DM dm,void* arrX,void* arrY,void* arrZ)
142: {
146: DMStagGetProductCoordinateArrays_Private(dm,arrX,arrY,arrZ,PETSC_TRUE);
147: return(0);
148: }
150: /*@C
151: DMStagGetProductCoordinateLocationSlot - get slot for use with local product coordinate arrays
153: Not Collective
155: High-level helper function to get slot indices for 1D coordinate DMs,
156: for use with DMStagGetProductCoordinateArrays() and related functions.
158: Input Parameters:
159: + dm - the DMStag object
160: - loc - the grid location
162: Output Parameter:
163: . slot - the index to use in local arrays
165: Notes:
166: Checks that the coordinates are actually set up so that using the
167: slots from the first 1d coordinate sub-DM is valid for all the 1D coordinate sub-DMs.
169: Level: intermediate
171: .seealso: DMSTAG, DMPRODUCT, DMStagGetProductCoordinateArrays(), DMStagGetProductCoordinateArraysRead(), DMStagSetUniformCoordinates()
172: @*/
173: PETSC_EXTERN PetscErrorCode DMStagGetProductCoordinateLocationSlot(DM dm,DMStagStencilLocation loc,PetscInt *slot)
174: {
176: DM dmCoord;
177: PetscInt dim,dofCheck[DMSTAG_MAX_STRATA],s,d;
181: DMGetDimension(dm,&dim);
182: DMGetCoordinateDM(dm,&dmCoord);
183: if (!dmCoord) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM does not have a coordinate DM");
184: {
185: PetscBool isProduct;
186: DMType dmType;
187: DMGetType(dmCoord,&dmType);
188: PetscStrcmp(DMPRODUCT,dmType,&isProduct);
189: if (!isProduct) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate DM is not of type DMPRODUCT");
190: }
191: for (s=0; s<DMSTAG_MAX_STRATA; ++s) dofCheck[s] = 0;
192: for (d=0; d<dim; ++d) {
193: DM subDM;
194: DMType dmType;
195: PetscBool isStag;
196: PetscInt dof[DMSTAG_MAX_STRATA],subDim;
197: DMProductGetDM(dmCoord,d,&subDM);
198: if (!subDM) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate DM is missing sub DM %D",d);
199: DMGetDimension(subDM,&subDim);
200: if (subDim != 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate sub-DM is not of dimension 1");
201: DMGetType(subDM,&dmType);
202: PetscStrcmp(DMSTAG,dmType,&isStag);
203: if (!isStag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate sub-DM is not of type DMSTAG");
204: DMStagGetDOF(subDM,&dof[0],&dof[1],&dof[2],&dof[3]);
205: if (d == 0) {
206: const PetscInt component = 0;
207: for (s=0; s<DMSTAG_MAX_STRATA; ++s) dofCheck[s] = dof[s];
208: DMStagGetLocationSlot(subDM,loc,component,slot);
209: } else {
210: for (s=0; s<DMSTAG_MAX_STRATA; ++s) {
211: if (dofCheck[s] != dof[s]) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate sub-DMs have different dofs");
212: }
213: }
214: }
215: return(0);
216: }
218: /*@C
219: DMStagGetCorners - return global element indices of the local region (excluding ghost points)
221: Not Collective
223: Input Parameter:
224: . dm - the DMStag object
226: Output Parameters:
227: + x,y,z - starting element indices in each direction
228: . m,n,p - element widths in each direction
229: - nExtrax,nExtray,nExtraz - number of extra partial elements in each direction.
231: Notes:
232: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
234: The number of extra partial elements is either 1 or 0.
235: The value is 1 on right, top, and front non-periodic domain ("physical") boundaries,
236: in the x, y, and z directions respectively, and otherwise 0.
238: Level: beginner
240: .seealso: DMSTAG, DMStagGetGhostCorners(), DMDAGetCorners()
241: @*/
242: PetscErrorCode DMStagGetCorners(DM dm,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p,PetscInt *nExtrax,PetscInt *nExtray,PetscInt *nExtraz)
243: {
244: const DM_Stag * const stag = (DM_Stag*)dm->data;
248: if (x) *x = stag->start[0];
249: if (y) *y = stag->start[1];
250: if (z) *z = stag->start[2];
251: if (m) *m = stag->n[0];
252: if (n) *n = stag->n[1];
253: if (p) *p = stag->n[2];
254: if (nExtrax) *nExtrax = stag->boundaryType[0] != DM_BOUNDARY_PERIODIC && stag->lastRank[0] ? 1 : 0;
255: if (nExtray) *nExtray = stag->boundaryType[1] != DM_BOUNDARY_PERIODIC && stag->lastRank[1] ? 1 : 0;
256: if (nExtraz) *nExtraz = stag->boundaryType[2] != DM_BOUNDARY_PERIODIC && stag->lastRank[2] ? 1 : 0;
257: return(0);
258: }
260: /*@C
261: DMStagGetDOF - get number of DOF associated with each stratum of the grid
263: Not Collective
265: Input Parameter:
266: . dm - the DMStag object
268: Output Parameters:
269: + dof0 - the number of points per 0-cell (vertex/node)
270: . dof1 - the number of points per 1-cell (element in 1D, edge in 2D and 3D)
271: . dof2 - the number of points per 2-cell (element in 2D, face in 3D)
272: - dof3 - the number of points per 3-cell (element in 3D)
274: Level: beginner
276: .seealso: DMSTAG, DMStagGetCorners(), DMStagGetGhostCorners(), DMStagGetGlobalSizes(), DMStagGetStencilWidth(), DMStagGetBoundaryTypes(), DMStagGetLocationDof(), DMDAGetDof()
277: @*/
278: PetscErrorCode DMStagGetDOF(DM dm,PetscInt *dof0,PetscInt *dof1,PetscInt *dof2,PetscInt *dof3)
279: {
280: const DM_Stag * const stag = (DM_Stag*)dm->data;
284: if (dof0) *dof0 = stag->dof[0];
285: if (dof1) *dof1 = stag->dof[1];
286: if (dof2) *dof2 = stag->dof[2];
287: if (dof3) *dof3 = stag->dof[3];
288: return(0);
289: }
291: /*@C
292: DMStagGetGhostCorners - return global element indices of the local region, including ghost points
294: Not Collective
296: Input Argument:
297: . dm - the DMStag object
299: Output Arguments:
300: + x,y,z - starting element indices in each direction
301: - m,n,p - element widths in each direction
303: Notes:
304: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
306: Level: beginner
308: .seealso: DMSTAG, DMStagGetCorners(), DMDAGetGhostCorners()
309: @*/
310: PetscErrorCode DMStagGetGhostCorners(DM dm,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p)
311: {
312: const DM_Stag * const stag = (DM_Stag*)dm->data;
316: if (x) *x = stag->startGhost[0];
317: if (y) *y = stag->startGhost[1];
318: if (z) *z = stag->startGhost[2];
319: if (m) *m = stag->nGhost[0];
320: if (n) *n = stag->nGhost[1];
321: if (p) *p = stag->nGhost[2];
322: return(0);
323: }
325: /*@C
326: DMStagGetGlobalSizes - get global element counts
328: Not Collective
330: Input Parameter:
331: . dm - the DMStag object
333: Output Parameters:
334: . M,N,P - global element counts in each direction
336: Notes:
337: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
339: Level: beginner
341: .seealso: DMSTAG, DMStagGetLocalSizes(), DMDAGetInfo()
342: @*/
343: PetscErrorCode DMStagGetGlobalSizes(DM dm,PetscInt* M,PetscInt* N,PetscInt* P)
344: {
345: const DM_Stag * const stag = (DM_Stag*)dm->data;
349: if (M) *M = stag->N[0];
350: if (N) *N = stag->N[1];
351: if (P) *P = stag->N[2];
352: return(0);
353: }
355: /*@C
356: DMStagGetIsFirstRank - get boolean value for whether this rank is first in each direction in the rank grid
358: Not Collective
360: Input Parameter:
361: . dm - the DMStag object
363: Output Parameters:
364: . isFirstRank0,isFirstRank1,isFirstRank2 - whether this rank is first in each direction
366: Level: intermediate
368: Notes:
369: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
371: .seealso: DMSTAG, DMStagGetIsLastRank()
372: @*/
373: PetscErrorCode DMStagGetIsFirstRank(DM dm,PetscBool *isFirstRank0,PetscBool *isFirstRank1,PetscBool *isFirstRank2)
374: {
375: const DM_Stag * const stag = (DM_Stag*)dm->data;
379: if (isFirstRank0) *isFirstRank0 = stag->firstRank[0];
380: if (isFirstRank1) *isFirstRank1 = stag->firstRank[1];
381: if (isFirstRank2) *isFirstRank2 = stag->firstRank[2];
382: return(0);
383: }
385: /*@C
386: DMStagGetIsLastRank - get boolean value for whether this rank is last in each direction in the rank grid
388: Not Collective
390: Input Parameter:
391: . dm - the DMStag object
393: Output Parameters:
394: . isLastRank0,isLastRank1,isLastRank2 - whether this rank is last in each direction
396: Level: intermediate
398: Notes:
399: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
400: Level: intermediate
402: .seealso: DMSTAG, DMStagGetIsFirstRank()
403: @*/
404: PetscErrorCode DMStagGetIsLastRank(DM dm,PetscBool *isLastRank0,PetscBool *isLastRank1,PetscBool *isLastRank2)
405: {
406: const DM_Stag * const stag = (DM_Stag*)dm->data;
410: if (isLastRank0) *isLastRank0 = stag->lastRank[0];
411: if (isLastRank1) *isLastRank1 = stag->lastRank[1];
412: if (isLastRank2) *isLastRank2 = stag->lastRank[2];
413: return(0);
414: }
416: /*@C
417: DMStagGetLocalSizes - get local elementwise sizes
419: Not Collective
421: Input Parameter:
422: . dm - the DMStag object
424: Output Parameters:
425: . m,n,p - local element counts (excluding ghosts) in each direction
427: Notes:
428: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
429: Level: intermediate
431: Level: beginner
433: .seealso: DMSTAG, DMStagGetGlobalSizes(), DMStagGetDOF(), DMStagGetNumRanks(), DMDAGetLocalInfo()
434: @*/
435: PetscErrorCode DMStagGetLocalSizes(DM dm,PetscInt* m,PetscInt* n,PetscInt* p)
436: {
437: const DM_Stag * const stag = (DM_Stag*)dm->data;
441: if (m) *m = stag->n[0];
442: if (n) *n = stag->n[1];
443: if (p) *p = stag->n[2];
444: return(0);
445: }
447: /*@C
448: DMStagGetNumRanks - get number of ranks in each direction in the global grid decomposition
450: Not Collective
452: Input Parameter:
453: . dm - the DMStag object
455: Output Parameters:
456: . nRanks0,nRanks1,nRanks2 - number of ranks in each direction in the grid decomposition
458: Notes:
459: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
460: Level: intermediate
462: Level: beginner
464: .seealso: DMSTAG, DMStagGetGlobalSizes(), DMStagGetLocalSize(), DMStagSetNumRanks(), DMDAGetInfo()
465: @*/
466: PetscErrorCode DMStagGetNumRanks(DM dm,PetscInt *nRanks0,PetscInt *nRanks1,PetscInt *nRanks2)
467: {
468: const DM_Stag * const stag = (DM_Stag*)dm->data;
472: if (nRanks0) *nRanks0 = stag->nRanks[0];
473: if (nRanks1) *nRanks1 = stag->nRanks[1];
474: if (nRanks2) *nRanks2 = stag->nRanks[2];
475: return(0);
476: }
478: /*@C
479: DMStagGetEntriesPerElement - get number of entries per element in the local representation
481: Not Collective
483: Input Parameter:
484: . dm - the DMStag object
486: Output Parameters:
487: . entriesPerElement - number of entries associated with each element in the local representation
489: Notes:
490: This is the natural block size for most local operations. In 1D it is equal to dof0 + dof1,
491: in 2D it is equal to dof0 + 2*dof1 + dof2, and in 3D it is equal to dof0 + 3*dof1 + 3*dof2 + dof3
493: Level: developer
495: .seealso: DMSTAG, DMStagGetDOF()
496: @*/
497: PetscErrorCode DMStagGetEntriesPerElement(DM dm,PetscInt *entriesPerElement)
498: {
499: const DM_Stag * const stag = (DM_Stag*)dm->data;
503: if (entriesPerElement) *entriesPerElement = stag->entriesPerElement;
504: return(0);
505: }
507: /*@C
508: DMStagGetStencilType - get elementwise ghost/halo stencil type
510: Not Collective
512: Input Parameter:
513: . dm - the DMStag object
515: Output Parameter:
516: . stencilType - the elementwise ghost stencil type: DMSTAG_STENCIL_BOX, DMSTAG_STENCIL_STAR, or DMSTAG_STENCIL_NONE
518: Level: beginner
520: .seealso: DMSTAG, DMStagSetStencilType(), DMStagStencilType, DMDAGetInfo()
521: @*/
522: PetscErrorCode DMStagGetStencilType(DM dm,DMStagStencilType *stencilType)
523: {
524: DM_Stag * const stag = (DM_Stag*)dm->data;
528: *stencilType = stag->stencilType;
529: return(0);
530: }
532: /*@C
533: DMStagGetStencilWidth - get elementwise stencil width
535: Not Collective
537: Input Parameter:
538: . dm - the DMStag object
540: Output Parameters:
541: . stencilWidth - stencil/halo/ghost width in elements
543: Level: beginner
545: .seealso: DMSTAG, DMDAGetStencilWidth(), DMDAGetInfo()
546: @*/
547: PetscErrorCode DMStagGetStencilWidth(DM dm,PetscInt *stencilWidth)
548: {
549: const DM_Stag * const stag = (DM_Stag*)dm->data;
553: if (stencilWidth) *stencilWidth = stag->stencilWidth;
554: return(0);
555: }
557: /*@C
558: DMStagGetOwnershipRanges - get elements per rank in each direction
560: Not Collective
562: Input Parameter:
563: . dm - the DMStag object
565: Output Parameters:
566: + lx - ownership along x direction (optional)
567: . ly - ownership along y direction (optional)
568: - lz - ownership along z direction (optional)
570: Notes:
571: These correspond to the optional final arguments passed to DMStagCreate1d(), DMStagCreate2d(), and DMStagCreate3d().
573: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
575: In C you should not free these arrays, nor change the values in them.
576: They will only have valid values while the DMStag they came from still exists (has not been destroyed).
578: Level: intermediate
580: .seealso: DMSTAG, DMStagSetGlobalSizes(), DMStagSetOwnershipRanges(), DMStagCreate1d(), DMStagCreate2d(), DMStagCreate3d(), DMDAGetOwnershipRanges()
581: @*/
582: PetscErrorCode DMStagGetOwnershipRanges(DM dm,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
583: {
584: const DM_Stag * const stag = (DM_Stag*)dm->data;
588: if (lx) *lx = stag->l[0];
589: if (ly) *ly = stag->l[1];
590: if (lz) *lz = stag->l[2];
591: return(0);
592: }
594: /*@C
595: DMStagCreateCompatibleDMStag - create a compatible DMStag with different dof/stratum
597: Collective
599: Input Parameters:
600: + dm - the DMStag object
601: - dof0,dof1,dof2,dof3 - number of dof on each stratum in the new DMStag
603: Output Parameters:
604: . newdm - the new, compatible DMStag
606: Notes:
607: Dof supplied for strata too big for the dimension are ignored; these may be set to 0.
608: In contrast to DMDACreateCompatibleDMDA(), coordinates are not reused.
610: Level: intermediate
612: .seealso: DMSTAG, DMDACreateCompatibleDMDA(), DMGetCompatibility(), DMStagMigrateVec()
613: @*/
614: PetscErrorCode DMStagCreateCompatibleDMStag(DM dm,PetscInt dof0,PetscInt dof1,PetscInt dof2,PetscInt dof3,DM *newdm)
615: {
616: PetscErrorCode ierr;
620: DMStagDuplicateWithoutSetup(dm,PetscObjectComm((PetscObject)dm),newdm);
621: DMStagSetDOF(*newdm,dof0,dof1,dof2,dof3);
622: DMSetUp(*newdm);
623: return(0);
624: }
626: /*@C
627: DMStagGetLocationSlot - get index to use in accessing raw local arrays
629: Not Collective
631: Input Parameters:
632: + dm - the DMStag object
633: . loc - location relative to an element
634: - c - component
636: Output Parameter:
637: . slot - index to use
639: Notes:
640: Provides an appropriate index to use with DMStagVecGetArray() and friends.
641: This is required so that the user doesn't need to know about the ordering of
642: dof associated with each local element.
644: Level: beginner
646: .seealso: DMSTAG, DMStagVecGetArray(), DMStagVecGetArrayRead(), DMStagGetDOF(), DMStagGetEntriesPerElement()
647: @*/
648: PetscErrorCode DMStagGetLocationSlot(DM dm,DMStagStencilLocation loc,PetscInt c,PetscInt *slot)
649: {
650: DM_Stag * const stag = (DM_Stag*)dm->data;
654: #if defined(PETSC_USE_DEBUG)
655: {
657: PetscInt dof;
658: DMStagGetLocationDOF(dm,loc,&dof);
659: if (dof < 1) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Location %s has no dof attached",DMStagStencilLocations[loc]);
660: if (c > dof-1) SETERRQ3(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Supplied component number (%D) for location %s is too big (maximum %D)",c,DMStagStencilLocations[loc],dof-1);
661: }
662: #endif
663: *slot = stag->locationOffsets[loc] + c;
664: return(0);
665: }
667: /*@C
668: DMStagMigrateVec - transfer a vector associated with a DMStag to a vector associated with a compatible DMStag
670: Collective
672: Input Parameters:
673: + dm - the source DMStag object
674: . vec - the source vector, compatible with dm
675: . dmTo - the compatible destination DMStag object
676: - vecTo - the destination vector, compatible with dmTo
678: Notes:
679: Extra dof are ignored, and unfilled dof are zeroed.
680: Currently only implemented to migrate global vectors to global vectors.
682: Level: advanced
684: .seealso: DMSTAG, DMStagCreateCompatibleDMStag(), DMGetCompatibility(), DMStagVecSplitToDMDA()
685: @*/
686: PetscErrorCode DMStagMigrateVec(DM dm,Vec vec,DM dmTo,Vec vecTo)
687: {
688: PetscErrorCode ierr;
689: DM_Stag * const stag = (DM_Stag*)dm->data;
690: DM_Stag * const stagTo = (DM_Stag*)dmTo->data;
691: PetscInt nLocalTo,nLocal,dim,i,j,k;
692: PetscInt start[DMSTAG_MAX_DIM],startGhost[DMSTAG_MAX_DIM],n[DMSTAG_MAX_DIM],nExtra[DMSTAG_MAX_DIM],offset[DMSTAG_MAX_DIM];
693: Vec vecToLocal,vecLocal;
694: PetscBool compatible,compatibleSet;
695: const PetscScalar *arr;
696: PetscScalar *arrTo;
697: const PetscInt epe = stag->entriesPerElement;
698: const PetscInt epeTo = stagTo->entriesPerElement;
705: DMGetCompatibility(dm,dmTo,&compatible,&compatibleSet);
706: if (!compatibleSet || !compatible) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"DMStag objects must be shown to be compatible");
707: DMGetDimension(dm,&dim);
708: VecGetLocalSize(vec,&nLocal);
709: VecGetLocalSize(vecTo,&nLocalTo);
710: if (nLocal != stag->entries|| nLocalTo !=stagTo->entries) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Vector migration only implemented for global vector to global vector.");
711: DMStagGetCorners(dm,&start[0],&start[1],&start[2],&n[0],&n[1],&n[2],&nExtra[0],&nExtra[1],&nExtra[2]);
712: DMStagGetGhostCorners(dm,&startGhost[0],&startGhost[1],&startGhost[2],NULL,NULL,NULL);
713: for (i=0; i<DMSTAG_MAX_DIM; ++i) offset[i] = start[i]-startGhost[i];
715: /* Proceed by transferring to a local vector, copying, and transferring back to a global vector */
716: DMGetLocalVector(dm,&vecLocal);
717: DMGetLocalVector(dmTo,&vecToLocal);
718: DMGlobalToLocalBegin(dm,vec,INSERT_VALUES,vecLocal);
719: DMGlobalToLocalEnd(dm,vec,INSERT_VALUES,vecLocal);
720: VecGetArrayRead(vecLocal,&arr);
721: VecGetArray(vecToLocal,&arrTo);
722: /* Note that some superfluous copying of entries on partial dummy elements is done */
723: if (dim == 1) {
724: for (i=offset[0]; i<offset[0] + n[0] + nExtra[0]; ++i) {
725: PetscInt d = 0,dTo = 0,b = 0,bTo = 0;
726: PetscInt si;
727: for (si=0; si<2; ++si) {
728: b += stag->dof[si];
729: bTo += stagTo->dof[si];
730: for (; d < b && dTo < bTo; ++d,++dTo) arrTo[i*epeTo + dTo] = arr[i*epe + d];
731: for (; dTo < bTo ; ++dTo) arrTo[i*epeTo + dTo] = 0.0;
732: d = b;
733: }
734: }
735: } else if (dim == 2) {
736: const PetscInt epr = stag->nGhost[0] * epe;
737: const PetscInt eprTo = stagTo->nGhost[0] * epeTo;
738: for (j=offset[1]; j<offset[1] + n[1] + nExtra[1]; ++j) {
739: for (i=offset[0]; i<offset[0] + n[0] + nExtra[0]; ++i) {
740: const PetscInt base = j*epr + i*epe;
741: const PetscInt baseTo = j*eprTo + i*epeTo;
742: PetscInt d = 0,dTo = 0,b = 0,bTo = 0;
743: const PetscInt s[4] = {0,1,1,2}; /* Dimensions of points, in order */
744: PetscInt si;
745: for (si=0; si<4; ++si) {
746: b += stag->dof[s[si]];
747: bTo += stagTo->dof[s[si]];
748: for (; d < b && dTo < bTo; ++d,++dTo) arrTo[baseTo + dTo] = arr[base + d];
749: for (; dTo < bTo; ++dTo) arrTo[baseTo + dTo] = 0.0;
750: d = b;
751: }
752: }
753: }
754: } else if (dim == 3) {
755: const PetscInt epr = stag->nGhost[0] * epe;
756: const PetscInt eprTo = stagTo->nGhost[0] * epeTo;
757: const PetscInt epl = stag->nGhost[1] * epr;
758: const PetscInt eplTo = stagTo->nGhost[1] * eprTo;
759: for (k=offset[2]; k<offset[2] + n[2] + nExtra[2]; ++k) {
760: for (j=offset[1]; j<offset[1] + n[1] + nExtra[1]; ++j) {
761: for (i=offset[0]; i<offset[0] + n[0] + nExtra[0]; ++i) {
762: PetscInt d = 0,dTo = 0,b = 0,bTo = 0;
763: const PetscInt base = k*epl + j*epr + i*epe;
764: const PetscInt baseTo = k*eplTo + j*eprTo + i*epeTo;
765: const PetscInt s[8] = {0,1,1,2,1,2,2,3}; /* dimensions of points, in order */
766: PetscInt is;
767: for (is=0; is<8; ++is) {
768: b += stag->dof[s[is]];
769: bTo += stagTo->dof[s[is]];
770: for (; d < b && dTo < bTo; ++d,++dTo) arrTo[baseTo + dTo] = arr[base + d];
771: for (; dTo < bTo; ++dTo) arrTo[baseTo + dTo] = 0.0;
772: d = b;
773: }
774: }
775: }
776: }
777: } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
778: VecRestoreArrayRead(vecLocal,&arr);
779: VecRestoreArray(vecToLocal,&arrTo);
780: DMRestoreLocalVector(dm,&vecLocal);
781: DMLocalToGlobalBegin(dmTo,vecToLocal,INSERT_VALUES,vecTo);
782: DMLocalToGlobalEnd(dmTo,vecToLocal,INSERT_VALUES,vecTo);
783: DMRestoreLocalVector(dmTo,&vecToLocal);
784: return(0);
785: }
787: /*@C
788: DMStagPopulateLocalToGlobalInjective - populate an internal 1-to-1 local-to-global map
790: Collective
792: Creates an internal object which explicitly maps a single local degree of
793: freedom to each global degree of freedom. This is used, if populated,
794: instead of SCATTER_REVERSE_LOCAL with the (1-to-many, in general)
795: global-to-local map, when DMLocalToGlobal() is called with INSERT_VALUES.
796: This allows usage, for example, even in the periodic, 1-rank case, where
797: the inverse of the global-to-local map, even when restricted to on-rank
798: communication, is non-injective. This is at the cost of storing an additional
799: VecScatter object inside each DMStag object.
801: Input Parameter:
802: . dm - the DMStag object
804: Notes:
805: In normal usage, library users shouldn't be concerned with this function,
806: as it is called during DMSetUp(), when required.
808: Returns immediately if the internal map is already populated.
810: Developer Notes:
811: This could, if desired, be moved up to a general DM routine. It would allow,
812: for example, DMDA to support DMLocalToGlobal() with INSERT_VALUES,
813: even in the single-rank periodic case.
815: Level: developer
817: .seealso: DMSTAG, DMLocalToGlobal(), VecScatter
818: @*/
819: PetscErrorCode DMStagPopulateLocalToGlobalInjective(DM dm)
820: {
821: PetscErrorCode ierr;
822: PetscInt dim;
823: DM_Stag * const stag = (DM_Stag*)dm->data;
827: if (stag->ltog_injective) return(0); /* Don't re-populate */
828: DMGetDimension(dm,&dim);
829: switch (dim) {
830: case 1: DMStagPopulateLocalToGlobalInjective_1d(dm); break;
831: case 2: DMStagPopulateLocalToGlobalInjective_2d(dm); break;
832: case 3: DMStagPopulateLocalToGlobalInjective_3d(dm); break;
833: default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
834: }
835: return(0);
836: }
838: static PetscErrorCode DMStagRestoreProductCoordinateArrays_Private(DM dm,void *arrX,void *arrY,void *arrZ,PetscBool read)
839: {
840: PetscErrorCode ierr;
841: PetscInt dim,d;
842: void* arr[DMSTAG_MAX_DIM];
843: DM dmCoord;
847: DMGetDimension(dm,&dim);
848: if (dim > DMSTAG_MAX_DIM) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not implemented for %D dimensions",dim);
849: arr[0] = arrX; arr[1] = arrY; arr[2] = arrZ;
850: DMGetCoordinateDM(dm,&dmCoord);
851: for (d=0; d<dim; ++d) {
852: DM subDM;
853: Vec coord1d_local;
855: DMProductGetDM(dmCoord,d,&subDM);
856: DMGetCoordinatesLocal(subDM,&coord1d_local);
857: if (read) {
858: DMStagVecRestoreArrayRead(subDM,coord1d_local,arr[d]);
859: } else {
860: DMStagVecRestoreArray(subDM,coord1d_local,arr[d]);
861: }
862: }
863: return(0);
864: }
866: /*@C
867: DMStagRestoreProductCoordinateArrays - restore local array access
869: Logically Collective
871: Input Parameter:
872: . dm - the DMStag object
874: Output Parameters:
875: . arrX,arrY,arrZ - local 1D coordinate arrays
877: Level: intermediate
879: Notes:
880: This function does not automatically perform a local->global scatter to populate global coordinates from the local coordinates. Thus, it may be required to explicitly perform these operations in some situations, as in the following partial example:
882: $ DMGetCoordinateDM(dm,&cdm);
883: $ for (d=0; d<3; ++d) {
884: $ DM subdm;
885: $ Vec coor,coor_local;
887: $ DMProductGetDM(cdm,d,&subdm);
888: $ DMGetCoordinates(subdm,&coor);
889: $ DMGetCoordinatesLocal(subdm,&coor_local);
890: $ DMLocalToGlobal(subdm,coor_local,INSERT_VALUES,coor);
891: $ PetscPrintf(PETSC_COMM_WORLD,"Coordinates dim %D:\n",d);
892: $ VecView(coor,PETSC_VIEWER_STDOUT_WORLD);
893: $ }
895: .seealso: DMSTAG, DMStagGetProductCoordinateArrays(), DMStagGetProductCoordinateArraysRead()
896: @*/
897: PetscErrorCode DMStagRestoreProductCoordinateArrays(DM dm,void *arrX,void *arrY,void *arrZ)
898: {
899: PetscErrorCode ierr;
902: DMStagRestoreProductCoordinateArrays_Private(dm,arrX,arrY,arrZ,PETSC_FALSE);
903: return(0);
904: }
906: /*@C
907: DMStagRestoreProductCoordinateArraysRead - restore local product array access, read-only
909: Logically Collective
911: Input Parameter:
912: . dm - the DMStag object
914: Output Parameters:
915: . arrX,arrY,arrZ - local 1D coordinate arrays
917: Level: intermediate
919: .seealso: DMSTAG, DMStagGetProductCoordinateArrays(), DMStagGetProductCoordinateArraysRead()
920: @*/
921: PetscErrorCode DMStagRestoreProductCoordinateArraysRead(DM dm,void *arrX,void *arrY,void *arrZ)
922: {
923: PetscErrorCode ierr;
926: DMStagRestoreProductCoordinateArrays_Private(dm,arrX,arrY,arrZ,PETSC_TRUE);
927: return(0);
928: }
930: /*@C
931: DMStagSetBoundaryTypes - set DMStag boundary types
933: Logically Collective; boundaryType0, boundaryType1, and boundaryType2 must contain common values
935: Input Parameters:
936: + dm - the DMStag object
937: - boundaryType0,boundaryType1,boundaryType2 - boundary types in each direction
939: Notes:
940: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
942: Level: advanced
944: .seealso: DMSTAG, DMBoundaryType, DMStagCreate1d(), DMStagCreate2d(), DMStagCreate3d(), DMDASetBoundaryType()
945: @*/
946: PetscErrorCode DMStagSetBoundaryTypes(DM dm,DMBoundaryType boundaryType0,DMBoundaryType boundaryType1,DMBoundaryType boundaryType2)
947: {
948: PetscErrorCode ierr;
949: DM_Stag * const stag = (DM_Stag*)dm->data;
950: PetscInt dim;
953: DMGetDimension(dm,&dim);
958: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
959: if (boundaryType0 ) stag->boundaryType[0] = boundaryType0;
960: if (boundaryType1 && dim > 1) stag->boundaryType[1] = boundaryType1;
961: if (boundaryType2 && dim > 2) stag->boundaryType[2] = boundaryType2;
962: return(0);
963: }
965: /*@C
966: DMStagSetCoordinateDMType - set DM type to store coordinates
968: Logically Collective; dmtype must contain common value
970: Input Parameters:
971: + dm - the DMStag object
972: - dmtype - DMtype for coordinates, either DMSTAG or DMPRODUCT
974: Level: advanced
976: .seealso: DMSTAG, DMPRODUCT, DMGetCoordinateDM(), DMStagSetUniformCoordinates(), DMStagSetUniformCoordinatesExplicit(), DMStagSetUniformCoordinatesProduct(), DMType
977: @*/
978: PetscErrorCode DMStagSetCoordinateDMType(DM dm,DMType dmtype)
979: {
980: PetscErrorCode ierr;
981: DM_Stag * const stag = (DM_Stag*)dm->data;
985: PetscFree(stag->coordinateDMType);
986: PetscStrallocpy(dmtype,(char**)&stag->coordinateDMType);
987: return(0);
988: }
990: /*@C
991: DMStagSetDOF - set dof/stratum
993: Logically Collective; dof0, dof1, dof2, and dof3 must contain common values
995: Input Parameters:
996: + dm - the DMStag object
997: - dof0,dof1,dof2,dof3 - dof per stratum
999: Notes:
1000: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1002: Level: advanced
1004: .seealso: DMSTAG, DMDASetDof()
1005: @*/
1006: PetscErrorCode DMStagSetDOF(DM dm,PetscInt dof0, PetscInt dof1,PetscInt dof2,PetscInt dof3)
1007: {
1008: PetscErrorCode ierr;
1009: DM_Stag * const stag = (DM_Stag*)dm->data;
1010: PetscInt dim;
1018: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1019: DMGetDimension(dm,&dim);
1020: if (dof0 < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"dof0 cannot be negative");
1021: if (dof1 < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"dof1 cannot be negative");
1022: if (dim > 1 && dof2 < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"dof2 cannot be negative");
1023: if (dim > 2 && dof3 < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"dof3 cannot be negative");
1024: stag->dof[0] = dof0;
1025: stag->dof[1] = dof1;
1026: if (dim > 1) stag->dof[2] = dof2;
1027: if (dim > 2) stag->dof[3] = dof3;
1028: return(0);
1029: }
1031: /*@C
1032: DMStagSetNumRanks - set ranks in each direction in the global rank grid
1034: Logically Collective; nRanks0, nRanks1, and nRanks2 must contain common values
1036: Input Parameters:
1037: + dm - the DMStag object
1038: - nRanks0,nRanks1,nRanks2 - number of ranks in each direction
1040: Notes:
1041: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1043: Level: developer
1045: .seealso: DMSTAG, DMDASetNumProcs()
1046: @*/
1047: PetscErrorCode DMStagSetNumRanks(DM dm,PetscInt nRanks0,PetscInt nRanks1,PetscInt nRanks2)
1048: {
1049: PetscErrorCode ierr;
1050: DM_Stag * const stag = (DM_Stag*)dm->data;
1051: PetscInt dim;
1058: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1059: DMGetDimension(dm,&dim);
1060: if (nRanks0 != PETSC_DECIDE && nRanks0 < 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"number of ranks in X direction cannot be less than 1");
1061: if (dim > 1 && nRanks1 != PETSC_DECIDE && nRanks1 < 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"number of ranks in Y direction cannot be less than 1");
1062: if (dim > 2 && nRanks2 != PETSC_DECIDE && nRanks2 < 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"number of ranks in Z direction cannot be less than 1");
1063: if (nRanks0) stag->nRanks[0] = nRanks0;
1064: if (dim > 1 && nRanks1) stag->nRanks[1] = nRanks1;
1065: if (dim > 2 && nRanks2) stag->nRanks[2] = nRanks2;
1066: return(0);
1067: }
1069: /*@C
1070: DMStagSetStencilType - set elementwise ghost/halo stencil type
1072: Logically Collective; stencilType must contain common value
1074: Input Parameters:
1075: + dm - the DMStag object
1076: - stencilType - the elementwise ghost stencil type: DMSTAG_STENCIL_BOX, DMSTAG_STENCIL_STAR, or DMSTAG_STENCIL_NONE
1078: Level: beginner
1080: .seealso: DMSTAG, DMStagGetStencilType(), DMStagStencilType, DMDASetStencilType()
1081: @*/
1082: PetscErrorCode DMStagSetStencilType(DM dm,DMStagStencilType stencilType)
1083: {
1084: DM_Stag * const stag = (DM_Stag*)dm->data;
1089: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1090: stag->stencilType = stencilType;
1091: return(0);
1092: }
1094: /*@C
1095: DMStagSetStencilWidth - set elementwise stencil width
1097: Logically Collective; stencilWidth must contain common value
1099: Input Parameters:
1100: + dm - the DMStag object
1101: - stencilWidth - stencil/halo/ghost width in elements
1103: Level: beginner
1105: .seealso: DMSTAG, DMDASetStencilWidth()
1106: @*/
1107: PetscErrorCode DMStagSetStencilWidth(DM dm,PetscInt stencilWidth)
1108: {
1109: DM_Stag * const stag = (DM_Stag*)dm->data;
1114: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1115: if (stencilWidth < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Stencil width must be non-negative");
1116: stag->stencilWidth = stencilWidth;
1117: return(0);
1118: }
1120: /*@C
1121: DMStagSetGlobalSizes - set global element counts in each direction
1123: Logically Collective; N0, N1, and N2 must contain common values
1125: Input Parameters:
1126: + dm - the DMStag object
1127: - N0,N1,N2 - global elementwise sizes
1129: Notes:
1130: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1132: Level: advanced
1134: .seealso: DMSTAG, DMStagGetGlobalSizes(), DMDASetSizes()
1135: @*/
1136: PetscErrorCode DMStagSetGlobalSizes(DM dm,PetscInt N0,PetscInt N1,PetscInt N2)
1137: {
1138: PetscErrorCode ierr;
1139: DM_Stag * const stag = (DM_Stag*)dm->data;
1140: PetscInt dim;
1144: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1145: DMGetDimension(dm,&dim);
1146: if (N0 < 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_SIZ,"Number of elements in X direction must be positive");
1147: if (dim > 1 && N1 < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_SIZ,"Number of elements in Y direction must be positive");
1148: if (dim > 2 && N2 < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_SIZ,"Number of elements in Z direction must be positive");
1149: if (N0) stag->N[0] = N0;
1150: if (N1) stag->N[1] = N1;
1151: if (N2) stag->N[2] = N2;
1152: return(0);
1153: }
1155: /*@C
1156: DMStagSetOwnershipRanges - set elements per rank in each direction
1158: Logically Collective; lx, ly, and lz must contain common values
1160: Input Parameters:
1161: + dm - the DMStag object
1162: - lx,ly,lz - element counts for each rank in each direction
1164: Notes:
1165: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
1167: Level: developer
1169: .seealso: DMSTAG, DMStagSetGlobalSizes(), DMStagGetOwnershipRanges(), DMDASetOwnershipRanges()
1170: @*/
1171: PetscErrorCode DMStagSetOwnershipRanges(DM dm,PetscInt const *lx,PetscInt const *ly,PetscInt const *lz)
1172: {
1173: PetscErrorCode ierr;
1174: DM_Stag * const stag = (DM_Stag*)dm->data;
1175: const PetscInt *lin[3];
1176: PetscInt d,dim;
1180: if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1181: lin[0] = lx; lin[1] = ly; lin[2] = lz;
1182: DMGetDimension(dm,&dim);
1183: for (d=0; d<dim; ++d) {
1184: if (lin[d]) {
1185: if (stag->nRanks[d] < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of ranks");
1186: if (!stag->l[d]) {
1187: PetscMalloc1(stag->nRanks[d], &stag->l[d]);
1188: }
1189: PetscArraycpy(stag->l[d], lin[d], stag->nRanks[d]);
1190: }
1191: }
1192: return(0);
1193: }
1195: /*@C
1196: DMStagSetUniformCoordinates - set DMStag coordinates to be a uniform grid
1198: Collective
1200: Input Parameters:
1201: + dm - the DMStag object
1202: - xmin,xmax,ymin,ymax,zmin,zmax - maximum and minimum global coordinate values
1204: Notes:
1205: DMStag supports 2 different types of coordinate DM: DMSTAG and DMPRODUCT.
1206: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1208: Level: advanced
1210: .seealso: DMSTAG, DMPRODUCT, DMStagSetUniformCoordinatesExplicit(), DMStagSetUniformCoordinatesProduct(), DMStagSetCoordinateDMType(), DMGetCoordinateDM(), DMGetCoordinates(), DMDASetUniformCoordinates()
1211: @*/
1212: PetscErrorCode DMStagSetUniformCoordinates(DM dm,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)
1213: {
1214: PetscErrorCode ierr;
1215: DM_Stag * const stag = (DM_Stag*)dm->data;
1216: PetscBool flg_stag,flg_product;
1220: if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called after DMSetUp()");
1221: if (!stag->coordinateDMType) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"You must first call DMStagSetCoordinateDMType()");
1222: PetscStrcmp(stag->coordinateDMType,DMSTAG,&flg_stag);
1223: PetscStrcmp(stag->coordinateDMType,DMPRODUCT,&flg_product);
1224: if (flg_stag) {
1225: DMStagSetUniformCoordinatesExplicit(dm,xmin,xmax,ymin,ymax,zmin,zmax);
1226: } else if (flg_product) {
1227: DMStagSetUniformCoordinatesProduct(dm,xmin,xmax,ymin,ymax,zmin,zmax);
1228: } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported DM Type %s",stag->coordinateDMType);
1229: return(0);
1230: }
1232: /*@C
1233: DMStagSetUniformCoordinatesExplicit - set DMStag coordinates to be a uniform grid, storing all values
1235: Collective
1237: Input Parameters:
1238: + dm - the DMStag object
1239: - xmin,xmax,ymin,ymax,zmin,zmax - maximum and minimum global coordinate values
1241: Notes:
1242: DMStag supports 2 different types of coordinate DM: either another DMStag, or a DMProduct.
1243: If the grid is orthogonal, using DMProduct should be more efficient.
1244: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1246: Level: beginner
1248: .seealso: DMSTAG, DMStagSetUniformCoordinates(), DMStagSetUniformCoordinatesProduct(), DMStagSetCoordinateDMType()
1249: @*/
1250: PetscErrorCode DMStagSetUniformCoordinatesExplicit(DM dm,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)
1251: {
1252: PetscErrorCode ierr;
1253: DM_Stag * const stag = (DM_Stag*)dm->data;
1254: PetscInt dim;
1255: PetscBool flg;
1259: if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called after DMSetUp()");
1260: PetscStrcmp(stag->coordinateDMType,DMSTAG,&flg);
1261: if (stag->coordinateDMType && !flg) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Refusing to change an already-set DM coordinate type");
1262: DMStagSetCoordinateDMType(dm,DMSTAG);
1263: DMGetDimension(dm,&dim);
1264: switch (dim) {
1265: case 1: DMStagSetUniformCoordinatesExplicit_1d(dm,xmin,xmax); break;
1266: case 2: DMStagSetUniformCoordinatesExplicit_2d(dm,xmin,xmax,ymin,ymax); break;
1267: case 3: DMStagSetUniformCoordinatesExplicit_3d(dm,xmin,xmax,ymin,ymax,zmin,zmax); break;
1268: default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
1269: }
1270: return(0);
1271: }
1273: /*@C
1274: DMStagSetUniformCoordinatesProduct - create uniform coordinates, as a product of 1D arrays
1276: 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
1277: coordinates.
1279: Collective
1281: Input Parameters:
1282: + dm - the DMStag object
1283: - xmin,xmax,ymin,ymax,zmin,zmax - maximum and minimum global coordinate values
1285: Notes:
1286: Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.
1288: Level: intermediate
1290: .seealso: DMSTAG, DMPRODUCT, DMStagSetUniformCoordinates(), DMStagSetUniformCoordinatesExplicit(), DMStagSetCoordinateDMType()
1291: @*/
1292: PetscErrorCode DMStagSetUniformCoordinatesProduct(DM dm,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)
1293: {
1294: PetscErrorCode ierr;
1295: DM_Stag * const stag = (DM_Stag*)dm->data;
1296: DM dmc;
1297: PetscInt dim,d,dof0,dof1;
1298: PetscBool flg;
1302: if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called after DMSetUp()");
1303: PetscStrcmp(stag->coordinateDMType,DMPRODUCT,&flg);
1304: if (stag->coordinateDMType && !flg) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Refusing to change an already-set DM coordinate type");
1305: DMStagSetCoordinateDMType(dm,DMPRODUCT);
1306: DMGetCoordinateDM(dm,&dmc);
1307: DMGetDimension(dm,&dim);
1309: /* Create 1D sub-DMs, living on subcommunicators */
1311: dof0 = 0;
1312: for (d=0; d<dim; ++d) { /* Include point dof in the sub-DMs if any non-element strata are active in the original DMStag */
1313: if (stag->dof[d]) {
1314: dof0 = 1;
1315: break;
1316: }
1317: }
1318: dof1 = stag->dof[dim] ? 1 : 0; /* Include element dof in the sub-DMs if the elements are active in the original DMStag */
1320: for (d=0; d<dim; ++d) {
1321: DM subdm;
1322: MPI_Comm subcomm;
1323: PetscMPIInt color;
1324: const PetscMPIInt key = 0; /* let existing rank break ties */
1326: /* Choose colors based on position in the plane orthogonal to this dim, and split */
1327: switch (d) {
1328: case 0: color = (dim > 1 ? stag->rank[1] : 0) + (dim > 2 ? stag->nRanks[1]*stag->rank[2] : 0); break;
1329: case 1: color = stag->rank[0] + (dim > 2 ? stag->nRanks[0]*stag->rank[2] : 0); break;
1330: case 2: color = stag->rank[0] + stag->nRanks[0]*stag->rank[1] ; break;
1331: default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP,"Unsupported dimension index %D",d);
1332: }
1333: MPI_Comm_split(PetscObjectComm((PetscObject)dm),color,key,&subcomm);
1335: /* Create sub-DMs living on these new communicators (which are destroyed by DMProduct) */
1336: DMStagCreate1d(subcomm,stag->boundaryType[d],stag->N[d],dof0,dof1,stag->stencilType,stag->stencilWidth,stag->l[d],&subdm);
1337: DMSetUp(subdm);
1338: switch (d) {
1339: case 0:
1340: DMStagSetUniformCoordinatesExplicit(subdm,xmin,xmax,0,0,0,0);
1341: break;
1342: case 1:
1343: DMStagSetUniformCoordinatesExplicit(subdm,ymin,ymax,0,0,0,0);
1344: break;
1345: case 2:
1346: DMStagSetUniformCoordinatesExplicit(subdm,zmin,zmax,0,0,0,0);
1347: break;
1348: default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP,"Unsupported dimension index %D",d);
1349: }
1350: DMProductSetDM(dmc,d,subdm);
1351: DMProductSetDimensionIndex(dmc,d,0);
1352: DMDestroy(&subdm);
1353: MPI_Comm_free(&subcomm);
1354: }
1355: return(0);
1356: }
1358: /*@C
1359: DMStagVecGetArray - get access to local array
1361: Logically Collective
1363: This function returns a (dim+1)-dimensional array for a dim-dimensional
1364: DMStag.
1366: The first 1-3 dimensions indicate an element in the global
1367: numbering, using the standard C ordering.
1369: The final dimension in this array corresponds to a degree
1370: of freedom with respect to this element, for example corresponding to
1371: the element or one of its neighboring faces, edges, or vertices.
1373: For example, for a 3D DMStag, indexing is array[k][j][i][idx], where k is the
1374: index in the z-direction, j is the index in the y-direction, and i is the
1375: index in the x-direction.
1377: "idx" is obtained with DMStagGetLocationSlot(), since the correct offset
1378: into the (dim+1)-dimensional C array depends on the grid size and the number
1379: of dof stored at each location.
1381: Input Parameters:
1382: + dm - the DMStag object
1383: - vec - the Vec object
1385: Output Parameters:
1386: . array - the array
1388: Notes:
1389: DMStagVecRestoreArray() must be called, once finished with the array
1391: Level: beginner
1393: .seealso: DMSTAG, DMStagVecGetArrayRead(), DMStagGetLocationSlot(), DMGetLocalVector(), DMCreateLocalVector(), DMGetGlobalVector(), DMCreateGlobalVector(), DMDAVecGetArray(), DMDAVecGetArrayDOF()
1394: @*/
1395: PetscErrorCode DMStagVecGetArray(DM dm,Vec vec,void *array)
1396: {
1397: PetscErrorCode ierr;
1398: DM_Stag * const stag = (DM_Stag*)dm->data;
1399: PetscInt dim;
1400: PetscInt nLocal;
1405: DMGetDimension(dm,&dim);
1406: VecGetLocalSize(vec,&nLocal);
1407: if (nLocal != stag->entriesGhost) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMStag local size %D\n",nLocal,stag->entriesGhost);
1408: switch (dim) {
1409: case 1:
1410: VecGetArray2d(vec,stag->nGhost[0],stag->entriesPerElement,stag->startGhost[0],0,(PetscScalar***)array);
1411: break;
1412: case 2:
1413: VecGetArray3d(vec,stag->nGhost[1],stag->nGhost[0],stag->entriesPerElement,stag->startGhost[1],stag->startGhost[0],0,(PetscScalar****)array);
1414: break;
1415: case 3:
1416: 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);
1417: break;
1418: default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
1419: }
1420: return(0);
1421: }
1423: /*@C
1424: DMStagVecGetArrayRead - get read-only access to a local array
1426: Logically Collective
1428: See the man page for DMStagVecGetArray() for more information.
1430: Input Parameters:
1431: + dm - the DMStag object
1432: - vec - the Vec object
1434: Output Parameters:
1435: . array - the read-only array
1437: Notes:
1438: DMStagVecRestoreArrayRead() must be called, once finished with the array
1440: Level: beginner
1442: .seealso: DMSTAG, DMStagVecGetArrayRead(), DMStagGetLocationSlot(), DMGetLocalVector(), DMCreateLocalVector(), DMGetGlobalVector(), DMCreateGlobalVector(), DMDAVecGetArrayRead(), DMDAVecGetArrayDOFRead()
1443: @*/
1444: PetscErrorCode DMStagVecGetArrayRead(DM dm,Vec vec,void *array)
1445: {
1446: PetscErrorCode ierr;
1447: DM_Stag * const stag = (DM_Stag*)dm->data;
1448: PetscInt dim;
1449: PetscInt nLocal;
1454: DMGetDimension(dm,&dim);
1455: VecGetLocalSize(vec,&nLocal);
1456: if (nLocal != stag->entriesGhost) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMStag local size %D\n",nLocal,stag->entriesGhost);
1457: switch (dim) {
1458: case 1:
1459: VecGetArray2dRead(vec,stag->nGhost[0],stag->entriesPerElement,stag->startGhost[0],0,(PetscScalar***)array);
1460: break;
1461: case 2:
1462: VecGetArray3dRead(vec,stag->nGhost[1],stag->nGhost[0],stag->entriesPerElement,stag->startGhost[1],stag->startGhost[0],0,(PetscScalar****)array);
1463: break;
1464: case 3:
1465: 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);
1466: break;
1467: default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
1468: }
1469: return(0);
1470: }
1472: /*@C
1473: DMStagVecRestoreArray - restore access to a raw array
1475: Logically Collective
1477: Input Parameters:
1478: + dm - the DMStag object
1479: - vec - the Vec object
1481: Output Parameters:
1482: . array - the array
1484: Level: beginner
1486: .seealso: DMSTAG, DMStagVecGetArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
1487: @*/
1488: PetscErrorCode DMStagVecRestoreArray(DM dm,Vec vec,void *array)
1489: {
1490: PetscErrorCode ierr;
1491: DM_Stag * const stag = (DM_Stag*)dm->data;
1492: PetscInt dim;
1493: PetscInt nLocal;
1498: DMGetDimension(dm,&dim);
1499: VecGetLocalSize(vec,&nLocal);
1500: if (nLocal != stag->entriesGhost) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMStag local size %D\n",nLocal,stag->entriesGhost);
1501: switch (dim) {
1502: case 1:
1503: VecRestoreArray2d(vec,stag->nGhost[0],stag->entriesPerElement,stag->startGhost[0],0,(PetscScalar***)array);
1504: break;
1505: case 2:
1506: VecRestoreArray3d(vec,stag->nGhost[1],stag->nGhost[0],stag->entriesPerElement,stag->startGhost[1],stag->startGhost[0],0,(PetscScalar****)array);
1507: break;
1508: case 3:
1509: 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);
1510: break;
1511: default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
1512: }
1513: return(0);
1514: }
1516: /*@C
1517: DMStagVecRestoreArrayRead - restore read-only access to a raw array
1519: Logically Collective
1521: Input Parameters:
1522: + dm - the DMStag object
1523: - vec - the Vec object
1525: Output Parameters:
1526: . array - the read-only array
1528: Level: beginner
1530: .seealso: DMSTAG, DMStagVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMDAVecRestoreArrayDOFRead()
1531: @*/
1532: PetscErrorCode DMStagVecRestoreArrayRead(DM dm,Vec vec,void *array)
1533: {
1534: PetscErrorCode ierr;
1535: DM_Stag * const stag = (DM_Stag*)dm->data;
1536: PetscInt dim;
1537: PetscInt nLocal;
1542: DMGetDimension(dm,&dim);
1543: VecGetLocalSize(vec,&nLocal);
1544: if (nLocal != stag->entriesGhost) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMStag local size %D\n",nLocal,stag->entriesGhost);
1545: switch (dim) {
1546: case 1:
1547: VecRestoreArray2dRead(vec,stag->nGhost[0],stag->entriesPerElement,stag->startGhost[0],0,(PetscScalar***)array);
1548: break;
1549: case 2:
1550: VecRestoreArray3dRead(vec,stag->nGhost[1],stag->nGhost[0],stag->entriesPerElement,stag->startGhost[1],stag->startGhost[0],0,(PetscScalar****)array);
1551: break;
1552: case 3:
1553: 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);
1554: break;
1555: default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
1556: }
1557: return(0);
1558: }