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