xref: /petsc/src/dm/impls/da/dadd.c (revision e30e807f52449ee56cb8e6e4bd3267d3371e2f11)
1*e30e807fSPeter Brune #include <petsc-private/daimpl.h>
2*e30e807fSPeter Brune 
3*e30e807fSPeter Brune #undef __FUNCT__
4*e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainDA_Private"
5*e30e807fSPeter Brune PetscErrorCode DMDASubDomainDA_Private(DM dm, DM *dddm) {
6*e30e807fSPeter Brune   DM               da;
7*e30e807fSPeter Brune   DM_DA            *dd;
8*e30e807fSPeter Brune   PetscErrorCode   ierr;
9*e30e807fSPeter Brune   DMDALocalInfo    info;
10*e30e807fSPeter Brune   PetscReal        lmin[3],lmax[3];
11*e30e807fSPeter Brune   PetscInt         xsize,ysize,zsize;
12*e30e807fSPeter Brune   PetscInt         xo,yo,zo;
13*e30e807fSPeter Brune 
14*e30e807fSPeter Brune   PetscFunctionBegin;
15*e30e807fSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
16*e30e807fSPeter Brune   ierr = DMDACreate(PETSC_COMM_SELF,&da);CHKERRQ(ierr);
17*e30e807fSPeter Brune   ierr = DMSetOptionsPrefix(da,"sub_");CHKERRQ(ierr);
18*e30e807fSPeter Brune 
19*e30e807fSPeter Brune   ierr = DMDASetDim(da, info.dim);CHKERRQ(ierr);
20*e30e807fSPeter Brune   ierr = DMDASetDof(da, info.dof);CHKERRQ(ierr);
21*e30e807fSPeter Brune 
22*e30e807fSPeter Brune   ierr = DMDASetStencilType(da,info.st);CHKERRQ(ierr);
23*e30e807fSPeter Brune   ierr = DMDASetStencilWidth(da,info.sw);CHKERRQ(ierr);
24*e30e807fSPeter Brune 
25*e30e807fSPeter Brune   dd = (DM_DA *)dm->data;
26*e30e807fSPeter Brune 
27*e30e807fSPeter Brune   /* local with overlap */
28*e30e807fSPeter Brune   xsize = info.xm;
29*e30e807fSPeter Brune   ysize = info.ym;
30*e30e807fSPeter Brune   zsize = info.zm;
31*e30e807fSPeter Brune   xo    = info.xs;
32*e30e807fSPeter Brune   yo    = info.ys;
33*e30e807fSPeter Brune   zo    = info.zs;
34*e30e807fSPeter Brune   if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs != 0)) {
35*e30e807fSPeter Brune     xsize += dd->overlap;
36*e30e807fSPeter Brune     xo -= dd->overlap;
37*e30e807fSPeter Brune   }
38*e30e807fSPeter Brune   if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys != 0)) {
39*e30e807fSPeter Brune     ysize += dd->overlap;
40*e30e807fSPeter Brune     yo    -= dd->overlap;
41*e30e807fSPeter Brune   }
42*e30e807fSPeter Brune   if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs != 0)) {
43*e30e807fSPeter Brune     zsize += dd->overlap;
44*e30e807fSPeter Brune     zo    -= dd->overlap;
45*e30e807fSPeter Brune   }
46*e30e807fSPeter Brune 
47*e30e807fSPeter Brune   if (info.bx == DMDA_BOUNDARY_PERIODIC || (info.xs+info.xm != info.mx)) {
48*e30e807fSPeter Brune     xsize += dd->overlap;
49*e30e807fSPeter Brune   }
50*e30e807fSPeter Brune   if (info.by == DMDA_BOUNDARY_PERIODIC || (info.ys+info.ym != info.my)) {
51*e30e807fSPeter Brune     ysize += dd->overlap;
52*e30e807fSPeter Brune   }
53*e30e807fSPeter Brune   if (info.bz == DMDA_BOUNDARY_PERIODIC || (info.zs+info.zm != info.mz)) {
54*e30e807fSPeter Brune     zsize += dd->overlap;
55*e30e807fSPeter Brune   }
56*e30e807fSPeter Brune 
57*e30e807fSPeter Brune   ierr = DMDASetSizes(da, xsize, ysize, zsize);CHKERRQ(ierr);
58*e30e807fSPeter Brune   ierr = DMDASetNumProcs(da, 1, 1, 1);CHKERRQ(ierr);
59*e30e807fSPeter Brune   ierr = DMDASetBoundaryType(da, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_GHOSTED);CHKERRQ(ierr);
60*e30e807fSPeter Brune 
61*e30e807fSPeter Brune   /* set up as a block instead */
62*e30e807fSPeter Brune   ierr = DMSetUp(da);CHKERRQ(ierr);
63*e30e807fSPeter Brune   ierr = DMDASetOffset(da,xo,yo,zo);CHKERRQ(ierr);
64*e30e807fSPeter Brune 
65*e30e807fSPeter Brune 
66*e30e807fSPeter Brune   /* todo - nonuniform coordinates */
67*e30e807fSPeter Brune   ierr = DMDAGetLocalBoundingBox(dm,lmin,lmax);CHKERRQ(ierr);
68*e30e807fSPeter Brune   ierr = DMDASetUniformCoordinates(da,lmin[0],lmax[0],lmin[1],lmax[1],lmin[2],lmax[2]);CHKERRQ(ierr);
69*e30e807fSPeter Brune 
70*e30e807fSPeter Brune   dd = (DM_DA *)da->data;
71*e30e807fSPeter Brune   dd->Mo = info.mx;
72*e30e807fSPeter Brune   dd->No = info.my;
73*e30e807fSPeter Brune   dd->Po = info.mz;
74*e30e807fSPeter Brune 
75*e30e807fSPeter Brune   *dddm = da;
76*e30e807fSPeter Brune 
77*e30e807fSPeter Brune   PetscFunctionReturn(0);
78*e30e807fSPeter Brune }
79*e30e807fSPeter Brune 
80*e30e807fSPeter Brune #undef __FUNCT__
81*e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecompositionScatters_DA"
82*e30e807fSPeter Brune /*
83*e30e807fSPeter Brune  Fills the local vector problem on the subdomain from the global problem.
84*e30e807fSPeter Brune 
85*e30e807fSPeter Brune  */
86*e30e807fSPeter Brune PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM dm,PetscInt nsubdms,DM *subdms,VecScatter **iscat,VecScatter **oscat, VecScatter **lscat) {
87*e30e807fSPeter Brune   PetscErrorCode   ierr;
88*e30e807fSPeter Brune   DMDALocalInfo    dinfo,sinfo;
89*e30e807fSPeter Brune   IS               isis,idis,osis,odis,gsis,gdis;
90*e30e807fSPeter Brune   PetscInt         *ididx,*isidx,*odidx,*osidx,*gdidx,*gsidx,*idx_global,n_global,*idx_sub,n_sub;
91*e30e807fSPeter Brune   PetscInt         l,i,j,k,d,n_i,n_o,n_g,sl,dl,di,dj,dk,si,sj,sk;
92*e30e807fSPeter Brune   Vec              dvec,svec,slvec;
93*e30e807fSPeter Brune   DM               subdm;
94*e30e807fSPeter Brune 
95*e30e807fSPeter Brune   PetscFunctionBegin;
96*e30e807fSPeter Brune 
97*e30e807fSPeter Brune   /* allocate the arrays of scatters */
98*e30e807fSPeter Brune   if (iscat) ierr = PetscMalloc(sizeof(VecScatter *),iscat);CHKERRQ(ierr);
99*e30e807fSPeter Brune   if (oscat) ierr = PetscMalloc(sizeof(VecScatter *),oscat);CHKERRQ(ierr);
100*e30e807fSPeter Brune   if (lscat) ierr = PetscMalloc(sizeof(VecScatter *),lscat);CHKERRQ(ierr);
101*e30e807fSPeter Brune 
102*e30e807fSPeter Brune   ierr = DMDAGetLocalInfo(dm,&dinfo);CHKERRQ(ierr);
103*e30e807fSPeter Brune   ierr = DMDAGetGlobalIndices(dm,&n_global,&idx_global);CHKERRQ(ierr);
104*e30e807fSPeter Brune   for (l = 0;l < nsubdms;l++) {
105*e30e807fSPeter Brune     n_i = 0;
106*e30e807fSPeter Brune     n_o = 0;
107*e30e807fSPeter Brune     n_g = 0;
108*e30e807fSPeter Brune     subdm = subdms[l];
109*e30e807fSPeter Brune     ierr = DMDAGetLocalInfo(subdm,&sinfo);CHKERRQ(ierr);
110*e30e807fSPeter Brune     ierr = DMDAGetGlobalIndices(subdm,&n_sub,&idx_sub);CHKERRQ(ierr);
111*e30e807fSPeter Brune     /* count the three region sizes */
112*e30e807fSPeter Brune     for (k=sinfo.gzs;k<sinfo.gzs+sinfo.gzm;k++) {
113*e30e807fSPeter Brune       for (j=sinfo.gys;j<sinfo.gys+sinfo.gym;j++) {
114*e30e807fSPeter Brune         for (i=sinfo.gxs;i<sinfo.gxs+sinfo.gxm;i++) {
115*e30e807fSPeter Brune           for (d=0;d<sinfo.dof;d++) {
116*e30e807fSPeter Brune             if (k >= sinfo.zs           && j >= sinfo.ys         && i >= sinfo.xs &&
117*e30e807fSPeter Brune                 k <  sinfo.zs+sinfo.zm  && j < sinfo.ys+sinfo.ym && i < sinfo.xs+sinfo.xm) {
118*e30e807fSPeter Brune 
119*e30e807fSPeter Brune               /* interior - subinterior overlap */
120*e30e807fSPeter Brune               if (k >= dinfo.zs            && j >= dinfo.ys          && i >= dinfo.xs &&
121*e30e807fSPeter Brune                   k <  dinfo.zs+dinfo.zm  && j < dinfo.ys+dinfo.ym && i < dinfo.xs+dinfo.xm) {
122*e30e807fSPeter Brune                 n_i++;
123*e30e807fSPeter Brune               }
124*e30e807fSPeter Brune               /* ghost - subinterior overlap */
125*e30e807fSPeter Brune               if (k >= dinfo.gzs            && j >= dinfo.gys          && i >= dinfo.gxs &&
126*e30e807fSPeter Brune                   k <  dinfo.gzs+dinfo.gzm  && j < dinfo.gys+dinfo.gym && i < dinfo.gxs+dinfo.gxm) {
127*e30e807fSPeter Brune                 n_o++;
128*e30e807fSPeter Brune               }
129*e30e807fSPeter Brune             }
130*e30e807fSPeter Brune 
131*e30e807fSPeter Brune             /* ghost - subghost overlap */
132*e30e807fSPeter Brune             if (k >= dinfo.gzs            && j >= dinfo.gys          && i >= dinfo.gxs &&
133*e30e807fSPeter Brune                 k <  dinfo.gzs+dinfo.gzm  && j < dinfo.gys+dinfo.gym && i < dinfo.gxs+dinfo.gxm) {
134*e30e807fSPeter Brune               n_g++;
135*e30e807fSPeter Brune             }
136*e30e807fSPeter Brune           }
137*e30e807fSPeter Brune         }
138*e30e807fSPeter Brune       }
139*e30e807fSPeter Brune     }
140*e30e807fSPeter Brune 
141*e30e807fSPeter Brune     if (n_g == 0) SETERRQ(((PetscObject)subdm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Processor-local domain and subdomain do not intersect!");
142*e30e807fSPeter Brune 
143*e30e807fSPeter Brune     /* local and subdomain local index set indices */
144*e30e807fSPeter Brune     ierr = PetscMalloc(sizeof(PetscInt)*n_i,&ididx);CHKERRQ(ierr);
145*e30e807fSPeter Brune     ierr = PetscMalloc(sizeof(PetscInt)*n_i,&isidx);CHKERRQ(ierr);
146*e30e807fSPeter Brune 
147*e30e807fSPeter Brune     ierr = PetscMalloc(sizeof(PetscInt)*n_o,&odidx);CHKERRQ(ierr);
148*e30e807fSPeter Brune     ierr = PetscMalloc(sizeof(PetscInt)*n_o,&osidx);CHKERRQ(ierr);
149*e30e807fSPeter Brune 
150*e30e807fSPeter Brune     ierr = PetscMalloc(sizeof(PetscInt)*n_g,&gdidx);CHKERRQ(ierr);
151*e30e807fSPeter Brune     ierr = PetscMalloc(sizeof(PetscInt)*n_g,&gsidx);CHKERRQ(ierr);
152*e30e807fSPeter Brune 
153*e30e807fSPeter Brune     n_i = 0; n_o = 0;n_g = 0;
154*e30e807fSPeter Brune     for (k=sinfo.gzs;k<sinfo.gzs+sinfo.gzm;k++) {
155*e30e807fSPeter Brune       for (j=sinfo.gys;j<sinfo.gys+sinfo.gym;j++) {
156*e30e807fSPeter Brune         for (i=sinfo.gxs;i<sinfo.gxs+sinfo.gxm;i++) {
157*e30e807fSPeter Brune           for (d=0;d<sinfo.dof;d++) {
158*e30e807fSPeter Brune             si = i - sinfo.gxs;
159*e30e807fSPeter Brune             sj = j - sinfo.gys;
160*e30e807fSPeter Brune             sk = k - sinfo.gzs;
161*e30e807fSPeter Brune             sl = d + sinfo.dof*(si + sinfo.gxm*(sj + sinfo.gym*sk));
162*e30e807fSPeter Brune             di = i - dinfo.gxs;
163*e30e807fSPeter Brune             dj = j - dinfo.gys;
164*e30e807fSPeter Brune             dk = k - dinfo.gzs;
165*e30e807fSPeter Brune             dl = d + dinfo.dof*(di + dinfo.gxm*(dj + dinfo.gym*dk));
166*e30e807fSPeter Brune 
167*e30e807fSPeter Brune             if (k >= sinfo.zs           && j >= sinfo.ys         && i >= sinfo.xs &&
168*e30e807fSPeter Brune                 k <  sinfo.zs+sinfo.zm  && j < sinfo.ys+sinfo.ym && i < sinfo.xs+sinfo.xm) {
169*e30e807fSPeter Brune 
170*e30e807fSPeter Brune               /* interior - subinterior overlap */
171*e30e807fSPeter Brune               if (k >= dinfo.zs            && j >= dinfo.ys          && i >= dinfo.xs &&
172*e30e807fSPeter Brune                   k <  dinfo.zs+dinfo.zm  && j < dinfo.ys+dinfo.ym && i < dinfo.xs+dinfo.xm) {
173*e30e807fSPeter Brune                 ididx[n_i] = idx_global[dl];
174*e30e807fSPeter Brune                 isidx[n_i] = idx_sub[sl];
175*e30e807fSPeter Brune                 n_i++;
176*e30e807fSPeter Brune               }
177*e30e807fSPeter Brune               /* ghost - subinterior overlap */
178*e30e807fSPeter Brune               if (k >= dinfo.gzs            && j >= dinfo.gys          && i >= dinfo.gxs &&
179*e30e807fSPeter Brune                   k <  dinfo.gzs+dinfo.gzm  && j < dinfo.gys+dinfo.gym && i < dinfo.gxs+dinfo.gxm) {
180*e30e807fSPeter Brune                 odidx[n_o] = idx_global[dl];
181*e30e807fSPeter Brune                 osidx[n_o] = idx_sub[sl];
182*e30e807fSPeter Brune                 n_o++;
183*e30e807fSPeter Brune               }
184*e30e807fSPeter Brune             }
185*e30e807fSPeter Brune 
186*e30e807fSPeter Brune             /* ghost - subghost overlap */
187*e30e807fSPeter Brune             if (k >= dinfo.gzs            && j >= dinfo.gys          && i >= dinfo.gxs &&
188*e30e807fSPeter Brune                 k <  dinfo.gzs+dinfo.gzm  && j < dinfo.gys+dinfo.gym && i < dinfo.gxs+dinfo.gxm) {
189*e30e807fSPeter Brune               gdidx[n_g] = idx_global[dl];
190*e30e807fSPeter Brune               gsidx[n_g] = sl;
191*e30e807fSPeter Brune               n_g++;
192*e30e807fSPeter Brune             }
193*e30e807fSPeter Brune           }
194*e30e807fSPeter Brune         }
195*e30e807fSPeter Brune       }
196*e30e807fSPeter Brune     }
197*e30e807fSPeter Brune 
198*e30e807fSPeter Brune     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_i,ididx,PETSC_OWN_POINTER,&idis);CHKERRQ(ierr);
199*e30e807fSPeter Brune     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_i,isidx,PETSC_OWN_POINTER,&isis);CHKERRQ(ierr);
200*e30e807fSPeter Brune 
201*e30e807fSPeter Brune     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_o,odidx,PETSC_OWN_POINTER,&odis);CHKERRQ(ierr);
202*e30e807fSPeter Brune     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_o,osidx,PETSC_OWN_POINTER,&osis);CHKERRQ(ierr);
203*e30e807fSPeter Brune 
204*e30e807fSPeter Brune     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_g,gdidx,PETSC_OWN_POINTER,&gdis);CHKERRQ(ierr);
205*e30e807fSPeter Brune     ierr = ISCreateGeneral(PETSC_COMM_SELF,n_g,gsidx,PETSC_OWN_POINTER,&gsis);CHKERRQ(ierr);
206*e30e807fSPeter Brune 
207*e30e807fSPeter Brune     /* form the scatter */
208*e30e807fSPeter Brune     ierr = DMGetGlobalVector(dm,&dvec);CHKERRQ(ierr);
209*e30e807fSPeter Brune     ierr = DMGetGlobalVector(subdm,&svec);CHKERRQ(ierr);
210*e30e807fSPeter Brune     ierr = DMGetLocalVector(subdm,&slvec);CHKERRQ(ierr);
211*e30e807fSPeter Brune 
212*e30e807fSPeter Brune     if (iscat) {ierr = VecScatterCreate(dvec,idis,svec,isis,&(*iscat)[l]);CHKERRQ(ierr);}
213*e30e807fSPeter Brune     if (oscat) {ierr = VecScatterCreate(dvec,odis,svec,osis,&(*oscat)[l]);CHKERRQ(ierr);}
214*e30e807fSPeter Brune     if (lscat) {ierr = VecScatterCreate(dvec,gdis,slvec,gsis,&(*lscat)[l]);CHKERRQ(ierr);}
215*e30e807fSPeter Brune 
216*e30e807fSPeter Brune     ierr = DMRestoreGlobalVector(dm,&dvec);CHKERRQ(ierr);
217*e30e807fSPeter Brune     ierr = DMRestoreGlobalVector(subdm,&svec);CHKERRQ(ierr);
218*e30e807fSPeter Brune     ierr = DMRestoreLocalVector(subdm,&slvec);CHKERRQ(ierr);
219*e30e807fSPeter Brune 
220*e30e807fSPeter Brune     ierr = ISDestroy(&idis);CHKERRQ(ierr);
221*e30e807fSPeter Brune     ierr = ISDestroy(&isis);CHKERRQ(ierr);
222*e30e807fSPeter Brune 
223*e30e807fSPeter Brune     ierr = ISDestroy(&odis);CHKERRQ(ierr);
224*e30e807fSPeter Brune     ierr = ISDestroy(&osis);CHKERRQ(ierr);
225*e30e807fSPeter Brune 
226*e30e807fSPeter Brune     ierr = ISDestroy(&gdis);CHKERRQ(ierr);
227*e30e807fSPeter Brune     ierr = ISDestroy(&gsis);CHKERRQ(ierr);
228*e30e807fSPeter Brune   }
229*e30e807fSPeter Brune   PetscFunctionReturn(0);
230*e30e807fSPeter Brune 
231*e30e807fSPeter Brune }
232*e30e807fSPeter Brune 
233*e30e807fSPeter Brune #undef __FUNCT__
234*e30e807fSPeter Brune #define __FUNCT__ "DMDASubDomainIS_Private"
235*e30e807fSPeter Brune /* We have that the interior regions are going to be the same, but the ghost regions might not match up
236*e30e807fSPeter Brune 
237*e30e807fSPeter Brune ----------
238*e30e807fSPeter Brune ----------
239*e30e807fSPeter Brune --++++++o=
240*e30e807fSPeter Brune --++++++o=
241*e30e807fSPeter Brune --++++++o=
242*e30e807fSPeter Brune --++++++o=
243*e30e807fSPeter Brune --ooooooo=
244*e30e807fSPeter Brune --========
245*e30e807fSPeter Brune 
246*e30e807fSPeter Brune Therefore, for each point in the overall, we must check if it's:
247*e30e807fSPeter Brune 
248*e30e807fSPeter Brune 1. +: In the interior of the global dm; it lines up
249*e30e807fSPeter Brune 2. o: In the overlap region -- for now the same as 1; no overlap
250*e30e807fSPeter Brune 3. =: In the shared ghost region -- handled by DMCreateDomainDecompositionLocalScatter()
251*e30e807fSPeter Brune 4. -: In the nonshared ghost region
252*e30e807fSPeter Brune  */
253*e30e807fSPeter Brune 
254*e30e807fSPeter Brune PetscErrorCode DMDASubDomainIS_Private(DM dm,DM subdm,IS *iis,IS *ois) {
255*e30e807fSPeter Brune   PetscErrorCode   ierr;
256*e30e807fSPeter Brune   DMDALocalInfo    info,subinfo;
257*e30e807fSPeter Brune   PetscInt         *iiidx,*oiidx,*gidx,gindx;
258*e30e807fSPeter Brune   PetscInt         i,j,k,d,n,nsub,nover,llindx,lindx,li,lj,lk,gi,gj,gk;
259*e30e807fSPeter Brune 
260*e30e807fSPeter Brune   PetscFunctionBegin;
261*e30e807fSPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
262*e30e807fSPeter Brune   ierr = DMDAGetLocalInfo(subdm,&subinfo);CHKERRQ(ierr);
263*e30e807fSPeter Brune   ierr = DMDAGetGlobalIndices(dm,&n,&gidx);CHKERRQ(ierr);
264*e30e807fSPeter Brune   /* todo -- overlap */
265*e30e807fSPeter Brune   nsub = info.xm*info.ym*info.zm*info.dof;
266*e30e807fSPeter Brune   nover = subinfo.xm*subinfo.ym*subinfo.zm*subinfo.dof;
267*e30e807fSPeter Brune   /* iis is going to have size of the local problem's global part but have a lot of fill-in */
268*e30e807fSPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*nsub,&iiidx);CHKERRQ(ierr);
269*e30e807fSPeter Brune   /* ois is going to have size of the local problem's global part */
270*e30e807fSPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*nover,&oiidx);CHKERRQ(ierr);
271*e30e807fSPeter Brune   /* loop over the ghost region of the subdm and fill in the indices */
272*e30e807fSPeter Brune   for (k=subinfo.gzs;k<subinfo.gzs+subinfo.gzm;k++) {
273*e30e807fSPeter Brune     for (j=subinfo.gys;j<subinfo.gys+subinfo.gym;j++) {
274*e30e807fSPeter Brune       for (i=subinfo.gxs;i<subinfo.gxs+subinfo.gxm;i++) {
275*e30e807fSPeter Brune         for (d=0;d<subinfo.dof;d++) {
276*e30e807fSPeter Brune           li = i - subinfo.xs;
277*e30e807fSPeter Brune           lj = j - subinfo.ys;
278*e30e807fSPeter Brune           lk = k - subinfo.zs;
279*e30e807fSPeter Brune           lindx = d + subinfo.dof*(li + subinfo.xm*(lj + subinfo.ym*lk));
280*e30e807fSPeter Brune           li = i - info.xs;
281*e30e807fSPeter Brune           lj = j - info.ys;
282*e30e807fSPeter Brune           lk = k - info.zs;
283*e30e807fSPeter Brune           llindx = d + info.dof*(li + info.xm*(lj + info.ym*lk));
284*e30e807fSPeter Brune           gi = i - info.gxs;
285*e30e807fSPeter Brune           gj = j - info.gys;
286*e30e807fSPeter Brune           gk = k - info.gzs;
287*e30e807fSPeter Brune           gindx = d + info.dof*(gi + info.gxm*(gj + info.gym*gk));
288*e30e807fSPeter Brune 
289*e30e807fSPeter Brune           /* check if the current point is inside the interior region */
290*e30e807fSPeter Brune           if (k >= info.zs          && j >= info.ys          && i >= info.xs &&
291*e30e807fSPeter Brune               k <  info.zs+info.zm  && j < info.ys+info.ym   && i < info.xs+info.xm) {
292*e30e807fSPeter Brune             iiidx[llindx] = gidx[gindx];
293*e30e807fSPeter Brune             oiidx[lindx] = gidx[gindx];
294*e30e807fSPeter Brune             /* overlap region */
295*e30e807fSPeter Brune           } else if (k >= subinfo.zs             && j >= subinfo.ys                && i >= subinfo.xs &&
296*e30e807fSPeter Brune                      k <  subinfo.zs+subinfo.zm  && j < subinfo.ys+subinfo.ym   && i < subinfo.xs+subinfo.xm) {
297*e30e807fSPeter Brune             oiidx[lindx] = gidx[gindx];
298*e30e807fSPeter Brune           }
299*e30e807fSPeter Brune         }
300*e30e807fSPeter Brune       }
301*e30e807fSPeter Brune     }
302*e30e807fSPeter Brune   }
303*e30e807fSPeter Brune 
304*e30e807fSPeter Brune   /* create the index sets */
305*e30e807fSPeter Brune   ierr = ISCreateGeneral(PETSC_COMM_SELF,nsub,iiidx,PETSC_OWN_POINTER,iis);CHKERRQ(ierr);
306*e30e807fSPeter Brune   ierr = ISCreateGeneral(PETSC_COMM_SELF,nover,oiidx,PETSC_OWN_POINTER,ois);CHKERRQ(ierr);
307*e30e807fSPeter Brune   PetscFunctionReturn(0);
308*e30e807fSPeter Brune }
309*e30e807fSPeter Brune 
310*e30e807fSPeter Brune #undef __FUNCT__
311*e30e807fSPeter Brune #define __FUNCT__ "DMCreateDomainDecomposition_DA"
312*e30e807fSPeter Brune PetscErrorCode DMCreateDomainDecomposition_DA(DM dm,PetscInt *len,char ***names,IS **iis,IS **ois,DM **subdm) {
313*e30e807fSPeter Brune   PetscErrorCode ierr;
314*e30e807fSPeter Brune   IS             iis0,ois0;
315*e30e807fSPeter Brune   DM             subdm0;
316*e30e807fSPeter Brune   PetscFunctionBegin;
317*e30e807fSPeter Brune   if (len)*len = 1;
318*e30e807fSPeter Brune 
319*e30e807fSPeter Brune   if (iis) {ierr = PetscMalloc(sizeof(IS *),iis);CHKERRQ(ierr);}
320*e30e807fSPeter Brune   if (ois) {ierr = PetscMalloc(sizeof(IS *),ois);CHKERRQ(ierr);}
321*e30e807fSPeter Brune   if (subdm) {ierr = PetscMalloc(sizeof(DM *),subdm);CHKERRQ(ierr);}
322*e30e807fSPeter Brune   if (names) {ierr = PetscMalloc(sizeof(char *),names);CHKERRQ(ierr);}
323*e30e807fSPeter Brune   ierr = DMDASubDomainDA_Private(dm,&subdm0);CHKERRQ(ierr);
324*e30e807fSPeter Brune   ierr = DMDASubDomainIS_Private(dm,subdm0,&iis0,&ois0);CHKERRQ(ierr);
325*e30e807fSPeter Brune   if (iis) {
326*e30e807fSPeter Brune     (*iis)[0] = iis0;
327*e30e807fSPeter Brune   } else {
328*e30e807fSPeter Brune     ierr = ISDestroy(&iis0);CHKERRQ(ierr);
329*e30e807fSPeter Brune   }
330*e30e807fSPeter Brune   if (ois) {
331*e30e807fSPeter Brune     (*ois)[0] = ois0;
332*e30e807fSPeter Brune   } else {
333*e30e807fSPeter Brune     ierr = ISDestroy(&ois0);CHKERRQ(ierr);
334*e30e807fSPeter Brune   }
335*e30e807fSPeter Brune   if (subdm) (*subdm)[0] = subdm0;
336*e30e807fSPeter Brune   if (names) (*names)[0] = 0;
337*e30e807fSPeter Brune   PetscFunctionReturn(0);
338*e30e807fSPeter Brune }
339