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