xref: /petsc/src/dm/impls/da/dacorn.c (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
347c6ae99SBarry Smith   Code for manipulating distributed regular arrays in parallel.
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"   I*/
7f19dbd58SToby Isaac #include <petscdmfield.h>
847c6ae99SBarry Smith 
9d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateCoordinateDM_DA(DM dm, DM *cdm)
10d71ae5a4SJacob Faibussowitsch {
11dd4c3f67SMatthew G. Knepley   const char *prefix;
12dd4c3f67SMatthew G. Knepley 
133f2d3b52SPeter Brune   PetscFunctionBegin;
149566063dSJacob Faibussowitsch   PetscCall(DMDACreateCompatibleDMDA(dm, dm->dim, cdm));
15dd4c3f67SMatthew G. Knepley   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
16dd4c3f67SMatthew G. Knepley   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*cdm, prefix));
17dd4c3f67SMatthew G. Knepley   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)*cdm, "cdm_"));
183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1947c6ae99SBarry Smith }
2047c6ae99SBarry Smith 
21d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateCoordinateField_DA(DM dm, DMField *field)
22d71ae5a4SJacob Faibussowitsch {
23f19dbd58SToby Isaac   PetscReal   gmin[3], gmax[3];
244d1a973fSToby Isaac   PetscScalar corners[24];
25f19dbd58SToby Isaac   PetscInt    dim;
26f19dbd58SToby Isaac   PetscInt    i, j;
27f19dbd58SToby Isaac   DM          cdm;
28f19dbd58SToby Isaac 
29f19dbd58SToby Isaac   PetscFunctionBegin;
309566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
31f19dbd58SToby Isaac   /* TODO: this is wrong if coordinates are not rectilinear */
329566063dSJacob Faibussowitsch   PetscCall(DMGetBoundingBox(dm, gmin, gmax));
33f19dbd58SToby Isaac   for (i = 0; i < (1 << dim); i++) {
34ad540459SPierre Jolivet     for (j = 0; j < dim; j++) corners[i * dim + j] = (i & (1 << j)) ? gmax[j] : gmin[j];
35f19dbd58SToby Isaac   }
369566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &cdm));
379566063dSJacob Faibussowitsch   PetscCall(DMFieldCreateDA(cdm, dim, corners, field));
389566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&cdm));
393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40f19dbd58SToby Isaac }
41f19dbd58SToby Isaac 
4247c6ae99SBarry Smith /*@C
43aa219208SBarry Smith   DMDASetFieldName - Sets the names of individual field components in multicomponent
44dce8aebaSBarry Smith   vectors associated with a `DMDA`.
4547c6ae99SBarry Smith 
4620f4b53cSBarry Smith   Logically Collective; name must contain a common value
4747c6ae99SBarry Smith 
4847c6ae99SBarry Smith   Input Parameters:
4947c6ae99SBarry Smith + da   - the distributed array
50dce8aebaSBarry Smith . nf   - field number for the `DMDA` (0, 1, ... dof-1), where dof indicates the
51dce8aebaSBarry Smith         number of degrees of freedom per node within the `DMDA`
5260225df5SJacob Faibussowitsch - name - the name of the field (component)
5347c6ae99SBarry Smith 
5447c6ae99SBarry Smith   Level: intermediate
5547c6ae99SBarry Smith 
56dce8aebaSBarry Smith   Note:
57dce8aebaSBarry Smith   It must be called after having called `DMSetUp()`.
58dce8aebaSBarry Smith 
59dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetFieldName()`, `DMDASetCoordinateName()`, `DMDAGetCoordinateName()`, `DMDASetFieldNames()`, `DMSetUp()`
6047c6ae99SBarry Smith @*/
61d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetFieldName(DM da, PetscInt nf, const char name[])
62d71ae5a4SJacob Faibussowitsch {
6347c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
6447c6ae99SBarry Smith 
6547c6ae99SBarry Smith   PetscFunctionBegin;
66a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
671dca8a05SBarry Smith   PetscCheck(nf >= 0 && nf < dd->w, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid field number: %" PetscInt_FMT, nf);
687a8be351SBarry Smith   PetscCheck(dd->fieldname, PetscObjectComm((PetscObject)da), PETSC_ERR_ORDER, "You should call DMSetUp() first");
699566063dSJacob Faibussowitsch   PetscCall(PetscFree(dd->fieldname[nf]));
709566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &dd->fieldname[nf]));
713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7247c6ae99SBarry Smith }
7347c6ae99SBarry Smith 
74c629b14aSBarry Smith /*@C
75dce8aebaSBarry Smith   DMDAGetFieldNames - Gets the name of each component in the vector associated with the `DMDA`
76c629b14aSBarry Smith 
7720f4b53cSBarry Smith   Not Collective; names will contain a common value; No Fortran Support
78c629b14aSBarry Smith 
79c629b14aSBarry Smith   Input Parameter:
8060225df5SJacob Faibussowitsch . da - the `DMDA` object
81c629b14aSBarry Smith 
82c629b14aSBarry Smith   Output Parameter:
8320f4b53cSBarry Smith . 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`
84c629b14aSBarry Smith 
85c629b14aSBarry Smith   Level: intermediate
86c629b14aSBarry Smith 
8760225df5SJacob Faibussowitsch   Fortran Notes:
8820f4b53cSBarry Smith   Use `DMDAGetFieldName()`
89f5f57ec0SBarry Smith 
90dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetFieldName()`, `DMDASetCoordinateName()`, `DMDAGetCoordinateName()`, `DMDASetFieldName()`, `DMDASetFieldNames()`
91c629b14aSBarry Smith @*/
92d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetFieldNames(DM da, const char *const **names)
93d71ae5a4SJacob Faibussowitsch {
94c629b14aSBarry Smith   DM_DA *dd = (DM_DA *)da->data;
95c629b14aSBarry Smith 
96c629b14aSBarry Smith   PetscFunctionBegin;
97c629b14aSBarry Smith   *names = (const char *const *)dd->fieldname;
983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99c629b14aSBarry Smith }
100c629b14aSBarry Smith 
101c629b14aSBarry Smith /*@C
102c629b14aSBarry Smith   DMDASetFieldNames - Sets the name of each component in the vector associated with the DMDA
103c629b14aSBarry Smith 
10420f4b53cSBarry Smith   Logically Collective; names must contain a common value; No Fortran Support
105c629b14aSBarry Smith 
106c629b14aSBarry Smith   Input Parameters:
10760225df5SJacob Faibussowitsch + da    - the `DMDA` object
108dce8aebaSBarry Smith - 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`
1093eac1ceeSStefano Zampini 
110c629b14aSBarry Smith   Level: intermediate
111c629b14aSBarry Smith 
112dce8aebaSBarry Smith   Note:
113dce8aebaSBarry Smith   It must be called after having called `DMSetUp()`.
114f5f57ec0SBarry Smith 
11560225df5SJacob Faibussowitsch   Fortran Notes:
11620f4b53cSBarry Smith   Use `DMDASetFieldName()`
117dce8aebaSBarry Smith 
118dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetFieldName()`, `DMDASetCoordinateName()`, `DMDAGetCoordinateName()`, `DMDASetFieldName()`, `DMSetUp()`
119c629b14aSBarry Smith @*/
120d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetFieldNames(DM da, const char *const *names)
121d71ae5a4SJacob Faibussowitsch {
122c629b14aSBarry Smith   DM_DA   *dd = (DM_DA *)da->data;
1233eac1ceeSStefano Zampini   char   **fieldname;
1243eac1ceeSStefano Zampini   PetscInt nf = 0;
125c629b14aSBarry Smith 
126c629b14aSBarry Smith   PetscFunctionBegin;
1277a8be351SBarry Smith   PetscCheck(dd->fieldname, PetscObjectComm((PetscObject)da), PETSC_ERR_ORDER, "You should call DMSetUp() first");
128a8f51744SPierre Jolivet   while (names[nf++]) { }
12963a3b9bcSJacob Faibussowitsch   PetscCheck(nf == dd->w + 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of fields %" PetscInt_FMT, nf - 1);
1309566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy(names, &fieldname));
1319566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&dd->fieldname));
1323eac1ceeSStefano Zampini   dd->fieldname = fieldname;
1333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
134c629b14aSBarry Smith }
135c629b14aSBarry Smith 
13647c6ae99SBarry Smith /*@C
137aa219208SBarry Smith   DMDAGetFieldName - Gets the names of individual field components in multicomponent
138dce8aebaSBarry Smith   vectors associated with a `DMDA`.
13947c6ae99SBarry Smith 
14020f4b53cSBarry Smith   Not Collective; name will contain a common value
14147c6ae99SBarry Smith 
142d8d19677SJose E. Roman   Input Parameters:
14347c6ae99SBarry Smith + da - the distributed array
144dce8aebaSBarry Smith - nf - field number for the `DMDA` (0, 1, ... dof-1), where dof indicates the
145dce8aebaSBarry Smith         number of degrees of freedom per node within the `DMDA`
14647c6ae99SBarry Smith 
14747c6ae99SBarry Smith   Output Parameter:
14860225df5SJacob Faibussowitsch . name - the name of the field (component)
14947c6ae99SBarry Smith 
15047c6ae99SBarry Smith   Level: intermediate
15147c6ae99SBarry Smith 
152dce8aebaSBarry Smith   Note:
153dce8aebaSBarry Smith   It must be called after having called `DMSetUp()`.
154dce8aebaSBarry Smith 
155dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDASetFieldName()`, `DMDASetCoordinateName()`, `DMDAGetCoordinateName()`, `DMSetUp()`
15647c6ae99SBarry Smith @*/
157d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetFieldName(DM da, PetscInt nf, const char **name)
158d71ae5a4SJacob Faibussowitsch {
15947c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
16047c6ae99SBarry Smith 
16147c6ae99SBarry Smith   PetscFunctionBegin;
162a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1634f572ea9SToby Isaac   PetscAssertPointer(name, 3);
1641dca8a05SBarry Smith   PetscCheck(nf >= 0 && nf < dd->w, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid field number: %" PetscInt_FMT, nf);
1657a8be351SBarry Smith   PetscCheck(dd->fieldname, PetscObjectComm((PetscObject)da), PETSC_ERR_ORDER, "You should call DMSetUp() first");
16647c6ae99SBarry Smith   *name = dd->fieldname[nf];
1673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16847c6ae99SBarry Smith }
16947c6ae99SBarry Smith 
170109c9344SBarry Smith /*@C
171dce8aebaSBarry Smith   DMDASetCoordinateName - Sets the name of the coordinate directions associated with a `DMDA`, for example "x" or "y"
172109c9344SBarry Smith 
17320f4b53cSBarry Smith   Logically Collective; name must contain a common value; No Fortran Support
174109c9344SBarry Smith 
175109c9344SBarry Smith   Input Parameters:
176dce8aebaSBarry Smith + dm   - the `DMDA`
177109c9344SBarry Smith . nf   - coordinate number for the DMDA (0, 1, ... dim-1),
178109c9344SBarry Smith - name - the name of the coordinate
179109c9344SBarry Smith 
180109c9344SBarry Smith   Level: intermediate
181109c9344SBarry Smith 
182dce8aebaSBarry Smith   Note:
18320f4b53cSBarry Smith   Must be called after having called `DMSetUp()`.
184f5f57ec0SBarry Smith 
185dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetCoordinateName()`, `DMDASetFieldName()`, `DMDAGetFieldName()`, `DMSetUp()`
186109c9344SBarry Smith @*/
187d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetCoordinateName(DM dm, PetscInt nf, const char name[])
188d71ae5a4SJacob Faibussowitsch {
189c73cfb54SMatthew G. Knepley   DM_DA *dd = (DM_DA *)dm->data;
190109c9344SBarry Smith 
191109c9344SBarry Smith   PetscFunctionBegin;
192a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
1931dca8a05SBarry Smith   PetscCheck(nf >= 0 && nf < dm->dim, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid coordinate number: %" PetscInt_FMT, nf);
1947a8be351SBarry Smith   PetscCheck(dd->coordinatename, PetscObjectComm((PetscObject)dm), PETSC_ERR_ORDER, "You should call DMSetUp() first");
1959566063dSJacob Faibussowitsch   PetscCall(PetscFree(dd->coordinatename[nf]));
1969566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &dd->coordinatename[nf]));
1973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
198109c9344SBarry Smith }
199109c9344SBarry Smith 
200109c9344SBarry Smith /*@C
201da81f932SPierre Jolivet   DMDAGetCoordinateName - Gets the name of a coordinate direction associated with a `DMDA`.
202109c9344SBarry Smith 
20320f4b53cSBarry Smith   Not Collective; name will contain a common value; No Fortran Support
204109c9344SBarry Smith 
205d8d19677SJose E. Roman   Input Parameters:
206dce8aebaSBarry Smith + dm - the `DMDA`
20720f4b53cSBarry Smith - nf - number for the `DMDA` (0, 1, ... dim-1)
208109c9344SBarry Smith 
209109c9344SBarry Smith   Output Parameter:
21060225df5SJacob Faibussowitsch . name - the name of the coordinate direction
211109c9344SBarry Smith 
212109c9344SBarry Smith   Level: intermediate
213109c9344SBarry Smith 
214dce8aebaSBarry Smith   Note:
215dce8aebaSBarry Smith   It must be called after having called `DMSetUp()`.
216dce8aebaSBarry Smith 
217dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDASetCoordinateName()`, `DMDASetFieldName()`, `DMDAGetFieldName()`, `DMSetUp()`
218109c9344SBarry Smith @*/
219d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetCoordinateName(DM dm, PetscInt nf, const char **name)
220d71ae5a4SJacob Faibussowitsch {
221c73cfb54SMatthew G. Knepley   DM_DA *dd = (DM_DA *)dm->data;
222109c9344SBarry Smith 
223109c9344SBarry Smith   PetscFunctionBegin;
224a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
2254f572ea9SToby Isaac   PetscAssertPointer(name, 3);
2261dca8a05SBarry Smith   PetscCheck(nf >= 0 && nf < dm->dim, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid coordinate number: %" PetscInt_FMT, nf);
2277a8be351SBarry Smith   PetscCheck(dd->coordinatename, PetscObjectComm((PetscObject)dm), PETSC_ERR_ORDER, "You should call DMSetUp() first");
228109c9344SBarry Smith   *name = dd->coordinatename[nf];
2293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
230109c9344SBarry Smith }
231109c9344SBarry Smith 
232a5d1443cSVincent Le Chenadec /*@C
233aa219208SBarry Smith   DMDAGetCorners - Returns the global (x,y,z) indices of the lower left
23459f3ab6dSMatthew G. Knepley   corner and size of the local region, excluding ghost points.
23547c6ae99SBarry Smith 
23620f4b53cSBarry Smith   Not Collective
23747c6ae99SBarry Smith 
23847c6ae99SBarry Smith   Input Parameter:
23947c6ae99SBarry Smith . da - the distributed array
24047c6ae99SBarry Smith 
24147c6ae99SBarry Smith   Output Parameters:
2426b867d5aSJose E. Roman + x - the corner index for the first dimension
2436b867d5aSJose E. Roman . y - the corner index for the second dimension (only used in 2D and 3D problems)
2446b867d5aSJose E. Roman . z - the corner index for the third dimension (only used in 3D problems)
2456b867d5aSJose E. Roman . m - the width in the first dimension
2466b867d5aSJose E. Roman . n - the width in the second dimension (only used in 2D and 3D problems)
2476b867d5aSJose E. Roman - p - the width in the third dimension (only used in 3D problems)
24847c6ae99SBarry Smith 
249dce8aebaSBarry Smith   Level: beginner
250dce8aebaSBarry Smith 
25147c6ae99SBarry Smith   Note:
25247c6ae99SBarry Smith   The corner information is independent of the number of degrees of
253dce8aebaSBarry Smith   freedom per node set with the `DMDACreateXX()` routine. Thus the x, y, z, and
25447c6ae99SBarry Smith   m, n, p can be thought of as coordinates on a logical grid, where each
25547c6ae99SBarry Smith   grid point has (potentially) several degrees of freedom.
2560298fd71SBarry Smith   Any of y, z, n, and p can be passed in as NULL if not needed.
25747c6ae99SBarry Smith 
258dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetOwnershipRanges()`, `DMStagGetCorners()`
25947c6ae99SBarry Smith @*/
260d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetCorners(DM da, PetscInt *x, PetscInt *y, PetscInt *z, PetscInt *m, PetscInt *n, PetscInt *p)
261d71ae5a4SJacob Faibussowitsch {
26247c6ae99SBarry Smith   PetscInt w;
26347c6ae99SBarry Smith   DM_DA   *dd = (DM_DA *)da->data;
26447c6ae99SBarry Smith 
26547c6ae99SBarry Smith   PetscFunctionBegin;
266a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
26747c6ae99SBarry Smith   /* since the xs, xe ... have all been multiplied by the number of degrees
26847c6ae99SBarry Smith      of freedom per cell, w = dd->w, we divide that out before returning.*/
26947c6ae99SBarry Smith   w = dd->w;
27059bc5b24SSatish Balay   if (x) *x = dd->xs / w + dd->xo;
27147c6ae99SBarry Smith   /* the y and z have NOT been multiplied by w */
27259bc5b24SSatish Balay   if (y) *y = dd->ys + dd->yo;
27359bc5b24SSatish Balay   if (z) *z = dd->zs + dd->zo;
27459bc5b24SSatish Balay   if (m) *m = (dd->xe - dd->xs) / w;
27559bc5b24SSatish Balay   if (n) *n = (dd->ye - dd->ys);
27659bc5b24SSatish Balay   if (p) *p = (dd->ze - dd->zs);
2773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27847c6ae99SBarry Smith }
27947c6ae99SBarry Smith 
280d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetLocalBoundingIndices_DMDA(DM dm, PetscReal lmin[], PetscReal lmax[])
281d71ae5a4SJacob Faibussowitsch {
2827324c66bSJed Brown   DMDALocalInfo info;
28347c6ae99SBarry Smith 
28447c6ae99SBarry Smith   PetscFunctionBegin;
2859566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
286b2e4378dSMatthew G. Knepley   lmin[0] = info.xs;
287b2e4378dSMatthew G. Knepley   lmin[1] = info.ys;
288b2e4378dSMatthew G. Knepley   lmin[2] = info.zs;
289b2e4378dSMatthew G. Knepley   lmax[0] = info.xs + info.xm - 1;
290b2e4378dSMatthew G. Knepley   lmax[1] = info.ys + info.ym - 1;
291b2e4378dSMatthew G. Knepley   lmax[2] = info.zs + info.zm - 1;
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29347c6ae99SBarry Smith }
294bc2bf880SBarry Smith 
295*a4e35b19SJacob Faibussowitsch // PetscClangLinter pragma ignore: -fdoc-*
296bc2bf880SBarry Smith /*@
29783d91586SPatrick Sanan   DMDAGetReducedDMDA - Deprecated; use DMDACreateCompatibleDMDA()
29883d91586SPatrick Sanan 
29983d91586SPatrick Sanan   Level: deprecated
30083d91586SPatrick Sanan @*/
301d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetReducedDMDA(DM da, PetscInt nfields, DM *nda)
302d71ae5a4SJacob Faibussowitsch {
30383d91586SPatrick Sanan   PetscFunctionBegin;
3049566063dSJacob Faibussowitsch   PetscCall(DMDACreateCompatibleDMDA(da, nfields, nda));
3053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30683d91586SPatrick Sanan }
30783d91586SPatrick Sanan 
30883d91586SPatrick Sanan /*@
309dce8aebaSBarry Smith   DMDACreateCompatibleDMDA - Creates a `DMDA` with the same layout but with fewer or more fields
310bc2bf880SBarry Smith 
311b1dd8793SDave May   Collective
312bc2bf880SBarry Smith 
313907376e6SBarry Smith   Input Parameters:
314bc2bf880SBarry Smith + da      - the distributed array
315dce8aebaSBarry Smith - nfields - number of fields in new `DMDA`
316bc2bf880SBarry Smith 
317bc2bf880SBarry Smith   Output Parameter:
318dce8aebaSBarry Smith . nda - the new `DMDA`
319bc2bf880SBarry Smith 
320bc2bf880SBarry Smith   Level: intermediate
321bc2bf880SBarry Smith 
322dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMSetCoordinates()`, `DMDASetUniformCoordinates()`, `DMGetCoordinates()`, `DMDAGetGhostedCoordinates()`, `DMStagCreateCompatibleDMStag()`
323bc2bf880SBarry Smith @*/
324d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDACreateCompatibleDMDA(DM da, PetscInt nfields, DM *nda)
325d71ae5a4SJacob Faibussowitsch {
326bc2bf880SBarry Smith   DM_DA          *dd = (DM_DA *)da->data;
32795c13181SPeter Brune   PetscInt        s, m, n, p, M, N, P, dim, Mo, No, Po;
328320964c4SBlaise Bourdin   const PetscInt *lx, *ly, *lz;
329bff4a2f0SMatthew G. Knepley   DMBoundaryType  bx, by, bz;
330320964c4SBlaise Bourdin   DMDAStencilType stencil_type;
3316858538eSMatthew G. Knepley   Vec             coords;
33295c13181SPeter Brune   PetscInt        ox, oy, oz;
33395c13181SPeter Brune   PetscInt        cl, rl;
334320964c4SBlaise Bourdin 
335320964c4SBlaise Bourdin   PetscFunctionBegin;
336c73cfb54SMatthew G. Knepley   dim = da->dim;
33795c13181SPeter Brune   M   = dd->M;
33895c13181SPeter Brune   N   = dd->N;
33995c13181SPeter Brune   P   = dd->P;
34095c13181SPeter Brune   m   = dd->m;
34195c13181SPeter Brune   n   = dd->n;
34295c13181SPeter Brune   p   = dd->p;
34395c13181SPeter Brune   s   = dd->s;
34495c13181SPeter Brune   bx  = dd->bx;
34595c13181SPeter Brune   by  = dd->by;
34695c13181SPeter Brune   bz  = dd->bz;
3478865f1eaSKarl Rupp 
34895c13181SPeter Brune   stencil_type = dd->stencil_type;
3498865f1eaSKarl Rupp 
3509566063dSJacob Faibussowitsch   PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, &lz));
351320964c4SBlaise Bourdin   if (dim == 1) {
3529566063dSJacob Faibussowitsch     PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da), bx, M, nfields, s, dd->lx, nda));
353320964c4SBlaise Bourdin   } else if (dim == 2) {
3549566063dSJacob Faibussowitsch     PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)da), bx, by, stencil_type, M, N, m, n, nfields, s, lx, ly, nda));
355320964c4SBlaise Bourdin   } else if (dim == 3) {
3569566063dSJacob Faibussowitsch     PetscCall(DMDACreate3d(PetscObjectComm((PetscObject)da), bx, by, bz, stencil_type, M, N, P, m, n, p, nfields, s, lx, ly, lz, nda));
357bc2bf880SBarry Smith   }
3589566063dSJacob Faibussowitsch   PetscCall(DMSetUp(*nda));
3596858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(da, &coords));
3606858538eSMatthew G. Knepley   PetscCall(DMSetCoordinates(*nda, coords));
36195c13181SPeter Brune 
36295c13181SPeter Brune   /* allow for getting a reduced DA corresponding to a domain decomposition */
3639566063dSJacob Faibussowitsch   PetscCall(DMDAGetOffset(da, &ox, &oy, &oz, &Mo, &No, &Po));
3649566063dSJacob Faibussowitsch   PetscCall(DMDASetOffset(*nda, ox, oy, oz, Mo, No, Po));
36595c13181SPeter Brune 
36695c13181SPeter Brune   /* allow for getting a reduced DA corresponding to a coarsened DA */
3679566063dSJacob Faibussowitsch   PetscCall(DMGetCoarsenLevel(da, &cl));
3689566063dSJacob Faibussowitsch   PetscCall(DMGetRefineLevel(da, &rl));
3698865f1eaSKarl Rupp 
37095c13181SPeter Brune   (*nda)->levelup   = rl;
37195c13181SPeter Brune   (*nda)->leveldown = cl;
3723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
373bc2bf880SBarry Smith }
374bc2bf880SBarry Smith 
375c593f006SBarry Smith /*@C
376dce8aebaSBarry Smith   DMDAGetCoordinateArray - Gets an array containing the coordinates of the `DMDA`
377c593f006SBarry Smith 
37820f4b53cSBarry Smith   Not Collective; No Fortran Support
379c593f006SBarry Smith 
380c593f006SBarry Smith   Input Parameter:
381dce8aebaSBarry Smith . dm - the `DMDA`
382c593f006SBarry Smith 
383c593f006SBarry Smith   Output Parameter:
384c593f006SBarry Smith . xc - the coordinates
385c593f006SBarry Smith 
386c593f006SBarry Smith   Level: intermediate
387c593f006SBarry Smith 
388dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDASetCoordinateName()`, `DMDASetFieldName()`, `DMDAGetFieldName()`, `DMDARestoreCoordinateArray()`
389c593f006SBarry Smith @*/
390d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetCoordinateArray(DM dm, void *xc)
391d71ae5a4SJacob Faibussowitsch {
392c593f006SBarry Smith   DM  cdm;
393c593f006SBarry Smith   Vec x;
394c593f006SBarry Smith 
395c593f006SBarry Smith   PetscFunctionBegin;
396c593f006SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3979566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(dm, &x));
3989566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
3999566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cdm, x, xc));
4003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
401c593f006SBarry Smith }
402c593f006SBarry Smith 
403c593f006SBarry Smith /*@C
404dce8aebaSBarry Smith   DMDARestoreCoordinateArray - Sets an array containing the coordinates of the `DMDA`
405c593f006SBarry Smith 
40620f4b53cSBarry Smith   Not Collective; No Fortran Support
407c593f006SBarry Smith 
408d8d19677SJose E. Roman   Input Parameters:
409dce8aebaSBarry Smith + dm - the `DMDA`
410c593f006SBarry Smith - xc - the coordinates
411c593f006SBarry Smith 
412c593f006SBarry Smith   Level: intermediate
413c593f006SBarry Smith 
414dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDASetCoordinateName()`, `DMDASetFieldName()`, `DMDAGetFieldName()`, `DMDAGetCoordinateArray()`
415c593f006SBarry Smith @*/
416d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreCoordinateArray(DM dm, void *xc)
417d71ae5a4SJacob Faibussowitsch {
418c593f006SBarry Smith   DM  cdm;
419c593f006SBarry Smith   Vec x;
420c593f006SBarry Smith 
421c593f006SBarry Smith   PetscFunctionBegin;
422c593f006SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4239566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(dm, &x));
4249566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
4259566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cdm, x, xc));
4263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
427c593f006SBarry Smith }
428