1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 2e30e807fSPeter Brune 3063654ddSPeter Brune /*@ 4dce8aebaSBarry Smith DMDACreatePatchIS - Creates an index set corresponding to a patch of the `DMDA`. 5ff9846d9SPeter Brune 63b5e53cdSSajid Ali Collective 7063654ddSPeter Brune 8063654ddSPeter Brune Input Parameters: 9dce8aebaSBarry Smith + da - the `DMDA` 10063654ddSPeter Brune . lower - a matstencil with i, j and k corresponding to the lower corner of the patch 113b5e53cdSSajid Ali . upper - a matstencil with i, j and k corresponding to the upper corner of the patch 123b5e53cdSSajid Ali - offproc - indicate whether the returned IS will contain off process indices 13063654ddSPeter Brune 14*2fe279fdSBarry Smith Output Parameter: 15dce8aebaSBarry Smith . is - the `IS` corresponding to the patch 16063654ddSPeter Brune 17063654ddSPeter Brune Level: developer 18063654ddSPeter Brune 193b5e53cdSSajid Ali Notes: 20dce8aebaSBarry Smith This routine always returns an `IS` on the `DMDA` comm, if offproc is set to `PETSC_TRUE`, 21dce8aebaSBarry Smith the routine returns an `IS` with all the indices requested regardless of whether these indices 223b5e53cdSSajid Ali are present on the requesting rank or not. Thus, it is upon the caller to ensure that 23dce8aebaSBarry Smith the indices returned in this mode are appropriate. If offproc is set to `PETSC_FALSE`, 24dce8aebaSBarry Smith the `IS` only returns the subset of indices that are present on the requesting rank and there 253b5e53cdSSajid Ali is no duplication of indices. 263b5e53cdSSajid Ali 27dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMCreateDomainDecompositionScatters()` 28063654ddSPeter Brune @*/ 29d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDACreatePatchIS(DM da, MatStencil *lower, MatStencil *upper, IS *is, PetscBool offproc) 30d71ae5a4SJacob Faibussowitsch { 31063654ddSPeter Brune PetscInt ms = 0, ns = 0, ps = 0; 323b5e53cdSSajid Ali PetscInt mw = 0, nw = 0, pw = 0; 33063654ddSPeter Brune PetscInt me = 1, ne = 1, pe = 1; 34063654ddSPeter Brune PetscInt mr = 0, nr = 0, pr = 0; 35ff9846d9SPeter Brune PetscInt ii, jj, kk; 36063654ddSPeter Brune PetscInt si, sj, sk; 373b5e53cdSSajid Ali PetscInt i, j, k, l, idx = 0; 38ff9846d9SPeter Brune PetscInt base; 39063654ddSPeter Brune PetscInt xm = 1, ym = 1, zm = 1; 40ff9846d9SPeter Brune PetscInt ox, oy, oz; 41063654ddSPeter Brune PetscInt m, n, p, M, N, P, dof; 423b5e53cdSSajid Ali const PetscInt *lx, *ly, *lz; 43063654ddSPeter Brune PetscInt nindices; 44063654ddSPeter Brune PetscInt *indices; 45063654ddSPeter Brune DM_DA *dd = (DM_DA *)da->data; 463b5e53cdSSajid Ali PetscBool skip_i = PETSC_TRUE, skip_j = PETSC_TRUE, skip_k = PETSC_TRUE; 473b5e53cdSSajid Ali PetscBool valid_j = PETSC_FALSE, valid_k = PETSC_FALSE; /* DMDA has at least 1 dimension */ 48ff9846d9SPeter Brune 490adebc6cSBarry Smith PetscFunctionBegin; 509371c9d4SSatish Balay M = dd->M; 519371c9d4SSatish Balay N = dd->N; 529371c9d4SSatish Balay P = dd->P; 539371c9d4SSatish Balay m = dd->m; 549371c9d4SSatish Balay n = dd->n; 559371c9d4SSatish Balay p = dd->p; 56063654ddSPeter Brune dof = dd->w; 573b5e53cdSSajid Ali 583b5e53cdSSajid Ali nindices = -1; 593b5e53cdSSajid Ali if (PetscLikely(upper->i - lower->i)) { 603b5e53cdSSajid Ali nindices = nindices * (upper->i - lower->i); 613b5e53cdSSajid Ali skip_i = PETSC_FALSE; 623b5e53cdSSajid Ali } 633b5e53cdSSajid Ali if (N > 1) { 643b5e53cdSSajid Ali valid_j = PETSC_TRUE; 653b5e53cdSSajid Ali if (PetscLikely(upper->j - lower->j)) { 663b5e53cdSSajid Ali nindices = nindices * (upper->j - lower->j); 673b5e53cdSSajid Ali skip_j = PETSC_FALSE; 683b5e53cdSSajid Ali } 693b5e53cdSSajid Ali } 703b5e53cdSSajid Ali if (P > 1) { 713b5e53cdSSajid Ali valid_k = PETSC_TRUE; 723b5e53cdSSajid Ali if (PetscLikely(upper->k - lower->k)) { 733b5e53cdSSajid Ali nindices = nindices * (upper->k - lower->k); 743b5e53cdSSajid Ali skip_k = PETSC_FALSE; 753b5e53cdSSajid Ali } 763b5e53cdSSajid Ali } 773b5e53cdSSajid Ali if (PetscLikely(nindices < 0)) { 783b5e53cdSSajid Ali if (PetscUnlikely(skip_i && skip_j && skip_k)) { 793b5e53cdSSajid Ali nindices = 0; 803b5e53cdSSajid Ali } else nindices = nindices * (-1); 813b5e53cdSSajid Ali } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Lower and Upper stencils are identical! Please check inputs."); 823b5e53cdSSajid Ali 839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nindices * dof, &indices)); 849566063dSJacob Faibussowitsch PetscCall(DMDAGetOffset(da, &ox, &oy, &oz, NULL, NULL, NULL)); 853b5e53cdSSajid Ali 863b5e53cdSSajid Ali if (!valid_k) { 873b5e53cdSSajid Ali k = 0; 883b5e53cdSSajid Ali upper->k = 0; 893b5e53cdSSajid Ali lower->k = 0; 903b5e53cdSSajid Ali } 913b5e53cdSSajid Ali if (!valid_j) { 923b5e53cdSSajid Ali j = 0; 933b5e53cdSSajid Ali upper->j = 0; 943b5e53cdSSajid Ali lower->j = 0; 953b5e53cdSSajid Ali } 963b5e53cdSSajid Ali 973b5e53cdSSajid Ali if (offproc) { 989566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, &lz)); 99063654ddSPeter Brune /* start at index 0 on processor 0 */ 100063654ddSPeter Brune mr = 0; 101063654ddSPeter Brune nr = 0; 102063654ddSPeter Brune pr = 0; 103063654ddSPeter Brune ms = 0; 104063654ddSPeter Brune ns = 0; 105063654ddSPeter Brune ps = 0; 106063654ddSPeter Brune if (lx) me = lx[0]; 107063654ddSPeter Brune if (ly) ne = ly[0]; 108063654ddSPeter Brune if (lz) pe = lz[0]; 1093b5e53cdSSajid Ali /* 1103b5e53cdSSajid Ali If no indices are to be returned, create an empty is, 1113b5e53cdSSajid Ali this prevents hanging in while loops 1123b5e53cdSSajid Ali */ 1133b5e53cdSSajid Ali if (skip_i && skip_j && skip_k) goto createis; 1143b5e53cdSSajid Ali /* 1153b5e53cdSSajid Ali do..while loops to ensure the block gets entered once, 1163b5e53cdSSajid Ali regardless of control condition being met, necessary for 1173b5e53cdSSajid Ali cases when a subset of skip_i/j/k is true 1183b5e53cdSSajid Ali */ 1199371c9d4SSatish Balay if (skip_k) k = upper->k - oz; 1209371c9d4SSatish Balay else k = lower->k - oz; 1213b5e53cdSSajid Ali do { 1229371c9d4SSatish Balay if (skip_j) j = upper->j - oy; 1239371c9d4SSatish Balay else j = lower->j - oy; 1243b5e53cdSSajid Ali do { 1259371c9d4SSatish Balay if (skip_i) i = upper->i - ox; 1269371c9d4SSatish Balay else i = lower->i - ox; 1273b5e53cdSSajid Ali do { 128063654ddSPeter Brune /* "actual" indices rather than ones outside of the domain */ 129063654ddSPeter Brune ii = i; 130063654ddSPeter Brune jj = j; 131063654ddSPeter Brune kk = k; 132063654ddSPeter Brune if (ii < 0) ii = M + ii; 133063654ddSPeter Brune if (jj < 0) jj = N + jj; 134063654ddSPeter Brune if (kk < 0) kk = P + kk; 135063654ddSPeter Brune if (ii > M - 1) ii = ii - M; 136063654ddSPeter Brune if (jj > N - 1) jj = jj - N; 137063654ddSPeter Brune if (kk > P - 1) kk = kk - P; 138063654ddSPeter Brune /* gone out of processor range on x axis */ 139063654ddSPeter Brune while (ii > me - 1 || ii < ms) { 140063654ddSPeter Brune if (mr == m - 1) { 141063654ddSPeter Brune ms = 0; 142063654ddSPeter Brune me = lx[0]; 143063654ddSPeter Brune mr = 0; 144063654ddSPeter Brune } else { 145063654ddSPeter Brune mr++; 146063654ddSPeter Brune ms = me; 147063654ddSPeter Brune me += lx[mr]; 148063654ddSPeter Brune } 149063654ddSPeter Brune } 150063654ddSPeter Brune /* gone out of processor range on y axis */ 151063654ddSPeter Brune while (jj > ne - 1 || jj < ns) { 152063654ddSPeter Brune if (nr == n - 1) { 153063654ddSPeter Brune ns = 0; 154063654ddSPeter Brune ne = ly[0]; 155063654ddSPeter Brune nr = 0; 156063654ddSPeter Brune } else { 157063654ddSPeter Brune nr++; 158063654ddSPeter Brune ns = ne; 159063654ddSPeter Brune ne += ly[nr]; 160063654ddSPeter Brune } 161063654ddSPeter Brune } 162063654ddSPeter Brune /* gone out of processor range on z axis */ 163063654ddSPeter Brune while (kk > pe - 1 || kk < ps) { 164063654ddSPeter Brune if (pr == p - 1) { 165063654ddSPeter Brune ps = 0; 166063654ddSPeter Brune pe = lz[0]; 167063654ddSPeter Brune pr = 0; 168063654ddSPeter Brune } else { 169063654ddSPeter Brune pr++; 170063654ddSPeter Brune ps = pe; 171063654ddSPeter Brune pe += lz[pr]; 172063654ddSPeter Brune } 173063654ddSPeter Brune } 174063654ddSPeter Brune /* compute the vector base on owning processor */ 175063654ddSPeter Brune xm = me - ms; 176063654ddSPeter Brune ym = ne - ns; 177063654ddSPeter Brune zm = pe - ps; 178b0bff951SPeter Brune base = ms * ym * zm + ns * M * zm + ps * M * N; 179063654ddSPeter Brune /* compute the local coordinates on owning processor */ 180063654ddSPeter Brune si = ii - ms; 181063654ddSPeter Brune sj = jj - ns; 182063654ddSPeter Brune sk = kk - ps; 183063654ddSPeter Brune for (l = 0; l < dof; l++) { 184063654ddSPeter Brune indices[idx] = l + dof * (base + si + xm * sj + xm * ym * sk); 185ff9846d9SPeter Brune idx++; 186ff9846d9SPeter Brune } 1873b5e53cdSSajid Ali i++; 1883b5e53cdSSajid Ali } while (i < upper->i - ox); 1893b5e53cdSSajid Ali j++; 1903b5e53cdSSajid Ali } while (j < upper->j - oy); 1913b5e53cdSSajid Ali k++; 1923b5e53cdSSajid Ali } while (k < upper->k - oz); 1933b5e53cdSSajid Ali } 1943b5e53cdSSajid Ali 1953b5e53cdSSajid Ali if (!offproc) { 1969566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &ms, &ns, &ps, &mw, &nw, &pw)); 1973b5e53cdSSajid Ali me = ms + mw; 1983b5e53cdSSajid Ali if (N > 1) ne = ns + nw; 1993b5e53cdSSajid Ali if (P > 1) pe = ps + pw; 2003b5e53cdSSajid Ali /* Account for DM offsets */ 2019371c9d4SSatish Balay ms = ms - ox; 2029371c9d4SSatish Balay me = me - ox; 2039371c9d4SSatish Balay ns = ns - oy; 2049371c9d4SSatish Balay ne = ne - oy; 2059371c9d4SSatish Balay ps = ps - oz; 2069371c9d4SSatish Balay pe = pe - oz; 2073b5e53cdSSajid Ali 2083b5e53cdSSajid Ali /* compute the vector base on owning processor */ 2093b5e53cdSSajid Ali xm = me - ms; 2103b5e53cdSSajid Ali ym = ne - ns; 2113b5e53cdSSajid Ali zm = pe - ps; 2123b5e53cdSSajid Ali base = ms * ym * zm + ns * M * zm + ps * M * N; 2133b5e53cdSSajid Ali /* 2143b5e53cdSSajid Ali if no indices are to be returned, create an empty is, 2153b5e53cdSSajid Ali this prevents hanging in while loops 2163b5e53cdSSajid Ali */ 2173b5e53cdSSajid Ali if (skip_i && skip_j && skip_k) goto createis; 2183b5e53cdSSajid Ali /* 2193b5e53cdSSajid Ali do..while loops to ensure the block gets entered once, 2203b5e53cdSSajid Ali regardless of control condition being met, necessary for 2213b5e53cdSSajid Ali cases when a subset of skip_i/j/k is true 2223b5e53cdSSajid Ali */ 2239371c9d4SSatish Balay if (skip_k) k = upper->k - oz; 2249371c9d4SSatish Balay else k = lower->k - oz; 2253b5e53cdSSajid Ali do { 2269371c9d4SSatish Balay if (skip_j) j = upper->j - oy; 2279371c9d4SSatish Balay else j = lower->j - oy; 2283b5e53cdSSajid Ali do { 2299371c9d4SSatish Balay if (skip_i) i = upper->i - ox; 2309371c9d4SSatish Balay else i = lower->i - ox; 2313b5e53cdSSajid Ali do { 2323b5e53cdSSajid Ali if (k >= ps && k <= pe - 1) { 2333b5e53cdSSajid Ali if (j >= ns && j <= ne - 1) { 2343b5e53cdSSajid Ali if (i >= ms && i <= me - 1) { 2353b5e53cdSSajid Ali /* compute the local coordinates on owning processor */ 2363b5e53cdSSajid Ali si = i - ms; 2373b5e53cdSSajid Ali sj = j - ns; 2383b5e53cdSSajid Ali sk = k - ps; 2393b5e53cdSSajid Ali for (l = 0; l < dof; l++) { 2403b5e53cdSSajid Ali indices[idx] = l + dof * (base + si + xm * sj + xm * ym * sk); 2413b5e53cdSSajid Ali idx++; 242ff9846d9SPeter Brune } 243ff9846d9SPeter Brune } 244063654ddSPeter Brune } 2453b5e53cdSSajid Ali } 2463b5e53cdSSajid Ali i++; 2473b5e53cdSSajid Ali } while (i < upper->i - ox); 2483b5e53cdSSajid Ali j++; 2493b5e53cdSSajid Ali } while (j < upper->j - oy); 2503b5e53cdSSajid Ali k++; 2513b5e53cdSSajid Ali } while (k < upper->k - oz); 2523b5e53cdSSajid Ali 2539566063dSJacob Faibussowitsch PetscCall(PetscRealloc((size_t)(idx * sizeof(PetscInt)), (void *)&indices)); 2543b5e53cdSSajid Ali } 2553b5e53cdSSajid Ali 2563b5e53cdSSajid Ali createis: 2579566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)da), idx, indices, PETSC_OWN_POINTER, is)); 2583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 259ff9846d9SPeter Brune } 260ff9846d9SPeter Brune 261d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASubDomainDA_Private(DM dm, PetscInt *nlocal, DM **sdm) 262d71ae5a4SJacob Faibussowitsch { 26390c77898SPeter Brune DM *da; 26490c77898SPeter Brune PetscInt dim, size, i, j, k, idx; 265e30e807fSPeter Brune DMDALocalInfo info; 266ff9846d9SPeter Brune PetscInt xsize, ysize, zsize; 267e30e807fSPeter Brune PetscInt xo, yo, zo; 26890c77898SPeter Brune PetscInt xs, ys, zs; 26990c77898SPeter Brune PetscInt xm = 1, ym = 1, zm = 1; 2707ddda789SPeter Brune PetscInt xol, yol, zol; 27190c77898SPeter Brune PetscInt m = 1, n = 1, p = 1; 27290c77898SPeter Brune PetscInt M, N, P; 27390c77898SPeter Brune PetscInt pm, mtmp; 274e30e807fSPeter Brune 275e30e807fSPeter Brune PetscFunctionBegin; 2769566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 2779566063dSJacob Faibussowitsch PetscCall(DMDAGetOverlap(dm, &xol, &yol, &zol)); 2789566063dSJacob Faibussowitsch PetscCall(DMDAGetNumLocalSubDomains(dm, &size)); 2799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &da)); 2807ddda789SPeter Brune 28190c77898SPeter Brune if (nlocal) *nlocal = size; 282e30e807fSPeter Brune 28390c77898SPeter Brune dim = info.dim; 284e30e807fSPeter Brune 28590c77898SPeter Brune M = info.xm; 28690c77898SPeter Brune N = info.ym; 28790c77898SPeter Brune P = info.zm; 28890c77898SPeter Brune 28990c77898SPeter Brune if (dim == 1) { 29090c77898SPeter Brune m = size; 29190c77898SPeter Brune } else if (dim == 2) { 29290c77898SPeter Brune m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)N))); 29390c77898SPeter Brune while (m > 0) { 29490c77898SPeter Brune n = size / m; 29590c77898SPeter Brune if (m * n * p == size) break; 29690c77898SPeter Brune m--; 29790c77898SPeter Brune } 29890c77898SPeter Brune } else if (dim == 3) { 2999371c9d4SSatish Balay n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N * N) * ((PetscReal)size) / ((PetscReal)P * M), (PetscReal)(1. / 3.))); 3009371c9d4SSatish Balay if (!n) n = 1; 30190c77898SPeter Brune while (n > 0) { 30290c77898SPeter Brune pm = size / n; 30390c77898SPeter Brune if (n * pm == size) break; 30490c77898SPeter Brune n--; 30590c77898SPeter Brune } 30690c77898SPeter Brune if (!n) n = 1; 30790c77898SPeter Brune m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)P * n))); 30890c77898SPeter Brune if (!m) m = 1; 30990c77898SPeter Brune while (m > 0) { 31090c77898SPeter Brune p = size / (m * n); 31190c77898SPeter Brune if (m * n * p == size) break; 31290c77898SPeter Brune m--; 31390c77898SPeter Brune } 3149371c9d4SSatish Balay if (M > P && m < p) { 3159371c9d4SSatish Balay mtmp = m; 3169371c9d4SSatish Balay m = p; 3179371c9d4SSatish Balay p = mtmp; 3189371c9d4SSatish Balay } 31990c77898SPeter Brune } 32090c77898SPeter Brune 32190c77898SPeter Brune zs = info.zs; 32290c77898SPeter Brune idx = 0; 32390c77898SPeter Brune for (k = 0; k < p; k++) { 32490c77898SPeter Brune ys = info.ys; 32590c77898SPeter Brune for (j = 0; j < n; j++) { 32690c77898SPeter Brune xs = info.xs; 32790c77898SPeter Brune for (i = 0; i < m; i++) { 32890c77898SPeter Brune if (dim == 1) { 32990c77898SPeter Brune xm = M / m + ((M % m) > i); 33090c77898SPeter Brune } else if (dim == 2) { 33190c77898SPeter Brune xm = M / m + ((M % m) > i); 33290c77898SPeter Brune ym = N / n + ((N % n) > j); 33390c77898SPeter Brune } else if (dim == 3) { 33490c77898SPeter Brune xm = M / m + ((M % m) > i); 33590c77898SPeter Brune ym = N / n + ((N % n) > j); 33690c77898SPeter Brune zm = P / p + ((P % p) > k); 33790c77898SPeter Brune } 33890c77898SPeter Brune 33990c77898SPeter Brune xsize = xm; 34090c77898SPeter Brune ysize = ym; 34190c77898SPeter Brune zsize = zm; 34290c77898SPeter Brune xo = xs; 34390c77898SPeter Brune yo = ys; 34490c77898SPeter Brune zo = zs; 34590c77898SPeter Brune 3469566063dSJacob Faibussowitsch PetscCall(DMDACreate(PETSC_COMM_SELF, &(da[idx]))); 3479566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(da[idx], "sub_")); 3489566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da[idx], info.dim)); 3499566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da[idx], info.dof)); 35090c77898SPeter Brune 3519566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da[idx], info.st)); 3529566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da[idx], info.sw)); 35390c77898SPeter Brune 354bff4a2f0SMatthew G. Knepley if (info.bx == DM_BOUNDARY_PERIODIC || (xs != 0)) { 3557ddda789SPeter Brune xsize += xol; 3567ddda789SPeter Brune xo -= xol; 357e30e807fSPeter Brune } 358bff4a2f0SMatthew G. Knepley if (info.by == DM_BOUNDARY_PERIODIC || (ys != 0)) { 3597ddda789SPeter Brune ysize += yol; 3607ddda789SPeter Brune yo -= yol; 361e30e807fSPeter Brune } 362bff4a2f0SMatthew G. Knepley if (info.bz == DM_BOUNDARY_PERIODIC || (zs != 0)) { 3637ddda789SPeter Brune zsize += zol; 3647ddda789SPeter Brune zo -= zol; 365e30e807fSPeter Brune } 366e30e807fSPeter Brune 367bff4a2f0SMatthew G. Knepley if (info.bx == DM_BOUNDARY_PERIODIC || (xs + xm != info.mx)) xsize += xol; 368bff4a2f0SMatthew G. Knepley if (info.by == DM_BOUNDARY_PERIODIC || (ys + ym != info.my)) ysize += yol; 369bff4a2f0SMatthew G. Knepley if (info.bz == DM_BOUNDARY_PERIODIC || (zs + zm != info.mz)) zsize += zol; 3708865f1eaSKarl Rupp 371bff4a2f0SMatthew G. Knepley if (info.bx != DM_BOUNDARY_PERIODIC) { 37290c77898SPeter Brune if (xo < 0) { 37390c77898SPeter Brune xsize += xo; 37490c77898SPeter Brune xo = 0; 37590c77898SPeter Brune } 376ad540459SPierre Jolivet if (xo + xsize > info.mx - 1) xsize -= xo + xsize - info.mx; 37790c77898SPeter Brune } 378bff4a2f0SMatthew G. Knepley if (info.by != DM_BOUNDARY_PERIODIC) { 37990c77898SPeter Brune if (yo < 0) { 38090c77898SPeter Brune ysize += yo; 38190c77898SPeter Brune yo = 0; 38290c77898SPeter Brune } 383ad540459SPierre Jolivet if (yo + ysize > info.my - 1) ysize -= yo + ysize - info.my; 38490c77898SPeter Brune } 385bff4a2f0SMatthew G. Knepley if (info.bz != DM_BOUNDARY_PERIODIC) { 38690c77898SPeter Brune if (zo < 0) { 38790c77898SPeter Brune zsize += zo; 38890c77898SPeter Brune zo = 0; 38990c77898SPeter Brune } 390ad540459SPierre Jolivet if (zo + zsize > info.mz - 1) zsize -= zo + zsize - info.mz; 39190c77898SPeter Brune } 39290c77898SPeter Brune 3939566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da[idx], xsize, ysize, zsize)); 3949566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(da[idx], 1, 1, 1)); 3959566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da[idx], DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED)); 396e30e807fSPeter Brune 397e30e807fSPeter Brune /* set up as a block instead */ 3989566063dSJacob Faibussowitsch PetscCall(DMSetUp(da[idx])); 399e30e807fSPeter Brune 40040234c92SPeter Brune /* nonoverlapping region */ 4019566063dSJacob Faibussowitsch PetscCall(DMDASetNonOverlappingRegion(da[idx], xs, ys, zs, xm, ym, zm)); 40240234c92SPeter Brune 40395c13181SPeter Brune /* this alters the behavior of DMDAGetInfo, DMDAGetLocalInfo, DMDAGetCorners, and DMDAGetGhostedCorners and should be used with care */ 4049566063dSJacob Faibussowitsch PetscCall(DMDASetOffset(da[idx], xo, yo, zo, info.mx, info.my, info.mz)); 40590c77898SPeter Brune xs += xm; 40690c77898SPeter Brune idx++; 40790c77898SPeter Brune } 40890c77898SPeter Brune ys += ym; 40990c77898SPeter Brune } 41090c77898SPeter Brune zs += zm; 41190c77898SPeter Brune } 41290c77898SPeter Brune *sdm = da; 4133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 414e30e807fSPeter Brune } 415e30e807fSPeter Brune 416e30e807fSPeter Brune /* 417e30e807fSPeter Brune Fills the local vector problem on the subdomain from the global problem. 418e30e807fSPeter Brune 419285ae305SPeter Brune Right now this assumes one subdomain per processor. 420285ae305SPeter Brune 421e30e807fSPeter Brune */ 422d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm, PetscInt nsubdms, DM *subdms, VecScatter **iscat, VecScatter **oscat, VecScatter **lscat) 423d71ae5a4SJacob Faibussowitsch { 424285ae305SPeter Brune DMDALocalInfo info, subinfo; 425e30e807fSPeter Brune DM subdm; 426285ae305SPeter Brune MatStencil upper, lower; 427285ae305SPeter Brune IS idis, isis, odis, osis, gdis; 428285ae305SPeter Brune Vec svec, dvec, slvec; 42940234c92SPeter Brune PetscInt xm, ym, zm, xs, ys, zs; 43090c77898SPeter Brune PetscInt i; 4313b5e53cdSSajid Ali PetscBool patchis_offproc = PETSC_TRUE; 432e30e807fSPeter Brune 433e30e807fSPeter Brune PetscFunctionBegin; 434e30e807fSPeter Brune /* allocate the arrays of scatters */ 4359566063dSJacob Faibussowitsch if (iscat) PetscCall(PetscMalloc1(nsubdms, iscat)); 4369566063dSJacob Faibussowitsch if (oscat) PetscCall(PetscMalloc1(nsubdms, oscat)); 4379566063dSJacob Faibussowitsch if (lscat) PetscCall(PetscMalloc1(nsubdms, lscat)); 438e30e807fSPeter Brune 4399566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 44090c77898SPeter Brune for (i = 0; i < nsubdms; i++) { 44190c77898SPeter Brune subdm = subdms[i]; 4429566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(subdm, &subinfo)); 4439566063dSJacob Faibussowitsch PetscCall(DMDAGetNonOverlappingRegion(subdm, &xs, &ys, &zs, &xm, &ym, &zm)); 444e30e807fSPeter Brune 445285ae305SPeter Brune /* create the global and subdomain index sets for the inner domain */ 44640234c92SPeter Brune lower.i = xs; 44740234c92SPeter Brune lower.j = ys; 44840234c92SPeter Brune lower.k = zs; 44940234c92SPeter Brune upper.i = xs + xm; 45040234c92SPeter Brune upper.j = ys + ym; 45140234c92SPeter Brune upper.k = zs + zm; 4529566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &idis, patchis_offproc)); 4539566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(subdm, &lower, &upper, &isis, patchis_offproc)); 454e30e807fSPeter Brune 455285ae305SPeter Brune /* create the global and subdomain index sets for the outer subdomain */ 456285ae305SPeter Brune lower.i = subinfo.xs; 457285ae305SPeter Brune lower.j = subinfo.ys; 458285ae305SPeter Brune lower.k = subinfo.zs; 459285ae305SPeter Brune upper.i = subinfo.xs + subinfo.xm; 460285ae305SPeter Brune upper.j = subinfo.ys + subinfo.ym; 461285ae305SPeter Brune upper.k = subinfo.zs + subinfo.zm; 4629566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &odis, patchis_offproc)); 4639566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(subdm, &lower, &upper, &osis, patchis_offproc)); 464e30e807fSPeter Brune 465285ae305SPeter Brune /* global and subdomain ISes for the local indices of the subdomain */ 466285ae305SPeter Brune /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */ 467285ae305SPeter Brune lower.i = subinfo.gxs; 468285ae305SPeter Brune lower.j = subinfo.gys; 469285ae305SPeter Brune lower.k = subinfo.gzs; 470285ae305SPeter Brune upper.i = subinfo.gxs + subinfo.gxm; 471285ae305SPeter Brune upper.j = subinfo.gys + subinfo.gym; 472285ae305SPeter Brune upper.k = subinfo.gzs + subinfo.gzm; 4739566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &gdis, patchis_offproc)); 474e30e807fSPeter Brune 475e30e807fSPeter Brune /* form the scatter */ 4769566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &dvec)); 4779566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(subdm, &svec)); 4789566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(subdm, &slvec)); 479e30e807fSPeter Brune 4809566063dSJacob Faibussowitsch if (iscat) PetscCall(VecScatterCreate(dvec, idis, svec, isis, &(*iscat)[i])); 4819566063dSJacob Faibussowitsch if (oscat) PetscCall(VecScatterCreate(dvec, odis, svec, osis, &(*oscat)[i])); 4829566063dSJacob Faibussowitsch if (lscat) PetscCall(VecScatterCreate(dvec, gdis, slvec, NULL, &(*lscat)[i])); 483e30e807fSPeter Brune 4849566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &dvec)); 4859566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(subdm, &svec)); 4869566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(subdm, &slvec)); 487e30e807fSPeter Brune 4889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&idis)); 4899566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isis)); 490e30e807fSPeter Brune 4919566063dSJacob Faibussowitsch PetscCall(ISDestroy(&odis)); 4929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&osis)); 493e30e807fSPeter Brune 4949566063dSJacob Faibussowitsch PetscCall(ISDestroy(&gdis)); 49590c77898SPeter Brune } 4963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 497e30e807fSPeter Brune } 498e30e807fSPeter Brune 499d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASubDomainIS_Private(DM dm, PetscInt n, DM *subdm, IS **iis, IS **ois) 500d71ae5a4SJacob Faibussowitsch { 50190c77898SPeter Brune PetscInt i; 502e30e807fSPeter Brune DMDALocalInfo info, subinfo; 503285ae305SPeter Brune MatStencil lower, upper; 5043b5e53cdSSajid Ali PetscBool patchis_offproc = PETSC_TRUE; 505e30e807fSPeter Brune 506e30e807fSPeter Brune PetscFunctionBegin; 5079566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 5089566063dSJacob Faibussowitsch if (iis) PetscCall(PetscMalloc1(n, iis)); 5099566063dSJacob Faibussowitsch if (ois) PetscCall(PetscMalloc1(n, ois)); 510e30e807fSPeter Brune 51190c77898SPeter Brune for (i = 0; i < n; i++) { 5129566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(subdm[i], &subinfo)); 51390c77898SPeter Brune if (iis) { 514285ae305SPeter Brune /* create the inner IS */ 515285ae305SPeter Brune lower.i = info.xs; 516285ae305SPeter Brune lower.j = info.ys; 517285ae305SPeter Brune lower.k = info.zs; 518285ae305SPeter Brune upper.i = info.xs + info.xm; 519285ae305SPeter Brune upper.j = info.ys + info.ym; 520285ae305SPeter Brune upper.k = info.zs + info.zm; 5219566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &(*iis)[i], patchis_offproc)); 52290c77898SPeter Brune } 523e30e807fSPeter Brune 52490c77898SPeter Brune if (ois) { 525285ae305SPeter Brune /* create the outer IS */ 526285ae305SPeter Brune lower.i = subinfo.xs; 527285ae305SPeter Brune lower.j = subinfo.ys; 528285ae305SPeter Brune lower.k = subinfo.zs; 529285ae305SPeter Brune upper.i = subinfo.xs + subinfo.xm; 530285ae305SPeter Brune upper.j = subinfo.ys + subinfo.ym; 531285ae305SPeter Brune upper.k = subinfo.zs + subinfo.zm; 5329566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &(*ois)[i], patchis_offproc)); 53390c77898SPeter Brune } 53490c77898SPeter Brune } 5353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 536e30e807fSPeter Brune } 537e30e807fSPeter Brune 538d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateDomainDecomposition_DA(DM dm, PetscInt *len, char ***names, IS **iis, IS **ois, DM **subdm) 539d71ae5a4SJacob Faibussowitsch { 54090c77898SPeter Brune DM *sdm; 54190c77898SPeter Brune PetscInt n, i; 5420a696f66SPeter Brune 5430adebc6cSBarry Smith PetscFunctionBegin; 5449566063dSJacob Faibussowitsch PetscCall(DMDASubDomainDA_Private(dm, &n, &sdm)); 54590c77898SPeter Brune if (names) { 5469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, names)); 547ea78f98cSLisandro Dalcin for (i = 0; i < n; i++) (*names)[i] = NULL; 548e30e807fSPeter Brune } 5499566063dSJacob Faibussowitsch PetscCall(DMDASubDomainIS_Private(dm, n, sdm, iis, ois)); 55090c77898SPeter Brune if (subdm) *subdm = sdm; 551e2c616c8SPeter Brune else { 55248a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(DMDestroy(&sdm[i])); 553e2c616c8SPeter Brune } 55490c77898SPeter Brune if (len) *len = n; 5553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 556e30e807fSPeter Brune } 557