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