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