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