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