1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 29804daf3SBarry Smith #include <petscdraw.h> 347c6ae99SBarry Smith 4d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMView_DA_2d(DM da, PetscViewer viewer) 5d71ae5a4SJacob Faibussowitsch { 647c6ae99SBarry Smith PetscMPIInt rank; 7c9493c35SStefano Zampini PetscBool iascii, isdraw, isglvis, isbinary; 847c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 9d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 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)); 20d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 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) { 28*1690c2aeSBarry Smith PetscInt i, nmax = 0, nmin = PETSC_INT_MAX, 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)); 443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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)); 683ba16761SJacob Faibussowitsch if (isnull) PetscFunctionReturn(PETSC_SUCCESS); 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; 7948a46eb9SPierre 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; 8248a46eb9SPierre 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)); 136d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 1379a42bb27SBarry Smith } else if (ismatlab) { 1389566063dSJacob Faibussowitsch PetscCall(DMView_DA_Matlab(da, viewer)); 1399a42bb27SBarry Smith #endif 14011aeaf0aSBarry Smith } 1413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14247c6ae99SBarry Smith } 14347c6ae99SBarry Smith 14447c6ae99SBarry Smith #if defined(new) 14547c6ae99SBarry Smith /* 14601c1178eSBarry 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 */ 154d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetDiagonal_MFFD(DM da, Vec U, Vec a) 155d71ae5a4SJacob Faibussowitsch { 15647c6ae99SBarry Smith PetscScalar h, *aa, *ww, v; 15747c6ae99SBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON, umin = 100.0 * PETSC_SQRT_MACHINE_EPSILON; 15847c6ae99SBarry Smith PetscInt gI, nI; 15947c6ae99SBarry Smith MatStencil stencil; 160aa219208SBarry Smith DMDALocalInfo info; 16147c6ae99SBarry Smith 16247c6ae99SBarry Smith PetscFunctionBegin; 1639566063dSJacob Faibussowitsch PetscCall((*ctx->func)(0, U, a, ctx->funcctx)); 1649566063dSJacob Faibussowitsch PetscCall((*ctx->funcisetbase)(U, ctx->funcctx)); 16547c6ae99SBarry Smith 1669566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &ww)); 1679566063dSJacob Faibussowitsch PetscCall(VecGetArray(a, &aa)); 16847c6ae99SBarry Smith 16947c6ae99SBarry Smith nI = 0; 17047c6ae99SBarry Smith h = ww[gI]; 17147c6ae99SBarry Smith if (h == 0.0) h = 1.0; 17247c6ae99SBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 17347c6ae99SBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 17447c6ae99SBarry Smith h *= epsilon; 17547c6ae99SBarry Smith 17647c6ae99SBarry Smith ww[gI] += h; 1779566063dSJacob Faibussowitsch PetscCall((*ctx->funci)(i, w, &v, ctx->funcctx)); 17847c6ae99SBarry Smith aa[nI] = (v - aa[nI]) / h; 17947c6ae99SBarry Smith ww[gI] -= h; 18047c6ae99SBarry Smith nI++; 1818865f1eaSKarl Rupp 1829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &ww)); 1839566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a, &aa)); 1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18547c6ae99SBarry Smith } 18647c6ae99SBarry Smith #endif 18747c6ae99SBarry Smith 188d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetUp_DA_2D(DM da) 189d71ae5a4SJacob Faibussowitsch { 19047c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 19147c6ae99SBarry Smith const PetscInt M = dd->M; 19247c6ae99SBarry Smith const PetscInt N = dd->N; 1936497c311SBarry Smith PetscMPIInt m, n; 19447c6ae99SBarry Smith const PetscInt dof = dd->w; 19547c6ae99SBarry Smith const PetscInt s = dd->s; 196bff4a2f0SMatthew G. Knepley DMBoundaryType bx = dd->bx; 197bff4a2f0SMatthew G. Knepley DMBoundaryType by = dd->by; 19819fd82e9SBarry Smith DMDAStencilType stencil_type = dd->stencil_type; 19947c6ae99SBarry Smith PetscInt *lx = dd->lx; 20047c6ae99SBarry Smith PetscInt *ly = dd->ly; 20147c6ae99SBarry Smith MPI_Comm comm; 2026497c311SBarry Smith PetscMPIInt rank, size, n0, n1, n2, n3, n5, n6, n7, n8; 203bd1fc5aeSBarry Smith PetscInt xs, xe, ys, ye, x, y, Xs, Xe, Ys, Ye, IXs, IXe, IYs, IYe; 2046497c311SBarry Smith PetscInt up, down, left, right, i, *idx, nn; 205db87c5ecSEthan Coon PetscInt xbase, *bases, *ldims, j, x_t, y_t, s_t, base, count; 20647c6ae99SBarry Smith PetscInt s_x, s_y; /* s proportionalized to w */ 2076497c311SBarry Smith PetscMPIInt sn0 = 0, sn2 = 0, sn6 = 0, sn8 = 0; 20847c6ae99SBarry Smith Vec local, global; 209bd1fc5aeSBarry Smith VecScatter gtol; 21045b6f7e9SBarry Smith IS to, from; 21147c6ae99SBarry Smith 21247c6ae99SBarry Smith PetscFunctionBegin; 2137a8be351SBarry 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"); 2149566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 2153855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES) 2167de69702SBarry Smith 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); 2173855c12bSBarry Smith #endif 2186497c311SBarry Smith PetscCall(PetscMPIIntCast(dd->m, &m)); 2196497c311SBarry Smith PetscCall(PetscMPIIntCast(dd->n, &n)); 22047c6ae99SBarry Smith 2219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 22347c6ae99SBarry Smith 2247d310018SBarry Smith dd->p = 1; 22547c6ae99SBarry Smith if (m != PETSC_DECIDE) { 2266497c311SBarry Smith PetscCheck(m >= 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in X direction: %d", m); 2276497c311SBarry Smith PetscCheck(m <= size, comm, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in X direction: %d %d", m, size); 22847c6ae99SBarry Smith } 22947c6ae99SBarry Smith if (n != PETSC_DECIDE) { 2306497c311SBarry Smith PetscCheck(n >= 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors in Y direction: %d", n); 2316497c311SBarry Smith PetscCheck(n <= size, comm, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in Y direction: %d %d", n, size); 23247c6ae99SBarry Smith } 23347c6ae99SBarry Smith 23447c6ae99SBarry Smith if (m == PETSC_DECIDE || n == PETSC_DECIDE) { 23547c6ae99SBarry Smith if (n != PETSC_DECIDE) { 23647c6ae99SBarry Smith m = size / n; 23747c6ae99SBarry Smith } else if (m != PETSC_DECIDE) { 23847c6ae99SBarry Smith n = size / m; 23947c6ae99SBarry Smith } else { 24047c6ae99SBarry Smith /* try for squarish distribution */ 2416497c311SBarry Smith m = (PetscMPIInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)N))); 24247c6ae99SBarry Smith if (!m) m = 1; 24347c6ae99SBarry Smith while (m > 0) { 24447c6ae99SBarry Smith n = size / m; 24547c6ae99SBarry Smith if (m * n == size) break; 24647c6ae99SBarry Smith m--; 24747c6ae99SBarry Smith } 2489371c9d4SSatish Balay if (M > N && m < n) { 2496497c311SBarry Smith PetscMPIInt _m = m; 2509371c9d4SSatish Balay m = n; 2519371c9d4SSatish Balay n = _m; 2529371c9d4SSatish Balay } 25347c6ae99SBarry Smith } 2547a8be351SBarry Smith PetscCheck(m * n == size, comm, PETSC_ERR_PLIB, "Unable to create partition, check the size of the communicator and input m and n "); 2557a8be351SBarry Smith } else PetscCheck(m * n == size, comm, PETSC_ERR_ARG_OUTOFRANGE, "Given Bad partition"); 25647c6ae99SBarry Smith 2576497c311SBarry Smith PetscCheck(M >= m, comm, PETSC_ERR_ARG_OUTOFRANGE, "Partition in x direction is too fine! %" PetscInt_FMT " %d", M, m); 2586497c311SBarry Smith PetscCheck(N >= n, comm, PETSC_ERR_ARG_OUTOFRANGE, "Partition in y direction is too fine! %" PetscInt_FMT " %d", N, n); 25947c6ae99SBarry Smith 26047c6ae99SBarry Smith /* 26147c6ae99SBarry Smith Determine locally owned region 26247c6ae99SBarry Smith xs is the first local node number, x is the number of local nodes 26347c6ae99SBarry Smith */ 26447c6ae99SBarry Smith if (!lx) { 2659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &dd->lx)); 26647c6ae99SBarry Smith lx = dd->lx; 267ad540459SPierre Jolivet for (i = 0; i < m; i++) lx[i] = M / m + ((M % m) > i); 26847c6ae99SBarry Smith } 26947c6ae99SBarry Smith x = lx[rank % m]; 27047c6ae99SBarry Smith xs = 0; 271ad540459SPierre Jolivet for (i = 0; i < (rank % m); i++) xs += lx[i]; 27276bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 27347c6ae99SBarry Smith left = xs; 274ad540459SPierre Jolivet for (i = (rank % m); i < m; i++) left += lx[i]; 27563a3b9bcSJacob 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); 27676bd3646SJed Brown } 27747c6ae99SBarry Smith 27847c6ae99SBarry Smith /* 27947c6ae99SBarry Smith Determine locally owned region 28047c6ae99SBarry Smith ys is the first local node number, y is the number of local nodes 28147c6ae99SBarry Smith */ 28247c6ae99SBarry Smith if (!ly) { 2839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dd->ly)); 28447c6ae99SBarry Smith ly = dd->ly; 285ad540459SPierre Jolivet for (i = 0; i < n; i++) ly[i] = N / n + ((N % n) > i); 28647c6ae99SBarry Smith } 28747c6ae99SBarry Smith y = ly[rank / m]; 28847c6ae99SBarry Smith ys = 0; 289ad540459SPierre Jolivet for (i = 0; i < (rank / m); i++) ys += ly[i]; 29076bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 29147c6ae99SBarry Smith left = ys; 292ad540459SPierre Jolivet for (i = (rank / m); i < n; i++) left += ly[i]; 29363a3b9bcSJacob 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); 29476bd3646SJed Brown } 29547c6ae99SBarry Smith 296bcea557cSEthan Coon /* 297bcea557cSEthan Coon check if the scatter requires more than one process neighbor or wraps around 298bcea557cSEthan Coon the domain more than once 299bcea557cSEthan Coon */ 30063a3b9bcSJacob 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); 30163a3b9bcSJacob 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); 3024ad8454bSPierre Jolivet PetscCheck((x > s) || (bx != DM_BOUNDARY_MIRROR), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local x-width of domain x %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT " with mirror", x, s); 3034ad8454bSPierre Jolivet PetscCheck((y > s) || (by != DM_BOUNDARY_MIRROR), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local y-width of domain y %" PetscInt_FMT " is smaller than stencil width s %" PetscInt_FMT " with mirror", y, s); 30447c6ae99SBarry Smith xe = xs + x; 30547c6ae99SBarry Smith ye = ys + y; 30647c6ae99SBarry Smith 307ce00eea3SSatish Balay /* determine ghost region (Xs) and region scattered into (IXs) */ 308d9c9ebe5SPeter Brune if (xs - s > 0) { 3099371c9d4SSatish Balay Xs = xs - s; 3109371c9d4SSatish Balay IXs = xs - s; 31188661749SPeter Brune } else { 31288661749SPeter Brune if (bx) { 31388661749SPeter Brune Xs = xs - s; 31488661749SPeter Brune } else { 31588661749SPeter Brune Xs = 0; 31688661749SPeter Brune } 31788661749SPeter Brune IXs = 0; 31888661749SPeter Brune } 319d9c9ebe5SPeter Brune if (xe + s <= M) { 3209371c9d4SSatish Balay Xe = xe + s; 3219371c9d4SSatish Balay IXe = xe + s; 32288661749SPeter Brune } else { 32388661749SPeter Brune if (bx) { 3249371c9d4SSatish Balay Xs = xs - s; 3259371c9d4SSatish Balay Xe = xe + s; 32688661749SPeter Brune } else { 32788661749SPeter Brune Xe = M; 32888661749SPeter Brune } 32988661749SPeter Brune IXe = M; 33088661749SPeter Brune } 33147c6ae99SBarry Smith 332bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || bx == DM_BOUNDARY_MIRROR) { 333d9c9ebe5SPeter Brune IXs = xs - s; 334d9c9ebe5SPeter Brune IXe = xe + s; 335d9c9ebe5SPeter Brune Xs = xs - s; 336d9c9ebe5SPeter Brune Xe = xe + s; 33788661749SPeter Brune } 33847c6ae99SBarry Smith 339d9c9ebe5SPeter Brune if (ys - s > 0) { 3409371c9d4SSatish Balay Ys = ys - s; 3419371c9d4SSatish Balay IYs = ys - s; 34288661749SPeter Brune } else { 34388661749SPeter Brune if (by) { 34488661749SPeter Brune Ys = ys - s; 34588661749SPeter Brune } else { 34688661749SPeter Brune Ys = 0; 34788661749SPeter Brune } 34888661749SPeter Brune IYs = 0; 34988661749SPeter Brune } 350d9c9ebe5SPeter Brune if (ye + s <= N) { 3519371c9d4SSatish Balay Ye = ye + s; 3529371c9d4SSatish Balay IYe = ye + s; 35388661749SPeter Brune } else { 35488661749SPeter Brune if (by) { 35588661749SPeter Brune Ye = ye + s; 35688661749SPeter Brune } else { 35788661749SPeter Brune Ye = N; 35888661749SPeter Brune } 35988661749SPeter Brune IYe = N; 36088661749SPeter Brune } 36188661749SPeter Brune 362bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_MIRROR) { 363d9c9ebe5SPeter Brune IYs = ys - s; 364d9c9ebe5SPeter Brune IYe = ye + s; 365d9c9ebe5SPeter Brune Ys = ys - s; 366d9c9ebe5SPeter Brune Ye = ye + s; 36788661749SPeter Brune } 36888661749SPeter Brune 36988661749SPeter Brune /* stencil length in each direction */ 370d9c9ebe5SPeter Brune s_x = s; 371d9c9ebe5SPeter Brune s_y = s; 37247c6ae99SBarry Smith 37347c6ae99SBarry Smith /* determine starting point of each processor */ 37447c6ae99SBarry Smith nn = x * y; 3759566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size + 1, &bases, size, &ldims)); 3769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nn, 1, MPIU_INT, ldims, 1, MPIU_INT, comm)); 37747c6ae99SBarry Smith bases[0] = 0; 378ad540459SPierre Jolivet for (i = 1; i <= size; i++) bases[i] = ldims[i - 1]; 379ad540459SPierre Jolivet for (i = 1; i <= size; i++) bases[i] += bases[i - 1]; 380ce00eea3SSatish Balay base = bases[rank] * dof; 38147c6ae99SBarry Smith 38247c6ae99SBarry Smith /* allocate the base parallel and sequential vectors */ 383ce00eea3SSatish Balay dd->Nlocal = x * y * dof; 3849566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, dof, dd->Nlocal, PETSC_DECIDE, NULL, &global)); 385ce00eea3SSatish Balay dd->nlocal = (Xe - Xs) * (Ye - Ys) * dof; 3869566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, dof, dd->nlocal, NULL, &local)); 38747c6ae99SBarry Smith 3884104a7a0SPatrick Sanan /* generate global to local vector scatter and local to global mapping*/ 38947c6ae99SBarry Smith 390ce00eea3SSatish Balay /* global to local must include ghost points within the domain, 391ce00eea3SSatish Balay but not ghost points outside the domain that aren't periodic */ 3929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((IXe - IXs) * (IYe - IYs), &idx)); 393d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_BOX) { 3949371c9d4SSatish Balay left = IXs - Xs; 3959371c9d4SSatish Balay right = left + (IXe - IXs); 3969371c9d4SSatish Balay down = IYs - Ys; 3979371c9d4SSatish Balay up = down + (IYe - IYs); 398ce00eea3SSatish Balay count = 0; 399ce00eea3SSatish Balay for (i = down; i < up; i++) { 400ad540459SPierre Jolivet for (j = left; j < right; j++) idx[count++] = j + i * (Xe - Xs); 401ce00eea3SSatish Balay } 4029566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm, dof, count, idx, PETSC_OWN_POINTER, &to)); 403ce00eea3SSatish Balay 40447c6ae99SBarry Smith } else { 40547c6ae99SBarry Smith /* must drop into cross shape region */ 40647c6ae99SBarry Smith /* ---------| 40747c6ae99SBarry Smith | top | 408ce00eea3SSatish Balay |--- ---| up 40947c6ae99SBarry Smith | middle | 41047c6ae99SBarry Smith | | 411ce00eea3SSatish Balay ---- ---- down 41247c6ae99SBarry Smith | bottom | 41347c6ae99SBarry Smith ----------- 41447c6ae99SBarry Smith Xs xs xe Xe */ 4159371c9d4SSatish Balay left = xs - Xs; 4169371c9d4SSatish Balay right = left + x; 4179371c9d4SSatish Balay down = ys - Ys; 4189371c9d4SSatish Balay up = down + y; 41947c6ae99SBarry Smith count = 0; 420ce00eea3SSatish Balay /* bottom */ 421ce00eea3SSatish Balay for (i = (IYs - Ys); i < down; i++) { 422ad540459SPierre Jolivet for (j = left; j < right; j++) idx[count++] = j + i * (Xe - Xs); 42347c6ae99SBarry Smith } 42447c6ae99SBarry Smith /* middle */ 42547c6ae99SBarry Smith for (i = down; i < up; i++) { 426ad540459SPierre Jolivet for (j = (IXs - Xs); j < (IXe - Xs); j++) idx[count++] = j + i * (Xe - Xs); 42747c6ae99SBarry Smith } 42847c6ae99SBarry Smith /* top */ 429ce00eea3SSatish Balay for (i = up; i < up + IYe - ye; i++) { 430ad540459SPierre Jolivet for (j = left; j < right; j++) idx[count++] = j + i * (Xe - Xs); 43147c6ae99SBarry Smith } 4329566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm, dof, count, idx, PETSC_OWN_POINTER, &to)); 43347c6ae99SBarry Smith } 43447c6ae99SBarry Smith 43547c6ae99SBarry Smith /* determine who lies on each side of us stored in n6 n7 n8 43647c6ae99SBarry Smith n3 n5 43747c6ae99SBarry Smith n0 n1 n2 43847c6ae99SBarry Smith */ 43947c6ae99SBarry Smith 44047c6ae99SBarry Smith /* Assume the Non-Periodic Case */ 44147c6ae99SBarry Smith n1 = rank - m; 44247c6ae99SBarry Smith if (rank % m) { 44347c6ae99SBarry Smith n0 = n1 - 1; 44447c6ae99SBarry Smith } else { 44547c6ae99SBarry Smith n0 = -1; 44647c6ae99SBarry Smith } 44747c6ae99SBarry Smith if ((rank + 1) % m) { 44847c6ae99SBarry Smith n2 = n1 + 1; 44947c6ae99SBarry Smith n5 = rank + 1; 4509371c9d4SSatish Balay n8 = rank + m + 1; 4519371c9d4SSatish Balay if (n8 >= m * n) n8 = -1; 45247c6ae99SBarry Smith } else { 4539371c9d4SSatish Balay n2 = -1; 4549371c9d4SSatish Balay n5 = -1; 4559371c9d4SSatish Balay n8 = -1; 45647c6ae99SBarry Smith } 45747c6ae99SBarry Smith if (rank % m) { 45847c6ae99SBarry Smith n3 = rank - 1; 4599371c9d4SSatish Balay n6 = n3 + m; 4609371c9d4SSatish Balay if (n6 >= m * n) n6 = -1; 46147c6ae99SBarry Smith } else { 4629371c9d4SSatish Balay n3 = -1; 4639371c9d4SSatish Balay n6 = -1; 46447c6ae99SBarry Smith } 4659371c9d4SSatish Balay n7 = rank + m; 4669371c9d4SSatish Balay if (n7 >= m * n) n7 = -1; 46747c6ae99SBarry Smith 468bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC && by == DM_BOUNDARY_PERIODIC) { 46947c6ae99SBarry Smith /* Modify for Periodic Cases */ 47047c6ae99SBarry Smith /* Handle all four corners */ 47147c6ae99SBarry Smith if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m - 1; 47247c6ae99SBarry Smith if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0; 47347c6ae99SBarry Smith if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size - m; 47447c6ae99SBarry Smith if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size - 1; 47547c6ae99SBarry Smith 47647c6ae99SBarry Smith /* Handle Top and Bottom Sides */ 47747c6ae99SBarry Smith if (n1 < 0) n1 = rank + m * (n - 1); 47847c6ae99SBarry Smith if (n7 < 0) n7 = rank - m * (n - 1); 47947c6ae99SBarry Smith if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 48047c6ae99SBarry Smith if ((n3 >= 0) && (n6 < 0)) n6 = (rank % m) - 1; 48147c6ae99SBarry Smith if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 48247c6ae99SBarry Smith if ((n5 >= 0) && (n8 < 0)) n8 = (rank % m) + 1; 48347c6ae99SBarry Smith 48447c6ae99SBarry Smith /* Handle Left and Right Sides */ 48547c6ae99SBarry Smith if (n3 < 0) n3 = rank + (m - 1); 48647c6ae99SBarry Smith if (n5 < 0) n5 = rank - (m - 1); 48747c6ae99SBarry Smith if ((n1 >= 0) && (n0 < 0)) n0 = rank - 1; 48847c6ae99SBarry Smith if ((n1 >= 0) && (n2 < 0)) n2 = rank - 2 * m + 1; 48947c6ae99SBarry Smith if ((n7 >= 0) && (n6 < 0)) n6 = rank + 2 * m - 1; 49047c6ae99SBarry Smith if ((n7 >= 0) && (n8 < 0)) n8 = rank + 1; 491bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_PERIODIC) { /* Handle Top and Bottom Sides */ 492ce00eea3SSatish Balay if (n1 < 0) n1 = rank + m * (n - 1); 493ce00eea3SSatish Balay if (n7 < 0) n7 = rank - m * (n - 1); 494ce00eea3SSatish Balay if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1; 495ce00eea3SSatish Balay if ((n3 >= 0) && (n6 < 0)) n6 = (rank % m) - 1; 496ce00eea3SSatish Balay if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1; 497ce00eea3SSatish Balay if ((n5 >= 0) && (n8 < 0)) n8 = (rank % m) + 1; 498bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */ 499ce00eea3SSatish Balay if (n3 < 0) n3 = rank + (m - 1); 500ce00eea3SSatish Balay if (n5 < 0) n5 = rank - (m - 1); 501ce00eea3SSatish Balay if ((n1 >= 0) && (n0 < 0)) n0 = rank - 1; 502ce00eea3SSatish Balay if ((n1 >= 0) && (n2 < 0)) n2 = rank - 2 * m + 1; 503ce00eea3SSatish Balay if ((n7 >= 0) && (n6 < 0)) n6 = rank + 2 * m - 1; 504ce00eea3SSatish Balay if ((n7 >= 0) && (n8 < 0)) n8 = rank + 1; 50547c6ae99SBarry Smith } 506ce00eea3SSatish Balay 5079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(9, &dd->neighbors)); 5088865f1eaSKarl Rupp 50947c6ae99SBarry Smith dd->neighbors[0] = n0; 51047c6ae99SBarry Smith dd->neighbors[1] = n1; 51147c6ae99SBarry Smith dd->neighbors[2] = n2; 51247c6ae99SBarry Smith dd->neighbors[3] = n3; 51347c6ae99SBarry Smith dd->neighbors[4] = rank; 51447c6ae99SBarry Smith dd->neighbors[5] = n5; 51547c6ae99SBarry Smith dd->neighbors[6] = n6; 51647c6ae99SBarry Smith dd->neighbors[7] = n7; 51747c6ae99SBarry Smith dd->neighbors[8] = n8; 51847c6ae99SBarry Smith 519d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 52047c6ae99SBarry Smith /* save corner processor numbers */ 5219371c9d4SSatish Balay sn0 = n0; 5229371c9d4SSatish Balay sn2 = n2; 5239371c9d4SSatish Balay sn6 = n6; 5249371c9d4SSatish Balay sn8 = n8; 52547c6ae99SBarry Smith n0 = n2 = n6 = n8 = -1; 52647c6ae99SBarry Smith } 52747c6ae99SBarry Smith 5289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((Xe - Xs) * (Ye - Ys), &idx)); 52947c6ae99SBarry Smith 530ce00eea3SSatish Balay nn = 0; 53147c6ae99SBarry Smith xbase = bases[rank]; 53247c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 53347c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 534ce00eea3SSatish Balay x_t = lx[n0 % m]; 5354ad8454bSPierre Jolivet y_t = ly[n0 / m]; 53647c6ae99SBarry Smith s_t = bases[n0] + x_t * y_t - (s_y - i) * x_t - s_x; 5378865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 53847c6ae99SBarry Smith } 539ac119b13SBarry Smith 54047c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 54147c6ae99SBarry Smith x_t = x; 5424ad8454bSPierre Jolivet y_t = ly[n1 / m]; 54347c6ae99SBarry Smith s_t = bases[n1] + x_t * y_t - (s_y + 1 - i) * x_t; 5448865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 545bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 5468865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = bases[rank] + x * (s_y - i + 1) + j; 54747c6ae99SBarry Smith } 548ac119b13SBarry Smith 54947c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 550ce00eea3SSatish Balay x_t = lx[n2 % m]; 5514ad8454bSPierre Jolivet y_t = ly[n2 / m]; 55247c6ae99SBarry Smith s_t = bases[n2] + x_t * y_t - (s_y + 1 - i) * x_t; 5538865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 55447c6ae99SBarry Smith } 55547c6ae99SBarry Smith } 55647c6ae99SBarry Smith 55747c6ae99SBarry Smith for (i = 0; i < y; i++) { 55847c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 559ce00eea3SSatish Balay x_t = lx[n3 % m]; 56047c6ae99SBarry Smith /* y_t = y; */ 56147c6ae99SBarry Smith s_t = bases[n3] + (i + 1) * x_t - s_x; 5628865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 563bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 5648865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x * i + s_x - j; 56547c6ae99SBarry Smith } 56647c6ae99SBarry Smith 5678865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = xbase++; /* interior */ 56847c6ae99SBarry Smith 56947c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 570ce00eea3SSatish Balay x_t = lx[n5 % m]; 57147c6ae99SBarry Smith /* y_t = y; */ 57247c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 5738865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 574bff4a2f0SMatthew G. Knepley } else if (bx == DM_BOUNDARY_MIRROR) { 5758865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x * (i + 1) - 2 - j; 57647c6ae99SBarry Smith } 57747c6ae99SBarry Smith } 57847c6ae99SBarry Smith 57947c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 58047c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 581ce00eea3SSatish Balay x_t = lx[n6 % m]; 5824ad8454bSPierre Jolivet /* y_t = ly[n6 / m]; */ 58347c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 5848865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 58547c6ae99SBarry Smith } 586ac119b13SBarry Smith 58747c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 58847c6ae99SBarry Smith x_t = x; 5894ad8454bSPierre Jolivet /* y_t = ly[n7 / m]; */ 59047c6ae99SBarry Smith s_t = bases[n7] + (i - 1) * x_t; 5918865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 592bff4a2f0SMatthew G. Knepley } else if (by == DM_BOUNDARY_MIRROR) { 5938865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = bases[rank] + x * (y - i - 1) + j; 59447c6ae99SBarry Smith } 595ac119b13SBarry Smith 59647c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 597ce00eea3SSatish Balay x_t = lx[n8 % m]; 5984ad8454bSPierre Jolivet /* y_t = ly[n8 / m]; */ 59947c6ae99SBarry Smith s_t = bases[n8] + (i - 1) * x_t; 6008865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 60147c6ae99SBarry Smith } 60247c6ae99SBarry Smith } 60347c6ae99SBarry Smith 6049566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(comm, dof, nn, idx, PETSC_USE_POINTER, &from)); 6059566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(global, from, local, to, >ol)); 6069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 6079566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 60847c6ae99SBarry Smith 609d9c9ebe5SPeter Brune if (stencil_type == DMDA_STENCIL_STAR) { 6109371c9d4SSatish Balay n0 = sn0; 6119371c9d4SSatish Balay n2 = sn2; 6129371c9d4SSatish Balay n6 = sn6; 6139371c9d4SSatish Balay n8 = sn8; 614ce00eea3SSatish Balay } 615ce00eea3SSatish Balay 616f4f49eeaSPierre Jolivet if ((stencil_type == DMDA_STENCIL_STAR) || (bx && bx != DM_BOUNDARY_PERIODIC) || (by && by != DM_BOUNDARY_PERIODIC)) { 61747c6ae99SBarry Smith /* 61847c6ae99SBarry Smith Recompute the local to global mappings, this time keeping the 619ce00eea3SSatish Balay information about the cross corner processor numbers and any ghosted 620ce00eea3SSatish Balay but not periodic indices. 62147c6ae99SBarry Smith */ 62247c6ae99SBarry Smith nn = 0; 62347c6ae99SBarry Smith xbase = bases[rank]; 62447c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 62547c6ae99SBarry Smith if (n0 >= 0) { /* left below */ 626ce00eea3SSatish Balay x_t = lx[n0 % m]; 6274ad8454bSPierre Jolivet y_t = ly[n0 / m]; 62847c6ae99SBarry Smith s_t = bases[n0] + x_t * y_t - (s_y - i) * x_t - s_x; 6298865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 630ce00eea3SSatish Balay } else if (xs - Xs > 0 && ys - Ys > 0) { 6318865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 63247c6ae99SBarry Smith } 63347c6ae99SBarry Smith if (n1 >= 0) { /* directly below */ 63447c6ae99SBarry Smith x_t = x; 6354ad8454bSPierre Jolivet y_t = ly[n1 / m]; 63647c6ae99SBarry Smith s_t = bases[n1] + x_t * y_t - (s_y + 1 - i) * x_t; 6378865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 638ce00eea3SSatish Balay } else if (ys - Ys > 0) { 639bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 6408865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = bases[rank] + x * (s_y - i + 1) + j; 641624904c4SBarry Smith } else { 6428865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = -1; 64347c6ae99SBarry Smith } 644624904c4SBarry Smith } 64547c6ae99SBarry Smith if (n2 >= 0) { /* right below */ 646ce00eea3SSatish Balay x_t = lx[n2 % m]; 6474ad8454bSPierre Jolivet y_t = ly[n2 / m]; 64847c6ae99SBarry Smith s_t = bases[n2] + x_t * y_t - (s_y + 1 - i) * x_t; 6498865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 650ce00eea3SSatish Balay } else if (Xe - xe > 0 && ys - Ys > 0) { 6518865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 65247c6ae99SBarry Smith } 65347c6ae99SBarry Smith } 65447c6ae99SBarry Smith 65547c6ae99SBarry Smith for (i = 0; i < y; i++) { 65647c6ae99SBarry Smith if (n3 >= 0) { /* directly left */ 657ce00eea3SSatish Balay x_t = lx[n3 % m]; 65847c6ae99SBarry Smith /* y_t = y; */ 65947c6ae99SBarry Smith s_t = bases[n3] + (i + 1) * x_t - s_x; 6608865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 661ce00eea3SSatish Balay } else if (xs - Xs > 0) { 662bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 6638865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x * i + s_x - j; 664624904c4SBarry Smith } else { 6658865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 66647c6ae99SBarry Smith } 667624904c4SBarry Smith } 66847c6ae99SBarry Smith 6698865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = xbase++; /* interior */ 67047c6ae99SBarry Smith 67147c6ae99SBarry Smith if (n5 >= 0) { /* directly right */ 672ce00eea3SSatish Balay x_t = lx[n5 % m]; 67347c6ae99SBarry Smith /* y_t = y; */ 67447c6ae99SBarry Smith s_t = bases[n5] + (i)*x_t; 6758865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 676ce00eea3SSatish Balay } else if (Xe - xe > 0) { 677bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_MIRROR) { 6788865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x * (i + 1) - 2 - j; 679624904c4SBarry Smith } else { 6808865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 68147c6ae99SBarry Smith } 68247c6ae99SBarry Smith } 683624904c4SBarry Smith } 68447c6ae99SBarry Smith 68547c6ae99SBarry Smith for (i = 1; i <= s_y; i++) { 68647c6ae99SBarry Smith if (n6 >= 0) { /* left above */ 687ce00eea3SSatish Balay x_t = lx[n6 % m]; 6884ad8454bSPierre Jolivet /* y_t = ly[n6 / m]; */ 68947c6ae99SBarry Smith s_t = bases[n6] + (i)*x_t - s_x; 6908865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 691ce00eea3SSatish Balay } else if (xs - Xs > 0 && Ye - ye > 0) { 6928865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 69347c6ae99SBarry Smith } 69447c6ae99SBarry Smith if (n7 >= 0) { /* directly above */ 69547c6ae99SBarry Smith x_t = x; 6964ad8454bSPierre Jolivet /* y_t = ly[n7 / m]; */ 69747c6ae99SBarry Smith s_t = bases[n7] + (i - 1) * x_t; 6988865f1eaSKarl Rupp for (j = 0; j < x_t; j++) idx[nn++] = s_t++; 699ce00eea3SSatish Balay } else if (Ye - ye > 0) { 700bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_MIRROR) { 7018865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = bases[rank] + x * (y - i - 1) + j; 702624904c4SBarry Smith } else { 7038865f1eaSKarl Rupp for (j = 0; j < x; j++) idx[nn++] = -1; 70447c6ae99SBarry Smith } 705624904c4SBarry Smith } 70647c6ae99SBarry Smith if (n8 >= 0) { /* right above */ 707ce00eea3SSatish Balay x_t = lx[n8 % m]; 7084ad8454bSPierre Jolivet /* y_t = ly[n8 / m]; */ 70947c6ae99SBarry Smith s_t = bases[n8] + (i - 1) * x_t; 7108865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = s_t++; 711ce00eea3SSatish Balay } else if (Xe - xe > 0 && Ye - ye > 0) { 7128865f1eaSKarl Rupp for (j = 0; j < s_x; j++) idx[nn++] = -1; 71347c6ae99SBarry Smith } 71447c6ae99SBarry Smith } 71547c6ae99SBarry Smith } 716ce00eea3SSatish Balay /* 717ce00eea3SSatish Balay Set the local to global ordering in the global vector, this allows use 718ce00eea3SSatish Balay of VecSetValuesLocal(). 719ce00eea3SSatish Balay */ 7209566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, dof, nn, idx, PETSC_OWN_POINTER, &da->ltogmap)); 72147c6ae99SBarry Smith 7229566063dSJacob Faibussowitsch PetscCall(PetscFree2(bases, ldims)); 7239371c9d4SSatish Balay dd->m = m; 7249371c9d4SSatish Balay dd->n = n; 725ce00eea3SSatish Balay /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */ 7269371c9d4SSatish Balay dd->xs = xs * dof; 7279371c9d4SSatish Balay dd->xe = xe * dof; 7289371c9d4SSatish Balay dd->ys = ys; 7299371c9d4SSatish Balay dd->ye = ye; 7309371c9d4SSatish Balay dd->zs = 0; 7319371c9d4SSatish Balay dd->ze = 1; 7329371c9d4SSatish Balay dd->Xs = Xs * dof; 7339371c9d4SSatish Balay dd->Xe = Xe * dof; 7349371c9d4SSatish Balay dd->Ys = Ys; 7359371c9d4SSatish Balay dd->Ye = Ye; 7369371c9d4SSatish Balay dd->Zs = 0; 7379371c9d4SSatish Balay dd->Ze = 1; 73847c6ae99SBarry Smith 7399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 7409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 74147c6ae99SBarry Smith 74247c6ae99SBarry Smith dd->gtol = gtol; 74347c6ae99SBarry Smith dd->base = base; 7449a42bb27SBarry Smith da->ops->view = DMView_DA_2d; 7450298fd71SBarry Smith dd->ltol = NULL; 7460298fd71SBarry Smith dd->ao = NULL; 7473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 74847c6ae99SBarry Smith } 74947c6ae99SBarry Smith 750cc4c1da9SBarry Smith /*@ 751aa219208SBarry Smith DMDACreate2d - Creates an object that will manage the communication of two-dimensional 75212b4a537SBarry Smith regular array data that is distributed across one or more MPI processes. 75347c6ae99SBarry Smith 754d083f849SBarry Smith Collective 75547c6ae99SBarry Smith 75647c6ae99SBarry Smith Input Parameters: 75747c6ae99SBarry Smith + comm - MPI communicator 758a4e35b19SJacob Faibussowitsch . bx - type of ghost nodes the x array have. Use one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`. 759a4e35b19SJacob Faibussowitsch . by - type of ghost nodes the y array have. Use one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`. 760dce8aebaSBarry Smith . stencil_type - stencil type. Use either `DMDA_STENCIL_BOX` or `DMDA_STENCIL_STAR`. 761a4e35b19SJacob Faibussowitsch . M - global dimension in x direction of the array 762a4e35b19SJacob Faibussowitsch . N - global dimension in y direction of the array 763a4e35b19SJacob Faibussowitsch . m - corresponding number of processors in x dimension (or `PETSC_DECIDE` to have calculated) 764a4e35b19SJacob Faibussowitsch . n - corresponding number of processors in y dimension (or `PETSC_DECIDE` to have calculated) 76547c6ae99SBarry Smith . dof - number of degrees of freedom per node 76647c6ae99SBarry Smith . s - stencil width 767a4e35b19SJacob Faibussowitsch . lx - arrays containing the number of nodes in each cell along the x coordinates, or `NULL`. 768a4e35b19SJacob Faibussowitsch - ly - arrays containing the number of nodes in each cell along the y coordinates, or `NULL`. 76947c6ae99SBarry Smith 77047c6ae99SBarry Smith Output Parameter: 77147c6ae99SBarry Smith . da - the resulting distributed array object 77247c6ae99SBarry Smith 773dce8aebaSBarry Smith Options Database Keys: 774dce8aebaSBarry Smith + -dm_view - Calls `DMView()` at the conclusion of `DMDACreate2d()` 775e43dc8daSBarry Smith . -da_grid_x <nx> - number of grid points in x direction 776e43dc8daSBarry Smith . -da_grid_y <ny> - number of grid points in y direction 77747c6ae99SBarry Smith . -da_processors_x <nx> - number of processors in x direction 77847c6ae99SBarry Smith . -da_processors_y <ny> - number of processors in y direction 779fe58071bSMatthew G. Knepley . -da_bd_x <bx> - boundary type in x direction 780fe58071bSMatthew G. Knepley . -da_bd_y <by> - boundary type in y direction 781fe58071bSMatthew G. Knepley . -da_bd_all <bt> - boundary type in all directions 782e0f5d30fSBarry Smith . -da_refine_x <rx> - refinement ratio in x direction 783e0f5d30fSBarry Smith . -da_refine_y <ry> - refinement ratio in y direction 78412b4a537SBarry Smith - -da_refine <n> - refine the `DMDA` n times before creating 785e0f5d30fSBarry Smith 78647c6ae99SBarry Smith Level: beginner 78747c6ae99SBarry Smith 78847c6ae99SBarry Smith Notes: 789a4e35b19SJacob Faibussowitsch If `lx` or `ly` are non-null, these must be of length as `m` and `n`, and the corresponding 790a4e35b19SJacob Faibussowitsch `m` and `n` cannot be `PETSC_DECIDE`. The sum of the `lx` entries must be `M`, and the sum of 791a4e35b19SJacob Faibussowitsch the `ly` entries must be `N`. 792a4e35b19SJacob Faibussowitsch 793dce8aebaSBarry Smith The stencil type `DMDA_STENCIL_STAR` with width 1 corresponds to the 794dce8aebaSBarry Smith standard 5-pt stencil, while `DMDA_STENCIL_BOX` with width 1 denotes 79547c6ae99SBarry Smith the standard 9-pt stencil. 79647c6ae99SBarry Smith 797dce8aebaSBarry Smith The array data itself is NOT stored in the `DMDA`, it is stored in `Vec` objects; 798dce8aebaSBarry Smith The appropriate vector objects can be obtained with calls to `DMCreateGlobalVector()` 799dce8aebaSBarry Smith and DMCreateLocalVector() and calls to `VecDuplicate()` if more are needed. 80047c6ae99SBarry Smith 801dce8aebaSBarry Smith You must call `DMSetUp()` after this call before using this `DM`. 802897f7067SBarry Smith 80312b4a537SBarry Smith To use the options database to change values in the `DMDA` call `DMSetFromOptions()` after this call 804dce8aebaSBarry Smith but before `DMSetUp()`. 805897f7067SBarry Smith 80612b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDestroy()`, `DMView()`, `DMDACreate1d()`, `DMDACreate3d()`, `DMGlobalToLocalBegin()`, `DMDAGetRefinementFactor()`, 807db781477SPatrick Sanan `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToLocalBegin()`, `DMLocalToLocalEnd()`, `DMDASetRefinementFactor()`, 808db781477SPatrick Sanan `DMDAGetInfo()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMDACreateNaturalVector()`, `DMLoad()`, `DMDAGetOwnershipRanges()`, 80912b4a537SBarry Smith `DMStagCreate2d()`, `DMBoundaryType` 81047c6ae99SBarry Smith @*/ 811d71ae5a4SJacob Faibussowitsch 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) 812d71ae5a4SJacob Faibussowitsch { 81347c6ae99SBarry Smith PetscFunctionBegin; 8149566063dSJacob Faibussowitsch PetscCall(DMDACreate(comm, da)); 8159566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*da, 2)); 8169566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(*da, M, N, 1)); 8179566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(*da, m, n, PETSC_DECIDE)); 8189566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(*da, bx, by, DM_BOUNDARY_NONE)); 8199566063dSJacob Faibussowitsch PetscCall(DMDASetDof(*da, dof)); 8209566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(*da, stencil_type)); 8219566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(*da, s)); 8229566063dSJacob Faibussowitsch PetscCall(DMDASetOwnershipRanges(*da, lx, ly, NULL)); 8233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82447c6ae99SBarry Smith } 825