xref: /petsc/src/dm/impls/da/dadd.c (revision 063654dd348a6ea583232b5381e465b2ed5ab0f0)
14035e84dSBarry Smith #include <petsc-private/dmdaimpl.h>
2e30e807fSPeter Brune 
3e30e807fSPeter Brune #undef __FUNCT__
4ff9846d9SPeter Brune #define __FUNCT__ "DMDACreatePatchIS"
5*063654ddSPeter Brune /*@
6*063654ddSPeter Brune   DMDACreatePatchIS - Creates an index set corresponding to a patch of the DA.
7ff9846d9SPeter Brune 
8*063654ddSPeter Brune   Not Collective
9*063654ddSPeter Brune 
10*063654ddSPeter Brune   Input Parameters:
11*063654ddSPeter Brune +  da - the DMDA
12*063654ddSPeter Brune .  lower - a matstencil with i, j and k corresponding to the lower corner of the patch
13*063654ddSPeter Brune -  upper - a matstencil with i, j and k corresponding to the upper corner of the patch
14*063654ddSPeter Brune 
15*063654ddSPeter Brune   Output Parameters:
16*063654ddSPeter Brune .  is - the IS corresponding to the patch
17*063654ddSPeter Brune 
18*063654ddSPeter Brune   Level: developer
19*063654ddSPeter Brune 
20*063654ddSPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDACreateDomainDecompositionScatters()
21*063654ddSPeter Brune @*/
22ff9846d9SPeter Brune PetscErrorCode DMDACreatePatchIS(DM da,MatStencil *lower,MatStencil *upper,IS *is)
23ff9846d9SPeter Brune {
24*063654ddSPeter Brune   PetscInt       ms=0,ns=0,ps=0;
25*063654ddSPeter Brune   PetscInt       me=1,ne=1,pe=1;
26*063654ddSPeter Brune   PetscInt       mr=0,nr=0,pr=0;
27ff9846d9SPeter Brune   PetscInt       ii,jj,kk;
28*063654ddSPeter Brune   PetscInt       si,sj,sk;
29*063654ddSPeter Brune   PetscInt       i,j,k,l,idx;
30ff9846d9SPeter Brune   PetscInt       base;
31*063654ddSPeter Brune   PetscInt       xm=1,ym=1,zm=1;
32*063654ddSPeter Brune   const PetscInt *lx,*ly,*lz;
33ff9846d9SPeter Brune   PetscInt       ox,oy,oz;
34*063654ddSPeter Brune   PetscInt       m,n,p,M,N,P,dof;
35*063654ddSPeter Brune   PetscInt       nindices;
36*063654ddSPeter Brune   PetscInt       *indices;
37*063654ddSPeter Brune   DM_DA          *dd = (DM_DA*)da->data;
38*063654ddSPeter Brune   PetscErrorCode ierr;
39ff9846d9SPeter Brune 
400adebc6cSBarry Smith   PetscFunctionBegin;
41*063654ddSPeter Brune   /* need to get the sizes of the actual DM rather than the "global" space of a subdomain DM */
42*063654ddSPeter Brune   M = dd->M;N = dd->N;P=dd->P;
43*063654ddSPeter Brune   m = dd->m;n = dd->n;p=dd->p;
44*063654ddSPeter Brune   dof = dd->w;
450298fd71SBarry Smith   ierr = DMDAGetOffset(da,&ox,&oy,&oz,NULL,NULL,NULL);CHKERRQ(ierr);
46*063654ddSPeter Brune   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr);
47*063654ddSPeter Brune   nindices = (upper->i - lower->i)*(upper->j - lower->j)*(upper->k - lower->k)*dof;
48*063654ddSPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*nindices,&indices);CHKERRQ(ierr);
49*063654ddSPeter Brune   /* start at index 0 on processor 0 */
50*063654ddSPeter Brune   mr = 0;
51*063654ddSPeter Brune   nr = 0;
52*063654ddSPeter Brune   pr = 0;
53*063654ddSPeter Brune   ms = 0;
54*063654ddSPeter Brune   ns = 0;
55*063654ddSPeter Brune   ps = 0;
56*063654ddSPeter Brune   if (lx) me = lx[0];
57*063654ddSPeter Brune   if (ly) ne = ly[0];
58*063654ddSPeter Brune   if (lz) pe = lz[0];
59ff9846d9SPeter Brune   idx = 0;
60ff9846d9SPeter Brune   for (k=lower->k-oz;k<upper->k-oz;k++) {
61ff9846d9SPeter Brune     for (j=lower->j-oy;j < upper->j-oy;j++) {
62*063654ddSPeter Brune       for (i=lower->i-ox;i < upper->i-ox;i++) {
63*063654ddSPeter Brune         /* "actual" indices rather than ones outside of the domain */
64*063654ddSPeter Brune         ii = i;
65*063654ddSPeter Brune         jj = j;
66*063654ddSPeter Brune         kk = k;
67*063654ddSPeter Brune         if (ii < 0) ii = M + ii;
68*063654ddSPeter Brune         if (jj < 0) jj = N + jj;
69*063654ddSPeter Brune         if (kk < 0) kk = P + kk;
70*063654ddSPeter Brune         if (ii > M-1) ii = ii - M;
71*063654ddSPeter Brune         if (jj > N-1) jj = jj - N;
72*063654ddSPeter Brune         if (kk > P-1) kk = kk - P;
73*063654ddSPeter Brune         /* gone out of processor range on x axis */
74*063654ddSPeter Brune         while(ii > me-1 || ii < ms) {
75*063654ddSPeter Brune           if (mr == m-1) {
76*063654ddSPeter Brune             ms = 0;
77*063654ddSPeter Brune             me = lx[0];
78*063654ddSPeter Brune             mr = 0;
79*063654ddSPeter Brune           } else {
80*063654ddSPeter Brune             mr++;
81*063654ddSPeter Brune             ms = me;
82*063654ddSPeter Brune             me += lx[mr];
83*063654ddSPeter Brune           }
84*063654ddSPeter Brune         }
85*063654ddSPeter Brune         /* gone out of processor range on y axis */
86*063654ddSPeter Brune         while(jj > ne-1 || jj < ns) {
87*063654ddSPeter Brune           if (nr == n-1) {
88*063654ddSPeter Brune             ns = 0;
89*063654ddSPeter Brune             ne = ly[0];
90*063654ddSPeter Brune             nr = 0;
91*063654ddSPeter Brune           } else {
92*063654ddSPeter Brune             nr++;
93*063654ddSPeter Brune             ns = ne;
94*063654ddSPeter Brune             ne += ly[nr];
95*063654ddSPeter Brune           }
96*063654ddSPeter Brune         }
97*063654ddSPeter Brune         /* gone out of processor range on z axis */
98*063654ddSPeter Brune         while(kk > pe-1 || kk < ps) {
99*063654ddSPeter Brune           if (pr == p-1) {
100*063654ddSPeter Brune             ps = 0;
101*063654ddSPeter Brune             pe = lz[0];
102*063654ddSPeter Brune             pr = 0;
103*063654ddSPeter Brune           } else {
104*063654ddSPeter Brune             pr++;
105*063654ddSPeter Brune             ps = pe;
106*063654ddSPeter Brune             pe += lz[pr];
107*063654ddSPeter Brune           }
108*063654ddSPeter Brune         }
109*063654ddSPeter Brune         /* compute the vector base on owning processor */
110*063654ddSPeter Brune         xm = me - ms;
111*063654ddSPeter Brune         ym = ne - ns;
112*063654ddSPeter Brune         zm = pe - ps;
113*063654ddSPeter Brune         base = ms*ym*zm + ns*M + ps*M*N;
114*063654ddSPeter Brune         /* compute the local coordinates on owning processor */
115*063654ddSPeter Brune         si = ii - ms;
116*063654ddSPeter Brune         sj = jj - ns;
117*063654ddSPeter Brune         sk = kk - ps;
118*063654ddSPeter Brune         for (l=0;l<dof;l++) {
119*063654ddSPeter Brune           indices[idx] = l + dof*(base + si + xm*sj + xm*ym*sk);
120ff9846d9SPeter Brune           idx++;
121ff9846d9SPeter Brune         }
122ff9846d9SPeter Brune       }
123ff9846d9SPeter Brune     }
124*063654ddSPeter Brune   }
125*063654ddSPeter Brune   ISCreateGeneral(PETSC_COMM_SELF,idx,indices,PETSC_OWN_POINTER,is);CHKERRQ(ierr);
126ff9846d9SPeter Brune   PetscFunctionReturn(0);
127ff9846d9SPeter Brune }
128ff9846d9SPeter Brune 
129ff9846d9SPeter Brune #undef __FUNCT__
130e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainDA_Private"
1310adebc6cSBarry Smith PetscErrorCode DMDASubDomainDA_Private(DM dm, DM *dddm)
1320adebc6cSBarry Smith {
133e30e807fSPeter Brune   DM             da;
134e30e807fSPeter Brune   PetscErrorCode ierr;
135e30e807fSPeter Brune   DMDALocalInfo  info;
136e30e807fSPeter Brune   PetscReal      lmin[3],lmax[3];
137ff9846d9SPeter Brune   PetscInt       xsize,ysize,zsize;
138e30e807fSPeter Brune   PetscInt       xo,yo,zo;
1397ddda789SPeter Brune   PetscInt       xol,yol,zol;
140e30e807fSPeter Brune 
141e30e807fSPeter Brune   PetscFunctionBegin;
142e30e807fSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
1437ddda789SPeter Brune   ierr = DMDAGetOverlap(dm,&xol,&yol,&zol);CHKERRQ(ierr);
1447ddda789SPeter Brune 
145e30e807fSPeter Brune   ierr = DMDACreate(PETSC_COMM_SELF,&da);CHKERRQ(ierr);
146e30e807fSPeter Brune   ierr = DMSetOptionsPrefix(da,"sub_");CHKERRQ(ierr);
147e30e807fSPeter Brune   ierr = DMDASetDim(da, info.dim);CHKERRQ(ierr);
148e30e807fSPeter Brune   ierr = DMDASetDof(da, info.dof);CHKERRQ(ierr);
149e30e807fSPeter Brune 
150e30e807fSPeter Brune   ierr = DMDASetStencilType(da,info.st);CHKERRQ(ierr);
151e30e807fSPeter Brune   ierr = DMDASetStencilWidth(da,info.sw);CHKERRQ(ierr);
152e30e807fSPeter Brune 
153e30e807fSPeter Brune   /* local with overlap */
154e30e807fSPeter Brune   xsize = info.xm;
155e30e807fSPeter Brune   ysize = info.ym;
156e30e807fSPeter Brune   zsize = info.zm;
157e30e807fSPeter Brune   xo    = info.xs;
158e30e807fSPeter Brune   yo    = info.ys;
159e30e807fSPeter Brune   zo    = info.zs;
160e30e807fSPeter Brune   if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs != 0)) {
1617ddda789SPeter Brune     xsize += xol;
1627ddda789SPeter Brune     xo    -= xol;
163e30e807fSPeter Brune   }
164e30e807fSPeter Brune   if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys != 0)) {
1657ddda789SPeter Brune     ysize += yol;
1667ddda789SPeter Brune     yo    -= yol;
167e30e807fSPeter Brune   }
168e30e807fSPeter Brune   if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs != 0)) {
1697ddda789SPeter Brune     zsize += zol;
1707ddda789SPeter Brune     zo    -= zol;
171e30e807fSPeter Brune   }
172e30e807fSPeter Brune 
1737ddda789SPeter Brune   if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs+info.xm != info.mx)) xsize += xol;
1747ddda789SPeter Brune   if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys+info.ym != info.my)) ysize += yol;
1757ddda789SPeter Brune   if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs+info.zm != info.mz)) zsize += zol;
1768865f1eaSKarl Rupp 
177e30e807fSPeter Brune   ierr = DMDASetSizes(da, xsize, ysize, zsize);CHKERRQ(ierr);
178e30e807fSPeter Brune   ierr = DMDASetNumProcs(da, 1, 1, 1);CHKERRQ(ierr);
179e30e807fSPeter Brune   ierr = DMDASetBoundaryType(da, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED);CHKERRQ(ierr);
180e30e807fSPeter Brune 
181e30e807fSPeter Brune   /* set up as a block instead */
182e30e807fSPeter Brune   ierr = DMSetUp(da);CHKERRQ(ierr);
183e30e807fSPeter Brune 
184e30e807fSPeter Brune   /* todo - nonuniform coordinates */
185e30e807fSPeter Brune   ierr = DMDAGetLocalBoundingBox(dm,lmin,lmax);CHKERRQ(ierr);
186e30e807fSPeter Brune   ierr = DMDASetUniformCoordinates(da,lmin[0],lmax[0],lmin[1],lmax[1],lmin[2],lmax[2]);CHKERRQ(ierr);
187e30e807fSPeter Brune 
18840234c92SPeter Brune   /* nonoverlapping region */
18940234c92SPeter Brune   ierr = DMDASetNonOverlappingRegion(da,info.xs,info.ys,info.zs,info.xm,info.ym,info.zm);CHKERRQ(ierr);
19040234c92SPeter Brune 
19195c13181SPeter Brune   /* this alters the behavior of DMDAGetInfo, DMDAGetLocalInfo, DMDAGetCorners, and DMDAGetGhostedCorners and should be used with care */
19295c13181SPeter Brune   ierr = DMDASetOffset(da,xo,yo,zo,info.mx,info.my,info.mz);CHKERRQ(ierr);
193e30e807fSPeter Brune   *dddm = da;
194e30e807fSPeter Brune   PetscFunctionReturn(0);
195e30e807fSPeter Brune }
196e30e807fSPeter Brune 
197e30e807fSPeter Brune #undef __FUNCT__
198e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecompositionScatters_DA"
199e30e807fSPeter Brune /*
200e30e807fSPeter Brune  Fills the local vector problem on the subdomain from the global problem.
201e30e807fSPeter Brune 
202285ae305SPeter Brune  Right now this assumes one subdomain per processor.
203285ae305SPeter Brune 
204e30e807fSPeter Brune  */
2050adebc6cSBarry Smith PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat)
2060adebc6cSBarry Smith {
207e30e807fSPeter Brune   PetscErrorCode ierr;
208285ae305SPeter Brune   DMDALocalInfo  info,subinfo;
209e30e807fSPeter Brune   DM             subdm;
210285ae305SPeter Brune   MatStencil     upper,lower;
211285ae305SPeter Brune   IS             idis,isis,odis,osis,gdis;
212285ae305SPeter Brune   Vec            svec,dvec,slvec;
21340234c92SPeter Brune   PetscInt       xm,ym,zm,xs,ys,zs;
214e30e807fSPeter Brune 
215e30e807fSPeter Brune   PetscFunctionBegin;
216ce94432eSBarry Smith   if (nsubdms != 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have more than one subdomain per processor (yet)");
217285ae305SPeter Brune 
218e30e807fSPeter Brune   /* allocate the arrays of scatters */
219050c5e14SJed Brown   if (iscat) {ierr = PetscMalloc(sizeof(VecScatter*),iscat);CHKERRQ(ierr);}
220050c5e14SJed Brown   if (oscat) {ierr = PetscMalloc(sizeof(VecScatter*),oscat);CHKERRQ(ierr);}
221050c5e14SJed Brown   if (lscat) {ierr = PetscMalloc(sizeof(VecScatter*),lscat);CHKERRQ(ierr);}
222e30e807fSPeter Brune 
223285ae305SPeter Brune   ierr  = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
224285ae305SPeter Brune   subdm = subdms[0];
225285ae305SPeter Brune   ierr  = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr);
226e30e807fSPeter Brune 
227285ae305SPeter Brune   /* create the global and subdomain index sets for the inner domain */
228285ae305SPeter Brune   /* TODO - make this actually support multiple subdomains -- subdomain needs to provide where it's nonoverlapping portion belongs */
22940234c92SPeter Brune   ierr = DMDAGetNonOverlappingRegion(subdm,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
23040234c92SPeter Brune   lower.i = xs;
23140234c92SPeter Brune   lower.j = ys;
23240234c92SPeter Brune   lower.k = zs;
23340234c92SPeter Brune   upper.i = xs+xm;
23440234c92SPeter Brune   upper.j = ys+ym;
23540234c92SPeter Brune   upper.k = zs+zm;
236285ae305SPeter Brune   ierr    = DMDACreatePatchIS(dm,&lower,&upper,&idis);CHKERRQ(ierr);
237285ae305SPeter Brune   ierr    = DMDACreatePatchIS(subdm,&lower,&upper,&isis);CHKERRQ(ierr);
238e30e807fSPeter Brune 
239285ae305SPeter Brune   /* create the global and subdomain index sets for the outer subdomain */
240285ae305SPeter Brune   lower.i = subinfo.xs;
241285ae305SPeter Brune   lower.j = subinfo.ys;
242285ae305SPeter Brune   lower.k = subinfo.zs;
243285ae305SPeter Brune   upper.i = subinfo.xs+subinfo.xm;
244285ae305SPeter Brune   upper.j = subinfo.ys+subinfo.ym;
245285ae305SPeter Brune   upper.k = subinfo.zs+subinfo.zm;
246285ae305SPeter Brune   ierr    = DMDACreatePatchIS(dm,&lower,&upper,&odis);CHKERRQ(ierr);
247285ae305SPeter Brune   ierr    = DMDACreatePatchIS(subdm,&lower,&upper,&osis);CHKERRQ(ierr);
248e30e807fSPeter Brune 
249285ae305SPeter Brune   /* global and subdomain ISes for the local indices of the subdomain */
250285ae305SPeter Brune   /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */
251285ae305SPeter Brune   lower.i = subinfo.gxs;
252285ae305SPeter Brune   lower.j = subinfo.gys;
253285ae305SPeter Brune   lower.k = subinfo.gzs;
254285ae305SPeter Brune   upper.i = subinfo.gxs+subinfo.gxm;
255285ae305SPeter Brune   upper.j = subinfo.gys+subinfo.gym;
256285ae305SPeter Brune   upper.k = subinfo.gzs+subinfo.gzm;
257e30e807fSPeter Brune 
258285ae305SPeter Brune   ierr = DMDACreatePatchIS(dm,&lower,&upper,&gdis);CHKERRQ(ierr);
259e30e807fSPeter Brune 
260e30e807fSPeter Brune   /* form the scatter */
261e30e807fSPeter Brune   ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr);
262e30e807fSPeter Brune   ierr = DMGetGlobalVector(subdm,&svec);CHKERRQ(ierr);
263e30e807fSPeter Brune   ierr = DMGetLocalVector(subdm,&slvec);CHKERRQ(ierr);
264e30e807fSPeter Brune 
265285ae305SPeter Brune   if (iscat) {ierr = VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[0]);CHKERRQ(ierr);}
266285ae305SPeter Brune   if (oscat) {ierr = VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[0]);CHKERRQ(ierr);}
2670298fd71SBarry Smith   if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,NULL,&(*lscat)[0]);CHKERRQ(ierr);}
268e30e807fSPeter Brune 
269e30e807fSPeter Brune   ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr);
270e30e807fSPeter Brune   ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr);
271e30e807fSPeter Brune   ierr = DMRestoreLocalVector(subdm,&slvec);CHKERRQ(ierr);
272e30e807fSPeter Brune 
273e30e807fSPeter Brune   ierr = ISDestroy(&idis);CHKERRQ(ierr);
274e30e807fSPeter Brune   ierr = ISDestroy(&isis);CHKERRQ(ierr);
275e30e807fSPeter Brune 
276e30e807fSPeter Brune   ierr = ISDestroy(&odis);CHKERRQ(ierr);
277e30e807fSPeter Brune   ierr = ISDestroy(&osis);CHKERRQ(ierr);
278e30e807fSPeter Brune 
279e30e807fSPeter Brune   ierr = ISDestroy(&gdis);CHKERRQ(ierr);
280e30e807fSPeter Brune   PetscFunctionReturn(0);
281e30e807fSPeter Brune }
282e30e807fSPeter Brune 
283e30e807fSPeter Brune #undef __FUNCT__
284e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainIS_Private"
285e30e807fSPeter Brune /* We have that the interior regions are going to be the same, but the ghost regions might not match up
286e30e807fSPeter Brune 
287e30e807fSPeter Brune ----------
288e30e807fSPeter Brune ----------
289e30e807fSPeter Brune --++++++o=
290e30e807fSPeter Brune --++++++o=
291e30e807fSPeter Brune --++++++o=
292e30e807fSPeter Brune --++++++o=
293e30e807fSPeter Brune --ooooooo=
294e30e807fSPeter Brune --========
295e30e807fSPeter Brune 
296e30e807fSPeter Brune Therefore, for each point in the overall, we must check if it's:
297e30e807fSPeter Brune 
298e30e807fSPeter Brune 1. +: In the interior of the global dm; it lines up
299e30e807fSPeter Brune 2. o: In the overlap region -- for now the same as 1; no overlap
300e30e807fSPeter Brune 3. =: In the shared ghost region -- handled by DMCreateDomainDecompositionLocalScatter()
301e30e807fSPeter Brune 4. -: In the nonshared ghost region
302e30e807fSPeter Brune  */
303e30e807fSPeter Brune 
3040adebc6cSBarry Smith PetscErrorCode DMDASubDomainIS_Private(DM dm,DM subdm,IS *iis,IS *ois)
3050adebc6cSBarry Smith {
306e30e807fSPeter Brune   PetscErrorCode ierr;
307e30e807fSPeter Brune   DMDALocalInfo  info,subinfo;
308285ae305SPeter Brune   MatStencil     lower,upper;
309e30e807fSPeter Brune 
310e30e807fSPeter Brune   PetscFunctionBegin;
311e30e807fSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
312e30e807fSPeter Brune   ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr);
313e30e807fSPeter Brune 
314285ae305SPeter Brune   /* create the inner IS */
315285ae305SPeter Brune   lower.i = info.xs;
316285ae305SPeter Brune   lower.j = info.ys;
317285ae305SPeter Brune   lower.k = info.zs;
318285ae305SPeter Brune   upper.i = info.xs+info.xm;
319285ae305SPeter Brune   upper.j = info.ys+info.ym;
320285ae305SPeter Brune   upper.k = info.zs+info.zm;
321e30e807fSPeter Brune 
322285ae305SPeter Brune   ierr = DMDACreatePatchIS(dm,&lower,&upper,iis);CHKERRQ(ierr);
323285ae305SPeter Brune 
324285ae305SPeter Brune   /* create the outer IS */
325285ae305SPeter Brune   lower.i = subinfo.xs;
326285ae305SPeter Brune   lower.j = subinfo.ys;
327285ae305SPeter Brune   lower.k = subinfo.zs;
328285ae305SPeter Brune   upper.i = subinfo.xs+subinfo.xm;
329285ae305SPeter Brune   upper.j = subinfo.ys+subinfo.ym;
330285ae305SPeter Brune   upper.k = subinfo.zs+subinfo.zm;
331285ae305SPeter Brune   ierr    = DMDACreatePatchIS(dm,&lower,&upper,ois);CHKERRQ(ierr);
332e30e807fSPeter Brune   PetscFunctionReturn(0);
333e30e807fSPeter Brune }
334e30e807fSPeter Brune 
335e30e807fSPeter Brune #undef __FUNCT__
336e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecomposition_DA"
3370adebc6cSBarry Smith PetscErrorCode DMCreateDomainDecomposition_DA(DM dm,PetscInt *len,char ***names,IS **iis,IS **ois,DM **subdm)
3380adebc6cSBarry Smith {
339e30e807fSPeter Brune   PetscErrorCode ierr;
340e30e807fSPeter Brune   IS             iis0,ois0;
341e30e807fSPeter Brune   DM             subdm0;
3420a696f66SPeter Brune 
3430adebc6cSBarry Smith   PetscFunctionBegin;
344e30e807fSPeter Brune   if (len) *len = 1;
345e30e807fSPeter Brune 
346e30e807fSPeter Brune   if (iis) {ierr = PetscMalloc(sizeof(IS*),iis);CHKERRQ(ierr);}
347e30e807fSPeter Brune   if (ois) {ierr = PetscMalloc(sizeof(IS*),ois);CHKERRQ(ierr);}
348e30e807fSPeter Brune   if (subdm) {ierr = PetscMalloc(sizeof(DM*),subdm);CHKERRQ(ierr);}
349e30e807fSPeter Brune   if (names) {ierr = PetscMalloc(sizeof(char*),names);CHKERRQ(ierr);}
350e30e807fSPeter Brune   ierr = DMDASubDomainDA_Private(dm,&subdm0);CHKERRQ(ierr);
351e30e807fSPeter Brune   ierr = DMDASubDomainIS_Private(dm,subdm0,&iis0,&ois0);CHKERRQ(ierr);
3528865f1eaSKarl Rupp   if (iis) (*iis)[0] = iis0;
3538865f1eaSKarl Rupp   else {
354e30e807fSPeter Brune     ierr = ISDestroy(&iis0);CHKERRQ(ierr);
355e30e807fSPeter Brune   }
3568865f1eaSKarl Rupp   if (ois) (*ois)[0] = ois0;
3578865f1eaSKarl Rupp   else {
358e30e807fSPeter Brune     ierr = ISDestroy(&ois0);CHKERRQ(ierr);
359e30e807fSPeter Brune   }
360e30e807fSPeter Brune   if (subdm) (*subdm)[0] = subdm0;
361e30e807fSPeter Brune   if (names) (*names)[0] = 0;
362e30e807fSPeter Brune   PetscFunctionReturn(0);
363e30e807fSPeter Brune }
364