19a42bb27SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 39804daf3SBarry Smith #include <petscdraw.h> 447c6ae99SBarry Smith 59371c9d4SSatish Balay static PetscErrorCode DMView_DA_2d(DM da, PetscViewer viewer) { 647c6ae99SBarry Smith PetscMPIInt rank; 7c9493c35SStefano Zampini PetscBool iascii, isdraw, isglvis, isbinary; 847c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 99a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 109a42bb27SBarry Smith PetscBool ismatlab; 119a42bb27SBarry Smith #endif 1247c6ae99SBarry Smith 1347c6ae99SBarry Smith PetscFunctionBegin; 149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)da), &rank)); 1547c6ae99SBarry Smith 169566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 179566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 189566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis)); 199566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 209a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 219566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERMATLAB, &ismatlab)); 229a42bb27SBarry Smith #endif 2347c6ae99SBarry Smith if (iascii) { 2447c6ae99SBarry Smith PetscViewerFormat format; 2547c6ae99SBarry Smith 269566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 2776a8abe0SBarry Smith if (format == PETSC_VIEWER_LOAD_BALANCE) { 2876a8abe0SBarry Smith PetscInt i, nmax = 0, nmin = PETSC_MAX_INT, navg = 0, *nz, nzlocal; 2976a8abe0SBarry Smith DMDALocalInfo info; 3076a8abe0SBarry Smith PetscMPIInt size; 319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da), &size)); 329566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 3376a8abe0SBarry Smith nzlocal = info.xm * info.ym; 349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &nz)); 359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nzlocal, 1, MPIU_INT, nz, 1, MPIU_INT, PetscObjectComm((PetscObject)da))); 3676a8abe0SBarry Smith for (i = 0; i < (PetscInt)size; i++) { 3776a8abe0SBarry Smith nmax = PetscMax(nmax, nz[i]); 3876a8abe0SBarry Smith nmin = PetscMin(nmin, nz[i]); 3976a8abe0SBarry Smith navg += nz[i]; 4076a8abe0SBarry Smith } 419566063dSJacob Faibussowitsch PetscCall(PetscFree(nz)); 4276a8abe0SBarry Smith navg = navg / size; 4363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Load Balance - Grid Points: Min %" PetscInt_FMT " avg %" PetscInt_FMT " max %" PetscInt_FMT "\n", nmin, navg, nmax)); 4476a8abe0SBarry Smith PetscFunctionReturn(0); 4576a8abe0SBarry Smith } 468ec8862eSJed Brown if (format != PETSC_VIEWER_ASCII_VTK_DEPRECATED && format != PETSC_VIEWER_ASCII_VTK_CELL_DEPRECATED && format != PETSC_VIEWER_ASCII_GLVIS) { 47aa219208SBarry Smith DMDALocalInfo info; 489566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 5063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Processor [%d] M %" PetscInt_FMT " N %" PetscInt_FMT " m %" PetscInt_FMT " n %" PetscInt_FMT " w %" PetscInt_FMT " s %" PetscInt_FMT "\n", rank, dd->M, dd->N, dd->m, dd->n, dd->w, dd->s)); 5163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "X range of indices: %" PetscInt_FMT " %" PetscInt_FMT ", Y range of indices: %" PetscInt_FMT " %" PetscInt_FMT "\n", info.xs, info.xs + info.xm, info.ys, info.ys + info.ym)); 529566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 541baa6e33SBarry Smith } else if (format == PETSC_VIEWER_ASCII_GLVIS) PetscCall(DMView_DA_GLVis(da, viewer)); 551baa6e33SBarry Smith else PetscCall(DMView_DA_VTK(da, viewer)); 5647c6ae99SBarry Smith } else if (isdraw) { 5747c6ae99SBarry Smith PetscDraw draw; 5847c6ae99SBarry Smith double ymin = -1 * dd->s - 1, ymax = dd->N + dd->s; 5947c6ae99SBarry Smith double xmin = -1 * dd->s - 1, xmax = dd->M + dd->s; 6047c6ae99SBarry Smith double x, y; 618ea3bf28SBarry Smith PetscInt base; 628ea3bf28SBarry Smith const PetscInt *idx; 6347c6ae99SBarry Smith char node[10]; 6447c6ae99SBarry Smith PetscBool isnull; 6547c6ae99SBarry Smith 669566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 679566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw, &isnull)); 685b399a63SLisandro Dalcin if (isnull) PetscFunctionReturn(0); 6947c6ae99SBarry Smith 709566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 719566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 729566063dSJacob Faibussowitsch PetscCall(PetscDrawSetCoordinates(draw, xmin, ymin, xmax, ymax)); 735b399a63SLisandro Dalcin 74d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 7547c6ae99SBarry Smith /* first processor draw all node lines */ 76dd400576SPatrick Sanan if (rank == 0) { 779371c9d4SSatish Balay ymin = 0.0; 789371c9d4SSatish Balay ymax = dd->N - 1; 79*48a46eb9SPierre Jolivet for (xmin = 0; xmin < dd->M; xmin++) PetscCall(PetscDrawLine(draw, xmin, ymin, xmin, ymax, PETSC_DRAW_BLACK)); 809371c9d4SSatish Balay xmin = 0.0; 819371c9d4SSatish Balay xmax = dd->M - 1; 82*48a46eb9SPierre Jolivet for (ymin = 0; ymin < dd->N; ymin++) PetscCall(PetscDrawLine(draw, xmin, ymin, xmax, ymin, PETSC_DRAW_BLACK)); 8347c6ae99SBarry Smith } 84d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 859566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 869566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 8747c6ae99SBarry Smith 88d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 8947c6ae99SBarry Smith /* draw my box */ 909371c9d4SSatish Balay xmin = dd->xs / dd->w; 919371c9d4SSatish Balay xmax = (dd->xe - 1) / dd->w; 929371c9d4SSatish Balay ymin = dd->ys; 939371c9d4SSatish Balay ymax = dd->ye - 1; 949566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymin, xmax, ymin, PETSC_DRAW_RED)); 959566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymin, xmin, ymax, PETSC_DRAW_RED)); 969566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmin, ymax, xmax, ymax, PETSC_DRAW_RED)); 979566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xmax, ymin, xmax, ymax, PETSC_DRAW_RED)); 9847c6ae99SBarry Smith /* put in numbers */ 9947c6ae99SBarry Smith base = (dd->base) / dd->w; 10047c6ae99SBarry Smith for (y = ymin; y <= ymax; y++) { 10147c6ae99SBarry Smith for (x = xmin; x <= xmax; x++) { 1029566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node, sizeof(node), "%d", (int)base++)); 1039566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, x, y, PETSC_DRAW_BLACK, node)); 10447c6ae99SBarry Smith } 10547c6ae99SBarry Smith } 106d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1079566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1089566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 10947c6ae99SBarry Smith 110d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 1115b399a63SLisandro Dalcin /* overlay ghost numbers, useful for error checking */ 1129566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockIndices(da->ltogmap, &idx)); 1139371c9d4SSatish Balay base = 0; 1149371c9d4SSatish Balay xmin = dd->Xs; 1159371c9d4SSatish Balay xmax = dd->Xe; 1169371c9d4SSatish Balay ymin = dd->Ys; 1179371c9d4SSatish Balay ymax = dd->Ye; 11847c6ae99SBarry Smith for (y = ymin; y < ymax; y++) { 11947c6ae99SBarry Smith for (x = xmin; x < xmax; x++) { 12047c6ae99SBarry Smith if ((base % dd->w) == 0) { 1219566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(node, sizeof(node), "%d", (int)(idx[base / dd->w]))); 1229566063dSJacob Faibussowitsch PetscCall(PetscDrawString(draw, x / dd->w, y, PETSC_DRAW_BLUE, node)); 12347c6ae99SBarry Smith } 12447c6ae99SBarry Smith base++; 12547c6ae99SBarry Smith } 12647c6ae99SBarry Smith } 1279566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap, &idx)); 128d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 1299566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1309566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 1319566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 132c9493c35SStefano Zampini } else if (isglvis) { 1339566063dSJacob Faibussowitsch PetscCall(DMView_DA_GLVis(da, viewer)); 1349a42bb27SBarry Smith } else if (isbinary) { 1359566063dSJacob Faibussowitsch PetscCall(DMView_DA_Binary(da, viewer)); 1369a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 1379a42bb27SBarry Smith } else if (ismatlab) { 1389566063dSJacob Faibussowitsch PetscCall(DMView_DA_Matlab(da, viewer)); 1399a42bb27SBarry Smith #endif 14011aeaf0aSBarry Smith } 14147c6ae99SBarry Smith PetscFunctionReturn(0); 14247c6ae99SBarry Smith } 14347c6ae99SBarry Smith 14447c6ae99SBarry Smith #if defined(new) 14547c6ae99SBarry Smith /* 146aa219208SBarry Smith DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local 147aa219208SBarry Smith function lives on a DMDA 14847c6ae99SBarry Smith 14947c6ae99SBarry Smith y ~= (F(u + ha) - F(u))/h, 15047c6ae99SBarry Smith where F = nonlinear function, as set by SNESSetFunction() 15147c6ae99SBarry Smith u = current iterate 15247c6ae99SBarry Smith h = difference interval 15347c6ae99SBarry Smith */ 1549371c9d4SSatish Balay PetscErrorCode DMDAGetDiagonal_MFFD(DM da, Vec U, Vec a) { 15547c6ae99SBarry Smith PetscScalar h, *aa, *ww, v; 15647c6ae99SBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON, umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON; 15747c6ae99SBarry Smith PetscInt gI, nI; 15847c6ae99SBarry Smith MatStencil stencil; 159aa219208SBarry Smith DMDALocalInfo info; 16047c6ae99SBarry Smith 16147c6ae99SBarry Smith PetscFunctionBegin; 1629566063dSJacob Faibussowitsch PetscCall((*ctx->func)(0, U, a, ctx->funcctx)); 1639566063dSJacob Faibussowitsch PetscCall((*ctx->funcisetbase)(U, ctx->funcctx)); 16447c6ae99SBarry Smith 1659566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &ww)); 1669566063dSJacob Faibussowitsch PetscCall(VecGetArray(a, &aa)); 16747c6ae99SBarry Smith 16847c6ae99SBarry Smith nI = 0; 16947c6ae99SBarry Smith h = ww[gI]; 17047c6ae99SBarry Smith if (h == 0.0) h = 1.0; 17147c6ae99SBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 17247c6ae99SBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 17347c6ae99SBarry Smith h *= epsilon; 17447c6ae99SBarry Smith 17547c6ae99SBarry Smith ww[gI] += h; 1769566063dSJacob Faibussowitsch PetscCall((*ctx->funci)(i, w, &v, ctx->funcctx)); 17747c6ae99SBarry Smith aa[nI] = (v - aa[nI]) / h; 17847c6ae99SBarry Smith ww[gI] -= h; 17947c6ae99SBarry Smith nI++; 1808865f1eaSKarl Rupp 1819566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &ww)); 1829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a, &aa)); 18347c6ae99SBarry Smith PetscFunctionReturn(0); 18447c6ae99SBarry Smith } 18547c6ae99SBarry Smith #endif 18647c6ae99SBarry Smith 1879371c9d4SSatish Balay PetscErrorCode DMSetUp_DA_2D(DM da) { 18847c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 18947c6ae99SBarry Smith const PetscInt M = dd->M; 19047c6ae99SBarry Smith const PetscInt N = dd->N; 19147c6ae99SBarry Smith PetscInt m = dd->m; 19247c6ae99SBarry Smith PetscInt n = dd->n; 19347c6ae99SBarry Smith const PetscInt dof = dd->w; 19447c6ae99SBarry Smith const PetscInt s = dd->s; 195bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 196bff4a2f0SMatthew G. Knepley DMBoundaryType by = dd->by; 19719fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 19847c6ae99SBarry Smith PetscInt *lx = dd->lx; 19947c6ae99SBarry Smith PetscInt *ly = dd->ly; 20047c6ae99SBarry Smith MPI_Comm comm; 20147c6ae99SBarry Smith PetscMPIInt rank, size; 202bd1fc5aeSBarry Smith PetscInt xs, xe, ys, ye, x, y, Xs, Xe, Ys, Ye, IXs, IXe, IYs, IYe; 2038ea3bf28SBarry Smith PetscInt up, down, left, right, i, n0, n1, n2, n3, n5, n6, n7, n8, *idx, nn; 204db87c5ecSEthan Coon PetscInt xbase, *bases, *ldims, j, x_t, y_t, s_t, base, count; 20547c6ae99SBarry Smith PetscInt s_x, s_y; /* s proportionalized to w */ 20647c6ae99SBarry Smith PetscInt sn0 = 0, sn2 = 0, sn6 = 0, sn8 = 0; 20747c6ae99SBarry Smith Vec local, global; 208bd1fc5aeSBarry Smith VecScatter gtol; 20945b6f7e9SBarry Smith IS to, from; 21047c6ae99SBarry Smith 21147c6ae99SBarry Smith PetscFunctionBegin; 2127a8be351SBarry Smith PetscCheck(stencil_type != DMDA_STENCIL_BOX || (bx != DM_BOUNDARY_MIRROR && by != DM_BOUNDARY_MIRROR), PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Mirror boundary and box stencil"); 2139566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 2143855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 21563a3b9bcSJacob Faibussowitsch PetscCheck(((PetscInt64)M) * ((PetscInt64)N) * ((PetscInt64)dof) <= (PetscInt64)PETSC_MPI_INT_MAX, comm, PETSC_ERR_INT_OVERFLOW, "Mesh of %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT " (dof) is too large for 32 bit indices", M, N, dof); 2163855c12bSBarry Smith #endif 21747c6ae99SBarry Smith 2189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 22047c6ae99SBarry Smith 2217d310018SBarry Smith dd->p = 1; 22247c6ae99SBarry Smith if (m != PETSC_DECIDE) { 22363a3b9bcSJacob Faibussowitsch PetscCheck(m >= 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in X direction: %" PetscInt_FMT, m); 224f7d195e4SLawrence Mitchell PetscCheck(m <= size, comm, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in X direction: %" PetscInt_FMT " %d", m, size); 22547c6ae99SBarry Smith } 22647c6ae99SBarry Smith if (n != PETSC_DECIDE) { 22763a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in Y direction: %" PetscInt_FMT, n); 228f7d195e4SLawrence Mitchell PetscCheck(n <= size, comm, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in Y direction: %" PetscInt_FMT " %d", n, size); 22947c6ae99SBarry Smith } 23047c6ae99SBarry Smith 23147c6ae99SBarry Smith if (m == PETSC_DECIDE || n == PETSC_DECIDE) { 23247c6ae99SBarry Smith if (n != PETSC_DECIDE) { 23347c6ae99SBarry Smith m = size / n; 23447c6ae99SBarry Smith } else if (m != PETSC_DECIDE) { 23547c6ae99SBarry Smith n = size / m; 23647c6ae99SBarry Smith } else { 23747c6ae99SBarry Smith /* try for squarish distribution */ 238369cc0aeSBarry Smith m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)N))); 23947c6ae99SBarry Smith if (!m) m = 1; 24047c6ae99SBarry Smith while (m > 0) { 24147c6ae99SBarry Smith n = size / m; 24247c6ae99SBarry Smith if (m * n == size) break; 24347c6ae99SBarry Smith m--; 24447c6ae99SBarry Smith } 2459371c9d4SSatish Balay if (M > N && m < n) { 2469371c9d4SSatish Balay PetscInt _m = m; 2479371c9d4SSatish Balay m = n; 2489371c9d4SSatish Balay n = _m; 2499371c9d4SSatish Balay } 25047c6ae99SBarry Smith } 2517a8be351SBarry Smith PetscCheck(m * n == size, comm, PETSC_ERR_PLIB, "Unable to create partition, check the size of the communicator and input m and n "); 2527a8be351SBarry Smith } else PetscCheck(m * n == size, comm, PETSC_ERR_ARG_OUTOFRANGE, "Given Bad partition"); 25347c6ae99SBarry Smith 25463a3b9bcSJacob Faibussowitsch PetscCheck(M >= m, comm, PETSC_ERR_ARG_OUTOFRANGE, "Partition in x direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT, M, m); 25563a3b9bcSJacob Faibussowitsch PetscCheck(N >= n, comm, PETSC_ERR_ARG_OUTOFRANGE, "Partition in y direction is too fine! %" PetscInt_FMT " %" PetscInt_FMT, N, n); 25647c6ae99SBarry Smith 25747c6ae99SBarry Smith /* 25847c6ae99SBarry Smith Determine locally owned region 25947c6ae99SBarry Smith xs is the first local node number, x is the number of local nodes 26047c6ae99SBarry Smith */ 26147c6ae99SBarry Smith if (!lx) { 2629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &dd->lx)); 26347c6ae99SBarry Smith lx = dd->lx; 2649371c9d4SSatish Balay for (i = 0; i < m; i++) { lx[i] = M / m + ((M % m) > i); } 26547c6ae99SBarry Smith } 26647c6ae99SBarry Smith x = lx[rank % m]; 26747c6ae99SBarry Smith xs = 0; 2689371c9d4SSatish Balay for (i = 0; i < (rank % m); i++) { xs += lx[i]; } 26976bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 27047c6ae99SBarry Smith left = xs; 2719371c9d4SSatish Balay for (i = (rank % m); i < m; i++) { left += lx[i]; } 27263a3b9bcSJacob Faibussowitsch PetscCheck(left == M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Sum of lx across processors not equal to M: %" PetscInt_FMT " %" PetscInt_FMT, left, M); 27376bd3646SJed Brown } 27447c6ae99SBarry Smith 27547c6ae99SBarry Smith /* 27647c6ae99SBarry Smith Determine locally owned region 27747c6ae99SBarry Smith ys is the first local node number, y is the number of local nodes 27847c6ae99SBarry Smith */ 27947c6ae99SBarry Smith if (!ly) { 2809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->ly)); 28147c6ae99SBarry Smith ly = dd->ly; 2829371c9d4SSatish Balay for (i = 0; i < n; i++) { ly[i] = N / n + ((N % n) > i); } 28347c6ae99SBarry Smith } 28447c6ae99SBarry Smith y = ly[rank / m]; 28547c6ae99SBarry Smith ys = 0; 2869371c9d4SSatish Balay for (i = 0; i < (rank / m); i++) { ys += ly[i]; } 28776bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 28847c6ae99SBarry Smith left = ys; 2899371c9d4SSatish Balay for (i = (rank / m); i < n; i++) { left += ly[i]; } 29063a3b9bcSJacob Faibussowitsch PetscCheck(left == N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Sum of ly across processors not equal to N: %" PetscInt_FMT " %" PetscInt_FMT, left, N); 29176bd3646SJed Brown } 29247c6ae99SBarry Smith 293bcea557cSEthan Coon /* 294bcea557cSEthan Coon check if the scatter requires more than one process neighbor or wraps around 295bcea557cSEthan Coon the domain more than once 296bcea557cSEthan Coon */ 29763a3b9bcSJacob Faibussowitsch 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); 29863a3b9bcSJacob Faibussowitsch 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); 29947c6ae99SBarry Smith xe = xs + x; 30047c6ae99SBarry Smith ye = ys + y; 30147c6ae99SBarry Smith 302ce00eea3SSatish Balay /* determine ghost region (Xs) and region scattered into (IXs) */ 303d9c9ebe5SPeter Brune if (xs - s > 0) { 3049371c9d4SSatish Balay Xs = xs - s; 3059371c9d4SSatish Balay IXs = xs - s; 30688661749SPeter Brune } else { 30788661749SPeter Brune if (bx) { 30888661749SPeter Brune Xs = xs - s; 30988661749SPeter Brune } else { 31088661749SPeter Brune Xs = 0; 31188661749SPeter Brune } 31288661749SPeter Brune IXs = 0; 31388661749SPeter Brune } 314d9c9ebe5SPeter Brune if (xe + s <= M) { 3159371c9d4SSatish Balay Xe = xe + s; 3169371c9d4SSatish Balay IXe = xe + s; 31788661749SPeter Brune } else { 31888661749SPeter Brune if (bx) { 3199371c9d4SSatish Balay Xs = xs - s; 3209371c9d4SSatish Balay Xe = xe + s; 32188661749SPeter Brune } else { 32288661749SPeter Brune Xe = M; 32388661749SPeter Brune } 32488661749SPeter Brune IXe = M; 32588661749SPeter Brune } 32647c6ae99SBarry Smith 327bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 328d9c9ebe5SPeter Brune IXs = xs - s; 329d9c9ebe5SPeter Brune IXe = xe + s; 330d9c9ebe5SPeter Brune Xs = xs - s; 331d9c9ebe5SPeter Brune Xe = xe + s; 33288661749SPeter Brune } 33347c6ae99SBarry Smith 334d9c9ebe5SPeter Brune if (ys - s > 0) { 3359371c9d4SSatish Balay Ys = ys - s; 3369371c9d4SSatish Balay IYs = ys - s; 33788661749SPeter Brune } else { 33888661749SPeter Brune if (by) { 33988661749SPeter Brune Ys = ys - s; 34088661749SPeter Brune } else { 34188661749SPeter Brune Ys = 0; 34288661749SPeter Brune } 34388661749SPeter Brune IYs = 0; 34488661749SPeter Brune } 345d9c9ebe5SPeter Brune if (ye + s <= N) { 3469371c9d4SSatish Balay Ye = ye + s; 3479371c9d4SSatish Balay IYe = ye + s; 34888661749SPeter Brune } else { 34988661749SPeter Brune if (by) { 35088661749SPeter Brune Ye = ye + s; 35188661749SPeter Brune } else { 35288661749SPeter Brune Ye = N; 35388661749SPeter Brune } 35488661749SPeter Brune IYe = N; 35588661749SPeter Brune } 35688661749SPeter Brune 357bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 358d9c9ebe5SPeter Brune IYs = ys - s; 359d9c9ebe5SPeter Brune IYe = ye + s; 360d9c9ebe5SPeter Brune Ys = ys - s; 361d9c9ebe5SPeter Brune Ye = ye + s; 36288661749SPeter Brune } 36388661749SPeter Brune 36488661749SPeter Brune /* stencil length in each direction */ 365d9c9ebe5SPeter Brune s_x = s; 366d9c9ebe5SPeter Brune s_y = s; 36747c6ae99SBarry Smith 36847c6ae99SBarry Smith /* determine starting point of each processor */ 36947c6ae99SBarry Smith nn = x * y; 3709566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size + 1, &bases, size, &ldims)); 3719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nn, 1, MPIU_INT, ldims, 1, MPIU_INT, comm)); 37247c6ae99SBarry Smith bases[0] = 0; 3739371c9d4SSatish Balay for (i = 1; i <= size; i++) { bases[i] = ldims[i - 1]; } 3749371c9d4SSatish Balay for (i = 1; i <= size; i++) { bases[i] += bases[i - 1]; } 375ce00eea3SSatish Balay base = bases[rank] * dof; 37647c6ae99SBarry Smith 37747c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 378ce00eea3SSatish Balay dd->Nlocal = x * y * dof; 3799566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, dof, dd->Nlocal, PETSC_DECIDE, NULL, &global)); 380ce00eea3SSatish Balay dd->nlocal = (Xe - Xs) * (Ye - Ys) * dof; 3819566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, dof, dd->nlocal, NULL, &local)); 38247c6ae99SBarry Smith 3834104a7a0SPatrick Sanan /* generate global to local vector scatter and local to global mapping*/ 38447c6ae99SBarry Smith 385ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 386ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 3879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((IXe - IXs) * (IYe - IYs), &idx)); 388d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 3899371c9d4SSatish Balay left = IXs - Xs; 3909371c9d4SSatish Balay right = left + (IXe - IXs); 3919371c9d4SSatish Balay down = IYs - Ys; 3929371c9d4SSatish Balay up = down + (IYe - IYs); 393ce00eea3SSatish Balay count = 0; 394ce00eea3SSatish Balay for (i = down; i < up; i++) { 3959371c9d4SSatish Balay for (j = left; j < right; j++) { idx[count++] = j + i * (Xe - Xs); } 396ce00eea3SSatish Balay } 3979566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm, dof, count, idx, PETSC_OWN_POINTER, &to)); 398ce00eea3SSatish Balay 39947c6ae99SBarry Smith } else { 40047c6ae99SBarry Smith /* must drop into cross shape region */ 40147c6ae99SBarry Smith /* ---------| 40247c6ae99SBarry Smith | top | 403ce00eea3SSatish Balay |--- ---| up 40447c6ae99SBarry Smith | middle | 40547c6ae99SBarry Smith | | 406ce00eea3SSatish Balay ---- ---- down 40747c6ae99SBarry Smith | bottom | 40847c6ae99SBarry Smith ----------- 40947c6ae99SBarry Smith Xs xs xe Xe */ 4109371c9d4SSatish Balay left = xs - Xs; 4119371c9d4SSatish Balay right = left + x; 4129371c9d4SSatish Balay down = ys - Ys; 4139371c9d4SSatish Balay up = down + y; 41447c6ae99SBarry Smith count = 0; 415ce00eea3SSatish Balay /* bottom */ 416ce00eea3SSatish Balay for (i = (IYs - Ys); i < down; i++) { 4179371c9d4SSatish Balay for (j = left; j < right; j++) { idx[count++] = j + i * (Xe - Xs); } 41847c6ae99SBarry Smith } 41947c6ae99SBarry Smith /* middle */ 42047c6ae99SBarry Smith for (i = down; i < up; i++) { 4219371c9d4SSatish Balay for (j = (IXs - Xs); j < (IXe - Xs); j++) { idx[count++] = j + i * (Xe - Xs); } 42247c6ae99SBarry Smith } 42347c6ae99SBarry Smith /* top */ 424ce00eea3SSatish Balay for (i = up; i < up + IYe - ye; i++) { 4259371c9d4SSatish Balay for (j = left; j < right; j++) { idx[count++] = j + i * (Xe - Xs); } 42647c6ae99SBarry Smith } 4279566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm, dof, count, idx, PETSC_OWN_POINTER, &to)); 42847c6ae99SBarry Smith } 42947c6ae99SBarry Smith 43047c6ae99SBarry Smith /* determine who lies on each side of us stored in n6 n7 n8 43147c6ae99SBarry Smith n3 n5 43247c6ae99SBarry Smith n0 n1 n2 43347c6ae99SBarry Smith */ 43447c6ae99SBarry Smith 43547c6ae99SBarry Smith /* Assume the Non-Periodic Case */ 43647c6ae99SBarry Smith n1 = rank - m; 43747c6ae99SBarry Smith if (rank % m) { 43847c6ae99SBarry Smith n0 = n1 - 1; 43947c6ae99SBarry Smith } else { 44047c6ae99SBarry Smith n0 = -1; 44147c6ae99SBarry Smith } 44247c6ae99SBarry Smith if ((rank + 1) % m) { 44347c6ae99SBarry Smith n2 = n1 + 1; 44447c6ae99SBarry Smith n5 = rank + 1; 4459371c9d4SSatish Balay n8 = rank + m + 1; 4469371c9d4SSatish Balay if (n8 >= m * n) n8 = -1; 44747c6ae99SBarry Smith } else { 4489371c9d4SSatish Balay n2 = -1; 4499371c9d4SSatish Balay n5 = -1; 4509371c9d4SSatish Balay n8 = -1; 45147c6ae99SBarry Smith } 45247c6ae99SBarry Smith if (rank % m) { 45347c6ae99SBarry Smith n3 = rank - 1; 4549371c9d4SSatish Balay n6 = n3 + m; 4559371c9d4SSatish Balay if (n6 >= m * n) n6 = -1; 45647c6ae99SBarry Smith } else { 4579371c9d4SSatish Balay n3 = -1; 4589371c9d4SSatish Balay n6 = -1; 45947c6ae99SBarry Smith } 4609371c9d4SSatish Balay n7 = rank + m; 4619371c9d4SSatish Balay if (n7 >= m * n) n7 = -1; 46247c6ae99SBarry Smith 463bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) { 46447c6ae99SBarry Smith /* Modify for Periodic Cases */ 46547c6ae99SBarry Smith /* Handle all four corners */ 46647c6ae99SBarry Smith if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m - 1; 46747c6ae99SBarry Smith if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0; 46847c6ae99SBarry Smith if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size - m; 46947c6ae99SBarry Smith if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size - 1; 47047c6ae99SBarry Smith 47147c6ae99SBarry Smith /* Handle Top and Bottom Sides */ 47247c6ae99SBarry Smith if (n1 < 0) n1 = rank + m * (n - 1); 47347c6ae99SBarry Smith if (n7 < 0) n7 = rank - m * (n - 1); 47447c6ae99SBarry Smith if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 47547c6ae99SBarry Smith if ((n3 >= 0) && (n6 < 0)) n6 = (rank % m) - 1; 47647c6ae99SBarry Smith if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 47747c6ae99SBarry Smith if ((n5 >= 0) && (n8 < 0)) n8 = (rank % m) + 1; 47847c6ae99SBarry Smith 47947c6ae99SBarry Smith /* Handle Left and Right Sides */ 48047c6ae99SBarry Smith if (n3 < 0) n3 = rank + (m - 1); 48147c6ae99SBarry Smith if (n5 < 0) n5 = rank - (m - 1); 48247c6ae99SBarry Smith if ((n1 >= 0) && (n0 < 0)) n0 = rank - 1; 48347c6ae99SBarry Smith if ((n1 >= 0) && (n2 < 0)) n2 = rank - 2 * m + 1; 48447c6ae99SBarry Smith if ((n7 >= 0) && (n6 < 0)) n6 = rank + 2 * m - 1; 48547c6ae99SBarry Smith if ((n7 >= 0) && (n8 < 0)) n8 = rank + 1; 486bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_PERIODIC) { /* Handle Top and Bottom Sides */ 487ce00eea3SSatish Balay if (n1 < 0) n1 = rank + m * (n - 1); 488ce00eea3SSatish Balay if (n7 < 0) n7 = rank - m * (n - 1); 489ce00eea3SSatish Balay if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 490ce00eea3SSatish Balay if ((n3 >= 0) && (n6 < 0)) n6 = (rank % m) - 1; 491ce00eea3SSatish Balay if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 492ce00eea3SSatish Balay if ((n5 >= 0) && (n8 < 0)) n8 = (rank % m) + 1; 493bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */ 494ce00eea3SSatish Balay if (n3 < 0) n3 = rank + (m - 1); 495ce00eea3SSatish Balay if (n5 < 0) n5 = rank - (m - 1); 496ce00eea3SSatish Balay if ((n1 >= 0) && (n0 < 0)) n0 = rank - 1; 497ce00eea3SSatish Balay if ((n1 >= 0) && (n2 < 0)) n2 = rank - 2 * m + 1; 498ce00eea3SSatish Balay if ((n7 >= 0) && (n6 < 0)) n6 = rank + 2 * m - 1; 499ce00eea3SSatish Balay if ((n7 >= 0) && (n8 < 0)) n8 = rank + 1; 50047c6ae99SBarry Smith } 501ce00eea3SSatish Balay 5029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(9, &dd->neighbors)); 5038865f1eaSKarl Rupp 50447c6ae99SBarry Smith dd->neighbors[0] = n0; 50547c6ae99SBarry Smith dd->neighbors[1] = n1; 50647c6ae99SBarry Smith dd->neighbors[2] = n2; 50747c6ae99SBarry Smith dd->neighbors[3] = n3; 50847c6ae99SBarry Smith dd->neighbors[4] = rank; 50947c6ae99SBarry Smith dd->neighbors[5] = n5; 51047c6ae99SBarry Smith dd->neighbors[6] = n6; 51147c6ae99SBarry Smith dd->neighbors[7] = n7; 51247c6ae99SBarry Smith dd->neighbors[8] = n8; 51347c6ae99SBarry Smith 514d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 51547c6ae99SBarry Smith /* save corner processor numbers */ 5169371c9d4SSatish Balay sn0 = n0; 5179371c9d4SSatish Balay sn2 = n2; 5189371c9d4SSatish Balay sn6 = n6; 5199371c9d4SSatish Balay sn8 = n8; 52047c6ae99SBarry Smith n0 = n2 = n6 = n8 = -1; 52147c6ae99SBarry Smith } 52247c6ae99SBarry Smith 5239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((Xe - Xs) * (Ye - Ys), &idx)); 52447c6ae99SBarry Smith 525ce00eea3SSatish Balay nn = 0; 52647c6ae99SBarry Smith xbase = bases[rank]; 52747c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 52847c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 529ce00eea3SSatish Balay x_t = lx[n0 % m]; 53047c6ae99SBarry Smith y_t = ly[(n0 / m)]; 53147c6ae99SBarry Smith s_t = bases[n0] + x_t * y_t - (s_y - i) * x_t - s_x; 5328865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 53347c6ae99SBarry Smith } 534ac119b13SBarry Smith 53547c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 53647c6ae99SBarry Smith x_t = x; 53747c6ae99SBarry Smith y_t = ly[(n1 / m)]; 53847c6ae99SBarry Smith s_t = bases[n1] + x_t * y_t - (s_y + 1 - i) * x_t; 5398865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 540bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 5418865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = bases[rank] + x * (s_y - i + 1) + j; 54247c6ae99SBarry Smith } 543ac119b13SBarry Smith 54447c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 545ce00eea3SSatish Balay x_t = lx[n2 % m]; 54647c6ae99SBarry Smith y_t = ly[(n2 / m)]; 54747c6ae99SBarry Smith s_t = bases[n2] + x_t * y_t - (s_y + 1 - i) * x_t; 5488865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 54947c6ae99SBarry Smith } 55047c6ae99SBarry Smith } 55147c6ae99SBarry Smith 55247c6ae99SBarry Smith for (i = 0; i < y; i++) { 55347c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 554ce00eea3SSatish Balay x_t = lx[n3 % m]; 55547c6ae99SBarry Smith /* y_t = y; */ 55647c6ae99SBarry Smith s_t = bases[n3] + (i + 1) * x_t - s_x; 5578865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 558bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 5598865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x * i + s_x - j; 56047c6ae99SBarry Smith } 56147c6ae99SBarry Smith 5628865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = xbase++; /* interior */ 56347c6ae99SBarry Smith 56447c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 565ce00eea3SSatish Balay x_t = lx[n5 % m]; 56647c6ae99SBarry Smith /* y_t = y; */ 56747c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 5688865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 569bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 5708865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x * (i + 1) - 2 - j; 57147c6ae99SBarry Smith } 57247c6ae99SBarry Smith } 57347c6ae99SBarry Smith 57447c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 57547c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 576ce00eea3SSatish Balay x_t = lx[n6 % m]; 57747c6ae99SBarry Smith /* y_t = ly[(n6/m)]; */ 57847c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 5798865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 58047c6ae99SBarry Smith } 581ac119b13SBarry Smith 58247c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 58347c6ae99SBarry Smith x_t = x; 58447c6ae99SBarry Smith /* y_t = ly[(n7/m)]; */ 58547c6ae99SBarry Smith s_t = bases[n7] + (i - 1) * x_t; 5868865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 587bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 5888865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = bases[rank] + x * (y - i - 1) + j; 58947c6ae99SBarry Smith } 590ac119b13SBarry Smith 59147c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 592ce00eea3SSatish Balay x_t = lx[n8 % m]; 59347c6ae99SBarry Smith /* y_t = ly[(n8/m)]; */ 59447c6ae99SBarry Smith s_t = bases[n8] + (i - 1) * x_t; 5958865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 59647c6ae99SBarry Smith } 59747c6ae99SBarry Smith } 59847c6ae99SBarry Smith 5999566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm, dof, nn, idx, PETSC_USE_POINTER, &from)); 6009566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(global, from, local, to, >ol)); 6019566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da, (PetscObject)gtol)); 6029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 6039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 60447c6ae99SBarry Smith 605d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 6069371c9d4SSatish Balay n0 = sn0; 6079371c9d4SSatish Balay n2 = sn2; 6089371c9d4SSatish Balay n6 = sn6; 6099371c9d4SSatish Balay n8 = sn8; 610ce00eea3SSatish Balay } 611ce00eea3SSatish Balay 612288e7d53SBarry Smith if (((stencil_type == DMDA_STENCIL_STAR) || (bx && bx != DM_BOUNDARY_PERIODIC) || (by && by != DM_BOUNDARY_PERIODIC))) { 61347c6ae99SBarry Smith /* 61447c6ae99SBarry Smith Recompute the local to global mappings, this time keeping the 615ce00eea3SSatish Balay information about the cross corner processor numbers and any ghosted 616ce00eea3SSatish Balay but not periodic indices. 61747c6ae99SBarry Smith */ 61847c6ae99SBarry Smith nn = 0; 61947c6ae99SBarry Smith xbase = bases[rank]; 62047c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 62147c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 622ce00eea3SSatish Balay x_t = lx[n0 % m]; 62347c6ae99SBarry Smith y_t = ly[(n0 / m)]; 62447c6ae99SBarry Smith s_t = bases[n0] + x_t * y_t - (s_y - i) * x_t - s_x; 6258865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 626ce00eea3SSatish Balay } else if (xs - Xs > 0 && ys - Ys > 0) { 6278865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 62847c6ae99SBarry Smith } 62947c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 63047c6ae99SBarry Smith x_t = x; 63147c6ae99SBarry Smith y_t = ly[(n1 / m)]; 63247c6ae99SBarry Smith s_t = bases[n1] + x_t * y_t - (s_y + 1 - i) * x_t; 6338865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 634ce00eea3SSatish Balay } else if (ys - Ys > 0) { 635bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 6368865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = bases[rank] + x * (s_y - i + 1) + j; 637624904c4SBarry Smith } else { 6388865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = -1; 63947c6ae99SBarry Smith } 640624904c4SBarry Smith } 64147c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 642ce00eea3SSatish Balay x_t = lx[n2 % m]; 64347c6ae99SBarry Smith y_t = ly[(n2 / m)]; 64447c6ae99SBarry Smith s_t = bases[n2] + x_t * y_t - (s_y + 1 - i) * x_t; 6458865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 646ce00eea3SSatish Balay } else if (Xe - xe > 0 && ys - Ys > 0) { 6478865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 64847c6ae99SBarry Smith } 64947c6ae99SBarry Smith } 65047c6ae99SBarry Smith 65147c6ae99SBarry Smith for (i = 0; i < y; i++) { 65247c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 653ce00eea3SSatish Balay x_t = lx[n3 % m]; 65447c6ae99SBarry Smith /* y_t = y; */ 65547c6ae99SBarry Smith s_t = bases[n3] + (i + 1) * x_t - s_x; 6568865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 657ce00eea3SSatish Balay } else if (xs - Xs > 0) { 658bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 6598865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x * i + s_x - j; 660624904c4SBarry Smith } else { 6618865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 66247c6ae99SBarry Smith } 663624904c4SBarry Smith } 66447c6ae99SBarry Smith 6658865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = xbase++; /* interior */ 66647c6ae99SBarry Smith 66747c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 668ce00eea3SSatish Balay x_t = lx[n5 % m]; 66947c6ae99SBarry Smith /* y_t = y; */ 67047c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 6718865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 672ce00eea3SSatish Balay } else if (Xe - xe > 0) { 673bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 6748865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x * (i + 1) - 2 - j; 675624904c4SBarry Smith } else { 6768865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 67747c6ae99SBarry Smith } 67847c6ae99SBarry Smith } 679624904c4SBarry Smith } 68047c6ae99SBarry Smith 68147c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 68247c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 683ce00eea3SSatish Balay x_t = lx[n6 % m]; 68447c6ae99SBarry Smith /* y_t = ly[(n6/m)]; */ 68547c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 6868865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 687ce00eea3SSatish Balay } else if (xs - Xs > 0 && Ye - ye > 0) { 6888865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 68947c6ae99SBarry Smith } 69047c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 69147c6ae99SBarry Smith x_t = x; 69247c6ae99SBarry Smith /* y_t = ly[(n7/m)]; */ 69347c6ae99SBarry Smith s_t = bases[n7] + (i - 1) * x_t; 6948865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 695ce00eea3SSatish Balay } else if (Ye - ye > 0) { 696bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 6978865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = bases[rank] + x * (y - i - 1) + j; 698624904c4SBarry Smith } else { 6998865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = -1; 70047c6ae99SBarry Smith } 701624904c4SBarry Smith } 70247c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 703ce00eea3SSatish Balay x_t = lx[n8 % m]; 70447c6ae99SBarry Smith /* y_t = ly[(n8/m)]; */ 70547c6ae99SBarry Smith s_t = bases[n8] + (i - 1) * x_t; 7068865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 707ce00eea3SSatish Balay } else if (Xe - xe > 0 && Ye - ye > 0) { 7088865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 70947c6ae99SBarry Smith } 71047c6ae99SBarry Smith } 71147c6ae99SBarry Smith } 712ce00eea3SSatish Balay /* 713ce00eea3SSatish Balay Set the local to global ordering in the global vector, this allows use 714ce00eea3SSatish Balay of VecSetValuesLocal(). 715ce00eea3SSatish Balay */ 7169566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, dof, nn, idx, PETSC_OWN_POINTER, &da->ltogmap)); 7179566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da, (PetscObject)da->ltogmap)); 71847c6ae99SBarry Smith 7199566063dSJacob Faibussowitsch PetscCall(PetscFree2(bases, ldims)); 7209371c9d4SSatish Balay dd->m = m; 7219371c9d4SSatish Balay dd->n = n; 722ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 7239371c9d4SSatish Balay dd->xs = xs * dof; 7249371c9d4SSatish Balay dd->xe = xe * dof; 7259371c9d4SSatish Balay dd->ys = ys; 7269371c9d4SSatish Balay dd->ye = ye; 7279371c9d4SSatish Balay dd->zs = 0; 7289371c9d4SSatish Balay dd->ze = 1; 7299371c9d4SSatish Balay dd->Xs = Xs * dof; 7309371c9d4SSatish Balay dd->Xe = Xe * dof; 7319371c9d4SSatish Balay dd->Ys = Ys; 7329371c9d4SSatish Balay dd->Ye = Ye; 7339371c9d4SSatish Balay dd->Zs = 0; 7349371c9d4SSatish Balay dd->Ze = 1; 73547c6ae99SBarry Smith 7369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 7379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 73847c6ae99SBarry Smith 73947c6ae99SBarry Smith dd->gtol = gtol; 74047c6ae99SBarry Smith dd->base = base; 7419a42bb27SBarry Smith da->ops->view = DMView_DA_2d; 7420298fd71SBarry Smith dd->ltol = NULL; 7430298fd71SBarry Smith dd->ao = NULL; 74447c6ae99SBarry Smith PetscFunctionReturn(0); 74547c6ae99SBarry Smith } 74647c6ae99SBarry Smith 74747c6ae99SBarry Smith /*@C 748aa219208SBarry Smith DMDACreate2d - Creates an object that will manage the communication of two-dimensional 74947c6ae99SBarry Smith regular array data that is distributed across some processors. 75047c6ae99SBarry Smith 751d083f849SBarry Smith Collective 75247c6ae99SBarry Smith 75347c6ae99SBarry Smith Input Parameters: 75447c6ae99SBarry Smith + comm - MPI communicator 7551321219cSEthan Coon . bx,by - type of ghost nodes the array have. 756bff4a2f0SMatthew G. Knepley Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC. 757aa219208SBarry Smith . stencil_type - stencil type. Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 758897f7067SBarry Smith . M,N - global dimension in each direction of the array 75947c6ae99SBarry Smith . m,n - corresponding number of processors in each dimension 76047c6ae99SBarry Smith (or PETSC_DECIDE to have calculated) 76147c6ae99SBarry Smith . dof - number of degrees of freedom per node 76247c6ae99SBarry Smith . s - stencil width 76347c6ae99SBarry Smith - lx, ly - arrays containing the number of nodes in each cell along 7640298fd71SBarry Smith the x and y coordinates, or NULL. If non-null, these 76547c6ae99SBarry Smith must be of length as m and n, and the corresponding 76647c6ae99SBarry Smith m and n cannot be PETSC_DECIDE. The sum of the lx[] entries 76747c6ae99SBarry Smith must be M, and the sum of the ly[] entries must be N. 76847c6ae99SBarry Smith 76947c6ae99SBarry Smith Output Parameter: 77047c6ae99SBarry Smith . da - the resulting distributed array object 77147c6ae99SBarry Smith 77247c6ae99SBarry Smith Options Database Key: 773706a11cbSBarry Smith + -dm_view - Calls DMView() at the conclusion of DMDACreate2d() 774e43dc8daSBarry Smith . -da_grid_x <nx> - number of grid points in x direction 775e43dc8daSBarry Smith . -da_grid_y <ny> - number of grid points in y direction 77647c6ae99SBarry Smith . -da_processors_x <nx> - number of processors in x direction 77747c6ae99SBarry Smith . -da_processors_y <ny> - number of processors in y direction 778e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 779e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 780e43dc8daSBarry Smith - -da_refine <n> - refine the DMDA n times before creating 781e0f5d30fSBarry Smith 78247c6ae99SBarry Smith Level: beginner 78347c6ae99SBarry Smith 78447c6ae99SBarry Smith Notes: 785aa219208SBarry Smith The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the 786aa219208SBarry Smith standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes 78747c6ae99SBarry Smith the standard 9-pt stencil. 78847c6ae99SBarry Smith 789aa219208SBarry Smith The array data itself is NOT stored in the DMDA, it is stored in Vec objects; 790564755cdSBarry Smith The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() 791564755cdSBarry Smith and DMCreateLocalVector() and calls to VecDuplicate() if more are needed. 79247c6ae99SBarry Smith 793897f7067SBarry Smith You must call DMSetUp() after this call before using this DM. 794897f7067SBarry Smith 795897f7067SBarry Smith If you wish to use the options database to change values in the DMDA call DMSetFromOptions() after this call 796897f7067SBarry Smith but before DMSetUp(). 797897f7067SBarry Smith 798db781477SPatrick Sanan .seealso: `DMDestroy()`, `DMView()`, `DMDACreate1d()`, `DMDACreate3d()`, `DMGlobalToLocalBegin()`, `DMDAGetRefinementFactor()`, 799db781477SPatrick Sanan `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDASetRefinementFactor()`, 800db781477SPatrick Sanan `DMDAGetInfo()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMDACreateNaturalVector()`, `DMLoad()`, `DMDAGetOwnershipRanges()`, 801db781477SPatrick Sanan `DMStagCreate2d()` 80247c6ae99SBarry Smith 80347c6ae99SBarry Smith @*/ 804fe16a2e9SBarry Smith 8059371c9d4SSatish Balay PetscErrorCode DMDACreate2d(MPI_Comm comm, DMBoundaryType bx, DMBoundaryType by, DMDAStencilType stencil_type, PetscInt M, PetscInt N, PetscInt m, PetscInt n, PetscInt dof, PetscInt s, const PetscInt lx[], const PetscInt ly[], DM *da) { 80647c6ae99SBarry Smith PetscFunctionBegin; 8079566063dSJacob Faibussowitsch PetscCall(DMDACreate(comm, da)); 8089566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*da, 2)); 8099566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*da, M, N, 1)); 8109566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*da, m, n, PETSC_DECIDE)); 8119566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE)); 8129566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*da, dof)); 8139566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*da, stencil_type)); 8149566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*da, s)); 8159566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*da, lx, ly, NULL)); 81647c6ae99SBarry Smith PetscFunctionReturn(0); 81747c6ae99SBarry Smith } 818