xref: /petsc/src/dm/impls/da/da3.c (revision 1690c2ae071c7584458d4e437df7b47bc4686b3c)
147c6ae99SBarry Smith /*
247c6ae99SBarry Smith    Code for manipulating distributed regular 3d arrays in parallel.
347c6ae99SBarry Smith    File created by Peter Mell  7/14/95
447c6ae99SBarry Smith  */
547c6ae99SBarry Smith 
6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"    I*/
747c6ae99SBarry Smith 
89804daf3SBarry Smith #include <petscdraw.h>
9d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMView_DA_3d(DM da, PetscViewer viewer)
10d71ae5a4SJacob Faibussowitsch {
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;
15d1e78c4fSBarry 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));
26d1e78c4fSBarry 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) {
35*1690c2aeSBarry Smith       PetscInt      i, nmax = 0, nmin = PETSC_INT_MAX, 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));
513ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
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));
863ba16761SJacob Faibussowitsch     if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
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) {
160a364092eSJacob Faibussowitsch             PetscCall(PetscSNPrintf(node, PETSC_STATIC_ARRAY_LENGTH(node), "%" PetscInt_FMT, 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));
183d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB)
1849a42bb27SBarry Smith   } else if (ismatlab) {
1859566063dSJacob Faibussowitsch     PetscCall(DMView_DA_Matlab(da, viewer));
1869a42bb27SBarry Smith #endif
18711aeaf0aSBarry Smith   }
1883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18947c6ae99SBarry Smith }
19047c6ae99SBarry Smith 
191d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetUp_DA_3D(DM da)
192d71ae5a4SJacob Faibussowitsch {
19347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA *)da->data;
19447c6ae99SBarry Smith   const PetscInt  M  = dd->M;
19547c6ae99SBarry Smith   const PetscInt  N  = dd->N;
19647c6ae99SBarry Smith   const PetscInt  P  = dd->P;
1976497c311SBarry Smith   PetscMPIInt     m, n, p;
19847c6ae99SBarry Smith   const PetscInt  dof          = dd->w;
19947c6ae99SBarry Smith   const PetscInt  s            = dd->s;
200bff4a2f0SMatthew G. Knepley   DMBoundaryType  bx           = dd->bx;
201bff4a2f0SMatthew G. Knepley   DMBoundaryType  by           = dd->by;
202bff4a2f0SMatthew G. Knepley   DMBoundaryType  bz           = dd->bz;
20319fd82e9SBarry Smith   DMDAStencilType stencil_type = dd->stencil_type;
20447c6ae99SBarry Smith   PetscInt       *lx           = dd->lx;
20547c6ae99SBarry Smith   PetscInt       *ly           = dd->ly;
20647c6ae99SBarry Smith   PetscInt       *lz           = dd->lz;
20747c6ae99SBarry Smith   MPI_Comm        comm;
20847c6ae99SBarry Smith   PetscMPIInt     rank, size;
209ce00eea3SSatish Balay   PetscInt        xs = 0, xe, ys = 0, ye, zs = 0, ze, x = 0, y = 0, z = 0;
210bd1fc5aeSBarry Smith   PetscInt        Xs, Xe, Ys, Ye, Zs, Ze, IXs, IXe, IYs, IYe, IZs, IZe, pm;
2118ea3bf28SBarry Smith   PetscInt        left, right, up, down, bottom, top, i, j, k, *idx, nn;
2126497c311SBarry Smith   PetscMPIInt     n0, n1, n2, n3, n4, n5, n6, n7, n8, n9, n10, n11, n12, n14;
2136497c311SBarry Smith   PetscMPIInt     n15, n16, n17, n18, n19, n20, n21, n22, n23, n24, n25, n26;
214db87c5ecSEthan Coon   PetscInt       *bases, *ldims, base, x_t, y_t, z_t, s_t, count, s_x, s_y, s_z;
2156497c311SBarry Smith   PetscMPIInt     sn0 = 0, sn1 = 0, sn2 = 0, sn3 = 0, sn5 = 0, sn6 = 0, sn7 = 0;
2166497c311SBarry Smith   PetscMPIInt     sn8 = 0, sn9 = 0, sn11 = 0, sn15 = 0, sn24 = 0, sn25 = 0, sn26 = 0;
2176497c311SBarry Smith   PetscMPIInt     sn17 = 0, sn18 = 0, sn19 = 0, sn20 = 0, sn21 = 0, sn23 = 0;
21847c6ae99SBarry Smith   Vec             local, global;
219bd1fc5aeSBarry Smith   VecScatter      gtol;
22045b6f7e9SBarry Smith   IS              to, from;
2216f951b95Secoon   PetscBool       twod;
22247c6ae99SBarry Smith 
22347c6ae99SBarry Smith   PetscFunctionBegin;
2247a8be351SBarry 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");
2259566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
2263855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
2277de69702SBarry Smith   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);
2283855c12bSBarry Smith #endif
2296497c311SBarry Smith   PetscCall(PetscMPIIntCast(dd->m, &m));
2306497c311SBarry Smith   PetscCall(PetscMPIIntCast(dd->n, &n));
2316497c311SBarry Smith   PetscCall(PetscMPIIntCast(dd->p, &p));
2323855c12bSBarry Smith 
2339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
23547c6ae99SBarry Smith 
23647c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
2376497c311SBarry Smith     PetscCheck(m >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in X direction: %d", m);
2386497c311SBarry Smith     PetscCheck(m <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in X direction: %d %d", m, size);
23947c6ae99SBarry Smith   }
24047c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
2416497c311SBarry Smith     PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in Y direction: %d", n);
2426497c311SBarry Smith     PetscCheck(n <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in Y direction: %d %d", n, size);
24347c6ae99SBarry Smith   }
24447c6ae99SBarry Smith   if (p != PETSC_DECIDE) {
2456497c311SBarry Smith     PetscCheck(p >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in Z direction: %d", p);
2466497c311SBarry Smith     PetscCheck(p <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in Z direction: %d %d", p, size);
24747c6ae99SBarry Smith   }
2486497c311SBarry Smith   PetscCheck(m <= 0 || n <= 0 || p <= 0 || m * n * p == size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "m %d * n %d * p %d != size %d", m, n, p, size);
24947c6ae99SBarry Smith 
25047c6ae99SBarry Smith   /* Partition the array among the processors */
25147c6ae99SBarry Smith   if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) {
25247c6ae99SBarry Smith     m = size / (n * p);
25347c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
25447c6ae99SBarry Smith     n = size / (m * p);
25547c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
25647c6ae99SBarry Smith     p = size / (m * n);
25747c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
25847c6ae99SBarry Smith     /* try for squarish distribution */
259369cc0aeSBarry Smith     m = (int)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)N * p)));
26047c6ae99SBarry Smith     if (!m) m = 1;
26147c6ae99SBarry Smith     while (m > 0) {
26247c6ae99SBarry Smith       n = size / (m * p);
26347c6ae99SBarry Smith       if (m * n * p == size) break;
26447c6ae99SBarry Smith       m--;
26547c6ae99SBarry Smith     }
2666497c311SBarry Smith     PetscCheck(m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad p value: p = %d", p);
2679371c9d4SSatish Balay     if (M > N && m < n) {
2686497c311SBarry Smith       PetscMPIInt _m = m;
2699371c9d4SSatish Balay       m              = n;
2709371c9d4SSatish Balay       n              = _m;
2719371c9d4SSatish Balay     }
27247c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
27347c6ae99SBarry Smith     /* try for squarish distribution */
274369cc0aeSBarry Smith     m = (int)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)P * n)));
27547c6ae99SBarry Smith     if (!m) m = 1;
27647c6ae99SBarry Smith     while (m > 0) {
27747c6ae99SBarry Smith       p = size / (m * n);
27847c6ae99SBarry Smith       if (m * n * p == size) break;
27947c6ae99SBarry Smith       m--;
28047c6ae99SBarry Smith     }
2816497c311SBarry Smith     PetscCheck(m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad n value: n = %d", n);
2829371c9d4SSatish Balay     if (M > P && m < p) {
2836497c311SBarry Smith       PetscMPIInt _m = m;
2849371c9d4SSatish Balay       m              = p;
2859371c9d4SSatish Balay       p              = _m;
2869371c9d4SSatish Balay     }
28747c6ae99SBarry Smith   } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
28847c6ae99SBarry Smith     /* try for squarish distribution */
289369cc0aeSBarry Smith     n = (int)(0.5 + PetscSqrtReal(((PetscReal)N) * ((PetscReal)size) / ((PetscReal)P * m)));
29047c6ae99SBarry Smith     if (!n) n = 1;
29147c6ae99SBarry Smith     while (n > 0) {
29247c6ae99SBarry Smith       p = size / (m * n);
29347c6ae99SBarry Smith       if (m * n * p == size) break;
29447c6ae99SBarry Smith       n--;
29547c6ae99SBarry Smith     }
2966497c311SBarry Smith     PetscCheck(n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad m value: m = %d", n);
2979371c9d4SSatish Balay     if (N > P && n < p) {
2986497c311SBarry Smith       PetscMPIInt _n = n;
2999371c9d4SSatish Balay       n              = p;
3009371c9d4SSatish Balay       p              = _n;
3019371c9d4SSatish Balay     }
30247c6ae99SBarry Smith   } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
30347c6ae99SBarry Smith     /* try for squarish distribution */
3046497c311SBarry Smith     n = (PetscMPIInt)(0.5 + PetscPowReal(((PetscReal)N * N) * ((PetscReal)size) / ((PetscReal)P * M), (PetscReal)(1. / 3.)));
30547c6ae99SBarry Smith     if (!n) n = 1;
30647c6ae99SBarry Smith     while (n > 0) {
30747c6ae99SBarry Smith       pm = size / n;
30847c6ae99SBarry Smith       if (n * pm == size) break;
30947c6ae99SBarry Smith       n--;
31047c6ae99SBarry Smith     }
31147c6ae99SBarry Smith     if (!n) n = 1;
3126497c311SBarry Smith     m = (PetscMPIInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)P * n)));
31347c6ae99SBarry Smith     if (!m) m = 1;
31447c6ae99SBarry Smith     while (m > 0) {
31547c6ae99SBarry Smith       p = size / (m * n);
31647c6ae99SBarry Smith       if (m * n * p == size) break;
31747c6ae99SBarry Smith       m--;
31847c6ae99SBarry Smith     }
3199371c9d4SSatish Balay     if (M > P && m < p) {
3206497c311SBarry Smith       PetscMPIInt _m = m;
3219371c9d4SSatish Balay       m              = p;
3229371c9d4SSatish Balay       p              = _m;
3239371c9d4SSatish Balay     }
32408401ef6SPierre Jolivet   } else PetscCheck(m * n * p == size, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Given Bad partition");
32547c6ae99SBarry Smith 
32608401ef6SPierre Jolivet   PetscCheck(m * n * p == size, PetscObjectComm((PetscObject)da), PETSC_ERR_PLIB, "Could not find good partition");
3276497c311SBarry Smith   PetscCheck(M >= m, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Partition in x direction is too fine! %" PetscInt_FMT " %d", M, m);
3286497c311SBarry Smith   PetscCheck(N >= n, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Partition in y direction is too fine! %" PetscInt_FMT " %d", N, n);
3296497c311SBarry Smith   PetscCheck(P >= p, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Partition in z direction is too fine! %" PetscInt_FMT " %d", P, p);
33047c6ae99SBarry Smith 
33147c6ae99SBarry Smith   /*
33247c6ae99SBarry Smith      Determine locally owned region
33347c6ae99SBarry Smith      [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes
33447c6ae99SBarry Smith   */
33547c6ae99SBarry Smith 
33647c6ae99SBarry Smith   if (!lx) {
3379566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &dd->lx));
33847c6ae99SBarry Smith     lx = dd->lx;
3398865f1eaSKarl Rupp     for (i = 0; i < m; i++) lx[i] = M / m + ((M % m) > (i % m));
34047c6ae99SBarry Smith   }
34147c6ae99SBarry Smith   x  = lx[rank % m];
34247c6ae99SBarry Smith   xs = 0;
3438865f1eaSKarl Rupp   for (i = 0; i < (rank % m); i++) xs += lx[i];
3441dca8a05SBarry 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);
34547c6ae99SBarry Smith 
34647c6ae99SBarry Smith   if (!ly) {
3479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &dd->ly));
34847c6ae99SBarry Smith     ly = dd->ly;
3498865f1eaSKarl Rupp     for (i = 0; i < n; i++) ly[i] = N / n + ((N % n) > (i % n));
35047c6ae99SBarry Smith   }
35147c6ae99SBarry Smith   y = ly[(rank % (m * n)) / m];
3521dca8a05SBarry 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);
35330729d88SBarry Smith 
35447c6ae99SBarry Smith   ys = 0;
3558865f1eaSKarl Rupp   for (i = 0; i < (rank % (m * n)) / m; i++) ys += ly[i];
35647c6ae99SBarry Smith 
35747c6ae99SBarry Smith   if (!lz) {
3589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(p, &dd->lz));
35947c6ae99SBarry Smith     lz = dd->lz;
3608865f1eaSKarl Rupp     for (i = 0; i < p; i++) lz[i] = P / p + ((P % p) > (i % p));
36147c6ae99SBarry Smith   }
36247c6ae99SBarry Smith   z = lz[rank / (m * n)];
363bcea557cSEthan Coon 
3644ad8454bSPierre Jolivet   PetscCheck((x > s) || (bx != DM_BOUNDARY_MIRROR), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local x-width of domain x %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT " with mirror", x, s);
3654ad8454bSPierre Jolivet   PetscCheck((y > s) || (by != DM_BOUNDARY_MIRROR), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local y-width of domain y %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT " with mirror", y, s);
3664ad8454bSPierre Jolivet   PetscCheck((z > s) || (bz != DM_BOUNDARY_MIRROR), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local z-width of domain z %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT " with mirror", z, s);
36718f812f4SBarry Smith 
368fdc81ce8SEthan Coon   /* note this is different than x- and y-, as we will handle as an important special
369bff4a2f0SMatthew G. Knepley    case when p=P=1 and DM_BOUNDARY_PERIODIC and s > z.  This is to deal with 2D problems
370fdc81ce8SEthan Coon    in a 3D code.  Additional code for this case is noted with "2d case" comments */
3716f951b95Secoon   twod = PETSC_FALSE;
3728865f1eaSKarl Rupp   if (P == 1) twod = PETSC_TRUE;
3731dca8a05SBarry 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);
37447c6ae99SBarry Smith   zs = 0;
3758865f1eaSKarl Rupp   for (i = 0; i < (rank / (m * n)); i++) zs += lz[i];
37647c6ae99SBarry Smith   ye = ys + y;
37747c6ae99SBarry Smith   xe = xs + x;
37847c6ae99SBarry Smith   ze = zs + z;
37947c6ae99SBarry Smith 
38088661749SPeter Brune   /* determine ghost region (Xs) and region scattered into (IXs)  */
381d9c9ebe5SPeter Brune   if (xs - s > 0) {
3829371c9d4SSatish Balay     Xs  = xs - s;
3839371c9d4SSatish Balay     IXs = xs - s;
38488661749SPeter Brune   } else {
3858865f1eaSKarl Rupp     if (bx) Xs = xs - s;
3868865f1eaSKarl Rupp     else Xs = 0;
38788661749SPeter Brune     IXs = 0;
38888661749SPeter Brune   }
389d9c9ebe5SPeter Brune   if (xe + s <= M) {
3909371c9d4SSatish Balay     Xe  = xe + s;
3919371c9d4SSatish Balay     IXe = xe + s;
39288661749SPeter Brune   } else {
39388661749SPeter Brune     if (bx) {
3949371c9d4SSatish Balay       Xs = xs - s;
3959371c9d4SSatish Balay       Xe = xe + s;
3968865f1eaSKarl Rupp     } else Xe = M;
39788661749SPeter Brune     IXe = M;
39888661749SPeter Brune   }
39947c6ae99SBarry Smith 
400bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) {
401d9c9ebe5SPeter Brune     IXs = xs - s;
402d9c9ebe5SPeter Brune     IXe = xe + s;
403d9c9ebe5SPeter Brune     Xs  = xs - s;
404d9c9ebe5SPeter Brune     Xe  = xe + s;
40588661749SPeter Brune   }
40688661749SPeter Brune 
407d9c9ebe5SPeter Brune   if (ys - s > 0) {
4089371c9d4SSatish Balay     Ys  = ys - s;
4099371c9d4SSatish Balay     IYs = ys - s;
41088661749SPeter Brune   } else {
4118865f1eaSKarl Rupp     if (by) Ys = ys - s;
4128865f1eaSKarl Rupp     else Ys = 0;
41388661749SPeter Brune     IYs = 0;
41488661749SPeter Brune   }
415d9c9ebe5SPeter Brune   if (ye + s <= N) {
4169371c9d4SSatish Balay     Ye  = ye + s;
4179371c9d4SSatish Balay     IYe = ye + s;
41888661749SPeter Brune   } else {
4198865f1eaSKarl Rupp     if (by) Ye = ye + s;
4208865f1eaSKarl Rupp     else Ye = N;
42188661749SPeter Brune     IYe = N;
42288661749SPeter Brune   }
42388661749SPeter Brune 
424bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) {
425d9c9ebe5SPeter Brune     IYs = ys - s;
426d9c9ebe5SPeter Brune     IYe = ye + s;
427d9c9ebe5SPeter Brune     Ys  = ys - s;
428d9c9ebe5SPeter Brune     Ye  = ye + s;
42988661749SPeter Brune   }
43088661749SPeter Brune 
431d9c9ebe5SPeter Brune   if (zs - s > 0) {
4329371c9d4SSatish Balay     Zs  = zs - s;
4339371c9d4SSatish Balay     IZs = zs - s;
43488661749SPeter Brune   } else {
4358865f1eaSKarl Rupp     if (bz) Zs = zs - s;
4368865f1eaSKarl Rupp     else Zs = 0;
43788661749SPeter Brune     IZs = 0;
43888661749SPeter Brune   }
439d9c9ebe5SPeter Brune   if (ze + s <= P) {
4409371c9d4SSatish Balay     Ze  = ze + s;
4419371c9d4SSatish Balay     IZe = ze + s;
44288661749SPeter Brune   } else {
4438865f1eaSKarl Rupp     if (bz) Ze = ze + s;
4448865f1eaSKarl Rupp     else Ze = P;
44588661749SPeter Brune     IZe = P;
44688661749SPeter Brune   }
44788661749SPeter Brune 
448bff4a2f0SMatthew G. Knepley   if (bz == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_MIRROR) {
449d9c9ebe5SPeter Brune     IZs = zs - s;
450d9c9ebe5SPeter Brune     IZe = ze + s;
451d9c9ebe5SPeter Brune     Zs  = zs - s;
452d9c9ebe5SPeter Brune     Ze  = ze + s;
45388661749SPeter Brune   }
45447c6ae99SBarry Smith 
45547c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
456d9c9ebe5SPeter Brune   s_x = s;
457d9c9ebe5SPeter Brune   s_y = s;
458d9c9ebe5SPeter Brune   s_z = s;
45947c6ae99SBarry Smith 
46047c6ae99SBarry Smith   /* determine starting point of each processor */
46147c6ae99SBarry Smith   nn = x * y * z;
4629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(size + 1, &bases, size, &ldims));
4639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&nn, 1, MPIU_INT, ldims, 1, MPIU_INT, comm));
46447c6ae99SBarry Smith   bases[0] = 0;
4658865f1eaSKarl Rupp   for (i = 1; i <= size; i++) bases[i] = ldims[i - 1];
4668865f1eaSKarl Rupp   for (i = 1; i <= size; i++) bases[i] += bases[i - 1];
467ce00eea3SSatish Balay   base = bases[rank] * dof;
46847c6ae99SBarry Smith 
46947c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
470ce00eea3SSatish Balay   dd->Nlocal = x * y * z * dof;
4719566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(comm, dof, dd->Nlocal, PETSC_DECIDE, NULL, &global));
472ce00eea3SSatish Balay   dd->nlocal = (Xe - Xs) * (Ye - Ys) * (Ze - Zs) * dof;
4739566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, dof, dd->nlocal, NULL, &local));
47447c6ae99SBarry Smith 
4754104a7a0SPatrick Sanan   /* generate global to local vector scatter and local to global mapping*/
47647c6ae99SBarry Smith 
477ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
478ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
4799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((IXe - IXs) * (IYe - IYs) * (IZe - IZs), &idx));
480d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_BOX) {
4819371c9d4SSatish Balay     left   = IXs - Xs;
4829371c9d4SSatish Balay     right  = left + (IXe - IXs);
4839371c9d4SSatish Balay     bottom = IYs - Ys;
4849371c9d4SSatish Balay     top    = bottom + (IYe - IYs);
4859371c9d4SSatish Balay     down   = IZs - Zs;
4869371c9d4SSatish Balay     up     = down + (IZe - IZs);
487ce00eea3SSatish Balay     count  = 0;
488ce00eea3SSatish Balay     for (i = down; i < up; i++) {
489ce00eea3SSatish Balay       for (j = bottom; j < top; j++) {
490ad540459SPierre Jolivet         for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k;
491ce00eea3SSatish Balay       }
492ce00eea3SSatish Balay     }
4939566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(comm, dof, count, idx, PETSC_OWN_POINTER, &to));
49447c6ae99SBarry Smith   } else {
49547c6ae99SBarry Smith     /* This is way ugly! We need to list the funny cross type region */
4969371c9d4SSatish Balay     left   = xs - Xs;
4979371c9d4SSatish Balay     right  = left + x;
4989371c9d4SSatish Balay     bottom = ys - Ys;
4999371c9d4SSatish Balay     top    = bottom + y;
5009371c9d4SSatish Balay     down   = zs - Zs;
5019371c9d4SSatish Balay     up     = down + z;
50247c6ae99SBarry Smith     count  = 0;
503a5b23f4aSJose E. Roman     /* the bottom chunk */
504ce00eea3SSatish Balay     for (i = (IZs - Zs); i < down; i++) {
50547c6ae99SBarry Smith       for (j = bottom; j < top; j++) {
506ce00eea3SSatish Balay         for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k;
50747c6ae99SBarry Smith       }
50847c6ae99SBarry Smith     }
50947c6ae99SBarry Smith     /* the middle piece */
51047c6ae99SBarry Smith     for (i = down; i < up; i++) {
51147c6ae99SBarry Smith       /* front */
512ce00eea3SSatish Balay       for (j = (IYs - Ys); j < bottom; j++) {
513ce00eea3SSatish Balay         for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k;
51447c6ae99SBarry Smith       }
51547c6ae99SBarry Smith       /* middle */
51647c6ae99SBarry Smith       for (j = bottom; j < top; j++) {
517ce00eea3SSatish Balay         for (k = IXs - Xs; k < IXe - Xs; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k;
51847c6ae99SBarry Smith       }
51947c6ae99SBarry Smith       /* back */
520ce00eea3SSatish Balay       for (j = top; j < top + IYe - ye; j++) {
521ce00eea3SSatish Balay         for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k;
52247c6ae99SBarry Smith       }
52347c6ae99SBarry Smith     }
52447c6ae99SBarry Smith     /* the top piece */
525ce00eea3SSatish Balay     for (i = up; i < up + IZe - ze; i++) {
52647c6ae99SBarry Smith       for (j = bottom; j < top; j++) {
527ce00eea3SSatish Balay         for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k;
52847c6ae99SBarry Smith       }
52947c6ae99SBarry Smith     }
5309566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(comm, dof, count, idx, PETSC_OWN_POINTER, &to));
53147c6ae99SBarry Smith   }
53247c6ae99SBarry Smith 
53347c6ae99SBarry Smith   /* determine who lies on each side of use stored in    n24 n25 n26
53447c6ae99SBarry Smith                                                          n21 n22 n23
53547c6ae99SBarry Smith                                                          n18 n19 n20
53647c6ae99SBarry Smith 
53747c6ae99SBarry Smith                                                          n15 n16 n17
53847c6ae99SBarry Smith                                                          n12     n14
53947c6ae99SBarry Smith                                                          n9  n10 n11
54047c6ae99SBarry Smith 
54147c6ae99SBarry Smith                                                          n6  n7  n8
54247c6ae99SBarry Smith                                                          n3  n4  n5
54347c6ae99SBarry Smith                                                          n0  n1  n2
54447c6ae99SBarry Smith   */
54547c6ae99SBarry Smith 
54647c6ae99SBarry Smith   /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
54747c6ae99SBarry Smith   /* Assume Nodes are Internal to the Cube */
54847c6ae99SBarry Smith   n0 = rank - m * n - m - 1;
54947c6ae99SBarry Smith   n1 = rank - m * n - m;
55047c6ae99SBarry Smith   n2 = rank - m * n - m + 1;
55147c6ae99SBarry Smith   n3 = rank - m * n - 1;
55247c6ae99SBarry Smith   n4 = rank - m * n;
55347c6ae99SBarry Smith   n5 = rank - m * n + 1;
55447c6ae99SBarry Smith   n6 = rank - m * n + m - 1;
55547c6ae99SBarry Smith   n7 = rank - m * n + m;
55647c6ae99SBarry Smith   n8 = rank - m * n + m + 1;
55747c6ae99SBarry Smith 
55847c6ae99SBarry Smith   n9  = rank - m - 1;
55947c6ae99SBarry Smith   n10 = rank - m;
56047c6ae99SBarry Smith   n11 = rank - m + 1;
56147c6ae99SBarry Smith   n12 = rank - 1;
56247c6ae99SBarry Smith   n14 = rank + 1;
56347c6ae99SBarry Smith   n15 = rank + m - 1;
56447c6ae99SBarry Smith   n16 = rank + m;
56547c6ae99SBarry Smith   n17 = rank + m + 1;
56647c6ae99SBarry Smith 
56747c6ae99SBarry Smith   n18 = rank + m * n - m - 1;
56847c6ae99SBarry Smith   n19 = rank + m * n - m;
56947c6ae99SBarry Smith   n20 = rank + m * n - m + 1;
57047c6ae99SBarry Smith   n21 = rank + m * n - 1;
57147c6ae99SBarry Smith   n22 = rank + m * n;
57247c6ae99SBarry Smith   n23 = rank + m * n + 1;
57347c6ae99SBarry Smith   n24 = rank + m * n + m - 1;
57447c6ae99SBarry Smith   n25 = rank + m * n + m;
57547c6ae99SBarry Smith   n26 = rank + m * n + m + 1;
57647c6ae99SBarry Smith 
57747c6ae99SBarry Smith   /* Assume Pieces are on Faces of Cube */
57847c6ae99SBarry Smith 
57947c6ae99SBarry Smith   if (xs == 0) { /* First assume not corner or edge */
58047c6ae99SBarry Smith     n0  = rank - 1 - (m * n);
58147c6ae99SBarry Smith     n3  = rank + m - 1 - (m * n);
58247c6ae99SBarry Smith     n6  = rank + 2 * m - 1 - (m * n);
58347c6ae99SBarry Smith     n9  = rank - 1;
58447c6ae99SBarry Smith     n12 = rank + m - 1;
58547c6ae99SBarry Smith     n15 = rank + 2 * m - 1;
58647c6ae99SBarry Smith     n18 = rank - 1 + (m * n);
58747c6ae99SBarry Smith     n21 = rank + m - 1 + (m * n);
58847c6ae99SBarry Smith     n24 = rank + 2 * m - 1 + (m * n);
58947c6ae99SBarry Smith   }
59047c6ae99SBarry Smith 
591ce00eea3SSatish Balay   if (xe == M) { /* First assume not corner or edge */
59247c6ae99SBarry Smith     n2  = rank - 2 * m + 1 - (m * n);
59347c6ae99SBarry Smith     n5  = rank - m + 1 - (m * n);
59447c6ae99SBarry Smith     n8  = rank + 1 - (m * n);
59547c6ae99SBarry Smith     n11 = rank - 2 * m + 1;
59647c6ae99SBarry Smith     n14 = rank - m + 1;
59747c6ae99SBarry Smith     n17 = rank + 1;
59847c6ae99SBarry Smith     n20 = rank - 2 * m + 1 + (m * n);
59947c6ae99SBarry Smith     n23 = rank - m + 1 + (m * n);
60047c6ae99SBarry Smith     n26 = rank + 1 + (m * n);
60147c6ae99SBarry Smith   }
60247c6ae99SBarry Smith 
60347c6ae99SBarry Smith   if (ys == 0) { /* First assume not corner or edge */
60447c6ae99SBarry Smith     n0  = rank + m * (n - 1) - 1 - (m * n);
60547c6ae99SBarry Smith     n1  = rank + m * (n - 1) - (m * n);
60647c6ae99SBarry Smith     n2  = rank + m * (n - 1) + 1 - (m * n);
60747c6ae99SBarry Smith     n9  = rank + m * (n - 1) - 1;
60847c6ae99SBarry Smith     n10 = rank + m * (n - 1);
60947c6ae99SBarry Smith     n11 = rank + m * (n - 1) + 1;
61047c6ae99SBarry Smith     n18 = rank + m * (n - 1) - 1 + (m * n);
61147c6ae99SBarry Smith     n19 = rank + m * (n - 1) + (m * n);
61247c6ae99SBarry Smith     n20 = rank + m * (n - 1) + 1 + (m * n);
61347c6ae99SBarry Smith   }
61447c6ae99SBarry Smith 
61547c6ae99SBarry Smith   if (ye == N) { /* First assume not corner or edge */
61647c6ae99SBarry Smith     n6  = rank - m * (n - 1) - 1 - (m * n);
61747c6ae99SBarry Smith     n7  = rank - m * (n - 1) - (m * n);
61847c6ae99SBarry Smith     n8  = rank - m * (n - 1) + 1 - (m * n);
61947c6ae99SBarry Smith     n15 = rank - m * (n - 1) - 1;
62047c6ae99SBarry Smith     n16 = rank - m * (n - 1);
62147c6ae99SBarry Smith     n17 = rank - m * (n - 1) + 1;
62247c6ae99SBarry Smith     n24 = rank - m * (n - 1) - 1 + (m * n);
62347c6ae99SBarry Smith     n25 = rank - m * (n - 1) + (m * n);
62447c6ae99SBarry Smith     n26 = rank - m * (n - 1) + 1 + (m * n);
62547c6ae99SBarry Smith   }
62647c6ae99SBarry Smith 
62747c6ae99SBarry Smith   if (zs == 0) { /* First assume not corner or edge */
62847c6ae99SBarry Smith     n0 = size - (m * n) + rank - m - 1;
62947c6ae99SBarry Smith     n1 = size - (m * n) + rank - m;
63047c6ae99SBarry Smith     n2 = size - (m * n) + rank - m + 1;
63147c6ae99SBarry Smith     n3 = size - (m * n) + rank - 1;
63247c6ae99SBarry Smith     n4 = size - (m * n) + rank;
63347c6ae99SBarry Smith     n5 = size - (m * n) + rank + 1;
63447c6ae99SBarry Smith     n6 = size - (m * n) + rank + m - 1;
63547c6ae99SBarry Smith     n7 = size - (m * n) + rank + m;
63647c6ae99SBarry Smith     n8 = size - (m * n) + rank + m + 1;
63747c6ae99SBarry Smith   }
63847c6ae99SBarry Smith 
63947c6ae99SBarry Smith   if (ze == P) { /* First assume not corner or edge */
64047c6ae99SBarry Smith     n18 = (m * n) - (size - rank) - m - 1;
64147c6ae99SBarry Smith     n19 = (m * n) - (size - rank) - m;
64247c6ae99SBarry Smith     n20 = (m * n) - (size - rank) - m + 1;
64347c6ae99SBarry Smith     n21 = (m * n) - (size - rank) - 1;
64447c6ae99SBarry Smith     n22 = (m * n) - (size - rank);
64547c6ae99SBarry Smith     n23 = (m * n) - (size - rank) + 1;
64647c6ae99SBarry Smith     n24 = (m * n) - (size - rank) + m - 1;
64747c6ae99SBarry Smith     n25 = (m * n) - (size - rank) + m;
64847c6ae99SBarry Smith     n26 = (m * n) - (size - rank) + m + 1;
64947c6ae99SBarry Smith   }
65047c6ae99SBarry Smith 
65147c6ae99SBarry Smith   if ((xs == 0) && (zs == 0)) { /* Assume an edge, not corner */
65247c6ae99SBarry Smith     n0 = size - m * n + rank + m - 1 - m;
65347c6ae99SBarry Smith     n3 = size - m * n + rank + m - 1;
65447c6ae99SBarry Smith     n6 = size - m * n + rank + m - 1 + m;
65547c6ae99SBarry Smith   }
65647c6ae99SBarry Smith 
65747c6ae99SBarry Smith   if ((xs == 0) && (ze == P)) { /* Assume an edge, not corner */
65847c6ae99SBarry Smith     n18 = m * n - (size - rank) + m - 1 - m;
65947c6ae99SBarry Smith     n21 = m * n - (size - rank) + m - 1;
66047c6ae99SBarry Smith     n24 = m * n - (size - rank) + m - 1 + m;
66147c6ae99SBarry Smith   }
66247c6ae99SBarry Smith 
66347c6ae99SBarry Smith   if ((xs == 0) && (ys == 0)) { /* Assume an edge, not corner */
66447c6ae99SBarry Smith     n0  = rank + m * n - 1 - m * n;
66547c6ae99SBarry Smith     n9  = rank + m * n - 1;
66647c6ae99SBarry Smith     n18 = rank + m * n - 1 + m * n;
66747c6ae99SBarry Smith   }
66847c6ae99SBarry Smith 
66947c6ae99SBarry Smith   if ((xs == 0) && (ye == N)) { /* Assume an edge, not corner */
67047c6ae99SBarry Smith     n6  = rank - m * (n - 1) + m - 1 - m * n;
67147c6ae99SBarry Smith     n15 = rank - m * (n - 1) + m - 1;
67247c6ae99SBarry Smith     n24 = rank - m * (n - 1) + m - 1 + m * n;
67347c6ae99SBarry Smith   }
67447c6ae99SBarry Smith 
675ce00eea3SSatish Balay   if ((xe == M) && (zs == 0)) { /* Assume an edge, not corner */
67647c6ae99SBarry Smith     n2 = size - (m * n - rank) - (m - 1) - m;
67747c6ae99SBarry Smith     n5 = size - (m * n - rank) - (m - 1);
67847c6ae99SBarry Smith     n8 = size - (m * n - rank) - (m - 1) + m;
67947c6ae99SBarry Smith   }
68047c6ae99SBarry Smith 
681ce00eea3SSatish Balay   if ((xe == M) && (ze == P)) { /* Assume an edge, not corner */
68247c6ae99SBarry Smith     n20 = m * n - (size - rank) - (m - 1) - m;
68347c6ae99SBarry Smith     n23 = m * n - (size - rank) - (m - 1);
68447c6ae99SBarry Smith     n26 = m * n - (size - rank) - (m - 1) + m;
68547c6ae99SBarry Smith   }
68647c6ae99SBarry Smith 
687ce00eea3SSatish Balay   if ((xe == M) && (ys == 0)) { /* Assume an edge, not corner */
68847c6ae99SBarry Smith     n2  = rank + m * (n - 1) - (m - 1) - m * n;
68947c6ae99SBarry Smith     n11 = rank + m * (n - 1) - (m - 1);
69047c6ae99SBarry Smith     n20 = rank + m * (n - 1) - (m - 1) + m * n;
69147c6ae99SBarry Smith   }
69247c6ae99SBarry Smith 
693ce00eea3SSatish Balay   if ((xe == M) && (ye == N)) { /* Assume an edge, not corner */
69447c6ae99SBarry Smith     n8  = rank - m * n + 1 - m * n;
69547c6ae99SBarry Smith     n17 = rank - m * n + 1;
69647c6ae99SBarry Smith     n26 = rank - m * n + 1 + m * n;
69747c6ae99SBarry Smith   }
69847c6ae99SBarry Smith 
69947c6ae99SBarry Smith   if ((ys == 0) && (zs == 0)) { /* Assume an edge, not corner */
70047c6ae99SBarry Smith     n0 = size - m + rank - 1;
70147c6ae99SBarry Smith     n1 = size - m + rank;
70247c6ae99SBarry Smith     n2 = size - m + rank + 1;
70347c6ae99SBarry Smith   }
70447c6ae99SBarry Smith 
70547c6ae99SBarry Smith   if ((ys == 0) && (ze == P)) { /* Assume an edge, not corner */
70647c6ae99SBarry Smith     n18 = m * n - (size - rank) + m * (n - 1) - 1;
70747c6ae99SBarry Smith     n19 = m * n - (size - rank) + m * (n - 1);
70847c6ae99SBarry Smith     n20 = m * n - (size - rank) + m * (n - 1) + 1;
70947c6ae99SBarry Smith   }
71047c6ae99SBarry Smith 
71147c6ae99SBarry Smith   if ((ye == N) && (zs == 0)) { /* Assume an edge, not corner */
71247c6ae99SBarry Smith     n6 = size - (m * n - rank) - m * (n - 1) - 1;
71347c6ae99SBarry Smith     n7 = size - (m * n - rank) - m * (n - 1);
71447c6ae99SBarry Smith     n8 = size - (m * n - rank) - m * (n - 1) + 1;
71547c6ae99SBarry Smith   }
71647c6ae99SBarry Smith 
71747c6ae99SBarry Smith   if ((ye == N) && (ze == P)) { /* Assume an edge, not corner */
71847c6ae99SBarry Smith     n24 = rank - (size - m) - 1;
71947c6ae99SBarry Smith     n25 = rank - (size - m);
72047c6ae99SBarry Smith     n26 = rank - (size - m) + 1;
72147c6ae99SBarry Smith   }
72247c6ae99SBarry Smith 
72347c6ae99SBarry Smith   /* Check for Corners */
7248865f1eaSKarl Rupp   if ((xs == 0) && (ys == 0) && (zs == 0)) n0 = size - 1;
7258865f1eaSKarl Rupp   if ((xs == 0) && (ys == 0) && (ze == P)) n18 = m * n - 1;
7268865f1eaSKarl Rupp   if ((xs == 0) && (ye == N) && (zs == 0)) n6 = (size - 1) - m * (n - 1);
7278865f1eaSKarl Rupp   if ((xs == 0) && (ye == N) && (ze == P)) n24 = m - 1;
7288865f1eaSKarl Rupp   if ((xe == M) && (ys == 0) && (zs == 0)) n2 = size - m;
7298865f1eaSKarl Rupp   if ((xe == M) && (ys == 0) && (ze == P)) n20 = m * n - m;
7308865f1eaSKarl Rupp   if ((xe == M) && (ye == N) && (zs == 0)) n8 = size - m * n;
7318865f1eaSKarl Rupp   if ((xe == M) && (ye == N) && (ze == P)) n26 = 0;
73247c6ae99SBarry Smith 
73347c6ae99SBarry Smith   /* Check for when not X,Y, and Z Periodic */
73447c6ae99SBarry Smith 
73547c6ae99SBarry Smith   /* If not X periodic */
736bff4a2f0SMatthew G. Knepley   if (bx != DM_BOUNDARY_PERIODIC) {
7378865f1eaSKarl Rupp     if (xs == 0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2;
7388865f1eaSKarl Rupp     if (xe == M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2;
73947c6ae99SBarry Smith   }
74047c6ae99SBarry Smith 
74147c6ae99SBarry Smith   /* If not Y periodic */
742bff4a2f0SMatthew G. Knepley   if (by != DM_BOUNDARY_PERIODIC) {
7438865f1eaSKarl Rupp     if (ys == 0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2;
7448865f1eaSKarl Rupp     if (ye == N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2;
74547c6ae99SBarry Smith   }
74647c6ae99SBarry Smith 
74747c6ae99SBarry Smith   /* If not Z periodic */
748bff4a2f0SMatthew G. Knepley   if (bz != DM_BOUNDARY_PERIODIC) {
7498865f1eaSKarl Rupp     if (zs == 0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2;
7508865f1eaSKarl Rupp     if (ze == P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;
75147c6ae99SBarry Smith   }
75247c6ae99SBarry Smith 
7539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(27, &dd->neighbors));
7548865f1eaSKarl Rupp 
75547c6ae99SBarry Smith   dd->neighbors[0]  = n0;
75647c6ae99SBarry Smith   dd->neighbors[1]  = n1;
75747c6ae99SBarry Smith   dd->neighbors[2]  = n2;
75847c6ae99SBarry Smith   dd->neighbors[3]  = n3;
75947c6ae99SBarry Smith   dd->neighbors[4]  = n4;
76047c6ae99SBarry Smith   dd->neighbors[5]  = n5;
76147c6ae99SBarry Smith   dd->neighbors[6]  = n6;
76247c6ae99SBarry Smith   dd->neighbors[7]  = n7;
76347c6ae99SBarry Smith   dd->neighbors[8]  = n8;
76447c6ae99SBarry Smith   dd->neighbors[9]  = n9;
76547c6ae99SBarry Smith   dd->neighbors[10] = n10;
76647c6ae99SBarry Smith   dd->neighbors[11] = n11;
76747c6ae99SBarry Smith   dd->neighbors[12] = n12;
76847c6ae99SBarry Smith   dd->neighbors[13] = rank;
76947c6ae99SBarry Smith   dd->neighbors[14] = n14;
77047c6ae99SBarry Smith   dd->neighbors[15] = n15;
77147c6ae99SBarry Smith   dd->neighbors[16] = n16;
77247c6ae99SBarry Smith   dd->neighbors[17] = n17;
77347c6ae99SBarry Smith   dd->neighbors[18] = n18;
77447c6ae99SBarry Smith   dd->neighbors[19] = n19;
77547c6ae99SBarry Smith   dd->neighbors[20] = n20;
77647c6ae99SBarry Smith   dd->neighbors[21] = n21;
77747c6ae99SBarry Smith   dd->neighbors[22] = n22;
77847c6ae99SBarry Smith   dd->neighbors[23] = n23;
77947c6ae99SBarry Smith   dd->neighbors[24] = n24;
78047c6ae99SBarry Smith   dd->neighbors[25] = n25;
78147c6ae99SBarry Smith   dd->neighbors[26] = n26;
78247c6ae99SBarry Smith 
78347c6ae99SBarry Smith   /* If star stencil then delete the corner neighbors */
784d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
78547c6ae99SBarry Smith     /* save information about corner neighbors */
7869371c9d4SSatish Balay     sn0  = n0;
7879371c9d4SSatish Balay     sn1  = n1;
7889371c9d4SSatish Balay     sn2  = n2;
7899371c9d4SSatish Balay     sn3  = n3;
7909371c9d4SSatish Balay     sn5  = n5;
7919371c9d4SSatish Balay     sn6  = n6;
7929371c9d4SSatish Balay     sn7  = n7;
7939371c9d4SSatish Balay     sn8  = n8;
7949371c9d4SSatish Balay     sn9  = n9;
7959371c9d4SSatish Balay     sn11 = n11;
7969371c9d4SSatish Balay     sn15 = n15;
7979371c9d4SSatish Balay     sn17 = n17;
7989371c9d4SSatish Balay     sn18 = n18;
7999371c9d4SSatish Balay     sn19 = n19;
8009371c9d4SSatish Balay     sn20 = n20;
8019371c9d4SSatish Balay     sn21 = n21;
8029371c9d4SSatish Balay     sn23 = n23;
8039371c9d4SSatish Balay     sn24 = n24;
8049371c9d4SSatish Balay     sn25 = n25;
80547c6ae99SBarry Smith     sn26 = n26;
8068865f1eaSKarl Rupp     n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
80747c6ae99SBarry Smith   }
80847c6ae99SBarry Smith 
8099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((Xe - Xs) * (Ye - Ys) * (Ze - Zs), &idx));
81047c6ae99SBarry Smith 
81147c6ae99SBarry Smith   nn = 0;
81247c6ae99SBarry Smith   /* Bottom Level */
81347c6ae99SBarry Smith   for (k = 0; k < s_z; k++) {
81447c6ae99SBarry Smith     for (i = 1; i <= s_y; i++) {
81547c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
816ce00eea3SSatish Balay         x_t = lx[n0 % m];
81747c6ae99SBarry Smith         y_t = ly[(n0 % (m * n)) / m];
81847c6ae99SBarry Smith         z_t = lz[n0 / (m * n)];
81947c6ae99SBarry 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;
8208865f1eaSKarl 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 */
8218865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
82247c6ae99SBarry Smith       }
82347c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
82447c6ae99SBarry Smith         x_t = x;
82547c6ae99SBarry Smith         y_t = ly[(n1 % (m * n)) / m];
82647c6ae99SBarry Smith         z_t = lz[n1 / (m * n)];
82747c6ae99SBarry 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;
8288865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n1] + x_t * y_t * z_t - (s_y + 1 - i) * x_t; /* 2D case */
8298865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
83047c6ae99SBarry Smith       }
83147c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
832ce00eea3SSatish Balay         x_t = lx[n2 % m];
83347c6ae99SBarry Smith         y_t = ly[(n2 % (m * n)) / m];
83447c6ae99SBarry Smith         z_t = lz[n2 / (m * n)];
83547c6ae99SBarry 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;
8368865f1eaSKarl Rupp         if (twod && (s_t < 0)) s_t = bases[n2] + x_t * y_t * z_t - (s_y + 1 - i) * x_t; /* 2D case */
8378865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
83847c6ae99SBarry Smith       }
83947c6ae99SBarry Smith     }
84047c6ae99SBarry Smith 
84147c6ae99SBarry Smith     for (i = 0; i < y; i++) {
84247c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
843ce00eea3SSatish Balay         x_t = lx[n3 % m];
84447c6ae99SBarry Smith         y_t = y;
84547c6ae99SBarry Smith         z_t = lz[n3 / (m * n)];
84647c6ae99SBarry Smith         s_t = bases[n3] + (i + 1) * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
8478865f1eaSKarl 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 */
8488865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
84947c6ae99SBarry Smith       }
85047c6ae99SBarry Smith 
85147c6ae99SBarry Smith       if (n4 >= 0) { /* middle */
85247c6ae99SBarry Smith         x_t = x;
85347c6ae99SBarry Smith         y_t = y;
85447c6ae99SBarry Smith         z_t = lz[n4 / (m * n)];
85547c6ae99SBarry Smith         s_t = bases[n4] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
8568865f1eaSKarl 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 */
8578865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
858bff4a2f0SMatthew G. Knepley       } else if (bz == DM_BOUNDARY_MIRROR) {
859a6fd6b0aSBarry Smith         for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + x * i + (s_z - k - 1) * x * y;
86047c6ae99SBarry Smith       }
86147c6ae99SBarry Smith 
86247c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
863ce00eea3SSatish Balay         x_t = lx[n5 % m];
86447c6ae99SBarry Smith         y_t = y;
86547c6ae99SBarry Smith         z_t = lz[n5 / (m * n)];
86647c6ae99SBarry Smith         s_t = bases[n5] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
8678865f1eaSKarl 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 */
8688865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
86947c6ae99SBarry Smith       }
87047c6ae99SBarry Smith     }
87147c6ae99SBarry Smith 
87247c6ae99SBarry Smith     for (i = 1; i <= s_y; i++) {
87347c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
874ce00eea3SSatish Balay         x_t = lx[n6 % m];
87547c6ae99SBarry Smith         y_t = ly[(n6 % (m * n)) / m];
87647c6ae99SBarry Smith         z_t = lz[n6 / (m * n)];
87747c6ae99SBarry Smith         s_t = bases[n6] + i * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
8788865f1eaSKarl 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 */
8798865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
88047c6ae99SBarry Smith       }
88147c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
88247c6ae99SBarry Smith         x_t = x;
88347c6ae99SBarry Smith         y_t = ly[(n7 % (m * n)) / m];
88447c6ae99SBarry Smith         z_t = lz[n7 / (m * n)];
88547c6ae99SBarry Smith         s_t = bases[n7] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
8868865f1eaSKarl 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 */
8878865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
88847c6ae99SBarry Smith       }
88947c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
890ce00eea3SSatish Balay         x_t = lx[n8 % m];
89147c6ae99SBarry Smith         y_t = ly[(n8 % (m * n)) / m];
89247c6ae99SBarry Smith         z_t = lz[n8 / (m * n)];
89347c6ae99SBarry Smith         s_t = bases[n8] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
8948865f1eaSKarl 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 */
8958865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
89647c6ae99SBarry Smith       }
89747c6ae99SBarry Smith     }
89847c6ae99SBarry Smith   }
89947c6ae99SBarry Smith 
90047c6ae99SBarry Smith   /* Middle Level */
90147c6ae99SBarry Smith   for (k = 0; k < z; k++) {
90247c6ae99SBarry Smith     for (i = 1; i <= s_y; i++) {
90347c6ae99SBarry Smith       if (n9 >= 0) { /* left below */
904ce00eea3SSatish Balay         x_t = lx[n9 % m];
90547c6ae99SBarry Smith         y_t = ly[(n9 % (m * n)) / m];
90647c6ae99SBarry Smith         /* z_t = z; */
90747c6ae99SBarry Smith         s_t = bases[n9] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t;
9088865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
90947c6ae99SBarry Smith       }
91047c6ae99SBarry Smith       if (n10 >= 0) { /* directly below */
91147c6ae99SBarry Smith         x_t = x;
91247c6ae99SBarry Smith         y_t = ly[(n10 % (m * n)) / m];
91347c6ae99SBarry Smith         /* z_t = z; */
91447c6ae99SBarry Smith         s_t = bases[n10] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
9158865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
916bff4a2f0SMatthew G. Knepley       } else if (by == DM_BOUNDARY_MIRROR) {
917a6fd6b0aSBarry Smith         for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (s_y - i) * x;
91847c6ae99SBarry Smith       }
91947c6ae99SBarry Smith       if (n11 >= 0) { /* right below */
920ce00eea3SSatish Balay         x_t = lx[n11 % m];
92147c6ae99SBarry Smith         y_t = ly[(n11 % (m * n)) / m];
92247c6ae99SBarry Smith         /* z_t = z; */
92347c6ae99SBarry Smith         s_t = bases[n11] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
9248865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
92547c6ae99SBarry Smith       }
92647c6ae99SBarry Smith     }
92747c6ae99SBarry Smith 
92847c6ae99SBarry Smith     for (i = 0; i < y; i++) {
92947c6ae99SBarry Smith       if (n12 >= 0) { /* directly left */
930ce00eea3SSatish Balay         x_t = lx[n12 % m];
93147c6ae99SBarry Smith         y_t = y;
93247c6ae99SBarry Smith         /* z_t = z; */
93347c6ae99SBarry Smith         s_t = bases[n12] + (i + 1) * x_t - s_x + k * x_t * y_t;
9348865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
935bff4a2f0SMatthew G. Knepley       } else if (bx == DM_BOUNDARY_MIRROR) {
936a6fd6b0aSBarry Smith         for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k * x * y + i * x;
93747c6ae99SBarry Smith       }
93847c6ae99SBarry Smith 
93947c6ae99SBarry Smith       /* Interior */
94047c6ae99SBarry Smith       s_t = bases[rank] + i * x + k * x * y;
9418865f1eaSKarl Rupp       for (j = 0; j < x; j++) idx[nn++] = s_t++;
94247c6ae99SBarry Smith 
94347c6ae99SBarry Smith       if (n14 >= 0) { /* directly right */
944ce00eea3SSatish Balay         x_t = lx[n14 % m];
94547c6ae99SBarry Smith         y_t = y;
94647c6ae99SBarry Smith         /* z_t = z; */
94747c6ae99SBarry Smith         s_t = bases[n14] + i * x_t + k * x_t * y_t;
9488865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
949bff4a2f0SMatthew G. Knepley       } else if (bx == DM_BOUNDARY_MIRROR) {
950a6fd6b0aSBarry Smith         for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k * x * y + i * x;
95147c6ae99SBarry Smith       }
95247c6ae99SBarry Smith     }
95347c6ae99SBarry Smith 
95447c6ae99SBarry Smith     for (i = 1; i <= s_y; i++) {
95547c6ae99SBarry Smith       if (n15 >= 0) { /* left above */
956ce00eea3SSatish Balay         x_t = lx[n15 % m];
95747c6ae99SBarry Smith         y_t = ly[(n15 % (m * n)) / m];
95847c6ae99SBarry Smith         /* z_t = z; */
95947c6ae99SBarry Smith         s_t = bases[n15] + i * x_t - s_x + k * x_t * y_t;
9608865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
96147c6ae99SBarry Smith       }
96247c6ae99SBarry Smith       if (n16 >= 0) { /* directly above */
96347c6ae99SBarry Smith         x_t = x;
96447c6ae99SBarry Smith         y_t = ly[(n16 % (m * n)) / m];
96547c6ae99SBarry Smith         /* z_t = z; */
96647c6ae99SBarry Smith         s_t = bases[n16] + (i - 1) * x_t + k * x_t * y_t;
9678865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
968bff4a2f0SMatthew G. Knepley       } else if (by == DM_BOUNDARY_MIRROR) {
969a6fd6b0aSBarry Smith         for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (y - i) * x;
97047c6ae99SBarry Smith       }
97147c6ae99SBarry Smith       if (n17 >= 0) { /* right above */
972ce00eea3SSatish Balay         x_t = lx[n17 % m];
97347c6ae99SBarry Smith         y_t = ly[(n17 % (m * n)) / m];
97447c6ae99SBarry Smith         /* z_t = z; */
97547c6ae99SBarry Smith         s_t = bases[n17] + (i - 1) * x_t + k * x_t * y_t;
9768865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
97747c6ae99SBarry Smith       }
97847c6ae99SBarry Smith     }
97947c6ae99SBarry Smith   }
98047c6ae99SBarry Smith 
98147c6ae99SBarry Smith   /* Upper Level */
98247c6ae99SBarry Smith   for (k = 0; k < s_z; k++) {
98347c6ae99SBarry Smith     for (i = 1; i <= s_y; i++) {
98447c6ae99SBarry Smith       if (n18 >= 0) { /* left below */
985ce00eea3SSatish Balay         x_t = lx[n18 % m];
98647c6ae99SBarry Smith         y_t = ly[(n18 % (m * n)) / m];
98747c6ae99SBarry Smith         /* z_t = lz[n18 / (m*n)]; */
98847c6ae99SBarry Smith         s_t = bases[n18] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t;
9898865f1eaSKarl 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 */
9908865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
99147c6ae99SBarry Smith       }
99247c6ae99SBarry Smith       if (n19 >= 0) { /* directly below */
99347c6ae99SBarry Smith         x_t = x;
99447c6ae99SBarry Smith         y_t = ly[(n19 % (m * n)) / m];
99547c6ae99SBarry Smith         /* z_t = lz[n19 / (m*n)]; */
99647c6ae99SBarry Smith         s_t = bases[n19] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
9978865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n19] - (s_y + 1 - i) * x_t + x_t * y_t; /* 2d case */
9988865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
99947c6ae99SBarry Smith       }
100047c6ae99SBarry Smith       if (n20 >= 0) { /* right below */
1001ce00eea3SSatish Balay         x_t = lx[n20 % m];
100247c6ae99SBarry Smith         y_t = ly[(n20 % (m * n)) / m];
100347c6ae99SBarry Smith         /* z_t = lz[n20 / (m*n)]; */
100447c6ae99SBarry Smith         s_t = bases[n20] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
10058865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n20] - (s_y + 1 - i) * x_t + x_t * y_t; /* 2d case */
10068865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
100747c6ae99SBarry Smith       }
100847c6ae99SBarry Smith     }
100947c6ae99SBarry Smith 
101047c6ae99SBarry Smith     for (i = 0; i < y; i++) {
101147c6ae99SBarry Smith       if (n21 >= 0) { /* directly left */
1012ce00eea3SSatish Balay         x_t = lx[n21 % m];
101347c6ae99SBarry Smith         y_t = y;
101447c6ae99SBarry Smith         /* z_t = lz[n21 / (m*n)]; */
101547c6ae99SBarry Smith         s_t = bases[n21] + (i + 1) * x_t - s_x + k * x_t * y_t;
10168865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n21] + (i + 1) * x_t - s_x; /* 2d case */
10178865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
101847c6ae99SBarry Smith       }
101947c6ae99SBarry Smith 
102047c6ae99SBarry Smith       if (n22 >= 0) { /* middle */
102147c6ae99SBarry Smith         x_t = x;
102247c6ae99SBarry Smith         y_t = y;
102347c6ae99SBarry Smith         /* z_t = lz[n22 / (m*n)]; */
102447c6ae99SBarry Smith         s_t = bases[n22] + i * x_t + k * x_t * y_t;
10258865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n22] + i * x_t; /* 2d case */
10268865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1027bff4a2f0SMatthew G. Knepley       } else if (bz == DM_BOUNDARY_MIRROR) {
1028a6fd6b0aSBarry Smith         for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + (z - k - 1) * x * y + i * x;
102947c6ae99SBarry Smith       }
103047c6ae99SBarry Smith 
103147c6ae99SBarry Smith       if (n23 >= 0) { /* directly right */
1032ce00eea3SSatish Balay         x_t = lx[n23 % m];
103347c6ae99SBarry Smith         y_t = y;
103447c6ae99SBarry Smith         /* z_t = lz[n23 / (m*n)]; */
103547c6ae99SBarry Smith         s_t = bases[n23] + i * x_t + k * x_t * y_t;
10368865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n23] + i * x_t; /* 2d case */
10378865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
103847c6ae99SBarry Smith       }
103947c6ae99SBarry Smith     }
104047c6ae99SBarry Smith 
104147c6ae99SBarry Smith     for (i = 1; i <= s_y; i++) {
104247c6ae99SBarry Smith       if (n24 >= 0) { /* left above */
1043ce00eea3SSatish Balay         x_t = lx[n24 % m];
104447c6ae99SBarry Smith         y_t = ly[(n24 % (m * n)) / m];
104547c6ae99SBarry Smith         /* z_t = lz[n24 / (m*n)]; */
104647c6ae99SBarry Smith         s_t = bases[n24] + i * x_t - s_x + k * x_t * y_t;
10478865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n24] + i * x_t - s_x; /* 2d case */
10488865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
104947c6ae99SBarry Smith       }
105047c6ae99SBarry Smith       if (n25 >= 0) { /* directly above */
105147c6ae99SBarry Smith         x_t = x;
105247c6ae99SBarry Smith         y_t = ly[(n25 % (m * n)) / m];
105347c6ae99SBarry Smith         /* z_t = lz[n25 / (m*n)]; */
105447c6ae99SBarry Smith         s_t = bases[n25] + (i - 1) * x_t + k * x_t * y_t;
10558865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n25] + (i - 1) * x_t; /* 2d case */
10568865f1eaSKarl Rupp         for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
105747c6ae99SBarry Smith       }
105847c6ae99SBarry Smith       if (n26 >= 0) { /* right above */
1059ce00eea3SSatish Balay         x_t = lx[n26 % m];
106047c6ae99SBarry Smith         y_t = ly[(n26 % (m * n)) / m];
106147c6ae99SBarry Smith         /* z_t = lz[n26 / (m*n)]; */
106247c6ae99SBarry Smith         s_t = bases[n26] + (i - 1) * x_t + k * x_t * y_t;
10638865f1eaSKarl Rupp         if (twod && (s_t >= M * N * P)) s_t = bases[n26] + (i - 1) * x_t; /* 2d case */
10648865f1eaSKarl Rupp         for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
106547c6ae99SBarry Smith       }
106647c6ae99SBarry Smith     }
106747c6ae99SBarry Smith   }
1068ce00eea3SSatish Balay 
10699566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(comm, dof, nn, idx, PETSC_USE_POINTER, &from));
10709566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(global, from, local, to, &gtol));
10719566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
10729566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
107347c6ae99SBarry Smith 
1074d9c9ebe5SPeter Brune   if (stencil_type == DMDA_STENCIL_STAR) {
10759371c9d4SSatish Balay     n0  = sn0;
10769371c9d4SSatish Balay     n1  = sn1;
10779371c9d4SSatish Balay     n2  = sn2;
10789371c9d4SSatish Balay     n3  = sn3;
10799371c9d4SSatish Balay     n5  = sn5;
10809371c9d4SSatish Balay     n6  = sn6;
10819371c9d4SSatish Balay     n7  = sn7;
10829371c9d4SSatish Balay     n8  = sn8;
10839371c9d4SSatish Balay     n9  = sn9;
10849371c9d4SSatish Balay     n11 = sn11;
10859371c9d4SSatish Balay     n15 = sn15;
10869371c9d4SSatish Balay     n17 = sn17;
10879371c9d4SSatish Balay     n18 = sn18;
10889371c9d4SSatish Balay     n19 = sn19;
10899371c9d4SSatish Balay     n20 = sn20;
10909371c9d4SSatish Balay     n21 = sn21;
10919371c9d4SSatish Balay     n23 = sn23;
10929371c9d4SSatish Balay     n24 = sn24;
10939371c9d4SSatish Balay     n25 = sn25;
109447c6ae99SBarry Smith     n26 = sn26;
1095ce00eea3SSatish Balay   }
109647c6ae99SBarry Smith 
1097f4f49eeaSPierre Jolivet   if ((stencil_type == DMDA_STENCIL_STAR) || (bx != DM_BOUNDARY_PERIODIC && bx) || (by != DM_BOUNDARY_PERIODIC && by) || (bz != DM_BOUNDARY_PERIODIC && bz)) {
1098ce00eea3SSatish Balay     /*
1099ce00eea3SSatish Balay         Recompute the local to global mappings, this time keeping the
1100ce00eea3SSatish Balay       information about the cross corner processor numbers.
1101ce00eea3SSatish Balay     */
110247c6ae99SBarry Smith     nn = 0;
110347c6ae99SBarry Smith     /* Bottom Level */
110447c6ae99SBarry Smith     for (k = 0; k < s_z; k++) {
110547c6ae99SBarry Smith       for (i = 1; i <= s_y; i++) {
110647c6ae99SBarry Smith         if (n0 >= 0) { /* left below */
1107ce00eea3SSatish Balay           x_t = lx[n0 % m];
110847c6ae99SBarry Smith           y_t = ly[(n0 % (m * n)) / m];
110947c6ae99SBarry Smith           z_t = lz[n0 / (m * n)];
111047c6ae99SBarry 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;
11118865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1112ce00eea3SSatish Balay         } else if (Xs - xs < 0 && Ys - ys < 0 && Zs - zs < 0) {
11138865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
111447c6ae99SBarry Smith         }
111547c6ae99SBarry Smith         if (n1 >= 0) { /* directly below */
111647c6ae99SBarry Smith           x_t = x;
111747c6ae99SBarry Smith           y_t = ly[(n1 % (m * n)) / m];
111847c6ae99SBarry Smith           z_t = lz[n1 / (m * n)];
111947c6ae99SBarry 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;
11208865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1121ce00eea3SSatish Balay         } else if (Ys - ys < 0 && Zs - zs < 0) {
11228865f1eaSKarl Rupp           for (j = 0; j < x; j++) idx[nn++] = -1;
112347c6ae99SBarry Smith         }
112447c6ae99SBarry Smith         if (n2 >= 0) { /* right below */
1125ce00eea3SSatish Balay           x_t = lx[n2 % m];
112647c6ae99SBarry Smith           y_t = ly[(n2 % (m * n)) / m];
112747c6ae99SBarry Smith           z_t = lz[n2 / (m * n)];
112847c6ae99SBarry 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;
11298865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1130ce00eea3SSatish Balay         } else if (xe - Xe < 0 && Ys - ys < 0 && Zs - zs < 0) {
11318865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
113247c6ae99SBarry Smith         }
113347c6ae99SBarry Smith       }
113447c6ae99SBarry Smith 
113547c6ae99SBarry Smith       for (i = 0; i < y; i++) {
113647c6ae99SBarry Smith         if (n3 >= 0) { /* directly left */
1137ce00eea3SSatish Balay           x_t = lx[n3 % m];
113847c6ae99SBarry Smith           y_t = y;
113947c6ae99SBarry Smith           z_t = lz[n3 / (m * n)];
114047c6ae99SBarry Smith           s_t = bases[n3] + (i + 1) * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
11418865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1142ce00eea3SSatish Balay         } else if (Xs - xs < 0 && Zs - zs < 0) {
11438865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
114447c6ae99SBarry Smith         }
114547c6ae99SBarry Smith 
114647c6ae99SBarry Smith         if (n4 >= 0) { /* middle */
114747c6ae99SBarry Smith           x_t = x;
114847c6ae99SBarry Smith           y_t = y;
114947c6ae99SBarry Smith           z_t = lz[n4 / (m * n)];
115047c6ae99SBarry Smith           s_t = bases[n4] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
11518865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1152ce00eea3SSatish Balay         } else if (Zs - zs < 0) {
1153bff4a2f0SMatthew G. Knepley           if (bz == DM_BOUNDARY_MIRROR) {
115418f812f4SBarry Smith             for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + x * i + (s_z - k) * x * y;
1155c2859e5eSBarry Smith           } else {
11568865f1eaSKarl Rupp             for (j = 0; j < x; j++) idx[nn++] = -1;
115747c6ae99SBarry Smith           }
1158c2859e5eSBarry Smith         }
115947c6ae99SBarry Smith 
116047c6ae99SBarry Smith         if (n5 >= 0) { /* directly right */
1161ce00eea3SSatish Balay           x_t = lx[n5 % m];
116247c6ae99SBarry Smith           y_t = y;
116347c6ae99SBarry Smith           z_t = lz[n5 / (m * n)];
116447c6ae99SBarry Smith           s_t = bases[n5] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
11658865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1166ce00eea3SSatish Balay         } else if (xe - Xe < 0 && Zs - zs < 0) {
11678865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
116847c6ae99SBarry Smith         }
116947c6ae99SBarry Smith       }
117047c6ae99SBarry Smith 
117147c6ae99SBarry Smith       for (i = 1; i <= s_y; i++) {
117247c6ae99SBarry Smith         if (n6 >= 0) { /* left above */
1173ce00eea3SSatish Balay           x_t = lx[n6 % m];
117447c6ae99SBarry Smith           y_t = ly[(n6 % (m * n)) / m];
117547c6ae99SBarry Smith           z_t = lz[n6 / (m * n)];
117647c6ae99SBarry Smith           s_t = bases[n6] + i * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
11778865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1178ce00eea3SSatish Balay         } else if (Xs - xs < 0 && ye - Ye < 0 && Zs - zs < 0) {
11798865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
118047c6ae99SBarry Smith         }
118147c6ae99SBarry Smith         if (n7 >= 0) { /* directly above */
118247c6ae99SBarry Smith           x_t = x;
118347c6ae99SBarry Smith           y_t = ly[(n7 % (m * n)) / m];
118447c6ae99SBarry Smith           z_t = lz[n7 / (m * n)];
118547c6ae99SBarry Smith           s_t = bases[n7] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
11868865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1187ce00eea3SSatish Balay         } else if (ye - Ye < 0 && Zs - zs < 0) {
11888865f1eaSKarl Rupp           for (j = 0; j < x; j++) idx[nn++] = -1;
118947c6ae99SBarry Smith         }
119047c6ae99SBarry Smith         if (n8 >= 0) { /* right above */
1191ce00eea3SSatish Balay           x_t = lx[n8 % m];
119247c6ae99SBarry Smith           y_t = ly[(n8 % (m * n)) / m];
119347c6ae99SBarry Smith           z_t = lz[n8 / (m * n)];
119447c6ae99SBarry Smith           s_t = bases[n8] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t;
11958865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1196ce00eea3SSatish Balay         } else if (xe - Xe < 0 && ye - Ye < 0 && Zs - zs < 0) {
11978865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
119847c6ae99SBarry Smith         }
119947c6ae99SBarry Smith       }
120047c6ae99SBarry Smith     }
120147c6ae99SBarry Smith 
120247c6ae99SBarry Smith     /* Middle Level */
120347c6ae99SBarry Smith     for (k = 0; k < z; k++) {
120447c6ae99SBarry Smith       for (i = 1; i <= s_y; i++) {
120547c6ae99SBarry Smith         if (n9 >= 0) { /* left below */
1206ce00eea3SSatish Balay           x_t = lx[n9 % m];
120747c6ae99SBarry Smith           y_t = ly[(n9 % (m * n)) / m];
120847c6ae99SBarry Smith           /* z_t = z; */
120947c6ae99SBarry Smith           s_t = bases[n9] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t;
12108865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1211ce00eea3SSatish Balay         } else if (Xs - xs < 0 && Ys - ys < 0) {
12128865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
121347c6ae99SBarry Smith         }
121447c6ae99SBarry Smith         if (n10 >= 0) { /* directly below */
121547c6ae99SBarry Smith           x_t = x;
121647c6ae99SBarry Smith           y_t = ly[(n10 % (m * n)) / m];
121747c6ae99SBarry Smith           /* z_t = z; */
121847c6ae99SBarry Smith           s_t = bases[n10] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
12198865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1220ce00eea3SSatish Balay         } else if (Ys - ys < 0) {
1221bff4a2f0SMatthew G. Knepley           if (by == DM_BOUNDARY_MIRROR) {
122218f812f4SBarry Smith             for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (s_y - i + 1) * x;
1223c2859e5eSBarry Smith           } else {
12248865f1eaSKarl Rupp             for (j = 0; j < x; j++) idx[nn++] = -1;
1225c2859e5eSBarry Smith           }
122647c6ae99SBarry Smith         }
122747c6ae99SBarry Smith         if (n11 >= 0) { /* right below */
1228ce00eea3SSatish Balay           x_t = lx[n11 % m];
122947c6ae99SBarry Smith           y_t = ly[(n11 % (m * n)) / m];
123047c6ae99SBarry Smith           /* z_t = z; */
123147c6ae99SBarry Smith           s_t = bases[n11] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
12328865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1233ce00eea3SSatish Balay         } else if (xe - Xe < 0 && Ys - ys < 0) {
12348865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
123547c6ae99SBarry Smith         }
123647c6ae99SBarry Smith       }
123747c6ae99SBarry Smith 
123847c6ae99SBarry Smith       for (i = 0; i < y; i++) {
123947c6ae99SBarry Smith         if (n12 >= 0) { /* directly left */
1240ce00eea3SSatish Balay           x_t = lx[n12 % m];
124147c6ae99SBarry Smith           y_t = y;
124247c6ae99SBarry Smith           /* z_t = z; */
124347c6ae99SBarry Smith           s_t = bases[n12] + (i + 1) * x_t - s_x + k * x_t * y_t;
12448865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1245ce00eea3SSatish Balay         } else if (Xs - xs < 0) {
1246bff4a2f0SMatthew G. Knepley           if (bx == DM_BOUNDARY_MIRROR) {
124718f812f4SBarry Smith             for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k * x * y + i * x + 1;
1248c2859e5eSBarry Smith           } else {
12498865f1eaSKarl Rupp             for (j = 0; j < s_x; j++) idx[nn++] = -1;
125047c6ae99SBarry Smith           }
1251c2859e5eSBarry Smith         }
125247c6ae99SBarry Smith 
125347c6ae99SBarry Smith         /* Interior */
125447c6ae99SBarry Smith         s_t = bases[rank] + i * x + k * x * y;
12558865f1eaSKarl Rupp         for (j = 0; j < x; j++) idx[nn++] = s_t++;
125647c6ae99SBarry Smith 
125747c6ae99SBarry Smith         if (n14 >= 0) { /* directly right */
1258ce00eea3SSatish Balay           x_t = lx[n14 % m];
125947c6ae99SBarry Smith           y_t = y;
126047c6ae99SBarry Smith           /* z_t = z; */
126147c6ae99SBarry Smith           s_t = bases[n14] + i * x_t + k * x_t * y_t;
12628865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1263ce00eea3SSatish Balay         } else if (xe - Xe < 0) {
1264bff4a2f0SMatthew G. Knepley           if (bx == DM_BOUNDARY_MIRROR) {
126518f812f4SBarry Smith             for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k * x * y + i * x - 1;
1266c2859e5eSBarry Smith           } else {
12678865f1eaSKarl Rupp             for (j = 0; j < s_x; j++) idx[nn++] = -1;
126847c6ae99SBarry Smith           }
126947c6ae99SBarry Smith         }
1270c2859e5eSBarry Smith       }
127147c6ae99SBarry Smith 
127247c6ae99SBarry Smith       for (i = 1; i <= s_y; i++) {
127347c6ae99SBarry Smith         if (n15 >= 0) { /* left above */
1274ce00eea3SSatish Balay           x_t = lx[n15 % m];
127547c6ae99SBarry Smith           y_t = ly[(n15 % (m * n)) / m];
127647c6ae99SBarry Smith           /* z_t = z; */
127747c6ae99SBarry Smith           s_t = bases[n15] + i * x_t - s_x + k * x_t * y_t;
12788865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1279ce00eea3SSatish Balay         } else if (Xs - xs < 0 && ye - Ye < 0) {
12808865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
128147c6ae99SBarry Smith         }
128247c6ae99SBarry Smith         if (n16 >= 0) { /* directly above */
128347c6ae99SBarry Smith           x_t = x;
128447c6ae99SBarry Smith           y_t = ly[(n16 % (m * n)) / m];
128547c6ae99SBarry Smith           /* z_t = z; */
128647c6ae99SBarry Smith           s_t = bases[n16] + (i - 1) * x_t + k * x_t * y_t;
12878865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1288ce00eea3SSatish Balay         } else if (ye - Ye < 0) {
1289bff4a2f0SMatthew G. Knepley           if (by == DM_BOUNDARY_MIRROR) {
129018f812f4SBarry Smith             for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (y - i - 1) * x;
1291c2859e5eSBarry Smith           } else {
12928865f1eaSKarl Rupp             for (j = 0; j < x; j++) idx[nn++] = -1;
129347c6ae99SBarry Smith           }
1294c2859e5eSBarry Smith         }
129547c6ae99SBarry Smith         if (n17 >= 0) { /* right above */
1296ce00eea3SSatish Balay           x_t = lx[n17 % m];
129747c6ae99SBarry Smith           y_t = ly[(n17 % (m * n)) / m];
129847c6ae99SBarry Smith           /* z_t = z; */
129947c6ae99SBarry Smith           s_t = bases[n17] + (i - 1) * x_t + k * x_t * y_t;
13008865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1301ce00eea3SSatish Balay         } else if (xe - Xe < 0 && ye - Ye < 0) {
13028865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
130347c6ae99SBarry Smith         }
130447c6ae99SBarry Smith       }
130547c6ae99SBarry Smith     }
130647c6ae99SBarry Smith 
130747c6ae99SBarry Smith     /* Upper Level */
130847c6ae99SBarry Smith     for (k = 0; k < s_z; k++) {
130947c6ae99SBarry Smith       for (i = 1; i <= s_y; i++) {
131047c6ae99SBarry Smith         if (n18 >= 0) { /* left below */
1311ce00eea3SSatish Balay           x_t = lx[n18 % m];
131247c6ae99SBarry Smith           y_t = ly[(n18 % (m * n)) / m];
131347c6ae99SBarry Smith           /* z_t = lz[n18 / (m*n)]; */
131447c6ae99SBarry Smith           s_t = bases[n18] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t;
13158865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1316ce00eea3SSatish Balay         } else if (Xs - xs < 0 && Ys - ys < 0 && ze - Ze < 0) {
13178865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
131847c6ae99SBarry Smith         }
131947c6ae99SBarry Smith         if (n19 >= 0) { /* directly below */
132047c6ae99SBarry Smith           x_t = x;
132147c6ae99SBarry Smith           y_t = ly[(n19 % (m * n)) / m];
132247c6ae99SBarry Smith           /* z_t = lz[n19 / (m*n)]; */
132347c6ae99SBarry Smith           s_t = bases[n19] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
13248865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1325ce00eea3SSatish Balay         } else if (Ys - ys < 0 && ze - Ze < 0) {
13268865f1eaSKarl Rupp           for (j = 0; j < x; j++) idx[nn++] = -1;
132747c6ae99SBarry Smith         }
132847c6ae99SBarry Smith         if (n20 >= 0) { /* right below */
1329ce00eea3SSatish Balay           x_t = lx[n20 % m];
133047c6ae99SBarry Smith           y_t = ly[(n20 % (m * n)) / m];
133147c6ae99SBarry Smith           /* z_t = lz[n20 / (m*n)]; */
133247c6ae99SBarry Smith           s_t = bases[n20] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t;
13338865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1334ce00eea3SSatish Balay         } else if (xe - Xe < 0 && Ys - ys < 0 && ze - Ze < 0) {
13358865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
133647c6ae99SBarry Smith         }
133747c6ae99SBarry Smith       }
133847c6ae99SBarry Smith 
133947c6ae99SBarry Smith       for (i = 0; i < y; i++) {
134047c6ae99SBarry Smith         if (n21 >= 0) { /* directly left */
1341ce00eea3SSatish Balay           x_t = lx[n21 % m];
134247c6ae99SBarry Smith           y_t = y;
134347c6ae99SBarry Smith           /* z_t = lz[n21 / (m*n)]; */
134447c6ae99SBarry Smith           s_t = bases[n21] + (i + 1) * x_t - s_x + k * x_t * y_t;
13458865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1346ce00eea3SSatish Balay         } else if (Xs - xs < 0 && ze - Ze < 0) {
13478865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
134847c6ae99SBarry Smith         }
134947c6ae99SBarry Smith 
135047c6ae99SBarry Smith         if (n22 >= 0) { /* middle */
135147c6ae99SBarry Smith           x_t = x;
135247c6ae99SBarry Smith           y_t = y;
135347c6ae99SBarry Smith           /* z_t = lz[n22 / (m*n)]; */
135447c6ae99SBarry Smith           s_t = bases[n22] + i * x_t + k * x_t * y_t;
13558865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1356ce00eea3SSatish Balay         } else if (ze - Ze < 0) {
1357bff4a2f0SMatthew G. Knepley           if (bz == DM_BOUNDARY_MIRROR) {
135818f812f4SBarry Smith             for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + (z - k - 2) * x * y + i * x;
1359c2859e5eSBarry Smith           } else {
13608865f1eaSKarl Rupp             for (j = 0; j < x; j++) idx[nn++] = -1;
136147c6ae99SBarry Smith           }
1362c2859e5eSBarry Smith         }
136347c6ae99SBarry Smith 
136447c6ae99SBarry Smith         if (n23 >= 0) { /* directly right */
1365ce00eea3SSatish Balay           x_t = lx[n23 % m];
136647c6ae99SBarry Smith           y_t = y;
136747c6ae99SBarry Smith           /* z_t = lz[n23 / (m*n)]; */
136847c6ae99SBarry Smith           s_t = bases[n23] + i * x_t + k * x_t * y_t;
13698865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1370ce00eea3SSatish Balay         } else if (xe - Xe < 0 && ze - Ze < 0) {
13718865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
137247c6ae99SBarry Smith         }
137347c6ae99SBarry Smith       }
137447c6ae99SBarry Smith 
137547c6ae99SBarry Smith       for (i = 1; i <= s_y; i++) {
137647c6ae99SBarry Smith         if (n24 >= 0) { /* left above */
1377ce00eea3SSatish Balay           x_t = lx[n24 % m];
137847c6ae99SBarry Smith           y_t = ly[(n24 % (m * n)) / m];
137947c6ae99SBarry Smith           /* z_t = lz[n24 / (m*n)]; */
138047c6ae99SBarry Smith           s_t = bases[n24] + i * x_t - s_x + k * x_t * y_t;
13818865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1382ce00eea3SSatish Balay         } else if (Xs - xs < 0 && ye - Ye < 0 && ze - Ze < 0) {
13838865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
138447c6ae99SBarry Smith         }
138547c6ae99SBarry Smith         if (n25 >= 0) { /* directly above */
138647c6ae99SBarry Smith           x_t = x;
138747c6ae99SBarry Smith           y_t = ly[(n25 % (m * n)) / m];
138847c6ae99SBarry Smith           /* z_t = lz[n25 / (m*n)]; */
138947c6ae99SBarry Smith           s_t = bases[n25] + (i - 1) * x_t + k * x_t * y_t;
13908865f1eaSKarl Rupp           for (j = 0; j < x_t; j++) idx[nn++] = s_t++;
1391ce00eea3SSatish Balay         } else if (ye - Ye < 0 && ze - Ze < 0) {
13928865f1eaSKarl Rupp           for (j = 0; j < x; j++) idx[nn++] = -1;
139347c6ae99SBarry Smith         }
139447c6ae99SBarry Smith         if (n26 >= 0) { /* right above */
1395ce00eea3SSatish Balay           x_t = lx[n26 % m];
139647c6ae99SBarry Smith           y_t = ly[(n26 % (m * n)) / m];
139747c6ae99SBarry Smith           /* z_t = lz[n26 / (m*n)]; */
139847c6ae99SBarry Smith           s_t = bases[n26] + (i - 1) * x_t + k * x_t * y_t;
13998865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = s_t++;
1400ce00eea3SSatish Balay         } else if (xe - Xe < 0 && ye - Ye < 0 && ze - Ze < 0) {
14018865f1eaSKarl Rupp           for (j = 0; j < s_x; j++) idx[nn++] = -1;
140247c6ae99SBarry Smith         }
140347c6ae99SBarry Smith       }
140447c6ae99SBarry Smith     }
140547c6ae99SBarry Smith   }
140647c6ae99SBarry Smith   /*
140747c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
140847c6ae99SBarry Smith      of VecSetValuesLocal().
140947c6ae99SBarry Smith   */
14109566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm, dof, nn, idx, PETSC_OWN_POINTER, &da->ltogmap));
141147c6ae99SBarry Smith 
14129566063dSJacob Faibussowitsch   PetscCall(PetscFree2(bases, ldims));
14139371c9d4SSatish Balay   dd->m = m;
14149371c9d4SSatish Balay   dd->n = n;
14159371c9d4SSatish Balay   dd->p = p;
1416ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
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;
14239371c9d4SSatish Balay   dd->Xs = Xs * dof;
14249371c9d4SSatish Balay   dd->Xe = Xe * dof;
14259371c9d4SSatish Balay   dd->Ys = Ys;
14269371c9d4SSatish Balay   dd->Ye = Ye;
14279371c9d4SSatish Balay   dd->Zs = Zs;
14289371c9d4SSatish Balay   dd->Ze = Ze;
1429ce00eea3SSatish Balay 
14309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&local));
14319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
1432ce00eea3SSatish Balay 
1433ce00eea3SSatish Balay   dd->gtol      = gtol;
1434ce00eea3SSatish Balay   dd->base      = base;
1435ce00eea3SSatish Balay   da->ops->view = DMView_DA_3d;
14360298fd71SBarry Smith   dd->ltol      = NULL;
14370298fd71SBarry Smith   dd->ao        = NULL;
14383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
143947c6ae99SBarry Smith }
1440564755cdSBarry Smith 
1441cc4c1da9SBarry Smith /*@
1442aa219208SBarry Smith   DMDACreate3d - Creates an object that will manage the communication of three-dimensional
144312b4a537SBarry Smith   regular array data that is distributed across one or more MPI processes.
144447c6ae99SBarry Smith 
1445d083f849SBarry Smith   Collective
144647c6ae99SBarry Smith 
144747c6ae99SBarry Smith   Input Parameters:
144847c6ae99SBarry Smith + comm         - MPI communicator
1449a4e35b19SJacob Faibussowitsch . bx           - type of x ghost nodes the array have.
1450a4e35b19SJacob Faibussowitsch                  Use one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`.
1451a4e35b19SJacob Faibussowitsch . by           - type of y ghost nodes the array have.
1452a4e35b19SJacob Faibussowitsch                  Use one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`.
1453a4e35b19SJacob Faibussowitsch . bz           - type of z ghost nodes the array have.
1454dce8aebaSBarry Smith                  Use one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`.
1455dce8aebaSBarry Smith . stencil_type - Type of stencil (`DMDA_STENCIL_STAR` or `DMDA_STENCIL_BOX`)
1456a4e35b19SJacob Faibussowitsch . M            - global dimension in x direction of the array
1457a4e35b19SJacob Faibussowitsch . N            - global dimension in y direction of the array
1458a4e35b19SJacob Faibussowitsch . P            - global dimension in z direction of the array
1459a4e35b19SJacob Faibussowitsch . m            - corresponding number of processors in x dimension (or `PETSC_DECIDE` to have calculated)
1460a4e35b19SJacob Faibussowitsch . n            - corresponding number of processors in y dimension (or `PETSC_DECIDE` to have calculated)
1461a4e35b19SJacob Faibussowitsch . p            - corresponding number of processors in z dimension (or `PETSC_DECIDE` to have calculated)
146247c6ae99SBarry Smith . dof          - number of degrees of freedom per node
146310d7c030SMatthew G Knepley . s            - stencil width
1464a4e35b19SJacob Faibussowitsch . lx           - arrays containing the number of nodes in each cell along the x  coordinates, or `NULL`.
1465a4e35b19SJacob Faibussowitsch . ly           - arrays containing the number of nodes in each cell along the y coordinates, or `NULL`.
1466a4e35b19SJacob Faibussowitsch - lz           - arrays containing the number of nodes in each cell along the z coordinates, or `NULL`.
146747c6ae99SBarry Smith 
146847c6ae99SBarry Smith   Output Parameter:
146947c6ae99SBarry Smith . da - the resulting distributed array object
147047c6ae99SBarry Smith 
1471dce8aebaSBarry Smith   Options Database Keys:
1472dce8aebaSBarry Smith + -dm_view              - Calls `DMView()` at the conclusion of `DMDACreate3d()`
1473e43dc8daSBarry Smith . -da_grid_x <nx>       - number of grid points in x direction
1474e43dc8daSBarry Smith . -da_grid_y <ny>       - number of grid points in y direction
1475e43dc8daSBarry Smith . -da_grid_z <nz>       - number of grid points in z direction
147647c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction
147747c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction
147847c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction
1479fe58071bSMatthew G. Knepley . -da_bd_x <bx>         - boundary type in x direction
1480fe58071bSMatthew G. Knepley . -da_bd_y <by>         - boundary type in y direction
1481fe58071bSMatthew G. Knepley . -da_bd_z <bz>         - boundary type in x direction
1482fe58071bSMatthew G. Knepley . -da_bd_all <bt>       - boundary type in all directions
1483e0f5d30fSBarry Smith . -da_refine_x <rx>     - refinement ratio in x direction
1484e0f5d30fSBarry Smith . -da_refine_y <ry>     - refinement ratio in y direction
1485e0f5d30fSBarry Smith . -da_refine_z <rz>     - refinement ratio in z directio
148612b4a537SBarry Smith - -da_refine <n>        - refine the `DMDA` n times before creating it
148747c6ae99SBarry Smith 
148847c6ae99SBarry Smith   Level: beginner
148947c6ae99SBarry Smith 
149047c6ae99SBarry Smith   Notes:
1491a4e35b19SJacob Faibussowitsch   If `lx`, `ly`, or `lz` are non-null, these must be of length as `m`, `n`, `p` and the
1492a4e35b19SJacob Faibussowitsch   corresponding `m`, `n`, or `p` cannot be `PETSC_DECIDE`. Sum of the `lx` entries must be `M`,
1493a4e35b19SJacob Faibussowitsch   sum of the `ly` must `N`, sum of the `lz` must be `P`.
1494a4e35b19SJacob Faibussowitsch 
1495dce8aebaSBarry Smith   The stencil type `DMDA_STENCIL_STAR` with width 1 corresponds to the
1496dce8aebaSBarry Smith   standard 7-pt stencil, while `DMDA_STENCIL_BOX` with width 1 denotes
149747c6ae99SBarry Smith   the standard 27-pt stencil.
149847c6ae99SBarry Smith 
1499dce8aebaSBarry Smith   The array data itself is NOT stored in the `DMDA`, it is stored in `Vec` objects;
1500dce8aebaSBarry Smith   The appropriate vector objects can be obtained with calls to `DMCreateGlobalVector()`
1501dce8aebaSBarry Smith   and `DMCreateLocalVector()` and calls to `VecDuplicate()` if more are needed.
150247c6ae99SBarry Smith 
1503dce8aebaSBarry Smith   You must call `DMSetUp()` after this call before using this `DM`.
1504897f7067SBarry Smith 
150512b4a537SBarry Smith   To use the options database to change values in the `DMDA` call `DMSetFromOptions()` after this call
1506dce8aebaSBarry Smith   but before `DMSetUp()`.
1507897f7067SBarry Smith 
150812b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDestroy()`, `DMView()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMGlobalToLocalBegin()`, `DMDAGetRefinementFactor()`,
1509db781477SPatrick Sanan           `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDASetRefinementFactor()`,
1510db781477SPatrick Sanan           `DMDAGetInfo()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMDACreateNaturalVector()`, `DMLoad()`, `DMDAGetOwnershipRanges()`,
151112b4a537SBarry Smith           `DMStagCreate3d()`, `DMBoundaryType`
151247c6ae99SBarry Smith @*/
1513d71ae5a4SJacob Faibussowitsch 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)
1514d71ae5a4SJacob Faibussowitsch {
151547c6ae99SBarry Smith   PetscFunctionBegin;
15169566063dSJacob Faibussowitsch   PetscCall(DMDACreate(comm, da));
15179566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*da, 3));
15189566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(*da, M, N, P));
15199566063dSJacob Faibussowitsch   PetscCall(DMDASetNumProcs(*da, m, n, p));
15209566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(*da, bx, by, bz));
15219566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(*da, dof));
15229566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilType(*da, stencil_type));
15239566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(*da, s));
15249566063dSJacob Faibussowitsch   PetscCall(DMDASetOwnershipRanges(*da, lx, ly, lz));
15253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
152647c6ae99SBarry Smith }
1527