xref: /petsc/src/dm/impls/da/dadd.c (revision 2b3cbbda7ef6bfb59aeed26278d0c91b4282b9fd)
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 
27*2b3cbbdaSStefano Zampini .seealso: `DMCreateDomainDecomposition()`, `DMCreateDomainDecompositionScatters()`
28063654ddSPeter Brune @*/
293b5e53cdSSajid Ali PetscErrorCode DMDACreatePatchIS(DM da,MatStencil *lower,MatStencil *upper,IS *is, PetscBool offproc)
30ff9846d9SPeter Brune {
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;
50063654ddSPeter Brune   M = dd->M; N = dd->N; P = dd->P;
51063654ddSPeter Brune   m = dd->m; n = dd->n; p = dd->p;
52063654ddSPeter Brune   dof = dd->w;
533b5e53cdSSajid Ali 
543b5e53cdSSajid Ali   nindices = -1;
553b5e53cdSSajid Ali   if (PetscLikely(upper->i - lower->i)) {
563b5e53cdSSajid Ali     nindices = nindices*(upper->i - lower->i);
573b5e53cdSSajid Ali     skip_i=PETSC_FALSE;
583b5e53cdSSajid Ali   }
593b5e53cdSSajid Ali   if (N>1) {
603b5e53cdSSajid Ali     valid_j = PETSC_TRUE;
613b5e53cdSSajid Ali     if (PetscLikely(upper->j - lower->j)) {
623b5e53cdSSajid Ali       nindices = nindices*(upper->j - lower->j);
633b5e53cdSSajid Ali       skip_j=PETSC_FALSE;
643b5e53cdSSajid Ali     }
653b5e53cdSSajid Ali   }
663b5e53cdSSajid Ali   if (P>1) {
673b5e53cdSSajid Ali     valid_k = PETSC_TRUE;
683b5e53cdSSajid Ali     if (PetscLikely(upper->k - lower->k)) {
693b5e53cdSSajid Ali       nindices = nindices*(upper->k - lower->k);
703b5e53cdSSajid Ali       skip_k=PETSC_FALSE;
713b5e53cdSSajid Ali     }
723b5e53cdSSajid Ali   }
733b5e53cdSSajid Ali   if (PetscLikely(nindices<0)) {
743b5e53cdSSajid Ali     if (PetscUnlikely(skip_i && skip_j && skip_k)) {
753b5e53cdSSajid Ali       nindices = 0;
763b5e53cdSSajid Ali     } else nindices = nindices*(-1);
773b5e53cdSSajid Ali   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Lower and Upper stencils are identical! Please check inputs.");
783b5e53cdSSajid Ali 
799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nindices*dof,&indices));
809566063dSJacob Faibussowitsch   PetscCall(DMDAGetOffset(da,&ox,&oy,&oz,NULL,NULL,NULL));
813b5e53cdSSajid Ali 
823b5e53cdSSajid Ali   if (!valid_k) {
833b5e53cdSSajid Ali     k = 0;
843b5e53cdSSajid Ali     upper->k=0;
853b5e53cdSSajid Ali     lower->k=0;
863b5e53cdSSajid Ali   }
873b5e53cdSSajid Ali   if (!valid_j) {
883b5e53cdSSajid Ali     j = 0;
893b5e53cdSSajid Ali     upper->j=0;
903b5e53cdSSajid Ali     lower->j=0;
913b5e53cdSSajid Ali   }
923b5e53cdSSajid Ali 
933b5e53cdSSajid Ali   if (offproc) {
949566063dSJacob Faibussowitsch     PetscCall(DMDAGetOwnershipRanges(da,&lx,&ly,&lz));
95063654ddSPeter Brune     /* start at index 0 on processor 0 */
96063654ddSPeter Brune     mr = 0;
97063654ddSPeter Brune     nr = 0;
98063654ddSPeter Brune     pr = 0;
99063654ddSPeter Brune     ms = 0;
100063654ddSPeter Brune     ns = 0;
101063654ddSPeter Brune     ps = 0;
102063654ddSPeter Brune     if (lx) me = lx[0];
103063654ddSPeter Brune     if (ly) ne = ly[0];
104063654ddSPeter Brune     if (lz) pe = lz[0];
1053b5e53cdSSajid Ali     /*
1063b5e53cdSSajid Ali        If no indices are to be returned, create an empty is,
1073b5e53cdSSajid Ali        this prevents hanging in while loops
1083b5e53cdSSajid Ali     */
1093b5e53cdSSajid Ali     if (skip_i && skip_j && skip_k) goto createis;
1103b5e53cdSSajid Ali     /*
1113b5e53cdSSajid Ali        do..while loops to ensure the block gets entered once,
1123b5e53cdSSajid Ali        regardless of control condition being met, necessary for
1133b5e53cdSSajid Ali        cases when a subset of skip_i/j/k is true
1143b5e53cdSSajid Ali     */
1153b5e53cdSSajid Ali     if (skip_k) k = upper->k-oz; else k = lower->k-oz;
1163b5e53cdSSajid Ali     do {
1173b5e53cdSSajid Ali       if (skip_j) j = upper->j-oy; else j = lower->j-oy;
1183b5e53cdSSajid Ali       do {
1193b5e53cdSSajid Ali         if (skip_i) i = upper->i-ox; else i = lower->i-ox;
1203b5e53cdSSajid Ali         do {
121063654ddSPeter Brune           /* "actual" indices rather than ones outside of the domain */
122063654ddSPeter Brune           ii = i;
123063654ddSPeter Brune           jj = j;
124063654ddSPeter Brune           kk = k;
125063654ddSPeter Brune           if (ii < 0) ii = M + ii;
126063654ddSPeter Brune           if (jj < 0) jj = N + jj;
127063654ddSPeter Brune           if (kk < 0) kk = P + kk;
128063654ddSPeter Brune           if (ii > M-1) ii = ii - M;
129063654ddSPeter Brune           if (jj > N-1) jj = jj - N;
130063654ddSPeter Brune           if (kk > P-1) kk = kk - P;
131063654ddSPeter Brune           /* gone out of processor range on x axis */
132063654ddSPeter Brune           while (ii > me-1 || ii < ms) {
133063654ddSPeter Brune             if (mr == m-1) {
134063654ddSPeter Brune               ms = 0;
135063654ddSPeter Brune               me = lx[0];
136063654ddSPeter Brune               mr = 0;
137063654ddSPeter Brune             } else {
138063654ddSPeter Brune               mr++;
139063654ddSPeter Brune               ms = me;
140063654ddSPeter Brune               me += lx[mr];
141063654ddSPeter Brune             }
142063654ddSPeter Brune           }
143063654ddSPeter Brune           /* gone out of processor range on y axis */
144063654ddSPeter Brune           while (jj > ne-1 || jj < ns) {
145063654ddSPeter Brune             if (nr == n-1) {
146063654ddSPeter Brune               ns = 0;
147063654ddSPeter Brune               ne = ly[0];
148063654ddSPeter Brune               nr = 0;
149063654ddSPeter Brune             } else {
150063654ddSPeter Brune               nr++;
151063654ddSPeter Brune               ns = ne;
152063654ddSPeter Brune               ne += ly[nr];
153063654ddSPeter Brune             }
154063654ddSPeter Brune           }
155063654ddSPeter Brune           /* gone out of processor range on z axis */
156063654ddSPeter Brune           while (kk > pe-1 || kk < ps) {
157063654ddSPeter Brune             if (pr == p-1) {
158063654ddSPeter Brune               ps = 0;
159063654ddSPeter Brune               pe = lz[0];
160063654ddSPeter Brune               pr = 0;
161063654ddSPeter Brune             } else {
162063654ddSPeter Brune               pr++;
163063654ddSPeter Brune               ps = pe;
164063654ddSPeter Brune               pe += lz[pr];
165063654ddSPeter Brune             }
166063654ddSPeter Brune           }
167063654ddSPeter Brune           /* compute the vector base on owning processor */
168063654ddSPeter Brune           xm = me - ms;
169063654ddSPeter Brune           ym = ne - ns;
170063654ddSPeter Brune           zm = pe - ps;
171b0bff951SPeter Brune           base = ms*ym*zm + ns*M*zm + ps*M*N;
172063654ddSPeter Brune           /* compute the local coordinates on owning processor */
173063654ddSPeter Brune           si = ii - ms;
174063654ddSPeter Brune           sj = jj - ns;
175063654ddSPeter Brune           sk = kk - ps;
176063654ddSPeter Brune           for (l=0;l<dof;l++) {
177063654ddSPeter Brune             indices[idx] = l + dof*(base + si + xm*sj + xm*ym*sk);
178ff9846d9SPeter Brune             idx++;
179ff9846d9SPeter Brune           }
1803b5e53cdSSajid Ali           i++;
1813b5e53cdSSajid Ali         } while (i<upper->i-ox);
1823b5e53cdSSajid Ali         j++;
1833b5e53cdSSajid Ali       } while (j<upper->j-oy);
1843b5e53cdSSajid Ali       k++;
1853b5e53cdSSajid Ali     } while (k<upper->k-oz);
1863b5e53cdSSajid Ali   }
1873b5e53cdSSajid Ali 
1883b5e53cdSSajid Ali   if (!offproc) {
1899566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(da, &ms, &ns, &ps, &mw, &nw, &pw));
1903b5e53cdSSajid Ali     me = ms + mw;
1913b5e53cdSSajid Ali     if (N>1) ne = ns + nw;
1923b5e53cdSSajid Ali     if (P>1) pe = ps + pw;
1933b5e53cdSSajid Ali     /* Account for DM offsets */
1943b5e53cdSSajid Ali     ms = ms - ox; me = me - ox;
1953b5e53cdSSajid Ali     ns = ns - oy; ne = ne - oy;
1963b5e53cdSSajid Ali     ps = ps - oz; pe = pe - oz;
1973b5e53cdSSajid Ali 
1983b5e53cdSSajid Ali     /* compute the vector base on owning processor */
1993b5e53cdSSajid Ali     xm = me - ms;
2003b5e53cdSSajid Ali     ym = ne - ns;
2013b5e53cdSSajid Ali     zm = pe - ps;
2023b5e53cdSSajid Ali     base = ms*ym*zm + ns*M*zm + ps*M*N;
2033b5e53cdSSajid Ali     /*
2043b5e53cdSSajid Ali        if no indices are to be returned, create an empty is,
2053b5e53cdSSajid Ali        this prevents hanging in while loops
2063b5e53cdSSajid Ali     */
2073b5e53cdSSajid Ali     if (skip_i && skip_j && skip_k) goto createis;
2083b5e53cdSSajid Ali     /*
2093b5e53cdSSajid Ali        do..while loops to ensure the block gets entered once,
2103b5e53cdSSajid Ali        regardless of control condition being met, necessary for
2113b5e53cdSSajid Ali        cases when a subset of skip_i/j/k is true
2123b5e53cdSSajid Ali     */
2133b5e53cdSSajid Ali     if (skip_k) k = upper->k-oz; else k = lower->k-oz;
2143b5e53cdSSajid Ali     do {
2153b5e53cdSSajid Ali       if (skip_j) j = upper->j-oy; else j = lower->j-oy;
2163b5e53cdSSajid Ali       do {
2173b5e53cdSSajid Ali         if (skip_i) i = upper->i-ox; else i = lower->i-ox;
2183b5e53cdSSajid Ali         do {
2193b5e53cdSSajid Ali           if (k>=ps && k<=pe-1) {
2203b5e53cdSSajid Ali             if (j>=ns && j<=ne-1) {
2213b5e53cdSSajid Ali               if (i>=ms && i<=me-1) {
2223b5e53cdSSajid Ali                 /* compute the local coordinates on owning processor */
2233b5e53cdSSajid Ali                 si = i - ms;
2243b5e53cdSSajid Ali                 sj = j - ns;
2253b5e53cdSSajid Ali                 sk = k - ps;
2263b5e53cdSSajid Ali                 for (l=0; l<dof; l++) {
2273b5e53cdSSajid Ali                   indices[idx] = l + dof*(base + si + xm*sj + xm*ym*sk);
2283b5e53cdSSajid Ali                   idx++;
229ff9846d9SPeter Brune                 }
230ff9846d9SPeter Brune               }
231063654ddSPeter Brune             }
2323b5e53cdSSajid Ali           }
2333b5e53cdSSajid Ali           i++;
2343b5e53cdSSajid Ali         } while (i<upper->i-ox);
2353b5e53cdSSajid Ali         j++;
2363b5e53cdSSajid Ali       } while (j<upper->j-oy);
2373b5e53cdSSajid Ali       k++;
2383b5e53cdSSajid Ali     } while (k<upper->k-oz);
2393b5e53cdSSajid Ali 
2409566063dSJacob Faibussowitsch     PetscCall(PetscRealloc((size_t)(idx*sizeof(PetscInt)), (void*)&indices));
2413b5e53cdSSajid Ali   }
2423b5e53cdSSajid Ali 
2433b5e53cdSSajid Ali   createis:
2449566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)da),idx,indices,PETSC_OWN_POINTER,is));
245ff9846d9SPeter Brune   PetscFunctionReturn(0);
246ff9846d9SPeter Brune }
247ff9846d9SPeter Brune 
24890c77898SPeter Brune PetscErrorCode DMDASubDomainDA_Private(DM dm, PetscInt *nlocal, DM **sdm)
2490adebc6cSBarry Smith {
25090c77898SPeter Brune   DM             *da;
25190c77898SPeter Brune   PetscInt       dim,size,i,j,k,idx;
252e30e807fSPeter Brune   DMDALocalInfo  info;
253ff9846d9SPeter Brune   PetscInt       xsize,ysize,zsize;
254e30e807fSPeter Brune   PetscInt       xo,yo,zo;
25590c77898SPeter Brune   PetscInt       xs,ys,zs;
25690c77898SPeter Brune   PetscInt       xm=1,ym=1,zm=1;
2577ddda789SPeter Brune   PetscInt       xol,yol,zol;
25890c77898SPeter Brune   PetscInt       m=1,n=1,p=1;
25990c77898SPeter Brune   PetscInt       M,N,P;
26090c77898SPeter Brune   PetscInt       pm,mtmp;
261e30e807fSPeter Brune 
262e30e807fSPeter Brune   PetscFunctionBegin;
2639566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
2649566063dSJacob Faibussowitsch   PetscCall(DMDAGetOverlap(dm,&xol,&yol,&zol));
2659566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumLocalSubDomains(dm,&size));
2669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&da));
2677ddda789SPeter Brune 
26890c77898SPeter Brune   if (nlocal) *nlocal = size;
269e30e807fSPeter Brune 
27090c77898SPeter Brune   dim = info.dim;
271e30e807fSPeter Brune 
27290c77898SPeter Brune   M = info.xm;
27390c77898SPeter Brune   N = info.ym;
27490c77898SPeter Brune   P = info.zm;
27590c77898SPeter Brune 
27690c77898SPeter Brune   if (dim == 1) {
27790c77898SPeter Brune     m = size;
27890c77898SPeter Brune   } else if (dim == 2) {
27990c77898SPeter Brune     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N)));
28090c77898SPeter Brune     while (m > 0) {
28190c77898SPeter Brune       n = size/m;
28290c77898SPeter Brune       if (m*n*p == size) break;
28390c77898SPeter Brune       m--;
28490c77898SPeter Brune     }
28590c77898SPeter Brune   } else if (dim == 3) {
28690c77898SPeter Brune     n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); if (!n) n = 1;
28790c77898SPeter Brune     while (n > 0) {
28890c77898SPeter Brune       pm = size/n;
28990c77898SPeter Brune       if (n*pm == size) break;
29090c77898SPeter Brune       n--;
29190c77898SPeter Brune     }
29290c77898SPeter Brune     if (!n) n = 1;
29390c77898SPeter Brune     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
29490c77898SPeter Brune     if (!m) m = 1;
29590c77898SPeter Brune     while (m > 0) {
29690c77898SPeter Brune       p = size/(m*n);
29790c77898SPeter Brune       if (m*n*p == size) break;
29890c77898SPeter Brune       m--;
29990c77898SPeter Brune     }
30090c77898SPeter Brune     if (M > P && m < p) {mtmp = m; m = p; p = mtmp;}
30190c77898SPeter Brune   }
30290c77898SPeter Brune 
30390c77898SPeter Brune   zs = info.zs;
30490c77898SPeter Brune   idx = 0;
30590c77898SPeter Brune   for (k = 0; k < p; k++) {
30690c77898SPeter Brune     ys = info.ys;
30790c77898SPeter Brune     for (j = 0; j < n; j++) {
30890c77898SPeter Brune       xs = info.xs;
30990c77898SPeter Brune       for (i = 0; i < m; i++) {
31090c77898SPeter Brune         if (dim == 1) {
31190c77898SPeter Brune           xm = M/m + ((M % m) > i);
31290c77898SPeter Brune         } else if (dim == 2) {
31390c77898SPeter Brune           xm = M/m + ((M % m) > i);
31490c77898SPeter Brune           ym = N/n + ((N % n) > j);
31590c77898SPeter Brune         } else if (dim == 3) {
31690c77898SPeter Brune           xm = M/m + ((M % m) > i);
31790c77898SPeter Brune           ym = N/n + ((N % n) > j);
31890c77898SPeter Brune           zm = P/p + ((P % p) > k);
31990c77898SPeter Brune         }
32090c77898SPeter Brune 
32190c77898SPeter Brune         xsize = xm;
32290c77898SPeter Brune         ysize = ym;
32390c77898SPeter Brune         zsize = zm;
32490c77898SPeter Brune         xo = xs;
32590c77898SPeter Brune         yo = ys;
32690c77898SPeter Brune         zo = zs;
32790c77898SPeter Brune 
3289566063dSJacob Faibussowitsch         PetscCall(DMDACreate(PETSC_COMM_SELF,&(da[idx])));
3299566063dSJacob Faibussowitsch         PetscCall(DMSetOptionsPrefix(da[idx],"sub_"));
3309566063dSJacob Faibussowitsch         PetscCall(DMSetDimension(da[idx], info.dim));
3319566063dSJacob Faibussowitsch         PetscCall(DMDASetDof(da[idx], info.dof));
33290c77898SPeter Brune 
3339566063dSJacob Faibussowitsch         PetscCall(DMDASetStencilType(da[idx],info.st));
3349566063dSJacob Faibussowitsch         PetscCall(DMDASetStencilWidth(da[idx],info.sw));
33590c77898SPeter Brune 
336bff4a2f0SMatthew G. Knepley         if (info.bx == DM_BOUNDARY_PERIODIC || (xs != 0)) {
3377ddda789SPeter Brune           xsize += xol;
3387ddda789SPeter Brune           xo    -= xol;
339e30e807fSPeter Brune         }
340bff4a2f0SMatthew G. Knepley         if (info.by == DM_BOUNDARY_PERIODIC || (ys != 0)) {
3417ddda789SPeter Brune           ysize += yol;
3427ddda789SPeter Brune           yo    -= yol;
343e30e807fSPeter Brune         }
344bff4a2f0SMatthew G. Knepley         if (info.bz == DM_BOUNDARY_PERIODIC || (zs != 0)) {
3457ddda789SPeter Brune           zsize += zol;
3467ddda789SPeter Brune           zo    -= zol;
347e30e807fSPeter Brune         }
348e30e807fSPeter Brune 
349bff4a2f0SMatthew G. Knepley         if (info.bx == DM_BOUNDARY_PERIODIC || (xs+xm != info.mx)) xsize += xol;
350bff4a2f0SMatthew G. Knepley         if (info.by == DM_BOUNDARY_PERIODIC || (ys+ym != info.my)) ysize += yol;
351bff4a2f0SMatthew G. Knepley         if (info.bz == DM_BOUNDARY_PERIODIC || (zs+zm != info.mz)) zsize += zol;
3528865f1eaSKarl Rupp 
353bff4a2f0SMatthew G. Knepley         if (info.bx != DM_BOUNDARY_PERIODIC) {
35490c77898SPeter Brune           if (xo < 0) {
35590c77898SPeter Brune             xsize += xo;
35690c77898SPeter Brune             xo = 0;
35790c77898SPeter Brune           }
35890c77898SPeter Brune           if (xo+xsize > info.mx-1) {
35990c77898SPeter Brune             xsize -= xo+xsize - info.mx;
36090c77898SPeter Brune           }
36190c77898SPeter Brune         }
362bff4a2f0SMatthew G. Knepley         if (info.by != DM_BOUNDARY_PERIODIC) {
36390c77898SPeter Brune           if (yo < 0) {
36490c77898SPeter Brune             ysize += yo;
36590c77898SPeter Brune             yo = 0;
36690c77898SPeter Brune           }
36790c77898SPeter Brune           if (yo+ysize > info.my-1) {
36890c77898SPeter Brune             ysize -= yo+ysize - info.my;
36990c77898SPeter Brune           }
37090c77898SPeter Brune         }
371bff4a2f0SMatthew G. Knepley         if (info.bz != DM_BOUNDARY_PERIODIC) {
37290c77898SPeter Brune           if (zo < 0) {
37390c77898SPeter Brune             zsize += zo;
37490c77898SPeter Brune             zo = 0;
37590c77898SPeter Brune           }
37690c77898SPeter Brune           if (zo+zsize > info.mz-1) {
37790c77898SPeter Brune             zsize -= zo+zsize - info.mz;
37890c77898SPeter Brune           }
37990c77898SPeter Brune         }
38090c77898SPeter Brune 
3819566063dSJacob Faibussowitsch         PetscCall(DMDASetSizes(da[idx], xsize, ysize, zsize));
3829566063dSJacob Faibussowitsch         PetscCall(DMDASetNumProcs(da[idx], 1, 1, 1));
3839566063dSJacob Faibussowitsch         PetscCall(DMDASetBoundaryType(da[idx], DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED));
384e30e807fSPeter Brune 
385e30e807fSPeter Brune         /* set up as a block instead */
3869566063dSJacob Faibussowitsch         PetscCall(DMSetUp(da[idx]));
387e30e807fSPeter Brune 
38840234c92SPeter Brune         /* nonoverlapping region */
3899566063dSJacob Faibussowitsch         PetscCall(DMDASetNonOverlappingRegion(da[idx],xs,ys,zs,xm,ym,zm));
39040234c92SPeter Brune 
39195c13181SPeter Brune         /* this alters the behavior of DMDAGetInfo, DMDAGetLocalInfo, DMDAGetCorners, and DMDAGetGhostedCorners and should be used with care */
3929566063dSJacob Faibussowitsch         PetscCall(DMDASetOffset(da[idx],xo,yo,zo,info.mx,info.my,info.mz));
39390c77898SPeter Brune         xs += xm;
39490c77898SPeter Brune         idx++;
39590c77898SPeter Brune       }
39690c77898SPeter Brune       ys += ym;
39790c77898SPeter Brune     }
39890c77898SPeter Brune     zs += zm;
39990c77898SPeter Brune   }
40090c77898SPeter Brune   *sdm = da;
401e30e807fSPeter Brune   PetscFunctionReturn(0);
402e30e807fSPeter Brune }
403e30e807fSPeter Brune 
404e30e807fSPeter Brune /*
405e30e807fSPeter Brune    Fills the local vector problem on the subdomain from the global problem.
406e30e807fSPeter Brune 
407285ae305SPeter Brune    Right now this assumes one subdomain per processor.
408285ae305SPeter Brune 
409e30e807fSPeter Brune */
4100adebc6cSBarry Smith PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat)
4110adebc6cSBarry Smith {
412285ae305SPeter Brune   DMDALocalInfo  info,subinfo;
413e30e807fSPeter Brune   DM             subdm;
414285ae305SPeter Brune   MatStencil     upper,lower;
415285ae305SPeter Brune   IS             idis,isis,odis,osis,gdis;
416285ae305SPeter Brune   Vec            svec,dvec,slvec;
41740234c92SPeter Brune   PetscInt       xm,ym,zm,xs,ys,zs;
41890c77898SPeter Brune   PetscInt       i;
4193b5e53cdSSajid Ali   PetscBool      patchis_offproc = PETSC_TRUE;
420e30e807fSPeter Brune 
421e30e807fSPeter Brune   PetscFunctionBegin;
422e30e807fSPeter Brune   /* allocate the arrays of scatters */
4239566063dSJacob Faibussowitsch   if (iscat) PetscCall(PetscMalloc1(nsubdms,iscat));
4249566063dSJacob Faibussowitsch   if (oscat) PetscCall(PetscMalloc1(nsubdms,oscat));
4259566063dSJacob Faibussowitsch   if (lscat) PetscCall(PetscMalloc1(nsubdms,lscat));
426e30e807fSPeter Brune 
4279566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
42890c77898SPeter Brune   for (i = 0; i < nsubdms; i++) {
42990c77898SPeter Brune     subdm = subdms[i];
4309566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(subdm,&subinfo));
4319566063dSJacob Faibussowitsch     PetscCall(DMDAGetNonOverlappingRegion(subdm,&xs,&ys,&zs,&xm,&ym,&zm));
432e30e807fSPeter Brune 
433285ae305SPeter Brune     /* create the global and subdomain index sets for the inner domain */
43440234c92SPeter Brune     lower.i = xs;
43540234c92SPeter Brune     lower.j = ys;
43640234c92SPeter Brune     lower.k = zs;
43740234c92SPeter Brune     upper.i = xs+xm;
43840234c92SPeter Brune     upper.j = ys+ym;
43940234c92SPeter Brune     upper.k = zs+zm;
4409566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(dm,&lower,&upper,&idis,patchis_offproc));
4419566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(subdm,&lower,&upper,&isis,patchis_offproc));
442e30e807fSPeter Brune 
443285ae305SPeter Brune     /* create the global and subdomain index sets for the outer subdomain */
444285ae305SPeter Brune     lower.i = subinfo.xs;
445285ae305SPeter Brune     lower.j = subinfo.ys;
446285ae305SPeter Brune     lower.k = subinfo.zs;
447285ae305SPeter Brune     upper.i = subinfo.xs+subinfo.xm;
448285ae305SPeter Brune     upper.j = subinfo.ys+subinfo.ym;
449285ae305SPeter Brune     upper.k = subinfo.zs+subinfo.zm;
4509566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(dm,&lower,&upper,&odis,patchis_offproc));
4519566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(subdm,&lower,&upper,&osis,patchis_offproc));
452e30e807fSPeter Brune 
453285ae305SPeter Brune     /* global and subdomain ISes for the local indices of the subdomain */
454285ae305SPeter Brune     /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */
455285ae305SPeter Brune     lower.i = subinfo.gxs;
456285ae305SPeter Brune     lower.j = subinfo.gys;
457285ae305SPeter Brune     lower.k = subinfo.gzs;
458285ae305SPeter Brune     upper.i = subinfo.gxs+subinfo.gxm;
459285ae305SPeter Brune     upper.j = subinfo.gys+subinfo.gym;
460285ae305SPeter Brune     upper.k = subinfo.gzs+subinfo.gzm;
4619566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(dm,&lower,&upper,&gdis,patchis_offproc));
462e30e807fSPeter Brune 
463e30e807fSPeter Brune     /* form the scatter */
4649566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm,&dvec));
4659566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(subdm,&svec));
4669566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(subdm,&slvec));
467e30e807fSPeter Brune 
4689566063dSJacob Faibussowitsch     if (iscat) PetscCall(VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[i]));
4699566063dSJacob Faibussowitsch     if (oscat) PetscCall(VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[i]));
4709566063dSJacob Faibussowitsch     if (lscat) PetscCall(VecScatterCreate(dvec,gdis,slvec,NULL,&(*lscat)[i]));
471e30e807fSPeter Brune 
4729566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm,&dvec));
4739566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(subdm,&svec));
4749566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(subdm,&slvec));
475e30e807fSPeter Brune 
4769566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&idis));
4779566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isis));
478e30e807fSPeter Brune 
4799566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&odis));
4809566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&osis));
481e30e807fSPeter Brune 
4829566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&gdis));
48390c77898SPeter Brune   }
484e30e807fSPeter Brune   PetscFunctionReturn(0);
485e30e807fSPeter Brune }
486e30e807fSPeter Brune 
48790c77898SPeter Brune PetscErrorCode DMDASubDomainIS_Private(DM dm,PetscInt n,DM *subdm,IS **iis,IS **ois)
4880adebc6cSBarry Smith {
48990c77898SPeter Brune   PetscInt       i;
490e30e807fSPeter Brune   DMDALocalInfo  info,subinfo;
491285ae305SPeter Brune   MatStencil     lower,upper;
4923b5e53cdSSajid Ali   PetscBool      patchis_offproc = PETSC_TRUE;
493e30e807fSPeter Brune 
494e30e807fSPeter Brune   PetscFunctionBegin;
4959566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm,&info));
4969566063dSJacob Faibussowitsch   if (iis) PetscCall(PetscMalloc1(n,iis));
4979566063dSJacob Faibussowitsch   if (ois) PetscCall(PetscMalloc1(n,ois));
498e30e807fSPeter Brune 
49990c77898SPeter Brune   for (i = 0;i < n; i++) {
5009566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(subdm[i],&subinfo));
50190c77898SPeter Brune     if (iis) {
502285ae305SPeter Brune       /* create the inner IS */
503285ae305SPeter Brune       lower.i = info.xs;
504285ae305SPeter Brune       lower.j = info.ys;
505285ae305SPeter Brune       lower.k = info.zs;
506285ae305SPeter Brune       upper.i = info.xs+info.xm;
507285ae305SPeter Brune       upper.j = info.ys+info.ym;
508285ae305SPeter Brune       upper.k = info.zs+info.zm;
5099566063dSJacob Faibussowitsch       PetscCall(DMDACreatePatchIS(dm,&lower,&upper,&(*iis)[i],patchis_offproc));
51090c77898SPeter Brune     }
511e30e807fSPeter Brune 
51290c77898SPeter Brune     if (ois) {
513285ae305SPeter Brune       /* create the outer IS */
514285ae305SPeter Brune       lower.i = subinfo.xs;
515285ae305SPeter Brune       lower.j = subinfo.ys;
516285ae305SPeter Brune       lower.k = subinfo.zs;
517285ae305SPeter Brune       upper.i = subinfo.xs+subinfo.xm;
518285ae305SPeter Brune       upper.j = subinfo.ys+subinfo.ym;
519285ae305SPeter Brune       upper.k = subinfo.zs+subinfo.zm;
5209566063dSJacob Faibussowitsch       PetscCall(DMDACreatePatchIS(dm,&lower,&upper,&(*ois)[i],patchis_offproc));
52190c77898SPeter Brune     }
52290c77898SPeter Brune   }
523e30e807fSPeter Brune   PetscFunctionReturn(0);
524e30e807fSPeter Brune }
525e30e807fSPeter Brune 
5260adebc6cSBarry Smith PetscErrorCode DMCreateDomainDecomposition_DA(DM dm,PetscInt *len,char ***names,IS **iis,IS **ois,DM **subdm)
5270adebc6cSBarry Smith {
52890c77898SPeter Brune   DM             *sdm;
52990c77898SPeter Brune   PetscInt       n,i;
5300a696f66SPeter Brune 
5310adebc6cSBarry Smith   PetscFunctionBegin;
5329566063dSJacob Faibussowitsch   PetscCall(DMDASubDomainDA_Private(dm,&n,&sdm));
53390c77898SPeter Brune   if (names) {
5349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,names));
535ea78f98cSLisandro Dalcin     for (i=0;i<n;i++) (*names)[i] = NULL;
536e30e807fSPeter Brune   }
5379566063dSJacob Faibussowitsch   PetscCall(DMDASubDomainIS_Private(dm,n,sdm,iis,ois));
53890c77898SPeter Brune   if (subdm) *subdm = sdm;
539e2c616c8SPeter Brune   else {
540e2c616c8SPeter Brune     for (i=0;i<n;i++) {
5419566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&sdm[i]));
542e2c616c8SPeter Brune     }
543e2c616c8SPeter Brune   }
54490c77898SPeter Brune   if (len) *len = n;
545e30e807fSPeter Brune   PetscFunctionReturn(0);
546e30e807fSPeter Brune }
547