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