Actual source code: dacreate.c

  1: #include <petsc/private/dmdaimpl.h>

  3: static PetscErrorCode DMSetFromOptions_DA(DM da, PetscOptionItems *PetscOptionsObject)
  4: {
  5:   DM_DA         *dd     = (DM_DA *)da->data;
  6:   PetscInt       refine = 0, dim = da->dim, maxnlevels = 100, refx[100], refy[100], refz[100], n, i;
  7:   DMBoundaryType bt = DM_BOUNDARY_NONE;
  8:   PetscBool      flg;

 10:   PetscFunctionBegin;
 11:   PetscCheck(dd->M >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Dimension must be non-negative, call DMSetFromOptions() if you want to change the value at runtime");
 12:   PetscCheck(dd->N >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Dimension must be non-negative, call DMSetFromOptions() if you want to change the value at runtime");
 13:   PetscCheck(dd->P >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Dimension must be non-negative, call DMSetFromOptions() if you want to change the value at runtime");

 15:   PetscOptionsHeadBegin(PetscOptionsObject, "DMDA Options");
 16:   PetscCall(PetscOptionsBoundedInt("-da_grid_x", "Number of grid points in x direction", "DMDASetSizes", dd->M, &dd->M, NULL, 1));
 17:   if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_grid_y", "Number of grid points in y direction", "DMDASetSizes", dd->N, &dd->N, NULL, 1));
 18:   if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_grid_z", "Number of grid points in z direction", "DMDASetSizes", dd->P, &dd->P, NULL, 1));

 20:   PetscCall(PetscOptionsBoundedInt("-da_overlap", "Decomposition overlap in all directions", "DMDASetOverlap", dd->xol, &dd->xol, &flg, 0));
 21:   if (flg) PetscCall(DMDASetOverlap(da, dd->xol, dd->xol, dd->xol));
 22:   PetscCall(PetscOptionsBoundedInt("-da_overlap_x", "Decomposition overlap in x direction", "DMDASetOverlap", dd->xol, &dd->xol, NULL, 0));
 23:   if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_overlap_y", "Decomposition overlap in y direction", "DMDASetOverlap", dd->yol, &dd->yol, NULL, 0));
 24:   if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_overlap_z", "Decomposition overlap in z direction", "DMDASetOverlap", dd->zol, &dd->zol, NULL, 0));

 26:   PetscCall(PetscOptionsBoundedInt("-da_local_subdomains", "", "DMDASetNumLocalSubdomains", dd->Nsub, &dd->Nsub, &flg, PETSC_DECIDE));
 27:   if (flg) PetscCall(DMDASetNumLocalSubDomains(da, dd->Nsub));

 29:   /* Handle DMDA parallel distribution */
 30:   PetscCall(PetscOptionsBoundedInt("-da_processors_x", "Number of processors in x direction", "DMDASetNumProcs", dd->m, &dd->m, NULL, PETSC_DECIDE));
 31:   if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_processors_y", "Number of processors in y direction", "DMDASetNumProcs", dd->n, &dd->n, NULL, PETSC_DECIDE));
 32:   if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_processors_z", "Number of processors in z direction", "DMDASetNumProcs", dd->p, &dd->p, NULL, PETSC_DECIDE));
 33:   // Handle boundaries
 34:   PetscCall(PetscOptionsEnum("-da_bd_x", "Boundary type for x direction", "DMDASetBoundaryType", DMBoundaryTypes, (PetscEnum)dd->bx, (PetscEnum *)&dd->bx, NULL));
 35:   if (dim > 1) PetscCall(PetscOptionsEnum("-da_bd_y", "Boundary type for y direction", "DMDASetBoundaryType", DMBoundaryTypes, (PetscEnum)dd->by, (PetscEnum *)&dd->by, NULL));
 36:   if (dim > 2) PetscCall(PetscOptionsEnum("-da_bd_z", "Boundary type for z direction", "DMDASetBoundaryType", DMBoundaryTypes, (PetscEnum)dd->bz, (PetscEnum *)&dd->bz, NULL));
 37:   PetscCall(PetscOptionsEnum("-da_bd_all", "Boundary type for every direction", "DMDASetBoundaryType", DMBoundaryTypes, (PetscEnum)bt, (PetscEnum *)&bt, &flg));
 38:   if (flg) PetscCall(DMDASetBoundaryType(da, bt, bt, bt));
 39:   /* Handle DMDA refinement */
 40:   PetscCall(PetscOptionsBoundedInt("-da_refine_x", "Refinement ratio in x direction", "DMDASetRefinementFactor", dd->refine_x, &dd->refine_x, NULL, 1));
 41:   if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_refine_y", "Refinement ratio in y direction", "DMDASetRefinementFactor", dd->refine_y, &dd->refine_y, NULL, 1));
 42:   if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_refine_z", "Refinement ratio in z direction", "DMDASetRefinementFactor", dd->refine_z, &dd->refine_z, NULL, 1));
 43:   dd->coarsen_x = dd->refine_x;
 44:   dd->coarsen_y = dd->refine_y;
 45:   dd->coarsen_z = dd->refine_z;

 47:   /* Get refinement factors, defaults taken from the coarse DMDA */
 48:   PetscCall(DMDAGetRefinementFactor(da, &refx[0], &refy[0], &refz[0]));
 49:   for (i = 1; i < maxnlevels; i++) {
 50:     refx[i] = refx[0];
 51:     refy[i] = refy[0];
 52:     refz[i] = refz[0];
 53:   }
 54:   n = maxnlevels;
 55:   PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_x", "Refinement factor for each level", "None", refx, &n, &flg));
 56:   if (flg) {
 57:     dd->refine_x        = refx[0];
 58:     dd->refine_x_hier_n = n;
 59:     PetscCall(PetscMalloc1(n, &dd->refine_x_hier));
 60:     PetscCall(PetscArraycpy(dd->refine_x_hier, refx, n));
 61:   }
 62:   if (dim > 1) {
 63:     n = maxnlevels;
 64:     PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_y", "Refinement factor for each level", "None", refy, &n, &flg));
 65:     if (flg) {
 66:       dd->refine_y        = refy[0];
 67:       dd->refine_y_hier_n = n;
 68:       PetscCall(PetscMalloc1(n, &dd->refine_y_hier));
 69:       PetscCall(PetscArraycpy(dd->refine_y_hier, refy, n));
 70:     }
 71:   }
 72:   if (dim > 2) {
 73:     n = maxnlevels;
 74:     PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_z", "Refinement factor for each level", "None", refz, &n, &flg));
 75:     if (flg) {
 76:       dd->refine_z        = refz[0];
 77:       dd->refine_z_hier_n = n;
 78:       PetscCall(PetscMalloc1(n, &dd->refine_z_hier));
 79:       PetscCall(PetscArraycpy(dd->refine_z_hier, refz, n));
 80:     }
 81:   }

 83:   PetscCall(PetscOptionsBoundedInt("-da_refine", "Uniformly refine DA one or more times", "None", refine, &refine, NULL, 0));
 84:   PetscOptionsHeadEnd();

 86:   while (refine--) {
 87:     if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
 88:       PetscCall(PetscIntMultError(dd->refine_x, dd->M, &dd->M));
 89:     } else {
 90:       PetscCall(PetscIntMultError(dd->refine_x, dd->M - 1, &dd->M));
 91:       dd->M += 1;
 92:     }
 93:     if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
 94:       PetscCall(PetscIntMultError(dd->refine_y, dd->N, &dd->N));
 95:     } else {
 96:       PetscCall(PetscIntMultError(dd->refine_y, dd->N - 1, &dd->N));
 97:       dd->N += 1;
 98:     }
 99:     if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
100:       PetscCall(PetscIntMultError(dd->refine_z, dd->P, &dd->P));
101:     } else {
102:       PetscCall(PetscIntMultError(dd->refine_z, dd->P - 1, &dd->P));
103:       dd->P += 1;
104:     }
105:     da->levelup++;
106:     if (da->levelup - da->leveldown >= 0) {
107:       dd->refine_x = refx[da->levelup - da->leveldown];
108:       dd->refine_y = refy[da->levelup - da->leveldown];
109:       dd->refine_z = refz[da->levelup - da->leveldown];
110:     }
111:     if (da->levelup - da->leveldown >= 1) {
112:       dd->coarsen_x = refx[da->levelup - da->leveldown - 1];
113:       dd->coarsen_y = refy[da->levelup - da->leveldown - 1];
114:       dd->coarsen_z = refz[da->levelup - da->leveldown - 1];
115:     }
116:   }
117:   PetscFunctionReturn(PETSC_SUCCESS);
118: }

120: static PetscErrorCode DMLoad_DA(DM da, PetscViewer viewer)
121: {
122:   PetscInt        dim, m, n, p, dof, swidth;
123:   DMDAStencilType stencil;
124:   DMBoundaryType  bx, by, bz;
125:   PetscBool       coors;
126:   DM              dac;
127:   Vec             c;

129:   PetscFunctionBegin;
130:   PetscCall(PetscViewerBinaryRead(viewer, &dim, 1, NULL, PETSC_INT));
131:   PetscCall(PetscViewerBinaryRead(viewer, &m, 1, NULL, PETSC_INT));
132:   PetscCall(PetscViewerBinaryRead(viewer, &n, 1, NULL, PETSC_INT));
133:   PetscCall(PetscViewerBinaryRead(viewer, &p, 1, NULL, PETSC_INT));
134:   PetscCall(PetscViewerBinaryRead(viewer, &dof, 1, NULL, PETSC_INT));
135:   PetscCall(PetscViewerBinaryRead(viewer, &swidth, 1, NULL, PETSC_INT));
136:   PetscCall(PetscViewerBinaryRead(viewer, &bx, 1, NULL, PETSC_ENUM));
137:   PetscCall(PetscViewerBinaryRead(viewer, &by, 1, NULL, PETSC_ENUM));
138:   PetscCall(PetscViewerBinaryRead(viewer, &bz, 1, NULL, PETSC_ENUM));
139:   PetscCall(PetscViewerBinaryRead(viewer, &stencil, 1, NULL, PETSC_ENUM));

141:   PetscCall(DMSetDimension(da, dim));
142:   PetscCall(DMDASetSizes(da, m, n, p));
143:   PetscCall(DMDASetBoundaryType(da, bx, by, bz));
144:   PetscCall(DMDASetDof(da, dof));
145:   PetscCall(DMDASetStencilType(da, stencil));
146:   PetscCall(DMDASetStencilWidth(da, swidth));
147:   PetscCall(DMSetUp(da));
148:   PetscCall(PetscViewerBinaryRead(viewer, &coors, 1, NULL, PETSC_ENUM));
149:   if (coors) {
150:     PetscCall(DMGetCoordinateDM(da, &dac));
151:     PetscCall(DMCreateGlobalVector(dac, &c));
152:     PetscCall(VecLoad(c, viewer));
153:     PetscCall(DMSetCoordinates(da, c));
154:     PetscCall(VecDestroy(&c));
155:   }
156:   PetscFunctionReturn(PETSC_SUCCESS);
157: }

159: static PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
160: {
161:   DM_DA *da = (DM_DA *)dm->data;

163:   PetscFunctionBegin;
164:   if (subdm) {
165:     PetscSF sf;
166:     Vec     coords;
167:     void   *ctx;
168:     /* Cannot use DMClone since the dof stuff is mixed in. Ugh
169:     PetscCall(DMClone(dm, subdm)); */
170:     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), subdm));
171:     PetscCall(DMGetPointSF(dm, &sf));
172:     PetscCall(DMSetPointSF(*subdm, sf));
173:     PetscCall(DMGetApplicationContext(dm, &ctx));
174:     PetscCall(DMSetApplicationContext(*subdm, ctx));
175:     PetscCall(DMGetCoordinatesLocal(dm, &coords));
176:     if (coords) {
177:       PetscCall(DMSetCoordinatesLocal(*subdm, coords));
178:     } else {
179:       PetscCall(DMGetCoordinates(dm, &coords));
180:       if (coords) PetscCall(DMSetCoordinates(*subdm, coords));
181:     }

183:     PetscCall(DMSetType(*subdm, DMDA));
184:     PetscCall(DMSetDimension(*subdm, dm->dim));
185:     PetscCall(DMDASetSizes(*subdm, da->M, da->N, da->P));
186:     PetscCall(DMDASetNumProcs(*subdm, da->m, da->n, da->p));
187:     PetscCall(DMDASetBoundaryType(*subdm, da->bx, da->by, da->bz));
188:     PetscCall(DMDASetDof(*subdm, numFields));
189:     PetscCall(DMDASetStencilType(*subdm, da->stencil_type));
190:     PetscCall(DMDASetStencilWidth(*subdm, da->s));
191:     PetscCall(DMDASetOwnershipRanges(*subdm, da->lx, da->ly, da->lz));
192:   }
193:   if (is) {
194:     PetscInt *indices, cnt = 0, dof = da->w, i, j;

196:     PetscCall(PetscMalloc1(da->Nlocal * numFields / dof, &indices));
197:     for (i = da->base / dof; i < (da->base + da->Nlocal) / dof; ++i) {
198:       for (j = 0; j < numFields; ++j) indices[cnt++] = dof * i + fields[j];
199:     }
200:     PetscCheck(cnt == da->Nlocal * numFields / dof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count %" PetscInt_FMT " does not equal expected value %" PetscInt_FMT, cnt, da->Nlocal * numFields / dof);
201:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), cnt, indices, PETSC_OWN_POINTER, is));
202:   }
203:   PetscFunctionReturn(PETSC_SUCCESS);
204: }

206: static PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
207: {
208:   PetscInt i;
209:   DM_DA   *dd  = (DM_DA *)dm->data;
210:   PetscInt dof = dd->w;

212:   PetscFunctionBegin;
213:   if (len) *len = dof;
214:   if (islist) {
215:     Vec      v;
216:     PetscInt rstart, n;

218:     PetscCall(DMGetGlobalVector(dm, &v));
219:     PetscCall(VecGetOwnershipRange(v, &rstart, NULL));
220:     PetscCall(VecGetLocalSize(v, &n));
221:     PetscCall(DMRestoreGlobalVector(dm, &v));
222:     PetscCall(PetscMalloc1(dof, islist));
223:     for (i = 0; i < dof; i++) PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), n / dof, rstart + i, dof, &(*islist)[i]));
224:   }
225:   if (namelist) {
226:     PetscCall(PetscMalloc1(dof, namelist));
227:     if (dd->fieldname) {
228:       for (i = 0; i < dof; i++) PetscCall(PetscStrallocpy(dd->fieldname[i], &(*namelist)[i]));
229:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Currently DMDA must have fieldnames");
230:   }
231:   if (dmlist) {
232:     DM da;

234:     PetscCall(DMDACreate(PetscObjectComm((PetscObject)dm), &da));
235:     PetscCall(DMSetDimension(da, dm->dim));
236:     PetscCall(DMDASetSizes(da, dd->M, dd->N, dd->P));
237:     PetscCall(DMDASetNumProcs(da, dd->m, dd->n, dd->p));
238:     PetscCall(DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz));
239:     PetscCall(DMDASetDof(da, 1));
240:     PetscCall(DMDASetStencilType(da, dd->stencil_type));
241:     PetscCall(DMDASetStencilWidth(da, dd->s));
242:     PetscCall(DMSetUp(da));
243:     PetscCall(PetscMalloc1(dof, dmlist));
244:     for (i = 0; i < dof - 1; i++) PetscCall(PetscObjectReference((PetscObject)da));
245:     for (i = 0; i < dof; i++) (*dmlist)[i] = da;
246:   }
247:   PetscFunctionReturn(PETSC_SUCCESS);
248: }

250: static PetscErrorCode DMClone_DA(DM dm, DM *newdm)
251: {
252:   DM_DA *da = (DM_DA *)dm->data;

254:   PetscFunctionBegin;
255:   PetscCall(DMSetType(*newdm, DMDA));
256:   PetscCall(DMSetDimension(*newdm, dm->dim));
257:   PetscCall(DMDASetSizes(*newdm, da->M, da->N, da->P));
258:   PetscCall(DMDASetNumProcs(*newdm, da->m, da->n, da->p));
259:   PetscCall(DMDASetBoundaryType(*newdm, da->bx, da->by, da->bz));
260:   PetscCall(DMDASetDof(*newdm, da->w));
261:   PetscCall(DMDASetStencilType(*newdm, da->stencil_type));
262:   PetscCall(DMDASetStencilWidth(*newdm, da->s));
263:   PetscCall(DMDASetOwnershipRanges(*newdm, da->lx, da->ly, da->lz));
264:   PetscCall(DMSetUp(*newdm));
265:   PetscFunctionReturn(PETSC_SUCCESS);
266: }

268: static PetscErrorCode DMHasCreateInjection_DA(DM dm, PetscBool *flg)
269: {
270:   DM_DA *da = (DM_DA *)dm->data;

272:   PetscFunctionBegin;
274:   PetscAssertPointer(flg, 2);
275:   *flg = da->interptype == DMDA_Q1 ? PETSC_TRUE : PETSC_FALSE;
276:   PetscFunctionReturn(PETSC_SUCCESS);
277: }

279: static PetscErrorCode DMGetDimPoints_DA(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
280: {
281:   PetscFunctionBegin;
282:   PetscCall(DMDAGetDepthStratum(dm, dim, pStart, pEnd));
283:   PetscFunctionReturn(PETSC_SUCCESS);
284: }

286: static PetscErrorCode DMGetNeighbors_DA(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
287: {
288:   PetscInt        dim;
289:   DMDAStencilType st;

291:   PetscFunctionBegin;
292:   PetscCall(DMDAGetNeighbors(dm, ranks));
293:   PetscCall(DMDAGetInfo(dm, &dim, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &st));

295:   switch (dim) {
296:   case 1:
297:     *nranks = 3;
298:     /* if (st == DMDA_STENCIL_STAR) *nranks = 3; */
299:     break;
300:   case 2:
301:     *nranks = 9;
302:     /* if (st == DMDA_STENCIL_STAR) *nranks = 5; */
303:     break;
304:   case 3:
305:     *nranks = 27;
306:     /* if (st == DMDA_STENCIL_STAR) *nranks = 7; */
307:     break;
308:   default:
309:     break;
310:   }
311:   PetscFunctionReturn(PETSC_SUCCESS);
312: }

314: /*MC
315:    DMDA = "da" - A `DM` object that is used to manage data for a structured grid in 1, 2, or 3 dimensions.
316:          In the global representation of the vector each process stores a non-overlapping rectangular (or slab in 3d) portion of the grid points.
317:          In the local representation these rectangular regions (slabs) are extended in all directions by a stencil width set with `DMDASetStencilWidth()`.

319:          The vectors can be thought of as either cell centered or vertex centered on the mesh. But some variables cannot be cell centered and others
320:          vertex centered; see the documentation for `DMSTAG`, a similar `DM` implementation which supports more general staggered grids.

322:   Level: intermediate

324: .seealso: [](sec_struct), `DMType`, `DMCOMPOSITE`, `DMSTAG`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMDASetStencilWidth()`, `DMDASetStencilType()`,
325:           `DMDAStencilType`
326: M*/

328: PETSC_EXTERN PetscErrorCode DMCreate_DA(DM da)
329: {
330:   DM_DA *dd;

332:   PetscFunctionBegin;
333:   PetscAssertPointer(da, 1);
334:   PetscCall(PetscNew(&dd));
335:   da->data = dd;

337:   da->dim        = -1;
338:   dd->interptype = DMDA_Q1;
339:   dd->refine_x   = 2;
340:   dd->refine_y   = 2;
341:   dd->refine_z   = 2;
342:   dd->coarsen_x  = 2;
343:   dd->coarsen_y  = 2;
344:   dd->coarsen_z  = 2;
345:   dd->fieldname  = NULL;
346:   dd->nlocal     = -1;
347:   dd->Nlocal     = -1;
348:   dd->M          = -1;
349:   dd->N          = -1;
350:   dd->P          = -1;
351:   dd->m          = -1;
352:   dd->n          = -1;
353:   dd->p          = -1;
354:   dd->w          = -1;
355:   dd->s          = -1;

357:   dd->xs = -1;
358:   dd->xe = -1;
359:   dd->ys = -1;
360:   dd->ye = -1;
361:   dd->zs = -1;
362:   dd->ze = -1;
363:   dd->Xs = -1;
364:   dd->Xe = -1;
365:   dd->Ys = -1;
366:   dd->Ye = -1;
367:   dd->Zs = -1;
368:   dd->Ze = -1;

370:   dd->Nsub = 1;
371:   dd->xol  = 0;
372:   dd->yol  = 0;
373:   dd->zol  = 0;
374:   dd->xo   = 0;
375:   dd->yo   = 0;
376:   dd->zo   = 0;
377:   dd->Mo   = -1;
378:   dd->No   = -1;
379:   dd->Po   = -1;

381:   dd->gtol = NULL;
382:   dd->ltol = NULL;
383:   dd->ao   = NULL;
384:   PetscCall(PetscStrallocpy(AOBASIC, (char **)&dd->aotype));
385:   dd->base         = -1;
386:   dd->bx           = DM_BOUNDARY_NONE;
387:   dd->by           = DM_BOUNDARY_NONE;
388:   dd->bz           = DM_BOUNDARY_NONE;
389:   dd->stencil_type = DMDA_STENCIL_BOX;
390:   dd->interptype   = DMDA_Q1;
391:   dd->lx           = NULL;
392:   dd->ly           = NULL;
393:   dd->lz           = NULL;

395:   dd->elementtype = DMDA_ELEMENT_Q1;

397:   da->ops->globaltolocalbegin        = DMGlobalToLocalBegin_DA;
398:   da->ops->globaltolocalend          = DMGlobalToLocalEnd_DA;
399:   da->ops->localtoglobalbegin        = DMLocalToGlobalBegin_DA;
400:   da->ops->localtoglobalend          = DMLocalToGlobalEnd_DA;
401:   da->ops->localtolocalbegin         = DMLocalToLocalBegin_DA;
402:   da->ops->localtolocalend           = DMLocalToLocalEnd_DA;
403:   da->ops->createglobalvector        = DMCreateGlobalVector_DA;
404:   da->ops->createlocalvector         = DMCreateLocalVector_DA;
405:   da->ops->createinterpolation       = DMCreateInterpolation_DA;
406:   da->ops->getcoloring               = DMCreateColoring_DA;
407:   da->ops->creatematrix              = DMCreateMatrix_DA;
408:   da->ops->refine                    = DMRefine_DA;
409:   da->ops->coarsen                   = DMCoarsen_DA;
410:   da->ops->refinehierarchy           = DMRefineHierarchy_DA;
411:   da->ops->coarsenhierarchy          = DMCoarsenHierarchy_DA;
412:   da->ops->createinjection           = DMCreateInjection_DA;
413:   da->ops->hascreateinjection        = DMHasCreateInjection_DA;
414:   da->ops->destroy                   = DMDestroy_DA;
415:   da->ops->view                      = NULL;
416:   da->ops->setfromoptions            = DMSetFromOptions_DA;
417:   da->ops->setup                     = DMSetUp_DA;
418:   da->ops->clone                     = DMClone_DA;
419:   da->ops->load                      = DMLoad_DA;
420:   da->ops->createcoordinatedm        = DMCreateCoordinateDM_DA;
421:   da->ops->createsubdm               = DMCreateSubDM_DA;
422:   da->ops->createfielddecomposition  = DMCreateFieldDecomposition_DA;
423:   da->ops->createdomaindecomposition = DMCreateDomainDecomposition_DA;
424:   da->ops->createddscatters          = DMCreateDomainDecompositionScatters_DA;
425:   da->ops->getdimpoints              = DMGetDimPoints_DA;
426:   da->ops->getneighbors              = DMGetNeighbors_DA;
427:   da->ops->getlocalboundingbox       = DMGetLocalBoundingBox_DA;
428:   da->ops->locatepoints              = DMLocatePoints_DA_Regular;
429:   da->ops->getcompatibility          = DMGetCompatibility_DA;
430:   PetscCall(PetscObjectComposeFunction((PetscObject)da, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_DMDA));
431:   PetscFunctionReturn(PETSC_SUCCESS);
432: }

434: /*@
435:   DMDACreate - Creates a `DMDA` object for managing structured grids.

437:   Collective

439:   Input Parameter:
440: . comm - The communicator for the `DMDA` object

442:   Output Parameter:
443: . da - the `DMDA` object

445:   Level: advanced

447:   Notes:
448:   See [](sec_struct) for details on the construction of a `DMDA`

450:   `DMDACreate1d()`, `DMDACreate2d()`, and `DMDACreate3d()` are convenience routines to quickly completely create a `DMDA`

452: .seealso: [](sec_struct), `DM`, `DMDA`, `DMSetUp()`, `DMDASetSizes()`, `DMClone()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
453: @*/
454: PetscErrorCode DMDACreate(MPI_Comm comm, DM *da)
455: {
456:   PetscFunctionBegin;
457:   PetscAssertPointer(da, 2);
458:   PetscCall(DMCreate(comm, da));
459:   PetscCall(DMSetType(*da, DMDA));
460:   PetscFunctionReturn(PETSC_SUCCESS);
461: }