xref: /petsc/src/dm/impls/da/da3.c (revision d1e78c4f2649daee4d31c45b7350b666f6b9ac81)
17d0a6c19SBarry Smith 
247c6ae99SBarry Smith /*
347c6ae99SBarry Smith    Code for manipulating distributed regular 3d arrays in parallel.
447c6ae99SBarry Smith    File created by Peter Mell  7/14/95
547c6ae99SBarry Smith  */
647c6ae99SBarry Smith 
7af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"    I*/
847c6ae99SBarry Smith 
99804daf3SBarry Smith #include <petscdraw.h>
109371c9d4SSatish Balay static PetscErrorCode DMView_DA_3d(DM da, PetscViewer viewer) {
1147c6ae99SBarry Smith   PetscMPIInt rank;
12c9493c35SStefano Zampini   PetscBool   iascii, isdraw, isglvis, isbinary;
1347c6ae99SBarry Smith   DM_DA      *dd = (DM_DA *)da->data;
146858538eSMatthew G. Knepley   Vec         coordinates;
15*d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB)
169a42bb27SBarry Smith   PetscBool ismatlab;
179a42bb27SBarry Smith #endif
1847c6ae99SBarry Smith 
1947c6ae99SBarry Smith   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)da), &rank));
2147c6ae99SBarry Smith 
229566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
239566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
249566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
259566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
26*d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB)
279566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERMATLAB, &ismatlab));
289a42bb27SBarry Smith #endif
2947c6ae99SBarry Smith   if (iascii) {
3047c6ae99SBarry Smith     PetscViewerFormat format;
3147c6ae99SBarry Smith 
329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
339566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
3476a8abe0SBarry Smith     if (format == PETSC_VIEWER_LOAD_BALANCE) {
3576a8abe0SBarry Smith       PetscInt      i, nmax = 0, nmin = PETSC_MAX_INT, navg = 0, *nz, nzlocal;
3676a8abe0SBarry Smith       DMDALocalInfo info;
3776a8abe0SBarry Smith       PetscMPIInt   size;
389566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da), &size));
399566063dSJacob Faibussowitsch       PetscCall(DMDAGetLocalInfo(da, &info));
4076a8abe0SBarry Smith       nzlocal = info.xm * info.ym * info.zm;
419566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(size, &nz));
429566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allgather(&nzlocal, 1, MPIU_INT, nz, 1, MPIU_INT, PetscObjectComm((PetscObject)da)));
4376a8abe0SBarry Smith       for (i = 0; i < (PetscInt)size; i++) {
4476a8abe0SBarry Smith         nmax = PetscMax(nmax, nz[i]);
4576a8abe0SBarry Smith         nmin = PetscMin(nmin, nz[i]);
4676a8abe0SBarry Smith         navg += nz[i];
4776a8abe0SBarry Smith       }
489566063dSJacob Faibussowitsch       PetscCall(PetscFree(nz));
4976a8abe0SBarry Smith       navg = navg / size;
5063a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Load Balance - Grid Points: Min %" PetscInt_FMT "  avg %" PetscInt_FMT "  max %" PetscInt_FMT "\n", nmin, navg, nmax));
5176a8abe0SBarry Smith       PetscFunctionReturn(0);
5276a8abe0SBarry Smith     }
538ec8862eSJed Brown     if (format != PETSC_VIEWER_ASCII_VTK_DEPRECATED && format != PETSC_VIEWER_ASCII_VTK_CELL_DEPRECATED && format != PETSC_VIEWER_ASCII_GLVIS) {
54aa219208SBarry Smith       DMDALocalInfo info;
559566063dSJacob Faibussowitsch       PetscCall(DMDAGetLocalInfo(da, &info));
5663a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Processor [%d] M %" PetscInt_FMT " N %" PetscInt_FMT " P %" PetscInt_FMT " m %" PetscInt_FMT " n %" PetscInt_FMT " p %" PetscInt_FMT " w %" PetscInt_FMT " s %" PetscInt_FMT "\n", rank, dd->M, dd->N, dd->P, dd->m, dd->n, dd->p, dd->w, dd->s));
579371c9d4SSatish Balay       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "X range of indices: %" PetscInt_FMT " %" PetscInt_FMT ", Y range of indices: %" PetscInt_FMT " %" PetscInt_FMT ", Z range of indices: %" PetscInt_FMT " %" PetscInt_FMT "\n", info.xs,
589371c9d4SSatish Balay                                                    info.xs + info.xm, info.ys, info.ys + info.ym, info.zs, info.zs + info.zm));
596858538eSMatthew G. Knepley       PetscCall(DMGetCoordinates(da, &coordinates));
6047c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
616858538eSMatthew G. Knepley       if (coordinates) {
6247c6ae99SBarry Smith         PetscInt         last;
6347c6ae99SBarry Smith         const PetscReal *coors;
646858538eSMatthew G. Knepley         PetscCall(VecGetArrayRead(coordinates, &coors));
656858538eSMatthew G. Knepley         PetscCall(VecGetLocalSize(coordinates, &last));
6647c6ae99SBarry Smith         last = last - 3;
679566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Lower left corner %g %g %g : Upper right %g %g %g\n", (double)coors[0], (double)coors[1], (double)coors[2], (double)coors[last], (double)coors[last + 1], (double)coors[last + 2]));
686858538eSMatthew G. Knepley         PetscCall(VecRestoreArrayRead(coordinates, &coors));
6947c6ae99SBarry Smith       }
7047c6ae99SBarry Smith #endif
719566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
729566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
731baa6e33SBarry Smith     } else if (format == PETSC_VIEWER_ASCII_GLVIS) PetscCall(DMView_DA_GLVis(da, viewer));
741baa6e33SBarry Smith     else PetscCall(DMView_DA_VTK(da, viewer));
7547c6ae99SBarry Smith   } else if (isdraw) {
7647c6ae99SBarry Smith     PetscDraw       draw;
7747c6ae99SBarry Smith     PetscReal       ymin = -1.0, ymax = (PetscReal)dd->N;
7847c6ae99SBarry Smith     PetscReal       xmin = -1.0, xmax = (PetscReal)((dd->M + 2) * dd->P), x, y, ycoord, xcoord;
798ea3bf28SBarry Smith     PetscInt        k, plane, base;
808ea3bf28SBarry Smith     const PetscInt *idx;
8147c6ae99SBarry Smith     char            node[10];
8247c6ae99SBarry Smith     PetscBool       isnull;
8347c6ae99SBarry Smith 
849566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
859566063dSJacob Faibussowitsch     PetscCall(PetscDrawIsNull(draw, &isnull));
865b399a63SLisandro Dalcin     if (isnull) PetscFunctionReturn(0);
8747c6ae99SBarry Smith 
889566063dSJacob Faibussowitsch     PetscCall(PetscDrawCheckResizedWindow(draw));
899566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
909566063dSJacob Faibussowitsch     PetscCall(PetscDrawSetCoordinates(draw, xmin, ymin, xmax, ymax));
915b399a63SLisandro Dalcin 
92d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
9347c6ae99SBarry Smith     /* first processor draw all node lines */
94dd400576SPatrick Sanan     if (rank == 0) {
9547c6ae99SBarry Smith       for (k = 0; k < dd->P; k++) {
969371c9d4SSatish Balay         ymin = 0.0;
979371c9d4SSatish Balay         ymax = (PetscReal)(dd->N - 1);
9848a46eb9SPierre Jolivet         for (xmin = (PetscReal)(k * (dd->M + 1)); xmin < (PetscReal)(dd->M + (k * (dd->M + 1))); xmin++) PetscCall(PetscDrawLine(draw, xmin, ymin, xmin, ymax, PETSC_DRAW_BLACK));
999371c9d4SSatish Balay         xmin = (PetscReal)(k * (dd->M + 1));
1009371c9d4SSatish Balay         xmax = xmin + (PetscReal)(dd->M - 1);
10148a46eb9SPierre Jolivet         for (ymin = 0; ymin < (PetscReal)dd->N; ymin++) PetscCall(PetscDrawLine(draw, xmin, ymin, xmax, ymin, PETSC_DRAW_BLACK));
10247c6ae99SBarry Smith       }
10347c6ae99SBarry Smith     }
104d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
1059566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
1069566063dSJacob Faibussowitsch     PetscCall(PetscDrawPause(draw));
10747c6ae99SBarry Smith 
108d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
1095b399a63SLisandro Dalcin     /*Go through and draw for each plane*/
1105b399a63SLisandro Dalcin     for (k = 0; k < dd->P; k++) {
11147c6ae99SBarry Smith       if ((k >= dd->zs) && (k < dd->ze)) {
11247c6ae99SBarry Smith         /* draw my box */
11347c6ae99SBarry Smith         ymin = dd->ys;
11447c6ae99SBarry Smith         ymax = dd->ye - 1;
11547c6ae99SBarry Smith         xmin = dd->xs / dd->w + (dd->M + 1) * k;
11647c6ae99SBarry Smith         xmax = (dd->xe - 1) / dd->w + (dd->M + 1) * k;
11747c6ae99SBarry Smith 
1189566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw, xmin, ymin, xmax, ymin, PETSC_DRAW_RED));
1199566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw, xmin, ymin, xmin, ymax, PETSC_DRAW_RED));
1209566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw, xmin, ymax, xmax, ymax, PETSC_DRAW_RED));
1219566063dSJacob Faibussowitsch         PetscCall(PetscDrawLine(draw, xmax, ymin, xmax, ymax, PETSC_DRAW_RED));
12247c6ae99SBarry Smith 
12347c6ae99SBarry Smith         xmin = dd->xs / dd->w;
12447c6ae99SBarry Smith         xmax = (dd->xe - 1) / dd->w;
12547c6ae99SBarry Smith 
126832b7cebSLisandro Dalcin         /* identify which processor owns the box */
1279566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(node, sizeof(node), "%d", (int)rank));
1289566063dSJacob Faibussowitsch         PetscCall(PetscDrawString(draw, xmin + (dd->M + 1) * k + .2, ymin + .3, PETSC_DRAW_RED, node));
12947c6ae99SBarry Smith         /* put in numbers*/
13047c6ae99SBarry Smith         base = (dd->base + (dd->xe - dd->xs) * (dd->ye - dd->ys) * (k - dd->zs)) / dd->w;
13147c6ae99SBarry Smith         for (y = ymin; y <= ymax; y++) {
13247c6ae99SBarry Smith           for (x = xmin + (dd->M + 1) * k; x <= xmax + (dd->M + 1) * k; x++) {
1339566063dSJacob Faibussowitsch             PetscCall(PetscSNPrintf(node, sizeof(node), "%d", (int)base++));
1349566063dSJacob Faibussowitsch             PetscCall(PetscDrawString(draw, x, y, PETSC_DRAW_BLACK, node));
13547c6ae99SBarry Smith           }
13647c6ae99SBarry Smith         }
13747c6ae99SBarry Smith       }
13847c6ae99SBarry Smith     }
139d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
1409566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
1419566063dSJacob Faibussowitsch     PetscCall(PetscDrawPause(draw));
14247c6ae99SBarry Smith 
143d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
14447c6ae99SBarry Smith     for (k = 0 - dd->s; k < dd->P + dd->s; k++) {
14547c6ae99SBarry Smith       /* Go through and draw for each plane */
14647c6ae99SBarry Smith       if ((k >= dd->Zs) && (k < dd->Ze)) {
14747c6ae99SBarry Smith         /* overlay ghost numbers, useful for error checking */
148565245c5SBarry Smith         base = (dd->Xe - dd->Xs) * (dd->Ye - dd->Ys) * (k - dd->Zs) / dd->w;
1499566063dSJacob Faibussowitsch         PetscCall(ISLocalToGlobalMappingGetBlockIndices(da->ltogmap, &idx));
15047c6ae99SBarry Smith         plane = k;
151565245c5SBarry Smith         /* Keep z wrap around points on the drawing */
1528865f1eaSKarl Rupp         if (k < 0) plane = dd->P + k;
1538865f1eaSKarl Rupp         if (k >= dd->P) plane = k - dd->P;
1549371c9d4SSatish Balay         ymin = dd->Ys;
1559371c9d4SSatish Balay         ymax = dd->Ye;
15647c6ae99SBarry Smith         xmin = (dd->M + 1) * plane * dd->w;
15747c6ae99SBarry Smith         xmax = (dd->M + 1) * plane * dd->w + dd->M * dd->w;
15847c6ae99SBarry Smith         for (y = ymin; y < ymax; y++) {
15947c6ae99SBarry Smith           for (x = xmin + dd->Xs; x < xmin + dd->Xe; x += dd->w) {
1608ea3bf28SBarry Smith             sprintf(node, "%d", (int)(idx[base]));
16147c6ae99SBarry Smith             ycoord = y;
16247c6ae99SBarry Smith             /*Keep y wrap around points on drawing */
1638865f1eaSKarl Rupp             if (y < 0) ycoord = dd->N + y;
1648865f1eaSKarl Rupp             if (y >= dd->N) ycoord = y - dd->N;
16547c6ae99SBarry Smith             xcoord = x; /* Keep x wrap points on drawing */
1668865f1eaSKarl Rupp             if (x < xmin) xcoord = xmax - (xmin - x);
1678865f1eaSKarl Rupp             if (x >= xmax) xcoord = xmin + (x - xmax);
1689566063dSJacob Faibussowitsch             PetscCall(PetscDrawString(draw, xcoord / dd->w, ycoord, PETSC_DRAW_BLUE, node));
169565245c5SBarry Smith             base++;
17047c6ae99SBarry Smith           }
17147c6ae99SBarry Smith         }
1729566063dSJacob Faibussowitsch         PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap, &idx));
17347c6ae99SBarry Smith       }
17447c6ae99SBarry Smith     }
175d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
1769566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
1779566063dSJacob Faibussowitsch     PetscCall(PetscDrawPause(draw));
1789566063dSJacob Faibussowitsch     PetscCall(PetscDrawSave(draw));
179c9493c35SStefano Zampini   } else if (isglvis) {
1809566063dSJacob Faibussowitsch     PetscCall(DMView_DA_GLVis(da, viewer));
1819a42bb27SBarry Smith   } else if (isbinary) {
1829566063dSJacob Faibussowitsch     PetscCall(DMView_DA_Binary(da, viewer));
183*d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB)
1849a42bb27SBarry Smith   } else if (ismatlab) {
1859566063dSJacob Faibussowitsch     PetscCall(DMView_DA_Matlab(da, viewer));
1869a42bb27SBarry Smith #endif
18711aeaf0aSBarry Smith   }
18847c6ae99SBarry Smith   PetscFunctionReturn(0);
18947c6ae99SBarry Smith }
19047c6ae99SBarry Smith 
1919371c9d4SSatish Balay PetscErrorCode DMSetUp_DA_3D(DM da) {
19247c6ae99SBarry Smith   DM_DA          *dd           = (DM_DA *)da->data;
19347c6ae99SBarry Smith   const PetscInt  M            = dd->M;
19447c6ae99SBarry Smith   const PetscInt  N            = dd->N;
19547c6ae99SBarry Smith   const PetscInt  P            = dd->P;
19647c6ae99SBarry Smith   PetscInt        m            = dd->m;
19747c6ae99SBarry Smith   PetscInt        n            = dd->n;
19847c6ae99SBarry Smith   PetscInt        p            = dd->p;
19947c6ae99SBarry Smith   const PetscInt  dof          = dd->w;
20047c6ae99SBarry Smith   const PetscInt  s            = dd->s;
201bff4a2f0SMatthew G. Knepley   DMBoundaryType  bx           = dd->bx;
202bff4a2f0SMatthew G. Knepley   DMBoundaryType  by           = dd->by;
203bff4a2f0SMatthew G. Knepley   DMBoundaryType  bz           = dd->bz;
20419fd82e9SBarry Smith   DMDAStencilType stencil_type = dd->stencil_type;
20547c6ae99SBarry Smith   PetscInt       *lx           = dd->lx;
20647c6ae99SBarry Smith   PetscInt       *ly           = dd->ly;
20747c6ae99SBarry Smith   PetscInt       *lz           = dd->lz;
20847c6ae99SBarry Smith   MPI_Comm        comm;
20947c6ae99SBarry Smith   PetscMPIInt     rank, size;
210ce00eea3SSatish Balay   PetscInt        xs = 0, xe, ys = 0, ye, zs = 0, ze, x = 0, y = 0, z = 0;
211bd1fc5aeSBarry Smith   PetscInt        Xs, Xe, Ys, Ye, Zs, Ze, IXs, IXe, IYs, IYe, IZs, IZe, pm;
2128ea3bf28SBarry Smith   PetscInt        left, right, up, down, bottom, top, i, j, k, *idx, nn;
21347c6ae99SBarry Smith   PetscInt        n0, n1, n2, n3, n4, n5, n6, n7, n8, n9, n10, n11, n12, n14;
21447c6ae99SBarry Smith   PetscInt        n15, n16, n17, n18, n19, n20, n21, n22, n23, n24, n25, n26;
215db87c5ecSEthan Coon   PetscInt       *bases, *ldims, base, x_t, y_t, z_t, s_t, count, s_x, s_y, s_z;
21647c6ae99SBarry Smith   PetscInt        sn0 = 0, sn1 = 0, sn2 = 0, sn3 = 0, sn5 = 0, sn6 = 0, sn7 = 0;
21747c6ae99SBarry Smith   PetscInt        sn8 = 0, sn9 = 0, sn11 = 0, sn15 = 0, sn24 = 0, sn25 = 0, sn26 = 0;
21847c6ae99SBarry Smith   PetscInt        sn17 = 0, sn18 = 0, sn19 = 0, sn20 = 0, sn21 = 0, sn23 = 0;
21947c6ae99SBarry Smith   Vec             local, global;
220bd1fc5aeSBarry Smith   VecScatter      gtol;
22145b6f7e9SBarry Smith   IS              to, from;
2226f951b95Secoon   PetscBool       twod;
22347c6ae99SBarry Smith 
22447c6ae99SBarry Smith   PetscFunctionBegin;
2257a8be351SBarry Smith   PetscCheck(stencil_type != DMDA_STENCIL_BOX || (bx != DM_BOUNDARY_MIRROR && by != DM_BOUNDARY_MIRROR && bz != DM_BOUNDARY_MIRROR), PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Mirror boundary and box stencil");
2269566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
2273855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
22863a3b9bcSJacob Faibussowitsch   PetscCheck(((PetscInt64)M) * ((PetscInt64)N) * ((PetscInt64)P) * ((PetscInt64)dof) <= (PetscInt64)PETSC_MPI_INT_MAX, comm, PETSC_ERR_INT_OVERFLOW, "Mesh of %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT " (dof) is too large for 32 bit indices", M, N, P, dof);
2293855c12bSBarry Smith #endif
2303855c12bSBarry Smith 
2319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
23347c6ae99SBarry Smith 
23447c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
23563a3b9bcSJacob Faibussowitsch     PetscCheck(m >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in X direction: %" PetscInt_FMT, m);
236f7d195e4SLawrence Mitchell     PetscCheck(m <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in X direction: %" PetscInt_FMT " %d", m, size);
23747c6ae99SBarry Smith   }
23847c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
23963a3b9bcSJacob Faibussowitsch     PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in Y direction: %" PetscInt_FMT, n);
240f7d195e4SLawrence Mitchell     PetscCheck(n <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in Y direction: %" PetscInt_FMT " %d", n, size);
24147c6ae99SBarry Smith   }
24247c6ae99SBarry Smith   if (p != PETSC_DECIDE) {
24363a3b9bcSJacob Faibussowitsch     PetscCheck(p >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in Z direction: %" PetscInt_FMT, p);
244f7d195e4SLawrence Mitchell     PetscCheck(p <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in Z direction: %" PetscInt_FMT " %d", p, size);
24547c6ae99SBarry Smith   }
2461dca8a05SBarry Smith   PetscCheck(m <= 0 || n <= 0 || p <= 0 || m * n * p == size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "m %" PetscInt_FMT " * n %" PetscInt_FMT " * p %" PetscInt_FMT " != size %d", m, n, p, size);
24747c6ae99SBarry Smith 
24847c6ae99SBarry Smith   /* Partition the array among the processors */
24947c6ae99SBarry Smith   if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) {
25047c6ae99SBarry Smith     m = size / (n * p);
25147c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
25247c6ae99SBarry Smith     n = size / (m * p);
25347c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
25447c6ae99SBarry Smith     p = size / (m * n);
25547c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
25647c6ae99SBarry Smith     /* try for squarish distribution */
257369cc0aeSBarry Smith     m = (int)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)N * p)));
25847c6ae99SBarry Smith     if (!m) m = 1;
25947c6ae99SBarry Smith     while (m > 0) {
26047c6ae99SBarry Smith       n = size / (m * p);
26147c6ae99SBarry Smith       if (m * n * p == size) break;
26247c6ae99SBarry Smith       m--;
26347c6ae99SBarry Smith     }
26463a3b9bcSJacob Faibussowitsch     PetscCheck(m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad p value: p = %" PetscInt_FMT, p);
2659371c9d4SSatish Balay     if (M > N && m < n) {
2669371c9d4SSatish Balay       PetscInt _m = m;
2679371c9d4SSatish Balay       m           = n;
2689371c9d4SSatish Balay       n           = _m;
2699371c9d4SSatish Balay     }
27047c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
27147c6ae99SBarry Smith     /* try for squarish distribution */
272369cc0aeSBarry Smith     m = (int)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)P * n)));
27347c6ae99SBarry Smith     if (!m) m = 1;
27447c6ae99SBarry Smith     while (m > 0) {
27547c6ae99SBarry Smith       p = size / (m * n);
27647c6ae99SBarry Smith       if (m * n * p == size) break;
27747c6ae99SBarry Smith       m--;
27847c6ae99SBarry Smith     }
27963a3b9bcSJacob Faibussowitsch     PetscCheck(m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad n value: n = %" PetscInt_FMT, n);
2809371c9d4SSatish Balay     if (M > P && m < p) {
2819371c9d4SSatish Balay       PetscInt _m = m;
2829371c9d4SSatish Balay       m           = p;
2839371c9d4SSatish Balay       p           = _m;
2849371c9d4SSatish Balay     }
28547c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
28647c6ae99SBarry Smith     /* try for squarish distribution */
287369cc0aeSBarry Smith     n = (int)(0.5 + PetscSqrtReal(((PetscReal)N) * ((PetscReal)size) / ((PetscReal)P * m)));
28847c6ae99SBarry Smith     if (!n) n = 1;
28947c6ae99SBarry Smith     while (n > 0) {
29047c6ae99SBarry Smith       p = size / (m * n);
29147c6ae99SBarry Smith       if (m * n * p == size) break;
29247c6ae99SBarry Smith       n--;
29347c6ae99SBarry Smith     }
29463a3b9bcSJacob Faibussowitsch     PetscCheck(n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad m value: m = %" PetscInt_FMT, n);
2959371c9d4SSatish Balay     if (N > P && n < p) {
2969371c9d4SSatish Balay       PetscInt _n = n;
2979371c9d4SSatish Balay       n           = p;
2989371c9d4SSatish Balay       p           = _n;
2999371c9d4SSatish Balay     }
30047c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
30147c6ae99SBarry Smith     /* try for squarish distribution */
302369cc0aeSBarry Smith     n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N * N) * ((PetscReal)size) / ((PetscReal)P * M), (PetscReal)(1. / 3.)));
30347c6ae99SBarry Smith     if (!n) n = 1;
30447c6ae99SBarry Smith     while (n > 0) {
30547c6ae99SBarry Smith       pm = size / n;
30647c6ae99SBarry Smith       if (n * pm == size) break;
30747c6ae99SBarry Smith       n--;
30847c6ae99SBarry Smith     }
30947c6ae99SBarry Smith     if (!n) n = 1;
310369cc0aeSBarry Smith     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)P * n)));
31147c6ae99SBarry Smith     if (!m) m = 1;
31247c6ae99SBarry Smith     while (m > 0) {
31347c6ae99SBarry Smith       p = size / (m * n);
31447c6ae99SBarry Smith       if (m * n * p == size) break;
31547c6ae99SBarry Smith       m--;
31647c6ae99SBarry Smith     }
3179371c9d4SSatish Balay     if (M > P && m < p) {
3189371c9d4SSatish Balay       PetscInt _m = m;
3199371c9d4SSatish Balay       m           = p;
3209371c9d4SSatish Balay       p           = _m;
3219371c9d4SSatish Balay     }
32208401ef6SPierre Jolivet   } else PetscCheck(m * n * p == size, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Given Bad partition");
32347c6ae99SBarry Smith 
32408401ef6SPierre Jolivet   PetscCheck(m * n * p == size, PetscObjectComm((PetscObject)da), PETSC_ERR_PLIB, "Could not find good partition");
32563a3b9bcSJacob Faibussowitsch   PetscCheck(M >= m, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Partition in x direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT, M, m);
32663a3b9bcSJacob Faibussowitsch   PetscCheck(N >= n, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Partition in y direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT, N, n);
32763a3b9bcSJacob Faibussowitsch   PetscCheck(P >= p, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Partition in z direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT, P, p);
32847c6ae99SBarry Smith 
32947c6ae99SBarry Smith   /*
33047c6ae99SBarry Smith      Determine locally owned region
33147c6ae99SBarry Smith      [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes
33247c6ae99SBarry Smith   */
33347c6ae99SBarry Smith 
33447c6ae99SBarry Smith   if (!lx) {
3359566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &dd->lx));
33647c6ae99SBarry Smith     lx = dd->lx;
3378865f1eaSKarl Rupp     for (i = 0; i < m; i++) lx[i] = M / m + ((M % m) > (i % m));
33847c6ae99SBarry Smith   }
33947c6ae99SBarry Smith   x  = lx[rank % m];
34047c6ae99SBarry Smith   xs = 0;
3418865f1eaSKarl Rupp   for (i = 0; i < (rank % m); i++) xs += lx[i];
3421dca8a05SBarry Smith   PetscCheck(x >= s || (m <= 1 && bx != DM_BOUNDARY_PERIODIC), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local x-width of domain x %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT, x, s);
34347c6ae99SBarry Smith 
34447c6ae99SBarry Smith   if (!ly) {
3459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &dd->ly));
34647c6ae99SBarry Smith     ly = dd->ly;
3478865f1eaSKarl Rupp     for (i = 0; i < n; i++) ly[i] = N / n + ((N % n) > (i % n));
34847c6ae99SBarry Smith   }
34947c6ae99SBarry Smith   y = ly[(rank % (m * n)) / m];
3501dca8a05SBarry Smith   PetscCheck(y >= s || (n <= 1 && by != DM_BOUNDARY_PERIODIC), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local y-width of domain y %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT, y, s);
35130729d88SBarry Smith 
35247c6ae99SBarry Smith   ys = 0;
3538865f1eaSKarl Rupp   for (i = 0; i < (rank % (m * n)) / m; i++) ys += ly[i];
35447c6ae99SBarry Smith 
35547c6ae99SBarry Smith   if (!lz) {
3569566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(p, &dd->lz));
35747c6ae99SBarry Smith     lz = dd->lz;
3588865f1eaSKarl Rupp     for (i = 0; i < p; i++) lz[i] = P / p + ((P % p) > (i % p));
35947c6ae99SBarry Smith   }
36047c6ae99SBarry Smith   z = lz[rank / (m * n)];
361bcea557cSEthan Coon 
362fdc81ce8SEthan Coon   /* note this is different than x- and y-, as we will handle as an important special
363bff4a2f0SMatthew G. Knepley    case when p=P=1 and DM_BOUNDARY_PERIODIC and s > z.  This is to deal with 2D problems
364fdc81ce8SEthan Coon    in a 3D code.  Additional code for this case is noted with "2d case" comments */
3656f951b95Secoon   twod = PETSC_FALSE;
3668865f1eaSKarl Rupp   if (P == 1) twod = PETSC_TRUE;
3671dca8a05SBarry Smith   else PetscCheck(z >= s || (p <= 1 && bz != DM_BOUNDARY_PERIODIC), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local z-width of domain z %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT, z, s);
36847c6ae99SBarry Smith   zs = 0;
3698865f1eaSKarl Rupp   for (i = 0; i < (rank / (m * n)); i++) zs += lz[i];
37047c6ae99SBarry Smith   ye = ys + y;
37147c6ae99SBarry Smith   xe = xs + x;
37247c6ae99SBarry Smith   ze = zs + z;
37347c6ae99SBarry Smith 
37488661749SPeter Brune   /* determine ghost region (Xs) and region scattered into (IXs)  */
375d9c9ebe5SPeter Brune   if (xs - s > 0) {
3769371c9d4SSatish Balay     Xs  = xs - s;
3779371c9d4SSatish Balay     IXs = xs - s;
37888661749SPeter Brune   } else {
3798865f1eaSKarl Rupp     if (bx) Xs = xs - s;
3808865f1eaSKarl Rupp     else Xs = 0;
38188661749SPeter Brune     IXs = 0;
38288661749SPeter Brune   }
383d9c9ebe5SPeter Brune   if (xe + s <= M) {
3849371c9d4SSatish Balay     Xe  = xe + s;
3859371c9d4SSatish Balay     IXe = xe + s;
38688661749SPeter Brune   } else {
38788661749SPeter Brune     if (bx) {
3889371c9d4SSatish Balay       Xs = xs - s;
3899371c9d4SSatish Balay       Xe = xe + s;
3908865f1eaSKarl Rupp     } else Xe = M;
39188661749SPeter Brune     IXe = M;
39288661749SPeter Brune   }
39347c6ae99SBarry Smith 
394bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
395d9c9ebe5SPeter Brune     IXs = xs - s;
396d9c9ebe5SPeter Brune     IXe = xe + s;
397d9c9ebe5SPeter Brune     Xs  = xs - s;
398d9c9ebe5SPeter Brune     Xe  = xe + s;
39988661749SPeter Brune   }
40088661749SPeter Brune 
401d9c9ebe5SPeter Brune   if (ys - s > 0) {
4029371c9d4SSatish Balay     Ys  = ys - s;
4039371c9d4SSatish Balay     IYs = ys - s;
40488661749SPeter Brune   } else {
4058865f1eaSKarl Rupp     if (by) Ys = ys - s;
4068865f1eaSKarl Rupp     else Ys = 0;
40788661749SPeter Brune     IYs = 0;
40888661749SPeter Brune   }
409d9c9ebe5SPeter Brune   if (ye + s <= N) {
4109371c9d4SSatish Balay     Ye  = ye + s;
4119371c9d4SSatish Balay     IYe = ye + s;
41288661749SPeter Brune   } else {
4138865f1eaSKarl Rupp     if (by) Ye = ye + s;
4148865f1eaSKarl Rupp     else Ye = N;
41588661749SPeter Brune     IYe = N;
41688661749SPeter Brune   }
41788661749SPeter Brune 
418bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
419d9c9ebe5SPeter Brune     IYs = ys - s;
420d9c9ebe5SPeter Brune     IYe = ye + s;
421d9c9ebe5SPeter Brune     Ys  = ys - s;
422d9c9ebe5SPeter Brune     Ye  = ye + s;
42388661749SPeter Brune   }
42488661749SPeter Brune 
425d9c9ebe5SPeter Brune   if (zs - s > 0) {
4269371c9d4SSatish Balay     Zs  = zs - s;
4279371c9d4SSatish Balay     IZs = zs - s;
42888661749SPeter Brune   } else {
4298865f1eaSKarl Rupp     if (bz) Zs = zs - s;
4308865f1eaSKarl Rupp     else Zs = 0;
43188661749SPeter Brune     IZs = 0;
43288661749SPeter Brune   }
433d9c9ebe5SPeter Brune   if (ze + s <= P) {
4349371c9d4SSatish Balay     Ze  = ze + s;
4359371c9d4SSatish Balay     IZe = ze + s;
43688661749SPeter Brune   } else {
4378865f1eaSKarl Rupp     if (bz) Ze = ze + s;
4388865f1eaSKarl Rupp     else Ze = P;
43988661749SPeter Brune     IZe = P;
44088661749SPeter Brune   }
44188661749SPeter Brune 
442bff4a2f0SMatthew G. Knepley   if (bz == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_MIRROR) {
443d9c9ebe5SPeter Brune     IZs = zs - s;
444d9c9ebe5SPeter Brune     IZe = ze + s;
445d9c9ebe5SPeter Brune     Zs  = zs - s;
446d9c9ebe5SPeter Brune     Ze  = ze + s;
44788661749SPeter Brune   }
44847c6ae99SBarry Smith 
44947c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
450d9c9ebe5SPeter Brune   s_x = s;
451d9c9ebe5SPeter Brune   s_y = s;
452d9c9ebe5SPeter Brune   s_z = s;
45347c6ae99SBarry Smith 
45447c6ae99SBarry Smith   /* determine starting point of each processor */
45547c6ae99SBarry Smith   nn = x * y * z;
4569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(size + 1, &bases, size, &ldims));
4579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&nn, 1, MPIU_INT, ldims, 1, MPIU_INT, comm));
45847c6ae99SBarry Smith   bases[0] = 0;
4598865f1eaSKarl Rupp   for (i = 1; i <= size; i++) bases[i] = ldims[i - 1];
4608865f1eaSKarl Rupp   for (i = 1; i <= size; i++) bases[i] += bases[i - 1];
461ce00eea3SSatish Balay   base = bases[rank] * dof;
46247c6ae99SBarry Smith 
46347c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
464ce00eea3SSatish Balay   dd->Nlocal = x * y * z * dof;
4659566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(comm, dof, dd->Nlocal, PETSC_DECIDE, NULL, &global));
466ce00eea3SSatish Balay   dd->nlocal = (Xe - Xs) * (Ye - Ys) * (Ze - Zs) * dof;
4679566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, dof, dd->nlocal, NULL, &local));
46847c6ae99SBarry Smith 
4694104a7a0SPatrick Sanan   /* generate global to local vector scatter and local to global mapping*/
47047c6ae99SBarry Smith 
471ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
472ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
4739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((IXe - IXs) * (IYe - IYs) * (IZe - IZs), &idx));
474d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
4759371c9d4SSatish Balay     left   = IXs - Xs;
4769371c9d4SSatish Balay     right  = left + (IXe - IXs);
4779371c9d4SSatish Balay     bottom = IYs - Ys;
4789371c9d4SSatish Balay     top    = bottom + (IYe - IYs);
4799371c9d4SSatish Balay     down   = IZs - Zs;
4809371c9d4SSatish Balay     up     = down + (IZe - IZs);
481ce00eea3SSatish Balay     count  = 0;
482ce00eea3SSatish Balay     for (i = down; i < up; i++) {
483ce00eea3SSatish Balay       for (j = bottom; j < top; j++) {
484ad540459SPierre Jolivet         for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k;
485ce00eea3SSatish Balay       }
486ce00eea3SSatish Balay     }
4879566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(comm, dof, count, idx, PETSC_OWN_POINTER, &to));
48847c6ae99SBarry Smith   } else {
48947c6ae99SBarry Smith     /* This is way ugly! We need to list the funny cross type region */
4909371c9d4SSatish Balay     left   = xs - Xs;
4919371c9d4SSatish Balay     right  = left + x;
4929371c9d4SSatish Balay     bottom = ys - Ys;
4939371c9d4SSatish Balay     top    = bottom + y;
4949371c9d4SSatish Balay     down   = zs - Zs;
4959371c9d4SSatish Balay     up     = down + z;
49647c6ae99SBarry Smith     count  = 0;
497a5b23f4aSJose E. Roman     /* the bottom chunk */
498ce00eea3SSatish Balay     for (i = (IZs - Zs); i < down; i++) {
49947c6ae99SBarry Smith       for (j = bottom; j < top; j++) {
500ce00eea3SSatish Balay         for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k;
50147c6ae99SBarry Smith       }
50247c6ae99SBarry Smith     }
50347c6ae99SBarry Smith     /* the middle piece */
50447c6ae99SBarry Smith     for (i = down; i < up; i++) {
50547c6ae99SBarry Smith       /* front */
506ce00eea3SSatish Balay       for (j = (IYs - Ys); j < bottom; j++) {
507ce00eea3SSatish Balay         for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k;
50847c6ae99SBarry Smith       }
50947c6ae99SBarry Smith       /* middle */
51047c6ae99SBarry Smith       for (j = bottom; j < top; j++) {
511ce00eea3SSatish Balay         for (k = IXs - Xs; k < IXe - Xs; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k;
51247c6ae99SBarry Smith       }
51347c6ae99SBarry Smith       /* back */
514ce00eea3SSatish Balay       for (j = top; j < top + IYe - ye; j++) {
515ce00eea3SSatish Balay         for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k;
51647c6ae99SBarry Smith       }
51747c6ae99SBarry Smith     }
51847c6ae99SBarry Smith     /* the top piece */
519ce00eea3SSatish Balay     for (i = up; i < up + IZe - ze; i++) {
52047c6ae99SBarry Smith       for (j = bottom; j < top; j++) {
521ce00eea3SSatish Balay         for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k;
52247c6ae99SBarry Smith       }
52347c6ae99SBarry Smith     }
5249566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(comm, dof, count, idx, PETSC_OWN_POINTER, &to));
52547c6ae99SBarry Smith   }
52647c6ae99SBarry Smith 
52747c6ae99SBarry Smith   /* determine who lies on each side of use stored in    n24 n25 n26
52847c6ae99SBarry Smith                                                          n21 n22 n23
52947c6ae99SBarry Smith                                                          n18 n19 n20
53047c6ae99SBarry Smith 
53147c6ae99SBarry Smith                                                          n15 n16 n17
53247c6ae99SBarry Smith                                                          n12     n14
53347c6ae99SBarry Smith                                                          n9  n10 n11
53447c6ae99SBarry Smith 
53547c6ae99SBarry Smith                                                          n6  n7  n8
53647c6ae99SBarry Smith                                                          n3  n4  n5
53747c6ae99SBarry Smith                                                          n0  n1  n2
53847c6ae99SBarry Smith   */
53947c6ae99SBarry Smith 
54047c6ae99SBarry Smith   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
54147c6ae99SBarry Smith   /* Assume Nodes are Internal to the Cube */
54247c6ae99SBarry Smith   n0 = rank - m * n - m - 1;
54347c6ae99SBarry Smith   n1 = rank - m * n - m;
54447c6ae99SBarry Smith   n2 = rank - m * n - m + 1;
54547c6ae99SBarry Smith   n3 = rank - m * n - 1;
54647c6ae99SBarry Smith   n4 = rank - m * n;
54747c6ae99SBarry Smith   n5 = rank - m * n + 1;
54847c6ae99SBarry Smith   n6 = rank - m * n + m - 1;
54947c6ae99SBarry Smith   n7 = rank - m * n + m;
55047c6ae99SBarry Smith   n8 = rank - m * n + m + 1;
55147c6ae99SBarry Smith 
55247c6ae99SBarry Smith   n9  = rank - m - 1;
55347c6ae99SBarry Smith   n10 = rank - m;
55447c6ae99SBarry Smith   n11 = rank - m + 1;
55547c6ae99SBarry Smith   n12 = rank - 1;
55647c6ae99SBarry Smith   n14 = rank + 1;
55747c6ae99SBarry Smith   n15 = rank + m - 1;
55847c6ae99SBarry Smith   n16 = rank + m;
55947c6ae99SBarry Smith   n17 = rank + m + 1;
56047c6ae99SBarry Smith 
56147c6ae99SBarry Smith   n18 = rank + m * n - m - 1;
56247c6ae99SBarry Smith   n19 = rank + m * n - m;
56347c6ae99SBarry Smith   n20 = rank + m * n - m + 1;
56447c6ae99SBarry Smith   n21 = rank + m * n - 1;
56547c6ae99SBarry Smith   n22 = rank + m * n;
56647c6ae99SBarry Smith   n23 = rank + m * n + 1;
56747c6ae99SBarry Smith   n24 = rank + m * n + m - 1;
56847c6ae99SBarry Smith   n25 = rank + m * n + m;
56947c6ae99SBarry Smith   n26 = rank + m * n + m + 1;
57047c6ae99SBarry Smith 
57147c6ae99SBarry Smith   /* Assume Pieces are on Faces of Cube */
57247c6ae99SBarry Smith 
57347c6ae99SBarry Smith   if (xs == 0) { /* First assume not corner or edge */
57447c6ae99SBarry Smith     n0  = rank - 1 - (m * n);
57547c6ae99SBarry Smith     n3  = rank + m - 1 - (m * n);
57647c6ae99SBarry Smith     n6  = rank + 2 * m - 1 - (m * n);
57747c6ae99SBarry Smith     n9  = rank - 1;
57847c6ae99SBarry Smith     n12 = rank + m - 1;
57947c6ae99SBarry Smith     n15 = rank + 2 * m - 1;
58047c6ae99SBarry Smith     n18 = rank - 1 + (m * n);
58147c6ae99SBarry Smith     n21 = rank + m - 1 + (m * n);
58247c6ae99SBarry Smith     n24 = rank + 2 * m - 1 + (m * n);
58347c6ae99SBarry Smith   }
58447c6ae99SBarry Smith 
585ce00eea3SSatish Balay   if (xe == M) { /* First assume not corner or edge */
58647c6ae99SBarry Smith     n2  = rank - 2 * m + 1 - (m * n);
58747c6ae99SBarry Smith     n5  = rank - m + 1 - (m * n);
58847c6ae99SBarry Smith     n8  = rank + 1 - (m * n);
58947c6ae99SBarry Smith     n11 = rank - 2 * m + 1;
59047c6ae99SBarry Smith     n14 = rank - m + 1;
59147c6ae99SBarry Smith     n17 = rank + 1;
59247c6ae99SBarry Smith     n20 = rank - 2 * m + 1 + (m * n);
59347c6ae99SBarry Smith     n23 = rank - m + 1 + (m * n);
59447c6ae99SBarry Smith     n26 = rank + 1 + (m * n);
59547c6ae99SBarry Smith   }
59647c6ae99SBarry Smith 
59747c6ae99SBarry Smith   if (ys == 0) { /* First assume not corner or edge */
59847c6ae99SBarry Smith     n0  = rank + m * (n - 1) - 1 - (m * n);
59947c6ae99SBarry Smith     n1  = rank + m * (n - 1) - (m * n);
60047c6ae99SBarry Smith     n2  = rank + m * (n - 1) + 1 - (m * n);
60147c6ae99SBarry Smith     n9  = rank + m * (n - 1) - 1;
60247c6ae99SBarry Smith     n10 = rank + m * (n - 1);
60347c6ae99SBarry Smith     n11 = rank + m * (n - 1) + 1;
60447c6ae99SBarry Smith     n18 = rank + m * (n - 1) - 1 + (m * n);
60547c6ae99SBarry Smith     n19 = rank + m * (n - 1) + (m * n);
60647c6ae99SBarry Smith     n20 = rank + m * (n - 1) + 1 + (m * n);
60747c6ae99SBarry Smith   }
60847c6ae99SBarry Smith 
60947c6ae99SBarry Smith   if (ye == N) { /* First assume not corner or edge */
61047c6ae99SBarry Smith     n6  = rank - m * (n - 1) - 1 - (m * n);
61147c6ae99SBarry Smith     n7  = rank - m * (n - 1) - (m * n);
61247c6ae99SBarry Smith     n8  = rank - m * (n - 1) + 1 - (m * n);
61347c6ae99SBarry Smith     n15 = rank - m * (n - 1) - 1;
61447c6ae99SBarry Smith     n16 = rank - m * (n - 1);
61547c6ae99SBarry Smith     n17 = rank - m * (n - 1) + 1;
61647c6ae99SBarry Smith     n24 = rank - m * (n - 1) - 1 + (m * n);
61747c6ae99SBarry Smith     n25 = rank - m * (n - 1) + (m * n);
61847c6ae99SBarry Smith     n26 = rank - m * (n - 1) + 1 + (m * n);
61947c6ae99SBarry Smith   }
62047c6ae99SBarry Smith 
62147c6ae99SBarry Smith   if (zs == 0) { /* First assume not corner or edge */
62247c6ae99SBarry Smith     n0 = size - (m * n) + rank - m - 1;
62347c6ae99SBarry Smith     n1 = size - (m * n) + rank - m;
62447c6ae99SBarry Smith     n2 = size - (m * n) + rank - m + 1;
62547c6ae99SBarry Smith     n3 = size - (m * n) + rank - 1;
62647c6ae99SBarry Smith     n4 = size - (m * n) + rank;
62747c6ae99SBarry Smith     n5 = size - (m * n) + rank + 1;
62847c6ae99SBarry Smith     n6 = size - (m * n) + rank + m - 1;
62947c6ae99SBarry Smith     n7 = size - (m * n) + rank + m;
63047c6ae99SBarry Smith     n8 = size - (m * n) + rank + m + 1;
63147c6ae99SBarry Smith   }
63247c6ae99SBarry Smith 
63347c6ae99SBarry Smith   if (ze == P) { /* First assume not corner or edge */
63447c6ae99SBarry Smith     n18 = (m * n) - (size - rank) - m - 1;
63547c6ae99SBarry Smith     n19 = (m * n) - (size - rank) - m;
63647c6ae99SBarry Smith     n20 = (m * n) - (size - rank) - m + 1;
63747c6ae99SBarry Smith     n21 = (m * n) - (size - rank) - 1;
63847c6ae99SBarry Smith     n22 = (m * n) - (size - rank);
63947c6ae99SBarry Smith     n23 = (m * n) - (size - rank) + 1;
64047c6ae99SBarry Smith     n24 = (m * n) - (size - rank) + m - 1;
64147c6ae99SBarry Smith     n25 = (m * n) - (size - rank) + m;
64247c6ae99SBarry Smith     n26 = (m * n) - (size - rank) + m + 1;
64347c6ae99SBarry Smith   }
64447c6ae99SBarry Smith 
64547c6ae99SBarry Smith   if ((xs == 0) && (zs == 0)) { /* Assume an edge, not corner */
64647c6ae99SBarry Smith     n0 = size - m * n + rank + m - 1 - m;
64747c6ae99SBarry Smith     n3 = size - m * n + rank + m - 1;
64847c6ae99SBarry Smith     n6 = size - m * n + rank + m - 1 + m;
64947c6ae99SBarry Smith   }
65047c6ae99SBarry Smith 
65147c6ae99SBarry Smith   if ((xs == 0) && (ze == P)) { /* Assume an edge, not corner */
65247c6ae99SBarry Smith     n18 = m * n - (size - rank) + m - 1 - m;
65347c6ae99SBarry Smith     n21 = m * n - (size - rank) + m - 1;
65447c6ae99SBarry Smith     n24 = m * n - (size - rank) + m - 1 + m;
65547c6ae99SBarry Smith   }
65647c6ae99SBarry Smith 
65747c6ae99SBarry Smith   if ((xs == 0) && (ys == 0)) { /* Assume an edge, not corner */
65847c6ae99SBarry Smith     n0  = rank + m * n - 1 - m * n;
65947c6ae99SBarry Smith     n9  = rank + m * n - 1;
66047c6ae99SBarry Smith     n18 = rank + m * n - 1 + m * n;
66147c6ae99SBarry Smith   }
66247c6ae99SBarry Smith 
66347c6ae99SBarry Smith   if ((xs == 0) && (ye == N)) { /* Assume an edge, not corner */
66447c6ae99SBarry Smith     n6  = rank - m * (n - 1) + m - 1 - m * n;
66547c6ae99SBarry Smith     n15 = rank - m * (n - 1) + m - 1;
66647c6ae99SBarry Smith     n24 = rank - m * (n - 1) + m - 1 + m * n;
66747c6ae99SBarry Smith   }
66847c6ae99SBarry Smith 
669ce00eea3SSatish Balay   if ((xe == M) && (zs == 0)) { /* Assume an edge, not corner */
67047c6ae99SBarry Smith     n2 = size - (m * n - rank) - (m - 1) - m;
67147c6ae99SBarry Smith     n5 = size - (m * n - rank) - (m - 1);
67247c6ae99SBarry Smith     n8 = size - (m * n - rank) - (m - 1) + m;
67347c6ae99SBarry Smith   }
67447c6ae99SBarry Smith 
675ce00eea3SSatish Balay   if ((xe == M) && (ze == P)) { /* Assume an edge, not corner */
67647c6ae99SBarry Smith     n20 = m * n - (size - rank) - (m - 1) - m;
67747c6ae99SBarry Smith     n23 = m * n - (size - rank) - (m - 1);
67847c6ae99SBarry Smith     n26 = m * n - (size - rank) - (m - 1) + m;
67947c6ae99SBarry Smith   }
68047c6ae99SBarry Smith 
681ce00eea3SSatish Balay   if ((xe == M) && (ys == 0)) { /* Assume an edge, not corner */
68247c6ae99SBarry Smith     n2  = rank + m * (n - 1) - (m - 1) - m * n;
68347c6ae99SBarry Smith     n11 = rank + m * (n - 1) - (m - 1);
68447c6ae99SBarry Smith     n20 = rank + m * (n - 1) - (m - 1) + m * n;
68547c6ae99SBarry Smith   }
68647c6ae99SBarry Smith 
687ce00eea3SSatish Balay   if ((xe == M) && (ye == N)) { /* Assume an edge, not corner */
68847c6ae99SBarry Smith     n8  = rank - m * n + 1 - m * n;
68947c6ae99SBarry Smith     n17 = rank - m * n + 1;
69047c6ae99SBarry Smith     n26 = rank - m * n + 1 + m * n;
69147c6ae99SBarry Smith   }
69247c6ae99SBarry Smith 
69347c6ae99SBarry Smith   if ((ys == 0) && (zs == 0)) { /* Assume an edge, not corner */
69447c6ae99SBarry Smith     n0 = size - m + rank - 1;
69547c6ae99SBarry Smith     n1 = size - m + rank;
69647c6ae99SBarry Smith     n2 = size - m + rank + 1;
69747c6ae99SBarry Smith   }
69847c6ae99SBarry Smith 
69947c6ae99SBarry Smith   if ((ys == 0) && (ze == P)) { /* Assume an edge, not corner */
70047c6ae99SBarry Smith     n18 = m * n - (size - rank) + m * (n - 1) - 1;
70147c6ae99SBarry Smith     n19 = m * n - (size - rank) + m * (n - 1);
70247c6ae99SBarry Smith     n20 = m * n - (size - rank) + m * (n - 1) + 1;
70347c6ae99SBarry Smith   }
70447c6ae99SBarry Smith 
70547c6ae99SBarry Smith   if ((ye == N) && (zs == 0)) { /* Assume an edge, not corner */
70647c6ae99SBarry Smith     n6 = size - (m * n - rank) - m * (n - 1) - 1;
70747c6ae99SBarry Smith     n7 = size - (m * n - rank) - m * (n - 1);
70847c6ae99SBarry Smith     n8 = size - (m * n - rank) - m * (n - 1) + 1;
70947c6ae99SBarry Smith   }
71047c6ae99SBarry Smith 
71147c6ae99SBarry Smith   if ((ye == N) && (ze == P)) { /* Assume an edge, not corner */
71247c6ae99SBarry Smith     n24 = rank - (size - m) - 1;
71347c6ae99SBarry Smith     n25 = rank - (size - m);
71447c6ae99SBarry Smith     n26 = rank - (size - m) + 1;
71547c6ae99SBarry Smith   }
71647c6ae99SBarry Smith 
71747c6ae99SBarry Smith   /* Check for Corners */
7188865f1eaSKarl Rupp   if ((xs == 0) && (ys == 0) && (zs == 0)) n0 = size - 1;
7198865f1eaSKarl Rupp   if ((xs == 0) && (ys == 0) && (ze == P)) n18 = m * n - 1;
7208865f1eaSKarl Rupp   if ((xs == 0) && (ye == N) && (zs == 0)) n6 = (size - 1) - m * (n - 1);
7218865f1eaSKarl Rupp   if ((xs == 0) && (ye == N) && (ze == P)) n24 = m - 1;
7228865f1eaSKarl Rupp   if ((xe == M) && (ys == 0) && (zs == 0)) n2 = size - m;
7238865f1eaSKarl Rupp   if ((xe == M) && (ys == 0) && (ze == P)) n20 = m * n - m;
7248865f1eaSKarl Rupp   if ((xe == M) && (ye == N) && (zs == 0)) n8 = size - m * n;
7258865f1eaSKarl Rupp   if ((xe == M) && (ye == N) && (ze == P)) n26 = 0;
72647c6ae99SBarry Smith 
72747c6ae99SBarry Smith   /* Check for when not X,Y, and Z Periodic */
72847c6ae99SBarry Smith 
72947c6ae99SBarry Smith   /* If not X periodic */
730bff4a2f0SMatthew G. Knepley   if (bx != DM_BOUNDARY_PERIODIC) {
7318865f1eaSKarl Rupp     if (xs == 0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2;
7328865f1eaSKarl Rupp     if (xe == M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2;
73347c6ae99SBarry Smith   }
73447c6ae99SBarry Smith 
73547c6ae99SBarry Smith   /* If not Y periodic */
736bff4a2f0SMatthew G. Knepley   if (by != DM_BOUNDARY_PERIODIC) {
7378865f1eaSKarl Rupp     if (ys == 0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2;
7388865f1eaSKarl Rupp     if (ye == N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2;
73947c6ae99SBarry Smith   }
74047c6ae99SBarry Smith 
74147c6ae99SBarry Smith   /* If not Z periodic */
742bff4a2f0SMatthew G. Knepley   if (bz != DM_BOUNDARY_PERIODIC) {
7438865f1eaSKarl Rupp     if (zs == 0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2;
7448865f1eaSKarl Rupp     if (ze == P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;
74547c6ae99SBarry Smith   }
74647c6ae99SBarry Smith 
7479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(27, &dd->neighbors));
7488865f1eaSKarl Rupp 
74947c6ae99SBarry Smith   dd->neighbors[0]  = n0;
75047c6ae99SBarry Smith   dd->neighbors[1]  = n1;
75147c6ae99SBarry Smith   dd->neighbors[2]  = n2;
75247c6ae99SBarry Smith   dd->neighbors[3]  = n3;
75347c6ae99SBarry Smith   dd->neighbors[4]  = n4;
75447c6ae99SBarry Smith   dd->neighbors[5]  = n5;
75547c6ae99SBarry Smith   dd->neighbors[6]  = n6;
75647c6ae99SBarry Smith   dd->neighbors[7]  = n7;
75747c6ae99SBarry Smith   dd->neighbors[8]  = n8;
75847c6ae99SBarry Smith   dd->neighbors[9]  = n9;
75947c6ae99SBarry Smith   dd->neighbors[10] = n10;
76047c6ae99SBarry Smith   dd->neighbors[11] = n11;
76147c6ae99SBarry Smith   dd->neighbors[12] = n12;
76247c6ae99SBarry Smith   dd->neighbors[13] = rank;
76347c6ae99SBarry Smith   dd->neighbors[14] = n14;
76447c6ae99SBarry Smith   dd->neighbors[15] = n15;
76547c6ae99SBarry Smith   dd->neighbors[16] = n16;
76647c6ae99SBarry Smith   dd->neighbors[17] = n17;
76747c6ae99SBarry Smith   dd->neighbors[18] = n18;
76847c6ae99SBarry Smith   dd->neighbors[19] = n19;
76947c6ae99SBarry Smith   dd->neighbors[20] = n20;
77047c6ae99SBarry Smith   dd->neighbors[21] = n21;
77147c6ae99SBarry Smith   dd->neighbors[22] = n22;
77247c6ae99SBarry Smith   dd->neighbors[23] = n23;
77347c6ae99SBarry Smith   dd->neighbors[24] = n24;
77447c6ae99SBarry Smith   dd->neighbors[25] = n25;
77547c6ae99SBarry Smith   dd->neighbors[26] = n26;
77647c6ae99SBarry Smith 
77747c6ae99SBarry Smith   /* If star stencil then delete the corner neighbors */
778d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
77947c6ae99SBarry Smith     /* save information about corner neighbors */
7809371c9d4SSatish Balay     sn0  = n0;
7819371c9d4SSatish Balay     sn1  = n1;
7829371c9d4SSatish Balay     sn2  = n2;
7839371c9d4SSatish Balay     sn3  = n3;
7849371c9d4SSatish Balay     sn5  = n5;
7859371c9d4SSatish Balay     sn6  = n6;
7869371c9d4SSatish Balay     sn7  = n7;
7879371c9d4SSatish Balay     sn8  = n8;
7889371c9d4SSatish Balay     sn9  = n9;
7899371c9d4SSatish Balay     sn11 = n11;
7909371c9d4SSatish Balay     sn15 = n15;
7919371c9d4SSatish Balay     sn17 = n17;
7929371c9d4SSatish Balay     sn18 = n18;
7939371c9d4SSatish Balay     sn19 = n19;
7949371c9d4SSatish Balay     sn20 = n20;
7959371c9d4SSatish Balay     sn21 = n21;
7969371c9d4SSatish Balay     sn23 = n23;
7979371c9d4SSatish Balay     sn24 = n24;
7989371c9d4SSatish Balay     sn25 = n25;
79947c6ae99SBarry Smith     sn26 = n26;
8008865f1eaSKarl Rupp     n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
80147c6ae99SBarry Smith   }
80247c6ae99SBarry Smith 
8039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((Xe - Xs) * (Ye - Ys) * (Ze - Zs), &idx));
80447c6ae99SBarry Smith 
80547c6ae99SBarry Smith   nn = 0;
80647c6ae99SBarry Smith   /* Bottom Level */
80747c6ae99SBarry Smith   for (k = 0; k < s_z; k++) {
80847c6ae99SBarry Smith     for (i = 1; i <= s_y; i++) {
80947c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
810ce00eea3SSatish Balay         x_t = lx[n0 % m];
81147c6ae99SBarry Smith         y_t = ly[(n0 % (m * n)) / m];
81247c6ae99SBarry Smith         z_t = lz[n0 / (m * n)];
81347c6ae99SBarry Smith         s_t = bases[n0] + x_t * y_t * z_t - (s_y - i) * x_t - s_x - (s_z - k - 1) * x_t * y_t;
8148865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n0] + x_t * y_t * z_t - (s_y - i) * x_t - s_x; /* 2D case */
8158865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
81647c6ae99SBarry Smith       }
81747c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
81847c6ae99SBarry Smith         x_t = x;
81947c6ae99SBarry Smith         y_t = ly[(n1 % (m * n)) / m];
82047c6ae99SBarry Smith         z_t = lz[n1 / (m * n)];
82147c6ae99SBarry Smith         s_t = bases[n1] + x_t * y_t * z_t - (s_y + 1 - i) * x_t - (s_z - k - 1) * x_t * y_t;
8228865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n1] + x_t * y_t * z_t - (s_y + 1 - i) * x_t; /* 2D case */
8238865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
82447c6ae99SBarry Smith       }
82547c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
826ce00eea3SSatish Balay         x_t = lx[n2 % m];
82747c6ae99SBarry Smith         y_t = ly[(n2 % (m * n)) / m];
82847c6ae99SBarry Smith         z_t = lz[n2 / (m * n)];
82947c6ae99SBarry Smith         s_t = bases[n2] + x_t * y_t * z_t - (s_y + 1 - i) * x_t - (s_z - k - 1) * x_t * y_t;
8308865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n2] + x_t * y_t * z_t - (s_y + 1 - i) * x_t; /* 2D case */
8318865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
83247c6ae99SBarry Smith       }
83347c6ae99SBarry Smith     }
83447c6ae99SBarry Smith 
83547c6ae99SBarry Smith     for (i = 0; i < y; i++) {
83647c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
837ce00eea3SSatish Balay         x_t = lx[n3 % m];
83847c6ae99SBarry Smith         y_t = y;
83947c6ae99SBarry Smith         z_t = lz[n3 / (m * n)];
84047c6ae99SBarry Smith         s_t = bases[n3] + (i + 1) * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
8418865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n3] + (i + 1) * x_t - s_x + x_t * y_t * z_t - x_t * y_t; /* 2D case */
8428865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
84347c6ae99SBarry Smith       }
84447c6ae99SBarry Smith 
84547c6ae99SBarry Smith       if (n4 >= 0) { /* middle */
84647c6ae99SBarry Smith         x_t = x;
84747c6ae99SBarry Smith         y_t = y;
84847c6ae99SBarry Smith         z_t = lz[n4 / (m * n)];
84947c6ae99SBarry Smith         s_t = bases[n4] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
8508865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n4] + i * x_t + x_t * y_t * z_t - x_t * y_t; /* 2D case */
8518865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
852bff4a2f0SMatthew G. Knepley       } else if (bz == DM_BOUNDARY_MIRROR) {
853a6fd6b0aSBarry Smith         for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + x * i + (s_z - k - 1) * x * y;
85447c6ae99SBarry Smith       }
85547c6ae99SBarry Smith 
85647c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
857ce00eea3SSatish Balay         x_t = lx[n5 % m];
85847c6ae99SBarry Smith         y_t = y;
85947c6ae99SBarry Smith         z_t = lz[n5 / (m * n)];
86047c6ae99SBarry Smith         s_t = bases[n5] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
8618865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n5] + i * x_t + x_t * y_t * z_t - x_t * y_t; /* 2D case */
8628865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
86347c6ae99SBarry Smith       }
86447c6ae99SBarry Smith     }
86547c6ae99SBarry Smith 
86647c6ae99SBarry Smith     for (i = 1; i <= s_y; i++) {
86747c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
868ce00eea3SSatish Balay         x_t = lx[n6 % m];
86947c6ae99SBarry Smith         y_t = ly[(n6 % (m * n)) / m];
87047c6ae99SBarry Smith         z_t = lz[n6 / (m * n)];
87147c6ae99SBarry Smith         s_t = bases[n6] + i * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
8728865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n6] + i * x_t - s_x + x_t * y_t * z_t - x_t * y_t; /* 2D case */
8738865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
87447c6ae99SBarry Smith       }
87547c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
87647c6ae99SBarry Smith         x_t = x;
87747c6ae99SBarry Smith         y_t = ly[(n7 % (m * n)) / m];
87847c6ae99SBarry Smith         z_t = lz[n7 / (m * n)];
87947c6ae99SBarry Smith         s_t = bases[n7] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
8808865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n7] + (i - 1) * x_t + x_t * y_t * z_t - x_t * y_t; /* 2D case */
8818865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
88247c6ae99SBarry Smith       }
88347c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
884ce00eea3SSatish Balay         x_t = lx[n8 % m];
88547c6ae99SBarry Smith         y_t = ly[(n8 % (m * n)) / m];
88647c6ae99SBarry Smith         z_t = lz[n8 / (m * n)];
88747c6ae99SBarry Smith         s_t = bases[n8] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
8888865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n8] + (i - 1) * x_t + x_t * y_t * z_t - x_t * y_t; /* 2D case */
8898865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
89047c6ae99SBarry Smith       }
89147c6ae99SBarry Smith     }
89247c6ae99SBarry Smith   }
89347c6ae99SBarry Smith 
89447c6ae99SBarry Smith   /* Middle Level */
89547c6ae99SBarry Smith   for (k = 0; k < z; k++) {
89647c6ae99SBarry Smith     for (i = 1; i <= s_y; i++) {
89747c6ae99SBarry Smith       if (n9 >= 0) { /* left below */
898ce00eea3SSatish Balay         x_t = lx[n9 % m];
89947c6ae99SBarry Smith         y_t = ly[(n9 % (m * n)) / m];
90047c6ae99SBarry Smith         /* z_t = z; */
90147c6ae99SBarry Smith         s_t = bases[n9] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t;
9028865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
90347c6ae99SBarry Smith       }
90447c6ae99SBarry Smith       if (n10 >= 0) { /* directly below */
90547c6ae99SBarry Smith         x_t = x;
90647c6ae99SBarry Smith         y_t = ly[(n10 % (m * n)) / m];
90747c6ae99SBarry Smith         /* z_t = z; */
90847c6ae99SBarry Smith         s_t = bases[n10] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
9098865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
910bff4a2f0SMatthew G. Knepley       } else if (by == DM_BOUNDARY_MIRROR) {
911a6fd6b0aSBarry Smith         for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (s_y - i) * x;
91247c6ae99SBarry Smith       }
91347c6ae99SBarry Smith       if (n11 >= 0) { /* right below */
914ce00eea3SSatish Balay         x_t = lx[n11 % m];
91547c6ae99SBarry Smith         y_t = ly[(n11 % (m * n)) / m];
91647c6ae99SBarry Smith         /* z_t = z; */
91747c6ae99SBarry Smith         s_t = bases[n11] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
9188865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
91947c6ae99SBarry Smith       }
92047c6ae99SBarry Smith     }
92147c6ae99SBarry Smith 
92247c6ae99SBarry Smith     for (i = 0; i < y; i++) {
92347c6ae99SBarry Smith       if (n12 >= 0) { /* directly left */
924ce00eea3SSatish Balay         x_t = lx[n12 % m];
92547c6ae99SBarry Smith         y_t = y;
92647c6ae99SBarry Smith         /* z_t = z; */
92747c6ae99SBarry Smith         s_t = bases[n12] + (i + 1) * x_t - s_x + k * x_t * y_t;
9288865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
929bff4a2f0SMatthew G. Knepley       } else if (bx == DM_BOUNDARY_MIRROR) {
930a6fd6b0aSBarry Smith         for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k * x * y + i * x;
93147c6ae99SBarry Smith       }
93247c6ae99SBarry Smith 
93347c6ae99SBarry Smith       /* Interior */
93447c6ae99SBarry Smith       s_t = bases[rank] + i * x + k * x * y;
9358865f1eaSKarl Rupp       for (j = 0; j < x; j++) idx[nn++] = s_t++;
93647c6ae99SBarry Smith 
93747c6ae99SBarry Smith       if (n14 >= 0) { /* directly right */
938ce00eea3SSatish Balay         x_t = lx[n14 % m];
93947c6ae99SBarry Smith         y_t = y;
94047c6ae99SBarry Smith         /* z_t = z; */
94147c6ae99SBarry Smith         s_t = bases[n14] + i * x_t + k * x_t * y_t;
9428865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
943bff4a2f0SMatthew G. Knepley       } else if (bx == DM_BOUNDARY_MIRROR) {
944a6fd6b0aSBarry Smith         for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k * x * y + i * x;
94547c6ae99SBarry Smith       }
94647c6ae99SBarry Smith     }
94747c6ae99SBarry Smith 
94847c6ae99SBarry Smith     for (i = 1; i <= s_y; i++) {
94947c6ae99SBarry Smith       if (n15 >= 0) { /* left above */
950ce00eea3SSatish Balay         x_t = lx[n15 % m];
95147c6ae99SBarry Smith         y_t = ly[(n15 % (m * n)) / m];
95247c6ae99SBarry Smith         /* z_t = z; */
95347c6ae99SBarry Smith         s_t = bases[n15] + i * x_t - s_x + k * x_t * y_t;
9548865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
95547c6ae99SBarry Smith       }
95647c6ae99SBarry Smith       if (n16 >= 0) { /* directly above */
95747c6ae99SBarry Smith         x_t = x;
95847c6ae99SBarry Smith         y_t = ly[(n16 % (m * n)) / m];
95947c6ae99SBarry Smith         /* z_t = z; */
96047c6ae99SBarry Smith         s_t = bases[n16] + (i - 1) * x_t + k * x_t * y_t;
9618865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
962bff4a2f0SMatthew G. Knepley       } else if (by == DM_BOUNDARY_MIRROR) {
963a6fd6b0aSBarry Smith         for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (y - i) * x;
96447c6ae99SBarry Smith       }
96547c6ae99SBarry Smith       if (n17 >= 0) { /* right above */
966ce00eea3SSatish Balay         x_t = lx[n17 % m];
96747c6ae99SBarry Smith         y_t = ly[(n17 % (m * n)) / m];
96847c6ae99SBarry Smith         /* z_t = z; */
96947c6ae99SBarry Smith         s_t = bases[n17] + (i - 1) * x_t + k * x_t * y_t;
9708865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
97147c6ae99SBarry Smith       }
97247c6ae99SBarry Smith     }
97347c6ae99SBarry Smith   }
97447c6ae99SBarry Smith 
97547c6ae99SBarry Smith   /* Upper Level */
97647c6ae99SBarry Smith   for (k = 0; k < s_z; k++) {
97747c6ae99SBarry Smith     for (i = 1; i <= s_y; i++) {
97847c6ae99SBarry Smith       if (n18 >= 0) { /* left below */
979ce00eea3SSatish Balay         x_t = lx[n18 % m];
98047c6ae99SBarry Smith         y_t = ly[(n18 % (m * n)) / m];
98147c6ae99SBarry Smith         /* z_t = lz[n18 / (m*n)]; */
98247c6ae99SBarry Smith         s_t = bases[n18] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t;
9838865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n18] - (s_y - i) * x_t - s_x + x_t * y_t; /* 2d case */
9848865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
98547c6ae99SBarry Smith       }
98647c6ae99SBarry Smith       if (n19 >= 0) { /* directly below */
98747c6ae99SBarry Smith         x_t = x;
98847c6ae99SBarry Smith         y_t = ly[(n19 % (m * n)) / m];
98947c6ae99SBarry Smith         /* z_t = lz[n19 / (m*n)]; */
99047c6ae99SBarry Smith         s_t = bases[n19] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
9918865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n19] - (s_y + 1 - i) * x_t + x_t * y_t; /* 2d case */
9928865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
99347c6ae99SBarry Smith       }
99447c6ae99SBarry Smith       if (n20 >= 0) { /* right below */
995ce00eea3SSatish Balay         x_t = lx[n20 % m];
99647c6ae99SBarry Smith         y_t = ly[(n20 % (m * n)) / m];
99747c6ae99SBarry Smith         /* z_t = lz[n20 / (m*n)]; */
99847c6ae99SBarry Smith         s_t = bases[n20] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
9998865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n20] - (s_y + 1 - i) * x_t + x_t * y_t; /* 2d case */
10008865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
100147c6ae99SBarry Smith       }
100247c6ae99SBarry Smith     }
100347c6ae99SBarry Smith 
100447c6ae99SBarry Smith     for (i = 0; i < y; i++) {
100547c6ae99SBarry Smith       if (n21 >= 0) { /* directly left */
1006ce00eea3SSatish Balay         x_t = lx[n21 % m];
100747c6ae99SBarry Smith         y_t = y;
100847c6ae99SBarry Smith         /* z_t = lz[n21 / (m*n)]; */
100947c6ae99SBarry Smith         s_t = bases[n21] + (i + 1) * x_t - s_x + k * x_t * y_t;
10108865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n21] + (i + 1) * x_t - s_x; /* 2d case */
10118865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
101247c6ae99SBarry Smith       }
101347c6ae99SBarry Smith 
101447c6ae99SBarry Smith       if (n22 >= 0) { /* middle */
101547c6ae99SBarry Smith         x_t = x;
101647c6ae99SBarry Smith         y_t = y;
101747c6ae99SBarry Smith         /* z_t = lz[n22 / (m*n)]; */
101847c6ae99SBarry Smith         s_t = bases[n22] + i * x_t + k * x_t * y_t;
10198865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n22] + i * x_t; /* 2d case */
10208865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1021bff4a2f0SMatthew G. Knepley       } else if (bz == DM_BOUNDARY_MIRROR) {
1022a6fd6b0aSBarry Smith         for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + (z - k - 1) * x * y + i * x;
102347c6ae99SBarry Smith       }
102447c6ae99SBarry Smith 
102547c6ae99SBarry Smith       if (n23 >= 0) { /* directly right */
1026ce00eea3SSatish Balay         x_t = lx[n23 % m];
102747c6ae99SBarry Smith         y_t = y;
102847c6ae99SBarry Smith         /* z_t = lz[n23 / (m*n)]; */
102947c6ae99SBarry Smith         s_t = bases[n23] + i * x_t + k * x_t * y_t;
10308865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n23] + i * x_t; /* 2d case */
10318865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
103247c6ae99SBarry Smith       }
103347c6ae99SBarry Smith     }
103447c6ae99SBarry Smith 
103547c6ae99SBarry Smith     for (i = 1; i <= s_y; i++) {
103647c6ae99SBarry Smith       if (n24 >= 0) { /* left above */
1037ce00eea3SSatish Balay         x_t = lx[n24 % m];
103847c6ae99SBarry Smith         y_t = ly[(n24 % (m * n)) / m];
103947c6ae99SBarry Smith         /* z_t = lz[n24 / (m*n)]; */
104047c6ae99SBarry Smith         s_t = bases[n24] + i * x_t - s_x + k * x_t * y_t;
10418865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n24] + i * x_t - s_x; /* 2d case */
10428865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
104347c6ae99SBarry Smith       }
104447c6ae99SBarry Smith       if (n25 >= 0) { /* directly above */
104547c6ae99SBarry Smith         x_t = x;
104647c6ae99SBarry Smith         y_t = ly[(n25 % (m * n)) / m];
104747c6ae99SBarry Smith         /* z_t = lz[n25 / (m*n)]; */
104847c6ae99SBarry Smith         s_t = bases[n25] + (i - 1) * x_t + k * x_t * y_t;
10498865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n25] + (i - 1) * x_t; /* 2d case */
10508865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
105147c6ae99SBarry Smith       }
105247c6ae99SBarry Smith       if (n26 >= 0) { /* right above */
1053ce00eea3SSatish Balay         x_t = lx[n26 % m];
105447c6ae99SBarry Smith         y_t = ly[(n26 % (m * n)) / m];
105547c6ae99SBarry Smith         /* z_t = lz[n26 / (m*n)]; */
105647c6ae99SBarry Smith         s_t = bases[n26] + (i - 1) * x_t + k * x_t * y_t;
10578865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n26] + (i - 1) * x_t; /* 2d case */
10588865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
105947c6ae99SBarry Smith       }
106047c6ae99SBarry Smith     }
106147c6ae99SBarry Smith   }
1062ce00eea3SSatish Balay 
10639566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(comm, dof, nn, idx, PETSC_USE_POINTER, &from));
10649566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(global, from, local, to, &gtol));
10659566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
10669566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
106747c6ae99SBarry Smith 
1068d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
10699371c9d4SSatish Balay     n0  = sn0;
10709371c9d4SSatish Balay     n1  = sn1;
10719371c9d4SSatish Balay     n2  = sn2;
10729371c9d4SSatish Balay     n3  = sn3;
10739371c9d4SSatish Balay     n5  = sn5;
10749371c9d4SSatish Balay     n6  = sn6;
10759371c9d4SSatish Balay     n7  = sn7;
10769371c9d4SSatish Balay     n8  = sn8;
10779371c9d4SSatish Balay     n9  = sn9;
10789371c9d4SSatish Balay     n11 = sn11;
10799371c9d4SSatish Balay     n15 = sn15;
10809371c9d4SSatish Balay     n17 = sn17;
10819371c9d4SSatish Balay     n18 = sn18;
10829371c9d4SSatish Balay     n19 = sn19;
10839371c9d4SSatish Balay     n20 = sn20;
10849371c9d4SSatish Balay     n21 = sn21;
10859371c9d4SSatish Balay     n23 = sn23;
10869371c9d4SSatish Balay     n24 = sn24;
10879371c9d4SSatish Balay     n25 = sn25;
108847c6ae99SBarry Smith     n26 = sn26;
1089ce00eea3SSatish Balay   }
109047c6ae99SBarry Smith 
1091288e7d53SBarry Smith   if (((stencil_type == DMDA_STENCIL_STAR) || (bx != DM_BOUNDARY_PERIODIC && bx) || (by != DM_BOUNDARY_PERIODIC && by) || (bz != DM_BOUNDARY_PERIODIC && bz))) {
1092ce00eea3SSatish Balay     /*
1093ce00eea3SSatish Balay         Recompute the local to global mappings, this time keeping the
1094ce00eea3SSatish Balay       information about the cross corner processor numbers.
1095ce00eea3SSatish Balay     */
109647c6ae99SBarry Smith     nn = 0;
109747c6ae99SBarry Smith     /* Bottom Level */
109847c6ae99SBarry Smith     for (k = 0; k < s_z; k++) {
109947c6ae99SBarry Smith       for (i = 1; i <= s_y; i++) {
110047c6ae99SBarry Smith         if (n0 >= 0) { /* left below */
1101ce00eea3SSatish Balay           x_t = lx[n0 % m];
110247c6ae99SBarry Smith           y_t = ly[(n0 % (m * n)) / m];
110347c6ae99SBarry Smith           z_t = lz[n0 / (m * n)];
110447c6ae99SBarry Smith           s_t = bases[n0] + x_t * y_t * z_t - (s_y - i) * x_t - s_x - (s_z - k - 1) * x_t * y_t;
11058865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1106ce00eea3SSatish Balay         } else if (Xs - xs < 0 && Ys - ys < 0 && Zs - zs < 0) {
11078865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
110847c6ae99SBarry Smith         }
110947c6ae99SBarry Smith         if (n1 >= 0) { /* directly below */
111047c6ae99SBarry Smith           x_t = x;
111147c6ae99SBarry Smith           y_t = ly[(n1 % (m * n)) / m];
111247c6ae99SBarry Smith           z_t = lz[n1 / (m * n)];
111347c6ae99SBarry Smith           s_t = bases[n1] + x_t * y_t * z_t - (s_y + 1 - i) * x_t - (s_z - k - 1) * x_t * y_t;
11148865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1115ce00eea3SSatish Balay         } else if (Ys - ys < 0 && Zs - zs < 0) {
11168865f1eaSKarl Rupp           for (j = 0; j < x; j++) idx[nn++] = -1;
111747c6ae99SBarry Smith         }
111847c6ae99SBarry Smith         if (n2 >= 0) { /* right below */
1119ce00eea3SSatish Balay           x_t = lx[n2 % m];
112047c6ae99SBarry Smith           y_t = ly[(n2 % (m * n)) / m];
112147c6ae99SBarry Smith           z_t = lz[n2 / (m * n)];
112247c6ae99SBarry Smith           s_t = bases[n2] + x_t * y_t * z_t - (s_y + 1 - i) * x_t - (s_z - k - 1) * x_t * y_t;
11238865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1124ce00eea3SSatish Balay         } else if (xe - Xe < 0 && Ys - ys < 0 && Zs - zs < 0) {
11258865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
112647c6ae99SBarry Smith         }
112747c6ae99SBarry Smith       }
112847c6ae99SBarry Smith 
112947c6ae99SBarry Smith       for (i = 0; i < y; i++) {
113047c6ae99SBarry Smith         if (n3 >= 0) { /* directly left */
1131ce00eea3SSatish Balay           x_t = lx[n3 % m];
113247c6ae99SBarry Smith           y_t = y;
113347c6ae99SBarry Smith           z_t = lz[n3 / (m * n)];
113447c6ae99SBarry Smith           s_t = bases[n3] + (i + 1) * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
11358865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1136ce00eea3SSatish Balay         } else if (Xs - xs < 0 && Zs - zs < 0) {
11378865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
113847c6ae99SBarry Smith         }
113947c6ae99SBarry Smith 
114047c6ae99SBarry Smith         if (n4 >= 0) { /* middle */
114147c6ae99SBarry Smith           x_t = x;
114247c6ae99SBarry Smith           y_t = y;
114347c6ae99SBarry Smith           z_t = lz[n4 / (m * n)];
114447c6ae99SBarry Smith           s_t = bases[n4] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
11458865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1146ce00eea3SSatish Balay         } else if (Zs - zs < 0) {
1147bff4a2f0SMatthew G. Knepley           if (bz == DM_BOUNDARY_MIRROR) {
1148a6fd6b0aSBarry Smith             for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + x * i + (s_z - k - 1) * x * y;
1149c2859e5eSBarry Smith           } else {
11508865f1eaSKarl Rupp             for (j = 0; j < x; j++) idx[nn++] = -1;
115147c6ae99SBarry Smith           }
1152c2859e5eSBarry Smith         }
115347c6ae99SBarry Smith 
115447c6ae99SBarry Smith         if (n5 >= 0) { /* directly right */
1155ce00eea3SSatish Balay           x_t = lx[n5 % m];
115647c6ae99SBarry Smith           y_t = y;
115747c6ae99SBarry Smith           z_t = lz[n5 / (m * n)];
115847c6ae99SBarry Smith           s_t = bases[n5] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
11598865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1160ce00eea3SSatish Balay         } else if (xe - Xe < 0 && Zs - zs < 0) {
11618865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
116247c6ae99SBarry Smith         }
116347c6ae99SBarry Smith       }
116447c6ae99SBarry Smith 
116547c6ae99SBarry Smith       for (i = 1; i <= s_y; i++) {
116647c6ae99SBarry Smith         if (n6 >= 0) { /* left above */
1167ce00eea3SSatish Balay           x_t = lx[n6 % m];
116847c6ae99SBarry Smith           y_t = ly[(n6 % (m * n)) / m];
116947c6ae99SBarry Smith           z_t = lz[n6 / (m * n)];
117047c6ae99SBarry Smith           s_t = bases[n6] + i * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
11718865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1172ce00eea3SSatish Balay         } else if (Xs - xs < 0 && ye - Ye < 0 && Zs - zs < 0) {
11738865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
117447c6ae99SBarry Smith         }
117547c6ae99SBarry Smith         if (n7 >= 0) { /* directly above */
117647c6ae99SBarry Smith           x_t = x;
117747c6ae99SBarry Smith           y_t = ly[(n7 % (m * n)) / m];
117847c6ae99SBarry Smith           z_t = lz[n7 / (m * n)];
117947c6ae99SBarry Smith           s_t = bases[n7] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
11808865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1181ce00eea3SSatish Balay         } else if (ye - Ye < 0 && Zs - zs < 0) {
11828865f1eaSKarl Rupp           for (j = 0; j < x; j++) idx[nn++] = -1;
118347c6ae99SBarry Smith         }
118447c6ae99SBarry Smith         if (n8 >= 0) { /* right above */
1185ce00eea3SSatish Balay           x_t = lx[n8 % m];
118647c6ae99SBarry Smith           y_t = ly[(n8 % (m * n)) / m];
118747c6ae99SBarry Smith           z_t = lz[n8 / (m * n)];
118847c6ae99SBarry Smith           s_t = bases[n8] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
11898865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1190ce00eea3SSatish Balay         } else if (xe - Xe < 0 && ye - Ye < 0 && Zs - zs < 0) {
11918865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
119247c6ae99SBarry Smith         }
119347c6ae99SBarry Smith       }
119447c6ae99SBarry Smith     }
119547c6ae99SBarry Smith 
119647c6ae99SBarry Smith     /* Middle Level */
119747c6ae99SBarry Smith     for (k = 0; k < z; k++) {
119847c6ae99SBarry Smith       for (i = 1; i <= s_y; i++) {
119947c6ae99SBarry Smith         if (n9 >= 0) { /* left below */
1200ce00eea3SSatish Balay           x_t = lx[n9 % m];
120147c6ae99SBarry Smith           y_t = ly[(n9 % (m * n)) / m];
120247c6ae99SBarry Smith           /* z_t = z; */
120347c6ae99SBarry Smith           s_t = bases[n9] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t;
12048865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1205ce00eea3SSatish Balay         } else if (Xs - xs < 0 && Ys - ys < 0) {
12068865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
120747c6ae99SBarry Smith         }
120847c6ae99SBarry Smith         if (n10 >= 0) { /* directly below */
120947c6ae99SBarry Smith           x_t = x;
121047c6ae99SBarry Smith           y_t = ly[(n10 % (m * n)) / m];
121147c6ae99SBarry Smith           /* z_t = z; */
121247c6ae99SBarry Smith           s_t = bases[n10] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
12138865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1214ce00eea3SSatish Balay         } else if (Ys - ys < 0) {
1215bff4a2f0SMatthew G. Knepley           if (by == DM_BOUNDARY_MIRROR) {
1216feb237baSPierre Jolivet             for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (s_y - i) * x;
1217c2859e5eSBarry Smith           } else {
12188865f1eaSKarl Rupp             for (j = 0; j < x; j++) idx[nn++] = -1;
1219c2859e5eSBarry Smith           }
122047c6ae99SBarry Smith         }
122147c6ae99SBarry Smith         if (n11 >= 0) { /* right below */
1222ce00eea3SSatish Balay           x_t = lx[n11 % m];
122347c6ae99SBarry Smith           y_t = ly[(n11 % (m * n)) / m];
122447c6ae99SBarry Smith           /* z_t = z; */
122547c6ae99SBarry Smith           s_t = bases[n11] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
12268865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1227ce00eea3SSatish Balay         } else if (xe - Xe < 0 && Ys - ys < 0) {
12288865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
122947c6ae99SBarry Smith         }
123047c6ae99SBarry Smith       }
123147c6ae99SBarry Smith 
123247c6ae99SBarry Smith       for (i = 0; i < y; i++) {
123347c6ae99SBarry Smith         if (n12 >= 0) { /* directly left */
1234ce00eea3SSatish Balay           x_t = lx[n12 % m];
123547c6ae99SBarry Smith           y_t = y;
123647c6ae99SBarry Smith           /* z_t = z; */
123747c6ae99SBarry Smith           s_t = bases[n12] + (i + 1) * x_t - s_x + k * x_t * y_t;
12388865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1239ce00eea3SSatish Balay         } else if (Xs - xs < 0) {
1240bff4a2f0SMatthew G. Knepley           if (bx == DM_BOUNDARY_MIRROR) {
1241a6fd6b0aSBarry Smith             for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k * x * y + i * x;
1242c2859e5eSBarry Smith           } else {
12438865f1eaSKarl Rupp             for (j = 0; j < s_x; j++) idx[nn++] = -1;
124447c6ae99SBarry Smith           }
1245c2859e5eSBarry Smith         }
124647c6ae99SBarry Smith 
124747c6ae99SBarry Smith         /* Interior */
124847c6ae99SBarry Smith         s_t = bases[rank] + i * x + k * x * y;
12498865f1eaSKarl Rupp         for (j = 0; j < x; j++) idx[nn++] = s_t++;
125047c6ae99SBarry Smith 
125147c6ae99SBarry Smith         if (n14 >= 0) { /* directly right */
1252ce00eea3SSatish Balay           x_t = lx[n14 % m];
125347c6ae99SBarry Smith           y_t = y;
125447c6ae99SBarry Smith           /* z_t = z; */
125547c6ae99SBarry Smith           s_t = bases[n14] + i * x_t + k * x_t * y_t;
12568865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1257ce00eea3SSatish Balay         } else if (xe - Xe < 0) {
1258bff4a2f0SMatthew G. Knepley           if (bx == DM_BOUNDARY_MIRROR) {
1259a6fd6b0aSBarry Smith             for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k * x * y + i * x;
1260c2859e5eSBarry Smith           } else {
12618865f1eaSKarl Rupp             for (j = 0; j < s_x; j++) idx[nn++] = -1;
126247c6ae99SBarry Smith           }
126347c6ae99SBarry Smith         }
1264c2859e5eSBarry Smith       }
126547c6ae99SBarry Smith 
126647c6ae99SBarry Smith       for (i = 1; i <= s_y; i++) {
126747c6ae99SBarry Smith         if (n15 >= 0) { /* left above */
1268ce00eea3SSatish Balay           x_t = lx[n15 % m];
126947c6ae99SBarry Smith           y_t = ly[(n15 % (m * n)) / m];
127047c6ae99SBarry Smith           /* z_t = z; */
127147c6ae99SBarry Smith           s_t = bases[n15] + i * x_t - s_x + k * x_t * y_t;
12728865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1273ce00eea3SSatish Balay         } else if (Xs - xs < 0 && ye - Ye < 0) {
12748865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
127547c6ae99SBarry Smith         }
127647c6ae99SBarry Smith         if (n16 >= 0) { /* directly above */
127747c6ae99SBarry Smith           x_t = x;
127847c6ae99SBarry Smith           y_t = ly[(n16 % (m * n)) / m];
127947c6ae99SBarry Smith           /* z_t = z; */
128047c6ae99SBarry Smith           s_t = bases[n16] + (i - 1) * x_t + k * x_t * y_t;
12818865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1282ce00eea3SSatish Balay         } else if (ye - Ye < 0) {
1283bff4a2f0SMatthew G. Knepley           if (by == DM_BOUNDARY_MIRROR) {
1284a6fd6b0aSBarry Smith             for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (y - i) * x;
1285c2859e5eSBarry Smith           } else {
12868865f1eaSKarl Rupp             for (j = 0; j < x; j++) idx[nn++] = -1;
128747c6ae99SBarry Smith           }
1288c2859e5eSBarry Smith         }
128947c6ae99SBarry Smith         if (n17 >= 0) { /* right above */
1290ce00eea3SSatish Balay           x_t = lx[n17 % m];
129147c6ae99SBarry Smith           y_t = ly[(n17 % (m * n)) / m];
129247c6ae99SBarry Smith           /* z_t = z; */
129347c6ae99SBarry Smith           s_t = bases[n17] + (i - 1) * x_t + k * x_t * y_t;
12948865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1295ce00eea3SSatish Balay         } else if (xe - Xe < 0 && ye - Ye < 0) {
12968865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
129747c6ae99SBarry Smith         }
129847c6ae99SBarry Smith       }
129947c6ae99SBarry Smith     }
130047c6ae99SBarry Smith 
130147c6ae99SBarry Smith     /* Upper Level */
130247c6ae99SBarry Smith     for (k = 0; k < s_z; k++) {
130347c6ae99SBarry Smith       for (i = 1; i <= s_y; i++) {
130447c6ae99SBarry Smith         if (n18 >= 0) { /* left below */
1305ce00eea3SSatish Balay           x_t = lx[n18 % m];
130647c6ae99SBarry Smith           y_t = ly[(n18 % (m * n)) / m];
130747c6ae99SBarry Smith           /* z_t = lz[n18 / (m*n)]; */
130847c6ae99SBarry Smith           s_t = bases[n18] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t;
13098865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1310ce00eea3SSatish Balay         } else if (Xs - xs < 0 && Ys - ys < 0 && ze - Ze < 0) {
13118865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
131247c6ae99SBarry Smith         }
131347c6ae99SBarry Smith         if (n19 >= 0) { /* directly below */
131447c6ae99SBarry Smith           x_t = x;
131547c6ae99SBarry Smith           y_t = ly[(n19 % (m * n)) / m];
131647c6ae99SBarry Smith           /* z_t = lz[n19 / (m*n)]; */
131747c6ae99SBarry Smith           s_t = bases[n19] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
13188865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1319ce00eea3SSatish Balay         } else if (Ys - ys < 0 && ze - Ze < 0) {
13208865f1eaSKarl Rupp           for (j = 0; j < x; j++) idx[nn++] = -1;
132147c6ae99SBarry Smith         }
132247c6ae99SBarry Smith         if (n20 >= 0) { /* right below */
1323ce00eea3SSatish Balay           x_t = lx[n20 % m];
132447c6ae99SBarry Smith           y_t = ly[(n20 % (m * n)) / m];
132547c6ae99SBarry Smith           /* z_t = lz[n20 / (m*n)]; */
132647c6ae99SBarry Smith           s_t = bases[n20] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
13278865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1328ce00eea3SSatish Balay         } else if (xe - Xe < 0 && Ys - ys < 0 && ze - Ze < 0) {
13298865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
133047c6ae99SBarry Smith         }
133147c6ae99SBarry Smith       }
133247c6ae99SBarry Smith 
133347c6ae99SBarry Smith       for (i = 0; i < y; i++) {
133447c6ae99SBarry Smith         if (n21 >= 0) { /* directly left */
1335ce00eea3SSatish Balay           x_t = lx[n21 % m];
133647c6ae99SBarry Smith           y_t = y;
133747c6ae99SBarry Smith           /* z_t = lz[n21 / (m*n)]; */
133847c6ae99SBarry Smith           s_t = bases[n21] + (i + 1) * x_t - s_x + k * x_t * y_t;
13398865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1340ce00eea3SSatish Balay         } else if (Xs - xs < 0 && ze - Ze < 0) {
13418865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
134247c6ae99SBarry Smith         }
134347c6ae99SBarry Smith 
134447c6ae99SBarry Smith         if (n22 >= 0) { /* middle */
134547c6ae99SBarry Smith           x_t = x;
134647c6ae99SBarry Smith           y_t = y;
134747c6ae99SBarry Smith           /* z_t = lz[n22 / (m*n)]; */
134847c6ae99SBarry Smith           s_t = bases[n22] + i * x_t + k * x_t * y_t;
13498865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1350ce00eea3SSatish Balay         } else if (ze - Ze < 0) {
1351bff4a2f0SMatthew G. Knepley           if (bz == DM_BOUNDARY_MIRROR) {
1352a6fd6b0aSBarry Smith             for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + (z - k - 1) * x * y + i * x;
1353c2859e5eSBarry Smith           } else {
13548865f1eaSKarl Rupp             for (j = 0; j < x; j++) idx[nn++] = -1;
135547c6ae99SBarry Smith           }
1356c2859e5eSBarry Smith         }
135747c6ae99SBarry Smith 
135847c6ae99SBarry Smith         if (n23 >= 0) { /* directly right */
1359ce00eea3SSatish Balay           x_t = lx[n23 % m];
136047c6ae99SBarry Smith           y_t = y;
136147c6ae99SBarry Smith           /* z_t = lz[n23 / (m*n)]; */
136247c6ae99SBarry Smith           s_t = bases[n23] + i * x_t + k * x_t * y_t;
13638865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1364ce00eea3SSatish Balay         } else if (xe - Xe < 0 && ze - Ze < 0) {
13658865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
136647c6ae99SBarry Smith         }
136747c6ae99SBarry Smith       }
136847c6ae99SBarry Smith 
136947c6ae99SBarry Smith       for (i = 1; i <= s_y; i++) {
137047c6ae99SBarry Smith         if (n24 >= 0) { /* left above */
1371ce00eea3SSatish Balay           x_t = lx[n24 % m];
137247c6ae99SBarry Smith           y_t = ly[(n24 % (m * n)) / m];
137347c6ae99SBarry Smith           /* z_t = lz[n24 / (m*n)]; */
137447c6ae99SBarry Smith           s_t = bases[n24] + i * x_t - s_x + k * x_t * y_t;
13758865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1376ce00eea3SSatish Balay         } else if (Xs - xs < 0 && ye - Ye < 0 && ze - Ze < 0) {
13778865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
137847c6ae99SBarry Smith         }
137947c6ae99SBarry Smith         if (n25 >= 0) { /* directly above */
138047c6ae99SBarry Smith           x_t = x;
138147c6ae99SBarry Smith           y_t = ly[(n25 % (m * n)) / m];
138247c6ae99SBarry Smith           /* z_t = lz[n25 / (m*n)]; */
138347c6ae99SBarry Smith           s_t = bases[n25] + (i - 1) * x_t + k * x_t * y_t;
13848865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1385ce00eea3SSatish Balay         } else if (ye - Ye < 0 && ze - Ze < 0) {
13868865f1eaSKarl Rupp           for (j = 0; j < x; j++) idx[nn++] = -1;
138747c6ae99SBarry Smith         }
138847c6ae99SBarry Smith         if (n26 >= 0) { /* right above */
1389ce00eea3SSatish Balay           x_t = lx[n26 % m];
139047c6ae99SBarry Smith           y_t = ly[(n26 % (m * n)) / m];
139147c6ae99SBarry Smith           /* z_t = lz[n26 / (m*n)]; */
139247c6ae99SBarry Smith           s_t = bases[n26] + (i - 1) * x_t + k * x_t * y_t;
13938865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1394ce00eea3SSatish Balay         } else if (xe - Xe < 0 && ye - Ye < 0 && ze - Ze < 0) {
13958865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
139647c6ae99SBarry Smith         }
139747c6ae99SBarry Smith       }
139847c6ae99SBarry Smith     }
139947c6ae99SBarry Smith   }
140047c6ae99SBarry Smith   /*
140147c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
140247c6ae99SBarry Smith      of VecSetValuesLocal().
140347c6ae99SBarry Smith   */
14049566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm, dof, nn, idx, PETSC_OWN_POINTER, &da->ltogmap));
140547c6ae99SBarry Smith 
14069566063dSJacob Faibussowitsch   PetscCall(PetscFree2(bases, ldims));
14079371c9d4SSatish Balay   dd->m  = m;
14089371c9d4SSatish Balay   dd->n  = n;
14099371c9d4SSatish Balay   dd->p  = p;
1410ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
14119371c9d4SSatish Balay   dd->xs = xs * dof;
14129371c9d4SSatish Balay   dd->xe = xe * dof;
14139371c9d4SSatish Balay   dd->ys = ys;
14149371c9d4SSatish Balay   dd->ye = ye;
14159371c9d4SSatish Balay   dd->zs = zs;
14169371c9d4SSatish Balay   dd->ze = ze;
14179371c9d4SSatish Balay   dd->Xs = Xs * dof;
14189371c9d4SSatish Balay   dd->Xe = Xe * dof;
14199371c9d4SSatish Balay   dd->Ys = Ys;
14209371c9d4SSatish Balay   dd->Ye = Ye;
14219371c9d4SSatish Balay   dd->Zs = Zs;
14229371c9d4SSatish Balay   dd->Ze = Ze;
1423ce00eea3SSatish Balay 
14249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&local));
14259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
1426ce00eea3SSatish Balay 
1427ce00eea3SSatish Balay   dd->gtol      = gtol;
1428ce00eea3SSatish Balay   dd->base      = base;
1429ce00eea3SSatish Balay   da->ops->view = DMView_DA_3d;
14300298fd71SBarry Smith   dd->ltol      = NULL;
14310298fd71SBarry Smith   dd->ao        = NULL;
143247c6ae99SBarry Smith   PetscFunctionReturn(0);
143347c6ae99SBarry Smith }
1434564755cdSBarry Smith 
143547c6ae99SBarry Smith /*@C
1436aa219208SBarry Smith    DMDACreate3d - Creates an object that will manage the communication of three-dimensional
143747c6ae99SBarry Smith    regular array data that is distributed across some processors.
143847c6ae99SBarry Smith 
1439d083f849SBarry Smith    Collective
144047c6ae99SBarry Smith 
144147c6ae99SBarry Smith    Input Parameters:
144247c6ae99SBarry Smith +  comm - MPI communicator
14431321219cSEthan Coon .  bx,by,bz - type of ghost nodes the array have.
1444bff4a2f0SMatthew G. Knepley          Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
1445aa219208SBarry Smith .  stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX)
1446897f7067SBarry Smith .  M,N,P - global dimension in each direction of the array
144747c6ae99SBarry Smith .  m,n,p - corresponding number of processors in each dimension
144847c6ae99SBarry Smith            (or PETSC_DECIDE to have calculated)
144947c6ae99SBarry Smith .  dof - number of degrees of freedom per node
145010d7c030SMatthew G Knepley .  s - stencil width
145110d7c030SMatthew G Knepley -  lx, ly, lz - arrays containing the number of nodes in each cell along
14520298fd71SBarry Smith           the x, y, and z coordinates, or NULL. If non-null, these
145347c6ae99SBarry Smith           must be of length as m,n,p and the corresponding
145447c6ae99SBarry Smith           m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of
145547c6ae99SBarry Smith           the ly[] must N, sum of the lz[] must be P
145647c6ae99SBarry Smith 
145747c6ae99SBarry Smith    Output Parameter:
145847c6ae99SBarry Smith .  da - the resulting distributed array object
145947c6ae99SBarry Smith 
146047c6ae99SBarry Smith    Options Database Key:
1461706a11cbSBarry Smith +  -dm_view - Calls DMView() at the conclusion of DMDACreate3d()
1462e43dc8daSBarry Smith .  -da_grid_x <nx> - number of grid points in x direction
1463e43dc8daSBarry Smith .  -da_grid_y <ny> - number of grid points in y direction
1464e43dc8daSBarry Smith .  -da_grid_z <nz> - number of grid points in z direction
146547c6ae99SBarry Smith .  -da_processors_x <MX> - number of processors in x direction
146647c6ae99SBarry Smith .  -da_processors_y <MY> - number of processors in y direction
146747c6ae99SBarry Smith .  -da_processors_z <MZ> - number of processors in z direction
1468e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
1469e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
1470e0f5d30fSBarry Smith .  -da_refine_z <rz>- refinement ratio in z directio
1471e43dc8daSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating it
147247c6ae99SBarry Smith 
147347c6ae99SBarry Smith    Level: beginner
147447c6ae99SBarry Smith 
147547c6ae99SBarry Smith    Notes:
1476aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1477aa219208SBarry Smith    standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
147847c6ae99SBarry Smith    the standard 27-pt stencil.
147947c6ae99SBarry Smith 
1480aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1481564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1482564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
148347c6ae99SBarry Smith 
1484897f7067SBarry Smith    You must call DMSetUp() after this call before using this DM.
1485897f7067SBarry Smith 
1486897f7067SBarry Smith    If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call
1487897f7067SBarry Smith    but before DMSetUp().
1488897f7067SBarry Smith 
1489db781477SPatrick Sanan .seealso: `DMDestroy()`, `DMView()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMGlobalToLocalBegin()`, `DMDAGetRefinementFactor()`,
1490db781477SPatrick Sanan           `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDASetRefinementFactor()`,
1491db781477SPatrick Sanan           `DMDAGetInfo()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMDACreateNaturalVector()`, `DMLoad()`, `DMDAGetOwnershipRanges()`,
1492db781477SPatrick Sanan           `DMStagCreate3d()`
149347c6ae99SBarry Smith 
149447c6ae99SBarry Smith @*/
14959371c9d4SSatish Balay PetscErrorCode DMDACreate3d(MPI_Comm comm, DMBoundaryType bx, DMBoundaryType by, DMBoundaryType bz, DMDAStencilType stencil_type, PetscInt M, PetscInt N, PetscInt P, PetscInt m, PetscInt n, PetscInt p, PetscInt dof, PetscInt s, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[], DM *da) {
149647c6ae99SBarry Smith   PetscFunctionBegin;
14979566063dSJacob Faibussowitsch   PetscCall(DMDACreate(comm, da));
14989566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*da, 3));
14999566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(*da, M, N, P));
15009566063dSJacob Faibussowitsch   PetscCall(DMDASetNumProcs(*da, m, n, p));
15019566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(*da, bx, by, bz));
15029566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(*da, dof));
15039566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilType(*da, stencil_type));
15049566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(*da, s));
15059566063dSJacob Faibussowitsch   PetscCall(DMDASetOwnershipRanges(*da, lx, ly, lz));
150647c6ae99SBarry Smith   PetscFunctionReturn(0);
150747c6ae99SBarry Smith }
1508