1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 2e30e807fSPeter Brune 3063654ddSPeter Brune /*@ 412b4a537SBarry Smith DMDACreatePatchIS - Creates an index set corresponding to a logically rectangular patch of the `DMDA`. 5ff9846d9SPeter Brune 63b5e53cdSSajid Ali Collective 7063654ddSPeter Brune 8063654ddSPeter Brune Input Parameters: 9dce8aebaSBarry Smith + da - the `DMDA` 1012b4a537SBarry Smith . lower - a `MatStencil` with i, j and k entries corresponding to the lower corner of the patch 1112b4a537SBarry Smith . upper - a `MatStencil` with i, j and k entries corresponding to the upper corner of the patch 1212b4a537SBarry Smith - offproc - indicate whether the returned `IS` will contain off process indices 13063654ddSPeter Brune 142fe279fdSBarry Smith Output Parameter: 15dce8aebaSBarry Smith . is - the `IS` corresponding to the patch 16063654ddSPeter Brune 17063654ddSPeter Brune Level: developer 18063654ddSPeter Brune 193b5e53cdSSajid Ali Notes: 2012b4a537SBarry Smith This routine always returns an `IS` on the `DMDA` communicator. 213b5e53cdSSajid Ali 2212b4a537SBarry Smith If `offproc` is set to `PETSC_TRUE`, 2312b4a537SBarry Smith the routine returns an `IS` with all the indices requested regardless of whether these indices 2412b4a537SBarry Smith are present on the requesting MPI process or not. Thus, it is upon the caller to ensure that 2512b4a537SBarry Smith the indices returned in this mode are appropriate. 2612b4a537SBarry Smith 2712b4a537SBarry Smith If `offproc` is set to `PETSC_FALSE`, 2812b4a537SBarry Smith the `IS` only returns the subset of indices that are present on the requesting MPI process and there 2912b4a537SBarry Smith is no duplication of indices between multiple MPI processes. 3012b4a537SBarry Smith 3112b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMCreateDomainDecompositionScatters()` 32063654ddSPeter Brune @*/ 33d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDACreatePatchIS(DM da, MatStencil *lower, MatStencil *upper, IS *is, PetscBool offproc) 34d71ae5a4SJacob Faibussowitsch { 35063654ddSPeter Brune PetscInt ms = 0, ns = 0, ps = 0; 363b5e53cdSSajid Ali PetscInt mw = 0, nw = 0, pw = 0; 37063654ddSPeter Brune PetscInt me = 1, ne = 1, pe = 1; 38063654ddSPeter Brune PetscInt mr = 0, nr = 0, pr = 0; 39ff9846d9SPeter Brune PetscInt ii, jj, kk; 40063654ddSPeter Brune PetscInt si, sj, sk; 41*63bfac88SBarry Smith PetscInt i, k = 0, l, idx = 0; 42ff9846d9SPeter Brune PetscInt base; 43063654ddSPeter Brune PetscInt xm = 1, ym = 1, zm = 1; 44ff9846d9SPeter Brune PetscInt ox, oy, oz; 45063654ddSPeter Brune PetscInt m, n, p, M, N, P, dof; 463b5e53cdSSajid Ali const PetscInt *lx, *ly, *lz; 47063654ddSPeter Brune PetscInt nindices; 48063654ddSPeter Brune PetscInt *indices; 49063654ddSPeter Brune DM_DA *dd = (DM_DA *)da->data; 503b5e53cdSSajid Ali PetscBool skip_i = PETSC_TRUE, skip_j = PETSC_TRUE, skip_k = PETSC_TRUE; 513b5e53cdSSajid Ali PetscBool valid_j = PETSC_FALSE, valid_k = PETSC_FALSE; /* DMDA has at least 1 dimension */ 52ff9846d9SPeter Brune 530adebc6cSBarry Smith PetscFunctionBegin; 549371c9d4SSatish Balay M = dd->M; 559371c9d4SSatish Balay N = dd->N; 569371c9d4SSatish Balay P = dd->P; 579371c9d4SSatish Balay m = dd->m; 589371c9d4SSatish Balay n = dd->n; 599371c9d4SSatish Balay p = dd->p; 60063654ddSPeter Brune dof = dd->w; 613b5e53cdSSajid Ali 623b5e53cdSSajid Ali nindices = -1; 633b5e53cdSSajid Ali if (PetscLikely(upper->i - lower->i)) { 643b5e53cdSSajid Ali nindices = nindices * (upper->i - lower->i); 653b5e53cdSSajid Ali skip_i = PETSC_FALSE; 663b5e53cdSSajid Ali } 673b5e53cdSSajid Ali if (N > 1) { 683b5e53cdSSajid Ali valid_j = PETSC_TRUE; 693b5e53cdSSajid Ali if (PetscLikely(upper->j - lower->j)) { 703b5e53cdSSajid Ali nindices = nindices * (upper->j - lower->j); 713b5e53cdSSajid Ali skip_j = PETSC_FALSE; 723b5e53cdSSajid Ali } 733b5e53cdSSajid Ali } 743b5e53cdSSajid Ali if (P > 1) { 753b5e53cdSSajid Ali valid_k = PETSC_TRUE; 763b5e53cdSSajid Ali if (PetscLikely(upper->k - lower->k)) { 773b5e53cdSSajid Ali nindices = nindices * (upper->k - lower->k); 783b5e53cdSSajid Ali skip_k = PETSC_FALSE; 793b5e53cdSSajid Ali } 803b5e53cdSSajid Ali } 813b5e53cdSSajid Ali if (PetscLikely(nindices < 0)) { 823b5e53cdSSajid Ali if (PetscUnlikely(skip_i && skip_j && skip_k)) { 833b5e53cdSSajid Ali nindices = 0; 843b5e53cdSSajid Ali } else nindices = nindices * (-1); 853b5e53cdSSajid Ali } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Lower and Upper stencils are identical! Please check inputs."); 863b5e53cdSSajid Ali 879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nindices * dof, &indices)); 889566063dSJacob Faibussowitsch PetscCall(DMDAGetOffset(da, &ox, &oy, &oz, NULL, NULL, NULL)); 893b5e53cdSSajid Ali 903b5e53cdSSajid Ali if (!valid_k) { 913b5e53cdSSajid Ali upper->k = 0; 923b5e53cdSSajid Ali lower->k = 0; 933b5e53cdSSajid Ali } 943b5e53cdSSajid Ali if (!valid_j) { 953b5e53cdSSajid Ali upper->j = 0; 963b5e53cdSSajid Ali lower->j = 0; 973b5e53cdSSajid Ali } 983b5e53cdSSajid Ali 993b5e53cdSSajid Ali if (offproc) { 1009566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, &lz)); 101063654ddSPeter Brune /* start at index 0 on processor 0 */ 102063654ddSPeter Brune mr = 0; 103063654ddSPeter Brune nr = 0; 104063654ddSPeter Brune pr = 0; 105063654ddSPeter Brune ms = 0; 106063654ddSPeter Brune ns = 0; 107063654ddSPeter Brune ps = 0; 108063654ddSPeter Brune if (lx) me = lx[0]; 109063654ddSPeter Brune if (ly) ne = ly[0]; 110063654ddSPeter Brune if (lz) pe = lz[0]; 1113b5e53cdSSajid Ali /* 1123b5e53cdSSajid Ali If no indices are to be returned, create an empty is, 1133b5e53cdSSajid Ali this prevents hanging in while loops 1143b5e53cdSSajid Ali */ 1153b5e53cdSSajid Ali if (skip_i && skip_j && skip_k) goto createis; 1163b5e53cdSSajid Ali /* 1173b5e53cdSSajid Ali do..while loops to ensure the block gets entered once, 1183b5e53cdSSajid Ali regardless of control condition being met, necessary for 1193b5e53cdSSajid Ali cases when a subset of skip_i/j/k is true 1203b5e53cdSSajid Ali */ 1219371c9d4SSatish Balay if (skip_k) k = upper->k - oz; 1229371c9d4SSatish Balay else k = lower->k - oz; 1233b5e53cdSSajid Ali do { 124*63bfac88SBarry Smith PetscInt j; 125*63bfac88SBarry Smith 1269371c9d4SSatish Balay if (skip_j) j = upper->j - oy; 1279371c9d4SSatish Balay else j = lower->j - oy; 1283b5e53cdSSajid Ali do { 1299371c9d4SSatish Balay if (skip_i) i = upper->i - ox; 1309371c9d4SSatish Balay else i = lower->i - ox; 1313b5e53cdSSajid Ali do { 132063654ddSPeter Brune /* "actual" indices rather than ones outside of the domain */ 133063654ddSPeter Brune ii = i; 134063654ddSPeter Brune jj = j; 135063654ddSPeter Brune kk = k; 136063654ddSPeter Brune if (ii < 0) ii = M + ii; 137063654ddSPeter Brune if (jj < 0) jj = N + jj; 138063654ddSPeter Brune if (kk < 0) kk = P + kk; 139063654ddSPeter Brune if (ii > M - 1) ii = ii - M; 140063654ddSPeter Brune if (jj > N - 1) jj = jj - N; 141063654ddSPeter Brune if (kk > P - 1) kk = kk - P; 142063654ddSPeter Brune /* gone out of processor range on x axis */ 143063654ddSPeter Brune while (ii > me - 1 || ii < ms) { 144063654ddSPeter Brune if (mr == m - 1) { 145063654ddSPeter Brune ms = 0; 146063654ddSPeter Brune me = lx[0]; 147063654ddSPeter Brune mr = 0; 148063654ddSPeter Brune } else { 149063654ddSPeter Brune mr++; 150063654ddSPeter Brune ms = me; 151063654ddSPeter Brune me += lx[mr]; 152063654ddSPeter Brune } 153063654ddSPeter Brune } 154063654ddSPeter Brune /* gone out of processor range on y axis */ 155063654ddSPeter Brune while (jj > ne - 1 || jj < ns) { 156063654ddSPeter Brune if (nr == n - 1) { 157063654ddSPeter Brune ns = 0; 158063654ddSPeter Brune ne = ly[0]; 159063654ddSPeter Brune nr = 0; 160063654ddSPeter Brune } else { 161063654ddSPeter Brune nr++; 162063654ddSPeter Brune ns = ne; 163063654ddSPeter Brune ne += ly[nr]; 164063654ddSPeter Brune } 165063654ddSPeter Brune } 166063654ddSPeter Brune /* gone out of processor range on z axis */ 167063654ddSPeter Brune while (kk > pe - 1 || kk < ps) { 168063654ddSPeter Brune if (pr == p - 1) { 169063654ddSPeter Brune ps = 0; 170063654ddSPeter Brune pe = lz[0]; 171063654ddSPeter Brune pr = 0; 172063654ddSPeter Brune } else { 173063654ddSPeter Brune pr++; 174063654ddSPeter Brune ps = pe; 175063654ddSPeter Brune pe += lz[pr]; 176063654ddSPeter Brune } 177063654ddSPeter Brune } 178063654ddSPeter Brune /* compute the vector base on owning processor */ 179063654ddSPeter Brune xm = me - ms; 180063654ddSPeter Brune ym = ne - ns; 181063654ddSPeter Brune zm = pe - ps; 182b0bff951SPeter Brune base = ms * ym * zm + ns * M * zm + ps * M * N; 183063654ddSPeter Brune /* compute the local coordinates on owning processor */ 184063654ddSPeter Brune si = ii - ms; 185063654ddSPeter Brune sj = jj - ns; 186063654ddSPeter Brune sk = kk - ps; 187063654ddSPeter Brune for (l = 0; l < dof; l++) { 188063654ddSPeter Brune indices[idx] = l + dof * (base + si + xm * sj + xm * ym * sk); 189ff9846d9SPeter Brune idx++; 190ff9846d9SPeter Brune } 1913b5e53cdSSajid Ali i++; 1923b5e53cdSSajid Ali } while (i < upper->i - ox); 1933b5e53cdSSajid Ali j++; 1943b5e53cdSSajid Ali } while (j < upper->j - oy); 1953b5e53cdSSajid Ali k++; 1963b5e53cdSSajid Ali } while (k < upper->k - oz); 1973b5e53cdSSajid Ali } 1983b5e53cdSSajid Ali 1993b5e53cdSSajid Ali if (!offproc) { 2009566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &ms, &ns, &ps, &mw, &nw, &pw)); 2013b5e53cdSSajid Ali me = ms + mw; 2023b5e53cdSSajid Ali if (N > 1) ne = ns + nw; 2033b5e53cdSSajid Ali if (P > 1) pe = ps + pw; 2043b5e53cdSSajid Ali /* Account for DM offsets */ 2059371c9d4SSatish Balay ms = ms - ox; 2069371c9d4SSatish Balay me = me - ox; 2079371c9d4SSatish Balay ns = ns - oy; 2089371c9d4SSatish Balay ne = ne - oy; 2099371c9d4SSatish Balay ps = ps - oz; 2109371c9d4SSatish Balay pe = pe - oz; 2113b5e53cdSSajid Ali 2123b5e53cdSSajid Ali /* compute the vector base on owning processor */ 2133b5e53cdSSajid Ali xm = me - ms; 2143b5e53cdSSajid Ali ym = ne - ns; 2153b5e53cdSSajid Ali zm = pe - ps; 2163b5e53cdSSajid Ali base = ms * ym * zm + ns * M * zm + ps * M * N; 2173b5e53cdSSajid Ali /* 2183b5e53cdSSajid Ali if no indices are to be returned, create an empty is, 2193b5e53cdSSajid Ali this prevents hanging in while loops 2203b5e53cdSSajid Ali */ 2213b5e53cdSSajid Ali if (skip_i && skip_j && skip_k) goto createis; 2223b5e53cdSSajid Ali /* 2233b5e53cdSSajid Ali do..while loops to ensure the block gets entered once, 2243b5e53cdSSajid Ali regardless of control condition being met, necessary for 2253b5e53cdSSajid Ali cases when a subset of skip_i/j/k is true 2263b5e53cdSSajid Ali */ 2279371c9d4SSatish Balay if (skip_k) k = upper->k - oz; 2289371c9d4SSatish Balay else k = lower->k - oz; 2293b5e53cdSSajid Ali do { 230*63bfac88SBarry Smith PetscInt j; 231*63bfac88SBarry Smith 2329371c9d4SSatish Balay if (skip_j) j = upper->j - oy; 2339371c9d4SSatish Balay else j = lower->j - oy; 2343b5e53cdSSajid Ali do { 2359371c9d4SSatish Balay if (skip_i) i = upper->i - ox; 2369371c9d4SSatish Balay else i = lower->i - ox; 2373b5e53cdSSajid Ali do { 2383b5e53cdSSajid Ali if (k >= ps && k <= pe - 1) { 2393b5e53cdSSajid Ali if (j >= ns && j <= ne - 1) { 2403b5e53cdSSajid Ali if (i >= ms && i <= me - 1) { 2413b5e53cdSSajid Ali /* compute the local coordinates on owning processor */ 2423b5e53cdSSajid Ali si = i - ms; 2433b5e53cdSSajid Ali sj = j - ns; 2443b5e53cdSSajid Ali sk = k - ps; 2453b5e53cdSSajid Ali for (l = 0; l < dof; l++) { 2463b5e53cdSSajid Ali indices[idx] = l + dof * (base + si + xm * sj + xm * ym * sk); 2473b5e53cdSSajid Ali idx++; 248ff9846d9SPeter Brune } 249ff9846d9SPeter Brune } 250063654ddSPeter Brune } 2513b5e53cdSSajid Ali } 2523b5e53cdSSajid Ali i++; 2533b5e53cdSSajid Ali } while (i < upper->i - ox); 2543b5e53cdSSajid Ali j++; 2553b5e53cdSSajid Ali } while (j < upper->j - oy); 2563b5e53cdSSajid Ali k++; 2573b5e53cdSSajid Ali } while (k < upper->k - oz); 2583b5e53cdSSajid Ali 2599566063dSJacob Faibussowitsch PetscCall(PetscRealloc((size_t)(idx * sizeof(PetscInt)), (void *)&indices)); 2603b5e53cdSSajid Ali } 2613b5e53cdSSajid Ali 2623b5e53cdSSajid Ali createis: 2639566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)da), idx, indices, PETSC_OWN_POINTER, is)); 2643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 265ff9846d9SPeter Brune } 266ff9846d9SPeter Brune 26766976f2fSJacob Faibussowitsch static PetscErrorCode DMDASubDomainDA_Private(DM dm, PetscInt *nlocal, DM **sdm) 268d71ae5a4SJacob Faibussowitsch { 26990c77898SPeter Brune DM *da; 27090c77898SPeter Brune PetscInt dim, size, i, j, k, idx; 271e30e807fSPeter Brune DMDALocalInfo info; 272ff9846d9SPeter Brune PetscInt xsize, ysize, zsize; 273e30e807fSPeter Brune PetscInt xo, yo, zo; 27490c77898SPeter Brune PetscInt xs, ys, zs; 27590c77898SPeter Brune PetscInt xm = 1, ym = 1, zm = 1; 2767ddda789SPeter Brune PetscInt xol, yol, zol; 27790c77898SPeter Brune PetscInt m = 1, n = 1, p = 1; 27890c77898SPeter Brune PetscInt M, N, P; 27990c77898SPeter Brune PetscInt pm, mtmp; 280e30e807fSPeter Brune 281e30e807fSPeter Brune PetscFunctionBegin; 2829566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 2839566063dSJacob Faibussowitsch PetscCall(DMDAGetOverlap(dm, &xol, &yol, &zol)); 2849566063dSJacob Faibussowitsch PetscCall(DMDAGetNumLocalSubDomains(dm, &size)); 2859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &da)); 2867ddda789SPeter Brune 28790c77898SPeter Brune if (nlocal) *nlocal = size; 288e30e807fSPeter Brune 28990c77898SPeter Brune dim = info.dim; 290e30e807fSPeter Brune 29190c77898SPeter Brune M = info.xm; 29290c77898SPeter Brune N = info.ym; 29390c77898SPeter Brune P = info.zm; 29490c77898SPeter Brune 29590c77898SPeter Brune if (dim == 1) { 29690c77898SPeter Brune m = size; 29790c77898SPeter Brune } else if (dim == 2) { 29890c77898SPeter Brune m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)N))); 29990c77898SPeter Brune while (m > 0) { 30090c77898SPeter Brune n = size / m; 30190c77898SPeter Brune if (m * n * p == size) break; 30290c77898SPeter Brune m--; 30390c77898SPeter Brune } 30490c77898SPeter Brune } else if (dim == 3) { 3059371c9d4SSatish Balay n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N * N) * ((PetscReal)size) / ((PetscReal)P * M), (PetscReal)(1. / 3.))); 3069371c9d4SSatish Balay if (!n) n = 1; 30790c77898SPeter Brune while (n > 0) { 30890c77898SPeter Brune pm = size / n; 30990c77898SPeter Brune if (n * pm == size) break; 31090c77898SPeter Brune n--; 31190c77898SPeter Brune } 31290c77898SPeter Brune if (!n) n = 1; 31390c77898SPeter Brune m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)P * n))); 31490c77898SPeter Brune if (!m) m = 1; 31590c77898SPeter Brune while (m > 0) { 31690c77898SPeter Brune p = size / (m * n); 31790c77898SPeter Brune if (m * n * p == size) break; 31890c77898SPeter Brune m--; 31990c77898SPeter Brune } 3209371c9d4SSatish Balay if (M > P && m < p) { 3219371c9d4SSatish Balay mtmp = m; 3229371c9d4SSatish Balay m = p; 3239371c9d4SSatish Balay p = mtmp; 3249371c9d4SSatish Balay } 32590c77898SPeter Brune } 32690c77898SPeter Brune 32790c77898SPeter Brune zs = info.zs; 32890c77898SPeter Brune idx = 0; 32990c77898SPeter Brune for (k = 0; k < p; k++) { 33090c77898SPeter Brune ys = info.ys; 33190c77898SPeter Brune for (j = 0; j < n; j++) { 33290c77898SPeter Brune xs = info.xs; 33390c77898SPeter Brune for (i = 0; i < m; i++) { 33490c77898SPeter Brune if (dim == 1) { 33590c77898SPeter Brune xm = M / m + ((M % m) > i); 33690c77898SPeter Brune } else if (dim == 2) { 33790c77898SPeter Brune xm = M / m + ((M % m) > i); 33890c77898SPeter Brune ym = N / n + ((N % n) > j); 33990c77898SPeter Brune } else if (dim == 3) { 34090c77898SPeter Brune xm = M / m + ((M % m) > i); 34190c77898SPeter Brune ym = N / n + ((N % n) > j); 34290c77898SPeter Brune zm = P / p + ((P % p) > k); 34390c77898SPeter Brune } 34490c77898SPeter Brune 34590c77898SPeter Brune xsize = xm; 34690c77898SPeter Brune ysize = ym; 34790c77898SPeter Brune zsize = zm; 34890c77898SPeter Brune xo = xs; 34990c77898SPeter Brune yo = ys; 35090c77898SPeter Brune zo = zs; 35190c77898SPeter Brune 352f4f49eeaSPierre Jolivet PetscCall(DMDACreate(PETSC_COMM_SELF, &da[idx])); 3539566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(da[idx], "sub_")); 3549566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da[idx], info.dim)); 3559566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da[idx], info.dof)); 35690c77898SPeter Brune 3579566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da[idx], info.st)); 3589566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da[idx], info.sw)); 35990c77898SPeter Brune 360bff4a2f0SMatthew G. Knepley if (info.bx == DM_BOUNDARY_PERIODIC || (xs != 0)) { 3617ddda789SPeter Brune xsize += xol; 3627ddda789SPeter Brune xo -= xol; 363e30e807fSPeter Brune } 364bff4a2f0SMatthew G. Knepley if (info.by == DM_BOUNDARY_PERIODIC || (ys != 0)) { 3657ddda789SPeter Brune ysize += yol; 3667ddda789SPeter Brune yo -= yol; 367e30e807fSPeter Brune } 368bff4a2f0SMatthew G. Knepley if (info.bz == DM_BOUNDARY_PERIODIC || (zs != 0)) { 3697ddda789SPeter Brune zsize += zol; 3707ddda789SPeter Brune zo -= zol; 371e30e807fSPeter Brune } 372e30e807fSPeter Brune 373bff4a2f0SMatthew G. Knepley if (info.bx == DM_BOUNDARY_PERIODIC || (xs + xm != info.mx)) xsize += xol; 374bff4a2f0SMatthew G. Knepley if (info.by == DM_BOUNDARY_PERIODIC || (ys + ym != info.my)) ysize += yol; 375bff4a2f0SMatthew G. Knepley if (info.bz == DM_BOUNDARY_PERIODIC || (zs + zm != info.mz)) zsize += zol; 3768865f1eaSKarl Rupp 377bff4a2f0SMatthew G. Knepley if (info.bx != DM_BOUNDARY_PERIODIC) { 37890c77898SPeter Brune if (xo < 0) { 37990c77898SPeter Brune xsize += xo; 38090c77898SPeter Brune xo = 0; 38190c77898SPeter Brune } 382ad540459SPierre Jolivet if (xo + xsize > info.mx - 1) xsize -= xo + xsize - info.mx; 38390c77898SPeter Brune } 384bff4a2f0SMatthew G. Knepley if (info.by != DM_BOUNDARY_PERIODIC) { 38590c77898SPeter Brune if (yo < 0) { 38690c77898SPeter Brune ysize += yo; 38790c77898SPeter Brune yo = 0; 38890c77898SPeter Brune } 389ad540459SPierre Jolivet if (yo + ysize > info.my - 1) ysize -= yo + ysize - info.my; 39090c77898SPeter Brune } 391bff4a2f0SMatthew G. Knepley if (info.bz != DM_BOUNDARY_PERIODIC) { 39290c77898SPeter Brune if (zo < 0) { 39390c77898SPeter Brune zsize += zo; 39490c77898SPeter Brune zo = 0; 39590c77898SPeter Brune } 396ad540459SPierre Jolivet if (zo + zsize > info.mz - 1) zsize -= zo + zsize - info.mz; 39790c77898SPeter Brune } 39890c77898SPeter Brune 3999566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da[idx], xsize, ysize, zsize)); 4009566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(da[idx], 1, 1, 1)); 4019566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da[idx], DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED)); 402e30e807fSPeter Brune 403e30e807fSPeter Brune /* set up as a block instead */ 4049566063dSJacob Faibussowitsch PetscCall(DMSetUp(da[idx])); 405e30e807fSPeter Brune 40640234c92SPeter Brune /* nonoverlapping region */ 4079566063dSJacob Faibussowitsch PetscCall(DMDASetNonOverlappingRegion(da[idx], xs, ys, zs, xm, ym, zm)); 40840234c92SPeter Brune 40995c13181SPeter Brune /* this alters the behavior of DMDAGetInfo, DMDAGetLocalInfo, DMDAGetCorners, and DMDAGetGhostedCorners and should be used with care */ 4109566063dSJacob Faibussowitsch PetscCall(DMDASetOffset(da[idx], xo, yo, zo, info.mx, info.my, info.mz)); 41190c77898SPeter Brune xs += xm; 41290c77898SPeter Brune idx++; 41390c77898SPeter Brune } 41490c77898SPeter Brune ys += ym; 41590c77898SPeter Brune } 41690c77898SPeter Brune zs += zm; 41790c77898SPeter Brune } 41890c77898SPeter Brune *sdm = da; 4193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 420e30e807fSPeter Brune } 421e30e807fSPeter Brune 422e30e807fSPeter Brune /* 423e30e807fSPeter Brune Fills the local vector problem on the subdomain from the global problem. 424e30e807fSPeter Brune 425285ae305SPeter Brune Right now this assumes one subdomain per processor. 426285ae305SPeter Brune 427e30e807fSPeter Brune */ 428d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm, PetscInt nsubdms, DM *subdms, VecScatter **iscat, VecScatter **oscat, VecScatter **lscat) 429d71ae5a4SJacob Faibussowitsch { 430285ae305SPeter Brune DMDALocalInfo info, subinfo; 431e30e807fSPeter Brune DM subdm; 432285ae305SPeter Brune MatStencil upper, lower; 433285ae305SPeter Brune IS idis, isis, odis, osis, gdis; 434285ae305SPeter Brune Vec svec, dvec, slvec; 43540234c92SPeter Brune PetscInt xm, ym, zm, xs, ys, zs; 43690c77898SPeter Brune PetscInt i; 4373b5e53cdSSajid Ali PetscBool patchis_offproc = PETSC_TRUE; 438e30e807fSPeter Brune 439e30e807fSPeter Brune PetscFunctionBegin; 440e30e807fSPeter Brune /* allocate the arrays of scatters */ 4419566063dSJacob Faibussowitsch if (iscat) PetscCall(PetscMalloc1(nsubdms, iscat)); 4429566063dSJacob Faibussowitsch if (oscat) PetscCall(PetscMalloc1(nsubdms, oscat)); 4439566063dSJacob Faibussowitsch if (lscat) PetscCall(PetscMalloc1(nsubdms, lscat)); 444e30e807fSPeter Brune 4459566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 44690c77898SPeter Brune for (i = 0; i < nsubdms; i++) { 44790c77898SPeter Brune subdm = subdms[i]; 4489566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(subdm, &subinfo)); 4499566063dSJacob Faibussowitsch PetscCall(DMDAGetNonOverlappingRegion(subdm, &xs, &ys, &zs, &xm, &ym, &zm)); 450e30e807fSPeter Brune 451285ae305SPeter Brune /* create the global and subdomain index sets for the inner domain */ 45240234c92SPeter Brune lower.i = xs; 45340234c92SPeter Brune lower.j = ys; 45440234c92SPeter Brune lower.k = zs; 45540234c92SPeter Brune upper.i = xs + xm; 45640234c92SPeter Brune upper.j = ys + ym; 45740234c92SPeter Brune upper.k = zs + zm; 4589566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &idis, patchis_offproc)); 4599566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(subdm, &lower, &upper, &isis, patchis_offproc)); 460e30e807fSPeter Brune 461285ae305SPeter Brune /* create the global and subdomain index sets for the outer subdomain */ 462285ae305SPeter Brune lower.i = subinfo.xs; 463285ae305SPeter Brune lower.j = subinfo.ys; 464285ae305SPeter Brune lower.k = subinfo.zs; 465285ae305SPeter Brune upper.i = subinfo.xs + subinfo.xm; 466285ae305SPeter Brune upper.j = subinfo.ys + subinfo.ym; 467285ae305SPeter Brune upper.k = subinfo.zs + subinfo.zm; 4689566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &odis, patchis_offproc)); 4699566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(subdm, &lower, &upper, &osis, patchis_offproc)); 470e30e807fSPeter Brune 471285ae305SPeter Brune /* global and subdomain ISes for the local indices of the subdomain */ 472285ae305SPeter Brune /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */ 473285ae305SPeter Brune lower.i = subinfo.gxs; 474285ae305SPeter Brune lower.j = subinfo.gys; 475285ae305SPeter Brune lower.k = subinfo.gzs; 476285ae305SPeter Brune upper.i = subinfo.gxs + subinfo.gxm; 477285ae305SPeter Brune upper.j = subinfo.gys + subinfo.gym; 478285ae305SPeter Brune upper.k = subinfo.gzs + subinfo.gzm; 4799566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &gdis, patchis_offproc)); 480e30e807fSPeter Brune 481e30e807fSPeter Brune /* form the scatter */ 4829566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &dvec)); 4839566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(subdm, &svec)); 4849566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(subdm, &slvec)); 485e30e807fSPeter Brune 4869566063dSJacob Faibussowitsch if (iscat) PetscCall(VecScatterCreate(dvec, idis, svec, isis, &(*iscat)[i])); 4879566063dSJacob Faibussowitsch if (oscat) PetscCall(VecScatterCreate(dvec, odis, svec, osis, &(*oscat)[i])); 4889566063dSJacob Faibussowitsch if (lscat) PetscCall(VecScatterCreate(dvec, gdis, slvec, NULL, &(*lscat)[i])); 489e30e807fSPeter Brune 4909566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &dvec)); 4919566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(subdm, &svec)); 4929566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(subdm, &slvec)); 493e30e807fSPeter Brune 4949566063dSJacob Faibussowitsch PetscCall(ISDestroy(&idis)); 4959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isis)); 496e30e807fSPeter Brune 4979566063dSJacob Faibussowitsch PetscCall(ISDestroy(&odis)); 4989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&osis)); 499e30e807fSPeter Brune 5009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&gdis)); 50190c77898SPeter Brune } 5023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 503e30e807fSPeter Brune } 504e30e807fSPeter Brune 50566976f2fSJacob Faibussowitsch static PetscErrorCode DMDASubDomainIS_Private(DM dm, PetscInt n, DM *subdm, IS **iis, IS **ois) 506d71ae5a4SJacob Faibussowitsch { 50790c77898SPeter Brune PetscInt i; 508e30e807fSPeter Brune DMDALocalInfo info, subinfo; 509285ae305SPeter Brune MatStencil lower, upper; 5103b5e53cdSSajid Ali PetscBool patchis_offproc = PETSC_TRUE; 511e30e807fSPeter Brune 512e30e807fSPeter Brune PetscFunctionBegin; 5139566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 5149566063dSJacob Faibussowitsch if (iis) PetscCall(PetscMalloc1(n, iis)); 5159566063dSJacob Faibussowitsch if (ois) PetscCall(PetscMalloc1(n, ois)); 516e30e807fSPeter Brune 51790c77898SPeter Brune for (i = 0; i < n; i++) { 5189566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(subdm[i], &subinfo)); 51990c77898SPeter Brune if (iis) { 520285ae305SPeter Brune /* create the inner IS */ 521285ae305SPeter Brune lower.i = info.xs; 522285ae305SPeter Brune lower.j = info.ys; 523285ae305SPeter Brune lower.k = info.zs; 524285ae305SPeter Brune upper.i = info.xs + info.xm; 525285ae305SPeter Brune upper.j = info.ys + info.ym; 526285ae305SPeter Brune upper.k = info.zs + info.zm; 5279566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &(*iis)[i], patchis_offproc)); 52890c77898SPeter Brune } 529e30e807fSPeter Brune 53090c77898SPeter Brune if (ois) { 531285ae305SPeter Brune /* create the outer IS */ 532285ae305SPeter Brune lower.i = subinfo.xs; 533285ae305SPeter Brune lower.j = subinfo.ys; 534285ae305SPeter Brune lower.k = subinfo.zs; 535285ae305SPeter Brune upper.i = subinfo.xs + subinfo.xm; 536285ae305SPeter Brune upper.j = subinfo.ys + subinfo.ym; 537285ae305SPeter Brune upper.k = subinfo.zs + subinfo.zm; 5389566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &(*ois)[i], patchis_offproc)); 53990c77898SPeter Brune } 54090c77898SPeter Brune } 5413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 542e30e807fSPeter Brune } 543e30e807fSPeter Brune 544d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateDomainDecomposition_DA(DM dm, PetscInt *len, char ***names, IS **iis, IS **ois, DM **subdm) 545d71ae5a4SJacob Faibussowitsch { 54666976f2fSJacob Faibussowitsch DM *sdm = NULL; 54766976f2fSJacob Faibussowitsch PetscInt n = 0, i; 5480a696f66SPeter Brune 5490adebc6cSBarry Smith PetscFunctionBegin; 5509566063dSJacob Faibussowitsch PetscCall(DMDASubDomainDA_Private(dm, &n, &sdm)); 55190c77898SPeter Brune if (names) { 5529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, names)); 553ea78f98cSLisandro Dalcin for (i = 0; i < n; i++) (*names)[i] = NULL; 554e30e807fSPeter Brune } 5559566063dSJacob Faibussowitsch PetscCall(DMDASubDomainIS_Private(dm, n, sdm, iis, ois)); 55690c77898SPeter Brune if (subdm) *subdm = sdm; 557e2c616c8SPeter Brune else { 55848a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(DMDestroy(&sdm[i])); 559e2c616c8SPeter Brune } 56090c77898SPeter Brune if (len) *len = n; 5613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 562e30e807fSPeter Brune } 563