xref: /petsc/src/dm/impls/da/dacorn.c (revision b41ce5d507ea9a58bfa83cf403107a702e77a67d)
1 
2 /*
3   Code for manipulating distributed regular arrays in parallel.
4 */
5 
6 #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
7 
8 PetscErrorCode DMCreateCoordinateDM_DA(DM dm, DM *cdm)
9 {
10   PetscErrorCode ierr;
11   PetscFunctionBegin;
12   ierr = DMDAGetReducedDMDA(dm,dm->dim,cdm);CHKERRQ(ierr);
13   PetscFunctionReturn(0);
14 }
15 
16 /*@C
17    DMDASetFieldName - Sets the names of individual field components in multicomponent
18    vectors associated with a DMDA.
19 
20    Not Collective
21 
22    Input Parameters:
23 +  da - the distributed array
24 .  nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the
25         number of degrees of freedom per node within the DMDA
26 -  names - the name of the field (component)
27 
28   Level: intermediate
29 
30 .keywords: distributed array, get, component name
31 
32 .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMDASetFieldNames()
33 @*/
34 PetscErrorCode  DMDASetFieldName(DM da,PetscInt nf,const char name[])
35 {
36   PetscErrorCode ierr;
37   DM_DA          *dd = (DM_DA*)da->data;
38 
39   PetscFunctionBegin;
40   PetscValidHeaderSpecific(da,DM_CLASSID,1);
41   if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf);
42   ierr = PetscFree(dd->fieldname[nf]);CHKERRQ(ierr);
43   ierr = PetscStrallocpy(name,&dd->fieldname[nf]);CHKERRQ(ierr);
44   PetscFunctionReturn(0);
45 }
46 
47 /*@C
48    DMDAGetFieldNames - Gets the name of each component in the vector associated with the DMDA
49 
50    Collective on TS
51 
52    Input Parameter:
53 .  dm - the DMDA object
54 
55    Output Parameter:
56 .  names - the names of the components, final string is NULL, will have the same number of entries as the dof used in creating the DMDA
57 
58    Level: intermediate
59 
60 .keywords: distributed array, get, component name
61 
62 .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMDASetFieldName(), DMDASetFieldNames()
63 @*/
64 PetscErrorCode  DMDAGetFieldNames(DM da,const char * const **names)
65 {
66   DM_DA             *dd = (DM_DA*)da->data;
67 
68   PetscFunctionBegin;
69   *names = (const char * const *) dd->fieldname;
70   PetscFunctionReturn(0);
71 }
72 
73 /*@C
74    DMDASetFieldNames - Sets the name of each component in the vector associated with the DMDA
75 
76    Collective on TS
77 
78    Input Parameters:
79 +  dm - the DMDA object
80 -  names - the names of the components, final string must be NULL, must have the same number of entries as the dof used in creating the DMDA
81 
82    Level: intermediate
83 
84 .keywords: distributed array, get, component name
85 
86 .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMDASetFieldName()
87 @*/
88 PetscErrorCode  DMDASetFieldNames(DM da,const char * const *names)
89 {
90   PetscErrorCode    ierr;
91   DM_DA             *dd = (DM_DA*)da->data;
92 
93   PetscFunctionBegin;
94   ierr = PetscStrArrayDestroy(&dd->fieldname);CHKERRQ(ierr);
95   ierr = PetscStrArrayallocpy(names,&dd->fieldname);CHKERRQ(ierr);
96   PetscFunctionReturn(0);
97 }
98 
99 /*@C
100    DMDAGetFieldName - Gets the names of individual field components in multicomponent
101    vectors associated with a DMDA.
102 
103    Not Collective
104 
105    Input Parameter:
106 +  da - the distributed array
107 -  nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the
108         number of degrees of freedom per node within the DMDA
109 
110    Output Parameter:
111 .  names - the name of the field (component)
112 
113   Level: intermediate
114 
115 .keywords: distributed array, get, component name
116 
117 .seealso: DMDASetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName()
118 @*/
119 PetscErrorCode  DMDAGetFieldName(DM da,PetscInt nf,const char **name)
120 {
121   DM_DA *dd = (DM_DA*)da->data;
122 
123   PetscFunctionBegin;
124   PetscValidHeaderSpecific(da,DM_CLASSID,1);
125   PetscValidPointer(name,3);
126   if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf);
127   *name = dd->fieldname[nf];
128   PetscFunctionReturn(0);
129 }
130 
131 /*@C
132    DMDASetCoordinateName - Sets the name of the coordinate directions associated with a DMDA, for example "x" or "y"
133 
134    Not Collective
135 
136    Input Parameters:
137 +  dm - the DM
138 .  nf - coordinate number for the DMDA (0, 1, ... dim-1),
139 -  name - the name of the coordinate
140 
141   Level: intermediate
142 
143 .keywords: distributed array, get, component name
144 
145 .seealso: DMDAGetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName()
146 @*/
147 PetscErrorCode DMDASetCoordinateName(DM dm,PetscInt nf,const char name[])
148 {
149   PetscErrorCode ierr;
150   DM_DA          *dd = (DM_DA*)dm->data;
151 
152   PetscFunctionBegin;
153   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
154   if (nf < 0 || nf >= dm->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf);
155   ierr = PetscFree(dd->coordinatename[nf]);CHKERRQ(ierr);
156   ierr = PetscStrallocpy(name,&dd->coordinatename[nf]);CHKERRQ(ierr);
157   PetscFunctionReturn(0);
158 }
159 
160 /*@C
161    DMDAGetCoordinateName - Gets the name of a coodinate direction associated with a DMDA.
162 
163    Not Collective
164 
165    Input Parameter:
166 +  dm - the DM
167 -  nf -  number for the DMDA (0, 1, ... dim-1)
168 
169    Output Parameter:
170 .  names - the name of the coordinate direction
171 
172   Level: intermediate
173 
174 .keywords: distributed array, get, component name
175 
176 .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName()
177 @*/
178 PetscErrorCode DMDAGetCoordinateName(DM dm,PetscInt nf,const char **name)
179 {
180   DM_DA *dd = (DM_DA*)dm->data;
181 
182   PetscFunctionBegin;
183   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
184   PetscValidPointer(name,3);
185   if (nf < 0 || nf >= dm->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf);
186   *name = dd->coordinatename[nf];
187   PetscFunctionReturn(0);
188 }
189 
190 /*@C
191    DMDAGetCorners - Returns the global (x,y,z) indices of the lower left
192    corner and size of the local region, excluding ghost points.
193 
194    Not Collective
195 
196    Input Parameter:
197 .  da - the distributed array
198 
199    Output Parameters:
200 +  x,y,z - the corner indices (where y and z are optional; these are used
201            for 2D and 3D problems)
202 -  m,n,p - widths in the corresponding directions (where n and p are optional;
203            these are used for 2D and 3D problems)
204 
205    Note:
206    The corner information is independent of the number of degrees of
207    freedom per node set with the DMDACreateXX() routine. Thus the x, y, z, and
208    m, n, p can be thought of as coordinates on a logical grid, where each
209    grid point has (potentially) several degrees of freedom.
210    Any of y, z, n, and p can be passed in as NULL if not needed.
211 
212   Level: beginner
213 
214 .keywords: distributed array, get, corners, nodes, local indices
215 
216 .seealso: DMDAGetGhostCorners(), DMDAGetOwnershipRanges()
217 @*/
218 PetscErrorCode  DMDAGetCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p)
219 {
220   PetscInt w;
221   DM_DA    *dd = (DM_DA*)da->data;
222 
223   PetscFunctionBegin;
224   PetscValidHeaderSpecific(da,DM_CLASSID,1);
225   /* since the xs, xe ... have all been multiplied by the number of degrees
226      of freedom per cell, w = dd->w, we divide that out before returning.*/
227   w = dd->w;
228   if (x) *x = dd->xs/w + dd->xo;
229   /* the y and z have NOT been multiplied by w */
230   if (y) *y = dd->ys + dd->yo;
231   if (z) *z = dd->zs + dd->zo;
232   if (m) *m = (dd->xe - dd->xs)/w;
233   if (n) *n = (dd->ye - dd->ys);
234   if (p) *p = (dd->ze - dd->zs);
235   PetscFunctionReturn(0);
236 }
237 
238 /*@
239    DMDAGetLocalBoundingBox - Returns the local bounding box for the DMDA.
240 
241    Not Collective
242 
243    Input Parameter:
244 .  dm - the DM
245 
246    Output Parameters:
247 +  lmin - local minimum coordinates (length dim, optional)
248 -  lmax - local maximim coordinates (length dim, optional)
249 
250   Level: beginner
251 
252 .keywords: distributed array, get, coordinates
253 
254 .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetBoundingBox()
255 @*/
256 PetscErrorCode DMDAGetLocalBoundingBox(DM dm,PetscReal lmin[],PetscReal lmax[])
257 {
258   PetscErrorCode    ierr;
259   Vec               coords = NULL;
260   PetscInt          dim,i,j;
261   const PetscScalar *local_coords;
262   PetscReal         min[3]={PETSC_MAX_REAL,PETSC_MAX_REAL,PETSC_MAX_REAL},max[3]={PETSC_MIN_REAL,PETSC_MIN_REAL,PETSC_MIN_REAL};
263   PetscInt          N,Ni;
264 
265   PetscFunctionBegin;
266   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
267   dim  = dm->dim;
268   ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr);
269   if (coords) {
270     ierr = VecGetArrayRead(coords,&local_coords);CHKERRQ(ierr);
271     ierr = VecGetLocalSize(coords,&N);CHKERRQ(ierr);
272     Ni   = N/dim;
273     for (i=0; i<Ni; i++) {
274       for (j=0; j<3; j++) {
275         min[j] = j < dim ? PetscMin(min[j],PetscRealPart(local_coords[i*dim+j])) : 0;
276         max[j] = j < dim ? PetscMax(max[j],PetscRealPart(local_coords[i*dim+j])) : 0;
277       }
278     }
279     ierr = VecRestoreArrayRead(coords,&local_coords);CHKERRQ(ierr);
280   } else {                      /* Just use grid indices */
281     DMDALocalInfo info;
282     ierr   = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
283     min[0] = info.xs;
284     min[1] = info.ys;
285     min[2] = info.zs;
286     max[0] = info.xs + info.xm-1;
287     max[1] = info.ys + info.ym-1;
288     max[2] = info.zs + info.zm-1;
289   }
290   if (lmin) {ierr = PetscMemcpy(lmin,min,dim*sizeof(PetscReal));CHKERRQ(ierr);}
291   if (lmax) {ierr = PetscMemcpy(lmax,max,dim*sizeof(PetscReal));CHKERRQ(ierr);}
292   PetscFunctionReturn(0);
293 }
294 
295 /*@
296    DMDAGetBoundingBox - Returns the global bounding box for the DMDA.
297 
298    Collective on DMDA
299 
300    Input Parameter:
301 .  dm - the DM
302 
303    Output Parameters:
304 +  gmin - global minimum coordinates (length dim, optional)
305 -  gmax - global maximim coordinates (length dim, optional)
306 
307   Level: beginner
308 
309 .keywords: distributed array, get, coordinates
310 
311 .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetLocalBoundingBox()
312 @*/
313 PetscErrorCode DMDAGetBoundingBox(DM dm,PetscReal gmin[],PetscReal gmax[])
314 {
315   PetscErrorCode ierr;
316   PetscMPIInt    count;
317   PetscReal      lmin[3],lmax[3];
318 
319   PetscFunctionBegin;
320   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
321   ierr = PetscMPIIntCast(dm->dim,&count);CHKERRQ(ierr);
322   ierr = DMDAGetLocalBoundingBox(dm,lmin,lmax);CHKERRQ(ierr);
323   if (gmin) {ierr = MPIU_Allreduce(lmin,gmin,count,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);}
324   if (gmax) {ierr = MPIU_Allreduce(lmax,gmax,count,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);}
325   PetscFunctionReturn(0);
326 }
327 
328 /*@
329    DMDAGetReducedDMDA - Gets the DMDA with the same layout but with fewer or more fields
330 
331    Collective on DMDA
332 
333    Input Parameters:
334 +  da - the distributed array
335 -  nfields - number of fields in new DMDA
336 
337    Output Parameter:
338 .  nda - the new DMDA
339 
340   Level: intermediate
341 
342 .keywords: distributed array, get, corners, nodes, local indices, coordinates
343 
344 .seealso: DMDAGetGhostCorners(), DMSetCoordinates(), DMDASetUniformCoordinates(), DMGetCoordinates(), DMDAGetGhostedCoordinates()
345 @*/
346 PetscErrorCode  DMDAGetReducedDMDA(DM da,PetscInt nfields,DM *nda)
347 {
348   PetscErrorCode   ierr;
349   DM_DA            *dd = (DM_DA*)da->data;
350   PetscInt         s,m,n,p,M,N,P,dim,Mo,No,Po;
351   const PetscInt   *lx,*ly,*lz;
352   DMBoundaryType   bx,by,bz;
353   DMDAStencilType  stencil_type;
354   PetscInt         ox,oy,oz;
355   PetscInt         cl,rl;
356 
357   PetscFunctionBegin;
358   dim = da->dim;
359   M   = dd->M;
360   N   = dd->N;
361   P   = dd->P;
362   m   = dd->m;
363   n   = dd->n;
364   p   = dd->p;
365   s   = dd->s;
366   bx  = dd->bx;
367   by  = dd->by;
368   bz  = dd->bz;
369 
370   stencil_type = dd->stencil_type;
371 
372   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr);
373   if (dim == 1) {
374     ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),bx,M,nfields,s,dd->lx,nda);CHKERRQ(ierr);
375   } else if (dim == 2) {
376     ierr = DMDACreate2d(PetscObjectComm((PetscObject)da),bx,by,stencil_type,M,N,m,n,nfields,s,lx,ly,nda);CHKERRQ(ierr);
377   } else if (dim == 3) {
378     ierr = DMDACreate3d(PetscObjectComm((PetscObject)da),bx,by,bz,stencil_type,M,N,P,m,n,p,nfields,s,lx,ly,lz,nda);CHKERRQ(ierr);
379   }
380   ierr = DMSetUp(*nda);CHKERRQ(ierr);
381   if (da->coordinates) {
382     ierr = PetscObjectReference((PetscObject)da->coordinates);CHKERRQ(ierr);
383     (*nda)->coordinates = da->coordinates;
384   }
385 
386   /* allow for getting a reduced DA corresponding to a domain decomposition */
387   ierr = DMDAGetOffset(da,&ox,&oy,&oz,&Mo,&No,&Po);CHKERRQ(ierr);
388   ierr = DMDASetOffset(*nda,ox,oy,oz,Mo,No,Po);CHKERRQ(ierr);
389 
390   /* allow for getting a reduced DA corresponding to a coarsened DA */
391   ierr = DMGetCoarsenLevel(da,&cl);CHKERRQ(ierr);
392   ierr = DMGetRefineLevel(da,&rl);CHKERRQ(ierr);
393 
394   (*nda)->levelup   = rl;
395   (*nda)->leveldown = cl;
396   PetscFunctionReturn(0);
397 }
398 
399 /*@C
400    DMDAGetCoordinateArray - Gets an array containing the coordinates of the DMDA
401 
402    Not Collective
403 
404    Input Parameter:
405 .  dm - the DM
406 
407    Output Parameter:
408 .  xc - the coordinates
409 
410   Level: intermediate
411 
412 .keywords: distributed array, get, component name
413 
414 .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName(), DMDARestoreCoordinateArray()
415 @*/
416 PetscErrorCode DMDAGetCoordinateArray(DM dm,void *xc)
417 {
418   PetscErrorCode ierr;
419   DM             cdm;
420   Vec            x;
421 
422   PetscFunctionBegin;
423   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
424   ierr = DMGetCoordinates(dm,&x);CHKERRQ(ierr);
425   ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr);
426   ierr = DMDAVecGetArray(cdm,x,xc);CHKERRQ(ierr);
427   PetscFunctionReturn(0);
428 }
429 
430 /*@C
431    DMDARestoreCoordinateArray - Sets an array containing the coordinates of the DMDA
432 
433    Not Collective
434 
435    Input Parameter:
436 +  dm - the DM
437 -  xc - the coordinates
438 
439   Level: intermediate
440 
441 .keywords: distributed array, get, component name
442 
443 .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName(), DMDAGetCoordinateArray()
444 @*/
445 PetscErrorCode DMDARestoreCoordinateArray(DM dm,void *xc)
446 {
447   PetscErrorCode ierr;
448   DM             cdm;
449   Vec            x;
450 
451   PetscFunctionBegin;
452   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
453   ierr = DMGetCoordinates(dm,&x);CHKERRQ(ierr);
454   ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr);
455   ierr = DMDAVecRestoreArray(cdm,x,xc);CHKERRQ(ierr);
456   PetscFunctionReturn(0);
457 }
458