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