1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 2e30e807fSPeter Brune 3063654ddSPeter Brune /*@ 4*12b4a537SBarry 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` 10*12b4a537SBarry Smith . lower - a `MatStencil` with i, j and k entries corresponding to the lower corner of the patch 11*12b4a537SBarry Smith . upper - a `MatStencil` with i, j and k entries corresponding to the upper corner of the patch 12*12b4a537SBarry 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: 20*12b4a537SBarry Smith This routine always returns an `IS` on the `DMDA` communicator. 213b5e53cdSSajid Ali 22*12b4a537SBarry Smith If `offproc` is set to `PETSC_TRUE`, 23*12b4a537SBarry Smith the routine returns an `IS` with all the indices requested regardless of whether these indices 24*12b4a537SBarry Smith are present on the requesting MPI process or not. Thus, it is upon the caller to ensure that 25*12b4a537SBarry Smith the indices returned in this mode are appropriate. 26*12b4a537SBarry Smith 27*12b4a537SBarry Smith If `offproc` is set to `PETSC_FALSE`, 28*12b4a537SBarry Smith the `IS` only returns the subset of indices that are present on the requesting MPI process and there 29*12b4a537SBarry Smith is no duplication of indices between multiple MPI processes. 30*12b4a537SBarry Smith 31*12b4a537SBarry 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; 413b5e53cdSSajid Ali PetscInt i, j, k, 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 k = 0; 923b5e53cdSSajid Ali upper->k = 0; 933b5e53cdSSajid Ali lower->k = 0; 943b5e53cdSSajid Ali } 953b5e53cdSSajid Ali if (!valid_j) { 963b5e53cdSSajid Ali j = 0; 973b5e53cdSSajid Ali upper->j = 0; 983b5e53cdSSajid Ali lower->j = 0; 993b5e53cdSSajid Ali } 1003b5e53cdSSajid Ali 1013b5e53cdSSajid Ali if (offproc) { 1029566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, &lz)); 103063654ddSPeter Brune /* start at index 0 on processor 0 */ 104063654ddSPeter Brune mr = 0; 105063654ddSPeter Brune nr = 0; 106063654ddSPeter Brune pr = 0; 107063654ddSPeter Brune ms = 0; 108063654ddSPeter Brune ns = 0; 109063654ddSPeter Brune ps = 0; 110063654ddSPeter Brune if (lx) me = lx[0]; 111063654ddSPeter Brune if (ly) ne = ly[0]; 112063654ddSPeter Brune if (lz) pe = lz[0]; 1133b5e53cdSSajid Ali /* 1143b5e53cdSSajid Ali If no indices are to be returned, create an empty is, 1153b5e53cdSSajid Ali this prevents hanging in while loops 1163b5e53cdSSajid Ali */ 1173b5e53cdSSajid Ali if (skip_i && skip_j && skip_k) goto createis; 1183b5e53cdSSajid Ali /* 1193b5e53cdSSajid Ali do..while loops to ensure the block gets entered once, 1203b5e53cdSSajid Ali regardless of control condition being met, necessary for 1213b5e53cdSSajid Ali cases when a subset of skip_i/j/k is true 1223b5e53cdSSajid Ali */ 1239371c9d4SSatish Balay if (skip_k) k = upper->k - oz; 1249371c9d4SSatish Balay else k = lower->k - oz; 1253b5e53cdSSajid Ali do { 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 { 2309371c9d4SSatish Balay if (skip_j) j = upper->j - oy; 2319371c9d4SSatish Balay else j = lower->j - oy; 2323b5e53cdSSajid Ali do { 2339371c9d4SSatish Balay if (skip_i) i = upper->i - ox; 2349371c9d4SSatish Balay else i = lower->i - ox; 2353b5e53cdSSajid Ali do { 2363b5e53cdSSajid Ali if (k >= ps && k <= pe - 1) { 2373b5e53cdSSajid Ali if (j >= ns && j <= ne - 1) { 2383b5e53cdSSajid Ali if (i >= ms && i <= me - 1) { 2393b5e53cdSSajid Ali /* compute the local coordinates on owning processor */ 2403b5e53cdSSajid Ali si = i - ms; 2413b5e53cdSSajid Ali sj = j - ns; 2423b5e53cdSSajid Ali sk = k - ps; 2433b5e53cdSSajid Ali for (l = 0; l < dof; l++) { 2443b5e53cdSSajid Ali indices[idx] = l + dof * (base + si + xm * sj + xm * ym * sk); 2453b5e53cdSSajid Ali idx++; 246ff9846d9SPeter Brune } 247ff9846d9SPeter Brune } 248063654ddSPeter Brune } 2493b5e53cdSSajid Ali } 2503b5e53cdSSajid Ali i++; 2513b5e53cdSSajid Ali } while (i < upper->i - ox); 2523b5e53cdSSajid Ali j++; 2533b5e53cdSSajid Ali } while (j < upper->j - oy); 2543b5e53cdSSajid Ali k++; 2553b5e53cdSSajid Ali } while (k < upper->k - oz); 2563b5e53cdSSajid Ali 2579566063dSJacob Faibussowitsch PetscCall(PetscRealloc((size_t)(idx * sizeof(PetscInt)), (void *)&indices)); 2583b5e53cdSSajid Ali } 2593b5e53cdSSajid Ali 2603b5e53cdSSajid Ali createis: 2619566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)da), idx, indices, PETSC_OWN_POINTER, is)); 2623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 263ff9846d9SPeter Brune } 264ff9846d9SPeter Brune 26566976f2fSJacob Faibussowitsch static PetscErrorCode DMDASubDomainDA_Private(DM dm, PetscInt *nlocal, DM **sdm) 266d71ae5a4SJacob Faibussowitsch { 26790c77898SPeter Brune DM *da; 26890c77898SPeter Brune PetscInt dim, size, i, j, k, idx; 269e30e807fSPeter Brune DMDALocalInfo info; 270ff9846d9SPeter Brune PetscInt xsize, ysize, zsize; 271e30e807fSPeter Brune PetscInt xo, yo, zo; 27290c77898SPeter Brune PetscInt xs, ys, zs; 27390c77898SPeter Brune PetscInt xm = 1, ym = 1, zm = 1; 2747ddda789SPeter Brune PetscInt xol, yol, zol; 27590c77898SPeter Brune PetscInt m = 1, n = 1, p = 1; 27690c77898SPeter Brune PetscInt M, N, P; 27790c77898SPeter Brune PetscInt pm, mtmp; 278e30e807fSPeter Brune 279e30e807fSPeter Brune PetscFunctionBegin; 2809566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 2819566063dSJacob Faibussowitsch PetscCall(DMDAGetOverlap(dm, &xol, &yol, &zol)); 2829566063dSJacob Faibussowitsch PetscCall(DMDAGetNumLocalSubDomains(dm, &size)); 2839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &da)); 2847ddda789SPeter Brune 28590c77898SPeter Brune if (nlocal) *nlocal = size; 286e30e807fSPeter Brune 28790c77898SPeter Brune dim = info.dim; 288e30e807fSPeter Brune 28990c77898SPeter Brune M = info.xm; 29090c77898SPeter Brune N = info.ym; 29190c77898SPeter Brune P = info.zm; 29290c77898SPeter Brune 29390c77898SPeter Brune if (dim == 1) { 29490c77898SPeter Brune m = size; 29590c77898SPeter Brune } else if (dim == 2) { 29690c77898SPeter Brune m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)N))); 29790c77898SPeter Brune while (m > 0) { 29890c77898SPeter Brune n = size / m; 29990c77898SPeter Brune if (m * n * p == size) break; 30090c77898SPeter Brune m--; 30190c77898SPeter Brune } 30290c77898SPeter Brune } else if (dim == 3) { 3039371c9d4SSatish Balay n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N * N) * ((PetscReal)size) / ((PetscReal)P * M), (PetscReal)(1. / 3.))); 3049371c9d4SSatish Balay if (!n) n = 1; 30590c77898SPeter Brune while (n > 0) { 30690c77898SPeter Brune pm = size / n; 30790c77898SPeter Brune if (n * pm == size) break; 30890c77898SPeter Brune n--; 30990c77898SPeter Brune } 31090c77898SPeter Brune if (!n) n = 1; 31190c77898SPeter Brune m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)P * n))); 31290c77898SPeter Brune if (!m) m = 1; 31390c77898SPeter Brune while (m > 0) { 31490c77898SPeter Brune p = size / (m * n); 31590c77898SPeter Brune if (m * n * p == size) break; 31690c77898SPeter Brune m--; 31790c77898SPeter Brune } 3189371c9d4SSatish Balay if (M > P && m < p) { 3199371c9d4SSatish Balay mtmp = m; 3209371c9d4SSatish Balay m = p; 3219371c9d4SSatish Balay p = mtmp; 3229371c9d4SSatish Balay } 32390c77898SPeter Brune } 32490c77898SPeter Brune 32590c77898SPeter Brune zs = info.zs; 32690c77898SPeter Brune idx = 0; 32790c77898SPeter Brune for (k = 0; k < p; k++) { 32890c77898SPeter Brune ys = info.ys; 32990c77898SPeter Brune for (j = 0; j < n; j++) { 33090c77898SPeter Brune xs = info.xs; 33190c77898SPeter Brune for (i = 0; i < m; i++) { 33290c77898SPeter Brune if (dim == 1) { 33390c77898SPeter Brune xm = M / m + ((M % m) > i); 33490c77898SPeter Brune } else if (dim == 2) { 33590c77898SPeter Brune xm = M / m + ((M % m) > i); 33690c77898SPeter Brune ym = N / n + ((N % n) > j); 33790c77898SPeter Brune } else if (dim == 3) { 33890c77898SPeter Brune xm = M / m + ((M % m) > i); 33990c77898SPeter Brune ym = N / n + ((N % n) > j); 34090c77898SPeter Brune zm = P / p + ((P % p) > k); 34190c77898SPeter Brune } 34290c77898SPeter Brune 34390c77898SPeter Brune xsize = xm; 34490c77898SPeter Brune ysize = ym; 34590c77898SPeter Brune zsize = zm; 34690c77898SPeter Brune xo = xs; 34790c77898SPeter Brune yo = ys; 34890c77898SPeter Brune zo = zs; 34990c77898SPeter Brune 3509566063dSJacob Faibussowitsch PetscCall(DMDACreate(PETSC_COMM_SELF, &(da[idx]))); 3519566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(da[idx], "sub_")); 3529566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da[idx], info.dim)); 3539566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da[idx], info.dof)); 35490c77898SPeter Brune 3559566063dSJacob Faibussowitsch PetscCall(DMDASetStencilType(da[idx], info.st)); 3569566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da[idx], info.sw)); 35790c77898SPeter Brune 358bff4a2f0SMatthew G. Knepley if (info.bx == DM_BOUNDARY_PERIODIC || (xs != 0)) { 3597ddda789SPeter Brune xsize += xol; 3607ddda789SPeter Brune xo -= xol; 361e30e807fSPeter Brune } 362bff4a2f0SMatthew G. Knepley if (info.by == DM_BOUNDARY_PERIODIC || (ys != 0)) { 3637ddda789SPeter Brune ysize += yol; 3647ddda789SPeter Brune yo -= yol; 365e30e807fSPeter Brune } 366bff4a2f0SMatthew G. Knepley if (info.bz == DM_BOUNDARY_PERIODIC || (zs != 0)) { 3677ddda789SPeter Brune zsize += zol; 3687ddda789SPeter Brune zo -= zol; 369e30e807fSPeter Brune } 370e30e807fSPeter Brune 371bff4a2f0SMatthew G. Knepley if (info.bx == DM_BOUNDARY_PERIODIC || (xs + xm != info.mx)) xsize += xol; 372bff4a2f0SMatthew G. Knepley if (info.by == DM_BOUNDARY_PERIODIC || (ys + ym != info.my)) ysize += yol; 373bff4a2f0SMatthew G. Knepley if (info.bz == DM_BOUNDARY_PERIODIC || (zs + zm != info.mz)) zsize += zol; 3748865f1eaSKarl Rupp 375bff4a2f0SMatthew G. Knepley if (info.bx != DM_BOUNDARY_PERIODIC) { 37690c77898SPeter Brune if (xo < 0) { 37790c77898SPeter Brune xsize += xo; 37890c77898SPeter Brune xo = 0; 37990c77898SPeter Brune } 380ad540459SPierre Jolivet if (xo + xsize > info.mx - 1) xsize -= xo + xsize - info.mx; 38190c77898SPeter Brune } 382bff4a2f0SMatthew G. Knepley if (info.by != DM_BOUNDARY_PERIODIC) { 38390c77898SPeter Brune if (yo < 0) { 38490c77898SPeter Brune ysize += yo; 38590c77898SPeter Brune yo = 0; 38690c77898SPeter Brune } 387ad540459SPierre Jolivet if (yo + ysize > info.my - 1) ysize -= yo + ysize - info.my; 38890c77898SPeter Brune } 389bff4a2f0SMatthew G. Knepley if (info.bz != DM_BOUNDARY_PERIODIC) { 39090c77898SPeter Brune if (zo < 0) { 39190c77898SPeter Brune zsize += zo; 39290c77898SPeter Brune zo = 0; 39390c77898SPeter Brune } 394ad540459SPierre Jolivet if (zo + zsize > info.mz - 1) zsize -= zo + zsize - info.mz; 39590c77898SPeter Brune } 39690c77898SPeter Brune 3979566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da[idx], xsize, ysize, zsize)); 3989566063dSJacob Faibussowitsch PetscCall(DMDASetNumProcs(da[idx], 1, 1, 1)); 3999566063dSJacob Faibussowitsch PetscCall(DMDASetBoundaryType(da[idx], DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED)); 400e30e807fSPeter Brune 401e30e807fSPeter Brune /* set up as a block instead */ 4029566063dSJacob Faibussowitsch PetscCall(DMSetUp(da[idx])); 403e30e807fSPeter Brune 40440234c92SPeter Brune /* nonoverlapping region */ 4059566063dSJacob Faibussowitsch PetscCall(DMDASetNonOverlappingRegion(da[idx], xs, ys, zs, xm, ym, zm)); 40640234c92SPeter Brune 40795c13181SPeter Brune /* this alters the behavior of DMDAGetInfo, DMDAGetLocalInfo, DMDAGetCorners, and DMDAGetGhostedCorners and should be used with care */ 4089566063dSJacob Faibussowitsch PetscCall(DMDASetOffset(da[idx], xo, yo, zo, info.mx, info.my, info.mz)); 40990c77898SPeter Brune xs += xm; 41090c77898SPeter Brune idx++; 41190c77898SPeter Brune } 41290c77898SPeter Brune ys += ym; 41390c77898SPeter Brune } 41490c77898SPeter Brune zs += zm; 41590c77898SPeter Brune } 41690c77898SPeter Brune *sdm = da; 4173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 418e30e807fSPeter Brune } 419e30e807fSPeter Brune 420e30e807fSPeter Brune /* 421e30e807fSPeter Brune Fills the local vector problem on the subdomain from the global problem. 422e30e807fSPeter Brune 423285ae305SPeter Brune Right now this assumes one subdomain per processor. 424285ae305SPeter Brune 425e30e807fSPeter Brune */ 426d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm, PetscInt nsubdms, DM *subdms, VecScatter **iscat, VecScatter **oscat, VecScatter **lscat) 427d71ae5a4SJacob Faibussowitsch { 428285ae305SPeter Brune DMDALocalInfo info, subinfo; 429e30e807fSPeter Brune DM subdm; 430285ae305SPeter Brune MatStencil upper, lower; 431285ae305SPeter Brune IS idis, isis, odis, osis, gdis; 432285ae305SPeter Brune Vec svec, dvec, slvec; 43340234c92SPeter Brune PetscInt xm, ym, zm, xs, ys, zs; 43490c77898SPeter Brune PetscInt i; 4353b5e53cdSSajid Ali PetscBool patchis_offproc = PETSC_TRUE; 436e30e807fSPeter Brune 437e30e807fSPeter Brune PetscFunctionBegin; 438e30e807fSPeter Brune /* allocate the arrays of scatters */ 4399566063dSJacob Faibussowitsch if (iscat) PetscCall(PetscMalloc1(nsubdms, iscat)); 4409566063dSJacob Faibussowitsch if (oscat) PetscCall(PetscMalloc1(nsubdms, oscat)); 4419566063dSJacob Faibussowitsch if (lscat) PetscCall(PetscMalloc1(nsubdms, lscat)); 442e30e807fSPeter Brune 4439566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 44490c77898SPeter Brune for (i = 0; i < nsubdms; i++) { 44590c77898SPeter Brune subdm = subdms[i]; 4469566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(subdm, &subinfo)); 4479566063dSJacob Faibussowitsch PetscCall(DMDAGetNonOverlappingRegion(subdm, &xs, &ys, &zs, &xm, &ym, &zm)); 448e30e807fSPeter Brune 449285ae305SPeter Brune /* create the global and subdomain index sets for the inner domain */ 45040234c92SPeter Brune lower.i = xs; 45140234c92SPeter Brune lower.j = ys; 45240234c92SPeter Brune lower.k = zs; 45340234c92SPeter Brune upper.i = xs + xm; 45440234c92SPeter Brune upper.j = ys + ym; 45540234c92SPeter Brune upper.k = zs + zm; 4569566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &idis, patchis_offproc)); 4579566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(subdm, &lower, &upper, &isis, patchis_offproc)); 458e30e807fSPeter Brune 459285ae305SPeter Brune /* create the global and subdomain index sets for the outer subdomain */ 460285ae305SPeter Brune lower.i = subinfo.xs; 461285ae305SPeter Brune lower.j = subinfo.ys; 462285ae305SPeter Brune lower.k = subinfo.zs; 463285ae305SPeter Brune upper.i = subinfo.xs + subinfo.xm; 464285ae305SPeter Brune upper.j = subinfo.ys + subinfo.ym; 465285ae305SPeter Brune upper.k = subinfo.zs + subinfo.zm; 4669566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &odis, patchis_offproc)); 4679566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(subdm, &lower, &upper, &osis, patchis_offproc)); 468e30e807fSPeter Brune 469285ae305SPeter Brune /* global and subdomain ISes for the local indices of the subdomain */ 470285ae305SPeter Brune /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */ 471285ae305SPeter Brune lower.i = subinfo.gxs; 472285ae305SPeter Brune lower.j = subinfo.gys; 473285ae305SPeter Brune lower.k = subinfo.gzs; 474285ae305SPeter Brune upper.i = subinfo.gxs + subinfo.gxm; 475285ae305SPeter Brune upper.j = subinfo.gys + subinfo.gym; 476285ae305SPeter Brune upper.k = subinfo.gzs + subinfo.gzm; 4779566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &gdis, patchis_offproc)); 478e30e807fSPeter Brune 479e30e807fSPeter Brune /* form the scatter */ 4809566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &dvec)); 4819566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(subdm, &svec)); 4829566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(subdm, &slvec)); 483e30e807fSPeter Brune 4849566063dSJacob Faibussowitsch if (iscat) PetscCall(VecScatterCreate(dvec, idis, svec, isis, &(*iscat)[i])); 4859566063dSJacob Faibussowitsch if (oscat) PetscCall(VecScatterCreate(dvec, odis, svec, osis, &(*oscat)[i])); 4869566063dSJacob Faibussowitsch if (lscat) PetscCall(VecScatterCreate(dvec, gdis, slvec, NULL, &(*lscat)[i])); 487e30e807fSPeter Brune 4889566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &dvec)); 4899566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(subdm, &svec)); 4909566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(subdm, &slvec)); 491e30e807fSPeter Brune 4929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&idis)); 4939566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isis)); 494e30e807fSPeter Brune 4959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&odis)); 4969566063dSJacob Faibussowitsch PetscCall(ISDestroy(&osis)); 497e30e807fSPeter Brune 4989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&gdis)); 49990c77898SPeter Brune } 5003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 501e30e807fSPeter Brune } 502e30e807fSPeter Brune 50366976f2fSJacob Faibussowitsch static PetscErrorCode DMDASubDomainIS_Private(DM dm, PetscInt n, DM *subdm, IS **iis, IS **ois) 504d71ae5a4SJacob Faibussowitsch { 50590c77898SPeter Brune PetscInt i; 506e30e807fSPeter Brune DMDALocalInfo info, subinfo; 507285ae305SPeter Brune MatStencil lower, upper; 5083b5e53cdSSajid Ali PetscBool patchis_offproc = PETSC_TRUE; 509e30e807fSPeter Brune 510e30e807fSPeter Brune PetscFunctionBegin; 5119566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 5129566063dSJacob Faibussowitsch if (iis) PetscCall(PetscMalloc1(n, iis)); 5139566063dSJacob Faibussowitsch if (ois) PetscCall(PetscMalloc1(n, ois)); 514e30e807fSPeter Brune 51590c77898SPeter Brune for (i = 0; i < n; i++) { 5169566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(subdm[i], &subinfo)); 51790c77898SPeter Brune if (iis) { 518285ae305SPeter Brune /* create the inner IS */ 519285ae305SPeter Brune lower.i = info.xs; 520285ae305SPeter Brune lower.j = info.ys; 521285ae305SPeter Brune lower.k = info.zs; 522285ae305SPeter Brune upper.i = info.xs + info.xm; 523285ae305SPeter Brune upper.j = info.ys + info.ym; 524285ae305SPeter Brune upper.k = info.zs + info.zm; 5259566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &(*iis)[i], patchis_offproc)); 52690c77898SPeter Brune } 527e30e807fSPeter Brune 52890c77898SPeter Brune if (ois) { 529285ae305SPeter Brune /* create the outer IS */ 530285ae305SPeter Brune lower.i = subinfo.xs; 531285ae305SPeter Brune lower.j = subinfo.ys; 532285ae305SPeter Brune lower.k = subinfo.zs; 533285ae305SPeter Brune upper.i = subinfo.xs + subinfo.xm; 534285ae305SPeter Brune upper.j = subinfo.ys + subinfo.ym; 535285ae305SPeter Brune upper.k = subinfo.zs + subinfo.zm; 5369566063dSJacob Faibussowitsch PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &(*ois)[i], patchis_offproc)); 53790c77898SPeter Brune } 53890c77898SPeter Brune } 5393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 540e30e807fSPeter Brune } 541e30e807fSPeter Brune 542d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateDomainDecomposition_DA(DM dm, PetscInt *len, char ***names, IS **iis, IS **ois, DM **subdm) 543d71ae5a4SJacob Faibussowitsch { 54466976f2fSJacob Faibussowitsch DM *sdm = NULL; 54566976f2fSJacob Faibussowitsch PetscInt n = 0, i; 5460a696f66SPeter Brune 5470adebc6cSBarry Smith PetscFunctionBegin; 5489566063dSJacob Faibussowitsch PetscCall(DMDASubDomainDA_Private(dm, &n, &sdm)); 54990c77898SPeter Brune if (names) { 5509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, names)); 551ea78f98cSLisandro Dalcin for (i = 0; i < n; i++) (*names)[i] = NULL; 552e30e807fSPeter Brune } 5539566063dSJacob Faibussowitsch PetscCall(DMDASubDomainIS_Private(dm, n, sdm, iis, ois)); 55490c77898SPeter Brune if (subdm) *subdm = sdm; 555e2c616c8SPeter Brune else { 55648a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(DMDestroy(&sdm[i])); 557e2c616c8SPeter Brune } 55890c77898SPeter Brune if (len) *len = n; 5593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 560e30e807fSPeter Brune } 561