xref: /petsc/src/dm/impls/da/dadd.c (revision 285ae30523e893cd524decdd3d5eaaab4faf2787) !
1e30e807fSPeter Brune #include <petsc-private/daimpl.h>
2e30e807fSPeter Brune 
3e30e807fSPeter Brune #undef __FUNCT__
4ff9846d9SPeter Brune #define __FUNCT__ "DMDACreatePatchIS"
5ff9846d9SPeter Brune 
6ff9846d9SPeter Brune PetscErrorCode DMDACreatePatchIS(DM da,MatStencil *lower,MatStencil *upper,IS *is)
7ff9846d9SPeter Brune {
8ff9846d9SPeter Brune   PetscErrorCode ierr;
9ff9846d9SPeter Brune   PetscInt       i,j,k,idx;
10ff9846d9SPeter Brune   PetscInt       ii,jj,kk;
11ff9846d9SPeter Brune   Vec            v;
12ff9846d9SPeter Brune   PetscInt       n,pn,bs;
13ff9846d9SPeter Brune   PetscMPIInt    rank;
14ff9846d9SPeter Brune   PetscSF        sf,psf;
15ff9846d9SPeter Brune   PetscLayout    map;
16ff9846d9SPeter Brune   MPI_Comm       comm;
17ff9846d9SPeter Brune   PetscInt       *natidx,*globidx,*leafidx;
18ff9846d9SPeter Brune   PetscInt       *pnatidx,*pleafidx;
19ff9846d9SPeter Brune   PetscInt       base;
20ff9846d9SPeter Brune   PetscInt       ox,oy,oz;
21ff9846d9SPeter Brune   DM_DA          *dd;
22ff9846d9SPeter Brune   PetscFunctionBegin;
23ff9846d9SPeter Brune 
24ff9846d9SPeter Brune   comm = ((PetscObject)da)->comm;
25ff9846d9SPeter Brune   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
26ff9846d9SPeter Brune   dd = (DM_DA *)da->data;
27ff9846d9SPeter Brune 
28ff9846d9SPeter Brune   /* construct the natural mapping */
29ff9846d9SPeter Brune   ierr = DMGetGlobalVector(da,&v);CHKERRQ(ierr);
30ff9846d9SPeter Brune   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
31ff9846d9SPeter Brune   ierr = VecGetOwnershipRange(v,&base,PETSC_NULL);CHKERRQ(ierr);
32ff9846d9SPeter Brune   ierr = VecGetBlockSize(v,&bs);CHKERRQ(ierr);
33ff9846d9SPeter Brune   ierr = DMRestoreGlobalVector(da,&v);CHKERRQ(ierr);
34ff9846d9SPeter Brune 
35ff9846d9SPeter Brune   /* construct the layout */
36ff9846d9SPeter Brune   ierr = PetscLayoutCreate(comm,&map);CHKERRQ(ierr);
37ff9846d9SPeter Brune   ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
38ff9846d9SPeter Brune   ierr = PetscLayoutSetLocalSize(map,n);CHKERRQ(ierr);
39ff9846d9SPeter Brune   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
40ff9846d9SPeter Brune 
41ff9846d9SPeter Brune   /* construct the list of natural indices on this process when PETSc ordering is considered */
42ff9846d9SPeter Brune   ierr = DMDAGetOffset(da,&ox,&oy,&oz);CHKERRQ(ierr);
43ff9846d9SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*n,&natidx);CHKERRQ(ierr);
44ff9846d9SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*n,&globidx);CHKERRQ(ierr);
45ff9846d9SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*n,&leafidx);CHKERRQ(ierr);
46ff9846d9SPeter Brune   idx = 0;
47ff9846d9SPeter Brune   for (k=dd->zs;k<dd->ze;k++) {
48ff9846d9SPeter Brune     for (j=dd->ys;j<dd->ye;j++) {
49ff9846d9SPeter Brune       for (i=dd->xs;i<dd->xe;i++) {
50ff9846d9SPeter Brune         natidx[idx] = i + dd->w*(j*dd->M + k*dd->M*dd->N);
51ff9846d9SPeter Brune         globidx[idx] = base + idx;
52ff9846d9SPeter Brune         leafidx[idx] = 0;
53ff9846d9SPeter Brune         idx++;
54ff9846d9SPeter Brune       }
55ff9846d9SPeter Brune     }
56ff9846d9SPeter Brune   }
57ff9846d9SPeter Brune 
58ff9846d9SPeter Brune   if (idx != n) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE, "for some reason the count is wrong.");
59ff9846d9SPeter Brune 
60ff9846d9SPeter Brune   /* construct the SF going from the natural indices to the local set of PETSc indices */
61ff9846d9SPeter Brune   ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr);
62ff9846d9SPeter Brune   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
63ff9846d9SPeter Brune   ierr = PetscSFSetGraphLayout(sf,map,n,PETSC_NULL,PETSC_OWN_POINTER,natidx);CHKERRQ(ierr);
64ff9846d9SPeter Brune 
65ff9846d9SPeter Brune   /* broadcast the global indices over to the corresponding natural indices */
66ff9846d9SPeter Brune   ierr = PetscSFGatherBegin(sf,MPIU_INT,globidx,leafidx);CHKERRQ(ierr);
67ff9846d9SPeter Brune   ierr = PetscSFGatherEnd(sf,MPIU_INT,globidx,leafidx);CHKERRQ(ierr);
68ff9846d9SPeter Brune   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
69ff9846d9SPeter Brune 
70ff9846d9SPeter Brune 
71ff9846d9SPeter Brune   pn = dd->w*(upper->k - lower->k)*(upper->j - lower->j)*(upper->i - lower->i);
72ff9846d9SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*pn,&pnatidx);CHKERRQ(ierr);
73ff9846d9SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*pn,&pleafidx);CHKERRQ(ierr);
74ff9846d9SPeter Brune   idx = 0;
75ff9846d9SPeter Brune   for (k=lower->k-oz;k<upper->k-oz;k++) {
76ff9846d9SPeter Brune     for (j=lower->j-oy;j<upper->j-oy;j++) {
77ff9846d9SPeter Brune       for (i=dd->w*(lower->i-ox);i<dd->w*(upper->i-ox);i++) {
78ff9846d9SPeter Brune         ii = i % (dd->w*dd->M);
79ff9846d9SPeter Brune         jj = j % dd->N;
80ff9846d9SPeter Brune         kk = k % dd->P;
81ff9846d9SPeter Brune         if (ii < 0) ii = dd->w*dd->M + ii;
82ff9846d9SPeter Brune         if (jj < 0) jj = dd->N + jj;
83ff9846d9SPeter Brune         if (kk < 0) kk = dd->P + kk;
84ff9846d9SPeter Brune         pnatidx[idx] = ii + dd->w*(jj*dd->M + kk*dd->M*dd->N);
85ff9846d9SPeter Brune         idx++;
86ff9846d9SPeter Brune       }
87ff9846d9SPeter Brune     }
88ff9846d9SPeter Brune   }
89ff9846d9SPeter Brune 
90ff9846d9SPeter Brune   if (idx != pn) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE, "for some reason the count is wrong");
91ff9846d9SPeter Brune 
92ff9846d9SPeter Brune   ierr = PetscSFCreate(comm,&psf);CHKERRQ(ierr);
93ff9846d9SPeter Brune   ierr = PetscSFSetFromOptions(psf);CHKERRQ(ierr);
94ff9846d9SPeter Brune   ierr = PetscSFSetGraphLayout(psf,map,pn,PETSC_NULL,PETSC_OWN_POINTER,pnatidx);CHKERRQ(ierr);
95ff9846d9SPeter Brune 
96ff9846d9SPeter Brune   /* broadcast the global indices through to the patch */
97ff9846d9SPeter Brune   ierr = PetscSFBcastBegin(psf,MPIU_INT,leafidx,pleafidx);CHKERRQ(ierr);
98ff9846d9SPeter Brune   ierr = PetscSFBcastEnd(psf,MPIU_INT,leafidx,pleafidx);CHKERRQ(ierr);
99ff9846d9SPeter Brune 
100ff9846d9SPeter Brune   /* create the IS */
101ff9846d9SPeter Brune   ierr = ISCreateGeneral(comm,pn,pleafidx,PETSC_OWN_POINTER,is);CHKERRQ(ierr);
102ff9846d9SPeter Brune 
103ff9846d9SPeter Brune   ierr = PetscSFDestroy(&psf);CHKERRQ(ierr);
104ff9846d9SPeter Brune 
105ff9846d9SPeter Brune   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
106ff9846d9SPeter Brune 
107ff9846d9SPeter Brune   ierr = PetscFree(globidx);CHKERRQ(ierr);
108ff9846d9SPeter Brune   ierr = PetscFree(leafidx);CHKERRQ(ierr);
109ff9846d9SPeter Brune   ierr = PetscFree(natidx);CHKERRQ(ierr);
110ff9846d9SPeter Brune   ierr = PetscFree(pnatidx);CHKERRQ(ierr);
111ff9846d9SPeter Brune 
112ff9846d9SPeter Brune   PetscFunctionReturn(0);
113ff9846d9SPeter Brune }
114ff9846d9SPeter Brune 
115ff9846d9SPeter Brune #undef __FUNCT__
116e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainDA_Private"
117e30e807fSPeter Brune PetscErrorCode DMDASubDomainDA_Private(DM dm, DM *dddm) {
118e30e807fSPeter Brune   DM               da;
119e30e807fSPeter Brune   DM_DA            *dd;
120e30e807fSPeter Brune   PetscErrorCode   ierr;
121e30e807fSPeter Brune   DMDALocalInfo    info;
122e30e807fSPeter Brune   PetscReal        lmin[3],lmax[3];
123ff9846d9SPeter Brune   PetscInt         xsize,ysize,zsize;
124e30e807fSPeter Brune   PetscInt         xo,yo,zo;
125e30e807fSPeter Brune 
126e30e807fSPeter Brune   PetscFunctionBegin;
127e30e807fSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
128e30e807fSPeter Brune   ierr = DMDACreate(PETSC_COMM_SELF,&da);CHKERRQ(ierr);
129e30e807fSPeter Brune   ierr = DMSetOptionsPrefix(da,"sub_");CHKERRQ(ierr);
130e30e807fSPeter Brune 
131e30e807fSPeter Brune   ierr = DMDASetDim(da, info.dim);CHKERRQ(ierr);
132e30e807fSPeter Brune   ierr = DMDASetDof(da, info.dof);CHKERRQ(ierr);
133e30e807fSPeter Brune 
134e30e807fSPeter Brune   ierr = DMDASetStencilType(da,info.st);CHKERRQ(ierr);
135e30e807fSPeter Brune   ierr = DMDASetStencilWidth(da,info.sw);CHKERRQ(ierr);
136e30e807fSPeter Brune 
137e30e807fSPeter Brune   dd = (DM_DA *)dm->data;
138e30e807fSPeter Brune 
139e30e807fSPeter Brune   /* local with overlap */
140e30e807fSPeter Brune   xsize = info.xm;
141e30e807fSPeter Brune   ysize = info.ym;
142e30e807fSPeter Brune   zsize = info.zm;
143e30e807fSPeter Brune   xo    = info.xs;
144e30e807fSPeter Brune   yo    = info.ys;
145e30e807fSPeter Brune   zo    = info.zs;
146e30e807fSPeter Brune   if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs != 0)) {
147e30e807fSPeter Brune     xsize += dd->overlap;
148e30e807fSPeter Brune     xo -= dd->overlap;
149e30e807fSPeter Brune   }
150e30e807fSPeter Brune   if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys != 0)) {
151e30e807fSPeter Brune     ysize += dd->overlap;
152e30e807fSPeter Brune     yo    -= dd->overlap;
153e30e807fSPeter Brune   }
154e30e807fSPeter Brune   if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs != 0)) {
155e30e807fSPeter Brune     zsize += dd->overlap;
156e30e807fSPeter Brune     zo    -= dd->overlap;
157e30e807fSPeter Brune   }
158e30e807fSPeter Brune 
159e30e807fSPeter Brune   if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs+info.xm != info.mx)) {
160e30e807fSPeter Brune     xsize += dd->overlap;
161e30e807fSPeter Brune   }
162e30e807fSPeter Brune   if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys+info.ym != info.my)) {
163e30e807fSPeter Brune     ysize += dd->overlap;
164e30e807fSPeter Brune   }
165e30e807fSPeter Brune   if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs+info.zm != info.mz)) {
166e30e807fSPeter Brune     zsize += dd->overlap;
167e30e807fSPeter Brune   }
168e30e807fSPeter Brune 
169e30e807fSPeter Brune   ierr = DMDASetSizes(da, xsize, ysize, zsize);CHKERRQ(ierr);
170e30e807fSPeter Brune   ierr = DMDASetNumProcs(da, 1, 1, 1);CHKERRQ(ierr);
171e30e807fSPeter Brune   ierr = DMDASetBoundaryType(da, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED);CHKERRQ(ierr);
172e30e807fSPeter Brune 
173e30e807fSPeter Brune   /* set up as a block instead */
174e30e807fSPeter Brune   ierr = DMSetUp(da);CHKERRQ(ierr);
175e30e807fSPeter Brune   ierr = DMDASetOffset(da,xo,yo,zo);CHKERRQ(ierr);
176e30e807fSPeter Brune 
177e30e807fSPeter Brune 
178e30e807fSPeter Brune   /* todo - nonuniform coordinates */
179e30e807fSPeter Brune   ierr = DMDAGetLocalBoundingBox(dm,lmin,lmax);CHKERRQ(ierr);
180e30e807fSPeter Brune   ierr = DMDASetUniformCoordinates(da,lmin[0],lmax[0],lmin[1],lmax[1],lmin[2],lmax[2]);CHKERRQ(ierr);
181e30e807fSPeter Brune 
182e30e807fSPeter Brune   dd = (DM_DA *)da->data;
183e30e807fSPeter Brune   dd->Mo = info.mx;
184e30e807fSPeter Brune   dd->No = info.my;
185e30e807fSPeter Brune   dd->Po = info.mz;
186e30e807fSPeter Brune 
187e30e807fSPeter Brune   *dddm = da;
188e30e807fSPeter Brune 
189e30e807fSPeter Brune   PetscFunctionReturn(0);
190e30e807fSPeter Brune }
191e30e807fSPeter Brune 
192e30e807fSPeter Brune #undef __FUNCT__
193e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecompositionScatters_DA"
194e30e807fSPeter Brune /*
195e30e807fSPeter Brune  Fills the local vector problem on the subdomain from the global problem.
196e30e807fSPeter Brune 
197*285ae305SPeter Brune  Right now this assumes one subdomain per processor.
198*285ae305SPeter Brune 
199e30e807fSPeter Brune  */
200e30e807fSPeter Brune PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat) {
201e30e807fSPeter Brune   PetscErrorCode   ierr;
202*285ae305SPeter Brune   DMDALocalInfo    info,subinfo;
203e30e807fSPeter Brune   DM               subdm;
204*285ae305SPeter Brune   MatStencil       upper,lower;
205*285ae305SPeter Brune   IS               idis,isis,odis,osis,gdis;
206*285ae305SPeter Brune   Vec              svec,dvec,slvec;
207e30e807fSPeter Brune 
208e30e807fSPeter Brune   PetscFunctionBegin;
209e30e807fSPeter Brune 
210*285ae305SPeter Brune   if (nsubdms != 1) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot have more than one subdomain per processor (yet)");
211*285ae305SPeter Brune 
212e30e807fSPeter Brune   /* allocate the arrays of scatters */
213050c5e14SJed Brown   if (iscat) {ierr = PetscMalloc(sizeof(VecScatter *),iscat);CHKERRQ(ierr);}
214050c5e14SJed Brown   if (oscat) {ierr = PetscMalloc(sizeof(VecScatter *),oscat);CHKERRQ(ierr);}
215050c5e14SJed Brown   if (lscat) {ierr = PetscMalloc(sizeof(VecScatter *),lscat);CHKERRQ(ierr);}
216e30e807fSPeter Brune 
217*285ae305SPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
218*285ae305SPeter Brune   subdm = subdms[0];
219*285ae305SPeter Brune   ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr);
220e30e807fSPeter Brune 
221*285ae305SPeter Brune   /* create the global and subdomain index sets for the inner domain */
222*285ae305SPeter Brune   /* TODO - make this actually support multiple subdomains -- subdomain needs to provide where it's nonoverlapping portion belongs */
223*285ae305SPeter Brune   lower.i = info.xs;
224*285ae305SPeter Brune   lower.j = info.ys;
225*285ae305SPeter Brune   lower.k = info.zs;
226*285ae305SPeter Brune   upper.i = info.xs+info.xm;
227*285ae305SPeter Brune   upper.j = info.ys+info.ym;
228*285ae305SPeter Brune   upper.k = info.zs+info.zm;
229*285ae305SPeter Brune   ierr = DMDACreatePatchIS(dm,&lower,&upper,&idis);CHKERRQ(ierr);
230*285ae305SPeter Brune   ierr = DMDACreatePatchIS(subdm,&lower,&upper,&isis);CHKERRQ(ierr);
231e30e807fSPeter Brune 
232*285ae305SPeter Brune   /* create the global and subdomain index sets for the outer subdomain */
233*285ae305SPeter Brune   lower.i = subinfo.xs;
234*285ae305SPeter Brune   lower.j = subinfo.ys;
235*285ae305SPeter Brune   lower.k = subinfo.zs;
236*285ae305SPeter Brune   upper.i = subinfo.xs+subinfo.xm;
237*285ae305SPeter Brune   upper.j = subinfo.ys+subinfo.ym;
238*285ae305SPeter Brune   upper.k = subinfo.zs+subinfo.zm;
239*285ae305SPeter Brune   ierr = DMDACreatePatchIS(dm,&lower,&upper,&odis);CHKERRQ(ierr);
240*285ae305SPeter Brune   ierr = DMDACreatePatchIS(subdm,&lower,&upper,&osis);CHKERRQ(ierr);
241e30e807fSPeter Brune 
242*285ae305SPeter Brune   /* global and subdomain ISes for the local indices of the subdomain */
243*285ae305SPeter Brune   /* todo - make this not loop over at nonperiodic boundaries, which will be more involved */
244*285ae305SPeter Brune   lower.i = subinfo.gxs;
245*285ae305SPeter Brune   lower.j = subinfo.gys;
246*285ae305SPeter Brune   lower.k = subinfo.gzs;
247*285ae305SPeter Brune   upper.i = subinfo.gxs+subinfo.gxm;
248*285ae305SPeter Brune   upper.j = subinfo.gys+subinfo.gym;
249*285ae305SPeter Brune   upper.k = subinfo.gzs+subinfo.gzm;
250e30e807fSPeter Brune 
251*285ae305SPeter Brune   ierr = DMDACreatePatchIS(dm,&lower,&upper,&gdis);CHKERRQ(ierr);
252e30e807fSPeter Brune 
253e30e807fSPeter Brune     /* form the scatter */
254e30e807fSPeter Brune     ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr);
255e30e807fSPeter Brune     ierr = DMGetGlobalVector(subdm,&svec);CHKERRQ(ierr);
256e30e807fSPeter Brune     ierr = DMGetLocalVector(subdm,&slvec);CHKERRQ(ierr);
257e30e807fSPeter Brune 
258*285ae305SPeter Brune     if (iscat) {ierr = VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[0]);CHKERRQ(ierr);}
259*285ae305SPeter Brune     if (oscat) {ierr = VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[0]);CHKERRQ(ierr);}
260*285ae305SPeter Brune     if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,PETSC_NULL,&(*lscat)[0]);CHKERRQ(ierr);}
261e30e807fSPeter Brune 
262e30e807fSPeter Brune     ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr);
263e30e807fSPeter Brune     ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr);
264e30e807fSPeter Brune     ierr = DMRestoreLocalVector(subdm,&slvec);CHKERRQ(ierr);
265e30e807fSPeter Brune 
266e30e807fSPeter Brune     ierr = ISDestroy(&idis);CHKERRQ(ierr);
267e30e807fSPeter Brune     ierr = ISDestroy(&isis);CHKERRQ(ierr);
268e30e807fSPeter Brune 
269e30e807fSPeter Brune     ierr = ISDestroy(&odis);CHKERRQ(ierr);
270e30e807fSPeter Brune     ierr = ISDestroy(&osis);CHKERRQ(ierr);
271e30e807fSPeter Brune 
272e30e807fSPeter Brune     ierr = ISDestroy(&gdis);CHKERRQ(ierr);
273e30e807fSPeter Brune     PetscFunctionReturn(0);
274e30e807fSPeter Brune }
275e30e807fSPeter Brune 
276e30e807fSPeter Brune #undef __FUNCT__
277e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainIS_Private"
278e30e807fSPeter Brune /* We have that the interior regions are going to be the same, but the ghost regions might not match up
279e30e807fSPeter Brune 
280e30e807fSPeter Brune ----------
281e30e807fSPeter Brune ----------
282e30e807fSPeter Brune --++++++o=
283e30e807fSPeter Brune --++++++o=
284e30e807fSPeter Brune --++++++o=
285e30e807fSPeter Brune --++++++o=
286e30e807fSPeter Brune --ooooooo=
287e30e807fSPeter Brune --========
288e30e807fSPeter Brune 
289e30e807fSPeter Brune Therefore, for each point in the overall, we must check if it's:
290e30e807fSPeter Brune 
291e30e807fSPeter Brune 1. +: In the interior of the global dm; it lines up
292e30e807fSPeter Brune 2. o: In the overlap region -- for now the same as 1; no overlap
293e30e807fSPeter Brune 3. =: In the shared ghost region -- handled by DMCreateDomainDecompositionLocalScatter()
294e30e807fSPeter Brune 4. -: In the nonshared ghost region
295e30e807fSPeter Brune  */
296e30e807fSPeter Brune 
297e30e807fSPeter Brune PetscErrorCode DMDASubDomainIS_Private(DM dm,DM subdm,IS *iis,IS *ois) {
298e30e807fSPeter Brune   PetscErrorCode   ierr;
299e30e807fSPeter Brune   DMDALocalInfo    info,subinfo;
300*285ae305SPeter Brune   MatStencil       lower,upper;
301e30e807fSPeter Brune 
302e30e807fSPeter Brune   PetscFunctionBegin;
303e30e807fSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
304e30e807fSPeter Brune   ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr);
305e30e807fSPeter Brune 
306*285ae305SPeter Brune   /* create the inner IS */
307*285ae305SPeter Brune   lower.i = info.xs;
308*285ae305SPeter Brune   lower.j = info.ys;
309*285ae305SPeter Brune   lower.k = info.zs;
310*285ae305SPeter Brune   upper.i = info.xs+info.xm;
311*285ae305SPeter Brune   upper.j = info.ys+info.ym;
312*285ae305SPeter Brune   upper.k = info.zs+info.zm;
313e30e807fSPeter Brune 
314*285ae305SPeter Brune   ierr = DMDACreatePatchIS(dm,&lower,&upper,iis);CHKERRQ(ierr);
315*285ae305SPeter Brune 
316*285ae305SPeter Brune   /* create the outer IS */
317*285ae305SPeter Brune   lower.i = subinfo.xs;
318*285ae305SPeter Brune   lower.j = subinfo.ys;
319*285ae305SPeter Brune   lower.k = subinfo.zs;
320*285ae305SPeter Brune   upper.i = subinfo.xs+subinfo.xm;
321*285ae305SPeter Brune   upper.j = subinfo.ys+subinfo.ym;
322*285ae305SPeter Brune   upper.k = subinfo.zs+subinfo.zm;
323*285ae305SPeter Brune   ierr = DMDACreatePatchIS(dm,&lower,&upper,ois);CHKERRQ(ierr);
324*285ae305SPeter Brune 
325e30e807fSPeter Brune   PetscFunctionReturn(0);
326e30e807fSPeter Brune }
327e30e807fSPeter Brune 
328e30e807fSPeter Brune #undef __FUNCT__
329e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecomposition_DA"
330e30e807fSPeter Brune PetscErrorCode DMCreateDomainDecomposition_DA(DM dm,PetscInt *len,char ***names,IS **iis,IS **ois,DM **subdm) {
331e30e807fSPeter Brune   PetscErrorCode ierr;
332e30e807fSPeter Brune   IS             iis0,ois0;
333e30e807fSPeter Brune   DM             subdm0;
334e30e807fSPeter Brune   PetscFunctionBegin;
335e30e807fSPeter Brune   if (len)*len = 1;
336e30e807fSPeter Brune 
337e30e807fSPeter Brune   if (iis) {ierr = PetscMalloc(sizeof(IS *),iis);CHKERRQ(ierr);}
338e30e807fSPeter Brune   if (ois) {ierr = PetscMalloc(sizeof(IS *),ois);CHKERRQ(ierr);}
339e30e807fSPeter Brune   if (subdm) {ierr = PetscMalloc(sizeof(DM *),subdm);CHKERRQ(ierr);}
340e30e807fSPeter Brune   if (names) {ierr = PetscMalloc(sizeof(char *),names);CHKERRQ(ierr);}
341e30e807fSPeter Brune   ierr = DMDASubDomainDA_Private(dm,&subdm0);CHKERRQ(ierr);
342e30e807fSPeter Brune   ierr = DMDASubDomainIS_Private(dm,subdm0,&iis0,&ois0);CHKERRQ(ierr);
343e30e807fSPeter Brune   if (iis) {
344e30e807fSPeter Brune     (*iis)[0] = iis0;
345e30e807fSPeter Brune   } else {
346e30e807fSPeter Brune     ierr = ISDestroy(&iis0);CHKERRQ(ierr);
347e30e807fSPeter Brune   }
348e30e807fSPeter Brune   if (ois) {
349e30e807fSPeter Brune     (*ois)[0] = ois0;
350e30e807fSPeter Brune   } else {
351e30e807fSPeter Brune     ierr = ISDestroy(&ois0);CHKERRQ(ierr);
352e30e807fSPeter Brune   }
353e30e807fSPeter Brune   if (subdm) (*subdm)[0] = subdm0;
354e30e807fSPeter Brune   if (names) (*names)[0] = 0;
355e30e807fSPeter Brune   PetscFunctionReturn(0);
356e30e807fSPeter Brune }
357