xref: /petsc/src/dm/impls/da/dadd.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"   I*/
2e30e807fSPeter Brune 
3063654ddSPeter Brune /*@
4063654ddSPeter Brune   DMDACreatePatchIS - Creates an index set corresponding to a patch of the DA.
5ff9846d9SPeter Brune 
63b5e53cdSSajid Ali   Collective
7063654ddSPeter Brune 
8063654ddSPeter Brune   Input Parameters:
9063654ddSPeter Brune +  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 
14063654ddSPeter Brune   Output Parameters:
15063654ddSPeter Brune .  is - the IS corresponding to the patch
16063654ddSPeter Brune 
17063654ddSPeter Brune   Level: developer
18063654ddSPeter Brune 
193b5e53cdSSajid Ali Notes:
203b5e53cdSSajid Ali This routine always returns an IS on the DMDA's comm, if offproc is set to PETSC_TRUE,
213b5e53cdSSajid Ali 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
233b5e53cdSSajid Ali the indices returned in this mode are appropriate. If offproc is set to PETSC_FALSE,
243b5e53cdSSajid Ali 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 
272b3cbbdaSStefano Zampini .seealso: `DMCreateDomainDecomposition()`, `DMCreateDomainDecompositionScatters()`
28063654ddSPeter Brune @*/
299371c9d4SSatish Balay PetscErrorCode DMDACreatePatchIS(DM da, MatStencil *lower, MatStencil *upper, IS *is, PetscBool offproc) {
30063654ddSPeter Brune   PetscInt        ms = 0, ns = 0, ps = 0;
313b5e53cdSSajid Ali   PetscInt        mw = 0, nw = 0, pw = 0;
32063654ddSPeter Brune   PetscInt        me = 1, ne = 1, pe = 1;
33063654ddSPeter Brune   PetscInt        mr = 0, nr = 0, pr = 0;
34ff9846d9SPeter Brune   PetscInt        ii, jj, kk;
35063654ddSPeter Brune   PetscInt        si, sj, sk;
363b5e53cdSSajid Ali   PetscInt        i, j, k, l, idx = 0;
37ff9846d9SPeter Brune   PetscInt        base;
38063654ddSPeter Brune   PetscInt        xm = 1, ym = 1, zm = 1;
39ff9846d9SPeter Brune   PetscInt        ox, oy, oz;
40063654ddSPeter Brune   PetscInt        m, n, p, M, N, P, dof;
413b5e53cdSSajid Ali   const PetscInt *lx, *ly, *lz;
42063654ddSPeter Brune   PetscInt        nindices;
43063654ddSPeter Brune   PetscInt       *indices;
44063654ddSPeter Brune   DM_DA          *dd     = (DM_DA *)da->data;
453b5e53cdSSajid Ali   PetscBool       skip_i = PETSC_TRUE, skip_j = PETSC_TRUE, skip_k = PETSC_TRUE;
463b5e53cdSSajid Ali   PetscBool       valid_j = PETSC_FALSE, valid_k = PETSC_FALSE; /* DMDA has at least 1 dimension */
47ff9846d9SPeter Brune 
480adebc6cSBarry Smith   PetscFunctionBegin;
499371c9d4SSatish Balay   M   = dd->M;
509371c9d4SSatish Balay   N   = dd->N;
519371c9d4SSatish Balay   P   = dd->P;
529371c9d4SSatish Balay   m   = dd->m;
539371c9d4SSatish Balay   n   = dd->n;
549371c9d4SSatish Balay   p   = dd->p;
55063654ddSPeter Brune   dof = dd->w;
563b5e53cdSSajid Ali 
573b5e53cdSSajid Ali   nindices = -1;
583b5e53cdSSajid Ali   if (PetscLikely(upper->i - lower->i)) {
593b5e53cdSSajid Ali     nindices = nindices * (upper->i - lower->i);
603b5e53cdSSajid Ali     skip_i   = PETSC_FALSE;
613b5e53cdSSajid Ali   }
623b5e53cdSSajid Ali   if (N > 1) {
633b5e53cdSSajid Ali     valid_j = PETSC_TRUE;
643b5e53cdSSajid Ali     if (PetscLikely(upper->j - lower->j)) {
653b5e53cdSSajid Ali       nindices = nindices * (upper->j - lower->j);
663b5e53cdSSajid Ali       skip_j   = PETSC_FALSE;
673b5e53cdSSajid Ali     }
683b5e53cdSSajid Ali   }
693b5e53cdSSajid Ali   if (P > 1) {
703b5e53cdSSajid Ali     valid_k = PETSC_TRUE;
713b5e53cdSSajid Ali     if (PetscLikely(upper->k - lower->k)) {
723b5e53cdSSajid Ali       nindices = nindices * (upper->k - lower->k);
733b5e53cdSSajid Ali       skip_k   = PETSC_FALSE;
743b5e53cdSSajid Ali     }
753b5e53cdSSajid Ali   }
763b5e53cdSSajid Ali   if (PetscLikely(nindices < 0)) {
773b5e53cdSSajid Ali     if (PetscUnlikely(skip_i && skip_j && skip_k)) {
783b5e53cdSSajid Ali       nindices = 0;
793b5e53cdSSajid Ali     } else nindices = nindices * (-1);
803b5e53cdSSajid Ali   } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Lower and Upper stencils are identical! Please check inputs.");
813b5e53cdSSajid Ali 
829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nindices * dof, &indices));
839566063dSJacob Faibussowitsch   PetscCall(DMDAGetOffset(da, &ox, &oy, &oz, NULL, NULL, NULL));
843b5e53cdSSajid Ali 
853b5e53cdSSajid Ali   if (!valid_k) {
863b5e53cdSSajid Ali     k        = 0;
873b5e53cdSSajid Ali     upper->k = 0;
883b5e53cdSSajid Ali     lower->k = 0;
893b5e53cdSSajid Ali   }
903b5e53cdSSajid Ali   if (!valid_j) {
913b5e53cdSSajid Ali     j        = 0;
923b5e53cdSSajid Ali     upper->j = 0;
933b5e53cdSSajid Ali     lower->j = 0;
943b5e53cdSSajid Ali   }
953b5e53cdSSajid Ali 
963b5e53cdSSajid Ali   if (offproc) {
979566063dSJacob Faibussowitsch     PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, &lz));
98063654ddSPeter Brune     /* start at index 0 on processor 0 */
99063654ddSPeter Brune     mr = 0;
100063654ddSPeter Brune     nr = 0;
101063654ddSPeter Brune     pr = 0;
102063654ddSPeter Brune     ms = 0;
103063654ddSPeter Brune     ns = 0;
104063654ddSPeter Brune     ps = 0;
105063654ddSPeter Brune     if (lx) me = lx[0];
106063654ddSPeter Brune     if (ly) ne = ly[0];
107063654ddSPeter Brune     if (lz) pe = lz[0];
1083b5e53cdSSajid Ali     /*
1093b5e53cdSSajid Ali        If no indices are to be returned, create an empty is,
1103b5e53cdSSajid Ali        this prevents hanging in while loops
1113b5e53cdSSajid Ali     */
1123b5e53cdSSajid Ali     if (skip_i && skip_j && skip_k) goto createis;
1133b5e53cdSSajid Ali     /*
1143b5e53cdSSajid Ali        do..while loops to ensure the block gets entered once,
1153b5e53cdSSajid Ali        regardless of control condition being met, necessary for
1163b5e53cdSSajid Ali        cases when a subset of skip_i/j/k is true
1173b5e53cdSSajid Ali     */
1189371c9d4SSatish Balay     if (skip_k) k = upper->k - oz;
1199371c9d4SSatish Balay     else k = lower->k - oz;
1203b5e53cdSSajid Ali     do {
1219371c9d4SSatish Balay       if (skip_j) j = upper->j - oy;
1229371c9d4SSatish Balay       else j = lower->j - oy;
1233b5e53cdSSajid Ali       do {
1249371c9d4SSatish Balay         if (skip_i) i = upper->i - ox;
1259371c9d4SSatish Balay         else i = lower->i - ox;
1263b5e53cdSSajid Ali         do {
127063654ddSPeter Brune           /* "actual" indices rather than ones outside of the domain */
128063654ddSPeter Brune           ii = i;
129063654ddSPeter Brune           jj = j;
130063654ddSPeter Brune           kk = k;
131063654ddSPeter Brune           if (ii < 0) ii = M + ii;
132063654ddSPeter Brune           if (jj < 0) jj = N + jj;
133063654ddSPeter Brune           if (kk < 0) kk = P + kk;
134063654ddSPeter Brune           if (ii > M - 1) ii = ii - M;
135063654ddSPeter Brune           if (jj > N - 1) jj = jj - N;
136063654ddSPeter Brune           if (kk > P - 1) kk = kk - P;
137063654ddSPeter Brune           /* gone out of processor range on x axis */
138063654ddSPeter Brune           while (ii > me - 1 || ii < ms) {
139063654ddSPeter Brune             if (mr == m - 1) {
140063654ddSPeter Brune               ms = 0;
141063654ddSPeter Brune               me = lx[0];
142063654ddSPeter Brune               mr = 0;
143063654ddSPeter Brune             } else {
144063654ddSPeter Brune               mr++;
145063654ddSPeter Brune               ms = me;
146063654ddSPeter Brune               me += lx[mr];
147063654ddSPeter Brune             }
148063654ddSPeter Brune           }
149063654ddSPeter Brune           /* gone out of processor range on y axis */
150063654ddSPeter Brune           while (jj > ne - 1 || jj < ns) {
151063654ddSPeter Brune             if (nr == n - 1) {
152063654ddSPeter Brune               ns = 0;
153063654ddSPeter Brune               ne = ly[0];
154063654ddSPeter Brune               nr = 0;
155063654ddSPeter Brune             } else {
156063654ddSPeter Brune               nr++;
157063654ddSPeter Brune               ns = ne;
158063654ddSPeter Brune               ne += ly[nr];
159063654ddSPeter Brune             }
160063654ddSPeter Brune           }
161063654ddSPeter Brune           /* gone out of processor range on z axis */
162063654ddSPeter Brune           while (kk > pe - 1 || kk < ps) {
163063654ddSPeter Brune             if (pr == p - 1) {
164063654ddSPeter Brune               ps = 0;
165063654ddSPeter Brune               pe = lz[0];
166063654ddSPeter Brune               pr = 0;
167063654ddSPeter Brune             } else {
168063654ddSPeter Brune               pr++;
169063654ddSPeter Brune               ps = pe;
170063654ddSPeter Brune               pe += lz[pr];
171063654ddSPeter Brune             }
172063654ddSPeter Brune           }
173063654ddSPeter Brune           /* compute the vector base on owning processor */
174063654ddSPeter Brune           xm   = me - ms;
175063654ddSPeter Brune           ym   = ne - ns;
176063654ddSPeter Brune           zm   = pe - ps;
177b0bff951SPeter Brune           base = ms * ym * zm + ns * M * zm + ps * M * N;
178063654ddSPeter Brune           /* compute the local coordinates on owning processor */
179063654ddSPeter Brune           si   = ii - ms;
180063654ddSPeter Brune           sj   = jj - ns;
181063654ddSPeter Brune           sk   = kk - ps;
182063654ddSPeter Brune           for (l = 0; l < dof; l++) {
183063654ddSPeter Brune             indices[idx] = l + dof * (base + si + xm * sj + xm * ym * sk);
184ff9846d9SPeter Brune             idx++;
185ff9846d9SPeter Brune           }
1863b5e53cdSSajid Ali           i++;
1873b5e53cdSSajid Ali         } while (i < upper->i - ox);
1883b5e53cdSSajid Ali         j++;
1893b5e53cdSSajid Ali       } while (j < upper->j - oy);
1903b5e53cdSSajid Ali       k++;
1913b5e53cdSSajid Ali     } while (k < upper->k - oz);
1923b5e53cdSSajid Ali   }
1933b5e53cdSSajid Ali 
1943b5e53cdSSajid Ali   if (!offproc) {
1959566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(da, &ms, &ns, &ps, &mw, &nw, &pw));
1963b5e53cdSSajid Ali     me = ms + mw;
1973b5e53cdSSajid Ali     if (N > 1) ne = ns + nw;
1983b5e53cdSSajid Ali     if (P > 1) pe = ps + pw;
1993b5e53cdSSajid Ali     /* Account for DM offsets */
2009371c9d4SSatish Balay     ms = ms - ox;
2019371c9d4SSatish Balay     me = me - ox;
2029371c9d4SSatish Balay     ns = ns - oy;
2039371c9d4SSatish Balay     ne = ne - oy;
2049371c9d4SSatish Balay     ps = ps - oz;
2059371c9d4SSatish Balay     pe = pe - oz;
2063b5e53cdSSajid Ali 
2073b5e53cdSSajid Ali     /* compute the vector base on owning processor */
2083b5e53cdSSajid Ali     xm   = me - ms;
2093b5e53cdSSajid Ali     ym   = ne - ns;
2103b5e53cdSSajid Ali     zm   = pe - ps;
2113b5e53cdSSajid Ali     base = ms * ym * zm + ns * M * zm + ps * M * N;
2123b5e53cdSSajid Ali     /*
2133b5e53cdSSajid Ali        if no indices are to be returned, create an empty is,
2143b5e53cdSSajid Ali        this prevents hanging in while loops
2153b5e53cdSSajid Ali     */
2163b5e53cdSSajid Ali     if (skip_i && skip_j && skip_k) goto createis;
2173b5e53cdSSajid Ali     /*
2183b5e53cdSSajid Ali        do..while loops to ensure the block gets entered once,
2193b5e53cdSSajid Ali        regardless of control condition being met, necessary for
2203b5e53cdSSajid Ali        cases when a subset of skip_i/j/k is true
2213b5e53cdSSajid Ali     */
2229371c9d4SSatish Balay     if (skip_k) k = upper->k - oz;
2239371c9d4SSatish Balay     else k = lower->k - oz;
2243b5e53cdSSajid Ali     do {
2259371c9d4SSatish Balay       if (skip_j) j = upper->j - oy;
2269371c9d4SSatish Balay       else j = lower->j - oy;
2273b5e53cdSSajid Ali       do {
2289371c9d4SSatish Balay         if (skip_i) i = upper->i - ox;
2299371c9d4SSatish Balay         else i = lower->i - ox;
2303b5e53cdSSajid Ali         do {
2313b5e53cdSSajid Ali           if (k >= ps && k <= pe - 1) {
2323b5e53cdSSajid Ali             if (j >= ns && j <= ne - 1) {
2333b5e53cdSSajid Ali               if (i >= ms && i <= me - 1) {
2343b5e53cdSSajid Ali                 /* compute the local coordinates on owning processor */
2353b5e53cdSSajid Ali                 si = i - ms;
2363b5e53cdSSajid Ali                 sj = j - ns;
2373b5e53cdSSajid Ali                 sk = k - ps;
2383b5e53cdSSajid Ali                 for (l = 0; l < dof; l++) {
2393b5e53cdSSajid Ali                   indices[idx] = l + dof * (base + si + xm * sj + xm * ym * sk);
2403b5e53cdSSajid Ali                   idx++;
241ff9846d9SPeter Brune                 }
242ff9846d9SPeter Brune               }
243063654ddSPeter Brune             }
2443b5e53cdSSajid Ali           }
2453b5e53cdSSajid Ali           i++;
2463b5e53cdSSajid Ali         } while (i < upper->i - ox);
2473b5e53cdSSajid Ali         j++;
2483b5e53cdSSajid Ali       } while (j < upper->j - oy);
2493b5e53cdSSajid Ali       k++;
2503b5e53cdSSajid Ali     } while (k < upper->k - oz);
2513b5e53cdSSajid Ali 
2529566063dSJacob Faibussowitsch     PetscCall(PetscRealloc((size_t)(idx * sizeof(PetscInt)), (void *)&indices));
2533b5e53cdSSajid Ali   }
2543b5e53cdSSajid Ali 
2553b5e53cdSSajid Ali createis:
2569566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)da), idx, indices, PETSC_OWN_POINTER, is));
257ff9846d9SPeter Brune   PetscFunctionReturn(0);
258ff9846d9SPeter Brune }
259ff9846d9SPeter Brune 
2609371c9d4SSatish Balay PetscErrorCode DMDASubDomainDA_Private(DM dm, PetscInt *nlocal, DM **sdm) {
26190c77898SPeter Brune   DM           *da;
26290c77898SPeter Brune   PetscInt      dim, size, i, j, k, idx;
263e30e807fSPeter Brune   DMDALocalInfo info;
264ff9846d9SPeter Brune   PetscInt      xsize, ysize, zsize;
265e30e807fSPeter Brune   PetscInt      xo, yo, zo;
26690c77898SPeter Brune   PetscInt      xs, ys, zs;
26790c77898SPeter Brune   PetscInt      xm = 1, ym = 1, zm = 1;
2687ddda789SPeter Brune   PetscInt      xol, yol, zol;
26990c77898SPeter Brune   PetscInt      m = 1, n = 1, p = 1;
27090c77898SPeter Brune   PetscInt      M, N, P;
27190c77898SPeter Brune   PetscInt      pm, mtmp;
272e30e807fSPeter Brune 
273e30e807fSPeter Brune   PetscFunctionBegin;
2749566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
2759566063dSJacob Faibussowitsch   PetscCall(DMDAGetOverlap(dm, &xol, &yol, &zol));
2769566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumLocalSubDomains(dm, &size));
2779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &da));
2787ddda789SPeter Brune 
27990c77898SPeter Brune   if (nlocal) *nlocal = size;
280e30e807fSPeter Brune 
28190c77898SPeter Brune   dim = info.dim;
282e30e807fSPeter Brune 
28390c77898SPeter Brune   M = info.xm;
28490c77898SPeter Brune   N = info.ym;
28590c77898SPeter Brune   P = info.zm;
28690c77898SPeter Brune 
28790c77898SPeter Brune   if (dim == 1) {
28890c77898SPeter Brune     m = size;
28990c77898SPeter Brune   } else if (dim == 2) {
29090c77898SPeter Brune     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)N)));
29190c77898SPeter Brune     while (m > 0) {
29290c77898SPeter Brune       n = size / m;
29390c77898SPeter Brune       if (m * n * p == size) break;
29490c77898SPeter Brune       m--;
29590c77898SPeter Brune     }
29690c77898SPeter Brune   } else if (dim == 3) {
2979371c9d4SSatish Balay     n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N * N) * ((PetscReal)size) / ((PetscReal)P * M), (PetscReal)(1. / 3.)));
2989371c9d4SSatish Balay     if (!n) n = 1;
29990c77898SPeter Brune     while (n > 0) {
30090c77898SPeter Brune       pm = size / n;
30190c77898SPeter Brune       if (n * pm == size) break;
30290c77898SPeter Brune       n--;
30390c77898SPeter Brune     }
30490c77898SPeter Brune     if (!n) n = 1;
30590c77898SPeter Brune     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)P * n)));
30690c77898SPeter Brune     if (!m) m = 1;
30790c77898SPeter Brune     while (m > 0) {
30890c77898SPeter Brune       p = size / (m * n);
30990c77898SPeter Brune       if (m * n * p == size) break;
31090c77898SPeter Brune       m--;
31190c77898SPeter Brune     }
3129371c9d4SSatish Balay     if (M > P && m < p) {
3139371c9d4SSatish Balay       mtmp = m;
3149371c9d4SSatish Balay       m    = p;
3159371c9d4SSatish Balay       p    = mtmp;
3169371c9d4SSatish Balay     }
31790c77898SPeter Brune   }
31890c77898SPeter Brune 
31990c77898SPeter Brune   zs  = info.zs;
32090c77898SPeter Brune   idx = 0;
32190c77898SPeter Brune   for (k = 0; k < p; k++) {
32290c77898SPeter Brune     ys = info.ys;
32390c77898SPeter Brune     for (j = 0; j < n; j++) {
32490c77898SPeter Brune       xs = info.xs;
32590c77898SPeter Brune       for (i = 0; i < m; i++) {
32690c77898SPeter Brune         if (dim == 1) {
32790c77898SPeter Brune           xm = M / m + ((M % m) > i);
32890c77898SPeter Brune         } else if (dim == 2) {
32990c77898SPeter Brune           xm = M / m + ((M % m) > i);
33090c77898SPeter Brune           ym = N / n + ((N % n) > j);
33190c77898SPeter Brune         } else if (dim == 3) {
33290c77898SPeter Brune           xm = M / m + ((M % m) > i);
33390c77898SPeter Brune           ym = N / n + ((N % n) > j);
33490c77898SPeter Brune           zm = P / p + ((P % p) > k);
33590c77898SPeter Brune         }
33690c77898SPeter Brune 
33790c77898SPeter Brune         xsize = xm;
33890c77898SPeter Brune         ysize = ym;
33990c77898SPeter Brune         zsize = zm;
34090c77898SPeter Brune         xo    = xs;
34190c77898SPeter Brune         yo    = ys;
34290c77898SPeter Brune         zo    = zs;
34390c77898SPeter Brune 
3449566063dSJacob Faibussowitsch         PetscCall(DMDACreate(PETSC_COMM_SELF, &(da[idx])));
3459566063dSJacob Faibussowitsch         PetscCall(DMSetOptionsPrefix(da[idx], "sub_"));
3469566063dSJacob Faibussowitsch         PetscCall(DMSetDimension(da[idx], info.dim));
3479566063dSJacob Faibussowitsch         PetscCall(DMDASetDof(da[idx], info.dof));
34890c77898SPeter Brune 
3499566063dSJacob Faibussowitsch         PetscCall(DMDASetStencilType(da[idx], info.st));
3509566063dSJacob Faibussowitsch         PetscCall(DMDASetStencilWidth(da[idx], info.sw));
35190c77898SPeter Brune 
352bff4a2f0SMatthew G. Knepley         if (info.bx == DM_BOUNDARY_PERIODIC || (xs != 0)) {
3537ddda789SPeter Brune           xsize += xol;
3547ddda789SPeter Brune           xo -= xol;
355e30e807fSPeter Brune         }
356bff4a2f0SMatthew G. Knepley         if (info.by == DM_BOUNDARY_PERIODIC || (ys != 0)) {
3577ddda789SPeter Brune           ysize += yol;
3587ddda789SPeter Brune           yo -= yol;
359e30e807fSPeter Brune         }
360bff4a2f0SMatthew G. Knepley         if (info.bz == DM_BOUNDARY_PERIODIC || (zs != 0)) {
3617ddda789SPeter Brune           zsize += zol;
3627ddda789SPeter Brune           zo -= zol;
363e30e807fSPeter Brune         }
364e30e807fSPeter Brune 
365bff4a2f0SMatthew G. Knepley         if (info.bx == DM_BOUNDARY_PERIODIC || (xs + xm != info.mx)) xsize += xol;
366bff4a2f0SMatthew G. Knepley         if (info.by == DM_BOUNDARY_PERIODIC || (ys + ym != info.my)) ysize += yol;
367bff4a2f0SMatthew G. Knepley         if (info.bz == DM_BOUNDARY_PERIODIC || (zs + zm != info.mz)) zsize += zol;
3688865f1eaSKarl Rupp 
369bff4a2f0SMatthew G. Knepley         if (info.bx != DM_BOUNDARY_PERIODIC) {
37090c77898SPeter Brune           if (xo < 0) {
37190c77898SPeter Brune             xsize += xo;
37290c77898SPeter Brune             xo = 0;
37390c77898SPeter Brune           }
3749371c9d4SSatish Balay           if (xo + xsize > info.mx - 1) { xsize -= xo + xsize - info.mx; }
37590c77898SPeter Brune         }
376bff4a2f0SMatthew G. Knepley         if (info.by != DM_BOUNDARY_PERIODIC) {
37790c77898SPeter Brune           if (yo < 0) {
37890c77898SPeter Brune             ysize += yo;
37990c77898SPeter Brune             yo = 0;
38090c77898SPeter Brune           }
3819371c9d4SSatish Balay           if (yo + ysize > info.my - 1) { ysize -= yo + ysize - info.my; }
38290c77898SPeter Brune         }
383bff4a2f0SMatthew G. Knepley         if (info.bz != DM_BOUNDARY_PERIODIC) {
38490c77898SPeter Brune           if (zo < 0) {
38590c77898SPeter Brune             zsize += zo;
38690c77898SPeter Brune             zo = 0;
38790c77898SPeter Brune           }
3889371c9d4SSatish Balay           if (zo + zsize > info.mz - 1) { zsize -= zo + zsize - info.mz; }
38990c77898SPeter Brune         }
39090c77898SPeter Brune 
3919566063dSJacob Faibussowitsch         PetscCall(DMDASetSizes(da[idx], xsize, ysize, zsize));
3929566063dSJacob Faibussowitsch         PetscCall(DMDASetNumProcs(da[idx], 1, 1, 1));
3939566063dSJacob Faibussowitsch         PetscCall(DMDASetBoundaryType(da[idx], DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED));
394e30e807fSPeter Brune 
395e30e807fSPeter Brune         /* set up as a block instead */
3969566063dSJacob Faibussowitsch         PetscCall(DMSetUp(da[idx]));
397e30e807fSPeter Brune 
39840234c92SPeter Brune         /* nonoverlapping region */
3999566063dSJacob Faibussowitsch         PetscCall(DMDASetNonOverlappingRegion(da[idx], xs, ys, zs, xm, ym, zm));
40040234c92SPeter Brune 
40195c13181SPeter Brune         /* this alters the behavior of DMDAGetInfo, DMDAGetLocalInfo, DMDAGetCorners, and DMDAGetGhostedCorners and should be used with care */
4029566063dSJacob Faibussowitsch         PetscCall(DMDASetOffset(da[idx], xo, yo, zo, info.mx, info.my, info.mz));
40390c77898SPeter Brune         xs += xm;
40490c77898SPeter Brune         idx++;
40590c77898SPeter Brune       }
40690c77898SPeter Brune       ys += ym;
40790c77898SPeter Brune     }
40890c77898SPeter Brune     zs += zm;
40990c77898SPeter Brune   }
41090c77898SPeter Brune   *sdm = da;
411e30e807fSPeter Brune   PetscFunctionReturn(0);
412e30e807fSPeter Brune }
413e30e807fSPeter Brune 
414e30e807fSPeter Brune /*
415e30e807fSPeter Brune    Fills the local vector problem on the subdomain from the global problem.
416e30e807fSPeter Brune 
417285ae305SPeter Brune    Right now this assumes one subdomain per processor.
418285ae305SPeter Brune 
419e30e807fSPeter Brune */
4209371c9d4SSatish Balay PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm, PetscInt nsubdms, DM *subdms, VecScatter **iscat, VecScatter **oscat, VecScatter **lscat) {
421285ae305SPeter Brune   DMDALocalInfo info, subinfo;
422e30e807fSPeter Brune   DM            subdm;
423285ae305SPeter Brune   MatStencil    upper, lower;
424285ae305SPeter Brune   IS            idis, isis, odis, osis, gdis;
425285ae305SPeter Brune   Vec           svec, dvec, slvec;
42640234c92SPeter Brune   PetscInt      xm, ym, zm, xs, ys, zs;
42790c77898SPeter Brune   PetscInt      i;
4283b5e53cdSSajid Ali   PetscBool     patchis_offproc = PETSC_TRUE;
429e30e807fSPeter Brune 
430e30e807fSPeter Brune   PetscFunctionBegin;
431e30e807fSPeter Brune   /* allocate the arrays of scatters */
4329566063dSJacob Faibussowitsch   if (iscat) PetscCall(PetscMalloc1(nsubdms, iscat));
4339566063dSJacob Faibussowitsch   if (oscat) PetscCall(PetscMalloc1(nsubdms, oscat));
4349566063dSJacob Faibussowitsch   if (lscat) PetscCall(PetscMalloc1(nsubdms, lscat));
435e30e807fSPeter Brune 
4369566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
43790c77898SPeter Brune   for (i = 0; i < nsubdms; i++) {
43890c77898SPeter Brune     subdm = subdms[i];
4399566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(subdm, &subinfo));
4409566063dSJacob Faibussowitsch     PetscCall(DMDAGetNonOverlappingRegion(subdm, &xs, &ys, &zs, &xm, &ym, &zm));
441e30e807fSPeter Brune 
442285ae305SPeter Brune     /* create the global and subdomain index sets for the inner domain */
44340234c92SPeter Brune     lower.i = xs;
44440234c92SPeter Brune     lower.j = ys;
44540234c92SPeter Brune     lower.k = zs;
44640234c92SPeter Brune     upper.i = xs + xm;
44740234c92SPeter Brune     upper.j = ys + ym;
44840234c92SPeter Brune     upper.k = zs + zm;
4499566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &idis, patchis_offproc));
4509566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(subdm, &lower, &upper, &isis, patchis_offproc));
451e30e807fSPeter Brune 
452285ae305SPeter Brune     /* create the global and subdomain index sets for the outer subdomain */
453285ae305SPeter Brune     lower.i = subinfo.xs;
454285ae305SPeter Brune     lower.j = subinfo.ys;
455285ae305SPeter Brune     lower.k = subinfo.zs;
456285ae305SPeter Brune     upper.i = subinfo.xs + subinfo.xm;
457285ae305SPeter Brune     upper.j = subinfo.ys + subinfo.ym;
458285ae305SPeter Brune     upper.k = subinfo.zs + subinfo.zm;
4599566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &odis, patchis_offproc));
4609566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(subdm, &lower, &upper, &osis, patchis_offproc));
461e30e807fSPeter Brune 
462285ae305SPeter Brune     /* global and subdomain ISes for the local indices of the subdomain */
463285ae305SPeter Brune     /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */
464285ae305SPeter Brune     lower.i = subinfo.gxs;
465285ae305SPeter Brune     lower.j = subinfo.gys;
466285ae305SPeter Brune     lower.k = subinfo.gzs;
467285ae305SPeter Brune     upper.i = subinfo.gxs + subinfo.gxm;
468285ae305SPeter Brune     upper.j = subinfo.gys + subinfo.gym;
469285ae305SPeter Brune     upper.k = subinfo.gzs + subinfo.gzm;
4709566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &gdis, patchis_offproc));
471e30e807fSPeter Brune 
472e30e807fSPeter Brune     /* form the scatter */
4739566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &dvec));
4749566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(subdm, &svec));
4759566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(subdm, &slvec));
476e30e807fSPeter Brune 
4779566063dSJacob Faibussowitsch     if (iscat) PetscCall(VecScatterCreate(dvec, idis, svec, isis, &(*iscat)[i]));
4789566063dSJacob Faibussowitsch     if (oscat) PetscCall(VecScatterCreate(dvec, odis, svec, osis, &(*oscat)[i]));
4799566063dSJacob Faibussowitsch     if (lscat) PetscCall(VecScatterCreate(dvec, gdis, slvec, NULL, &(*lscat)[i]));
480e30e807fSPeter Brune 
4819566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &dvec));
4829566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(subdm, &svec));
4839566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(subdm, &slvec));
484e30e807fSPeter Brune 
4859566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&idis));
4869566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isis));
487e30e807fSPeter Brune 
4889566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&odis));
4899566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&osis));
490e30e807fSPeter Brune 
4919566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&gdis));
49290c77898SPeter Brune   }
493e30e807fSPeter Brune   PetscFunctionReturn(0);
494e30e807fSPeter Brune }
495e30e807fSPeter Brune 
4969371c9d4SSatish Balay PetscErrorCode DMDASubDomainIS_Private(DM dm, PetscInt n, DM *subdm, IS **iis, IS **ois) {
49790c77898SPeter Brune   PetscInt      i;
498e30e807fSPeter Brune   DMDALocalInfo info, subinfo;
499285ae305SPeter Brune   MatStencil    lower, upper;
5003b5e53cdSSajid Ali   PetscBool     patchis_offproc = PETSC_TRUE;
501e30e807fSPeter Brune 
502e30e807fSPeter Brune   PetscFunctionBegin;
5039566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
5049566063dSJacob Faibussowitsch   if (iis) PetscCall(PetscMalloc1(n, iis));
5059566063dSJacob Faibussowitsch   if (ois) PetscCall(PetscMalloc1(n, ois));
506e30e807fSPeter Brune 
50790c77898SPeter Brune   for (i = 0; i < n; i++) {
5089566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(subdm[i], &subinfo));
50990c77898SPeter Brune     if (iis) {
510285ae305SPeter Brune       /* create the inner IS */
511285ae305SPeter Brune       lower.i = info.xs;
512285ae305SPeter Brune       lower.j = info.ys;
513285ae305SPeter Brune       lower.k = info.zs;
514285ae305SPeter Brune       upper.i = info.xs + info.xm;
515285ae305SPeter Brune       upper.j = info.ys + info.ym;
516285ae305SPeter Brune       upper.k = info.zs + info.zm;
5179566063dSJacob Faibussowitsch       PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &(*iis)[i], patchis_offproc));
51890c77898SPeter Brune     }
519e30e807fSPeter Brune 
52090c77898SPeter Brune     if (ois) {
521285ae305SPeter Brune       /* create the outer IS */
522285ae305SPeter Brune       lower.i = subinfo.xs;
523285ae305SPeter Brune       lower.j = subinfo.ys;
524285ae305SPeter Brune       lower.k = subinfo.zs;
525285ae305SPeter Brune       upper.i = subinfo.xs + subinfo.xm;
526285ae305SPeter Brune       upper.j = subinfo.ys + subinfo.ym;
527285ae305SPeter Brune       upper.k = subinfo.zs + subinfo.zm;
5289566063dSJacob Faibussowitsch       PetscCall(DMDACreatePatchIS(dm, &lower, &upper, &(*ois)[i], patchis_offproc));
52990c77898SPeter Brune     }
53090c77898SPeter Brune   }
531e30e807fSPeter Brune   PetscFunctionReturn(0);
532e30e807fSPeter Brune }
533e30e807fSPeter Brune 
5349371c9d4SSatish Balay PetscErrorCode DMCreateDomainDecomposition_DA(DM dm, PetscInt *len, char ***names, IS **iis, IS **ois, DM **subdm) {
53590c77898SPeter Brune   DM      *sdm;
53690c77898SPeter Brune   PetscInt n, i;
5370a696f66SPeter Brune 
5380adebc6cSBarry Smith   PetscFunctionBegin;
5399566063dSJacob Faibussowitsch   PetscCall(DMDASubDomainDA_Private(dm, &n, &sdm));
54090c77898SPeter Brune   if (names) {
5419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, names));
542ea78f98cSLisandro Dalcin     for (i = 0; i < n; i++) (*names)[i] = NULL;
543e30e807fSPeter Brune   }
5449566063dSJacob Faibussowitsch   PetscCall(DMDASubDomainIS_Private(dm, n, sdm, iis, ois));
54590c77898SPeter Brune   if (subdm) *subdm = sdm;
546e2c616c8SPeter Brune   else {
547*48a46eb9SPierre Jolivet     for (i = 0; i < n; i++) PetscCall(DMDestroy(&sdm[i]));
548e2c616c8SPeter Brune   }
54990c77898SPeter Brune   if (len) *len = n;
550e30e807fSPeter Brune   PetscFunctionReturn(0);
551e30e807fSPeter Brune }
552