xref: /petsc/src/dm/impls/da/dacreate.c (revision 72fd0fbdbfed0156e00d46940d6a00aa32e1f56c)
1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"   I*/
247c6ae99SBarry Smith 
3ba38deedSJacob Faibussowitsch static PetscErrorCode DMSetFromOptions_DA(DM da, PetscOptionItems *PetscOptionsObject)
4d71ae5a4SJacob Faibussowitsch {
547c6ae99SBarry Smith   DM_DA         *dd     = (DM_DA *)da->data;
6c73cfb54SMatthew G. Knepley   PetscInt       refine = 0, dim = da->dim, maxnlevels = 100, refx[100], refy[100], refz[100], n, i;
7fe58071bSMatthew G. Knepley   DMBoundaryType bt = DM_BOUNDARY_NONE;
8897f7067SBarry Smith   PetscBool      flg;
947c6ae99SBarry Smith 
1047c6ae99SBarry Smith   PetscFunctionBegin;
1108401ef6SPierre Jolivet   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");
1208401ef6SPierre Jolivet   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");
1308401ef6SPierre Jolivet   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");
14235683edSBarry Smith 
15d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "DMDA Options");
169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-da_grid_x", "Number of grid points in x direction", "DMDASetSizes", dd->M, &dd->M, NULL, 1));
179566063dSJacob Faibussowitsch   if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_grid_y", "Number of grid points in y direction", "DMDASetSizes", dd->N, &dd->N, NULL, 1));
189566063dSJacob Faibussowitsch   if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_grid_z", "Number of grid points in z direction", "DMDASetSizes", dd->P, &dd->P, NULL, 1));
197ddda789SPeter Brune 
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-da_overlap", "Decomposition overlap in all directions", "DMDASetOverlap", dd->xol, &dd->xol, &flg, 0));
219566063dSJacob Faibussowitsch   if (flg) PetscCall(DMDASetOverlap(da, dd->xol, dd->xol, dd->xol));
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-da_overlap_x", "Decomposition overlap in x direction", "DMDASetOverlap", dd->xol, &dd->xol, NULL, 0));
239566063dSJacob Faibussowitsch   if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_overlap_y", "Decomposition overlap in y direction", "DMDASetOverlap", dd->yol, &dd->yol, NULL, 0));
249566063dSJacob Faibussowitsch   if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_overlap_z", "Decomposition overlap in z direction", "DMDASetOverlap", dd->zol, &dd->zol, NULL, 0));
253e7870d2SPeter Brune 
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-da_local_subdomains", "", "DMDASetNumLocalSubdomains", dd->Nsub, &dd->Nsub, &flg, PETSC_DECIDE));
279566063dSJacob Faibussowitsch   if (flg) PetscCall(DMDASetNumLocalSubDomains(da, dd->Nsub));
283e7870d2SPeter Brune 
298d07fd27SPatrick Sanan   /* Handle DMDA parallel distribution */
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-da_processors_x", "Number of processors in x direction", "DMDASetNumProcs", dd->m, &dd->m, NULL, PETSC_DECIDE));
319566063dSJacob Faibussowitsch   if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_processors_y", "Number of processors in y direction", "DMDASetNumProcs", dd->n, &dd->n, NULL, PETSC_DECIDE));
329566063dSJacob Faibussowitsch   if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_processors_z", "Number of processors in z direction", "DMDASetNumProcs", dd->p, &dd->p, NULL, PETSC_DECIDE));
33fe58071bSMatthew G. Knepley   // Handle boundaries
34fe58071bSMatthew G. Knepley   PetscCall(PetscOptionsEnum("-da_bd_x", "Boundary type for x direction", "DMDASetBoundaryType", DMBoundaryTypes, (PetscEnum)dd->bx, (PetscEnum *)&dd->bx, NULL));
35fe58071bSMatthew G. Knepley   if (dim > 1) PetscCall(PetscOptionsEnum("-da_bd_y", "Boundary type for y direction", "DMDASetBoundaryType", DMBoundaryTypes, (PetscEnum)dd->by, (PetscEnum *)&dd->by, NULL));
36fe58071bSMatthew G. Knepley   if (dim > 2) PetscCall(PetscOptionsEnum("-da_bd_z", "Boundary type for z direction", "DMDASetBoundaryType", DMBoundaryTypes, (PetscEnum)dd->bz, (PetscEnum *)&dd->bz, NULL));
37fe58071bSMatthew G. Knepley   PetscCall(PetscOptionsEnum("-da_bd_all", "Boundary type for every direction", "DMDASetBoundaryType", DMBoundaryTypes, (PetscEnum)bt, (PetscEnum *)&bt, &flg));
38fe58071bSMatthew G. Knepley   if (flg) PetscCall(DMDASetBoundaryType(da, bt, bt, bt));
39aa219208SBarry Smith   /* Handle DMDA refinement */
409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-da_refine_x", "Refinement ratio in x direction", "DMDASetRefinementFactor", dd->refine_x, &dd->refine_x, NULL, 1));
419566063dSJacob Faibussowitsch   if (dim > 1) PetscCall(PetscOptionsBoundedInt("-da_refine_y", "Refinement ratio in y direction", "DMDASetRefinementFactor", dd->refine_y, &dd->refine_y, NULL, 1));
429566063dSJacob Faibussowitsch   if (dim > 2) PetscCall(PetscOptionsBoundedInt("-da_refine_z", "Refinement ratio in z direction", "DMDASetRefinementFactor", dd->refine_z, &dd->refine_z, NULL, 1));
439371c9d4SSatish Balay   dd->coarsen_x = dd->refine_x;
449371c9d4SSatish Balay   dd->coarsen_y = dd->refine_y;
459371c9d4SSatish Balay   dd->coarsen_z = dd->refine_z;
4647c6ae99SBarry Smith 
47397b6216SJed Brown   /* Get refinement factors, defaults taken from the coarse DMDA */
489566063dSJacob Faibussowitsch   PetscCall(DMDAGetRefinementFactor(da, &refx[0], &refy[0], &refz[0]));
49397b6216SJed Brown   for (i = 1; i < maxnlevels; i++) {
50397b6216SJed Brown     refx[i] = refx[0];
51397b6216SJed Brown     refy[i] = refy[0];
52397b6216SJed Brown     refz[i] = refz[0];
53397b6216SJed Brown   }
54397b6216SJed Brown   n = maxnlevels;
559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_x", "Refinement factor for each level", "None", refx, &n, &flg));
56897f7067SBarry Smith   if (flg) {
57897f7067SBarry Smith     dd->refine_x        = refx[0];
58897f7067SBarry Smith     dd->refine_x_hier_n = n;
599566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &dd->refine_x_hier));
609566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd->refine_x_hier, refx, n));
61897f7067SBarry Smith   }
62c73cfb54SMatthew G. Knepley   if (dim > 1) {
63397b6216SJed Brown     n = maxnlevels;
649566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_y", "Refinement factor for each level", "None", refy, &n, &flg));
65897f7067SBarry Smith     if (flg) {
66897f7067SBarry Smith       dd->refine_y        = refy[0];
67897f7067SBarry Smith       dd->refine_y_hier_n = n;
689566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &dd->refine_y_hier));
699566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(dd->refine_y_hier, refy, n));
70897f7067SBarry Smith     }
71397b6216SJed Brown   }
72c73cfb54SMatthew G. Knepley   if (dim > 2) {
73397b6216SJed Brown     n = maxnlevels;
749566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-da_refine_hierarchy_z", "Refinement factor for each level", "None", refz, &n, &flg));
75897f7067SBarry Smith     if (flg) {
76897f7067SBarry Smith       dd->refine_z        = refz[0];
77897f7067SBarry Smith       dd->refine_z_hier_n = n;
789566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &dd->refine_z_hier));
799566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(dd->refine_z_hier, refz, n));
80897f7067SBarry Smith     }
81397b6216SJed Brown   }
82397b6216SJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-da_refine", "Uniformly refine DA one or more times", "None", refine, &refine, NULL, 0));
84d0609cedSBarry Smith   PetscOptionsHeadEnd();
85235683edSBarry Smith 
86e0f5d30fSBarry Smith   while (refine--) {
87bff4a2f0SMatthew G. Knepley     if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
887a3b7fd5SBarry Smith       PetscCall(PetscIntMultError(dd->refine_x, dd->M, &dd->M));
89e0f5d30fSBarry Smith     } else {
907a3b7fd5SBarry Smith       PetscCall(PetscIntMultError(dd->refine_x, dd->M - 1, &dd->M));
917a3b7fd5SBarry Smith       dd->M += 1;
92e0f5d30fSBarry Smith     }
93bff4a2f0SMatthew G. Knepley     if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
947a3b7fd5SBarry Smith       PetscCall(PetscIntMultError(dd->refine_y, dd->N, &dd->N));
95e0f5d30fSBarry Smith     } else {
967a3b7fd5SBarry Smith       PetscCall(PetscIntMultError(dd->refine_y, dd->N - 1, &dd->N));
977a3b7fd5SBarry Smith       dd->N += 1;
98e0f5d30fSBarry Smith     }
99bff4a2f0SMatthew G. Knepley     if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1007a3b7fd5SBarry Smith       PetscCall(PetscIntMultError(dd->refine_z, dd->P, &dd->P));
101e0f5d30fSBarry Smith     } else {
1027a3b7fd5SBarry Smith       PetscCall(PetscIntMultError(dd->refine_z, dd->P - 1, &dd->P));
1037a3b7fd5SBarry Smith       dd->P += 1;
104e0f5d30fSBarry Smith     }
105e0f5d30fSBarry Smith     da->levelup++;
10685fc4b34SJed Brown     if (da->levelup - da->leveldown >= 0) {
10785fc4b34SJed Brown       dd->refine_x = refx[da->levelup - da->leveldown];
10885fc4b34SJed Brown       dd->refine_y = refy[da->levelup - da->leveldown];
10985fc4b34SJed Brown       dd->refine_z = refz[da->levelup - da->leveldown];
11085fc4b34SJed Brown     }
11185fc4b34SJed Brown     if (da->levelup - da->leveldown >= 1) {
11285fc4b34SJed Brown       dd->coarsen_x = refx[da->levelup - da->leveldown - 1];
11385fc4b34SJed Brown       dd->coarsen_y = refy[da->levelup - da->leveldown - 1];
11485fc4b34SJed Brown       dd->coarsen_z = refz[da->levelup - da->leveldown - 1];
11585fc4b34SJed Brown     }
116e0f5d30fSBarry Smith   }
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11847c6ae99SBarry Smith }
11947c6ae99SBarry Smith 
120ba38deedSJacob Faibussowitsch static PetscErrorCode DMLoad_DA(DM da, PetscViewer viewer)
121d71ae5a4SJacob Faibussowitsch {
122060da220SMatthew G. Knepley   PetscInt        dim, m, n, p, dof, swidth;
123b859378eSBarry Smith   DMDAStencilType stencil;
124bff4a2f0SMatthew G. Knepley   DMBoundaryType  bx, by, bz;
125bc2bf880SBarry Smith   PetscBool       coors;
126bc2bf880SBarry Smith   DM              dac;
127bc2bf880SBarry Smith   Vec             c;
128b859378eSBarry Smith 
129b859378eSBarry Smith   PetscFunctionBegin;
1309566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &dim, 1, NULL, PETSC_INT));
1319566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &m, 1, NULL, PETSC_INT));
1329566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &n, 1, NULL, PETSC_INT));
1339566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &p, 1, NULL, PETSC_INT));
1349566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &dof, 1, NULL, PETSC_INT));
1359566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &swidth, 1, NULL, PETSC_INT));
1369566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &bx, 1, NULL, PETSC_ENUM));
1379566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &by, 1, NULL, PETSC_ENUM));
1389566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &bz, 1, NULL, PETSC_ENUM));
1399566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &stencil, 1, NULL, PETSC_ENUM));
140b859378eSBarry Smith 
1419566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(da, dim));
1429566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(da, m, n, p));
1439566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(da, bx, by, bz));
1449566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(da, dof));
1459566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilType(da, stencil));
1469566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(da, swidth));
1479566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
1489566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &coors, 1, NULL, PETSC_ENUM));
149bc2bf880SBarry Smith   if (coors) {
1509566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(da, &dac));
1519566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dac, &c));
1529566063dSJacob Faibussowitsch     PetscCall(VecLoad(c, viewer));
1539566063dSJacob Faibussowitsch     PetscCall(DMSetCoordinates(da, c));
1549566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&c));
155bc2bf880SBarry Smith   }
1563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
157b859378eSBarry Smith }
158b859378eSBarry Smith 
159ba38deedSJacob Faibussowitsch static PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
160d71ae5a4SJacob Faibussowitsch {
161d724dfffSBarry Smith   DM_DA *da = (DM_DA *)dm->data;
162d724dfffSBarry Smith 
163d724dfffSBarry Smith   PetscFunctionBegin;
164d724dfffSBarry Smith   if (subdm) {
165c38c1269SMatthew G. Knepley     PetscSF sf;
166c38c1269SMatthew G. Knepley     Vec     coords;
167c38c1269SMatthew G. Knepley     void   *ctx;
168c38c1269SMatthew G. Knepley     /* Cannot use DMClone since the dof stuff is mixed in. Ugh
1699566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, subdm)); */
1709566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), subdm));
1719566063dSJacob Faibussowitsch     PetscCall(DMGetPointSF(dm, &sf));
1729566063dSJacob Faibussowitsch     PetscCall(DMSetPointSF(*subdm, sf));
1739566063dSJacob Faibussowitsch     PetscCall(DMGetApplicationContext(dm, &ctx));
1749566063dSJacob Faibussowitsch     PetscCall(DMSetApplicationContext(*subdm, ctx));
1759566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coords));
176c38c1269SMatthew G. Knepley     if (coords) {
1779566063dSJacob Faibussowitsch       PetscCall(DMSetCoordinatesLocal(*subdm, coords));
178c38c1269SMatthew G. Knepley     } else {
1799566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinates(dm, &coords));
1809566063dSJacob Faibussowitsch       if (coords) PetscCall(DMSetCoordinates(*subdm, coords));
181c38c1269SMatthew G. Knepley     }
182c38c1269SMatthew G. Knepley 
1839566063dSJacob Faibussowitsch     PetscCall(DMSetType(*subdm, DMDA));
1849566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(*subdm, dm->dim));
1859566063dSJacob Faibussowitsch     PetscCall(DMDASetSizes(*subdm, da->M, da->N, da->P));
1869566063dSJacob Faibussowitsch     PetscCall(DMDASetNumProcs(*subdm, da->m, da->n, da->p));
1879566063dSJacob Faibussowitsch     PetscCall(DMDASetBoundaryType(*subdm, da->bx, da->by, da->bz));
1889566063dSJacob Faibussowitsch     PetscCall(DMDASetDof(*subdm, numFields));
1899566063dSJacob Faibussowitsch     PetscCall(DMDASetStencilType(*subdm, da->stencil_type));
1909566063dSJacob Faibussowitsch     PetscCall(DMDASetStencilWidth(*subdm, da->s));
1919566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(*subdm, da->lx, da->ly, da->lz));
192d724dfffSBarry Smith   }
193d724dfffSBarry Smith   if (is) {
194d724dfffSBarry Smith     PetscInt *indices, cnt = 0, dof = da->w, i, j;
19538221697SMatthew G. Knepley 
1969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(da->Nlocal * numFields / dof, &indices));
19738221697SMatthew G. Knepley     for (i = da->base / dof; i < (da->base + da->Nlocal) / dof; ++i) {
198ad540459SPierre Jolivet       for (j = 0; j < numFields; ++j) indices[cnt++] = dof * i + fields[j];
199d724dfffSBarry Smith     }
20063a3b9bcSJacob Faibussowitsch     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);
2019566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), cnt, indices, PETSC_OWN_POINTER, is));
202d724dfffSBarry Smith   }
2033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
204d724dfffSBarry Smith }
205d724dfffSBarry Smith 
206ba38deedSJacob Faibussowitsch static PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
207d71ae5a4SJacob Faibussowitsch {
2081c3fb106SBarry Smith   PetscInt i;
2091c3fb106SBarry Smith   DM_DA   *dd  = (DM_DA *)dm->data;
2101c3fb106SBarry Smith   PetscInt dof = dd->w;
2111c3fb106SBarry Smith 
2121c3fb106SBarry Smith   PetscFunctionBegin;
213731c8d9eSDmitry Karpeev   if (len) *len = dof;
2141c3fb106SBarry Smith   if (islist) {
2151c3fb106SBarry Smith     Vec      v;
2161c3fb106SBarry Smith     PetscInt rstart, n;
2171c3fb106SBarry Smith 
2189566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &v));
2199566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(v, &rstart, NULL));
2209566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(v, &n));
2219566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &v));
2229566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dof, islist));
22348a46eb9SPierre Jolivet     for (i = 0; i < dof; i++) PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), n / dof, rstart + i, dof, &(*islist)[i]));
2241c3fb106SBarry Smith   }
2251c3fb106SBarry Smith   if (namelist) {
2269566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dof, namelist));
2271c3fb106SBarry Smith     if (dd->fieldname) {
22848a46eb9SPierre Jolivet       for (i = 0; i < dof; i++) PetscCall(PetscStrallocpy(dd->fieldname[i], &(*namelist)[i]));
2291c3fb106SBarry Smith     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Currently DMDA must have fieldnames");
2301c3fb106SBarry Smith   }
2311c3fb106SBarry Smith   if (dmlist) {
2321c3fb106SBarry Smith     DM da;
2331c3fb106SBarry Smith 
2349566063dSJacob Faibussowitsch     PetscCall(DMDACreate(PetscObjectComm((PetscObject)dm), &da));
2359566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(da, dm->dim));
2369566063dSJacob Faibussowitsch     PetscCall(DMDASetSizes(da, dd->M, dd->N, dd->P));
2379566063dSJacob Faibussowitsch     PetscCall(DMDASetNumProcs(da, dd->m, dd->n, dd->p));
2389566063dSJacob Faibussowitsch     PetscCall(DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz));
2399566063dSJacob Faibussowitsch     PetscCall(DMDASetDof(da, 1));
2409566063dSJacob Faibussowitsch     PetscCall(DMDASetStencilType(da, dd->stencil_type));
2419566063dSJacob Faibussowitsch     PetscCall(DMDASetStencilWidth(da, dd->s));
2429566063dSJacob Faibussowitsch     PetscCall(DMSetUp(da));
2439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dof, dmlist));
2449566063dSJacob Faibussowitsch     for (i = 0; i < dof - 1; i++) PetscCall(PetscObjectReference((PetscObject)da));
2451c3fb106SBarry Smith     for (i = 0; i < dof; i++) (*dmlist)[i] = da;
2461c3fb106SBarry Smith   }
2473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2481c3fb106SBarry Smith }
2491c3fb106SBarry Smith 
250ba38deedSJacob Faibussowitsch static PetscErrorCode DMClone_DA(DM dm, DM *newdm)
251d71ae5a4SJacob Faibussowitsch {
25238221697SMatthew G. Knepley   DM_DA *da = (DM_DA *)dm->data;
25338221697SMatthew G. Knepley 
25438221697SMatthew G. Knepley   PetscFunctionBegin;
2559566063dSJacob Faibussowitsch   PetscCall(DMSetType(*newdm, DMDA));
2569566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*newdm, dm->dim));
2579566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(*newdm, da->M, da->N, da->P));
2589566063dSJacob Faibussowitsch   PetscCall(DMDASetNumProcs(*newdm, da->m, da->n, da->p));
2599566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(*newdm, da->bx, da->by, da->bz));
2609566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(*newdm, da->w));
2619566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilType(*newdm, da->stencil_type));
2629566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(*newdm, da->s));
2639566063dSJacob Faibussowitsch   PetscCall(DMDASetOwnershipRanges(*newdm, da->lx, da->ly, da->lz));
2649566063dSJacob Faibussowitsch   PetscCall(DMSetUp(*newdm));
2653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26638221697SMatthew G. Knepley }
26738221697SMatthew G. Knepley 
268d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMHasCreateInjection_DA(DM dm, PetscBool *flg)
269d71ae5a4SJacob Faibussowitsch {
2704a7a4c06SLawrence Mitchell   DM_DA *da = (DM_DA *)dm->data;
2714a7a4c06SLawrence Mitchell 
2724a7a4c06SLawrence Mitchell   PetscFunctionBegin;
2734a7a4c06SLawrence Mitchell   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2744f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
275131f56ebSLawrence Mitchell   *flg = da->interptype == DMDA_Q1 ? PETSC_TRUE : PETSC_FALSE;
2763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2774a7a4c06SLawrence Mitchell }
2784a7a4c06SLawrence Mitchell 
279d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetDimPoints_DA(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
280d71ae5a4SJacob Faibussowitsch {
281793f3fe5SMatthew G. Knepley   PetscFunctionBegin;
2829566063dSJacob Faibussowitsch   PetscCall(DMDAGetDepthStratum(dm, dim, pStart, pEnd));
2833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
284793f3fe5SMatthew G. Knepley }
285793f3fe5SMatthew G. Knepley 
286d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetNeighbors_DA(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
287d71ae5a4SJacob Faibussowitsch {
288502a2867SDave May   PetscInt        dim;
289502a2867SDave May   DMDAStencilType st;
290502a2867SDave May 
291502a2867SDave May   PetscFunctionBegin;
2929566063dSJacob Faibussowitsch   PetscCall(DMDAGetNeighbors(dm, ranks));
2939566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, &dim, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &st));
294502a2867SDave May 
295502a2867SDave May   switch (dim) {
296502a2867SDave May   case 1:
297502a2867SDave May     *nranks = 3;
298ad540459SPierre Jolivet     /* if (st == DMDA_STENCIL_STAR) *nranks = 3; */
299502a2867SDave May     break;
300502a2867SDave May   case 2:
301502a2867SDave May     *nranks = 9;
302ad540459SPierre Jolivet     /* if (st == DMDA_STENCIL_STAR) *nranks = 5; */
303502a2867SDave May     break;
304502a2867SDave May   case 3:
305502a2867SDave May     *nranks = 27;
306ad540459SPierre Jolivet     /* if (st == DMDA_STENCIL_STAR) *nranks = 7; */
307502a2867SDave May     break;
308d71ae5a4SJacob Faibussowitsch   default:
309d71ae5a4SJacob Faibussowitsch     break;
310502a2867SDave May   }
3113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
312502a2867SDave May }
313502a2867SDave May 
3143efe6655SBarry Smith /*MC
315dce8aebaSBarry Smith    DMDA = "da" - A `DM` object that is used to manage data for a structured grid in 1, 2, or 3 dimensions.
3163efe6655SBarry Smith          In the global representation of the vector each process stores a non-overlapping rectangular (or slab in 3d) portion of the grid points.
31712b4a537SBarry Smith          In the local representation these rectangular regions (slabs) are extended in all directions by a stencil width set with `DMDASetStencilWidth()`.
3183efe6655SBarry Smith 
3193efe6655SBarry Smith          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
32012b4a537SBarry Smith          vertex centered; see the documentation for `DMSTAG`, a similar `DM` implementation which supports more general staggered grids.
3213efe6655SBarry Smith 
3223efe6655SBarry Smith   Level: intermediate
3233efe6655SBarry Smith 
32412b4a537SBarry Smith .seealso: [](sec_struct), `DMType`, `DMCOMPOSITE`, `DMSTAG`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMDASetStencilWidth()`, `DMDASetStencilType()`,
32512b4a537SBarry Smith           `DMDAStencilType`
3263efe6655SBarry Smith M*/
3273efe6655SBarry Smith 
328d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_DA(DM da)
329d71ae5a4SJacob Faibussowitsch {
33047c6ae99SBarry Smith   DM_DA *dd;
33147c6ae99SBarry Smith 
33247c6ae99SBarry Smith   PetscFunctionBegin;
3334f572ea9SToby Isaac   PetscAssertPointer(da, 1);
3344dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&dd));
335a4121054SBarry Smith   da->data = dd;
33647c6ae99SBarry Smith 
337c73cfb54SMatthew G. Knepley   da->dim        = -1;
338aa219208SBarry Smith   dd->interptype = DMDA_Q1;
33947c6ae99SBarry Smith   dd->refine_x   = 2;
34047c6ae99SBarry Smith   dd->refine_y   = 2;
34147c6ae99SBarry Smith   dd->refine_z   = 2;
34281c108dcSJed Brown   dd->coarsen_x  = 2;
34381c108dcSJed Brown   dd->coarsen_y  = 2;
34481c108dcSJed Brown   dd->coarsen_z  = 2;
3450298fd71SBarry Smith   dd->fieldname  = NULL;
34647c6ae99SBarry Smith   dd->nlocal     = -1;
34747c6ae99SBarry Smith   dd->Nlocal     = -1;
34847c6ae99SBarry Smith   dd->M          = -1;
34947c6ae99SBarry Smith   dd->N          = -1;
35047c6ae99SBarry Smith   dd->P          = -1;
35147c6ae99SBarry Smith   dd->m          = -1;
35247c6ae99SBarry Smith   dd->n          = -1;
35347c6ae99SBarry Smith   dd->p          = -1;
35447c6ae99SBarry Smith   dd->w          = -1;
35547c6ae99SBarry Smith   dd->s          = -1;
3568865f1eaSKarl Rupp 
3579371c9d4SSatish Balay   dd->xs = -1;
3589371c9d4SSatish Balay   dd->xe = -1;
3599371c9d4SSatish Balay   dd->ys = -1;
3609371c9d4SSatish Balay   dd->ye = -1;
3619371c9d4SSatish Balay   dd->zs = -1;
3629371c9d4SSatish Balay   dd->ze = -1;
3639371c9d4SSatish Balay   dd->Xs = -1;
3649371c9d4SSatish Balay   dd->Xe = -1;
3659371c9d4SSatish Balay   dd->Ys = -1;
3669371c9d4SSatish Balay   dd->Ye = -1;
3679371c9d4SSatish Balay   dd->Zs = -1;
3689371c9d4SSatish Balay   dd->Ze = -1;
36947c6ae99SBarry Smith 
3703e7870d2SPeter Brune   dd->Nsub = 1;
3717ddda789SPeter Brune   dd->xol  = 0;
3727ddda789SPeter Brune   dd->yol  = 0;
3737ddda789SPeter Brune   dd->zol  = 0;
374d886c4f4SPeter Brune   dd->xo   = 0;
375d886c4f4SPeter Brune   dd->yo   = 0;
376d886c4f4SPeter Brune   dd->zo   = 0;
377e30e807fSPeter Brune   dd->Mo   = -1;
378e30e807fSPeter Brune   dd->No   = -1;
379e30e807fSPeter Brune   dd->Po   = -1;
38088661749SPeter Brune 
3810298fd71SBarry Smith   dd->gtol = NULL;
3820298fd71SBarry Smith   dd->ltol = NULL;
3830298fd71SBarry Smith   dd->ao   = NULL;
3843ba16761SJacob Faibussowitsch   PetscCall(PetscStrallocpy(AOBASIC, (char **)&dd->aotype));
38547c6ae99SBarry Smith   dd->base         = -1;
386bff4a2f0SMatthew G. Knepley   dd->bx           = DM_BOUNDARY_NONE;
387bff4a2f0SMatthew G. Knepley   dd->by           = DM_BOUNDARY_NONE;
388bff4a2f0SMatthew G. Knepley   dd->bz           = DM_BOUNDARY_NONE;
389aa219208SBarry Smith   dd->stencil_type = DMDA_STENCIL_BOX;
390aa219208SBarry Smith   dd->interptype   = DMDA_Q1;
3910298fd71SBarry Smith   dd->lx           = NULL;
3920298fd71SBarry Smith   dd->ly           = NULL;
3930298fd71SBarry Smith   dd->lz           = NULL;
39447c6ae99SBarry Smith 
395454e267fSLisandro Dalcin   dd->elementtype = DMDA_ELEMENT_Q1;
396454e267fSLisandro Dalcin 
397a4121054SBarry Smith   da->ops->globaltolocalbegin        = DMGlobalToLocalBegin_DA;
398a4121054SBarry Smith   da->ops->globaltolocalend          = DMGlobalToLocalEnd_DA;
399a4121054SBarry Smith   da->ops->localtoglobalbegin        = DMLocalToGlobalBegin_DA;
400a4121054SBarry Smith   da->ops->localtoglobalend          = DMLocalToGlobalEnd_DA;
401d78e899eSRichard Tran Mills   da->ops->localtolocalbegin         = DMLocalToLocalBegin_DA;
402d78e899eSRichard Tran Mills   da->ops->localtolocalend           = DMLocalToLocalEnd_DA;
403a4121054SBarry Smith   da->ops->createglobalvector        = DMCreateGlobalVector_DA;
404a4121054SBarry Smith   da->ops->createlocalvector         = DMCreateLocalVector_DA;
40525296bd5SBarry Smith   da->ops->createinterpolation       = DMCreateInterpolation_DA;
406e727c939SJed Brown   da->ops->getcoloring               = DMCreateColoring_DA;
40725296bd5SBarry Smith   da->ops->creatematrix              = DMCreateMatrix_DA;
408a4121054SBarry Smith   da->ops->refine                    = DMRefine_DA;
409a4121054SBarry Smith   da->ops->coarsen                   = DMCoarsen_DA;
410a4121054SBarry Smith   da->ops->refinehierarchy           = DMRefineHierarchy_DA;
411a4121054SBarry Smith   da->ops->coarsenhierarchy          = DMCoarsenHierarchy_DA;
4125a84ad33SLisandro Dalcin   da->ops->createinjection           = DMCreateInjection_DA;
4134a7a4c06SLawrence Mitchell   da->ops->hascreateinjection        = DMHasCreateInjection_DA;
414a4121054SBarry Smith   da->ops->destroy                   = DMDestroy_DA;
415ea78f98cSLisandro Dalcin   da->ops->view                      = NULL;
416a4121054SBarry Smith   da->ops->setfromoptions            = DMSetFromOptions_DA;
417a4121054SBarry Smith   da->ops->setup                     = DMSetUp_DA;
41838221697SMatthew G. Knepley   da->ops->clone                     = DMClone_DA;
419b859378eSBarry Smith   da->ops->load                      = DMLoad_DA;
4206636e97aSMatthew G Knepley   da->ops->createcoordinatedm        = DMCreateCoordinateDM_DA;
421d724dfffSBarry Smith   da->ops->createsubdm               = DMCreateSubDM_DA;
42216621825SDmitry Karpeev   da->ops->createfielddecomposition  = DMCreateFieldDecomposition_DA;
423e30e807fSPeter Brune   da->ops->createdomaindecomposition = DMCreateDomainDecomposition_DA;
424e30e807fSPeter Brune   da->ops->createddscatters          = DMCreateDomainDecompositionScatters_DA;
425793f3fe5SMatthew G. Knepley   da->ops->getdimpoints              = DMGetDimPoints_DA;
426502a2867SDave May   da->ops->getneighbors              = DMGetNeighbors_DA;
4276c6a6b79SMatthew G. Knepley   da->ops->getlocalboundingbox       = DMGetLocalBoundingBox_DA;
42805264a50SDave May   da->ops->locatepoints              = DMLocatePoints_DA_Regular;
4297c3cd84eSPatrick Sanan   da->ops->getcompatibility          = DMGetCompatibility_DA;
4309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)da, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_DMDA));
4313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
432a4121054SBarry Smith }
43347c6ae99SBarry Smith 
434a4121054SBarry Smith /*@
43512b4a537SBarry Smith   DMDACreate - Creates a `DMDA` object for managing structured grids.
436a4121054SBarry Smith 
4375edebc4bSPatrick Sanan   Collective
438a4121054SBarry Smith 
439a4121054SBarry Smith   Input Parameter:
44012b4a537SBarry Smith . comm - The communicator for the `DMDA` object
441a4121054SBarry Smith 
442a4121054SBarry Smith   Output Parameter:
443*72fd0fbdSBarry Smith . da - the `DMDA` object
444a4121054SBarry Smith 
445e0f5d30fSBarry Smith   Level: advanced
446e0f5d30fSBarry Smith 
447*72fd0fbdSBarry Smith   Notes:
448*72fd0fbdSBarry Smith   See [](sec_struct) for details on the construction of a `DMDA`
449a4121054SBarry Smith 
450*72fd0fbdSBarry Smith   `DMDACreate1d()`, `DMDACreate2d()`, and `DMDACreate3d()` are convenience routines to quickly completely create a `DMDA`
451*72fd0fbdSBarry Smith 
452*72fd0fbdSBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMSetUp()`, `DMDASetSizes()`, `DMClone()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
453a4121054SBarry Smith @*/
454d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDACreate(MPI_Comm comm, DM *da)
455d71ae5a4SJacob Faibussowitsch {
456a4121054SBarry Smith   PetscFunctionBegin;
4574f572ea9SToby Isaac   PetscAssertPointer(da, 2);
4589566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, da));
4599566063dSJacob Faibussowitsch   PetscCall(DMSetType(*da, DMDA));
4603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46147c6ae99SBarry Smith }
462