xref: /petsc/src/dm/impls/da/dadd.c (revision 3b5e53cdb169e94712b800a2d1c9e10766bebf6e)
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 
6*3b5e53cdSSajid 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
11*3b5e53cdSSajid Ali .  upper - a matstencil with i, j and k corresponding to the upper corner of the patch
12*3b5e53cdSSajid 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 
19*3b5e53cdSSajid Ali Notes:
20*3b5e53cdSSajid Ali This routine always returns an IS on the DMDA's comm, if offproc is set to PETSC_TRUE,
21*3b5e53cdSSajid Ali the routine returns an IS with all the indices requested regardless of whether these indices
22*3b5e53cdSSajid Ali are present on the requesting rank or not. Thus, it is upon the caller to ensure that
23*3b5e53cdSSajid Ali the indices returned in this mode are appropriate. If offproc is set to PETSC_FALSE,
24*3b5e53cdSSajid Ali the IS only returns the subset of indices that are present on the requesting rank and there
25*3b5e53cdSSajid Ali is no duplication of indices.
26*3b5e53cdSSajid Ali 
27063654ddSPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDACreateDomainDecompositionScatters()
28063654ddSPeter Brune @*/
29*3b5e53cdSSajid Ali PetscErrorCode DMDACreatePatchIS(DM da,MatStencil *lower,MatStencil *upper,IS *is, PetscBool offproc)
30ff9846d9SPeter Brune {
31063654ddSPeter Brune   PetscInt       ms=0,ns=0,ps=0;
32*3b5e53cdSSajid 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;
37*3b5e53cdSSajid 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;
42*3b5e53cdSSajid Ali   const PetscInt *lx,*ly,*lz;
43063654ddSPeter Brune   PetscInt       nindices;
44063654ddSPeter Brune   PetscInt       *indices;
45063654ddSPeter Brune   DM_DA          *dd = (DM_DA*)da->data;
46*3b5e53cdSSajid Ali   PetscBool      skip_i=PETSC_TRUE, skip_j=PETSC_TRUE, skip_k=PETSC_TRUE;
47*3b5e53cdSSajid Ali   PetscBool      valid_j=PETSC_FALSE, valid_k=PETSC_FALSE; /* DMDA has at least 1 dimension */
48063654ddSPeter Brune   PetscErrorCode ierr;
49ff9846d9SPeter Brune 
500adebc6cSBarry Smith   PetscFunctionBegin;
51063654ddSPeter Brune   M = dd->M; N = dd->N; P = dd->P;
52063654ddSPeter Brune   m = dd->m; n = dd->n; p = dd->p;
53063654ddSPeter Brune   dof = dd->w;
54*3b5e53cdSSajid Ali 
55*3b5e53cdSSajid Ali   nindices = -1;
56*3b5e53cdSSajid Ali   if (PetscLikely(upper->i - lower->i)) {
57*3b5e53cdSSajid Ali     nindices = nindices*(upper->i - lower->i);
58*3b5e53cdSSajid Ali     skip_i=PETSC_FALSE;
59*3b5e53cdSSajid Ali   }
60*3b5e53cdSSajid Ali   if (N>1) {
61*3b5e53cdSSajid Ali     valid_j = PETSC_TRUE;
62*3b5e53cdSSajid Ali     if (PetscLikely(upper->j - lower->j)) {
63*3b5e53cdSSajid Ali       nindices = nindices*(upper->j - lower->j);
64*3b5e53cdSSajid Ali       skip_j=PETSC_FALSE;
65*3b5e53cdSSajid Ali     }
66*3b5e53cdSSajid Ali   }
67*3b5e53cdSSajid Ali   if (P>1) {
68*3b5e53cdSSajid Ali     valid_k = PETSC_TRUE;
69*3b5e53cdSSajid Ali     if (PetscLikely(upper->k - lower->k)) {
70*3b5e53cdSSajid Ali       nindices = nindices*(upper->k - lower->k);
71*3b5e53cdSSajid Ali       skip_k=PETSC_FALSE;
72*3b5e53cdSSajid Ali     }
73*3b5e53cdSSajid Ali   }
74*3b5e53cdSSajid Ali   if (PetscLikely(nindices<0)) {
75*3b5e53cdSSajid Ali     if (PetscUnlikely(skip_i && skip_j && skip_k)) {
76*3b5e53cdSSajid Ali       nindices = 0;
77*3b5e53cdSSajid Ali     } else nindices = nindices*(-1);
78*3b5e53cdSSajid Ali   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Lower and Upper stencils are identical! Please check inputs.");
79*3b5e53cdSSajid Ali 
80*3b5e53cdSSajid Ali   ierr = PetscMalloc1(nindices*dof,&indices);CHKERRQ(ierr);
810298fd71SBarry Smith   ierr = DMDAGetOffset(da,&ox,&oy,&oz,NULL,NULL,NULL);CHKERRQ(ierr);
82*3b5e53cdSSajid Ali 
83*3b5e53cdSSajid Ali   if (!valid_k) {
84*3b5e53cdSSajid Ali     k = 0;
85*3b5e53cdSSajid Ali     upper->k=0;
86*3b5e53cdSSajid Ali     lower->k=0;
87*3b5e53cdSSajid Ali   }
88*3b5e53cdSSajid Ali   if (!valid_j) {
89*3b5e53cdSSajid Ali     j = 0;
90*3b5e53cdSSajid Ali     upper->j=0;
91*3b5e53cdSSajid Ali     lower->j=0;
92*3b5e53cdSSajid Ali   }
93*3b5e53cdSSajid Ali 
94*3b5e53cdSSajid Ali   if (offproc) {
95063654ddSPeter Brune     ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr);
96063654ddSPeter Brune     /* start at index 0 on processor 0 */
97063654ddSPeter Brune     mr = 0;
98063654ddSPeter Brune     nr = 0;
99063654ddSPeter Brune     pr = 0;
100063654ddSPeter Brune     ms = 0;
101063654ddSPeter Brune     ns = 0;
102063654ddSPeter Brune     ps = 0;
103063654ddSPeter Brune     if (lx) me = lx[0];
104063654ddSPeter Brune     if (ly) ne = ly[0];
105063654ddSPeter Brune     if (lz) pe = lz[0];
106*3b5e53cdSSajid Ali     /*
107*3b5e53cdSSajid Ali        If no indices are to be returned, create an empty is,
108*3b5e53cdSSajid Ali        this prevents hanging in while loops
109*3b5e53cdSSajid Ali     */
110*3b5e53cdSSajid Ali     if (skip_i && skip_j && skip_k) goto createis;
111*3b5e53cdSSajid Ali     /*
112*3b5e53cdSSajid Ali        do..while loops to ensure the block gets entered once,
113*3b5e53cdSSajid Ali        regardless of control condition being met, necessary for
114*3b5e53cdSSajid Ali        cases when a subset of skip_i/j/k is true
115*3b5e53cdSSajid Ali     */
116*3b5e53cdSSajid Ali     if (skip_k) k = upper->k-oz; else k = lower->k-oz;
117*3b5e53cdSSajid Ali     do {
118*3b5e53cdSSajid Ali       if (skip_j) j = upper->j-oy; else j = lower->j-oy;
119*3b5e53cdSSajid Ali       do {
120*3b5e53cdSSajid Ali         if (skip_i) i = upper->i-ox; else i = lower->i-ox;
121*3b5e53cdSSajid Ali         do {
122063654ddSPeter Brune           /* "actual" indices rather than ones outside of the domain */
123063654ddSPeter Brune           ii = i;
124063654ddSPeter Brune           jj = j;
125063654ddSPeter Brune           kk = k;
126063654ddSPeter Brune           if (ii < 0) ii = M + ii;
127063654ddSPeter Brune           if (jj < 0) jj = N + jj;
128063654ddSPeter Brune           if (kk < 0) kk = P + kk;
129063654ddSPeter Brune           if (ii > M-1) ii = ii - M;
130063654ddSPeter Brune           if (jj > N-1) jj = jj - N;
131063654ddSPeter Brune           if (kk > P-1) kk = kk - P;
132063654ddSPeter Brune           /* gone out of processor range on x axis */
133063654ddSPeter Brune           while (ii > me-1 || ii < ms) {
134063654ddSPeter Brune             if (mr == m-1) {
135063654ddSPeter Brune               ms = 0;
136063654ddSPeter Brune               me = lx[0];
137063654ddSPeter Brune               mr = 0;
138063654ddSPeter Brune             } else {
139063654ddSPeter Brune               mr++;
140063654ddSPeter Brune               ms = me;
141063654ddSPeter Brune               me += lx[mr];
142063654ddSPeter Brune             }
143063654ddSPeter Brune           }
144063654ddSPeter Brune           /* gone out of processor range on y axis */
145063654ddSPeter Brune           while (jj > ne-1 || jj < ns) {
146063654ddSPeter Brune             if (nr == n-1) {
147063654ddSPeter Brune               ns = 0;
148063654ddSPeter Brune               ne = ly[0];
149063654ddSPeter Brune               nr = 0;
150063654ddSPeter Brune             } else {
151063654ddSPeter Brune               nr++;
152063654ddSPeter Brune               ns = ne;
153063654ddSPeter Brune               ne += ly[nr];
154063654ddSPeter Brune             }
155063654ddSPeter Brune           }
156063654ddSPeter Brune           /* gone out of processor range on z axis */
157063654ddSPeter Brune           while (kk > pe-1 || kk < ps) {
158063654ddSPeter Brune             if (pr == p-1) {
159063654ddSPeter Brune               ps = 0;
160063654ddSPeter Brune               pe = lz[0];
161063654ddSPeter Brune               pr = 0;
162063654ddSPeter Brune             } else {
163063654ddSPeter Brune               pr++;
164063654ddSPeter Brune               ps = pe;
165063654ddSPeter Brune               pe += lz[pr];
166063654ddSPeter Brune             }
167063654ddSPeter Brune           }
168063654ddSPeter Brune           /* compute the vector base on owning processor */
169063654ddSPeter Brune           xm = me - ms;
170063654ddSPeter Brune           ym = ne - ns;
171063654ddSPeter Brune           zm = pe - ps;
172b0bff951SPeter Brune           base = ms*ym*zm + ns*M*zm + ps*M*N;
173063654ddSPeter Brune           /* compute the local coordinates on owning processor */
174063654ddSPeter Brune           si = ii - ms;
175063654ddSPeter Brune           sj = jj - ns;
176063654ddSPeter Brune           sk = kk - ps;
177063654ddSPeter Brune           for (l=0;l<dof;l++) {
178063654ddSPeter Brune             indices[idx] = l + dof*(base + si + xm*sj + xm*ym*sk);
179ff9846d9SPeter Brune             idx++;
180ff9846d9SPeter Brune           }
181*3b5e53cdSSajid Ali           i++;
182*3b5e53cdSSajid Ali         } while (i<upper->i-ox);
183*3b5e53cdSSajid Ali         j++;
184*3b5e53cdSSajid Ali       } while (j<upper->j-oy);
185*3b5e53cdSSajid Ali       k++;
186*3b5e53cdSSajid Ali     } while (k<upper->k-oz);
187*3b5e53cdSSajid Ali   }
188*3b5e53cdSSajid Ali 
189*3b5e53cdSSajid Ali   if (!offproc){
190*3b5e53cdSSajid Ali     ierr = DMDAGetCorners(da, &ms, &ns, &ps, &mw, &nw, &pw);CHKERRQ(ierr);
191*3b5e53cdSSajid Ali     me = ms + mw;
192*3b5e53cdSSajid Ali     if (N>1) ne = ns + nw;
193*3b5e53cdSSajid Ali     if (P>1) pe = ps + pw;
194*3b5e53cdSSajid Ali     /* Account for DM offsets */
195*3b5e53cdSSajid Ali     ms = ms - ox; me = me - ox;
196*3b5e53cdSSajid Ali     ns = ns - oy; ne = ne - oy;
197*3b5e53cdSSajid Ali     ps = ps - oz; pe = pe - oz;
198*3b5e53cdSSajid Ali 
199*3b5e53cdSSajid Ali     /* compute the vector base on owning processor */
200*3b5e53cdSSajid Ali     xm = me - ms;
201*3b5e53cdSSajid Ali     ym = ne - ns;
202*3b5e53cdSSajid Ali     zm = pe - ps;
203*3b5e53cdSSajid Ali     base = ms*ym*zm + ns*M*zm + ps*M*N;
204*3b5e53cdSSajid Ali     /*
205*3b5e53cdSSajid Ali        if no indices are to be returned, create an empty is,
206*3b5e53cdSSajid Ali        this prevents hanging in while loops
207*3b5e53cdSSajid Ali     */
208*3b5e53cdSSajid Ali     if (skip_i && skip_j && skip_k) goto createis;
209*3b5e53cdSSajid Ali     /*
210*3b5e53cdSSajid Ali        do..while loops to ensure the block gets entered once,
211*3b5e53cdSSajid Ali        regardless of control condition being met, necessary for
212*3b5e53cdSSajid Ali        cases when a subset of skip_i/j/k is true
213*3b5e53cdSSajid Ali     */
214*3b5e53cdSSajid Ali     if (skip_k) k = upper->k-oz; else k = lower->k-oz;
215*3b5e53cdSSajid Ali     do {
216*3b5e53cdSSajid Ali       if (skip_j) j = upper->j-oy; else j = lower->j-oy;
217*3b5e53cdSSajid Ali       do {
218*3b5e53cdSSajid Ali         if (skip_i) i = upper->i-ox; else i = lower->i-ox;
219*3b5e53cdSSajid Ali         do {
220*3b5e53cdSSajid Ali           if (k>=ps && k<=pe-1) {
221*3b5e53cdSSajid Ali             if (j>=ns && j<=ne-1) {
222*3b5e53cdSSajid Ali               if (i>=ms && i<=me-1) {
223*3b5e53cdSSajid Ali                 /* compute the local coordinates on owning processor */
224*3b5e53cdSSajid Ali                 si = i - ms;
225*3b5e53cdSSajid Ali                 sj = j - ns;
226*3b5e53cdSSajid Ali                 sk = k - ps;
227*3b5e53cdSSajid Ali                 for (l=0; l<dof; l++) {
228*3b5e53cdSSajid Ali                   indices[idx] = l + dof*(base + si + xm*sj + xm*ym*sk);
229*3b5e53cdSSajid Ali                   idx++;
230ff9846d9SPeter Brune                 }
231ff9846d9SPeter Brune               }
232063654ddSPeter Brune             }
233*3b5e53cdSSajid Ali           }
234*3b5e53cdSSajid Ali           i++;
235*3b5e53cdSSajid Ali         } while (i<upper->i-ox);
236*3b5e53cdSSajid Ali         j++;
237*3b5e53cdSSajid Ali       } while (j<upper->j-oy);
238*3b5e53cdSSajid Ali       k++;
239*3b5e53cdSSajid Ali     } while (k<upper->k-oz);
240*3b5e53cdSSajid Ali 
241*3b5e53cdSSajid Ali     ierr = PetscRealloc((size_t)(idx*sizeof(PetscInt)), (void*)&indices);CHKERRQ(ierr);
242*3b5e53cdSSajid Ali   }
243*3b5e53cdSSajid Ali 
244*3b5e53cdSSajid Ali   createis:
245*3b5e53cdSSajid Ali   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)da),idx,indices,PETSC_OWN_POINTER,is);CHKERRQ(ierr);
246ff9846d9SPeter Brune   PetscFunctionReturn(0);
247ff9846d9SPeter Brune }
248ff9846d9SPeter Brune 
24990c77898SPeter Brune PetscErrorCode DMDASubDomainDA_Private(DM dm, PetscInt *nlocal, DM **sdm)
2500adebc6cSBarry Smith {
25190c77898SPeter Brune   DM             *da;
25290c77898SPeter Brune   PetscInt       dim,size,i,j,k,idx;
253e30e807fSPeter Brune   PetscErrorCode ierr;
254e30e807fSPeter Brune   DMDALocalInfo  info;
255ff9846d9SPeter Brune   PetscInt       xsize,ysize,zsize;
256e30e807fSPeter Brune   PetscInt       xo,yo,zo;
25790c77898SPeter Brune   PetscInt       xs,ys,zs;
25890c77898SPeter Brune   PetscInt       xm=1,ym=1,zm=1;
2597ddda789SPeter Brune   PetscInt       xol,yol,zol;
26090c77898SPeter Brune   PetscInt       m=1,n=1,p=1;
26190c77898SPeter Brune   PetscInt       M,N,P;
26290c77898SPeter Brune   PetscInt       pm,mtmp;
263e30e807fSPeter Brune 
264e30e807fSPeter Brune   PetscFunctionBegin;
265e30e807fSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
2667ddda789SPeter Brune   ierr = DMDAGetOverlap(dm,&xol,&yol,&zol);CHKERRQ(ierr);
26790c77898SPeter Brune   ierr = DMDAGetNumLocalSubDomains(dm,&size);CHKERRQ(ierr);
268785e854fSJed Brown   ierr = PetscMalloc1(size,&da);CHKERRQ(ierr);
2697ddda789SPeter Brune 
27090c77898SPeter Brune   if (nlocal) *nlocal = size;
271e30e807fSPeter Brune 
27290c77898SPeter Brune   dim = info.dim;
273e30e807fSPeter Brune 
27490c77898SPeter Brune   M = info.xm;
27590c77898SPeter Brune   N = info.ym;
27690c77898SPeter Brune   P = info.zm;
27790c77898SPeter Brune 
27890c77898SPeter Brune   if (dim == 1) {
27990c77898SPeter Brune     m = size;
28090c77898SPeter Brune   } else if (dim == 2) {
28190c77898SPeter Brune     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)N)));
28290c77898SPeter Brune     while (m > 0) {
28390c77898SPeter Brune       n = size/m;
28490c77898SPeter Brune       if (m*n*p == size) break;
28590c77898SPeter Brune       m--;
28690c77898SPeter Brune     }
28790c77898SPeter Brune   } else if (dim == 3) {
28890c77898SPeter Brune     n = (PetscInt)(0.5 + PetscPowReal(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),(PetscReal)(1./3.))); if (!n) n = 1;
28990c77898SPeter Brune     while (n > 0) {
29090c77898SPeter Brune       pm = size/n;
29190c77898SPeter Brune       if (n*pm == size) break;
29290c77898SPeter Brune       n--;
29390c77898SPeter Brune     }
29490c77898SPeter Brune     if (!n) n = 1;
29590c77898SPeter Brune     m = (PetscInt)(0.5 + PetscSqrtReal(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
29690c77898SPeter Brune     if (!m) m = 1;
29790c77898SPeter Brune     while (m > 0) {
29890c77898SPeter Brune       p = size/(m*n);
29990c77898SPeter Brune       if (m*n*p == size) break;
30090c77898SPeter Brune       m--;
30190c77898SPeter Brune     }
30290c77898SPeter Brune     if (M > P && m < p) {mtmp = m; m = p; p = mtmp;}
30390c77898SPeter Brune   }
30490c77898SPeter Brune 
30590c77898SPeter Brune   zs = info.zs;
30690c77898SPeter Brune   idx = 0;
30790c77898SPeter Brune   for (k = 0; k < p; k++) {
30890c77898SPeter Brune     ys = info.ys;
30990c77898SPeter Brune     for (j = 0; j < n; j++) {
31090c77898SPeter Brune       xs = info.xs;
31190c77898SPeter Brune       for (i = 0; i < m; i++) {
31290c77898SPeter Brune         if (dim == 1) {
31390c77898SPeter Brune           xm = M/m + ((M % m) > i);
31490c77898SPeter Brune         } else if (dim == 2) {
31590c77898SPeter Brune           xm = M/m + ((M % m) > i);
31690c77898SPeter Brune           ym = N/n + ((N % n) > j);
31790c77898SPeter Brune         } else if (dim == 3) {
31890c77898SPeter Brune           xm = M/m + ((M % m) > i);
31990c77898SPeter Brune           ym = N/n + ((N % n) > j);
32090c77898SPeter Brune           zm = P/p + ((P % p) > k);
32190c77898SPeter Brune         }
32290c77898SPeter Brune 
32390c77898SPeter Brune         xsize = xm;
32490c77898SPeter Brune         ysize = ym;
32590c77898SPeter Brune         zsize = zm;
32690c77898SPeter Brune         xo = xs;
32790c77898SPeter Brune         yo = ys;
32890c77898SPeter Brune         zo = zs;
32990c77898SPeter Brune 
33090c77898SPeter Brune         ierr = DMDACreate(PETSC_COMM_SELF,&(da[idx]));CHKERRQ(ierr);
33190c77898SPeter Brune         ierr = DMSetOptionsPrefix(da[idx],"sub_");CHKERRQ(ierr);
332c73cfb54SMatthew G. Knepley         ierr = DMSetDimension(da[idx], info.dim);CHKERRQ(ierr);
33390c77898SPeter Brune         ierr = DMDASetDof(da[idx], info.dof);CHKERRQ(ierr);
33490c77898SPeter Brune 
33590c77898SPeter Brune         ierr = DMDASetStencilType(da[idx],info.st);CHKERRQ(ierr);
33690c77898SPeter Brune         ierr = DMDASetStencilWidth(da[idx],info.sw);CHKERRQ(ierr);
33790c77898SPeter Brune 
338bff4a2f0SMatthew G. Knepley         if (info.bx == DM_BOUNDARY_PERIODIC || (xs != 0)) {
3397ddda789SPeter Brune           xsize += xol;
3407ddda789SPeter Brune           xo    -= xol;
341e30e807fSPeter Brune         }
342bff4a2f0SMatthew G. Knepley         if (info.by == DM_BOUNDARY_PERIODIC || (ys != 0)) {
3437ddda789SPeter Brune           ysize += yol;
3447ddda789SPeter Brune           yo    -= yol;
345e30e807fSPeter Brune         }
346bff4a2f0SMatthew G. Knepley         if (info.bz == DM_BOUNDARY_PERIODIC || (zs != 0)) {
3477ddda789SPeter Brune           zsize += zol;
3487ddda789SPeter Brune           zo    -= zol;
349e30e807fSPeter Brune         }
350e30e807fSPeter Brune 
351bff4a2f0SMatthew G. Knepley         if (info.bx == DM_BOUNDARY_PERIODIC || (xs+xm != info.mx)) xsize += xol;
352bff4a2f0SMatthew G. Knepley         if (info.by == DM_BOUNDARY_PERIODIC || (ys+ym != info.my)) ysize += yol;
353bff4a2f0SMatthew G. Knepley         if (info.bz == DM_BOUNDARY_PERIODIC || (zs+zm != info.mz)) zsize += zol;
3548865f1eaSKarl Rupp 
355bff4a2f0SMatthew G. Knepley         if (info.bx != DM_BOUNDARY_PERIODIC) {
35690c77898SPeter Brune           if (xo < 0) {
35790c77898SPeter Brune             xsize += xo;
35890c77898SPeter Brune             xo = 0;
35990c77898SPeter Brune           }
36090c77898SPeter Brune           if (xo+xsize > info.mx-1) {
36190c77898SPeter Brune             xsize -= xo+xsize - info.mx;
36290c77898SPeter Brune           }
36390c77898SPeter Brune         }
364bff4a2f0SMatthew G. Knepley         if (info.by != DM_BOUNDARY_PERIODIC) {
36590c77898SPeter Brune           if (yo < 0) {
36690c77898SPeter Brune             ysize += yo;
36790c77898SPeter Brune             yo = 0;
36890c77898SPeter Brune           }
36990c77898SPeter Brune           if (yo+ysize > info.my-1) {
37090c77898SPeter Brune             ysize -= yo+ysize - info.my;
37190c77898SPeter Brune           }
37290c77898SPeter Brune         }
373bff4a2f0SMatthew G. Knepley         if (info.bz != DM_BOUNDARY_PERIODIC) {
37490c77898SPeter Brune           if (zo < 0) {
37590c77898SPeter Brune             zsize += zo;
37690c77898SPeter Brune             zo = 0;
37790c77898SPeter Brune           }
37890c77898SPeter Brune           if (zo+zsize > info.mz-1) {
37990c77898SPeter Brune             zsize -= zo+zsize - info.mz;
38090c77898SPeter Brune           }
38190c77898SPeter Brune         }
38290c77898SPeter Brune 
38390c77898SPeter Brune         ierr = DMDASetSizes(da[idx], xsize, ysize, zsize);CHKERRQ(ierr);
38490c77898SPeter Brune         ierr = DMDASetNumProcs(da[idx], 1, 1, 1);CHKERRQ(ierr);
385bff4a2f0SMatthew G. Knepley         ierr = DMDASetBoundaryType(da[idx], DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED);CHKERRQ(ierr);
386e30e807fSPeter Brune 
387e30e807fSPeter Brune         /* set up as a block instead */
38890c77898SPeter Brune         ierr = DMSetUp(da[idx]);CHKERRQ(ierr);
389e30e807fSPeter Brune 
39040234c92SPeter Brune         /* nonoverlapping region */
39190c77898SPeter Brune         ierr = DMDASetNonOverlappingRegion(da[idx],xs,ys,zs,xm,ym,zm);CHKERRQ(ierr);
39240234c92SPeter Brune 
39395c13181SPeter Brune         /* this alters the behavior of DMDAGetInfo, DMDAGetLocalInfo, DMDAGetCorners, and DMDAGetGhostedCorners and should be used with care */
39490c77898SPeter Brune         ierr = DMDASetOffset(da[idx],xo,yo,zo,info.mx,info.my,info.mz);CHKERRQ(ierr);
39590c77898SPeter Brune         xs += xm;
39690c77898SPeter Brune         idx++;
39790c77898SPeter Brune       }
39890c77898SPeter Brune       ys += ym;
39990c77898SPeter Brune     }
40090c77898SPeter Brune     zs += zm;
40190c77898SPeter Brune   }
40290c77898SPeter Brune   *sdm = da;
403e30e807fSPeter Brune   PetscFunctionReturn(0);
404e30e807fSPeter Brune }
405e30e807fSPeter Brune 
406e30e807fSPeter Brune /*
407e30e807fSPeter Brune    Fills the local vector problem on the subdomain from the global problem.
408e30e807fSPeter Brune 
409285ae305SPeter Brune    Right now this assumes one subdomain per processor.
410285ae305SPeter Brune 
411e30e807fSPeter Brune */
4120adebc6cSBarry Smith PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat)
4130adebc6cSBarry Smith {
414e30e807fSPeter Brune   PetscErrorCode ierr;
415285ae305SPeter Brune   DMDALocalInfo  info,subinfo;
416e30e807fSPeter Brune   DM             subdm;
417285ae305SPeter Brune   MatStencil     upper,lower;
418285ae305SPeter Brune   IS             idis,isis,odis,osis,gdis;
419285ae305SPeter Brune   Vec            svec,dvec,slvec;
42040234c92SPeter Brune   PetscInt       xm,ym,zm,xs,ys,zs;
42190c77898SPeter Brune   PetscInt       i;
422*3b5e53cdSSajid Ali   PetscBool      patchis_offproc = PETSC_TRUE;
423e30e807fSPeter Brune 
424e30e807fSPeter Brune   PetscFunctionBegin;
425e30e807fSPeter Brune   /* allocate the arrays of scatters */
426785e854fSJed Brown   if (iscat) {ierr = PetscMalloc1(nsubdms,iscat);CHKERRQ(ierr);}
427785e854fSJed Brown   if (oscat) {ierr = PetscMalloc1(nsubdms,oscat);CHKERRQ(ierr);}
428785e854fSJed Brown   if (lscat) {ierr = PetscMalloc1(nsubdms,lscat);CHKERRQ(ierr);}
429e30e807fSPeter Brune 
430285ae305SPeter Brune   ierr  = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
43190c77898SPeter Brune   for (i = 0; i < nsubdms; i++) {
43290c77898SPeter Brune     subdm = subdms[i];
433285ae305SPeter Brune     ierr  = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr);
43490c77898SPeter Brune     ierr  = DMDAGetNonOverlappingRegion(subdm,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
435e30e807fSPeter Brune 
436285ae305SPeter Brune     /* create the global and subdomain index sets for the inner domain */
43740234c92SPeter Brune     lower.i = xs;
43840234c92SPeter Brune     lower.j = ys;
43940234c92SPeter Brune     lower.k = zs;
44040234c92SPeter Brune     upper.i = xs+xm;
44140234c92SPeter Brune     upper.j = ys+ym;
44240234c92SPeter Brune     upper.k = zs+zm;
443*3b5e53cdSSajid Ali     ierr    = DMDACreatePatchIS(dm,&lower,&upper,&idis,patchis_offproc);CHKERRQ(ierr);
444*3b5e53cdSSajid Ali     ierr    = DMDACreatePatchIS(subdm,&lower,&upper,&isis,patchis_offproc);CHKERRQ(ierr);
445e30e807fSPeter Brune 
446285ae305SPeter Brune     /* create the global and subdomain index sets for the outer subdomain */
447285ae305SPeter Brune     lower.i = subinfo.xs;
448285ae305SPeter Brune     lower.j = subinfo.ys;
449285ae305SPeter Brune     lower.k = subinfo.zs;
450285ae305SPeter Brune     upper.i = subinfo.xs+subinfo.xm;
451285ae305SPeter Brune     upper.j = subinfo.ys+subinfo.ym;
452285ae305SPeter Brune     upper.k = subinfo.zs+subinfo.zm;
453*3b5e53cdSSajid Ali     ierr    = DMDACreatePatchIS(dm,&lower,&upper,&odis,patchis_offproc);CHKERRQ(ierr);
454*3b5e53cdSSajid Ali     ierr    = DMDACreatePatchIS(subdm,&lower,&upper,&osis,patchis_offproc);CHKERRQ(ierr);
455e30e807fSPeter Brune 
456285ae305SPeter Brune     /* global and subdomain ISes for the local indices of the subdomain */
457285ae305SPeter Brune     /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */
458285ae305SPeter Brune     lower.i = subinfo.gxs;
459285ae305SPeter Brune     lower.j = subinfo.gys;
460285ae305SPeter Brune     lower.k = subinfo.gzs;
461285ae305SPeter Brune     upper.i = subinfo.gxs+subinfo.gxm;
462285ae305SPeter Brune     upper.j = subinfo.gys+subinfo.gym;
463285ae305SPeter Brune     upper.k = subinfo.gzs+subinfo.gzm;
464*3b5e53cdSSajid Ali     ierr    = DMDACreatePatchIS(dm,&lower,&upper,&gdis,patchis_offproc);CHKERRQ(ierr);
465e30e807fSPeter Brune 
466e30e807fSPeter Brune     /* form the scatter */
467e30e807fSPeter Brune     ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr);
468e30e807fSPeter Brune     ierr = DMGetGlobalVector(subdm,&svec);CHKERRQ(ierr);
469e30e807fSPeter Brune     ierr = DMGetLocalVector(subdm,&slvec);CHKERRQ(ierr);
470e30e807fSPeter Brune 
4719448b7f1SJunchao Zhang     if (iscat) {ierr = VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[i]);CHKERRQ(ierr);}
4729448b7f1SJunchao Zhang     if (oscat) {ierr = VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[i]);CHKERRQ(ierr);}
4739448b7f1SJunchao Zhang     if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,NULL,&(*lscat)[i]);CHKERRQ(ierr);}
474e30e807fSPeter Brune 
475e30e807fSPeter Brune     ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr);
476e30e807fSPeter Brune     ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr);
477e30e807fSPeter Brune     ierr = DMRestoreLocalVector(subdm,&slvec);CHKERRQ(ierr);
478e30e807fSPeter Brune 
479e30e807fSPeter Brune     ierr = ISDestroy(&idis);CHKERRQ(ierr);
480e30e807fSPeter Brune     ierr = ISDestroy(&isis);CHKERRQ(ierr);
481e30e807fSPeter Brune 
482e30e807fSPeter Brune     ierr = ISDestroy(&odis);CHKERRQ(ierr);
483e30e807fSPeter Brune     ierr = ISDestroy(&osis);CHKERRQ(ierr);
484e30e807fSPeter Brune 
485e30e807fSPeter Brune     ierr = ISDestroy(&gdis);CHKERRQ(ierr);
48690c77898SPeter Brune   }
487e30e807fSPeter Brune   PetscFunctionReturn(0);
488e30e807fSPeter Brune }
489e30e807fSPeter Brune 
49090c77898SPeter Brune PetscErrorCode DMDASubDomainIS_Private(DM dm,PetscInt n,DM *subdm,IS **iis,IS **ois)
4910adebc6cSBarry Smith {
492e30e807fSPeter Brune   PetscErrorCode ierr;
49390c77898SPeter Brune   PetscInt       i;
494e30e807fSPeter Brune   DMDALocalInfo  info,subinfo;
495285ae305SPeter Brune   MatStencil     lower,upper;
496*3b5e53cdSSajid Ali   PetscBool      patchis_offproc = PETSC_TRUE;
497e30e807fSPeter Brune 
498e30e807fSPeter Brune   PetscFunctionBegin;
499e30e807fSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
500785e854fSJed Brown   if (iis) {ierr = PetscMalloc1(n,iis);CHKERRQ(ierr);}
501785e854fSJed Brown   if (ois) {ierr = PetscMalloc1(n,ois);CHKERRQ(ierr);}
502e30e807fSPeter Brune 
50390c77898SPeter Brune   for (i = 0;i < n; i++) {
50490c77898SPeter Brune     ierr = DMDAGetLocalInfo(subdm[i],&subinfo);CHKERRQ(ierr);
50590c77898SPeter Brune     if (iis) {
506285ae305SPeter Brune       /* create the inner IS */
507285ae305SPeter Brune       lower.i = info.xs;
508285ae305SPeter Brune       lower.j = info.ys;
509285ae305SPeter Brune       lower.k = info.zs;
510285ae305SPeter Brune       upper.i = info.xs+info.xm;
511285ae305SPeter Brune       upper.j = info.ys+info.ym;
512285ae305SPeter Brune       upper.k = info.zs+info.zm;
513*3b5e53cdSSajid Ali       ierr = DMDACreatePatchIS(dm,&lower,&upper,&(*iis)[i],patchis_offproc);CHKERRQ(ierr);
51490c77898SPeter Brune     }
515e30e807fSPeter Brune 
51690c77898SPeter Brune     if (ois) {
517285ae305SPeter Brune       /* create the outer IS */
518285ae305SPeter Brune       lower.i = subinfo.xs;
519285ae305SPeter Brune       lower.j = subinfo.ys;
520285ae305SPeter Brune       lower.k = subinfo.zs;
521285ae305SPeter Brune       upper.i = subinfo.xs+subinfo.xm;
522285ae305SPeter Brune       upper.j = subinfo.ys+subinfo.ym;
523285ae305SPeter Brune       upper.k = subinfo.zs+subinfo.zm;
524*3b5e53cdSSajid Ali       ierr    = DMDACreatePatchIS(dm,&lower,&upper,&(*ois)[i],patchis_offproc);CHKERRQ(ierr);
52590c77898SPeter Brune     }
52690c77898SPeter Brune   }
527e30e807fSPeter Brune   PetscFunctionReturn(0);
528e30e807fSPeter Brune }
529e30e807fSPeter Brune 
5300adebc6cSBarry Smith PetscErrorCode DMCreateDomainDecomposition_DA(DM dm,PetscInt *len,char ***names,IS **iis,IS **ois,DM **subdm)
5310adebc6cSBarry Smith {
532e30e807fSPeter Brune   PetscErrorCode ierr;
53390c77898SPeter Brune   DM             *sdm;
53490c77898SPeter Brune   PetscInt       n,i;
5350a696f66SPeter Brune 
5360adebc6cSBarry Smith   PetscFunctionBegin;
53790c77898SPeter Brune   ierr = DMDASubDomainDA_Private(dm,&n,&sdm);CHKERRQ(ierr);
53890c77898SPeter Brune   if (names) {
539785e854fSJed Brown     ierr = PetscMalloc1(n,names);CHKERRQ(ierr);
540ea78f98cSLisandro Dalcin     for (i=0;i<n;i++) (*names)[i] = NULL;
541e30e807fSPeter Brune   }
54290c77898SPeter Brune   ierr = DMDASubDomainIS_Private(dm,n,sdm,iis,ois);CHKERRQ(ierr);
54390c77898SPeter Brune   if (subdm) *subdm = sdm;
544e2c616c8SPeter Brune   else {
545e2c616c8SPeter Brune     for (i=0;i<n;i++) {
546e2c616c8SPeter Brune       ierr = DMDestroy(&sdm[i]);CHKERRQ(ierr);
547e2c616c8SPeter Brune     }
548e2c616c8SPeter Brune   }
54990c77898SPeter Brune   if (len) *len = n;
550e30e807fSPeter Brune   PetscFunctionReturn(0);
551e30e807fSPeter Brune }
552