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