xref: /petsc/src/dm/impls/da/da.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"   I*/
247c6ae99SBarry Smith 
347c6ae99SBarry Smith /*@
4e43dc8daSBarry Smith   DMDASetSizes - Sets the number of grid points in the three dimensional directions
547c6ae99SBarry Smith 
6d083f849SBarry Smith   Logically Collective on da
747c6ae99SBarry Smith 
847c6ae99SBarry Smith   Input Parameters:
9aa219208SBarry Smith + da - the DMDA
10e43dc8daSBarry Smith . M - the global X size
11e43dc8daSBarry Smith . N - the global Y size
12e43dc8daSBarry Smith - P - the global Z size
1347c6ae99SBarry Smith 
1447c6ae99SBarry Smith   Level: intermediate
1547c6ae99SBarry Smith 
16628e512eSPatrick Sanan   Developer Notes:
17628e512eSPatrick Sanan   Since the dimension may not yet have been set the code cannot error check for non-positive Y and Z number of grid points
18e43dc8daSBarry Smith 
19db781477SPatrick Sanan .seealso: `PetscSplitOwnership()`
2047c6ae99SBarry Smith @*/
219371c9d4SSatish Balay PetscErrorCode DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P) {
2247c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   PetscFunctionBegin;
25a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
2647c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, M, 2);
2747c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, N, 3);
2847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, P, 4);
297a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
307a8be351SBarry Smith   PetscCheck(M >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Number of grid points in X direction must be positive");
317a8be351SBarry Smith   PetscCheck(N >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Number of grid points in Y direction must be positive");
327a8be351SBarry Smith   PetscCheck(P >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Number of grid points in Z direction must be positive");
3347c6ae99SBarry Smith 
3447c6ae99SBarry Smith   dd->M = M;
3547c6ae99SBarry Smith   dd->N = N;
3647c6ae99SBarry Smith   dd->P = P;
3747c6ae99SBarry Smith   PetscFunctionReturn(0);
3847c6ae99SBarry Smith }
3947c6ae99SBarry Smith 
4047c6ae99SBarry Smith /*@
41aa219208SBarry Smith   DMDASetNumProcs - Sets the number of processes in each dimension
4247c6ae99SBarry Smith 
43d083f849SBarry Smith   Logically Collective on da
4447c6ae99SBarry Smith 
4547c6ae99SBarry Smith   Input Parameters:
46aa219208SBarry Smith + da - the DMDA
4747c6ae99SBarry Smith . m - the number of X procs (or PETSC_DECIDE)
4847c6ae99SBarry Smith . n - the number of Y procs (or PETSC_DECIDE)
4947c6ae99SBarry Smith - p - the number of Z procs (or PETSC_DECIDE)
5047c6ae99SBarry Smith 
5147c6ae99SBarry Smith   Level: intermediate
5247c6ae99SBarry Smith 
53db781477SPatrick Sanan .seealso: `DMDASetSizes()`, `PetscSplitOwnership()`
5447c6ae99SBarry Smith @*/
559371c9d4SSatish Balay PetscErrorCode DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p) {
5647c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
5747c6ae99SBarry Smith 
5847c6ae99SBarry Smith   PetscFunctionBegin;
59a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
6047c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, m, 2);
6147c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, n, 3);
6247c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, p, 4);
637a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
6447c6ae99SBarry Smith   dd->m = m;
6547c6ae99SBarry Smith   dd->n = n;
6647c6ae99SBarry Smith   dd->p = p;
67c73cfb54SMatthew G. Knepley   if (da->dim == 2) {
68d3be247eSBarry Smith     PetscMPIInt size;
699566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da), &size));
70e3f3e4b6SBarry Smith     if ((dd->m > 0) && (dd->n < 0)) {
71e3f3e4b6SBarry Smith       dd->n = size / dd->m;
7263a3b9bcSJacob Faibussowitsch       PetscCheck(dd->n * dd->m == size, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " processes in X direction not divisible into comm size %d", m, size);
73e3f3e4b6SBarry Smith     }
74e3f3e4b6SBarry Smith     if ((dd->n > 0) && (dd->m < 0)) {
75e3f3e4b6SBarry Smith       dd->m = size / dd->n;
7663a3b9bcSJacob Faibussowitsch       PetscCheck(dd->n * dd->m == size, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " processes in Y direction not divisible into comm size %d", n, size);
77e3f3e4b6SBarry Smith     }
78e3f3e4b6SBarry Smith   }
7947c6ae99SBarry Smith   PetscFunctionReturn(0);
8047c6ae99SBarry Smith }
8147c6ae99SBarry Smith 
8247c6ae99SBarry Smith /*@
831321219cSEthan Coon   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries.
8447c6ae99SBarry Smith 
8547c6ae99SBarry Smith   Not collective
8647c6ae99SBarry Smith 
87d8d19677SJose E. Roman   Input Parameters:
88aa219208SBarry Smith + da    - The DMDA
89bff4a2f0SMatthew G. Knepley - bx,by,bz - One of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC
9047c6ae99SBarry Smith 
9147c6ae99SBarry Smith   Level: intermediate
9247c6ae99SBarry Smith 
93db781477SPatrick Sanan .seealso: `DMDACreate()`, `DMDestroy()`, `DMDA`, `DMBoundaryType`
9447c6ae99SBarry Smith @*/
959371c9d4SSatish Balay PetscErrorCode DMDASetBoundaryType(DM da, DMBoundaryType bx, DMBoundaryType by, DMBoundaryType bz) {
9647c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith   PetscFunctionBegin;
99a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1005a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da, bx, 2);
1015a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da, by, 3);
1025a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da, bz, 4);
1037a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
1041321219cSEthan Coon   dd->bx = bx;
1051321219cSEthan Coon   dd->by = by;
1061321219cSEthan Coon   dd->bz = bz;
10747c6ae99SBarry Smith   PetscFunctionReturn(0);
10847c6ae99SBarry Smith }
10947c6ae99SBarry Smith 
11047c6ae99SBarry Smith /*@
111aa219208SBarry Smith   DMDASetDof - Sets the number of degrees of freedom per vertex
11247c6ae99SBarry Smith 
11347c6ae99SBarry Smith   Not collective
11447c6ae99SBarry Smith 
11559f3ab6dSMatthew G. Knepley   Input Parameters:
116aa219208SBarry Smith + da  - The DMDA
11747c6ae99SBarry Smith - dof - Number of degrees of freedom
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith   Level: intermediate
12047c6ae99SBarry Smith 
121db781477SPatrick Sanan .seealso: `DMDAGetDof()`, `DMDACreate()`, `DMDestroy()`, `DMDA`
12247c6ae99SBarry Smith @*/
1239371c9d4SSatish Balay PetscErrorCode DMDASetDof(DM da, PetscInt dof) {
12447c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith   PetscFunctionBegin;
127a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
12854cfb0beSLisandro Dalcin   PetscValidLogicalCollectiveInt(da, dof, 2);
1297a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
13047c6ae99SBarry Smith   dd->w  = dof;
1311411c6eeSJed Brown   da->bs = dof;
13247c6ae99SBarry Smith   PetscFunctionReturn(0);
13347c6ae99SBarry Smith }
13447c6ae99SBarry Smith 
135fb6725baSMatthew G. Knepley /*@
136fb6725baSMatthew G. Knepley   DMDAGetDof - Gets the number of degrees of freedom per vertex
137fb6725baSMatthew G. Knepley 
138fb6725baSMatthew G. Knepley   Not collective
139fb6725baSMatthew G. Knepley 
140fb6725baSMatthew G. Knepley   Input Parameter:
141fb6725baSMatthew G. Knepley . da  - The DMDA
142fb6725baSMatthew G. Knepley 
143fb6725baSMatthew G. Knepley   Output Parameter:
144fb6725baSMatthew G. Knepley . dof - Number of degrees of freedom
145fb6725baSMatthew G. Knepley 
146fb6725baSMatthew G. Knepley   Level: intermediate
147fb6725baSMatthew G. Knepley 
148db781477SPatrick Sanan .seealso: `DMDASetDof()`, `DMDACreate()`, `DMDestroy()`, `DMDA`
149fb6725baSMatthew G. Knepley @*/
1509371c9d4SSatish Balay PetscErrorCode DMDAGetDof(DM da, PetscInt *dof) {
151fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *)da->data;
152fb6725baSMatthew G. Knepley 
153fb6725baSMatthew G. Knepley   PetscFunctionBegin;
154a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
155dadcf809SJacob Faibussowitsch   PetscValidIntPointer(dof, 2);
156fb6725baSMatthew G. Knepley   *dof = dd->w;
157fb6725baSMatthew G. Knepley   PetscFunctionReturn(0);
158fb6725baSMatthew G. Knepley }
159fb6725baSMatthew G. Knepley 
1607ddda789SPeter Brune /*@
1617ddda789SPeter Brune   DMDAGetOverlap - Gets the size of the per-processor overlap.
1627ddda789SPeter Brune 
1637ddda789SPeter Brune   Not collective
1647ddda789SPeter Brune 
165f899ff85SJose E. Roman   Input Parameter:
1667ddda789SPeter Brune . da  - The DMDA
1677ddda789SPeter Brune 
1687ddda789SPeter Brune   Output Parameters:
1697ddda789SPeter Brune + x   - Overlap in the x direction
1707ddda789SPeter Brune . y   - Overlap in the y direction
1717ddda789SPeter Brune - z   - Overlap in the z direction
1727ddda789SPeter Brune 
1737ddda789SPeter Brune   Level: intermediate
1747ddda789SPeter Brune 
1752b3cbbdaSStefano Zampini .seealso: `DMCreateDomainDecomposition()`, `DMDASetOverlap()`, `DMDA`
1767ddda789SPeter Brune @*/
1779371c9d4SSatish Balay PetscErrorCode DMDAGetOverlap(DM da, PetscInt *x, PetscInt *y, PetscInt *z) {
1787ddda789SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
1797ddda789SPeter Brune 
1807ddda789SPeter Brune   PetscFunctionBegin;
181a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1827ddda789SPeter Brune   if (x) *x = dd->xol;
1837ddda789SPeter Brune   if (y) *y = dd->yol;
1847ddda789SPeter Brune   if (z) *z = dd->zol;
1857ddda789SPeter Brune   PetscFunctionReturn(0);
1867ddda789SPeter Brune }
1877ddda789SPeter Brune 
18888661749SPeter Brune /*@
18988661749SPeter Brune   DMDASetOverlap - Sets the size of the per-processor overlap.
19088661749SPeter Brune 
19188661749SPeter Brune   Not collective
19288661749SPeter Brune 
1937ddda789SPeter Brune   Input Parameters:
19488661749SPeter Brune + da  - The DMDA
1957ddda789SPeter Brune . x   - Overlap in the x direction
1967ddda789SPeter Brune . y   - Overlap in the y direction
1977ddda789SPeter Brune - z   - Overlap in the z direction
19888661749SPeter Brune 
19988661749SPeter Brune   Level: intermediate
20088661749SPeter Brune 
2012b3cbbdaSStefano Zampini .seealso: `DMCreateDomainDecomposition()`, `DMDAGetOverlap()`, `DMDA`
20288661749SPeter Brune @*/
2039371c9d4SSatish Balay PetscErrorCode DMDASetOverlap(DM da, PetscInt x, PetscInt y, PetscInt z) {
20488661749SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
20588661749SPeter Brune 
20688661749SPeter Brune   PetscFunctionBegin;
207a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
2087ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da, x, 2);
2097ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da, y, 3);
2107ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da, z, 4);
2117ddda789SPeter Brune   dd->xol = x;
2127ddda789SPeter Brune   dd->yol = y;
2137ddda789SPeter Brune   dd->zol = z;
21488661749SPeter Brune   PetscFunctionReturn(0);
21588661749SPeter Brune }
21688661749SPeter Brune 
2173e7870d2SPeter Brune /*@
2183e7870d2SPeter Brune   DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition.
2193e7870d2SPeter Brune 
2203e7870d2SPeter Brune   Not collective
2213e7870d2SPeter Brune 
2223e7870d2SPeter Brune   Input Parameters:
2233e7870d2SPeter Brune . da  - The DMDA
2243e7870d2SPeter Brune 
2253e7870d2SPeter Brune   Output Parameters:
226a2b725a8SWilliam Gropp . Nsub   - Number of local subdomains created upon decomposition
2273e7870d2SPeter Brune 
2283e7870d2SPeter Brune   Level: intermediate
2293e7870d2SPeter Brune 
2302b3cbbdaSStefano Zampini .seealso: `DMCreateDomainDecomposition()`, `DMDASetNumLocalSubDomains()`, `DMDA`
2313e7870d2SPeter Brune @*/
2329371c9d4SSatish Balay PetscErrorCode DMDAGetNumLocalSubDomains(DM da, PetscInt *Nsub) {
2333e7870d2SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
2343e7870d2SPeter Brune 
2353e7870d2SPeter Brune   PetscFunctionBegin;
236a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
2373e7870d2SPeter Brune   if (Nsub) *Nsub = dd->Nsub;
2383e7870d2SPeter Brune   PetscFunctionReturn(0);
2393e7870d2SPeter Brune }
2403e7870d2SPeter Brune 
2413e7870d2SPeter Brune /*@
2423e7870d2SPeter Brune   DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition.
2433e7870d2SPeter Brune 
2443e7870d2SPeter Brune   Not collective
2453e7870d2SPeter Brune 
2463e7870d2SPeter Brune   Input Parameters:
2473e7870d2SPeter Brune + da  - The DMDA
2483e7870d2SPeter Brune - Nsub - The number of local subdomains requested
2493e7870d2SPeter Brune 
2503e7870d2SPeter Brune   Level: intermediate
2513e7870d2SPeter Brune 
2522b3cbbdaSStefano Zampini .seealso: `DMCreateDomainDecomposition()`, `DMDAGetNumLocalSubDomains()`, `DMDA`
2533e7870d2SPeter Brune @*/
2549371c9d4SSatish Balay PetscErrorCode DMDASetNumLocalSubDomains(DM da, PetscInt Nsub) {
2553e7870d2SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
2563e7870d2SPeter Brune 
2573e7870d2SPeter Brune   PetscFunctionBegin;
258a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
2593e7870d2SPeter Brune   PetscValidLogicalCollectiveInt(da, Nsub, 2);
2603e7870d2SPeter Brune   dd->Nsub = Nsub;
2613e7870d2SPeter Brune   PetscFunctionReturn(0);
2623e7870d2SPeter Brune }
2633e7870d2SPeter Brune 
264d886c4f4SPeter Brune /*@
265d886c4f4SPeter Brune   DMDASetOffset - Sets the index offset of the DA.
266d886c4f4SPeter Brune 
267d886c4f4SPeter Brune   Collective on DA
268d886c4f4SPeter Brune 
269d8d19677SJose E. Roman   Input Parameters:
270d886c4f4SPeter Brune + da  - The DMDA
271d886c4f4SPeter Brune . xo  - The offset in the x direction
272d886c4f4SPeter Brune . yo  - The offset in the y direction
273d886c4f4SPeter Brune - zo  - The offset in the z direction
274d886c4f4SPeter Brune 
275d886c4f4SPeter Brune   Level: intermediate
276d886c4f4SPeter Brune 
27795452b02SPatrick Sanan   Notes:
27895452b02SPatrick Sanan     This is used primarily to overlap a computation on a local DA with that on a global DA without
279d886c4f4SPeter Brune   changing boundary conditions or subdomain features that depend upon the global offsets.
280d886c4f4SPeter Brune 
281db781477SPatrick Sanan .seealso: `DMDAGetOffset()`, `DMDAVecGetArray()`
282d886c4f4SPeter Brune @*/
2839371c9d4SSatish Balay PetscErrorCode DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po) {
284d886c4f4SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
285d886c4f4SPeter Brune 
286d886c4f4SPeter Brune   PetscFunctionBegin;
287a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
288d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da, xo, 2);
28995c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, yo, 3);
29095c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, zo, 4);
29195c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, Mo, 5);
29295c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, No, 6);
29395c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, Po, 7);
294d886c4f4SPeter Brune   dd->xo = xo;
295d886c4f4SPeter Brune   dd->yo = yo;
296d886c4f4SPeter Brune   dd->zo = zo;
29795c13181SPeter Brune   dd->Mo = Mo;
29895c13181SPeter Brune   dd->No = No;
29995c13181SPeter Brune   dd->Po = Po;
30095c13181SPeter Brune 
3016858538eSMatthew G. Knepley   if (da->coordinates[0].dm) PetscCall(DMDASetOffset(da->coordinates[0].dm, xo, yo, zo, Mo, No, Po));
302d886c4f4SPeter Brune   PetscFunctionReturn(0);
303d886c4f4SPeter Brune }
304d886c4f4SPeter Brune 
305d886c4f4SPeter Brune /*@
306d886c4f4SPeter Brune   DMDAGetOffset - Gets the index offset of the DA.
307d886c4f4SPeter Brune 
308d886c4f4SPeter Brune   Not collective
309d886c4f4SPeter Brune 
310d886c4f4SPeter Brune   Input Parameter:
311d886c4f4SPeter Brune . da  - The DMDA
312d886c4f4SPeter Brune 
313d886c4f4SPeter Brune   Output Parameters:
314d886c4f4SPeter Brune + xo  - The offset in the x direction
315d886c4f4SPeter Brune . yo  - The offset in the y direction
31695c13181SPeter Brune . zo  - The offset in the z direction
31795c13181SPeter Brune . Mo  - The global size in the x direction
31895c13181SPeter Brune . No  - The global size in the y direction
31995c13181SPeter Brune - Po  - The global size in the z direction
320d886c4f4SPeter Brune 
321d886c4f4SPeter Brune   Level: intermediate
322d886c4f4SPeter Brune 
323db781477SPatrick Sanan .seealso: `DMDASetOffset()`, `DMDAVecGetArray()`
324d886c4f4SPeter Brune @*/
3259371c9d4SSatish Balay PetscErrorCode DMDAGetOffset(DM da, PetscInt *xo, PetscInt *yo, PetscInt *zo, PetscInt *Mo, PetscInt *No, PetscInt *Po) {
326d886c4f4SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
327d886c4f4SPeter Brune 
328d886c4f4SPeter Brune   PetscFunctionBegin;
329a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
330d886c4f4SPeter Brune   if (xo) *xo = dd->xo;
331d886c4f4SPeter Brune   if (yo) *yo = dd->yo;
332d886c4f4SPeter Brune   if (zo) *zo = dd->zo;
33395c13181SPeter Brune   if (Mo) *Mo = dd->Mo;
33495c13181SPeter Brune   if (No) *No = dd->No;
33595c13181SPeter Brune   if (Po) *Po = dd->Po;
336d886c4f4SPeter Brune   PetscFunctionReturn(0);
337d886c4f4SPeter Brune }
338d886c4f4SPeter Brune 
33940234c92SPeter Brune /*@
34040234c92SPeter Brune   DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain DM.
34140234c92SPeter Brune 
34240234c92SPeter Brune   Not collective
34340234c92SPeter Brune 
34440234c92SPeter Brune   Input Parameter:
34540234c92SPeter Brune . da  - The DMDA
34640234c92SPeter Brune 
34740234c92SPeter Brune   Output Parameters:
34840234c92SPeter Brune + xs  - The start of the region in x
34940234c92SPeter Brune . ys  - The start of the region in y
35040234c92SPeter Brune . zs  - The start of the region in z
35140234c92SPeter Brune . xs  - The size of the region in x
35240234c92SPeter Brune . ys  - The size of the region in y
353a2b725a8SWilliam Gropp - zs  - The size of the region in z
35440234c92SPeter Brune 
35540234c92SPeter Brune   Level: intermediate
35640234c92SPeter Brune 
357db781477SPatrick Sanan .seealso: `DMDAGetOffset()`, `DMDAVecGetArray()`
35840234c92SPeter Brune @*/
3599371c9d4SSatish Balay PetscErrorCode DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm) {
36040234c92SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
36140234c92SPeter Brune 
36240234c92SPeter Brune   PetscFunctionBegin;
363a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
36440234c92SPeter Brune   if (xs) *xs = dd->nonxs;
36540234c92SPeter Brune   if (ys) *ys = dd->nonys;
36640234c92SPeter Brune   if (zs) *zs = dd->nonzs;
36740234c92SPeter Brune   if (xm) *xm = dd->nonxm;
36840234c92SPeter Brune   if (ym) *ym = dd->nonym;
36940234c92SPeter Brune   if (zm) *zm = dd->nonzm;
37040234c92SPeter Brune   PetscFunctionReturn(0);
37140234c92SPeter Brune }
37240234c92SPeter Brune 
37340234c92SPeter Brune /*@
37440234c92SPeter Brune   DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain DM.
37540234c92SPeter Brune 
37640234c92SPeter Brune   Collective on DA
37740234c92SPeter Brune 
378d8d19677SJose E. Roman   Input Parameters:
37940234c92SPeter Brune + da  - The DMDA
38040234c92SPeter Brune . xs  - The start of the region in x
38140234c92SPeter Brune . ys  - The start of the region in y
38240234c92SPeter Brune . zs  - The start of the region in z
38340234c92SPeter Brune . xs  - The size of the region in x
38440234c92SPeter Brune . ys  - The size of the region in y
385a2b725a8SWilliam Gropp - zs  - The size of the region in z
38640234c92SPeter Brune 
38740234c92SPeter Brune   Level: intermediate
38840234c92SPeter Brune 
389db781477SPatrick Sanan .seealso: `DMDAGetOffset()`, `DMDAVecGetArray()`
39040234c92SPeter Brune @*/
3919371c9d4SSatish Balay PetscErrorCode DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm) {
39240234c92SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
39340234c92SPeter Brune 
39440234c92SPeter Brune   PetscFunctionBegin;
395a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
39640234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, xs, 2);
39740234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, ys, 3);
39840234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, zs, 4);
39940234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, xm, 5);
40040234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, ym, 6);
40140234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, zm, 7);
40240234c92SPeter Brune   dd->nonxs = xs;
40340234c92SPeter Brune   dd->nonys = ys;
40440234c92SPeter Brune   dd->nonzs = zs;
40540234c92SPeter Brune   dd->nonxm = xm;
40640234c92SPeter Brune   dd->nonym = ym;
40740234c92SPeter Brune   dd->nonzm = zm;
40840234c92SPeter Brune 
40940234c92SPeter Brune   PetscFunctionReturn(0);
41040234c92SPeter Brune }
41188661749SPeter Brune 
41247c6ae99SBarry Smith /*@
413aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
41447c6ae99SBarry Smith 
415d083f849SBarry Smith   Logically Collective on da
41647c6ae99SBarry Smith 
417d8d19677SJose E. Roman   Input Parameters:
418aa219208SBarry Smith + da    - The DMDA
419aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
42047c6ae99SBarry Smith 
42147c6ae99SBarry Smith   Level: intermediate
42247c6ae99SBarry Smith 
423db781477SPatrick Sanan .seealso: `DMDACreate()`, `DMDestroy()`, `DMDA`
42447c6ae99SBarry Smith @*/
4259371c9d4SSatish Balay PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) {
42647c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
42747c6ae99SBarry Smith 
42847c6ae99SBarry Smith   PetscFunctionBegin;
429a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
43047c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da, stype, 2);
4317a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
43247c6ae99SBarry Smith   dd->stencil_type = stype;
43347c6ae99SBarry Smith   PetscFunctionReturn(0);
43447c6ae99SBarry Smith }
43547c6ae99SBarry Smith 
436fb6725baSMatthew G. Knepley /*@
437fb6725baSMatthew G. Knepley   DMDAGetStencilType - Gets the type of the communication stencil
438fb6725baSMatthew G. Knepley 
439fb6725baSMatthew G. Knepley   Not collective
440fb6725baSMatthew G. Knepley 
441fb6725baSMatthew G. Knepley   Input Parameter:
442fb6725baSMatthew G. Knepley . da    - The DMDA
443fb6725baSMatthew G. Knepley 
444fb6725baSMatthew G. Knepley   Output Parameter:
445fb6725baSMatthew G. Knepley . stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
446fb6725baSMatthew G. Knepley 
447fb6725baSMatthew G. Knepley   Level: intermediate
448fb6725baSMatthew G. Knepley 
449db781477SPatrick Sanan .seealso: `DMDACreate()`, `DMDestroy()`, `DMDA`
450fb6725baSMatthew G. Knepley @*/
4519371c9d4SSatish Balay PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype) {
452fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *)da->data;
453fb6725baSMatthew G. Knepley 
454fb6725baSMatthew G. Knepley   PetscFunctionBegin;
455a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
456fb6725baSMatthew G. Knepley   PetscValidPointer(stype, 2);
457fb6725baSMatthew G. Knepley   *stype = dd->stencil_type;
458fb6725baSMatthew G. Knepley   PetscFunctionReturn(0);
459fb6725baSMatthew G. Knepley }
460fb6725baSMatthew G. Knepley 
46147c6ae99SBarry Smith /*@
462aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
46347c6ae99SBarry Smith 
464d083f849SBarry Smith   Logically Collective on da
46547c6ae99SBarry Smith 
466d8d19677SJose E. Roman   Input Parameters:
467aa219208SBarry Smith + da    - The DMDA
46847c6ae99SBarry Smith - width - The stencil width
46947c6ae99SBarry Smith 
47047c6ae99SBarry Smith   Level: intermediate
47147c6ae99SBarry Smith 
472db781477SPatrick Sanan .seealso: `DMDACreate()`, `DMDestroy()`, `DMDA`
47347c6ae99SBarry Smith @*/
4749371c9d4SSatish Balay PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) {
47547c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
47647c6ae99SBarry Smith 
47747c6ae99SBarry Smith   PetscFunctionBegin;
478a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
47947c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, width, 2);
4807a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
48147c6ae99SBarry Smith   dd->s = width;
48247c6ae99SBarry Smith   PetscFunctionReturn(0);
48347c6ae99SBarry Smith }
48447c6ae99SBarry Smith 
485fb6725baSMatthew G. Knepley /*@
486fb6725baSMatthew G. Knepley   DMDAGetStencilWidth - Gets the width of the communication stencil
487fb6725baSMatthew G. Knepley 
488fb6725baSMatthew G. Knepley   Not collective
489fb6725baSMatthew G. Knepley 
490fb6725baSMatthew G. Knepley   Input Parameter:
491fb6725baSMatthew G. Knepley . da    - The DMDA
492fb6725baSMatthew G. Knepley 
493fb6725baSMatthew G. Knepley   Output Parameter:
494fb6725baSMatthew G. Knepley . width - The stencil width
495fb6725baSMatthew G. Knepley 
496fb6725baSMatthew G. Knepley   Level: intermediate
497fb6725baSMatthew G. Knepley 
498db781477SPatrick Sanan .seealso: `DMDACreate()`, `DMDestroy()`, `DMDA`
499fb6725baSMatthew G. Knepley @*/
5009371c9d4SSatish Balay PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width) {
501fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *)da->data;
502fb6725baSMatthew G. Knepley 
503fb6725baSMatthew G. Knepley   PetscFunctionBegin;
504a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
505dadcf809SJacob Faibussowitsch   PetscValidIntPointer(width, 2);
506fb6725baSMatthew G. Knepley   *width = dd->s;
507fb6725baSMatthew G. Knepley   PetscFunctionReturn(0);
508fb6725baSMatthew G. Knepley }
509fb6725baSMatthew G. Knepley 
5109371c9d4SSatish Balay static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da, PetscInt M, PetscInt m, const PetscInt lx[]) {
51147c6ae99SBarry Smith   PetscInt i, sum;
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith   PetscFunctionBegin;
5147a8be351SBarry Smith   PetscCheck(M >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Global dimension not set");
51547c6ae99SBarry Smith   for (i = sum = 0; i < m; i++) sum += lx[i];
51663a3b9bcSJacob Faibussowitsch   PetscCheck(sum == M, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "Ownership ranges sum to %" PetscInt_FMT " but global dimension is %" PetscInt_FMT, sum, M);
51747c6ae99SBarry Smith   PetscFunctionReturn(0);
51847c6ae99SBarry Smith }
51947c6ae99SBarry Smith 
52047c6ae99SBarry Smith /*@
521aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
52247c6ae99SBarry Smith 
523d083f849SBarry Smith   Logically Collective on da
52447c6ae99SBarry Smith 
525d8d19677SJose E. Roman   Input Parameters:
526aa219208SBarry Smith + da - The DMDA
5270298fd71SBarry Smith . lx - array containing number of nodes in the X direction on each process, or NULL. If non-null, must be of length da->m
5280298fd71SBarry Smith . ly - array containing number of nodes in the Y direction on each process, or NULL. If non-null, must be of length da->n
5290298fd71SBarry Smith - lz - array containing number of nodes in the Z direction on each process, or NULL. If non-null, must be of length da->p.
53047c6ae99SBarry Smith 
53147c6ae99SBarry Smith   Level: intermediate
53247c6ae99SBarry Smith 
533e3f3e4b6SBarry Smith   Note: these numbers are NOT multiplied by the number of dof per node.
534e3f3e4b6SBarry Smith 
535db781477SPatrick Sanan .seealso: `DMDACreate()`, `DMDestroy()`, `DMDA`
53647c6ae99SBarry Smith @*/
5379371c9d4SSatish Balay PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) {
53847c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
53947c6ae99SBarry Smith 
54047c6ae99SBarry Smith   PetscFunctionBegin;
541a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
5427a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
54347c6ae99SBarry Smith   if (lx) {
5447a8be351SBarry Smith     PetscCheck(dd->m >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of procs");
5459566063dSJacob Faibussowitsch     PetscCall(DMDACheckOwnershipRanges_Private(da, dd->M, dd->m, lx));
546*48a46eb9SPierre Jolivet     if (!dd->lx) PetscCall(PetscMalloc1(dd->m, &dd->lx));
5479566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd->lx, lx, dd->m));
54847c6ae99SBarry Smith   }
54947c6ae99SBarry Smith   if (ly) {
5507a8be351SBarry Smith     PetscCheck(dd->n >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of procs");
5519566063dSJacob Faibussowitsch     PetscCall(DMDACheckOwnershipRanges_Private(da, dd->N, dd->n, ly));
552*48a46eb9SPierre Jolivet     if (!dd->ly) PetscCall(PetscMalloc1(dd->n, &dd->ly));
5539566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd->ly, ly, dd->n));
55447c6ae99SBarry Smith   }
55547c6ae99SBarry Smith   if (lz) {
5567a8be351SBarry Smith     PetscCheck(dd->p >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of procs");
5579566063dSJacob Faibussowitsch     PetscCall(DMDACheckOwnershipRanges_Private(da, dd->P, dd->p, lz));
558*48a46eb9SPierre Jolivet     if (!dd->lz) PetscCall(PetscMalloc1(dd->p, &dd->lz));
5599566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd->lz, lz, dd->p));
56047c6ae99SBarry Smith   }
56147c6ae99SBarry Smith   PetscFunctionReturn(0);
56247c6ae99SBarry Smith }
56347c6ae99SBarry Smith 
56447c6ae99SBarry Smith /*@
565aa219208SBarry Smith        DMDASetInterpolationType - Sets the type of interpolation that will be
566e727c939SJed Brown           returned by DMCreateInterpolation()
56747c6ae99SBarry Smith 
568d083f849SBarry Smith    Logically Collective on da
56947c6ae99SBarry Smith 
570d8d19677SJose E. Roman    Input Parameters:
57147c6ae99SBarry Smith +  da - initial distributed array
572a2b725a8SWilliam Gropp -  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
57347c6ae99SBarry Smith 
57447c6ae99SBarry Smith    Level: intermediate
57547c6ae99SBarry Smith 
57695452b02SPatrick Sanan    Notes:
57795452b02SPatrick Sanan     you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation()
57847c6ae99SBarry Smith 
579db781477SPatrick Sanan .seealso: `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDestroy()`, `DMDA`, `DMDAInterpolationType`
58047c6ae99SBarry Smith @*/
5819371c9d4SSatish Balay PetscErrorCode DMDASetInterpolationType(DM da, DMDAInterpolationType ctype) {
58247c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
58347c6ae99SBarry Smith 
58447c6ae99SBarry Smith   PetscFunctionBegin;
585a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
58647c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da, ctype, 2);
58747c6ae99SBarry Smith   dd->interptype = ctype;
58847c6ae99SBarry Smith   PetscFunctionReturn(0);
58947c6ae99SBarry Smith }
59047c6ae99SBarry Smith 
5912dde6fd4SLisandro Dalcin /*@
5922dde6fd4SLisandro Dalcin        DMDAGetInterpolationType - Gets the type of interpolation that will be
593e727c939SJed Brown           used by DMCreateInterpolation()
5942dde6fd4SLisandro Dalcin 
5952dde6fd4SLisandro Dalcin    Not Collective
5962dde6fd4SLisandro Dalcin 
5972dde6fd4SLisandro Dalcin    Input Parameter:
5982dde6fd4SLisandro Dalcin .  da - distributed array
5992dde6fd4SLisandro Dalcin 
6002dde6fd4SLisandro Dalcin    Output Parameter:
6012dde6fd4SLisandro Dalcin .  ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms)
6022dde6fd4SLisandro Dalcin 
6032dde6fd4SLisandro Dalcin    Level: intermediate
6042dde6fd4SLisandro Dalcin 
605db781477SPatrick Sanan .seealso: `DMDA`, `DMDAInterpolationType`, `DMDASetInterpolationType()`, `DMCreateInterpolation()`
6062dde6fd4SLisandro Dalcin @*/
6079371c9d4SSatish Balay PetscErrorCode DMDAGetInterpolationType(DM da, DMDAInterpolationType *ctype) {
6082dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA *)da->data;
6092dde6fd4SLisandro Dalcin 
6102dde6fd4SLisandro Dalcin   PetscFunctionBegin;
611a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
6122dde6fd4SLisandro Dalcin   PetscValidPointer(ctype, 2);
6132dde6fd4SLisandro Dalcin   *ctype = dd->interptype;
6142dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
6152dde6fd4SLisandro Dalcin }
61647c6ae99SBarry Smith 
6176a119db4SBarry Smith /*@C
618aa219208SBarry Smith       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
61947c6ae99SBarry Smith         processes neighbors.
62047c6ae99SBarry Smith 
62147c6ae99SBarry Smith     Not Collective
62247c6ae99SBarry Smith 
62347c6ae99SBarry Smith    Input Parameter:
624aa219208SBarry Smith .     da - the DMDA object
62547c6ae99SBarry Smith 
62647c6ae99SBarry Smith    Output Parameters:
62747c6ae99SBarry Smith .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
62847c6ae99SBarry Smith               this process itself is in the list
62947c6ae99SBarry Smith 
63095452b02SPatrick Sanan    Notes:
63195452b02SPatrick Sanan     In 2d the array is of length 9, in 3d of length 27
63247c6ae99SBarry Smith           Not supported in 1d
633aa219208SBarry Smith           Do not free the array, it is freed when the DMDA is destroyed.
63447c6ae99SBarry Smith 
63595452b02SPatrick Sanan    Fortran Notes:
63695452b02SPatrick Sanan     In fortran you must pass in an array of the appropriate length.
63747c6ae99SBarry Smith 
63847c6ae99SBarry Smith    Level: intermediate
63947c6ae99SBarry Smith 
64047c6ae99SBarry Smith @*/
6419371c9d4SSatish Balay PetscErrorCode DMDAGetNeighbors(DM da, const PetscMPIInt *ranks[]) {
64247c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
6435fd66863SKarl Rupp 
64447c6ae99SBarry Smith   PetscFunctionBegin;
645a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
64647c6ae99SBarry Smith   *ranks = dd->neighbors;
64747c6ae99SBarry Smith   PetscFunctionReturn(0);
64847c6ae99SBarry Smith }
64947c6ae99SBarry Smith 
65047c6ae99SBarry Smith /*@C
651aa219208SBarry Smith       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
65247c6ae99SBarry Smith 
65347c6ae99SBarry Smith     Not Collective
65447c6ae99SBarry Smith 
65547c6ae99SBarry Smith    Input Parameter:
656aa219208SBarry Smith .     da - the DMDA object
65747c6ae99SBarry Smith 
658d8d19677SJose E. Roman    Output Parameters:
65947c6ae99SBarry Smith +     lx - ownership along x direction (optional)
66047c6ae99SBarry Smith .     ly - ownership along y direction (optional)
66147c6ae99SBarry Smith -     lz - ownership along z direction (optional)
66247c6ae99SBarry Smith 
66347c6ae99SBarry Smith    Level: intermediate
66447c6ae99SBarry Smith 
665aa219208SBarry Smith     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
668aa219208SBarry Smith     eighth arguments from DMDAGetInfo()
66947c6ae99SBarry Smith 
67047c6ae99SBarry Smith      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
671aa219208SBarry Smith     DMDA they came from still exists (has not been destroyed).
67247c6ae99SBarry Smith 
673e3f3e4b6SBarry Smith     These numbers are NOT multiplied by the number of dof per node.
674e3f3e4b6SBarry Smith 
675db781477SPatrick Sanan .seealso: `DMDAGetCorners()`, `DMDAGetGhostCorners()`, `DMDACreate()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `VecGetOwnershipRanges()`
67647c6ae99SBarry Smith @*/
6779371c9d4SSatish Balay PetscErrorCode DMDAGetOwnershipRanges(DM da, const PetscInt *lx[], const PetscInt *ly[], const PetscInt *lz[]) {
67847c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
67947c6ae99SBarry Smith 
68047c6ae99SBarry Smith   PetscFunctionBegin;
681a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
68247c6ae99SBarry Smith   if (lx) *lx = dd->lx;
68347c6ae99SBarry Smith   if (ly) *ly = dd->ly;
68447c6ae99SBarry Smith   if (lz) *lz = dd->lz;
68547c6ae99SBarry Smith   PetscFunctionReturn(0);
68647c6ae99SBarry Smith }
68747c6ae99SBarry Smith 
68847c6ae99SBarry Smith /*@
689aa219208SBarry Smith      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
69047c6ae99SBarry Smith 
691d083f849SBarry Smith     Logically Collective on da
69247c6ae99SBarry Smith 
69347c6ae99SBarry Smith   Input Parameters:
694aa219208SBarry Smith +    da - the DMDA object
69547c6ae99SBarry Smith .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
69647c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
69747c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
69847c6ae99SBarry Smith 
69947c6ae99SBarry Smith   Options Database:
70048eeb7c8SBarry Smith +  -da_refine_x refine_x - refinement ratio in x direction
70148eeb7c8SBarry Smith .  -da_refine_y rafine_y - refinement ratio in y direction
70248eeb7c8SBarry Smith .  -da_refine_z refine_z - refinement ratio in z direction
70348eeb7c8SBarry Smith -  -da_refine <n> - refine the DMDA object n times when it is created.
70447c6ae99SBarry Smith 
70547c6ae99SBarry Smith   Level: intermediate
70647c6ae99SBarry Smith 
70795452b02SPatrick Sanan     Notes:
70895452b02SPatrick Sanan     Pass PETSC_IGNORE to leave a value unchanged
70947c6ae99SBarry Smith 
710db781477SPatrick Sanan .seealso: `DMRefine()`, `DMDAGetRefinementFactor()`
71147c6ae99SBarry Smith @*/
7129371c9d4SSatish Balay PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y, PetscInt refine_z) {
71347c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
71447c6ae99SBarry Smith 
71547c6ae99SBarry Smith   PetscFunctionBegin;
716a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
71747c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, refine_x, 2);
71847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, refine_y, 3);
71947c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, refine_z, 4);
72047c6ae99SBarry Smith 
72147c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
72247c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
72347c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
72447c6ae99SBarry Smith   PetscFunctionReturn(0);
72547c6ae99SBarry Smith }
72647c6ae99SBarry Smith 
72747c6ae99SBarry Smith /*@C
728aa219208SBarry Smith      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
72947c6ae99SBarry Smith 
73047c6ae99SBarry Smith     Not Collective
73147c6ae99SBarry Smith 
73247c6ae99SBarry Smith   Input Parameter:
733aa219208SBarry Smith .    da - the DMDA object
73447c6ae99SBarry Smith 
73547c6ae99SBarry Smith   Output Parameters:
73647c6ae99SBarry Smith +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
73747c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
73847c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
73947c6ae99SBarry Smith 
74047c6ae99SBarry Smith   Level: intermediate
74147c6ae99SBarry Smith 
74295452b02SPatrick Sanan     Notes:
74395452b02SPatrick Sanan     Pass NULL for values you do not need
74447c6ae99SBarry Smith 
745db781477SPatrick Sanan .seealso: `DMRefine()`, `DMDASetRefinementFactor()`
74647c6ae99SBarry Smith @*/
7479371c9d4SSatish Balay PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y, PetscInt *refine_z) {
74847c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
74947c6ae99SBarry Smith 
75047c6ae99SBarry Smith   PetscFunctionBegin;
751a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
75247c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
75347c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
75447c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
75547c6ae99SBarry Smith   PetscFunctionReturn(0);
75647c6ae99SBarry Smith }
75747c6ae99SBarry Smith 
75847c6ae99SBarry Smith /*@C
759aa219208SBarry Smith      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
76047c6ae99SBarry Smith 
761d083f849SBarry Smith     Logically Collective on da
76247c6ae99SBarry Smith 
76347c6ae99SBarry Smith   Input Parameters:
764aa219208SBarry Smith +    da - the DMDA object
765aa219208SBarry Smith -    f - the function that allocates the matrix for that specific DMDA
76647c6ae99SBarry Smith 
76747c6ae99SBarry Smith   Level: developer
76847c6ae99SBarry Smith 
76995452b02SPatrick Sanan    Notes:
77095452b02SPatrick Sanan     See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
77147c6ae99SBarry Smith        the diagonal and off-diagonal blocks of the matrix
77247c6ae99SBarry Smith 
7731fd49c25SBarry Smith    Not supported from Fortran
7741fd49c25SBarry Smith 
775db781477SPatrick Sanan .seealso: `DMCreateMatrix()`, `DMDASetBlockFills()`
77647c6ae99SBarry Smith @*/
7779371c9d4SSatish Balay PetscErrorCode DMDASetGetMatrix(DM da, PetscErrorCode (*f)(DM, Mat *)) {
77847c6ae99SBarry Smith   PetscFunctionBegin;
779a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
78025296bd5SBarry Smith   da->ops->creatematrix = f;
78147c6ae99SBarry Smith   PetscFunctionReturn(0);
78247c6ae99SBarry Smith }
78347c6ae99SBarry Smith 
78438fb4e8eSJunchao Zhang /*@
78538fb4e8eSJunchao Zhang    DMDAMapMatStencilToGlobal - Map a list of MatStencils on a grid to global indices.
78638fb4e8eSJunchao Zhang 
78738fb4e8eSJunchao Zhang    Not Collective
78838fb4e8eSJunchao Zhang 
78938fb4e8eSJunchao Zhang    Input Parameters:
79038fb4e8eSJunchao Zhang +  da - the DMDA object
79138fb4e8eSJunchao Zhang .  m - number of MatStencils
79238fb4e8eSJunchao Zhang -  idxm - grid points (and component number when dof > 1)
79338fb4e8eSJunchao Zhang 
7941179163eSBarry Smith    Output Parameter:
7951179163eSBarry Smith .  gidxm - global row indices
7961179163eSBarry Smith 
7971179163eSBarry Smith    Level: intermediate
79838fb4e8eSJunchao Zhang 
79938fb4e8eSJunchao Zhang .seealso: `MatStencil`
80038fb4e8eSJunchao Zhang @*/
8019371c9d4SSatish Balay PetscErrorCode DMDAMapMatStencilToGlobal(DM da, PetscInt m, const MatStencil idxm[], PetscInt gidxm[]) {
80238fb4e8eSJunchao Zhang   const DM_DA           *dd  = (const DM_DA *)da->data;
80338fb4e8eSJunchao Zhang   const PetscInt        *dxm = (const PetscInt *)idxm;
80438fb4e8eSJunchao Zhang   PetscInt               i, j, sdim, tmp, dim;
80538fb4e8eSJunchao Zhang   PetscInt               dims[4], starts[4], dims2[3], starts2[3], dof = dd->w;
80638fb4e8eSJunchao Zhang   ISLocalToGlobalMapping ltog;
80738fb4e8eSJunchao Zhang 
80838fb4e8eSJunchao Zhang   PetscFunctionBegin;
80938fb4e8eSJunchao Zhang   if (m <= 0) PetscFunctionReturn(0);
81038fb4e8eSJunchao Zhang 
8112f27f4d1SJunchao Zhang   /* Code adapted from DMDAGetGhostCorners() */
81238fb4e8eSJunchao Zhang   starts2[0] = dd->Xs / dof + dd->xo;
81338fb4e8eSJunchao Zhang   starts2[1] = dd->Ys + dd->yo;
81438fb4e8eSJunchao Zhang   starts2[2] = dd->Zs + dd->zo;
81538fb4e8eSJunchao Zhang   dims2[0]   = (dd->Xe - dd->Xs) / dof;
81638fb4e8eSJunchao Zhang   dims2[1]   = (dd->Ye - dd->Ys);
81738fb4e8eSJunchao Zhang   dims2[2]   = (dd->Ze - dd->Zs);
81838fb4e8eSJunchao Zhang 
8192f27f4d1SJunchao Zhang   /* As if we do MatSetStencil() to get dims[]/starts[] of mat->stencil */
8202f27f4d1SJunchao Zhang   dim  = da->dim;                   /* DA dim: 1 to 3 */
8212f27f4d1SJunchao Zhang   sdim = dim + (dof > 1 ? 1 : 0);   /* Dimensions in MatStencil's (k,j,i,c) view */
8222f27f4d1SJunchao Zhang   for (i = 0; i < dim; i++) {       /* Reverse the order and also skip the unused dimensions */
8232f27f4d1SJunchao Zhang     dims[i]   = dims2[dim - i - 1]; /* ex. dims/starts[] are in order of {i} for 1D, {j,i} for 2D and {k,j,i} for 3D */
82438fb4e8eSJunchao Zhang     starts[i] = starts2[dim - i - 1];
82538fb4e8eSJunchao Zhang   }
8262f27f4d1SJunchao Zhang   starts[dim] = 0; /* Append the extra dim for dof (won't be used below if dof=1) */
82738fb4e8eSJunchao Zhang   dims[dim]   = dof;
82838fb4e8eSJunchao Zhang 
82938fb4e8eSJunchao Zhang   /* Map stencils to local indices (code adapted from MatSetValuesStencil()) */
83038fb4e8eSJunchao Zhang   for (i = 0; i < m; i++) {
8312f27f4d1SJunchao Zhang     dxm += 3 - dim; /* Input is {k,j,i,c}; move the pointer to the first used index, e.g., j in 2D */
8322f27f4d1SJunchao Zhang     tmp = 0;
8332f27f4d1SJunchao Zhang     for (j = 0; j < sdim; j++) {                                                      /* Iter over, ex. j,i or j,i,c in 2D */
8342f27f4d1SJunchao Zhang       if (tmp < 0 || dxm[j] < starts[j] || dxm[j] >= (starts[j] + dims[j])) tmp = -1; /* Beyond the ghost region, therefore ignored with negative indices */
8352f27f4d1SJunchao Zhang       else tmp = tmp * dims[j] + (dxm[j] - starts[j]);
83638fb4e8eSJunchao Zhang     }
83738fb4e8eSJunchao Zhang     gidxm[i] = tmp;
8382f27f4d1SJunchao Zhang     /* Move to the next MatStencil point */
8392f27f4d1SJunchao Zhang     if (dof > 1) dxm += sdim; /* c is already counted in sdim */
8402f27f4d1SJunchao Zhang     else dxm += sdim + 1;     /* skip the unused c */
84138fb4e8eSJunchao Zhang   }
84238fb4e8eSJunchao Zhang 
84338fb4e8eSJunchao Zhang   /* Map local indices to global indices */
84438fb4e8eSJunchao Zhang   PetscCall(DMGetLocalToGlobalMapping(da, &ltog));
84538fb4e8eSJunchao Zhang   PetscCall(ISLocalToGlobalMappingApply(ltog, m, gidxm, gidxm));
84638fb4e8eSJunchao Zhang   PetscFunctionReturn(0);
84738fb4e8eSJunchao Zhang }
84838fb4e8eSJunchao Zhang 
84947c6ae99SBarry Smith /*
85047c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
85147c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
85247c6ae99SBarry Smith 
85347c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
85447c6ae99SBarry Smith */
8559371c9d4SSatish Balay static PetscErrorCode DMDARefineOwnershipRanges(DM da, PetscBool periodic, PetscInt stencil_width, PetscInt ratio, PetscInt m, const PetscInt lc[], PetscInt lf[]) {
85647c6ae99SBarry Smith   PetscInt i, totalc = 0, remaining, startc = 0, startf = 0;
85747c6ae99SBarry Smith 
85847c6ae99SBarry Smith   PetscFunctionBegin;
85963a3b9bcSJacob Faibussowitsch   PetscCheck(ratio >= 1, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Requested refinement ratio %" PetscInt_FMT " must be at least 1", ratio);
86047c6ae99SBarry Smith   if (ratio == 1) {
8619566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(lf, lc, m));
86247c6ae99SBarry Smith     PetscFunctionReturn(0);
86347c6ae99SBarry Smith   }
86447c6ae99SBarry Smith   for (i = 0; i < m; i++) totalc += lc[i];
86547c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
86647c6ae99SBarry Smith   for (i = 0; i < m; i++) {
86747c6ae99SBarry Smith     PetscInt want = remaining / (m - i) + !!(remaining % (m - i));
86847c6ae99SBarry Smith     if (i == m - 1) lf[i] = want;
86947c6ae99SBarry Smith     else {
8707aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
8717aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
8727aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
8737aca7175SJed Brown       while ((startf + want) / ratio < nextc - stencil_width) want++;
8747aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
8757aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
8767aca7175SJed Brown       while ((startf + want - 1 + ratio - 1) / ratio > nextc - 1 + stencil_width) want--;
8777aca7175SJed Brown       /* Make sure all constraints are satisfied */
8789371c9d4SSatish Balay       if (want < 0 || want > remaining || ((startf + want) / ratio < nextc - stencil_width) || ((startf + want - 1 + ratio - 1) / ratio > nextc - 1 + stencil_width))
8799371c9d4SSatish Balay         SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Could not find a compatible refined ownership range");
88047c6ae99SBarry Smith     }
88147c6ae99SBarry Smith     lf[i] = want;
88247c6ae99SBarry Smith     startc += lc[i];
88347c6ae99SBarry Smith     startf += lf[i];
88447c6ae99SBarry Smith     remaining -= lf[i];
88547c6ae99SBarry Smith   }
88647c6ae99SBarry Smith   PetscFunctionReturn(0);
88747c6ae99SBarry Smith }
88847c6ae99SBarry Smith 
8892be375d4SJed Brown /*
8902be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
8912be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
8922be375d4SJed Brown 
8932be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
8942be375d4SJed Brown */
8959371c9d4SSatish Balay static PetscErrorCode DMDACoarsenOwnershipRanges(DM da, PetscBool periodic, PetscInt stencil_width, PetscInt ratio, PetscInt m, const PetscInt lf[], PetscInt lc[]) {
8962be375d4SJed Brown   PetscInt i, totalf, remaining, startc, startf;
8972be375d4SJed Brown 
8982be375d4SJed Brown   PetscFunctionBegin;
89963a3b9bcSJacob Faibussowitsch   PetscCheck(ratio >= 1, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Requested refinement ratio %" PetscInt_FMT " must be at least 1", ratio);
9002be375d4SJed Brown   if (ratio == 1) {
9019566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(lc, lf, m));
9022be375d4SJed Brown     PetscFunctionReturn(0);
9032be375d4SJed Brown   }
9042be375d4SJed Brown   for (i = 0, totalf = 0; i < m; i++) totalf += lf[i];
9052be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
9062be375d4SJed Brown   for (i = 0, startc = 0, startf = 0; i < m; i++) {
9072be375d4SJed Brown     PetscInt want = remaining / (m - i) + !!(remaining % (m - i));
9082be375d4SJed Brown     if (i == m - 1) lc[i] = want;
9092be375d4SJed Brown     else {
9102be375d4SJed Brown       const PetscInt nextf = startf + lf[i];
9112be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
9122be375d4SJed Brown        * node is within one stencil width. */
9132be375d4SJed Brown       while (nextf / ratio < startc + want - stencil_width) want--;
9142be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
9152be375d4SJed Brown        * fine node is within one stencil width. */
9162be375d4SJed Brown       while ((nextf - 1 + ratio - 1) / ratio > startc + want - 1 + stencil_width) want++;
9179371c9d4SSatish Balay       if (want < 0 || want > remaining || (nextf / ratio < startc + want - stencil_width) || ((nextf - 1 + ratio - 1) / ratio > startc + want - 1 + stencil_width))
9189371c9d4SSatish Balay         SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Could not find a compatible coarsened ownership range");
9192be375d4SJed Brown     }
9202be375d4SJed Brown     lc[i] = want;
9212be375d4SJed Brown     startc += lc[i];
9222be375d4SJed Brown     startf += lf[i];
9232be375d4SJed Brown     remaining -= lc[i];
9242be375d4SJed Brown   }
9252be375d4SJed Brown   PetscFunctionReturn(0);
9262be375d4SJed Brown }
9272be375d4SJed Brown 
9289371c9d4SSatish Balay PetscErrorCode DMRefine_DA(DM da, MPI_Comm comm, DM *daref) {
929c73cfb54SMatthew G. Knepley   PetscInt M, N, P, i, dim;
9306858538eSMatthew G. Knepley   Vec      coordsc, coordsf;
9319a42bb27SBarry Smith   DM       da2;
93247c6ae99SBarry Smith   DM_DA   *dd = (DM_DA *)da->data, *dd2;
93347c6ae99SBarry Smith 
93447c6ae99SBarry Smith   PetscFunctionBegin;
935a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
93647c6ae99SBarry Smith   PetscValidPointer(daref, 3);
93747c6ae99SBarry Smith 
9389566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(da, &dim));
939bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
94047c6ae99SBarry Smith     M = dd->refine_x * dd->M;
94147c6ae99SBarry Smith   } else {
94247c6ae99SBarry Smith     M = 1 + dd->refine_x * (dd->M - 1);
94347c6ae99SBarry Smith   }
944bff4a2f0SMatthew G. Knepley   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
945c73cfb54SMatthew G. Knepley     if (dim > 1) {
94647c6ae99SBarry Smith       N = dd->refine_y * dd->N;
94747c6ae99SBarry Smith     } else {
9481860e6e9SBarry Smith       N = 1;
9491860e6e9SBarry Smith     }
9501860e6e9SBarry Smith   } else {
95147c6ae99SBarry Smith     N = 1 + dd->refine_y * (dd->N - 1);
95247c6ae99SBarry Smith   }
953bff4a2f0SMatthew G. Knepley   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
954c73cfb54SMatthew G. Knepley     if (dim > 2) {
95547c6ae99SBarry Smith       P = dd->refine_z * dd->P;
95647c6ae99SBarry Smith     } else {
9571860e6e9SBarry Smith       P = 1;
9581860e6e9SBarry Smith     }
9591860e6e9SBarry Smith   } else {
96047c6ae99SBarry Smith     P = 1 + dd->refine_z * (dd->P - 1);
96147c6ae99SBarry Smith   }
9629566063dSJacob Faibussowitsch   PetscCall(DMDACreate(PetscObjectComm((PetscObject)da), &da2));
9639566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(da2, ((PetscObject)da)->prefix));
9649566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(da2, dim));
9659566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(da2, M, N, P));
9669566063dSJacob Faibussowitsch   PetscCall(DMDASetNumProcs(da2, dd->m, dd->n, dd->p));
9679566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(da2, dd->bx, dd->by, dd->bz));
9689566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(da2, dd->w));
9699566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilType(da2, dd->stencil_type));
9709566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(da2, dd->s));
971c73cfb54SMatthew G. Knepley   if (dim == 3) {
97247c6ae99SBarry Smith     PetscInt *lx, *ly, *lz;
9739566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(dd->m, &lx, dd->n, &ly, dd->p, &lz));
9749566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx));
9759566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_y, dd->n, dd->ly, ly));
9769566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_z, dd->p, dd->lz, lz));
9779566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(da2, lx, ly, lz));
9789566063dSJacob Faibussowitsch     PetscCall(PetscFree3(lx, ly, lz));
979c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
98047c6ae99SBarry Smith     PetscInt *lx, *ly;
9819566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(dd->m, &lx, dd->n, &ly));
9829566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx));
9839566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_y, dd->n, dd->ly, ly));
9849566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(da2, lx, ly, NULL));
9859566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lx, ly));
986c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
98747c6ae99SBarry Smith     PetscInt *lx;
9889566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->m, &lx));
9899566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx));
9909566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(da2, lx, NULL, NULL));
9919566063dSJacob Faibussowitsch     PetscCall(PetscFree(lx));
99247c6ae99SBarry Smith   }
99347c6ae99SBarry Smith   dd2 = (DM_DA *)da2->data;
99447c6ae99SBarry Smith 
99547c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
99625296bd5SBarry Smith   da2->ops->creatematrix = da->ops->creatematrix;
99725296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
99847c6ae99SBarry Smith   da2->ops->getcoloring  = da->ops->getcoloring;
99947c6ae99SBarry Smith   dd2->interptype        = dd->interptype;
100047c6ae99SBarry Smith 
100147c6ae99SBarry Smith   /* copy fill information if given */
100247c6ae99SBarry Smith   if (dd->dfill) {
10039566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->dfill[dd->w] + dd->w + 1, &dd2->dfill));
10049566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->dfill, dd->dfill, dd->dfill[dd->w] + dd->w + 1));
100547c6ae99SBarry Smith   }
100647c6ae99SBarry Smith   if (dd->ofill) {
10079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->ofill[dd->w] + dd->w + 1, &dd2->ofill));
10089566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->ofill, dd->ofill, dd->ofill[dd->w] + dd->w + 1));
100947c6ae99SBarry Smith   }
101047c6ae99SBarry Smith   /* copy the refine information */
1011397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
1012397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
1013397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
101447c6ae99SBarry Smith 
1015897f7067SBarry Smith   if (dd->refine_z_hier) {
10169371c9d4SSatish Balay     if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_z_hier_n) { dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown + 1]; }
10179371c9d4SSatish Balay     if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_z_hier_n) { dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown]; }
1018897f7067SBarry Smith     dd2->refine_z_hier_n = dd->refine_z_hier_n;
10199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_z_hier_n, &dd2->refine_z_hier));
10209566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_z_hier, dd->refine_z_hier, dd2->refine_z_hier_n));
1021897f7067SBarry Smith   }
1022897f7067SBarry Smith   if (dd->refine_y_hier) {
10239371c9d4SSatish Balay     if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_y_hier_n) { dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown + 1]; }
10249371c9d4SSatish Balay     if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_y_hier_n) { dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown]; }
1025897f7067SBarry Smith     dd2->refine_y_hier_n = dd->refine_y_hier_n;
10269566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_y_hier_n, &dd2->refine_y_hier));
10279566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_y_hier, dd->refine_y_hier, dd2->refine_y_hier_n));
1028897f7067SBarry Smith   }
1029897f7067SBarry Smith   if (dd->refine_x_hier) {
10309371c9d4SSatish Balay     if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_x_hier_n) { dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown + 1]; }
10319371c9d4SSatish Balay     if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_x_hier_n) { dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown]; }
1032897f7067SBarry Smith     dd2->refine_x_hier_n = dd->refine_x_hier_n;
10339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_x_hier_n, &dd2->refine_x_hier));
10349566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_x_hier, dd->refine_x_hier, dd2->refine_x_hier_n));
1035897f7067SBarry Smith   }
1036897f7067SBarry Smith 
103747c6ae99SBarry Smith   /* copy vector type information */
10389566063dSJacob Faibussowitsch   PetscCall(DMSetVecType(da2, da->vectype));
1039ddcf8b74SDave May 
1040efd51863SBarry Smith   dd2->lf = dd->lf;
1041efd51863SBarry Smith   dd2->lj = dd->lj;
1042efd51863SBarry Smith 
10436e87535bSJed Brown   da2->leveldown = da->leveldown;
10446e87535bSJed Brown   da2->levelup   = da->levelup + 1;
10458865f1eaSKarl Rupp 
10469566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da2));
10476e87535bSJed Brown 
1048ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
10496858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(da, &coordsc));
10506858538eSMatthew G. Knepley   if (coordsc) {
1051ddcf8b74SDave May     DM  cdaf, cdac;
1052ddcf8b74SDave May     Mat II;
1053ddcf8b74SDave May 
10549566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(da, &cdac));
10559566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(da2, &cdaf));
1056b61d3410SDave May     /* force creation of the coordinate vector */
10579566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(da2, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
10589566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(da2, &coordsf));
10599566063dSJacob Faibussowitsch     PetscCall(DMCreateInterpolation(cdac, cdaf, &II, NULL));
10609566063dSJacob Faibussowitsch     PetscCall(MatInterpolate(II, coordsc, coordsf));
10619566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&II));
1062ddcf8b74SDave May   }
1063397b6216SJed Brown 
1064f3141302SJed Brown   for (i = 0; i < da->bs; i++) {
1065f3141302SJed Brown     const char *fieldname;
10669566063dSJacob Faibussowitsch     PetscCall(DMDAGetFieldName(da, i, &fieldname));
10679566063dSJacob Faibussowitsch     PetscCall(DMDASetFieldName(da2, i, fieldname));
1068f3141302SJed Brown   }
1069397b6216SJed Brown 
107047c6ae99SBarry Smith   *daref = da2;
107147c6ae99SBarry Smith   PetscFunctionReturn(0);
107247c6ae99SBarry Smith }
107347c6ae99SBarry Smith 
10749371c9d4SSatish Balay PetscErrorCode DMCoarsen_DA(DM dmf, MPI_Comm comm, DM *dmc) {
1075c73cfb54SMatthew G. Knepley   PetscInt M, N, P, i, dim;
10766858538eSMatthew G. Knepley   Vec      coordsc, coordsf;
1077a5bc1bf3SBarry Smith   DM       dmc2;
1078a5bc1bf3SBarry Smith   DM_DA   *dd = (DM_DA *)dmf->data, *dd2;
107947c6ae99SBarry Smith 
108047c6ae99SBarry Smith   PetscFunctionBegin;
1081a5bc1bf3SBarry Smith   PetscValidHeaderSpecificType(dmf, DM_CLASSID, 1, DMDA);
1082a5bc1bf3SBarry Smith   PetscValidPointer(dmc, 3);
108347c6ae99SBarry Smith 
10849566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dmf, &dim));
1085bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1086397b6216SJed Brown     M = dd->M / dd->coarsen_x;
108747c6ae99SBarry Smith   } else {
1088397b6216SJed Brown     M = 1 + (dd->M - 1) / dd->coarsen_x;
108947c6ae99SBarry Smith   }
1090bff4a2f0SMatthew G. Knepley   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1091c73cfb54SMatthew G. Knepley     if (dim > 1) {
1092397b6216SJed Brown       N = dd->N / dd->coarsen_y;
109347c6ae99SBarry Smith     } else {
10941860e6e9SBarry Smith       N = 1;
10951860e6e9SBarry Smith     }
10961860e6e9SBarry Smith   } else {
1097397b6216SJed Brown     N = 1 + (dd->N - 1) / dd->coarsen_y;
109847c6ae99SBarry Smith   }
1099bff4a2f0SMatthew G. Knepley   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1100c73cfb54SMatthew G. Knepley     if (dim > 2) {
1101397b6216SJed Brown       P = dd->P / dd->coarsen_z;
110247c6ae99SBarry Smith     } else {
11031860e6e9SBarry Smith       P = 1;
11041860e6e9SBarry Smith     }
11051860e6e9SBarry Smith   } else {
1106397b6216SJed Brown     P = 1 + (dd->P - 1) / dd->coarsen_z;
110747c6ae99SBarry Smith   }
11089566063dSJacob Faibussowitsch   PetscCall(DMDACreate(PetscObjectComm((PetscObject)dmf), &dmc2));
11099566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(dmc2, ((PetscObject)dmf)->prefix));
11109566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dmc2, dim));
11119566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(dmc2, M, N, P));
11129566063dSJacob Faibussowitsch   PetscCall(DMDASetNumProcs(dmc2, dd->m, dd->n, dd->p));
11139566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(dmc2, dd->bx, dd->by, dd->bz));
11149566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(dmc2, dd->w));
11159566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilType(dmc2, dd->stencil_type));
11169566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(dmc2, dd->s));
1117c73cfb54SMatthew G. Knepley   if (dim == 3) {
11182be375d4SJed Brown     PetscInt *lx, *ly, *lz;
11199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(dd->m, &lx, dd->n, &ly, dd->p, &lz));
11209566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx));
11219566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_y, dd->n, dd->ly, ly));
11229566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_z, dd->p, dd->lz, lz));
11239566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(dmc2, lx, ly, lz));
11249566063dSJacob Faibussowitsch     PetscCall(PetscFree3(lx, ly, lz));
1125c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
11262be375d4SJed Brown     PetscInt *lx, *ly;
11279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(dd->m, &lx, dd->n, &ly));
11289566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx));
11299566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_y, dd->n, dd->ly, ly));
11309566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(dmc2, lx, ly, NULL));
11319566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lx, ly));
1132c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
11332be375d4SJed Brown     PetscInt *lx;
11349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->m, &lx));
11359566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx));
11369566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(dmc2, lx, NULL, NULL));
11379566063dSJacob Faibussowitsch     PetscCall(PetscFree(lx));
113847c6ae99SBarry Smith   }
1139a5bc1bf3SBarry Smith   dd2 = (DM_DA *)dmc2->data;
114047c6ae99SBarry Smith 
11414dcab191SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
1142a5bc1bf3SBarry Smith   /* dmc2->ops->createinterpolation = dmf->ops->createinterpolation; copying this one causes trouble for DMSetVI */
1143a5bc1bf3SBarry Smith   dmc2->ops->creatematrix = dmf->ops->creatematrix;
1144a5bc1bf3SBarry Smith   dmc2->ops->getcoloring  = dmf->ops->getcoloring;
114547c6ae99SBarry Smith   dd2->interptype         = dd->interptype;
114647c6ae99SBarry Smith 
114747c6ae99SBarry Smith   /* copy fill information if given */
114847c6ae99SBarry Smith   if (dd->dfill) {
11499566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->dfill[dd->w] + dd->w + 1, &dd2->dfill));
11509566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->dfill, dd->dfill, dd->dfill[dd->w] + dd->w + 1));
115147c6ae99SBarry Smith   }
115247c6ae99SBarry Smith   if (dd->ofill) {
11539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->ofill[dd->w] + dd->w + 1, &dd2->ofill));
11549566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->ofill, dd->ofill, dd->ofill[dd->w] + dd->w + 1));
115547c6ae99SBarry Smith   }
115647c6ae99SBarry Smith   /* copy the refine information */
1157397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
1158397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
1159397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
116047c6ae99SBarry Smith 
1161897f7067SBarry Smith   if (dd->refine_z_hier) {
11629371c9d4SSatish Balay     if (dmf->levelup - dmf->leveldown - 1 > -1 && dmf->levelup - dmf->leveldown - 1 < dd->refine_z_hier_n) { dd2->refine_z = dd->refine_z_hier[dmf->levelup - dmf->leveldown - 1]; }
11639371c9d4SSatish Balay     if (dmf->levelup - dmf->leveldown - 2 > -1 && dmf->levelup - dmf->leveldown - 2 < dd->refine_z_hier_n) { dd2->coarsen_z = dd->refine_z_hier[dmf->levelup - dmf->leveldown - 2]; }
1164897f7067SBarry Smith     dd2->refine_z_hier_n = dd->refine_z_hier_n;
11659566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_z_hier_n, &dd2->refine_z_hier));
11669566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_z_hier, dd->refine_z_hier, dd2->refine_z_hier_n));
1167897f7067SBarry Smith   }
1168897f7067SBarry Smith   if (dd->refine_y_hier) {
11699371c9d4SSatish Balay     if (dmf->levelup - dmf->leveldown - 1 > -1 && dmf->levelup - dmf->leveldown - 1 < dd->refine_y_hier_n) { dd2->refine_y = dd->refine_y_hier[dmf->levelup - dmf->leveldown - 1]; }
11709371c9d4SSatish Balay     if (dmf->levelup - dmf->leveldown - 2 > -1 && dmf->levelup - dmf->leveldown - 2 < dd->refine_y_hier_n) { dd2->coarsen_y = dd->refine_y_hier[dmf->levelup - dmf->leveldown - 2]; }
1171897f7067SBarry Smith     dd2->refine_y_hier_n = dd->refine_y_hier_n;
11729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_y_hier_n, &dd2->refine_y_hier));
11739566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_y_hier, dd->refine_y_hier, dd2->refine_y_hier_n));
1174897f7067SBarry Smith   }
1175897f7067SBarry Smith   if (dd->refine_x_hier) {
11769371c9d4SSatish Balay     if (dmf->levelup - dmf->leveldown - 1 > -1 && dmf->levelup - dmf->leveldown - 1 < dd->refine_x_hier_n) { dd2->refine_x = dd->refine_x_hier[dmf->levelup - dmf->leveldown - 1]; }
11779371c9d4SSatish Balay     if (dmf->levelup - dmf->leveldown - 2 > -1 && dmf->levelup - dmf->leveldown - 2 < dd->refine_x_hier_n) { dd2->coarsen_x = dd->refine_x_hier[dmf->levelup - dmf->leveldown - 2]; }
1178897f7067SBarry Smith     dd2->refine_x_hier_n = dd->refine_x_hier_n;
11799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_x_hier_n, &dd2->refine_x_hier));
11809566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_x_hier, dd->refine_x_hier, dd2->refine_x_hier_n));
1181897f7067SBarry Smith   }
1182897f7067SBarry Smith 
118347c6ae99SBarry Smith   /* copy vector type information */
11849566063dSJacob Faibussowitsch   PetscCall(DMSetVecType(dmc2, dmf->vectype));
118547c6ae99SBarry Smith 
1186644e2e5bSBarry Smith   dd2->lf = dd->lf;
1187644e2e5bSBarry Smith   dd2->lj = dd->lj;
1188644e2e5bSBarry Smith 
1189a5bc1bf3SBarry Smith   dmc2->leveldown = dmf->leveldown + 1;
1190a5bc1bf3SBarry Smith   dmc2->levelup   = dmf->levelup;
11918865f1eaSKarl Rupp 
11929566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmc2));
11936e87535bSJed Brown 
1194ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
11956858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dmf, &coordsf));
11966858538eSMatthew G. Knepley   if (coordsf) {
1197ddcf8b74SDave May     DM         cdaf, cdac;
11986dbf9973SLawrence Mitchell     Mat        inject;
11996dbf9973SLawrence Mitchell     VecScatter vscat;
1200ddcf8b74SDave May 
12019566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmf, &cdaf));
12029566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmc2, &cdac));
1203b61d3410SDave May     /* force creation of the coordinate vector */
12049566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(dmc2, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
12059566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmc2, &coordsc));
1206ddcf8b74SDave May 
12079566063dSJacob Faibussowitsch     PetscCall(DMCreateInjection(cdac, cdaf, &inject));
12089566063dSJacob Faibussowitsch     PetscCall(MatScatterGetVecScatter(inject, &vscat));
12099566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(vscat, coordsf, coordsc, INSERT_VALUES, SCATTER_FORWARD));
12109566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(vscat, coordsf, coordsc, INSERT_VALUES, SCATTER_FORWARD));
12119566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&inject));
1212ddcf8b74SDave May   }
1213f98405f7SJed Brown 
1214a5bc1bf3SBarry Smith   for (i = 0; i < dmf->bs; i++) {
1215f98405f7SJed Brown     const char *fieldname;
12169566063dSJacob Faibussowitsch     PetscCall(DMDAGetFieldName(dmf, i, &fieldname));
12179566063dSJacob Faibussowitsch     PetscCall(DMDASetFieldName(dmc2, i, fieldname));
1218f98405f7SJed Brown   }
1219f98405f7SJed Brown 
1220a5bc1bf3SBarry Smith   *dmc = dmc2;
122147c6ae99SBarry Smith   PetscFunctionReturn(0);
122247c6ae99SBarry Smith }
122347c6ae99SBarry Smith 
12249371c9d4SSatish Balay PetscErrorCode DMRefineHierarchy_DA(DM da, PetscInt nlevels, DM daf[]) {
122547c6ae99SBarry Smith   PetscInt i, n, *refx, *refy, *refz;
122647c6ae99SBarry Smith 
122747c6ae99SBarry Smith   PetscFunctionBegin;
122847c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
12297a8be351SBarry Smith   PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative");
123047c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
123147c6ae99SBarry Smith   PetscValidPointer(daf, 3);
123247c6ae99SBarry Smith 
1233aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
12349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(nlevels, &refx, nlevels, &refy, nlevels, &refz));
1235*48a46eb9SPierre Jolivet   for (i = 0; i < nlevels; i++) PetscCall(DMDAGetRefinementFactor(da, &refx[i], &refy[i], &refz[i]));
123647c6ae99SBarry Smith   n = nlevels;
12379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_x", refx, &n, NULL));
123847c6ae99SBarry Smith   n = nlevels;
12399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_y", refy, &n, NULL));
124047c6ae99SBarry Smith   n = nlevels;
12419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_z", refz, &n, NULL));
124247c6ae99SBarry Smith 
12439566063dSJacob Faibussowitsch   PetscCall(DMDASetRefinementFactor(da, refx[0], refy[0], refz[0]));
12449566063dSJacob Faibussowitsch   PetscCall(DMRefine(da, PetscObjectComm((PetscObject)da), &daf[0]));
124547c6ae99SBarry Smith   for (i = 1; i < nlevels; i++) {
12469566063dSJacob Faibussowitsch     PetscCall(DMDASetRefinementFactor(daf[i - 1], refx[i], refy[i], refz[i]));
12479566063dSJacob Faibussowitsch     PetscCall(DMRefine(daf[i - 1], PetscObjectComm((PetscObject)da), &daf[i]));
124847c6ae99SBarry Smith   }
12499566063dSJacob Faibussowitsch   PetscCall(PetscFree3(refx, refy, refz));
125047c6ae99SBarry Smith   PetscFunctionReturn(0);
125147c6ae99SBarry Smith }
125247c6ae99SBarry Smith 
12539371c9d4SSatish Balay PetscErrorCode DMCoarsenHierarchy_DA(DM da, PetscInt nlevels, DM dac[]) {
125447c6ae99SBarry Smith   PetscInt i;
125547c6ae99SBarry Smith 
125647c6ae99SBarry Smith   PetscFunctionBegin;
125747c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
12587a8be351SBarry Smith   PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative");
125947c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
126047c6ae99SBarry Smith   PetscValidPointer(dac, 3);
12619566063dSJacob Faibussowitsch   PetscCall(DMCoarsen(da, PetscObjectComm((PetscObject)da), &dac[0]));
1262*48a46eb9SPierre Jolivet   for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dac[i - 1], PetscObjectComm((PetscObject)da), &dac[i]));
126347c6ae99SBarry Smith   PetscFunctionReturn(0);
126447c6ae99SBarry Smith }
126562197512SBarry Smith 
12669371c9d4SSatish Balay PetscErrorCode DMDASetGLLCoordinates_1d(DM dm, PetscInt n, PetscReal *nodes) {
12678272889dSSatish Balay   PetscInt     i, j, xs, xn, q;
126862197512SBarry Smith   PetscScalar *xx;
126962197512SBarry Smith   PetscReal    h;
127062197512SBarry Smith   Vec          x;
127162197512SBarry Smith   DM_DA       *da = (DM_DA *)dm->data;
127262197512SBarry Smith 
127362197512SBarry Smith   PetscFunctionBegin;
127462197512SBarry Smith   if (da->bx != DM_BOUNDARY_PERIODIC) {
12759566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(dm, NULL, &q, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
127662197512SBarry Smith     q = (q - 1) / (n - 1); /* number of spectral elements */
127762197512SBarry Smith     h = 2.0 / q;
12789566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xn, NULL, NULL));
127962197512SBarry Smith     xs = xs / (n - 1);
128062197512SBarry Smith     xn = xn / (n - 1);
12819566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(dm, -1., 1., 0., 0., 0., 0.));
12829566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dm, &x));
12839566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, x, &xx));
128462197512SBarry Smith 
128562197512SBarry Smith     /* loop over local spectral elements */
128662197512SBarry Smith     for (j = xs; j < xs + xn; j++) {
128762197512SBarry Smith       /*
128862197512SBarry Smith        Except for the first process, each process starts on the second GLL point of the first element on that process
128962197512SBarry Smith        */
12909371c9d4SSatish Balay       for (i = (j == xs && xs > 0) ? 1 : 0; i < n; i++) { xx[j * (n - 1) + i] = -1.0 + h * j + h * (nodes[i] + 1.0) / 2.; }
129162197512SBarry Smith     }
12929566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, x, &xx));
129362197512SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Not yet implemented for periodic");
129462197512SBarry Smith   PetscFunctionReturn(0);
129562197512SBarry Smith }
129662197512SBarry Smith 
12971fd49c25SBarry Smith /*@
129862197512SBarry Smith 
129962197512SBarry Smith      DMDASetGLLCoordinates - Sets the global coordinates from -1 to 1 to the GLL points of as many GLL elements that fit the number of grid points
130062197512SBarry Smith 
1301d083f849SBarry Smith    Collective on da
130262197512SBarry Smith 
130362197512SBarry Smith    Input Parameters:
130462197512SBarry Smith +   da - the DMDA object
13058272889dSSatish Balay -   n - the number of GLL nodes
13068272889dSSatish Balay -   nodes - the GLL nodes
130762197512SBarry Smith 
130895452b02SPatrick Sanan    Notes:
130995452b02SPatrick Sanan     the parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points
131062197512SBarry Smith           on each process much be divisible by the number of GLL elements needed per process. This depends on whether the DM is
131162197512SBarry Smith           periodic or not.
131262197512SBarry Smith 
1313edc382c3SSatish Balay    Level: advanced
1314edc382c3SSatish Balay 
1315db781477SPatrick Sanan .seealso: `DMDACreate()`, `PetscDTGaussLobattoLegendreQuadrature()`, `DMGetCoordinates()`
131662197512SBarry Smith @*/
13179371c9d4SSatish Balay PetscErrorCode DMDASetGLLCoordinates(DM da, PetscInt n, PetscReal *nodes) {
131862197512SBarry Smith   PetscFunctionBegin;
131962197512SBarry Smith   if (da->dim == 1) {
13209566063dSJacob Faibussowitsch     PetscCall(DMDASetGLLCoordinates_1d(da, n, nodes));
132162197512SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Not yet implemented for 2 or 3d");
132262197512SBarry Smith   PetscFunctionReturn(0);
132362197512SBarry Smith }
13247c3cd84eSPatrick Sanan 
13259371c9d4SSatish Balay PETSC_INTERN PetscErrorCode DMGetCompatibility_DA(DM da1, DM dm2, PetscBool *compatible, PetscBool *set) {
13267c3cd84eSPatrick Sanan   DM_DA    *dd1 = (DM_DA *)da1->data, *dd2;
13277c3cd84eSPatrick Sanan   DM        da2;
13287c3cd84eSPatrick Sanan   DMType    dmtype2;
13297c3cd84eSPatrick Sanan   PetscBool isda, compatibleLocal;
13307c3cd84eSPatrick Sanan   PetscInt  i;
13317c3cd84eSPatrick Sanan 
13327c3cd84eSPatrick Sanan   PetscFunctionBegin;
13337a8be351SBarry Smith   PetscCheck(da1->setupcalled, PetscObjectComm((PetscObject)da1), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() must be called on first DM before DMGetCompatibility()");
13349566063dSJacob Faibussowitsch   PetscCall(DMGetType(dm2, &dmtype2));
13359566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(dmtype2, DMDA, &isda));
13367c3cd84eSPatrick Sanan   if (isda) {
13377c3cd84eSPatrick Sanan     da2 = dm2;
13387c3cd84eSPatrick Sanan     dd2 = (DM_DA *)da2->data;
13397a8be351SBarry Smith     PetscCheck(da2->setupcalled, PetscObjectComm((PetscObject)da2), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() must be called on second DM before DMGetCompatibility()");
1340110623a0SKarl Rupp     compatibleLocal = (PetscBool)(da1->dim == da2->dim);
1341c790a739SKarl Rupp     if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->s == dd2->s)); /* Stencil width */
1342110623a0SKarl Rupp     /*                                                                           Global size              ranks               Boundary type */
1343c790a739SKarl Rupp     if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->M == dd2->M) && (dd1->m == dd2->m) && (dd1->bx == dd2->bx));
1344c790a739SKarl Rupp     if (compatibleLocal && da1->dim > 1) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->N == dd2->N) && (dd1->n == dd2->n) && (dd1->by == dd2->by));
1345c790a739SKarl Rupp     if (compatibleLocal && da1->dim > 2) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->P == dd2->P) && (dd1->p == dd2->p) && (dd1->bz == dd2->bz));
13467c3cd84eSPatrick Sanan     if (compatibleLocal) {
13479371c9d4SSatish Balay       for (i = 0; i < dd1->m; ++i) { compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lx[i] == dd2->lx[i])); /* Local size     */ }
13487c3cd84eSPatrick Sanan     }
13497c3cd84eSPatrick Sanan     if (compatibleLocal && da1->dim > 1) {
13509371c9d4SSatish Balay       for (i = 0; i < dd1->n; ++i) { compatibleLocal = (PetscBool)(compatibleLocal && (dd1->ly[i] == dd2->ly[i])); }
13517c3cd84eSPatrick Sanan     }
13527c3cd84eSPatrick Sanan     if (compatibleLocal && da1->dim > 2) {
13539371c9d4SSatish Balay       for (i = 0; i < dd1->p; ++i) { compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lz[i] == dd2->lz[i])); }
13547c3cd84eSPatrick Sanan     }
13557c3cd84eSPatrick Sanan     *compatible = compatibleLocal;
13567c3cd84eSPatrick Sanan     *set        = PETSC_TRUE;
13577c3cd84eSPatrick Sanan   } else {
13587c3cd84eSPatrick Sanan     /* Decline to determine compatibility with other DM types */
13597c3cd84eSPatrick Sanan     *set = PETSC_FALSE;
13607c3cd84eSPatrick Sanan   }
13617c3cd84eSPatrick Sanan   PetscFunctionReturn(0);
13627c3cd84eSPatrick Sanan }
1363