Actual source code: stagutils.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  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];
 40:   PetscBool      checkDof;

 44:   DMGetDimension(dm,&dim);
 45:   if (dim > DMSTAG_MAX_DIM) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not implemented for %D dimensions",dim);
 46:   arr[0] = arrX; arr[1] = arrY; arr[2] = arrZ;
 47:   DMGetCoordinateDM(dm,&dmCoord);
 48:   if (!dmCoord) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM does not have a coordinate DM");
 49:   {
 50:     PetscBool isProduct;
 51:     DMType    dmType;
 52:     DMGetType(dmCoord,&dmType);
 53:     PetscStrcmp(DMPRODUCT,dmType,&isProduct);
 54:     if (!isProduct) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate DM is not of type DMPRODUCT");
 55:   }
 56:   for (s=0; s<DMSTAG_MAX_STRATA; ++s) dofCheck[s] = 0;
 57:   checkDof = PETSC_FALSE;
 58:   for (d=0; d<dim; ++d) {
 59:     DM        subDM;
 60:     DMType    dmType;
 61:     PetscBool isStag;
 62:     PetscInt  dof[DMSTAG_MAX_STRATA],subDim;
 63:     Vec       coord1d_local;

 65:     /* Ignore unrequested arrays */
 66:     if (!arr[d]) continue;

 68:     DMProductGetDM(dmCoord,d,&subDM);
 69:     if (!subDM) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate DM is missing sub DM %D",d);
 70:     DMGetDimension(subDM,&subDim);
 71:     if (subDim != 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate sub-DM is not of dimension 1");
 72:     DMGetType(subDM,&dmType);
 73:     PetscStrcmp(DMSTAG,dmType,&isStag);
 74:     if (!isStag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate sub-DM is not of type DMSTAG");
 75:     DMStagGetDOF(subDM,&dof[0],&dof[1],&dof[2],&dof[3]);
 76:     if (!checkDof) {
 77:       for (s=0; s<DMSTAG_MAX_STRATA; ++s) dofCheck[s] = dof[s];
 78:       checkDof = PETSC_TRUE;
 79:     } else {
 80:       for (s=0; s<DMSTAG_MAX_STRATA; ++s) {
 81:         if (dofCheck[s] != dof[s]) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate sub-DMs have different dofs");
 82:       }
 83:     }
 84:     DMGetCoordinatesLocal(subDM,&coord1d_local);
 85:     if (read) {
 86:       DMStagVecGetArrayRead(subDM,coord1d_local,arr[d]);
 87:     } else {
 88:       DMStagVecGetArray(subDM,coord1d_local,arr[d]);
 89:     }
 90:   }
 91:   return(0);
 92: }

 94: /*@C
 95:   DMStagGetProductCoordinateArrays - extract local product coordinate arrays, one per dimension

 97:   Logically Collective

 99:   A high-level helper function to quickly extract local coordinate arrays.

101:   Note that 2-dimensional arrays are returned. See
102:   DMStagVecGetArray(), which is called internally to produce these arrays
103:   representing coordinates on elements and vertices (element boundaries)
104:   for a 1-dimensional DMStag in each coordinate direction.

106:   One should use DMStagGetProductCoordinateSlot() to determine appropriate
107:   indices for the second dimension in these returned arrays. This function
108:   checks that the coordinate array is a suitable product of 1-dimensional
109:   DMStag objects.

111:   Input Parameter:
112: . dm - the DMStag object

114:   Output Parameters:
115: . arrX,arrY,arrZ - local 1D coordinate arrays

117:   Level: intermediate

119: .seealso: DMSTAG, DMPRODUCT, DMStagGetProductCoordinateArraysRead(), DMStagSetUniformCoordinates(), DMStagSetUniformCoordinatesProduct(), DMStagGetProductCoordinateLocationSlot()
120: @*/
121: PetscErrorCode DMStagGetProductCoordinateArrays(DM dm,void* arrX,void* arrY,void* arrZ)
122: {

126:   DMStagGetProductCoordinateArrays_Private(dm,arrX,arrY,arrZ,PETSC_FALSE);
127:   return(0);
128: }

130: /*@C
131:   DMStagGetProductCoordinateArraysRead - extract product coordinate arrays, read-only

133:   Logically Collective

135:   See the man page for DMStagGetProductCoordinateArrays() for more information.

137:   Input Parameter:
138: . dm - the DMStag object

140:   Output Parameters:
141: . arrX,arrY,arrZ - local 1D coordinate arrays

143:   Level: intermediate

145: .seealso: DMSTAG, DMPRODUCT, DMStagGetProductCoordinateArrays(), DMStagSetUniformCoordinates(), DMStagSetUniformCoordinatesProduct(), DMStagGetProductCoordinateLocationSlot()
146: @*/
147: PetscErrorCode DMStagGetProductCoordinateArraysRead(DM dm,void* arrX,void* arrY,void* arrZ)
148: {

152:   DMStagGetProductCoordinateArrays_Private(dm,arrX,arrY,arrZ,PETSC_TRUE);
153:   return(0);
154: }

156: /*@C
157:   DMStagGetProductCoordinateLocationSlot - get slot for use with local product coordinate arrays

159:   Not Collective

161:   High-level helper function to get slot indices for 1D coordinate DMs,
162:   for use with DMStagGetProductCoordinateArrays() and related functions.

164:   Input Parameters:
165: + dm - the DMStag object
166: - loc - the grid location

168:   Output Parameter:
169: . slot - the index to use in local arrays

171:   Notes:
172:   Checks that the coordinates are actually set up so that using the
173:   slots from the first 1d coordinate sub-DM is valid for all the 1D coordinate sub-DMs.

175:   Level: intermediate

177: .seealso: DMSTAG, DMPRODUCT, DMStagGetProductCoordinateArrays(), DMStagGetProductCoordinateArraysRead(), DMStagSetUniformCoordinates()
178: @*/
179: PETSC_EXTERN PetscErrorCode DMStagGetProductCoordinateLocationSlot(DM dm,DMStagStencilLocation loc,PetscInt *slot)
180: {
182:   DM             dmCoord;
183:   PetscInt       dim,dofCheck[DMSTAG_MAX_STRATA],s,d;

187:   DMGetDimension(dm,&dim);
188:   DMGetCoordinateDM(dm,&dmCoord);
189:   if (!dmCoord) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM does not have a coordinate DM");
190:   {
191:     PetscBool isProduct;
192:     DMType    dmType;
193:     DMGetType(dmCoord,&dmType);
194:     PetscStrcmp(DMPRODUCT,dmType,&isProduct);
195:     if (!isProduct) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate DM is not of type DMPRODUCT");
196:   }
197:   for (s=0; s<DMSTAG_MAX_STRATA; ++s) dofCheck[s] = 0;
198:   for (d=0; d<dim; ++d) {
199:     DM        subDM;
200:     DMType    dmType;
201:     PetscBool isStag;
202:     PetscInt  dof[DMSTAG_MAX_STRATA],subDim;
203:     DMProductGetDM(dmCoord,d,&subDM);
204:     if (!subDM) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate DM is missing sub DM %D",d);
205:     DMGetDimension(subDM,&subDim);
206:     if (subDim != 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate sub-DM is not of dimension 1");
207:     DMGetType(subDM,&dmType);
208:     PetscStrcmp(DMSTAG,dmType,&isStag);
209:     if (!isStag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate sub-DM is not of type DMSTAG");
210:     DMStagGetDOF(subDM,&dof[0],&dof[1],&dof[2],&dof[3]);
211:     if (d == 0) {
212:       const PetscInt component = 0;
213:       for (s=0; s<DMSTAG_MAX_STRATA; ++s) dofCheck[s] = dof[s];
214:       DMStagGetLocationSlot(subDM,loc,component,slot);
215:     } else {
216:       for (s=0; s<DMSTAG_MAX_STRATA; ++s) {
217:         if (dofCheck[s] != dof[s]) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Coordinate sub-DMs have different dofs");
218:       }
219:     }
220:   }
221:   return(0);
222: }

224: /*@C
225:   DMStagGetCorners - return global element indices of the local region (excluding ghost points)

227:   Not Collective

229:   Input Parameter:
230: . dm - the DMStag object

232:   Output Parameters:
233: + x,y,z - starting element indices in each direction
234: . m,n,p - element widths in each direction
235: - nExtrax,nExtray,nExtraz - number of extra partial elements in each direction.

237:   Notes:
238:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.

240:   The number of extra partial elements is either 1 or 0.
241:   The value is 1 on right, top, and front non-periodic domain ("physical") boundaries,
242:   in the x, y, and z directions respectively, and otherwise 0.

244:   Level: beginner

246: .seealso: DMSTAG, DMStagGetGhostCorners(), DMDAGetCorners()
247: @*/
248: PetscErrorCode DMStagGetCorners(DM dm,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p,PetscInt *nExtrax,PetscInt *nExtray,PetscInt *nExtraz)
249: {
250:   const DM_Stag * const stag = (DM_Stag*)dm->data;

254:   if (x) *x = stag->start[0];
255:   if (y) *y = stag->start[1];
256:   if (z) *z = stag->start[2];
257:   if (m) *m = stag->n[0];
258:   if (n) *n = stag->n[1];
259:   if (p) *p = stag->n[2];
260:   if (nExtrax) *nExtrax = stag->boundaryType[0] != DM_BOUNDARY_PERIODIC && stag->lastRank[0] ? 1 : 0;
261:   if (nExtray) *nExtray = stag->boundaryType[1] != DM_BOUNDARY_PERIODIC && stag->lastRank[1] ? 1 : 0;
262:   if (nExtraz) *nExtraz = stag->boundaryType[2] != DM_BOUNDARY_PERIODIC && stag->lastRank[2] ? 1 : 0;
263:   return(0);
264: }

266: /*@C
267:   DMStagGetDOF - get number of DOF associated with each stratum of the grid

269:   Not Collective

271:   Input Parameter:
272: . dm - the DMStag object

274:   Output Parameters:
275: + dof0 - the number of points per 0-cell (vertex/node)
276: . dof1 - the number of points per 1-cell (element in 1D, edge in 2D and 3D)
277: . dof2 - the number of points per 2-cell (element in 2D, face in 3D)
278: - dof3 - the number of points per 3-cell (element in 3D)

280:   Level: beginner

282: .seealso: DMSTAG, DMStagGetCorners(), DMStagGetGhostCorners(), DMStagGetGlobalSizes(), DMStagGetStencilWidth(), DMStagGetBoundaryTypes(), DMStagGetLocationDof(), DMDAGetDof()
283: @*/
284: PetscErrorCode DMStagGetDOF(DM dm,PetscInt *dof0,PetscInt *dof1,PetscInt *dof2,PetscInt *dof3)
285: {
286:   const DM_Stag * const stag = (DM_Stag*)dm->data;

290:   if (dof0) *dof0 = stag->dof[0];
291:   if (dof1) *dof1 = stag->dof[1];
292:   if (dof2) *dof2 = stag->dof[2];
293:   if (dof3) *dof3 = stag->dof[3];
294:   return(0);
295: }

297: /*@C
298:   DMStagGetGhostCorners - return global element indices of the local region, including ghost points

300:   Not Collective

302:   Input Argument:
303: . dm - the DMStag object

305:   Output Arguments:
306: + x,y,z - starting element indices in each direction
307: - m,n,p - element widths in each direction

309:   Notes:
310:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.

312:   Level: beginner

314: .seealso: DMSTAG, DMStagGetCorners(), DMDAGetGhostCorners()
315: @*/
316: PetscErrorCode DMStagGetGhostCorners(DM dm,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p)
317: {
318:   const DM_Stag * const stag = (DM_Stag*)dm->data;

322:   if (x) *x = stag->startGhost[0];
323:   if (y) *y = stag->startGhost[1];
324:   if (z) *z = stag->startGhost[2];
325:   if (m) *m = stag->nGhost[0];
326:   if (n) *n = stag->nGhost[1];
327:   if (p) *p = stag->nGhost[2];
328:   return(0);
329: }

331: /*@C
332:   DMStagGetGlobalSizes - get global element counts

334:   Not Collective

336:   Input Parameter:
337: . dm - the DMStag object

339:   Output Parameters:
340: . M,N,P - global element counts in each direction

342:   Notes:
343:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.

345:   Level: beginner

347: .seealso: DMSTAG, DMStagGetLocalSizes(), DMDAGetInfo()
348: @*/
349: PetscErrorCode DMStagGetGlobalSizes(DM dm,PetscInt* M,PetscInt* N,PetscInt* P)
350: {
351:   const DM_Stag * const stag = (DM_Stag*)dm->data;

355:   if (M) *M = stag->N[0];
356:   if (N) *N = stag->N[1];
357:   if (P) *P = stag->N[2];
358:   return(0);
359: }

361: /*@C
362:   DMStagGetIsFirstRank - get boolean value for whether this rank is first in each direction in the rank grid

364:   Not Collective

366:   Input Parameter:
367: . dm - the DMStag object

369:   Output Parameters:
370: . isFirstRank0,isFirstRank1,isFirstRank2 - whether this rank is first in each direction

372:   Level: intermediate

374:   Notes:
375:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.

377: .seealso: DMSTAG, DMStagGetIsLastRank()
378: @*/
379: PetscErrorCode DMStagGetIsFirstRank(DM dm,PetscBool *isFirstRank0,PetscBool *isFirstRank1,PetscBool *isFirstRank2)
380: {
381:   const DM_Stag * const stag = (DM_Stag*)dm->data;

385:   if (isFirstRank0) *isFirstRank0 = stag->firstRank[0];
386:   if (isFirstRank1) *isFirstRank1 = stag->firstRank[1];
387:   if (isFirstRank2) *isFirstRank2 = stag->firstRank[2];
388:   return(0);
389: }

391: /*@C
392:   DMStagGetIsLastRank - get boolean value for whether this rank is last in each direction in the rank grid

394:   Not Collective

396:   Input Parameter:
397: . dm - the DMStag object

399:   Output Parameters:
400: . isLastRank0,isLastRank1,isLastRank2 - whether this rank is last in each direction

402:   Level: intermediate

404:   Notes:
405:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
406:   Level: intermediate

408: .seealso: DMSTAG, DMStagGetIsFirstRank()
409: @*/
410: PetscErrorCode DMStagGetIsLastRank(DM dm,PetscBool *isLastRank0,PetscBool *isLastRank1,PetscBool *isLastRank2)
411: {
412:   const DM_Stag * const stag = (DM_Stag*)dm->data;

416:   if (isLastRank0) *isLastRank0 = stag->lastRank[0];
417:   if (isLastRank1) *isLastRank1 = stag->lastRank[1];
418:   if (isLastRank2) *isLastRank2 = stag->lastRank[2];
419:   return(0);
420: }

422: /*@C
423:   DMStagGetLocalSizes - get local elementwise sizes

425:   Not Collective

427:   Input Parameter:
428: . dm - the DMStag object

430:   Output Parameters:
431: . m,n,p - local element counts (excluding ghosts) in each direction

433:   Notes:
434:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
435:   Level: intermediate

437:   Level: beginner

439: .seealso: DMSTAG, DMStagGetGlobalSizes(), DMStagGetDOF(), DMStagGetNumRanks(), DMDAGetLocalInfo()
440: @*/
441: PetscErrorCode DMStagGetLocalSizes(DM dm,PetscInt* m,PetscInt* n,PetscInt* p)
442: {
443:   const DM_Stag * const stag = (DM_Stag*)dm->data;

447:   if (m) *m = stag->n[0];
448:   if (n) *n = stag->n[1];
449:   if (p) *p = stag->n[2];
450:   return(0);
451: }

453: /*@C
454:   DMStagGetNumRanks - get number of ranks in each direction in the global grid decomposition

456:   Not Collective

458:   Input Parameter:
459: . dm - the DMStag object

461:   Output Parameters:
462: . nRanks0,nRanks1,nRanks2 - number of ranks in each direction in the grid decomposition

464:   Notes:
465:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.
466:   Level: intermediate

468:   Level: beginner

470: .seealso: DMSTAG, DMStagGetGlobalSizes(), DMStagGetLocalSize(), DMStagSetNumRanks(), DMDAGetInfo()
471: @*/
472: PetscErrorCode DMStagGetNumRanks(DM dm,PetscInt *nRanks0,PetscInt *nRanks1,PetscInt *nRanks2)
473: {
474:   const DM_Stag * const stag = (DM_Stag*)dm->data;

478:   if (nRanks0) *nRanks0 = stag->nRanks[0];
479:   if (nRanks1) *nRanks1 = stag->nRanks[1];
480:   if (nRanks2) *nRanks2 = stag->nRanks[2];
481:   return(0);
482: }

484: /*@C
485:   DMStagGetEntriesPerElement - get number of entries per element in the local representation

487:   Not Collective

489:   Input Parameter:
490: . dm - the DMStag object

492:   Output Parameters:
493: . entriesPerElement - number of entries associated with each element in the local representation

495:   Notes:
496:   This is the natural block size for most local operations. In 1D it is equal to dof0 + dof1,
497:   in 2D it is equal to dof0 + 2*dof1 + dof2, and in 3D it is equal to dof0 + 3*dof1 + 3*dof2 + dof3

499:   Level: developer

501: .seealso: DMSTAG, DMStagGetDOF()
502: @*/
503: PetscErrorCode DMStagGetEntriesPerElement(DM dm,PetscInt *entriesPerElement)
504: {
505:   const DM_Stag * const stag = (DM_Stag*)dm->data;

509:   if (entriesPerElement) *entriesPerElement = stag->entriesPerElement;
510:   return(0);
511: }

513: /*@C
514:   DMStagGetStencilType - get elementwise ghost/halo stencil type

516:   Not Collective

518:   Input Parameter:
519: . dm - the DMStag object

521:   Output Parameter:
522: . stencilType - the elementwise ghost stencil type: DMSTAG_STENCIL_BOX, DMSTAG_STENCIL_STAR, or DMSTAG_STENCIL_NONE

524:   Level: beginner

526: .seealso: DMSTAG, DMStagSetStencilType(), DMStagStencilType, DMDAGetInfo()
527: @*/
528: PetscErrorCode DMStagGetStencilType(DM dm,DMStagStencilType *stencilType)
529: {
530:   DM_Stag * const stag = (DM_Stag*)dm->data;

534:   *stencilType = stag->stencilType;
535:   return(0);
536: }

538: /*@C
539:   DMStagGetStencilWidth - get elementwise stencil width

541:   Not Collective

543:   Input Parameter:
544: . dm - the DMStag object

546:   Output Parameters:
547: . stencilWidth - stencil/halo/ghost width in elements

549:   Level: beginner

551: .seealso: DMSTAG, DMDAGetStencilWidth(), DMDAGetInfo()
552: @*/
553: PetscErrorCode DMStagGetStencilWidth(DM dm,PetscInt *stencilWidth)
554: {
555:   const DM_Stag * const stag = (DM_Stag*)dm->data;

559:   if (stencilWidth) *stencilWidth = stag->stencilWidth;
560:   return(0);
561: }

563: /*@C
564:   DMStagGetOwnershipRanges - get elements per rank in each direction

566:   Not Collective

568:   Input Parameter:
569: .     dm - the DMStag object

571:   Output Parameters:
572: +     lx - ownership along x direction (optional)
573: .     ly - ownership along y direction (optional)
574: -     lz - ownership along z direction (optional)

576:   Notes:
577:   These correspond to the optional final arguments passed to DMStagCreate1d(), DMStagCreate2d(), and DMStagCreate3d().

579:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.

581:   In C you should not free these arrays, nor change the values in them.
582:   They will only have valid values while the DMStag they came from still exists (has not been destroyed).

584:   Level: intermediate

586: .seealso: DMSTAG, DMStagSetGlobalSizes(), DMStagSetOwnershipRanges(), DMStagCreate1d(), DMStagCreate2d(), DMStagCreate3d(), DMDAGetOwnershipRanges()
587: @*/
588: PetscErrorCode DMStagGetOwnershipRanges(DM dm,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
589: {
590:   const DM_Stag * const stag = (DM_Stag*)dm->data;

594:   if (lx) *lx = stag->l[0];
595:   if (ly) *ly = stag->l[1];
596:   if (lz) *lz = stag->l[2];
597:   return(0);
598: }

600: /*@C
601:   DMStagCreateCompatibleDMStag - create a compatible DMStag with different dof/stratum

603:   Collective

605:   Input Parameters:
606: + dm - the DMStag object
607: - dof0,dof1,dof2,dof3 - number of dof on each stratum in the new DMStag

609:   Output Parameters:
610: . newdm - the new, compatible DMStag

612:   Notes:
613:   Dof supplied for strata too big for the dimension are ignored; these may be set to 0.
614:   For example, for a 2-dimensional DMStag, dof2 sets the number of dof per element,
615:   and dof3 is unused. For a 3-dimensional DMStag, dof3 sets the number of dof per element.
616:   
617:   In contrast to DMDACreateCompatibleDMDA(), coordinates are not reused.

619:   Level: intermediate

621: .seealso: DMSTAG, DMDACreateCompatibleDMDA(), DMGetCompatibility(), DMStagMigrateVec()
622: @*/
623: PetscErrorCode DMStagCreateCompatibleDMStag(DM dm,PetscInt dof0,PetscInt dof1,PetscInt dof2,PetscInt dof3,DM *newdm)
624: {
625:   PetscErrorCode  ierr;

629:   DMStagDuplicateWithoutSetup(dm,PetscObjectComm((PetscObject)dm),newdm);
630:   DMStagSetDOF(*newdm,dof0,dof1,dof2,dof3);
631:   DMSetUp(*newdm);
632:   return(0);
633: }

635: /*@C
636:   DMStagGetLocationSlot - get index to use in accessing raw local arrays

638:   Not Collective

640:   Input Parameters:
641: + dm - the DMStag object
642: . loc - location relative to an element
643: - c - component

645:   Output Parameter:
646: . slot - index to use

648:   Notes:
649:   Provides an appropriate index to use with DMStagVecGetArray() and friends.
650:   This is required so that the user doesn't need to know about the ordering of
651:   dof associated with each local element.

653:   Level: beginner

655: .seealso: DMSTAG, DMStagVecGetArray(), DMStagVecGetArrayRead(), DMStagGetDOF(), DMStagGetEntriesPerElement()
656: @*/
657: PetscErrorCode DMStagGetLocationSlot(DM dm,DMStagStencilLocation loc,PetscInt c,PetscInt *slot)
658: {
659:   DM_Stag * const stag = (DM_Stag*)dm->data;

663: #if defined(PETSC_USE_DEBUG)
664:   {
666:     PetscInt       dof;
667:     DMStagGetLocationDOF(dm,loc,&dof);
668:     if (dof < 1) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Location %s has no dof attached",DMStagStencilLocations[loc]);
669:     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);
670:   }
671: #endif
672:   *slot = stag->locationOffsets[loc] + c;
673:   return(0);
674: }

676: /*@C
677:   DMStagMigrateVec - transfer a vector associated with a DMStag to a vector associated with a compatible DMStag

679:   Collective

681:   Input Parameters:
682: + dm - the source DMStag object
683: . vec - the source vector, compatible with dm
684: . dmTo - the compatible destination DMStag object
685: - vecTo - the destination vector, compatible with dmTo

687:   Notes:
688:   Extra dof are ignored, and unfilled dof are zeroed.
689:   Currently only implemented to migrate global vectors to global vectors.

691:   Level: advanced

693: .seealso: DMSTAG, DMStagCreateCompatibleDMStag(), DMGetCompatibility(), DMStagVecSplitToDMDA()
694: @*/
695: PetscErrorCode DMStagMigrateVec(DM dm,Vec vec,DM dmTo,Vec vecTo)
696: {
697:   PetscErrorCode    ierr;
698:   DM_Stag * const   stag = (DM_Stag*)dm->data;
699:   DM_Stag * const   stagTo = (DM_Stag*)dmTo->data;
700:   PetscInt          nLocalTo,nLocal,dim,i,j,k;
701:   PetscInt          start[DMSTAG_MAX_DIM],startGhost[DMSTAG_MAX_DIM],n[DMSTAG_MAX_DIM],nExtra[DMSTAG_MAX_DIM],offset[DMSTAG_MAX_DIM];
702:   Vec               vecToLocal,vecLocal;
703:   PetscBool         compatible,compatibleSet;
704:   const PetscScalar *arr;
705:   PetscScalar       *arrTo;
706:   const PetscInt    epe   = stag->entriesPerElement;
707:   const PetscInt    epeTo = stagTo->entriesPerElement;

714:   DMGetCompatibility(dm,dmTo,&compatible,&compatibleSet);
715:   if (!compatibleSet || !compatible) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"DMStag objects must be shown to be compatible");
716:   DMGetDimension(dm,&dim);
717:   VecGetLocalSize(vec,&nLocal);
718:   VecGetLocalSize(vecTo,&nLocalTo);
719:   if (nLocal != stag->entries|| nLocalTo !=stagTo->entries) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Vector migration only implemented for global vector to global vector.");
720:   DMStagGetCorners(dm,&start[0],&start[1],&start[2],&n[0],&n[1],&n[2],&nExtra[0],&nExtra[1],&nExtra[2]);
721:   DMStagGetGhostCorners(dm,&startGhost[0],&startGhost[1],&startGhost[2],NULL,NULL,NULL);
722:   for (i=0; i<DMSTAG_MAX_DIM; ++i) offset[i] = start[i]-startGhost[i];

724:   /* Proceed by transferring to a local vector, copying, and transferring back to a global vector */
725:   DMGetLocalVector(dm,&vecLocal);
726:   DMGetLocalVector(dmTo,&vecToLocal);
727:   DMGlobalToLocalBegin(dm,vec,INSERT_VALUES,vecLocal);
728:   DMGlobalToLocalEnd(dm,vec,INSERT_VALUES,vecLocal);
729:   VecGetArrayRead(vecLocal,&arr);
730:   VecGetArray(vecToLocal,&arrTo);
731:   /* Note that some superfluous copying of entries on partial dummy elements is done */
732:   if (dim == 1) {
733:     for (i=offset[0]; i<offset[0] + n[0] + nExtra[0]; ++i) {
734:       PetscInt d = 0,dTo = 0,b = 0,bTo = 0;
735:       PetscInt si;
736:       for (si=0; si<2; ++si) {
737:         b   += stag->dof[si];
738:         bTo += stagTo->dof[si];
739:         for (; d < b && dTo < bTo; ++d,++dTo) arrTo[i*epeTo + dTo] = arr[i*epe + d];
740:         for (; dTo < bTo         ;     ++dTo) arrTo[i*epeTo + dTo] = 0.0;
741:         d = b;
742:       }
743:     }
744:   } else if (dim == 2) {
745:     const PetscInt epr   = stag->nGhost[0] * epe;
746:     const PetscInt eprTo = stagTo->nGhost[0] * epeTo;
747:     for (j=offset[1]; j<offset[1] + n[1] + nExtra[1]; ++j) {
748:       for (i=offset[0]; i<offset[0] + n[0] + nExtra[0]; ++i) {
749:         const PetscInt base   = j*epr   + i*epe;
750:         const PetscInt baseTo = j*eprTo + i*epeTo;
751:         PetscInt d = 0,dTo = 0,b = 0,bTo = 0;
752:         const PetscInt s[4] = {0,1,1,2}; /* Dimensions of points, in order */
753:         PetscInt si;
754:         for (si=0; si<4; ++si) {
755:             b   += stag->dof[s[si]];
756:             bTo += stagTo->dof[s[si]];
757:             for (; d < b && dTo < bTo; ++d,++dTo) arrTo[baseTo + dTo] = arr[base + d];
758:             for (;          dTo < bTo;     ++dTo) arrTo[baseTo + dTo] = 0.0;
759:             d = b;
760:         }
761:       }
762:     }
763:   } else if (dim == 3) {
764:     const PetscInt epr   = stag->nGhost[0]   * epe;
765:     const PetscInt eprTo = stagTo->nGhost[0] * epeTo;
766:     const PetscInt epl   = stag->nGhost[1]   * epr;
767:     const PetscInt eplTo = stagTo->nGhost[1] * eprTo;
768:     for (k=offset[2]; k<offset[2] + n[2] + nExtra[2]; ++k) {
769:       for (j=offset[1]; j<offset[1] + n[1] + nExtra[1]; ++j) {
770:         for (i=offset[0]; i<offset[0] + n[0] + nExtra[0]; ++i) {
771:           PetscInt d = 0,dTo = 0,b = 0,bTo = 0;
772:           const PetscInt base   = k*epl   + j*epr   + i*epe;
773:           const PetscInt baseTo = k*eplTo + j*eprTo + i*epeTo;
774:           const PetscInt s[8] = {0,1,1,2,1,2,2,3}; /* dimensions of points, in order */
775:           PetscInt is;
776:           for (is=0; is<8; ++is) {
777:             b   += stag->dof[s[is]];
778:             bTo += stagTo->dof[s[is]];
779:             for (; d < b && dTo < bTo; ++d,++dTo) arrTo[baseTo + dTo] = arr[base + d];
780:             for (;          dTo < bTo;     ++dTo) arrTo[baseTo + dTo] = 0.0;
781:             d = b;
782:           }
783:         }
784:       }
785:     }
786:   } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
787:   VecRestoreArrayRead(vecLocal,&arr);
788:   VecRestoreArray(vecToLocal,&arrTo);
789:   DMRestoreLocalVector(dm,&vecLocal);
790:   DMLocalToGlobalBegin(dmTo,vecToLocal,INSERT_VALUES,vecTo);
791:   DMLocalToGlobalEnd(dmTo,vecToLocal,INSERT_VALUES,vecTo);
792:   DMRestoreLocalVector(dmTo,&vecToLocal);
793:   return(0);
794: }

796: /*@C
797:   DMStagPopulateLocalToGlobalInjective - populate an internal 1-to-1 local-to-global map

799:   Collective

801:   Creates an internal object which explicitly maps a single local degree of
802:   freedom to each global degree of freedom. This is used, if populated,
803:   instead of SCATTER_REVERSE_LOCAL with the (1-to-many, in general)
804:   global-to-local map, when DMLocalToGlobal() is called with INSERT_VALUES.
805:   This allows usage, for example, even in the periodic, 1-rank case, where
806:   the inverse of the global-to-local map, even when restricted to on-rank
807:   communication, is non-injective. This is at the cost of storing an additional
808:   VecScatter object inside each DMStag object.

810:   Input Parameter:
811: . dm - the DMStag object

813:   Notes:
814:   In normal usage, library users shouldn't be concerned with this function,
815:   as it is called during DMSetUp(), when required.

817:   Returns immediately if the internal map is already populated.

819:   Developer Notes:
820:   This could, if desired, be moved up to a general DM routine. It would allow,
821:   for example, DMDA to support DMLocalToGlobal() with INSERT_VALUES,
822:   even in the single-rank periodic case.

824:   Level: developer

826: .seealso: DMSTAG, DMLocalToGlobal(), VecScatter
827: @*/
828: PetscErrorCode DMStagPopulateLocalToGlobalInjective(DM dm)
829: {
830:   PetscErrorCode  ierr;
831:   PetscInt        dim;
832:   DM_Stag * const stag  = (DM_Stag*)dm->data;

836:   if (stag->ltog_injective) return(0); /* Don't re-populate */
837:   DMGetDimension(dm,&dim);
838:   switch (dim) {
839:     case 1: DMStagPopulateLocalToGlobalInjective_1d(dm); break;
840:     case 2: DMStagPopulateLocalToGlobalInjective_2d(dm); break;
841:     case 3: DMStagPopulateLocalToGlobalInjective_3d(dm); break;
842:     default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
843:   }
844:   return(0);
845: }

847: static PetscErrorCode DMStagRestoreProductCoordinateArrays_Private(DM dm,void *arrX,void *arrY,void *arrZ,PetscBool read)
848: {
849:   PetscErrorCode  ierr;
850:   PetscInt        dim,d;
851:   void*           arr[DMSTAG_MAX_DIM];
852:   DM              dmCoord;

856:   DMGetDimension(dm,&dim);
857:   if (dim > DMSTAG_MAX_DIM) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not implemented for %D dimensions",dim);
858:   arr[0] = arrX; arr[1] = arrY; arr[2] = arrZ;
859:   DMGetCoordinateDM(dm,&dmCoord);
860:   for (d=0; d<dim; ++d) {
861:     DM  subDM;
862:     Vec coord1d_local;

864:     /* Ignore unrequested arrays */
865:     if (!arr[d]) continue;

867:     DMProductGetDM(dmCoord,d,&subDM);
868:     DMGetCoordinatesLocal(subDM,&coord1d_local);
869:     if (read) {
870:       DMStagVecRestoreArrayRead(subDM,coord1d_local,arr[d]);
871:     } else {
872:       DMStagVecRestoreArray(subDM,coord1d_local,arr[d]);
873:     }
874:   }
875:   return(0);
876: }

878: /*@C
879:   DMStagRestoreProductCoordinateArrays - restore local array access

881:   Logically Collective

883:   Input Parameter:
884: . dm - the DMStag object

886:   Output Parameters:
887: . arrX,arrY,arrZ - local 1D coordinate arrays

889:   Level: intermediate

891:   Notes:
892:   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:

894: $   DMGetCoordinateDM(dm,&cdm);
895: $   for (d=0; d<3; ++d) {
896: $     DM  subdm;
897: $     Vec coor,coor_local;

899: $     DMProductGetDM(cdm,d,&subdm);
900: $     DMGetCoordinates(subdm,&coor);
901: $     DMGetCoordinatesLocal(subdm,&coor_local);
902: $     DMLocalToGlobal(subdm,coor_local,INSERT_VALUES,coor);
903: $     PetscPrintf(PETSC_COMM_WORLD,"Coordinates dim %D:\n",d);
904: $     VecView(coor,PETSC_VIEWER_STDOUT_WORLD);
905: $   }

907: .seealso: DMSTAG, DMStagGetProductCoordinateArrays(), DMStagGetProductCoordinateArraysRead()
908: @*/
909: PetscErrorCode DMStagRestoreProductCoordinateArrays(DM dm,void *arrX,void *arrY,void *arrZ)
910: {
911:   PetscErrorCode  ierr;

914:   DMStagRestoreProductCoordinateArrays_Private(dm,arrX,arrY,arrZ,PETSC_FALSE);
915:   return(0);
916: }

918: /*@C
919:   DMStagRestoreProductCoordinateArraysRead - restore local product array access, read-only

921:   Logically Collective

923:   Input Parameter:
924: . dm - the DMStag object

926:   Output Parameters:
927: . arrX,arrY,arrZ - local 1D coordinate arrays

929:   Level: intermediate

931: .seealso: DMSTAG, DMStagGetProductCoordinateArrays(), DMStagGetProductCoordinateArraysRead()
932: @*/
933: PetscErrorCode DMStagRestoreProductCoordinateArraysRead(DM dm,void *arrX,void *arrY,void *arrZ)
934: {
935:   PetscErrorCode  ierr;

938:   DMStagRestoreProductCoordinateArrays_Private(dm,arrX,arrY,arrZ,PETSC_TRUE);
939:   return(0);
940: }

942: /*@C
943:   DMStagSetBoundaryTypes - set DMStag boundary types

945:   Logically Collective; boundaryType0, boundaryType1, and boundaryType2 must contain common values

947:   Input Parameters:
948: + dm - the DMStag object
949: - boundaryType0,boundaryType1,boundaryType2 - boundary types in each direction

951:   Notes:
952:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.

954:   Level: advanced

956: .seealso: DMSTAG, DMBoundaryType, DMStagCreate1d(), DMStagCreate2d(), DMStagCreate3d(), DMDASetBoundaryType()
957: @*/
958: PetscErrorCode DMStagSetBoundaryTypes(DM dm,DMBoundaryType boundaryType0,DMBoundaryType boundaryType1,DMBoundaryType boundaryType2)
959: {
960:   PetscErrorCode  ierr;
961:   DM_Stag * const stag  = (DM_Stag*)dm->data;
962:   PetscInt        dim;

965:   DMGetDimension(dm,&dim);
970:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
971:                stag->boundaryType[0] = boundaryType0;
972:   if (dim > 1) stag->boundaryType[1] = boundaryType1;
973:   if (dim > 2) stag->boundaryType[2] = boundaryType2;
974:   return(0);
975: }

977: /*@C
978:   DMStagSetCoordinateDMType - set DM type to store coordinates

980:   Logically Collective; dmtype must contain common value

982:   Input Parameters:
983: + dm - the DMStag object
984: - dmtype - DMtype for coordinates, either DMSTAG or DMPRODUCT

986:   Level: advanced

988: .seealso: DMSTAG, DMPRODUCT, DMGetCoordinateDM(), DMStagSetUniformCoordinates(), DMStagSetUniformCoordinatesExplicit(), DMStagSetUniformCoordinatesProduct(), DMType
989: @*/
990: PetscErrorCode DMStagSetCoordinateDMType(DM dm,DMType dmtype)
991: {
992:   PetscErrorCode  ierr;
993:   DM_Stag * const stag = (DM_Stag*)dm->data;

997:   PetscFree(stag->coordinateDMType);
998:   PetscStrallocpy(dmtype,(char**)&stag->coordinateDMType);
999:   return(0);
1000: }

1002: /*@C
1003:   DMStagSetDOF - set dof/stratum

1005:   Logically Collective; dof0, dof1, dof2, and dof3 must contain common values

1007:   Input Parameters:
1008: + dm - the DMStag object
1009: - dof0,dof1,dof2,dof3 - dof per stratum

1011:   Notes:
1012:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.

1014:   Level: advanced

1016: .seealso: DMSTAG, DMDASetDof()
1017: @*/
1018: PetscErrorCode DMStagSetDOF(DM dm,PetscInt dof0, PetscInt dof1,PetscInt dof2,PetscInt dof3)
1019: {
1020:   PetscErrorCode  ierr;
1021:   DM_Stag * const stag = (DM_Stag*)dm->data;
1022:   PetscInt        dim;

1030:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1031:   DMGetDimension(dm,&dim);
1032:   if (dof0 < 0)            SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"dof0 cannot be negative");
1033:   if (dof1 < 0)            SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"dof1 cannot be negative");
1034:   if (dim > 1 && dof2 < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"dof2 cannot be negative");
1035:   if (dim > 2 && dof3 < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"dof3 cannot be negative");
1036:                stag->dof[0] = dof0;
1037:                stag->dof[1] = dof1;
1038:   if (dim > 1) stag->dof[2] = dof2;
1039:   if (dim > 2) stag->dof[3] = dof3;
1040:   return(0);
1041: }

1043: /*@C
1044:   DMStagSetNumRanks - set ranks in each direction in the global rank grid

1046:   Logically Collective; nRanks0, nRanks1, and nRanks2 must contain common values

1048:   Input Parameters:
1049: + dm - the DMStag object
1050: - nRanks0,nRanks1,nRanks2 - number of ranks in each direction

1052:   Notes:
1053:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.

1055:   Level: developer

1057: .seealso: DMSTAG, DMDASetNumProcs()
1058: @*/
1059: PetscErrorCode DMStagSetNumRanks(DM dm,PetscInt nRanks0,PetscInt nRanks1,PetscInt nRanks2)
1060: {
1061:   PetscErrorCode  ierr;
1062:   DM_Stag * const stag = (DM_Stag*)dm->data;
1063:   PetscInt        dim;

1070:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1071:   DMGetDimension(dm,&dim);
1072:   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");
1073:   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");
1074:   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");
1075:   if (nRanks0) stag->nRanks[0] = nRanks0;
1076:   if (dim > 1 && nRanks1) stag->nRanks[1] = nRanks1;
1077:   if (dim > 2 && nRanks2) stag->nRanks[2] = nRanks2;
1078:   return(0);
1079: }

1081: /*@C
1082:   DMStagSetStencilType - set elementwise ghost/halo stencil type

1084:   Logically Collective; stencilType must contain common value

1086:   Input Parameters:
1087: + dm - the DMStag object
1088: - stencilType - the elementwise ghost stencil type: DMSTAG_STENCIL_BOX, DMSTAG_STENCIL_STAR, or DMSTAG_STENCIL_NONE

1090:   Level: beginner

1092: .seealso: DMSTAG, DMStagGetStencilType(), DMStagStencilType, DMDASetStencilType()
1093: @*/
1094: PetscErrorCode DMStagSetStencilType(DM dm,DMStagStencilType stencilType)
1095: {
1096:   DM_Stag * const stag = (DM_Stag*)dm->data;

1101:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1102:   stag->stencilType = stencilType;
1103:   return(0);
1104: }

1106: /*@C
1107:   DMStagSetStencilWidth - set elementwise stencil width

1109:   Logically Collective; stencilWidth must contain common value

1111:   Input Parameters:
1112: + dm - the DMStag object
1113: - stencilWidth - stencil/halo/ghost width in elements

1115:   Level: beginner

1117: .seealso: DMSTAG, DMDASetStencilWidth()
1118: @*/
1119: PetscErrorCode DMStagSetStencilWidth(DM dm,PetscInt stencilWidth)
1120: {
1121:   DM_Stag * const stag = (DM_Stag*)dm->data;

1126:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1127:   if (stencilWidth < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Stencil width must be non-negative");
1128:   stag->stencilWidth = stencilWidth;
1129:   return(0);
1130: }

1132: /*@C
1133:   DMStagSetGlobalSizes - set global element counts in each direction

1135:   Logically Collective; N0, N1, and N2 must contain common values

1137:   Input Parameters:
1138: + dm - the DMStag object
1139: - N0,N1,N2 - global elementwise sizes

1141:   Notes:
1142:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.

1144:   Level: advanced

1146: .seealso: DMSTAG, DMStagGetGlobalSizes(), DMDASetSizes()
1147: @*/
1148: PetscErrorCode DMStagSetGlobalSizes(DM dm,PetscInt N0,PetscInt N1,PetscInt N2)
1149: {
1150:   PetscErrorCode  ierr;
1151:   DM_Stag * const stag = (DM_Stag*)dm->data;
1152:   PetscInt        dim;

1156:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1157:   DMGetDimension(dm,&dim);
1158:   if (N0 < 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_SIZ,"Number of elements in X direction must be positive");
1159:   if (dim > 1 && N1 < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_SIZ,"Number of elements in Y direction must be positive");
1160:   if (dim > 2 && N2 < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_SIZ,"Number of elements in Z direction must be positive");
1161:   if (N0) stag->N[0] = N0;
1162:   if (N1) stag->N[1] = N1;
1163:   if (N2) stag->N[2] = N2;
1164:   return(0);
1165: }

1167: /*@C
1168:   DMStagSetOwnershipRanges - set elements per rank in each direction

1170:   Logically Collective; lx, ly, and lz must contain common values

1172:   Input Parameters:
1173: + dm - the DMStag object
1174: - lx,ly,lz - element counts for each rank in each direction

1176:   Notes:
1177:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids. These arguments may be set to NULL in this case.

1179:   Level: developer

1181: .seealso: DMSTAG, DMStagSetGlobalSizes(), DMStagGetOwnershipRanges(), DMDASetOwnershipRanges()
1182: @*/
1183: PetscErrorCode DMStagSetOwnershipRanges(DM dm,PetscInt const *lx,PetscInt const *ly,PetscInt const *lz)
1184: {
1185:   PetscErrorCode  ierr;
1186:   DM_Stag * const stag = (DM_Stag*)dm->data;
1187:   const PetscInt  *lin[3];
1188:   PetscInt        d,dim;

1192:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1193:   lin[0] = lx; lin[1] = ly; lin[2] = lz;
1194:   DMGetDimension(dm,&dim);
1195:   for (d=0; d<dim; ++d) {
1196:     if (lin[d]) {
1197:       if (stag->nRanks[d] < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of ranks");
1198:       if (!stag->l[d]) {
1199:         PetscMalloc1(stag->nRanks[d], &stag->l[d]);
1200:       }
1201:       PetscArraycpy(stag->l[d], lin[d], stag->nRanks[d]);
1202:     }
1203:   }
1204:   return(0);
1205: }

1207: /*@C
1208:   DMStagSetUniformCoordinates - set DMStag coordinates to be a uniform grid

1210:   Collective

1212:   Input Parameters:
1213: + dm - the DMStag object
1214: - xmin,xmax,ymin,ymax,zmin,zmax - maximum and minimum global coordinate values

1216:   Notes:
1217:   DMStag supports 2 different types of coordinate DM: DMSTAG and DMPRODUCT.
1218:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.

1220:   Level: advanced

1222: .seealso: DMSTAG, DMPRODUCT, DMStagSetUniformCoordinatesExplicit(), DMStagSetUniformCoordinatesProduct(), DMStagSetCoordinateDMType(), DMGetCoordinateDM(), DMGetCoordinates(), DMDASetUniformCoordinates()
1223: @*/
1224: PetscErrorCode DMStagSetUniformCoordinates(DM dm,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)
1225: {
1226:   PetscErrorCode  ierr;
1227:   DM_Stag * const stag = (DM_Stag*)dm->data;
1228:   PetscBool       flg_stag,flg_product;

1232:   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called after DMSetUp()");
1233:   if (!stag->coordinateDMType) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"You must first call DMStagSetCoordinateDMType()");
1234:   PetscStrcmp(stag->coordinateDMType,DMSTAG,&flg_stag);
1235:   PetscStrcmp(stag->coordinateDMType,DMPRODUCT,&flg_product);
1236:   if (flg_stag) {
1237:     DMStagSetUniformCoordinatesExplicit(dm,xmin,xmax,ymin,ymax,zmin,zmax);
1238:   } else if (flg_product) {
1239:     DMStagSetUniformCoordinatesProduct(dm,xmin,xmax,ymin,ymax,zmin,zmax);
1240:   } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported DM Type %s",stag->coordinateDMType);
1241:   return(0);
1242: }

1244: /*@C
1245:   DMStagSetUniformCoordinatesExplicit - set DMStag coordinates to be a uniform grid, storing all values

1247:   Collective

1249:   Input Parameters:
1250: + dm - the DMStag object
1251: - xmin,xmax,ymin,ymax,zmin,zmax - maximum and minimum global coordinate values

1253:   Notes:
1254:   DMStag supports 2 different types of coordinate DM: either another DMStag, or a DMProduct.
1255:   If the grid is orthogonal, using DMProduct should be more efficient.
1256:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.

1258:   Level: beginner

1260: .seealso: DMSTAG, DMStagSetUniformCoordinates(), DMStagSetUniformCoordinatesProduct(), DMStagSetCoordinateDMType()
1261: @*/
1262: PetscErrorCode DMStagSetUniformCoordinatesExplicit(DM dm,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)
1263: {
1264:   PetscErrorCode  ierr;
1265:   DM_Stag * const stag = (DM_Stag*)dm->data;
1266:   PetscInt        dim;
1267:   PetscBool       flg;

1271:   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called after DMSetUp()");
1272:   PetscStrcmp(stag->coordinateDMType,DMSTAG,&flg);
1273:   if (stag->coordinateDMType && !flg) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Refusing to change an already-set DM coordinate type");
1274:   DMStagSetCoordinateDMType(dm,DMSTAG);
1275:   DMGetDimension(dm,&dim);
1276:   switch (dim) {
1277:     case 1: DMStagSetUniformCoordinatesExplicit_1d(dm,xmin,xmax);                     break;
1278:     case 2: DMStagSetUniformCoordinatesExplicit_2d(dm,xmin,xmax,ymin,ymax);           break;
1279:     case 3: DMStagSetUniformCoordinatesExplicit_3d(dm,xmin,xmax,ymin,ymax,zmin,zmax); break;
1280:     default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
1281:   }
1282:   return(0);
1283: }

1285: /*@C
1286:   DMStagSetUniformCoordinatesProduct - create uniform coordinates, as a product of 1D arrays

1288:   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
1289:   coordinates.

1291:   Collective

1293:   Input Parameters:
1294: + dm - the DMStag object
1295: - xmin,xmax,ymin,ymax,zmin,zmax - maximum and minimum global coordinate values

1297:   Notes:
1298:   Arguments corresponding to higher dimensions are ignored for 1D and 2D grids.

1300:   Level: intermediate

1302: .seealso: DMSTAG, DMPRODUCT, DMStagSetUniformCoordinates(), DMStagSetUniformCoordinatesExplicit(), DMStagSetCoordinateDMType()
1303: @*/
1304: PetscErrorCode DMStagSetUniformCoordinatesProduct(DM dm,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)
1305: {
1306:   PetscErrorCode  ierr;
1307:   DM_Stag * const stag = (DM_Stag*)dm->data;
1308:   DM              dmc;
1309:   PetscInt        dim,d,dof0,dof1;
1310:   PetscBool       flg;

1314:   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"This function must be called after DMSetUp()");
1315:   PetscStrcmp(stag->coordinateDMType,DMPRODUCT,&flg);
1316:   if (stag->coordinateDMType && !flg) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Refusing to change an already-set DM coordinate type");
1317:   DMStagSetCoordinateDMType(dm,DMPRODUCT);
1318:   DMGetCoordinateDM(dm,&dmc);
1319:   DMGetDimension(dm,&dim);

1321:   /* Create 1D sub-DMs, living on subcommunicators */

1323:   dof0 = 0;
1324:   for (d=0; d<dim; ++d) {          /* Include point dof in the sub-DMs if any non-element strata are active in the original DMStag */
1325:     if (stag->dof[d]) {
1326:       dof0 = 1;
1327:       break;
1328:     }
1329:   }
1330:   dof1 = stag->dof[dim] ? 1 : 0; /*  Include element dof in the sub-DMs if the elements are active in the original DMStag */

1332:   for (d=0; d<dim; ++d) {
1333:     DM                subdm;
1334:     MPI_Comm          subcomm;
1335:     PetscMPIInt       color;
1336:     const PetscMPIInt key = 0; /* let existing rank break ties */

1338:     /* Choose colors based on position in the plane orthogonal to this dim, and split */
1339:     switch (d) {
1340:       case 0: color = (dim > 1 ? stag->rank[1] : 0)  + (dim > 2 ? stag->nRanks[1]*stag->rank[2] : 0); break;
1341:       case 1: color =            stag->rank[0]       + (dim > 2 ? stag->nRanks[0]*stag->rank[2] : 0); break;
1342:       case 2: color =            stag->rank[0]       +            stag->nRanks[0]*stag->rank[1]     ; break;
1343:       default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP,"Unsupported dimension index %D",d);
1344:     }
1345:     MPI_Comm_split(PetscObjectComm((PetscObject)dm),color,key,&subcomm);

1347:     /* Create sub-DMs living on these new communicators (which are destroyed by DMProduct) */
1348:     DMStagCreate1d(subcomm,stag->boundaryType[d],stag->N[d],dof0,dof1,stag->stencilType,stag->stencilWidth,stag->l[d],&subdm);
1349:     DMSetUp(subdm);
1350:     switch (d) {
1351:       case 0:
1352:         DMStagSetUniformCoordinatesExplicit(subdm,xmin,xmax,0,0,0,0);
1353:         break;
1354:       case 1:
1355:         DMStagSetUniformCoordinatesExplicit(subdm,ymin,ymax,0,0,0,0);
1356:         break;
1357:       case 2:
1358:         DMStagSetUniformCoordinatesExplicit(subdm,zmin,zmax,0,0,0,0);
1359:         break;
1360:       default: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP,"Unsupported dimension index %D",d);
1361:     }
1362:     DMProductSetDM(dmc,d,subdm);
1363:     DMProductSetDimensionIndex(dmc,d,0);
1364:     DMDestroy(&subdm);
1365:     MPI_Comm_free(&subcomm);
1366:   }
1367:   return(0);
1368: }

1370: /*@C
1371:   DMStagVecGetArray - get access to local array

1373:   Logically Collective

1375:   This function returns a (dim+1)-dimensional array for a dim-dimensional
1376:   DMStag.

1378:   The first 1-3 dimensions indicate an element in the global
1379:   numbering, using the standard C ordering.

1381:   The final dimension in this array corresponds to a degree
1382:   of freedom with respect to this element, for example corresponding to
1383:   the element or one of its neighboring faces, edges, or vertices.

1385:   For example, for a 3D DMStag, indexing is array[k][j][i][idx], where k is the
1386:   index in the z-direction, j is the index in the y-direction, and i is the
1387:   index in the x-direction.

1389:   "idx" is obtained with DMStagGetLocationSlot(), since the correct offset
1390:   into the (dim+1)-dimensional C array depends on the grid size and the number
1391:   of dof stored at each location.

1393:   Input Parameters:
1394: + dm - the DMStag object
1395: - vec - the Vec object

1397:   Output Parameters:
1398: . array - the array

1400:   Notes:
1401:   DMStagVecRestoreArray() must be called, once finished with the array

1403:   Level: beginner

1405: .seealso: DMSTAG, DMStagVecGetArrayRead(), DMStagGetLocationSlot(), DMGetLocalVector(), DMCreateLocalVector(), DMGetGlobalVector(), DMCreateGlobalVector(), DMDAVecGetArray(), DMDAVecGetArrayDOF()
1406: @*/
1407: PetscErrorCode DMStagVecGetArray(DM dm,Vec vec,void *array)
1408: {
1409:   PetscErrorCode  ierr;
1410:   DM_Stag * const stag = (DM_Stag*)dm->data;
1411:   PetscInt        dim;
1412:   PetscInt        nLocal;

1417:   DMGetDimension(dm,&dim);
1418:   VecGetLocalSize(vec,&nLocal);
1419:   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);
1420:   switch (dim) {
1421:     case 1:
1422:       VecGetArray2d(vec,stag->nGhost[0],stag->entriesPerElement,stag->startGhost[0],0,(PetscScalar***)array);
1423:       break;
1424:     case 2:
1425:       VecGetArray3d(vec,stag->nGhost[1],stag->nGhost[0],stag->entriesPerElement,stag->startGhost[1],stag->startGhost[0],0,(PetscScalar****)array);
1426:       break;
1427:     case 3:
1428:       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);
1429:       break;
1430:     default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
1431:   }
1432:   return(0);
1433: }

1435: /*@C
1436:   DMStagVecGetArrayRead - get read-only access to a local array

1438:   Logically Collective

1440:   See the man page for DMStagVecGetArray() for more information.

1442:   Input Parameters:
1443: + dm - the DMStag object
1444: - vec - the Vec object

1446:   Output Parameters:
1447: . array - the read-only array

1449:   Notes:
1450:   DMStagVecRestoreArrayRead() must be called, once finished with the array

1452:   Level: beginner

1454: .seealso: DMSTAG, DMStagVecGetArrayRead(), DMStagGetLocationSlot(), DMGetLocalVector(), DMCreateLocalVector(), DMGetGlobalVector(), DMCreateGlobalVector(), DMDAVecGetArrayRead(), DMDAVecGetArrayDOFRead()
1455: @*/
1456: PetscErrorCode DMStagVecGetArrayRead(DM dm,Vec vec,void *array)
1457: {
1458:   PetscErrorCode  ierr;
1459:   DM_Stag * const stag = (DM_Stag*)dm->data;
1460:   PetscInt        dim;
1461:   PetscInt        nLocal;

1466:   DMGetDimension(dm,&dim);
1467:   VecGetLocalSize(vec,&nLocal);
1468:   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);
1469:   switch (dim) {
1470:     case 1:
1471:       VecGetArray2dRead(vec,stag->nGhost[0],stag->entriesPerElement,stag->startGhost[0],0,(PetscScalar***)array);
1472:       break;
1473:     case 2:
1474:       VecGetArray3dRead(vec,stag->nGhost[1],stag->nGhost[0],stag->entriesPerElement,stag->startGhost[1],stag->startGhost[0],0,(PetscScalar****)array);
1475:       break;
1476:     case 3:
1477:       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);
1478:       break;
1479:     default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
1480:   }
1481:   return(0);
1482: }

1484: /*@C
1485:   DMStagVecRestoreArray - restore access to a raw array

1487:   Logically Collective

1489:   Input Parameters:
1490: + dm - the DMStag object
1491: - vec - the Vec object

1493:   Output Parameters:
1494: . array - the array

1496:   Level: beginner

1498: .seealso: DMSTAG, DMStagVecGetArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
1499: @*/
1500: PetscErrorCode DMStagVecRestoreArray(DM dm,Vec vec,void *array)
1501: {
1502:   PetscErrorCode  ierr;
1503:   DM_Stag * const stag = (DM_Stag*)dm->data;
1504:   PetscInt        dim;
1505:   PetscInt        nLocal;

1510:   DMGetDimension(dm,&dim);
1511:   VecGetLocalSize(vec,&nLocal);
1512:   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);
1513:   switch (dim) {
1514:     case 1:
1515:       VecRestoreArray2d(vec,stag->nGhost[0],stag->entriesPerElement,stag->startGhost[0],0,(PetscScalar***)array);
1516:       break;
1517:     case 2:
1518:       VecRestoreArray3d(vec,stag->nGhost[1],stag->nGhost[0],stag->entriesPerElement,stag->startGhost[1],stag->startGhost[0],0,(PetscScalar****)array);
1519:       break;
1520:     case 3:
1521:       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);
1522:       break;
1523:     default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
1524:   }
1525:   return(0);
1526: }

1528: /*@C
1529:   DMStagVecRestoreArrayRead - restore read-only access to a raw array

1531:   Logically Collective

1533:   Input Parameters:
1534: + dm - the DMStag object
1535: - vec - the Vec object

1537:   Output Parameters:
1538: . array - the read-only array

1540:   Level: beginner

1542: .seealso: DMSTAG, DMStagVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMDAVecRestoreArrayDOFRead()
1543: @*/
1544: PetscErrorCode DMStagVecRestoreArrayRead(DM dm,Vec vec,void *array)
1545: {
1546:   PetscErrorCode  ierr;
1547:   DM_Stag * const stag = (DM_Stag*)dm->data;
1548:   PetscInt        dim;
1549:   PetscInt        nLocal;

1554:   DMGetDimension(dm,&dim);
1555:   VecGetLocalSize(vec,&nLocal);
1556:   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);
1557:   switch (dim) {
1558:     case 1:
1559:       VecRestoreArray2dRead(vec,stag->nGhost[0],stag->entriesPerElement,stag->startGhost[0],0,(PetscScalar***)array);
1560:       break;
1561:     case 2:
1562:       VecRestoreArray3dRead(vec,stag->nGhost[1],stag->nGhost[0],stag->entriesPerElement,stag->startGhost[1],stag->startGhost[0],0,(PetscScalar****)array);
1563:       break;
1564:     case 3:
1565:       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);
1566:       break;
1567:     default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unsupported dimension %D",dim);
1568:   }
1569:   return(0);
1570: }