17d0a6c19SBarry Smith 247c6ae99SBarry Smith /* 347c6ae99SBarry Smith Code for manipulating distributed regular 3d arrays in parallel. 447c6ae99SBarry Smith File created by Peter Mell 7/14/95 547c6ae99SBarry Smith */ 647c6ae99SBarry Smith 7af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 847c6ae99SBarry Smith 99804daf3SBarry Smith #include <petscdraw.h> 10d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMView_DA_3d(DM da, PetscViewer viewer) 11d71ae5a4SJacob Faibussowitsch { 1247c6ae99SBarry Smith PetscMPIInt rank; 13c9493c35SStefano Zampini PetscBool iascii, isdraw, isglvis, isbinary; 1447c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 156858538eSMatthew G. Knepley Vec coordinates; 16d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 179a42bb27SBarry Smith PetscBool ismatlab; 189a42bb27SBarry Smith #endif 1947c6ae99SBarry Smith 2047c6ae99SBarry Smith PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)da), &rank)); 2247c6ae99SBarry Smith 239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis)); 269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 27d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERMATLAB, &ismatlab)); 299a42bb27SBarry Smith #endif 3047c6ae99SBarry Smith if (iascii) { 3147c6ae99SBarry Smith PetscViewerFormat format; 3247c6ae99SBarry Smith 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 349566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 3576a8abe0SBarry Smith if (format == PETSC_VIEWER_LOAD_BALANCE) { 3676a8abe0SBarry Smith PetscInt i, nmax = 0, nmin = PETSC_MAX_INT, navg = 0, *nz, nzlocal; 3776a8abe0SBarry Smith DMDALocalInfo info; 3876a8abe0SBarry Smith PetscMPIInt size; 399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da), &size)); 409566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 4176a8abe0SBarry Smith nzlocal = info.xm * info.ym * info.zm; 429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &nz)); 439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nzlocal, 1, MPIU_INT, nz, 1, MPIU_INT, PetscObjectComm((PetscObject)da))); 4476a8abe0SBarry Smith for (i = 0; i < (PetscInt)size; i++) { 4576a8abe0SBarry Smith nmax = PetscMax(nmax, nz[i]); 4676a8abe0SBarry Smith nmin = PetscMin(nmin, nz[i]); 4776a8abe0SBarry Smith navg += nz[i]; 4876a8abe0SBarry Smith } 499566063dSJacob Faibussowitsch PetscCall(PetscFree(nz)); 5076a8abe0SBarry Smith navg = navg / size; 5163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Load Balance - Grid Points: Min %" PetscInt_FMT " avg %" PetscInt_FMT " max %" PetscInt_FMT "\n", nmin, navg, nmax)); 523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5376a8abe0SBarry Smith } 548ec8862eSJed Brown if (format != PETSC_VIEWER_ASCII_VTK_DEPRECATED && format != PETSC_VIEWER_ASCII_VTK_CELL_DEPRECATED && format != PETSC_VIEWER_ASCII_GLVIS) { 55aa219208SBarry Smith DMDALocalInfo info; 569566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 5763a3b9bcSJacob 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)); 589371c9d4SSatish 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, 599371c9d4SSatish Balay info.xs + info.xm, info.ys, info.ys + info.ym, info.zs, info.zs + info.zm)); 606858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(da, &coordinates)); 6147c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 626858538eSMatthew G. Knepley if (coordinates) { 6347c6ae99SBarry Smith PetscInt last; 6447c6ae99SBarry Smith const PetscReal *coors; 656858538eSMatthew G. Knepley PetscCall(VecGetArrayRead(coordinates, &coors)); 666858538eSMatthew G. Knepley PetscCall(VecGetLocalSize(coordinates, &last)); 6747c6ae99SBarry Smith last = last - 3; 689566063dSJacob 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])); 696858538eSMatthew G. Knepley PetscCall(VecRestoreArrayRead(coordinates, &coors)); 7047c6ae99SBarry Smith } 7147c6ae99SBarry Smith #endif 729566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 741baa6e33SBarry Smith } else if (format == PETSC_VIEWER_ASCII_GLVIS) PetscCall(DMView_DA_GLVis(da, viewer)); 751baa6e33SBarry Smith else PetscCall(DMView_DA_VTK(da, viewer)); 7647c6ae99SBarry Smith } else if (isdraw) { 7747c6ae99SBarry Smith PetscDraw draw; 7847c6ae99SBarry Smith PetscReal ymin = -1.0, ymax = (PetscReal)dd->N; 7947c6ae99SBarry Smith PetscReal xmin = -1.0, xmax = (PetscReal)((dd->M + 2) * dd->P), x, y, ycoord, xcoord; 808ea3bf28SBarry Smith PetscInt k, plane, base; 818ea3bf28SBarry Smith const PetscInt *idx; 8247c6ae99SBarry Smith char node[10]; 8347c6ae99SBarry Smith PetscBool isnull; 8447c6ae99SBarry Smith 859566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 869566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw, &isnull)); 873ba16761SJacob Faibussowitsch if (isnull) PetscFunctionReturn(PETSC_SUCCESS); 8847c6ae99SBarry Smith 899566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 909566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 919566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw, xmin, ymin, xmax, ymax)); 925b399a63SLisandro Dalcin 93d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 9447c6ae99SBarry Smith /* first processor draw all node lines */ 95dd400576SPatrick Sanan if (rank == 0) { 9647c6ae99SBarry Smith for (k = 0; k < dd->P; k++) { 979371c9d4SSatish Balay ymin = 0.0; 989371c9d4SSatish Balay ymax = (PetscReal)(dd->N - 1); 9948a46eb9SPierre 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)); 1009371c9d4SSatish Balay xmin = (PetscReal)(k * (dd->M + 1)); 1019371c9d4SSatish Balay xmax = xmin + (PetscReal)(dd->M - 1); 10248a46eb9SPierre Jolivet for (ymin = 0; ymin < (PetscReal)dd->N; ymin++) PetscCall(PetscDrawLine(draw, xmin, ymin, xmax, ymin, PETSC_DRAW_BLACK)); 10347c6ae99SBarry Smith } 10447c6ae99SBarry Smith } 105d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1069566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1079566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 10847c6ae99SBarry Smith 109d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 1105b399a63SLisandro Dalcin /*Go through and draw for each plane*/ 1115b399a63SLisandro Dalcin for (k = 0; k < dd->P; k++) { 11247c6ae99SBarry Smith if ((k >= dd->zs) && (k < dd->ze)) { 11347c6ae99SBarry Smith /* draw my box */ 11447c6ae99SBarry Smith ymin = dd->ys; 11547c6ae99SBarry Smith ymax = dd->ye - 1; 11647c6ae99SBarry Smith xmin = dd->xs / dd->w + (dd->M + 1) * k; 11747c6ae99SBarry Smith xmax = (dd->xe - 1) / dd->w + (dd->M + 1) * k; 11847c6ae99SBarry Smith 1199566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymin, xmax, ymin, PETSC_DRAW_RED)); 1209566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymin, xmin, ymax, PETSC_DRAW_RED)); 1219566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymax, xmax, ymax, PETSC_DRAW_RED)); 1229566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmax, ymin, xmax, ymax, PETSC_DRAW_RED)); 12347c6ae99SBarry Smith 12447c6ae99SBarry Smith xmin = dd->xs / dd->w; 12547c6ae99SBarry Smith xmax = (dd->xe - 1) / dd->w; 12647c6ae99SBarry Smith 127832b7cebSLisandro Dalcin /* identify which processor owns the box */ 1289566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node, sizeof(node), "%d", (int)rank)); 1299566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, xmin + (dd->M + 1) * k + .2, ymin + .3, PETSC_DRAW_RED, node)); 13047c6ae99SBarry Smith /* put in numbers*/ 13147c6ae99SBarry Smith base = (dd->base + (dd->xe - dd->xs) * (dd->ye - dd->ys) * (k - dd->zs)) / dd->w; 13247c6ae99SBarry Smith for (y = ymin; y <= ymax; y++) { 13347c6ae99SBarry Smith for (x = xmin + (dd->M + 1) * k; x <= xmax + (dd->M + 1) * k; x++) { 1349566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node, sizeof(node), "%d", (int)base++)); 1359566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, x, y, PETSC_DRAW_BLACK, node)); 13647c6ae99SBarry Smith } 13747c6ae99SBarry Smith } 13847c6ae99SBarry Smith } 13947c6ae99SBarry Smith } 140d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1419566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1429566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 14347c6ae99SBarry Smith 144d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 14547c6ae99SBarry Smith for (k = 0 - dd->s; k < dd->P + dd->s; k++) { 14647c6ae99SBarry Smith /* Go through and draw for each plane */ 14747c6ae99SBarry Smith if ((k >= dd->Zs) && (k < dd->Ze)) { 14847c6ae99SBarry Smith /* overlay ghost numbers, useful for error checking */ 149565245c5SBarry Smith base = (dd->Xe - dd->Xs) * (dd->Ye - dd->Ys) * (k - dd->Zs) / dd->w; 1509566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(da->ltogmap, &idx)); 15147c6ae99SBarry Smith plane = k; 152565245c5SBarry Smith /* Keep z wrap around points on the drawing */ 1538865f1eaSKarl Rupp if (k < 0) plane = dd->P + k; 1548865f1eaSKarl Rupp if (k >= dd->P) plane = k - dd->P; 1559371c9d4SSatish Balay ymin = dd->Ys; 1569371c9d4SSatish Balay ymax = dd->Ye; 15747c6ae99SBarry Smith xmin = (dd->M + 1) * plane * dd->w; 15847c6ae99SBarry Smith xmax = (dd->M + 1) * plane * dd->w + dd->M * dd->w; 15947c6ae99SBarry Smith for (y = ymin; y < ymax; y++) { 16047c6ae99SBarry Smith for (x = xmin + dd->Xs; x < xmin + dd->Xe; x += dd->w) { 161a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(node, PETSC_STATIC_ARRAY_LENGTH(node), "%" PetscInt_FMT, idx[base])); 16247c6ae99SBarry Smith ycoord = y; 16347c6ae99SBarry Smith /*Keep y wrap around points on drawing */ 1648865f1eaSKarl Rupp if (y < 0) ycoord = dd->N + y; 1658865f1eaSKarl Rupp if (y >= dd->N) ycoord = y - dd->N; 16647c6ae99SBarry Smith xcoord = x; /* Keep x wrap points on drawing */ 1678865f1eaSKarl Rupp if (x < xmin) xcoord = xmax - (xmin - x); 1688865f1eaSKarl Rupp if (x >= xmax) xcoord = xmin + (x - xmax); 1699566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, xcoord / dd->w, ycoord, PETSC_DRAW_BLUE, node)); 170565245c5SBarry Smith base++; 17147c6ae99SBarry Smith } 17247c6ae99SBarry Smith } 1739566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap, &idx)); 17447c6ae99SBarry Smith } 17547c6ae99SBarry Smith } 176d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1779566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1789566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 1799566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 180c9493c35SStefano Zampini } else if (isglvis) { 1819566063dSJacob Faibussowitsch PetscCall(DMView_DA_GLVis(da, viewer)); 1829a42bb27SBarry Smith } else if (isbinary) { 1839566063dSJacob Faibussowitsch PetscCall(DMView_DA_Binary(da, viewer)); 184d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 1859a42bb27SBarry Smith } else if (ismatlab) { 1869566063dSJacob Faibussowitsch PetscCall(DMView_DA_Matlab(da, viewer)); 1879a42bb27SBarry Smith #endif 18811aeaf0aSBarry Smith } 1893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19047c6ae99SBarry Smith } 19147c6ae99SBarry Smith 192d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetUp_DA_3D(DM da) 193d71ae5a4SJacob Faibussowitsch { 19447c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 19547c6ae99SBarry Smith const PetscInt M = dd->M; 19647c6ae99SBarry Smith const PetscInt N = dd->N; 19747c6ae99SBarry Smith const PetscInt P = dd->P; 19847c6ae99SBarry Smith PetscInt m = dd->m; 19947c6ae99SBarry Smith PetscInt n = dd->n; 20047c6ae99SBarry Smith PetscInt p = dd->p; 20147c6ae99SBarry Smith const PetscInt dof = dd->w; 20247c6ae99SBarry Smith const PetscInt s = dd->s; 203bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 204bff4a2f0SMatthew G. Knepley DMBoundaryType by = dd->by; 205bff4a2f0SMatthew G. Knepley DMBoundaryType bz = dd->bz; 20619fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 20747c6ae99SBarry Smith PetscInt *lx = dd->lx; 20847c6ae99SBarry Smith PetscInt *ly = dd->ly; 20947c6ae99SBarry Smith PetscInt *lz = dd->lz; 21047c6ae99SBarry Smith MPI_Comm comm; 21147c6ae99SBarry Smith PetscMPIInt rank, size; 212ce00eea3SSatish Balay PetscInt xs = 0, xe, ys = 0, ye, zs = 0, ze, x = 0, y = 0, z = 0; 213bd1fc5aeSBarry Smith PetscInt Xs, Xe, Ys, Ye, Zs, Ze, IXs, IXe, IYs, IYe, IZs, IZe, pm; 2148ea3bf28SBarry Smith PetscInt left, right, up, down, bottom, top, i, j, k, *idx, nn; 21547c6ae99SBarry Smith PetscInt n0, n1, n2, n3, n4, n5, n6, n7, n8, n9, n10, n11, n12, n14; 21647c6ae99SBarry Smith PetscInt n15, n16, n17, n18, n19, n20, n21, n22, n23, n24, n25, n26; 217db87c5ecSEthan Coon PetscInt *bases, *ldims, base, x_t, y_t, z_t, s_t, count, s_x, s_y, s_z; 21847c6ae99SBarry Smith PetscInt sn0 = 0, sn1 = 0, sn2 = 0, sn3 = 0, sn5 = 0, sn6 = 0, sn7 = 0; 21947c6ae99SBarry Smith PetscInt sn8 = 0, sn9 = 0, sn11 = 0, sn15 = 0, sn24 = 0, sn25 = 0, sn26 = 0; 22047c6ae99SBarry Smith PetscInt sn17 = 0, sn18 = 0, sn19 = 0, sn20 = 0, sn21 = 0, sn23 = 0; 22147c6ae99SBarry Smith Vec local, global; 222bd1fc5aeSBarry Smith VecScatter gtol; 22345b6f7e9SBarry Smith IS to, from; 2246f951b95Secoon PetscBool twod; 22547c6ae99SBarry Smith 22647c6ae99SBarry Smith PetscFunctionBegin; 2277a8be351SBarry 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"); 2289566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 2293855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 2307de69702SBarry 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); 2313855c12bSBarry Smith #endif 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) { 23763a3b9bcSJacob Faibussowitsch PetscCheck(m >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in X direction: %" PetscInt_FMT, m); 238f7d195e4SLawrence Mitchell PetscCheck(m <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in X direction: %" PetscInt_FMT " %d", m, size); 23947c6ae99SBarry Smith } 24047c6ae99SBarry Smith if (n != PETSC_DECIDE) { 24163a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in Y direction: %" PetscInt_FMT, n); 242f7d195e4SLawrence Mitchell PetscCheck(n <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in Y direction: %" PetscInt_FMT " %d", n, size); 24347c6ae99SBarry Smith } 24447c6ae99SBarry Smith if (p != PETSC_DECIDE) { 24563a3b9bcSJacob Faibussowitsch PetscCheck(p >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in Z direction: %" PetscInt_FMT, p); 246f7d195e4SLawrence Mitchell PetscCheck(p <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in Z direction: %" PetscInt_FMT " %d", p, size); 24747c6ae99SBarry Smith } 2481dca8a05SBarry Smith PetscCheck(m <= 0 || n <= 0 || p <= 0 || m * n * p == size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "m %" PetscInt_FMT " * n %" PetscInt_FMT " * p %" PetscInt_FMT " != size %d", m, n, p, size); 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 } 26663a3b9bcSJacob Faibussowitsch PetscCheck(m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad p value: p = %" PetscInt_FMT, p); 2679371c9d4SSatish Balay if (M > N && m < n) { 2689371c9d4SSatish Balay PetscInt _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 } 28163a3b9bcSJacob Faibussowitsch PetscCheck(m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad n value: n = %" PetscInt_FMT, n); 2829371c9d4SSatish Balay if (M > P && m < p) { 2839371c9d4SSatish Balay PetscInt _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 } 29663a3b9bcSJacob Faibussowitsch PetscCheck(n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad m value: m = %" PetscInt_FMT, n); 2979371c9d4SSatish Balay if (N > P && n < p) { 2989371c9d4SSatish Balay PetscInt _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 */ 304369cc0aeSBarry Smith n = (PetscInt)(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; 312369cc0aeSBarry Smith m = (PetscInt)(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) { 3209371c9d4SSatish Balay PetscInt _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"); 32763a3b9bcSJacob Faibussowitsch PetscCheck(M >= m, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Partition in x direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT, M, m); 32863a3b9bcSJacob Faibussowitsch PetscCheck(N >= n, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Partition in y direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT, N, n); 32963a3b9bcSJacob Faibussowitsch PetscCheck(P >= p, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Partition in z direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT, 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 364fdc81ce8SEthan Coon /* note this is different than x- and y-, as we will handle as an important special 365bff4a2f0SMatthew G. Knepley case when p=P=1 and DM_BOUNDARY_PERIODIC and s > z. This is to deal with 2D problems 366fdc81ce8SEthan Coon in a 3D code. Additional code for this case is noted with "2d case" comments */ 3676f951b95Secoon twod = PETSC_FALSE; 3688865f1eaSKarl Rupp if (P == 1) twod = PETSC_TRUE; 3691dca8a05SBarry 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); 37047c6ae99SBarry Smith zs = 0; 3718865f1eaSKarl Rupp for (i = 0; i < (rank / (m * n)); i++) zs += lz[i]; 37247c6ae99SBarry Smith ye = ys + y; 37347c6ae99SBarry Smith xe = xs + x; 37447c6ae99SBarry Smith ze = zs + z; 37547c6ae99SBarry Smith 37688661749SPeter Brune /* determine ghost region (Xs) and region scattered into (IXs) */ 377d9c9ebe5SPeter Brune if (xs - s > 0) { 3789371c9d4SSatish Balay Xs = xs - s; 3799371c9d4SSatish Balay IXs = xs - s; 38088661749SPeter Brune } else { 3818865f1eaSKarl Rupp if (bx) Xs = xs - s; 3828865f1eaSKarl Rupp else Xs = 0; 38388661749SPeter Brune IXs = 0; 38488661749SPeter Brune } 385d9c9ebe5SPeter Brune if (xe + s <= M) { 3869371c9d4SSatish Balay Xe = xe + s; 3879371c9d4SSatish Balay IXe = xe + s; 38888661749SPeter Brune } else { 38988661749SPeter Brune if (bx) { 3909371c9d4SSatish Balay Xs = xs - s; 3919371c9d4SSatish Balay Xe = xe + s; 3928865f1eaSKarl Rupp } else Xe = M; 39388661749SPeter Brune IXe = M; 39488661749SPeter Brune } 39547c6ae99SBarry Smith 396bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 397d9c9ebe5SPeter Brune IXs = xs - s; 398d9c9ebe5SPeter Brune IXe = xe + s; 399d9c9ebe5SPeter Brune Xs = xs - s; 400d9c9ebe5SPeter Brune Xe = xe + s; 40188661749SPeter Brune } 40288661749SPeter Brune 403d9c9ebe5SPeter Brune if (ys - s > 0) { 4049371c9d4SSatish Balay Ys = ys - s; 4059371c9d4SSatish Balay IYs = ys - s; 40688661749SPeter Brune } else { 4078865f1eaSKarl Rupp if (by) Ys = ys - s; 4088865f1eaSKarl Rupp else Ys = 0; 40988661749SPeter Brune IYs = 0; 41088661749SPeter Brune } 411d9c9ebe5SPeter Brune if (ye + s <= N) { 4129371c9d4SSatish Balay Ye = ye + s; 4139371c9d4SSatish Balay IYe = ye + s; 41488661749SPeter Brune } else { 4158865f1eaSKarl Rupp if (by) Ye = ye + s; 4168865f1eaSKarl Rupp else Ye = N; 41788661749SPeter Brune IYe = N; 41888661749SPeter Brune } 41988661749SPeter Brune 420bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 421d9c9ebe5SPeter Brune IYs = ys - s; 422d9c9ebe5SPeter Brune IYe = ye + s; 423d9c9ebe5SPeter Brune Ys = ys - s; 424d9c9ebe5SPeter Brune Ye = ye + s; 42588661749SPeter Brune } 42688661749SPeter Brune 427d9c9ebe5SPeter Brune if (zs - s > 0) { 4289371c9d4SSatish Balay Zs = zs - s; 4299371c9d4SSatish Balay IZs = zs - s; 43088661749SPeter Brune } else { 4318865f1eaSKarl Rupp if (bz) Zs = zs - s; 4328865f1eaSKarl Rupp else Zs = 0; 43388661749SPeter Brune IZs = 0; 43488661749SPeter Brune } 435d9c9ebe5SPeter Brune if (ze + s <= P) { 4369371c9d4SSatish Balay Ze = ze + s; 4379371c9d4SSatish Balay IZe = ze + s; 43888661749SPeter Brune } else { 4398865f1eaSKarl Rupp if (bz) Ze = ze + s; 4408865f1eaSKarl Rupp else Ze = P; 44188661749SPeter Brune IZe = P; 44288661749SPeter Brune } 44388661749SPeter Brune 444bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_MIRROR) { 445d9c9ebe5SPeter Brune IZs = zs - s; 446d9c9ebe5SPeter Brune IZe = ze + s; 447d9c9ebe5SPeter Brune Zs = zs - s; 448d9c9ebe5SPeter Brune Ze = ze + s; 44988661749SPeter Brune } 45047c6ae99SBarry Smith 45147c6ae99SBarry Smith /* Resize all X parameters to reflect w */ 452d9c9ebe5SPeter Brune s_x = s; 453d9c9ebe5SPeter Brune s_y = s; 454d9c9ebe5SPeter Brune s_z = s; 45547c6ae99SBarry Smith 45647c6ae99SBarry Smith /* determine starting point of each processor */ 45747c6ae99SBarry Smith nn = x * y * z; 4589566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size + 1, &bases, size, &ldims)); 4599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nn, 1, MPIU_INT, ldims, 1, MPIU_INT, comm)); 46047c6ae99SBarry Smith bases[0] = 0; 4618865f1eaSKarl Rupp for (i = 1; i <= size; i++) bases[i] = ldims[i - 1]; 4628865f1eaSKarl Rupp for (i = 1; i <= size; i++) bases[i] += bases[i - 1]; 463ce00eea3SSatish Balay base = bases[rank] * dof; 46447c6ae99SBarry Smith 46547c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 466ce00eea3SSatish Balay dd->Nlocal = x * y * z * dof; 4679566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, dof, dd->Nlocal, PETSC_DECIDE, NULL, &global)); 468ce00eea3SSatish Balay dd->nlocal = (Xe - Xs) * (Ye - Ys) * (Ze - Zs) * dof; 4699566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, dof, dd->nlocal, NULL, &local)); 47047c6ae99SBarry Smith 4714104a7a0SPatrick Sanan /* generate global to local vector scatter and local to global mapping*/ 47247c6ae99SBarry Smith 473ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 474ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 4759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((IXe - IXs) * (IYe - IYs) * (IZe - IZs), &idx)); 476d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 4779371c9d4SSatish Balay left = IXs - Xs; 4789371c9d4SSatish Balay right = left + (IXe - IXs); 4799371c9d4SSatish Balay bottom = IYs - Ys; 4809371c9d4SSatish Balay top = bottom + (IYe - IYs); 4819371c9d4SSatish Balay down = IZs - Zs; 4829371c9d4SSatish Balay up = down + (IZe - IZs); 483ce00eea3SSatish Balay count = 0; 484ce00eea3SSatish Balay for (i = down; i < up; i++) { 485ce00eea3SSatish Balay for (j = bottom; j < top; j++) { 486ad540459SPierre Jolivet for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k; 487ce00eea3SSatish Balay } 488ce00eea3SSatish Balay } 4899566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm, dof, count, idx, PETSC_OWN_POINTER, &to)); 49047c6ae99SBarry Smith } else { 49147c6ae99SBarry Smith /* This is way ugly! We need to list the funny cross type region */ 4929371c9d4SSatish Balay left = xs - Xs; 4939371c9d4SSatish Balay right = left + x; 4949371c9d4SSatish Balay bottom = ys - Ys; 4959371c9d4SSatish Balay top = bottom + y; 4969371c9d4SSatish Balay down = zs - Zs; 4979371c9d4SSatish Balay up = down + z; 49847c6ae99SBarry Smith count = 0; 499a5b23f4aSJose E. Roman /* the bottom chunk */ 500ce00eea3SSatish Balay for (i = (IZs - Zs); i < down; i++) { 50147c6ae99SBarry Smith for (j = bottom; j < top; j++) { 502ce00eea3SSatish Balay for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k; 50347c6ae99SBarry Smith } 50447c6ae99SBarry Smith } 50547c6ae99SBarry Smith /* the middle piece */ 50647c6ae99SBarry Smith for (i = down; i < up; i++) { 50747c6ae99SBarry Smith /* front */ 508ce00eea3SSatish Balay for (j = (IYs - Ys); j < bottom; j++) { 509ce00eea3SSatish Balay for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k; 51047c6ae99SBarry Smith } 51147c6ae99SBarry Smith /* middle */ 51247c6ae99SBarry Smith for (j = bottom; j < top; j++) { 513ce00eea3SSatish Balay for (k = IXs - Xs; k < IXe - Xs; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k; 51447c6ae99SBarry Smith } 51547c6ae99SBarry Smith /* back */ 516ce00eea3SSatish Balay for (j = top; j < top + IYe - ye; j++) { 517ce00eea3SSatish Balay for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k; 51847c6ae99SBarry Smith } 51947c6ae99SBarry Smith } 52047c6ae99SBarry Smith /* the top piece */ 521ce00eea3SSatish Balay for (i = up; i < up + IZe - ze; i++) { 52247c6ae99SBarry Smith for (j = bottom; j < top; j++) { 523ce00eea3SSatish Balay for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k; 52447c6ae99SBarry Smith } 52547c6ae99SBarry Smith } 5269566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm, dof, count, idx, PETSC_OWN_POINTER, &to)); 52747c6ae99SBarry Smith } 52847c6ae99SBarry Smith 52947c6ae99SBarry Smith /* determine who lies on each side of use stored in n24 n25 n26 53047c6ae99SBarry Smith n21 n22 n23 53147c6ae99SBarry Smith n18 n19 n20 53247c6ae99SBarry Smith 53347c6ae99SBarry Smith n15 n16 n17 53447c6ae99SBarry Smith n12 n14 53547c6ae99SBarry Smith n9 n10 n11 53647c6ae99SBarry Smith 53747c6ae99SBarry Smith n6 n7 n8 53847c6ae99SBarry Smith n3 n4 n5 53947c6ae99SBarry Smith n0 n1 n2 54047c6ae99SBarry Smith */ 54147c6ae99SBarry Smith 54247c6ae99SBarry Smith /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */ 54347c6ae99SBarry Smith /* Assume Nodes are Internal to the Cube */ 54447c6ae99SBarry Smith n0 = rank - m * n - m - 1; 54547c6ae99SBarry Smith n1 = rank - m * n - m; 54647c6ae99SBarry Smith n2 = rank - m * n - m + 1; 54747c6ae99SBarry Smith n3 = rank - m * n - 1; 54847c6ae99SBarry Smith n4 = rank - m * n; 54947c6ae99SBarry Smith n5 = rank - m * n + 1; 55047c6ae99SBarry Smith n6 = rank - m * n + m - 1; 55147c6ae99SBarry Smith n7 = rank - m * n + m; 55247c6ae99SBarry Smith n8 = rank - m * n + m + 1; 55347c6ae99SBarry Smith 55447c6ae99SBarry Smith n9 = rank - m - 1; 55547c6ae99SBarry Smith n10 = rank - m; 55647c6ae99SBarry Smith n11 = rank - m + 1; 55747c6ae99SBarry Smith n12 = rank - 1; 55847c6ae99SBarry Smith n14 = rank + 1; 55947c6ae99SBarry Smith n15 = rank + m - 1; 56047c6ae99SBarry Smith n16 = rank + m; 56147c6ae99SBarry Smith n17 = rank + m + 1; 56247c6ae99SBarry Smith 56347c6ae99SBarry Smith n18 = rank + m * n - m - 1; 56447c6ae99SBarry Smith n19 = rank + m * n - m; 56547c6ae99SBarry Smith n20 = rank + m * n - m + 1; 56647c6ae99SBarry Smith n21 = rank + m * n - 1; 56747c6ae99SBarry Smith n22 = rank + m * n; 56847c6ae99SBarry Smith n23 = rank + m * n + 1; 56947c6ae99SBarry Smith n24 = rank + m * n + m - 1; 57047c6ae99SBarry Smith n25 = rank + m * n + m; 57147c6ae99SBarry Smith n26 = rank + m * n + m + 1; 57247c6ae99SBarry Smith 57347c6ae99SBarry Smith /* Assume Pieces are on Faces of Cube */ 57447c6ae99SBarry Smith 57547c6ae99SBarry Smith if (xs == 0) { /* First assume not corner or edge */ 57647c6ae99SBarry Smith n0 = rank - 1 - (m * n); 57747c6ae99SBarry Smith n3 = rank + m - 1 - (m * n); 57847c6ae99SBarry Smith n6 = rank + 2 * m - 1 - (m * n); 57947c6ae99SBarry Smith n9 = rank - 1; 58047c6ae99SBarry Smith n12 = rank + m - 1; 58147c6ae99SBarry Smith n15 = rank + 2 * m - 1; 58247c6ae99SBarry Smith n18 = rank - 1 + (m * n); 58347c6ae99SBarry Smith n21 = rank + m - 1 + (m * n); 58447c6ae99SBarry Smith n24 = rank + 2 * m - 1 + (m * n); 58547c6ae99SBarry Smith } 58647c6ae99SBarry Smith 587ce00eea3SSatish Balay if (xe == M) { /* First assume not corner or edge */ 58847c6ae99SBarry Smith n2 = rank - 2 * m + 1 - (m * n); 58947c6ae99SBarry Smith n5 = rank - m + 1 - (m * n); 59047c6ae99SBarry Smith n8 = rank + 1 - (m * n); 59147c6ae99SBarry Smith n11 = rank - 2 * m + 1; 59247c6ae99SBarry Smith n14 = rank - m + 1; 59347c6ae99SBarry Smith n17 = rank + 1; 59447c6ae99SBarry Smith n20 = rank - 2 * m + 1 + (m * n); 59547c6ae99SBarry Smith n23 = rank - m + 1 + (m * n); 59647c6ae99SBarry Smith n26 = rank + 1 + (m * n); 59747c6ae99SBarry Smith } 59847c6ae99SBarry Smith 59947c6ae99SBarry Smith if (ys == 0) { /* First assume not corner or edge */ 60047c6ae99SBarry Smith n0 = rank + m * (n - 1) - 1 - (m * n); 60147c6ae99SBarry Smith n1 = rank + m * (n - 1) - (m * n); 60247c6ae99SBarry Smith n2 = rank + m * (n - 1) + 1 - (m * n); 60347c6ae99SBarry Smith n9 = rank + m * (n - 1) - 1; 60447c6ae99SBarry Smith n10 = rank + m * (n - 1); 60547c6ae99SBarry Smith n11 = rank + m * (n - 1) + 1; 60647c6ae99SBarry Smith n18 = rank + m * (n - 1) - 1 + (m * n); 60747c6ae99SBarry Smith n19 = rank + m * (n - 1) + (m * n); 60847c6ae99SBarry Smith n20 = rank + m * (n - 1) + 1 + (m * n); 60947c6ae99SBarry Smith } 61047c6ae99SBarry Smith 61147c6ae99SBarry Smith if (ye == N) { /* First assume not corner or edge */ 61247c6ae99SBarry Smith n6 = rank - m * (n - 1) - 1 - (m * n); 61347c6ae99SBarry Smith n7 = rank - m * (n - 1) - (m * n); 61447c6ae99SBarry Smith n8 = rank - m * (n - 1) + 1 - (m * n); 61547c6ae99SBarry Smith n15 = rank - m * (n - 1) - 1; 61647c6ae99SBarry Smith n16 = rank - m * (n - 1); 61747c6ae99SBarry Smith n17 = rank - m * (n - 1) + 1; 61847c6ae99SBarry Smith n24 = rank - m * (n - 1) - 1 + (m * n); 61947c6ae99SBarry Smith n25 = rank - m * (n - 1) + (m * n); 62047c6ae99SBarry Smith n26 = rank - m * (n - 1) + 1 + (m * n); 62147c6ae99SBarry Smith } 62247c6ae99SBarry Smith 62347c6ae99SBarry Smith if (zs == 0) { /* First assume not corner or edge */ 62447c6ae99SBarry Smith n0 = size - (m * n) + rank - m - 1; 62547c6ae99SBarry Smith n1 = size - (m * n) + rank - m; 62647c6ae99SBarry Smith n2 = size - (m * n) + rank - m + 1; 62747c6ae99SBarry Smith n3 = size - (m * n) + rank - 1; 62847c6ae99SBarry Smith n4 = size - (m * n) + rank; 62947c6ae99SBarry Smith n5 = size - (m * n) + rank + 1; 63047c6ae99SBarry Smith n6 = size - (m * n) + rank + m - 1; 63147c6ae99SBarry Smith n7 = size - (m * n) + rank + m; 63247c6ae99SBarry Smith n8 = size - (m * n) + rank + m + 1; 63347c6ae99SBarry Smith } 63447c6ae99SBarry Smith 63547c6ae99SBarry Smith if (ze == P) { /* First assume not corner or edge */ 63647c6ae99SBarry Smith n18 = (m * n) - (size - rank) - m - 1; 63747c6ae99SBarry Smith n19 = (m * n) - (size - rank) - m; 63847c6ae99SBarry Smith n20 = (m * n) - (size - rank) - m + 1; 63947c6ae99SBarry Smith n21 = (m * n) - (size - rank) - 1; 64047c6ae99SBarry Smith n22 = (m * n) - (size - rank); 64147c6ae99SBarry Smith n23 = (m * n) - (size - rank) + 1; 64247c6ae99SBarry Smith n24 = (m * n) - (size - rank) + m - 1; 64347c6ae99SBarry Smith n25 = (m * n) - (size - rank) + m; 64447c6ae99SBarry Smith n26 = (m * n) - (size - rank) + m + 1; 64547c6ae99SBarry Smith } 64647c6ae99SBarry Smith 64747c6ae99SBarry Smith if ((xs == 0) && (zs == 0)) { /* Assume an edge, not corner */ 64847c6ae99SBarry Smith n0 = size - m * n + rank + m - 1 - m; 64947c6ae99SBarry Smith n3 = size - m * n + rank + m - 1; 65047c6ae99SBarry Smith n6 = size - m * n + rank + m - 1 + m; 65147c6ae99SBarry Smith } 65247c6ae99SBarry Smith 65347c6ae99SBarry Smith if ((xs == 0) && (ze == P)) { /* Assume an edge, not corner */ 65447c6ae99SBarry Smith n18 = m * n - (size - rank) + m - 1 - m; 65547c6ae99SBarry Smith n21 = m * n - (size - rank) + m - 1; 65647c6ae99SBarry Smith n24 = m * n - (size - rank) + m - 1 + m; 65747c6ae99SBarry Smith } 65847c6ae99SBarry Smith 65947c6ae99SBarry Smith if ((xs == 0) && (ys == 0)) { /* Assume an edge, not corner */ 66047c6ae99SBarry Smith n0 = rank + m * n - 1 - m * n; 66147c6ae99SBarry Smith n9 = rank + m * n - 1; 66247c6ae99SBarry Smith n18 = rank + m * n - 1 + m * n; 66347c6ae99SBarry Smith } 66447c6ae99SBarry Smith 66547c6ae99SBarry Smith if ((xs == 0) && (ye == N)) { /* Assume an edge, not corner */ 66647c6ae99SBarry Smith n6 = rank - m * (n - 1) + m - 1 - m * n; 66747c6ae99SBarry Smith n15 = rank - m * (n - 1) + m - 1; 66847c6ae99SBarry Smith n24 = rank - m * (n - 1) + m - 1 + m * n; 66947c6ae99SBarry Smith } 67047c6ae99SBarry Smith 671ce00eea3SSatish Balay if ((xe == M) && (zs == 0)) { /* Assume an edge, not corner */ 67247c6ae99SBarry Smith n2 = size - (m * n - rank) - (m - 1) - m; 67347c6ae99SBarry Smith n5 = size - (m * n - rank) - (m - 1); 67447c6ae99SBarry Smith n8 = size - (m * n - rank) - (m - 1) + m; 67547c6ae99SBarry Smith } 67647c6ae99SBarry Smith 677ce00eea3SSatish Balay if ((xe == M) && (ze == P)) { /* Assume an edge, not corner */ 67847c6ae99SBarry Smith n20 = m * n - (size - rank) - (m - 1) - m; 67947c6ae99SBarry Smith n23 = m * n - (size - rank) - (m - 1); 68047c6ae99SBarry Smith n26 = m * n - (size - rank) - (m - 1) + m; 68147c6ae99SBarry Smith } 68247c6ae99SBarry Smith 683ce00eea3SSatish Balay if ((xe == M) && (ys == 0)) { /* Assume an edge, not corner */ 68447c6ae99SBarry Smith n2 = rank + m * (n - 1) - (m - 1) - m * n; 68547c6ae99SBarry Smith n11 = rank + m * (n - 1) - (m - 1); 68647c6ae99SBarry Smith n20 = rank + m * (n - 1) - (m - 1) + m * n; 68747c6ae99SBarry Smith } 68847c6ae99SBarry Smith 689ce00eea3SSatish Balay if ((xe == M) && (ye == N)) { /* Assume an edge, not corner */ 69047c6ae99SBarry Smith n8 = rank - m * n + 1 - m * n; 69147c6ae99SBarry Smith n17 = rank - m * n + 1; 69247c6ae99SBarry Smith n26 = rank - m * n + 1 + m * n; 69347c6ae99SBarry Smith } 69447c6ae99SBarry Smith 69547c6ae99SBarry Smith if ((ys == 0) && (zs == 0)) { /* Assume an edge, not corner */ 69647c6ae99SBarry Smith n0 = size - m + rank - 1; 69747c6ae99SBarry Smith n1 = size - m + rank; 69847c6ae99SBarry Smith n2 = size - m + rank + 1; 69947c6ae99SBarry Smith } 70047c6ae99SBarry Smith 70147c6ae99SBarry Smith if ((ys == 0) && (ze == P)) { /* Assume an edge, not corner */ 70247c6ae99SBarry Smith n18 = m * n - (size - rank) + m * (n - 1) - 1; 70347c6ae99SBarry Smith n19 = m * n - (size - rank) + m * (n - 1); 70447c6ae99SBarry Smith n20 = m * n - (size - rank) + m * (n - 1) + 1; 70547c6ae99SBarry Smith } 70647c6ae99SBarry Smith 70747c6ae99SBarry Smith if ((ye == N) && (zs == 0)) { /* Assume an edge, not corner */ 70847c6ae99SBarry Smith n6 = size - (m * n - rank) - m * (n - 1) - 1; 70947c6ae99SBarry Smith n7 = size - (m * n - rank) - m * (n - 1); 71047c6ae99SBarry Smith n8 = size - (m * n - rank) - m * (n - 1) + 1; 71147c6ae99SBarry Smith } 71247c6ae99SBarry Smith 71347c6ae99SBarry Smith if ((ye == N) && (ze == P)) { /* Assume an edge, not corner */ 71447c6ae99SBarry Smith n24 = rank - (size - m) - 1; 71547c6ae99SBarry Smith n25 = rank - (size - m); 71647c6ae99SBarry Smith n26 = rank - (size - m) + 1; 71747c6ae99SBarry Smith } 71847c6ae99SBarry Smith 71947c6ae99SBarry Smith /* Check for Corners */ 7208865f1eaSKarl Rupp if ((xs == 0) && (ys == 0) && (zs == 0)) n0 = size - 1; 7218865f1eaSKarl Rupp if ((xs == 0) && (ys == 0) && (ze == P)) n18 = m * n - 1; 7228865f1eaSKarl Rupp if ((xs == 0) && (ye == N) && (zs == 0)) n6 = (size - 1) - m * (n - 1); 7238865f1eaSKarl Rupp if ((xs == 0) && (ye == N) && (ze == P)) n24 = m - 1; 7248865f1eaSKarl Rupp if ((xe == M) && (ys == 0) && (zs == 0)) n2 = size - m; 7258865f1eaSKarl Rupp if ((xe == M) && (ys == 0) && (ze == P)) n20 = m * n - m; 7268865f1eaSKarl Rupp if ((xe == M) && (ye == N) && (zs == 0)) n8 = size - m * n; 7278865f1eaSKarl Rupp if ((xe == M) && (ye == N) && (ze == P)) n26 = 0; 72847c6ae99SBarry Smith 72947c6ae99SBarry Smith /* Check for when not X,Y, and Z Periodic */ 73047c6ae99SBarry Smith 73147c6ae99SBarry Smith /* If not X periodic */ 732bff4a2f0SMatthew G. Knepley if (bx != DM_BOUNDARY_PERIODIC) { 7338865f1eaSKarl Rupp if (xs == 0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2; 7348865f1eaSKarl Rupp if (xe == M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2; 73547c6ae99SBarry Smith } 73647c6ae99SBarry Smith 73747c6ae99SBarry Smith /* If not Y periodic */ 738bff4a2f0SMatthew G. Knepley if (by != DM_BOUNDARY_PERIODIC) { 7398865f1eaSKarl Rupp if (ys == 0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2; 7408865f1eaSKarl Rupp if (ye == N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2; 74147c6ae99SBarry Smith } 74247c6ae99SBarry Smith 74347c6ae99SBarry Smith /* If not Z periodic */ 744bff4a2f0SMatthew G. Knepley if (bz != DM_BOUNDARY_PERIODIC) { 7458865f1eaSKarl Rupp if (zs == 0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2; 7468865f1eaSKarl Rupp if (ze == P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2; 74747c6ae99SBarry Smith } 74847c6ae99SBarry Smith 7499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(27, &dd->neighbors)); 7508865f1eaSKarl Rupp 75147c6ae99SBarry Smith dd->neighbors[0] = n0; 75247c6ae99SBarry Smith dd->neighbors[1] = n1; 75347c6ae99SBarry Smith dd->neighbors[2] = n2; 75447c6ae99SBarry Smith dd->neighbors[3] = n3; 75547c6ae99SBarry Smith dd->neighbors[4] = n4; 75647c6ae99SBarry Smith dd->neighbors[5] = n5; 75747c6ae99SBarry Smith dd->neighbors[6] = n6; 75847c6ae99SBarry Smith dd->neighbors[7] = n7; 75947c6ae99SBarry Smith dd->neighbors[8] = n8; 76047c6ae99SBarry Smith dd->neighbors[9] = n9; 76147c6ae99SBarry Smith dd->neighbors[10] = n10; 76247c6ae99SBarry Smith dd->neighbors[11] = n11; 76347c6ae99SBarry Smith dd->neighbors[12] = n12; 76447c6ae99SBarry Smith dd->neighbors[13] = rank; 76547c6ae99SBarry Smith dd->neighbors[14] = n14; 76647c6ae99SBarry Smith dd->neighbors[15] = n15; 76747c6ae99SBarry Smith dd->neighbors[16] = n16; 76847c6ae99SBarry Smith dd->neighbors[17] = n17; 76947c6ae99SBarry Smith dd->neighbors[18] = n18; 77047c6ae99SBarry Smith dd->neighbors[19] = n19; 77147c6ae99SBarry Smith dd->neighbors[20] = n20; 77247c6ae99SBarry Smith dd->neighbors[21] = n21; 77347c6ae99SBarry Smith dd->neighbors[22] = n22; 77447c6ae99SBarry Smith dd->neighbors[23] = n23; 77547c6ae99SBarry Smith dd->neighbors[24] = n24; 77647c6ae99SBarry Smith dd->neighbors[25] = n25; 77747c6ae99SBarry Smith dd->neighbors[26] = n26; 77847c6ae99SBarry Smith 77947c6ae99SBarry Smith /* If star stencil then delete the corner neighbors */ 780d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 78147c6ae99SBarry Smith /* save information about corner neighbors */ 7829371c9d4SSatish Balay sn0 = n0; 7839371c9d4SSatish Balay sn1 = n1; 7849371c9d4SSatish Balay sn2 = n2; 7859371c9d4SSatish Balay sn3 = n3; 7869371c9d4SSatish Balay sn5 = n5; 7879371c9d4SSatish Balay sn6 = n6; 7889371c9d4SSatish Balay sn7 = n7; 7899371c9d4SSatish Balay sn8 = n8; 7909371c9d4SSatish Balay sn9 = n9; 7919371c9d4SSatish Balay sn11 = n11; 7929371c9d4SSatish Balay sn15 = n15; 7939371c9d4SSatish Balay sn17 = n17; 7949371c9d4SSatish Balay sn18 = n18; 7959371c9d4SSatish Balay sn19 = n19; 7969371c9d4SSatish Balay sn20 = n20; 7979371c9d4SSatish Balay sn21 = n21; 7989371c9d4SSatish Balay sn23 = n23; 7999371c9d4SSatish Balay sn24 = n24; 8009371c9d4SSatish Balay sn25 = n25; 80147c6ae99SBarry Smith sn26 = n26; 8028865f1eaSKarl Rupp n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; 80347c6ae99SBarry Smith } 80447c6ae99SBarry Smith 8059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((Xe - Xs) * (Ye - Ys) * (Ze - Zs), &idx)); 80647c6ae99SBarry Smith 80747c6ae99SBarry Smith nn = 0; 80847c6ae99SBarry Smith /* Bottom Level */ 80947c6ae99SBarry Smith for (k = 0; k < s_z; k++) { 81047c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 81147c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 812ce00eea3SSatish Balay x_t = lx[n0 % m]; 81347c6ae99SBarry Smith y_t = ly[(n0 % (m * n)) / m]; 81447c6ae99SBarry Smith z_t = lz[n0 / (m * n)]; 81547c6ae99SBarry 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; 8168865f1eaSKarl 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 */ 8178865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 81847c6ae99SBarry Smith } 81947c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 82047c6ae99SBarry Smith x_t = x; 82147c6ae99SBarry Smith y_t = ly[(n1 % (m * n)) / m]; 82247c6ae99SBarry Smith z_t = lz[n1 / (m * n)]; 82347c6ae99SBarry 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; 8248865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n1] + x_t * y_t * z_t - (s_y + 1 - i) * x_t; /* 2D case */ 8258865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 82647c6ae99SBarry Smith } 82747c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 828ce00eea3SSatish Balay x_t = lx[n2 % m]; 82947c6ae99SBarry Smith y_t = ly[(n2 % (m * n)) / m]; 83047c6ae99SBarry Smith z_t = lz[n2 / (m * n)]; 83147c6ae99SBarry 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; 8328865f1eaSKarl Rupp if (twod && (s_t < 0)) s_t = bases[n2] + x_t * y_t * z_t - (s_y + 1 - i) * x_t; /* 2D case */ 8338865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 83447c6ae99SBarry Smith } 83547c6ae99SBarry Smith } 83647c6ae99SBarry Smith 83747c6ae99SBarry Smith for (i = 0; i < y; i++) { 83847c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 839ce00eea3SSatish Balay x_t = lx[n3 % m]; 84047c6ae99SBarry Smith y_t = y; 84147c6ae99SBarry Smith z_t = lz[n3 / (m * n)]; 84247c6ae99SBarry Smith s_t = bases[n3] + (i + 1) * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t; 8438865f1eaSKarl 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 */ 8448865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 84547c6ae99SBarry Smith } 84647c6ae99SBarry Smith 84747c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 84847c6ae99SBarry Smith x_t = x; 84947c6ae99SBarry Smith y_t = y; 85047c6ae99SBarry Smith z_t = lz[n4 / (m * n)]; 85147c6ae99SBarry Smith s_t = bases[n4] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; 8528865f1eaSKarl 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 */ 8538865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 854bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 855a6fd6b0aSBarry Smith for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + x * i + (s_z - k - 1) * x * y; 85647c6ae99SBarry Smith } 85747c6ae99SBarry Smith 85847c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 859ce00eea3SSatish Balay x_t = lx[n5 % m]; 86047c6ae99SBarry Smith y_t = y; 86147c6ae99SBarry Smith z_t = lz[n5 / (m * n)]; 86247c6ae99SBarry Smith s_t = bases[n5] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; 8638865f1eaSKarl 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 */ 8648865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 86547c6ae99SBarry Smith } 86647c6ae99SBarry Smith } 86747c6ae99SBarry Smith 86847c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 86947c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 870ce00eea3SSatish Balay x_t = lx[n6 % m]; 87147c6ae99SBarry Smith y_t = ly[(n6 % (m * n)) / m]; 87247c6ae99SBarry Smith z_t = lz[n6 / (m * n)]; 87347c6ae99SBarry Smith s_t = bases[n6] + i * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t; 8748865f1eaSKarl 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 */ 8758865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 87647c6ae99SBarry Smith } 87747c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 87847c6ae99SBarry Smith x_t = x; 87947c6ae99SBarry Smith y_t = ly[(n7 % (m * n)) / m]; 88047c6ae99SBarry Smith z_t = lz[n7 / (m * n)]; 88147c6ae99SBarry Smith s_t = bases[n7] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; 8828865f1eaSKarl 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 */ 8838865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 88447c6ae99SBarry Smith } 88547c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 886ce00eea3SSatish Balay x_t = lx[n8 % m]; 88747c6ae99SBarry Smith y_t = ly[(n8 % (m * n)) / m]; 88847c6ae99SBarry Smith z_t = lz[n8 / (m * n)]; 88947c6ae99SBarry Smith s_t = bases[n8] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; 8908865f1eaSKarl 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 */ 8918865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 89247c6ae99SBarry Smith } 89347c6ae99SBarry Smith } 89447c6ae99SBarry Smith } 89547c6ae99SBarry Smith 89647c6ae99SBarry Smith /* Middle Level */ 89747c6ae99SBarry Smith for (k = 0; k < z; k++) { 89847c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 89947c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 900ce00eea3SSatish Balay x_t = lx[n9 % m]; 90147c6ae99SBarry Smith y_t = ly[(n9 % (m * n)) / m]; 90247c6ae99SBarry Smith /* z_t = z; */ 90347c6ae99SBarry Smith s_t = bases[n9] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t; 9048865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 90547c6ae99SBarry Smith } 90647c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 90747c6ae99SBarry Smith x_t = x; 90847c6ae99SBarry Smith y_t = ly[(n10 % (m * n)) / m]; 90947c6ae99SBarry Smith /* z_t = z; */ 91047c6ae99SBarry Smith s_t = bases[n10] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; 9118865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 912bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 913a6fd6b0aSBarry Smith for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (s_y - i) * x; 91447c6ae99SBarry Smith } 91547c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 916ce00eea3SSatish Balay x_t = lx[n11 % m]; 91747c6ae99SBarry Smith y_t = ly[(n11 % (m * n)) / m]; 91847c6ae99SBarry Smith /* z_t = z; */ 91947c6ae99SBarry Smith s_t = bases[n11] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; 9208865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 92147c6ae99SBarry Smith } 92247c6ae99SBarry Smith } 92347c6ae99SBarry Smith 92447c6ae99SBarry Smith for (i = 0; i < y; i++) { 92547c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 926ce00eea3SSatish Balay x_t = lx[n12 % m]; 92747c6ae99SBarry Smith y_t = y; 92847c6ae99SBarry Smith /* z_t = z; */ 92947c6ae99SBarry Smith s_t = bases[n12] + (i + 1) * x_t - s_x + k * x_t * y_t; 9308865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 931bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 932a6fd6b0aSBarry Smith for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k * x * y + i * x; 93347c6ae99SBarry Smith } 93447c6ae99SBarry Smith 93547c6ae99SBarry Smith /* Interior */ 93647c6ae99SBarry Smith s_t = bases[rank] + i * x + k * x * y; 9378865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = s_t++; 93847c6ae99SBarry Smith 93947c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 940ce00eea3SSatish Balay x_t = lx[n14 % m]; 94147c6ae99SBarry Smith y_t = y; 94247c6ae99SBarry Smith /* z_t = z; */ 94347c6ae99SBarry Smith s_t = bases[n14] + i * x_t + k * x_t * y_t; 9448865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 945bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 946a6fd6b0aSBarry Smith for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k * x * y + i * x; 94747c6ae99SBarry Smith } 94847c6ae99SBarry Smith } 94947c6ae99SBarry Smith 95047c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 95147c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 952ce00eea3SSatish Balay x_t = lx[n15 % m]; 95347c6ae99SBarry Smith y_t = ly[(n15 % (m * n)) / m]; 95447c6ae99SBarry Smith /* z_t = z; */ 95547c6ae99SBarry Smith s_t = bases[n15] + i * x_t - s_x + k * x_t * y_t; 9568865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 95747c6ae99SBarry Smith } 95847c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 95947c6ae99SBarry Smith x_t = x; 96047c6ae99SBarry Smith y_t = ly[(n16 % (m * n)) / m]; 96147c6ae99SBarry Smith /* z_t = z; */ 96247c6ae99SBarry Smith s_t = bases[n16] + (i - 1) * x_t + k * x_t * y_t; 9638865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 964bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 965a6fd6b0aSBarry Smith for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (y - i) * x; 96647c6ae99SBarry Smith } 96747c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 968ce00eea3SSatish Balay x_t = lx[n17 % m]; 96947c6ae99SBarry Smith y_t = ly[(n17 % (m * n)) / m]; 97047c6ae99SBarry Smith /* z_t = z; */ 97147c6ae99SBarry Smith s_t = bases[n17] + (i - 1) * x_t + k * x_t * y_t; 9728865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 97347c6ae99SBarry Smith } 97447c6ae99SBarry Smith } 97547c6ae99SBarry Smith } 97647c6ae99SBarry Smith 97747c6ae99SBarry Smith /* Upper Level */ 97847c6ae99SBarry Smith for (k = 0; k < s_z; k++) { 97947c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 98047c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 981ce00eea3SSatish Balay x_t = lx[n18 % m]; 98247c6ae99SBarry Smith y_t = ly[(n18 % (m * n)) / m]; 98347c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 98447c6ae99SBarry Smith s_t = bases[n18] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t; 9858865f1eaSKarl 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 */ 9868865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 98747c6ae99SBarry Smith } 98847c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 98947c6ae99SBarry Smith x_t = x; 99047c6ae99SBarry Smith y_t = ly[(n19 % (m * n)) / m]; 99147c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 99247c6ae99SBarry Smith s_t = bases[n19] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; 9938865f1eaSKarl Rupp if (twod && (s_t >= M * N * P)) s_t = bases[n19] - (s_y + 1 - i) * x_t + x_t * y_t; /* 2d case */ 9948865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 99547c6ae99SBarry Smith } 99647c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 997ce00eea3SSatish Balay x_t = lx[n20 % m]; 99847c6ae99SBarry Smith y_t = ly[(n20 % (m * n)) / m]; 99947c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 100047c6ae99SBarry Smith s_t = bases[n20] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; 10018865f1eaSKarl Rupp if (twod && (s_t >= M * N * P)) s_t = bases[n20] - (s_y + 1 - i) * x_t + x_t * y_t; /* 2d case */ 10028865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 100347c6ae99SBarry Smith } 100447c6ae99SBarry Smith } 100547c6ae99SBarry Smith 100647c6ae99SBarry Smith for (i = 0; i < y; i++) { 100747c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1008ce00eea3SSatish Balay x_t = lx[n21 % m]; 100947c6ae99SBarry Smith y_t = y; 101047c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 101147c6ae99SBarry Smith s_t = bases[n21] + (i + 1) * x_t - s_x + k * x_t * y_t; 10128865f1eaSKarl Rupp if (twod && (s_t >= M * N * P)) s_t = bases[n21] + (i + 1) * x_t - s_x; /* 2d case */ 10138865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 101447c6ae99SBarry Smith } 101547c6ae99SBarry Smith 101647c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 101747c6ae99SBarry Smith x_t = x; 101847c6ae99SBarry Smith y_t = y; 101947c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 102047c6ae99SBarry Smith s_t = bases[n22] + i * x_t + k * x_t * y_t; 10218865f1eaSKarl Rupp if (twod && (s_t >= M * N * P)) s_t = bases[n22] + i * x_t; /* 2d case */ 10228865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 1023bff4a2f0SMatthew G. Knepley } else if (bz == DM_BOUNDARY_MIRROR) { 1024a6fd6b0aSBarry Smith for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + (z - k - 1) * x * y + i * x; 102547c6ae99SBarry Smith } 102647c6ae99SBarry Smith 102747c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1028ce00eea3SSatish Balay x_t = lx[n23 % m]; 102947c6ae99SBarry Smith y_t = y; 103047c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 103147c6ae99SBarry Smith s_t = bases[n23] + i * x_t + k * x_t * y_t; 10328865f1eaSKarl Rupp if (twod && (s_t >= M * N * P)) s_t = bases[n23] + i * x_t; /* 2d case */ 10338865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 103447c6ae99SBarry Smith } 103547c6ae99SBarry Smith } 103647c6ae99SBarry Smith 103747c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 103847c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1039ce00eea3SSatish Balay x_t = lx[n24 % m]; 104047c6ae99SBarry Smith y_t = ly[(n24 % (m * n)) / m]; 104147c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 104247c6ae99SBarry Smith s_t = bases[n24] + i * x_t - s_x + k * x_t * y_t; 10438865f1eaSKarl Rupp if (twod && (s_t >= M * N * P)) s_t = bases[n24] + i * x_t - s_x; /* 2d case */ 10448865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 104547c6ae99SBarry Smith } 104647c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 104747c6ae99SBarry Smith x_t = x; 104847c6ae99SBarry Smith y_t = ly[(n25 % (m * n)) / m]; 104947c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 105047c6ae99SBarry Smith s_t = bases[n25] + (i - 1) * x_t + k * x_t * y_t; 10518865f1eaSKarl Rupp if (twod && (s_t >= M * N * P)) s_t = bases[n25] + (i - 1) * x_t; /* 2d case */ 10528865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 105347c6ae99SBarry Smith } 105447c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1055ce00eea3SSatish Balay x_t = lx[n26 % m]; 105647c6ae99SBarry Smith y_t = ly[(n26 % (m * n)) / m]; 105747c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 105847c6ae99SBarry Smith s_t = bases[n26] + (i - 1) * x_t + k * x_t * y_t; 10598865f1eaSKarl Rupp if (twod && (s_t >= M * N * P)) s_t = bases[n26] + (i - 1) * x_t; /* 2d case */ 10608865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 106147c6ae99SBarry Smith } 106247c6ae99SBarry Smith } 106347c6ae99SBarry Smith } 1064ce00eea3SSatish Balay 10659566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm, dof, nn, idx, PETSC_USE_POINTER, &from)); 10669566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(global, from, local, to, >ol)); 10679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 10689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 106947c6ae99SBarry Smith 1070d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 10719371c9d4SSatish Balay n0 = sn0; 10729371c9d4SSatish Balay n1 = sn1; 10739371c9d4SSatish Balay n2 = sn2; 10749371c9d4SSatish Balay n3 = sn3; 10759371c9d4SSatish Balay n5 = sn5; 10769371c9d4SSatish Balay n6 = sn6; 10779371c9d4SSatish Balay n7 = sn7; 10789371c9d4SSatish Balay n8 = sn8; 10799371c9d4SSatish Balay n9 = sn9; 10809371c9d4SSatish Balay n11 = sn11; 10819371c9d4SSatish Balay n15 = sn15; 10829371c9d4SSatish Balay n17 = sn17; 10839371c9d4SSatish Balay n18 = sn18; 10849371c9d4SSatish Balay n19 = sn19; 10859371c9d4SSatish Balay n20 = sn20; 10869371c9d4SSatish Balay n21 = sn21; 10879371c9d4SSatish Balay n23 = sn23; 10889371c9d4SSatish Balay n24 = sn24; 10899371c9d4SSatish Balay n25 = sn25; 109047c6ae99SBarry Smith n26 = sn26; 1091ce00eea3SSatish Balay } 109247c6ae99SBarry Smith 1093288e7d53SBarry Smith if (((stencil_type == DMDA_STENCIL_STAR) || (bx != DM_BOUNDARY_PERIODIC && bx) || (by != DM_BOUNDARY_PERIODIC && by) || (bz != DM_BOUNDARY_PERIODIC && bz))) { 1094ce00eea3SSatish Balay /* 1095ce00eea3SSatish Balay Recompute the local to global mappings, this time keeping the 1096ce00eea3SSatish Balay information about the cross corner processor numbers. 1097ce00eea3SSatish Balay */ 109847c6ae99SBarry Smith nn = 0; 109947c6ae99SBarry Smith /* Bottom Level */ 110047c6ae99SBarry Smith for (k = 0; k < s_z; k++) { 110147c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 110247c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 1103ce00eea3SSatish Balay x_t = lx[n0 % m]; 110447c6ae99SBarry Smith y_t = ly[(n0 % (m * n)) / m]; 110547c6ae99SBarry Smith z_t = lz[n0 / (m * n)]; 110647c6ae99SBarry 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; 11078865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1108ce00eea3SSatish Balay } else if (Xs - xs < 0 && Ys - ys < 0 && Zs - zs < 0) { 11098865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 111047c6ae99SBarry Smith } 111147c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 111247c6ae99SBarry Smith x_t = x; 111347c6ae99SBarry Smith y_t = ly[(n1 % (m * n)) / m]; 111447c6ae99SBarry Smith z_t = lz[n1 / (m * n)]; 111547c6ae99SBarry 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; 11168865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 1117ce00eea3SSatish Balay } else if (Ys - ys < 0 && Zs - zs < 0) { 11188865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = -1; 111947c6ae99SBarry Smith } 112047c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 1121ce00eea3SSatish Balay x_t = lx[n2 % m]; 112247c6ae99SBarry Smith y_t = ly[(n2 % (m * n)) / m]; 112347c6ae99SBarry Smith z_t = lz[n2 / (m * n)]; 112447c6ae99SBarry 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; 11258865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1126ce00eea3SSatish Balay } else if (xe - Xe < 0 && Ys - ys < 0 && Zs - zs < 0) { 11278865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 112847c6ae99SBarry Smith } 112947c6ae99SBarry Smith } 113047c6ae99SBarry Smith 113147c6ae99SBarry Smith for (i = 0; i < y; i++) { 113247c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 1133ce00eea3SSatish Balay x_t = lx[n3 % m]; 113447c6ae99SBarry Smith y_t = y; 113547c6ae99SBarry Smith z_t = lz[n3 / (m * n)]; 113647c6ae99SBarry Smith s_t = bases[n3] + (i + 1) * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t; 11378865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1138ce00eea3SSatish Balay } else if (Xs - xs < 0 && Zs - zs < 0) { 11398865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 114047c6ae99SBarry Smith } 114147c6ae99SBarry Smith 114247c6ae99SBarry Smith if (n4 >= 0) { /* middle */ 114347c6ae99SBarry Smith x_t = x; 114447c6ae99SBarry Smith y_t = y; 114547c6ae99SBarry Smith z_t = lz[n4 / (m * n)]; 114647c6ae99SBarry Smith s_t = bases[n4] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; 11478865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 1148ce00eea3SSatish Balay } else if (Zs - zs < 0) { 1149bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 1150a6fd6b0aSBarry Smith for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + x * i + (s_z - k - 1) * x * y; 1151c2859e5eSBarry Smith } else { 11528865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = -1; 115347c6ae99SBarry Smith } 1154c2859e5eSBarry Smith } 115547c6ae99SBarry Smith 115647c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 1157ce00eea3SSatish Balay x_t = lx[n5 % m]; 115847c6ae99SBarry Smith y_t = y; 115947c6ae99SBarry Smith z_t = lz[n5 / (m * n)]; 116047c6ae99SBarry Smith s_t = bases[n5] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; 11618865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1162ce00eea3SSatish Balay } else if (xe - Xe < 0 && Zs - zs < 0) { 11638865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 116447c6ae99SBarry Smith } 116547c6ae99SBarry Smith } 116647c6ae99SBarry Smith 116747c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 116847c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 1169ce00eea3SSatish Balay x_t = lx[n6 % m]; 117047c6ae99SBarry Smith y_t = ly[(n6 % (m * n)) / m]; 117147c6ae99SBarry Smith z_t = lz[n6 / (m * n)]; 117247c6ae99SBarry Smith s_t = bases[n6] + i * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t; 11738865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1174ce00eea3SSatish Balay } else if (Xs - xs < 0 && ye - Ye < 0 && Zs - zs < 0) { 11758865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 117647c6ae99SBarry Smith } 117747c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 117847c6ae99SBarry Smith x_t = x; 117947c6ae99SBarry Smith y_t = ly[(n7 % (m * n)) / m]; 118047c6ae99SBarry Smith z_t = lz[n7 / (m * n)]; 118147c6ae99SBarry Smith s_t = bases[n7] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; 11828865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 1183ce00eea3SSatish Balay } else if (ye - Ye < 0 && Zs - zs < 0) { 11848865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = -1; 118547c6ae99SBarry Smith } 118647c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 1187ce00eea3SSatish Balay x_t = lx[n8 % m]; 118847c6ae99SBarry Smith y_t = ly[(n8 % (m * n)) / m]; 118947c6ae99SBarry Smith z_t = lz[n8 / (m * n)]; 119047c6ae99SBarry Smith s_t = bases[n8] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; 11918865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1192ce00eea3SSatish Balay } else if (xe - Xe < 0 && ye - Ye < 0 && Zs - zs < 0) { 11938865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 119447c6ae99SBarry Smith } 119547c6ae99SBarry Smith } 119647c6ae99SBarry Smith } 119747c6ae99SBarry Smith 119847c6ae99SBarry Smith /* Middle Level */ 119947c6ae99SBarry Smith for (k = 0; k < z; k++) { 120047c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 120147c6ae99SBarry Smith if (n9 >= 0) { /* left below */ 1202ce00eea3SSatish Balay x_t = lx[n9 % m]; 120347c6ae99SBarry Smith y_t = ly[(n9 % (m * n)) / m]; 120447c6ae99SBarry Smith /* z_t = z; */ 120547c6ae99SBarry Smith s_t = bases[n9] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t; 12068865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1207ce00eea3SSatish Balay } else if (Xs - xs < 0 && Ys - ys < 0) { 12088865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 120947c6ae99SBarry Smith } 121047c6ae99SBarry Smith if (n10 >= 0) { /* directly below */ 121147c6ae99SBarry Smith x_t = x; 121247c6ae99SBarry Smith y_t = ly[(n10 % (m * n)) / m]; 121347c6ae99SBarry Smith /* z_t = z; */ 121447c6ae99SBarry Smith s_t = bases[n10] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; 12158865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 1216ce00eea3SSatish Balay } else if (Ys - ys < 0) { 1217bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 1218feb237baSPierre Jolivet for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (s_y - i) * x; 1219c2859e5eSBarry Smith } else { 12208865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = -1; 1221c2859e5eSBarry Smith } 122247c6ae99SBarry Smith } 122347c6ae99SBarry Smith if (n11 >= 0) { /* right below */ 1224ce00eea3SSatish Balay x_t = lx[n11 % m]; 122547c6ae99SBarry Smith y_t = ly[(n11 % (m * n)) / m]; 122647c6ae99SBarry Smith /* z_t = z; */ 122747c6ae99SBarry Smith s_t = bases[n11] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; 12288865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1229ce00eea3SSatish Balay } else if (xe - Xe < 0 && Ys - ys < 0) { 12308865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 123147c6ae99SBarry Smith } 123247c6ae99SBarry Smith } 123347c6ae99SBarry Smith 123447c6ae99SBarry Smith for (i = 0; i < y; i++) { 123547c6ae99SBarry Smith if (n12 >= 0) { /* directly left */ 1236ce00eea3SSatish Balay x_t = lx[n12 % m]; 123747c6ae99SBarry Smith y_t = y; 123847c6ae99SBarry Smith /* z_t = z; */ 123947c6ae99SBarry Smith s_t = bases[n12] + (i + 1) * x_t - s_x + k * x_t * y_t; 12408865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1241ce00eea3SSatish Balay } else if (Xs - xs < 0) { 1242bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 1243a6fd6b0aSBarry Smith for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k * x * y + i * x; 1244c2859e5eSBarry Smith } else { 12458865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 124647c6ae99SBarry Smith } 1247c2859e5eSBarry Smith } 124847c6ae99SBarry Smith 124947c6ae99SBarry Smith /* Interior */ 125047c6ae99SBarry Smith s_t = bases[rank] + i * x + k * x * y; 12518865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = s_t++; 125247c6ae99SBarry Smith 125347c6ae99SBarry Smith if (n14 >= 0) { /* directly right */ 1254ce00eea3SSatish Balay x_t = lx[n14 % m]; 125547c6ae99SBarry Smith y_t = y; 125647c6ae99SBarry Smith /* z_t = z; */ 125747c6ae99SBarry Smith s_t = bases[n14] + i * x_t + k * x_t * y_t; 12588865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1259ce00eea3SSatish Balay } else if (xe - Xe < 0) { 1260bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 1261a6fd6b0aSBarry Smith for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k * x * y + i * x; 1262c2859e5eSBarry Smith } else { 12638865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 126447c6ae99SBarry Smith } 126547c6ae99SBarry Smith } 1266c2859e5eSBarry Smith } 126747c6ae99SBarry Smith 126847c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 126947c6ae99SBarry Smith if (n15 >= 0) { /* left above */ 1270ce00eea3SSatish Balay x_t = lx[n15 % m]; 127147c6ae99SBarry Smith y_t = ly[(n15 % (m * n)) / m]; 127247c6ae99SBarry Smith /* z_t = z; */ 127347c6ae99SBarry Smith s_t = bases[n15] + i * x_t - s_x + k * x_t * y_t; 12748865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1275ce00eea3SSatish Balay } else if (Xs - xs < 0 && ye - Ye < 0) { 12768865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 127747c6ae99SBarry Smith } 127847c6ae99SBarry Smith if (n16 >= 0) { /* directly above */ 127947c6ae99SBarry Smith x_t = x; 128047c6ae99SBarry Smith y_t = ly[(n16 % (m * n)) / m]; 128147c6ae99SBarry Smith /* z_t = z; */ 128247c6ae99SBarry Smith s_t = bases[n16] + (i - 1) * x_t + k * x_t * y_t; 12838865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 1284ce00eea3SSatish Balay } else if (ye - Ye < 0) { 1285bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 1286a6fd6b0aSBarry Smith for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (y - i) * x; 1287c2859e5eSBarry Smith } else { 12888865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = -1; 128947c6ae99SBarry Smith } 1290c2859e5eSBarry Smith } 129147c6ae99SBarry Smith if (n17 >= 0) { /* right above */ 1292ce00eea3SSatish Balay x_t = lx[n17 % m]; 129347c6ae99SBarry Smith y_t = ly[(n17 % (m * n)) / m]; 129447c6ae99SBarry Smith /* z_t = z; */ 129547c6ae99SBarry Smith s_t = bases[n17] + (i - 1) * x_t + k * x_t * y_t; 12968865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1297ce00eea3SSatish Balay } else if (xe - Xe < 0 && ye - Ye < 0) { 12988865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 129947c6ae99SBarry Smith } 130047c6ae99SBarry Smith } 130147c6ae99SBarry Smith } 130247c6ae99SBarry Smith 130347c6ae99SBarry Smith /* Upper Level */ 130447c6ae99SBarry Smith for (k = 0; k < s_z; k++) { 130547c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 130647c6ae99SBarry Smith if (n18 >= 0) { /* left below */ 1307ce00eea3SSatish Balay x_t = lx[n18 % m]; 130847c6ae99SBarry Smith y_t = ly[(n18 % (m * n)) / m]; 130947c6ae99SBarry Smith /* z_t = lz[n18 / (m*n)]; */ 131047c6ae99SBarry Smith s_t = bases[n18] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t; 13118865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1312ce00eea3SSatish Balay } else if (Xs - xs < 0 && Ys - ys < 0 && ze - Ze < 0) { 13138865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 131447c6ae99SBarry Smith } 131547c6ae99SBarry Smith if (n19 >= 0) { /* directly below */ 131647c6ae99SBarry Smith x_t = x; 131747c6ae99SBarry Smith y_t = ly[(n19 % (m * n)) / m]; 131847c6ae99SBarry Smith /* z_t = lz[n19 / (m*n)]; */ 131947c6ae99SBarry Smith s_t = bases[n19] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; 13208865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 1321ce00eea3SSatish Balay } else if (Ys - ys < 0 && ze - Ze < 0) { 13228865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = -1; 132347c6ae99SBarry Smith } 132447c6ae99SBarry Smith if (n20 >= 0) { /* right below */ 1325ce00eea3SSatish Balay x_t = lx[n20 % m]; 132647c6ae99SBarry Smith y_t = ly[(n20 % (m * n)) / m]; 132747c6ae99SBarry Smith /* z_t = lz[n20 / (m*n)]; */ 132847c6ae99SBarry Smith s_t = bases[n20] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; 13298865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1330ce00eea3SSatish Balay } else if (xe - Xe < 0 && Ys - ys < 0 && ze - Ze < 0) { 13318865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 133247c6ae99SBarry Smith } 133347c6ae99SBarry Smith } 133447c6ae99SBarry Smith 133547c6ae99SBarry Smith for (i = 0; i < y; i++) { 133647c6ae99SBarry Smith if (n21 >= 0) { /* directly left */ 1337ce00eea3SSatish Balay x_t = lx[n21 % m]; 133847c6ae99SBarry Smith y_t = y; 133947c6ae99SBarry Smith /* z_t = lz[n21 / (m*n)]; */ 134047c6ae99SBarry Smith s_t = bases[n21] + (i + 1) * x_t - s_x + k * x_t * y_t; 13418865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1342ce00eea3SSatish Balay } else if (Xs - xs < 0 && ze - Ze < 0) { 13438865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 134447c6ae99SBarry Smith } 134547c6ae99SBarry Smith 134647c6ae99SBarry Smith if (n22 >= 0) { /* middle */ 134747c6ae99SBarry Smith x_t = x; 134847c6ae99SBarry Smith y_t = y; 134947c6ae99SBarry Smith /* z_t = lz[n22 / (m*n)]; */ 135047c6ae99SBarry Smith s_t = bases[n22] + i * x_t + k * x_t * y_t; 13518865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 1352ce00eea3SSatish Balay } else if (ze - Ze < 0) { 1353bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_MIRROR) { 1354a6fd6b0aSBarry Smith for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + (z - k - 1) * x * y + i * x; 1355c2859e5eSBarry Smith } else { 13568865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = -1; 135747c6ae99SBarry Smith } 1358c2859e5eSBarry Smith } 135947c6ae99SBarry Smith 136047c6ae99SBarry Smith if (n23 >= 0) { /* directly right */ 1361ce00eea3SSatish Balay x_t = lx[n23 % m]; 136247c6ae99SBarry Smith y_t = y; 136347c6ae99SBarry Smith /* z_t = lz[n23 / (m*n)]; */ 136447c6ae99SBarry Smith s_t = bases[n23] + i * x_t + k * x_t * y_t; 13658865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1366ce00eea3SSatish Balay } else if (xe - Xe < 0 && ze - Ze < 0) { 13678865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 136847c6ae99SBarry Smith } 136947c6ae99SBarry Smith } 137047c6ae99SBarry Smith 137147c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 137247c6ae99SBarry Smith if (n24 >= 0) { /* left above */ 1373ce00eea3SSatish Balay x_t = lx[n24 % m]; 137447c6ae99SBarry Smith y_t = ly[(n24 % (m * n)) / m]; 137547c6ae99SBarry Smith /* z_t = lz[n24 / (m*n)]; */ 137647c6ae99SBarry Smith s_t = bases[n24] + i * x_t - s_x + k * x_t * y_t; 13778865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1378ce00eea3SSatish Balay } else if (Xs - xs < 0 && ye - Ye < 0 && ze - Ze < 0) { 13798865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 138047c6ae99SBarry Smith } 138147c6ae99SBarry Smith if (n25 >= 0) { /* directly above */ 138247c6ae99SBarry Smith x_t = x; 138347c6ae99SBarry Smith y_t = ly[(n25 % (m * n)) / m]; 138447c6ae99SBarry Smith /* z_t = lz[n25 / (m*n)]; */ 138547c6ae99SBarry Smith s_t = bases[n25] + (i - 1) * x_t + k * x_t * y_t; 13868865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 1387ce00eea3SSatish Balay } else if (ye - Ye < 0 && ze - Ze < 0) { 13888865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = -1; 138947c6ae99SBarry Smith } 139047c6ae99SBarry Smith if (n26 >= 0) { /* right above */ 1391ce00eea3SSatish Balay x_t = lx[n26 % m]; 139247c6ae99SBarry Smith y_t = ly[(n26 % (m * n)) / m]; 139347c6ae99SBarry Smith /* z_t = lz[n26 / (m*n)]; */ 139447c6ae99SBarry Smith s_t = bases[n26] + (i - 1) * x_t + k * x_t * y_t; 13958865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 1396ce00eea3SSatish Balay } else if (xe - Xe < 0 && ye - Ye < 0 && ze - Ze < 0) { 13978865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 139847c6ae99SBarry Smith } 139947c6ae99SBarry Smith } 140047c6ae99SBarry Smith } 140147c6ae99SBarry Smith } 140247c6ae99SBarry Smith /* 140347c6ae99SBarry Smith Set the local to global ordering in the global vector, this allows use 140447c6ae99SBarry Smith of VecSetValuesLocal(). 140547c6ae99SBarry Smith */ 14069566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, dof, nn, idx, PETSC_OWN_POINTER, &da->ltogmap)); 140747c6ae99SBarry Smith 14089566063dSJacob Faibussowitsch PetscCall(PetscFree2(bases, ldims)); 14099371c9d4SSatish Balay dd->m = m; 14109371c9d4SSatish Balay dd->n = n; 14119371c9d4SSatish Balay dd->p = p; 1412ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 14139371c9d4SSatish Balay dd->xs = xs * dof; 14149371c9d4SSatish Balay dd->xe = xe * dof; 14159371c9d4SSatish Balay dd->ys = ys; 14169371c9d4SSatish Balay dd->ye = ye; 14179371c9d4SSatish Balay dd->zs = zs; 14189371c9d4SSatish Balay dd->ze = ze; 14199371c9d4SSatish Balay dd->Xs = Xs * dof; 14209371c9d4SSatish Balay dd->Xe = Xe * dof; 14219371c9d4SSatish Balay dd->Ys = Ys; 14229371c9d4SSatish Balay dd->Ye = Ye; 14239371c9d4SSatish Balay dd->Zs = Zs; 14249371c9d4SSatish Balay dd->Ze = Ze; 1425ce00eea3SSatish Balay 14269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 14279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 1428ce00eea3SSatish Balay 1429ce00eea3SSatish Balay dd->gtol = gtol; 1430ce00eea3SSatish Balay dd->base = base; 1431ce00eea3SSatish Balay da->ops->view = DMView_DA_3d; 14320298fd71SBarry Smith dd->ltol = NULL; 14330298fd71SBarry Smith dd->ao = NULL; 14343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 143547c6ae99SBarry Smith } 1436564755cdSBarry Smith 143747c6ae99SBarry Smith /*@C 1438aa219208SBarry Smith DMDACreate3d - Creates an object that will manage the communication of three-dimensional 143947c6ae99SBarry Smith regular array data that is distributed across some processors. 144047c6ae99SBarry Smith 1441d083f849SBarry Smith Collective 144247c6ae99SBarry Smith 144347c6ae99SBarry Smith Input Parameters: 144447c6ae99SBarry Smith + comm - MPI communicator 1445*a4e35b19SJacob Faibussowitsch . bx - type of x ghost nodes the array have. 1446*a4e35b19SJacob Faibussowitsch Use one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`. 1447*a4e35b19SJacob Faibussowitsch . by - type of y ghost nodes the array have. 1448*a4e35b19SJacob Faibussowitsch Use one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`. 1449*a4e35b19SJacob Faibussowitsch . bz - type of z ghost nodes the array have. 1450dce8aebaSBarry Smith Use one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`. 1451dce8aebaSBarry Smith . stencil_type - Type of stencil (`DMDA_STENCIL_STAR` or `DMDA_STENCIL_BOX`) 1452*a4e35b19SJacob Faibussowitsch . M - global dimension in x direction of the array 1453*a4e35b19SJacob Faibussowitsch . N - global dimension in y direction of the array 1454*a4e35b19SJacob Faibussowitsch . P - global dimension in z direction of the array 1455*a4e35b19SJacob Faibussowitsch . m - corresponding number of processors in x dimension (or `PETSC_DECIDE` to have calculated) 1456*a4e35b19SJacob Faibussowitsch . n - corresponding number of processors in y dimension (or `PETSC_DECIDE` to have calculated) 1457*a4e35b19SJacob Faibussowitsch . p - corresponding number of processors in z dimension (or `PETSC_DECIDE` to have calculated) 145847c6ae99SBarry Smith . dof - number of degrees of freedom per node 145910d7c030SMatthew G Knepley . s - stencil width 1460*a4e35b19SJacob Faibussowitsch . lx - arrays containing the number of nodes in each cell along the x coordinates, or `NULL`. 1461*a4e35b19SJacob Faibussowitsch . ly - arrays containing the number of nodes in each cell along the y coordinates, or `NULL`. 1462*a4e35b19SJacob Faibussowitsch - lz - arrays containing the number of nodes in each cell along the z coordinates, or `NULL`. 146347c6ae99SBarry Smith 146447c6ae99SBarry Smith Output Parameter: 146547c6ae99SBarry Smith . da - the resulting distributed array object 146647c6ae99SBarry Smith 1467dce8aebaSBarry Smith Options Database Keys: 1468dce8aebaSBarry Smith + -dm_view - Calls `DMView()` at the conclusion of `DMDACreate3d()` 1469e43dc8daSBarry Smith . -da_grid_x <nx> - number of grid points in x direction 1470e43dc8daSBarry Smith . -da_grid_y <ny> - number of grid points in y direction 1471e43dc8daSBarry Smith . -da_grid_z <nz> - number of grid points in z direction 147247c6ae99SBarry Smith . -da_processors_x <MX> - number of processors in x direction 147347c6ae99SBarry Smith . -da_processors_y <MY> - number of processors in y direction 147447c6ae99SBarry Smith . -da_processors_z <MZ> - number of processors in z direction 1475e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 1476e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 1477e0f5d30fSBarry Smith . -da_refine_z <rz> - refinement ratio in z directio 1478e43dc8daSBarry Smith - -da_refine <n> - refine the DMDA n times before creating it 147947c6ae99SBarry Smith 148047c6ae99SBarry Smith Level: beginner 148147c6ae99SBarry Smith 148247c6ae99SBarry Smith Notes: 1483*a4e35b19SJacob Faibussowitsch If `lx`, `ly`, or `lz` are non-null, these must be of length as `m`, `n`, `p` and the 1484*a4e35b19SJacob Faibussowitsch corresponding `m`, `n`, or `p` cannot be `PETSC_DECIDE`. Sum of the `lx` entries must be `M`, 1485*a4e35b19SJacob Faibussowitsch sum of the `ly` must `N`, sum of the `lz` must be `P`. 1486*a4e35b19SJacob Faibussowitsch 1487dce8aebaSBarry Smith The stencil type `DMDA_STENCIL_STAR` with width 1 corresponds to the 1488dce8aebaSBarry Smith standard 7-pt stencil, while `DMDA_STENCIL_BOX` with width 1 denotes 148947c6ae99SBarry Smith the standard 27-pt stencil. 149047c6ae99SBarry Smith 1491dce8aebaSBarry Smith The array data itself is NOT stored in the `DMDA`, it is stored in `Vec` objects; 1492dce8aebaSBarry Smith The appropriate vector objects can be obtained with calls to `DMCreateGlobalVector()` 1493dce8aebaSBarry Smith and `DMCreateLocalVector()` and calls to `VecDuplicate()` if more are needed. 149447c6ae99SBarry Smith 1495dce8aebaSBarry Smith You must call `DMSetUp()` after this call before using this `DM`. 1496897f7067SBarry Smith 1497dce8aebaSBarry Smith If you wish to use the options database to change values in the `DMDA` call `DMSetFromOptions()` after this call 1498dce8aebaSBarry Smith but before `DMSetUp()`. 1499897f7067SBarry Smith 1500dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDestroy()`, `DMView()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMGlobalToLocalBegin()`, `DMDAGetRefinementFactor()`, 1501db781477SPatrick Sanan `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDASetRefinementFactor()`, 1502db781477SPatrick Sanan `DMDAGetInfo()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMDACreateNaturalVector()`, `DMLoad()`, `DMDAGetOwnershipRanges()`, 1503db781477SPatrick Sanan `DMStagCreate3d()` 150447c6ae99SBarry Smith @*/ 1505d71ae5a4SJacob 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) 1506d71ae5a4SJacob Faibussowitsch { 150747c6ae99SBarry Smith PetscFunctionBegin; 15089566063dSJacob Faibussowitsch PetscCall(DMDACreate(comm, da)); 15099566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*da, 3)); 15109566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*da, M, N, P)); 15119566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*da, m, n, p)); 15129566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*da, bx, by, bz)); 15139566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*da, dof)); 15149566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*da, stencil_type)); 15159566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*da, s)); 15169566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*da, lx, ly, lz)); 15173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 151847c6ae99SBarry Smith } 1519