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