xref: /petsc/src/dm/impls/da/da.c (revision 12b4a53753ecbae42c98ba33876a303b79054923)
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 
620f4b53cSBarry Smith   Logically Collective
747c6ae99SBarry Smith 
847c6ae99SBarry Smith   Input Parameters:
9dce8aebaSBarry 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 
16*12b4a537SBarry Smith   Developer Note:
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 
19*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `PetscSplitOwnership()`
2047c6ae99SBarry Smith @*/
21d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P)
22d71ae5a4SJacob Faibussowitsch {
2347c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
2447c6ae99SBarry Smith 
2547c6ae99SBarry Smith   PetscFunctionBegin;
26a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
2747c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, M, 2);
2847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, N, 3);
2947c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, P, 4);
307a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
317a8be351SBarry Smith   PetscCheck(M >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Number of grid points in X direction must be positive");
327a8be351SBarry Smith   PetscCheck(N >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Number of grid points in Y direction must be positive");
337a8be351SBarry Smith   PetscCheck(P >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Number of grid points in Z direction must be positive");
3447c6ae99SBarry Smith 
3547c6ae99SBarry Smith   dd->M = M;
3647c6ae99SBarry Smith   dd->N = N;
3747c6ae99SBarry Smith   dd->P = P;
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3947c6ae99SBarry Smith }
4047c6ae99SBarry Smith 
4147c6ae99SBarry Smith /*@
42aa219208SBarry Smith   DMDASetNumProcs - Sets the number of processes in each dimension
4347c6ae99SBarry Smith 
4420f4b53cSBarry Smith   Logically Collective
4547c6ae99SBarry Smith 
4647c6ae99SBarry Smith   Input Parameters:
47dce8aebaSBarry Smith + da - the `DMDA`
48*12b4a537SBarry Smith . m  - the number of X processes (or `PETSC_DECIDE`)
49*12b4a537SBarry Smith . n  - the number of Y processes (or `PETSC_DECIDE`)
50*12b4a537SBarry Smith - p  - the number of Z processes (or `PETSC_DECIDE`)
5147c6ae99SBarry Smith 
5247c6ae99SBarry Smith   Level: intermediate
5347c6ae99SBarry Smith 
54*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDASetSizes()`, `PetscSplitOwnership()`
5547c6ae99SBarry Smith @*/
56d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p)
57d71ae5a4SJacob Faibussowitsch {
5847c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
5947c6ae99SBarry Smith 
6047c6ae99SBarry Smith   PetscFunctionBegin;
61a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
6247c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, m, 2);
6347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, n, 3);
6447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, p, 4);
657a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
6647c6ae99SBarry Smith   dd->m = m;
6747c6ae99SBarry Smith   dd->n = n;
6847c6ae99SBarry Smith   dd->p = p;
69c73cfb54SMatthew G. Knepley   if (da->dim == 2) {
70d3be247eSBarry Smith     PetscMPIInt size;
719566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da), &size));
72e3f3e4b6SBarry Smith     if ((dd->m > 0) && (dd->n < 0)) {
73e3f3e4b6SBarry Smith       dd->n = size / dd->m;
7463a3b9bcSJacob 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);
75e3f3e4b6SBarry Smith     }
76e3f3e4b6SBarry Smith     if ((dd->n > 0) && (dd->m < 0)) {
77e3f3e4b6SBarry Smith       dd->m = size / dd->n;
7863a3b9bcSJacob 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);
79e3f3e4b6SBarry Smith     }
80e3f3e4b6SBarry Smith   }
813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8247c6ae99SBarry Smith }
8347c6ae99SBarry Smith 
8447c6ae99SBarry Smith /*@
851321219cSEthan Coon   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries.
8647c6ae99SBarry Smith 
8720f4b53cSBarry Smith   Not Collective
8847c6ae99SBarry Smith 
89d8d19677SJose E. Roman   Input Parameters:
90dce8aebaSBarry Smith + da - The `DMDA`
91a4e35b19SJacob Faibussowitsch . bx - x boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`
92a4e35b19SJacob Faibussowitsch . by - y boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`
93a4e35b19SJacob Faibussowitsch - bz - z boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`
9447c6ae99SBarry Smith 
9547c6ae99SBarry Smith   Level: intermediate
9647c6ae99SBarry Smith 
97*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMBoundaryType`, `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`
9847c6ae99SBarry Smith @*/
99d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetBoundaryType(DM da, DMBoundaryType bx, DMBoundaryType by, DMBoundaryType bz)
100d71ae5a4SJacob Faibussowitsch {
10147c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
10247c6ae99SBarry Smith 
10347c6ae99SBarry Smith   PetscFunctionBegin;
104a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1055a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da, bx, 2);
1065a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da, by, 3);
1075a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da, bz, 4);
1087a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
1091321219cSEthan Coon   dd->bx = bx;
1101321219cSEthan Coon   dd->by = by;
1111321219cSEthan Coon   dd->bz = bz;
1123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11347c6ae99SBarry Smith }
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith /*@
116aa219208SBarry Smith   DMDASetDof - Sets the number of degrees of freedom per vertex
11747c6ae99SBarry Smith 
11820f4b53cSBarry Smith   Not Collective
11947c6ae99SBarry Smith 
12059f3ab6dSMatthew G. Knepley   Input Parameters:
121dce8aebaSBarry Smith + da  - The `DMDA`
122*12b4a537SBarry Smith - dof - Number of degrees of freedom per vertex
12347c6ae99SBarry Smith 
12447c6ae99SBarry Smith   Level: intermediate
12547c6ae99SBarry Smith 
126*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetDof()`, `DMDACreate()`, `DMDestroy()`
12747c6ae99SBarry Smith @*/
128d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetDof(DM da, PetscInt dof)
129d71ae5a4SJacob Faibussowitsch {
13047c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith   PetscFunctionBegin;
133a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
13454cfb0beSLisandro Dalcin   PetscValidLogicalCollectiveInt(da, dof, 2);
1357a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
13647c6ae99SBarry Smith   dd->w  = dof;
1371411c6eeSJed Brown   da->bs = dof;
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13947c6ae99SBarry Smith }
14047c6ae99SBarry Smith 
141fb6725baSMatthew G. Knepley /*@
142fb6725baSMatthew G. Knepley   DMDAGetDof - Gets the number of degrees of freedom per vertex
143fb6725baSMatthew G. Knepley 
14420f4b53cSBarry Smith   Not Collective
145fb6725baSMatthew G. Knepley 
146fb6725baSMatthew G. Knepley   Input Parameter:
147dce8aebaSBarry Smith . da - The `DMDA`
148fb6725baSMatthew G. Knepley 
149fb6725baSMatthew G. Knepley   Output Parameter:
150*12b4a537SBarry Smith . dof - Number of degrees of freedom per vertex
151fb6725baSMatthew G. Knepley 
152fb6725baSMatthew G. Knepley   Level: intermediate
153fb6725baSMatthew G. Knepley 
154*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDASetDof()`, `DMDACreate()`, `DMDestroy()`
155fb6725baSMatthew G. Knepley @*/
156d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetDof(DM da, PetscInt *dof)
157d71ae5a4SJacob Faibussowitsch {
158fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *)da->data;
159fb6725baSMatthew G. Knepley 
160fb6725baSMatthew G. Knepley   PetscFunctionBegin;
161a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1624f572ea9SToby Isaac   PetscAssertPointer(dof, 2);
163fb6725baSMatthew G. Knepley   *dof = dd->w;
1643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165fb6725baSMatthew G. Knepley }
166fb6725baSMatthew G. Knepley 
1677ddda789SPeter Brune /*@
1687ddda789SPeter Brune   DMDAGetOverlap - Gets the size of the per-processor overlap.
1697ddda789SPeter Brune 
17020f4b53cSBarry Smith   Not Collective
1717ddda789SPeter Brune 
172f899ff85SJose E. Roman   Input Parameter:
173dce8aebaSBarry Smith . da - The `DMDA`
1747ddda789SPeter Brune 
1757ddda789SPeter Brune   Output Parameters:
1767ddda789SPeter Brune + x - Overlap in the x direction
1777ddda789SPeter Brune . y - Overlap in the y direction
1787ddda789SPeter Brune - z - Overlap in the z direction
1797ddda789SPeter Brune 
1807ddda789SPeter Brune   Level: intermediate
1817ddda789SPeter Brune 
182*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDASetOverlap()`
1837ddda789SPeter Brune @*/
184d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetOverlap(DM da, PetscInt *x, PetscInt *y, PetscInt *z)
185d71ae5a4SJacob Faibussowitsch {
1867ddda789SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
1877ddda789SPeter Brune 
1887ddda789SPeter Brune   PetscFunctionBegin;
189a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1907ddda789SPeter Brune   if (x) *x = dd->xol;
1917ddda789SPeter Brune   if (y) *y = dd->yol;
1927ddda789SPeter Brune   if (z) *z = dd->zol;
1933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1947ddda789SPeter Brune }
1957ddda789SPeter Brune 
19688661749SPeter Brune /*@
19788661749SPeter Brune   DMDASetOverlap - Sets the size of the per-processor overlap.
19888661749SPeter Brune 
19920f4b53cSBarry Smith   Not Collective
20088661749SPeter Brune 
2017ddda789SPeter Brune   Input Parameters:
202dce8aebaSBarry Smith + da - The `DMDA`
2037ddda789SPeter Brune . x  - Overlap in the x direction
2047ddda789SPeter Brune . y  - Overlap in the y direction
2057ddda789SPeter Brune - z  - Overlap in the z direction
20688661749SPeter Brune 
20788661749SPeter Brune   Level: intermediate
20888661749SPeter Brune 
209*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDAGetOverlap()`
21088661749SPeter Brune @*/
211d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOverlap(DM da, PetscInt x, PetscInt y, PetscInt z)
212d71ae5a4SJacob Faibussowitsch {
21388661749SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
21488661749SPeter Brune 
21588661749SPeter Brune   PetscFunctionBegin;
216a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
2177ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da, x, 2);
2187ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da, y, 3);
2197ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da, z, 4);
2207ddda789SPeter Brune   dd->xol = x;
2217ddda789SPeter Brune   dd->yol = y;
2227ddda789SPeter Brune   dd->zol = z;
2233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22488661749SPeter Brune }
22588661749SPeter Brune 
2263e7870d2SPeter Brune /*@
227*12b4a537SBarry Smith   DMDAGetNumLocalSubDomains - Gets the number of local subdomains that would be created upon decomposition.
2283e7870d2SPeter Brune 
22920f4b53cSBarry Smith   Not Collective
2303e7870d2SPeter Brune 
2312fe279fdSBarry Smith   Input Parameter:
232dce8aebaSBarry Smith . da - The `DMDA`
2333e7870d2SPeter Brune 
2342fe279fdSBarry Smith   Output Parameter:
235a2b725a8SWilliam Gropp . Nsub - Number of local subdomains created upon decomposition
2363e7870d2SPeter Brune 
2373e7870d2SPeter Brune   Level: intermediate
2383e7870d2SPeter Brune 
239*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDASetNumLocalSubDomains()`
2403e7870d2SPeter Brune @*/
241d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumLocalSubDomains(DM da, PetscInt *Nsub)
242d71ae5a4SJacob Faibussowitsch {
2433e7870d2SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
2443e7870d2SPeter Brune 
2453e7870d2SPeter Brune   PetscFunctionBegin;
246a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
2473e7870d2SPeter Brune   if (Nsub) *Nsub = dd->Nsub;
2483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2493e7870d2SPeter Brune }
2503e7870d2SPeter Brune 
2513e7870d2SPeter Brune /*@
252*12b4a537SBarry Smith   DMDASetNumLocalSubDomains - Sets the number of local subdomains to create when decomposing with `DMCreateDomainDecomposition()`
2533e7870d2SPeter Brune 
25420f4b53cSBarry Smith   Not Collective
2553e7870d2SPeter Brune 
2563e7870d2SPeter Brune   Input Parameters:
257dce8aebaSBarry Smith + da   - The `DMDA`
2583e7870d2SPeter Brune - Nsub - The number of local subdomains requested
2593e7870d2SPeter Brune 
2603e7870d2SPeter Brune   Level: intermediate
2613e7870d2SPeter Brune 
262*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDAGetNumLocalSubDomains()`
2633e7870d2SPeter Brune @*/
264d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetNumLocalSubDomains(DM da, PetscInt Nsub)
265d71ae5a4SJacob Faibussowitsch {
2663e7870d2SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
2673e7870d2SPeter Brune 
2683e7870d2SPeter Brune   PetscFunctionBegin;
269a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
2703e7870d2SPeter Brune   PetscValidLogicalCollectiveInt(da, Nsub, 2);
2713e7870d2SPeter Brune   dd->Nsub = Nsub;
2723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2733e7870d2SPeter Brune }
2743e7870d2SPeter Brune 
275d886c4f4SPeter Brune /*@
276*12b4a537SBarry Smith   DMDASetOffset - Sets the index offset of the `DMDA`.
277d886c4f4SPeter Brune 
27820f4b53cSBarry Smith   Collective
279d886c4f4SPeter Brune 
280d8d19677SJose E. Roman   Input Parameters:
281dce8aebaSBarry Smith + da - The `DMDA`
282d886c4f4SPeter Brune . xo - The offset in the x direction
283d886c4f4SPeter Brune . yo - The offset in the y direction
284a4e35b19SJacob Faibussowitsch . zo - The offset in the z direction
285a4e35b19SJacob Faibussowitsch . Mo - The problem offset in the x direction
286a4e35b19SJacob Faibussowitsch . No - The problem offset in the y direction
287a4e35b19SJacob Faibussowitsch - Po - The problem offset in the z direction
288d886c4f4SPeter Brune 
289*12b4a537SBarry Smith   Level: developer
290d886c4f4SPeter Brune 
291dce8aebaSBarry Smith   Note:
292dce8aebaSBarry Smith   This is used primarily to overlap a computation on a local `DMDA` with that on a global `DMDA` without
293d886c4f4SPeter Brune   changing boundary conditions or subdomain features that depend upon the global offsets.
294d886c4f4SPeter Brune 
295*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()`
296d886c4f4SPeter Brune @*/
297d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po)
298d71ae5a4SJacob Faibussowitsch {
299d886c4f4SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
300d886c4f4SPeter Brune 
301d886c4f4SPeter Brune   PetscFunctionBegin;
302a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
303d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da, xo, 2);
30495c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, yo, 3);
30595c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, zo, 4);
30695c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, Mo, 5);
30795c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, No, 6);
30895c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, Po, 7);
309d886c4f4SPeter Brune   dd->xo = xo;
310d886c4f4SPeter Brune   dd->yo = yo;
311d886c4f4SPeter Brune   dd->zo = zo;
31295c13181SPeter Brune   dd->Mo = Mo;
31395c13181SPeter Brune   dd->No = No;
31495c13181SPeter Brune   dd->Po = Po;
31595c13181SPeter Brune 
3166858538eSMatthew G. Knepley   if (da->coordinates[0].dm) PetscCall(DMDASetOffset(da->coordinates[0].dm, xo, yo, zo, Mo, No, Po));
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
318d886c4f4SPeter Brune }
319d886c4f4SPeter Brune 
320d886c4f4SPeter Brune /*@
321dce8aebaSBarry Smith   DMDAGetOffset - Gets the index offset of the `DMDA`.
322d886c4f4SPeter Brune 
32320f4b53cSBarry Smith   Not Collective
324d886c4f4SPeter Brune 
325d886c4f4SPeter Brune   Input Parameter:
326dce8aebaSBarry Smith . da - The `DMDA`
327d886c4f4SPeter Brune 
328d886c4f4SPeter Brune   Output Parameters:
329d886c4f4SPeter Brune + xo - The offset in the x direction
330d886c4f4SPeter Brune . yo - The offset in the y direction
33195c13181SPeter Brune . zo - The offset in the z direction
33295c13181SPeter Brune . Mo - The global size in the x direction
33395c13181SPeter Brune . No - The global size in the y direction
33495c13181SPeter Brune - Po - The global size in the z direction
335d886c4f4SPeter Brune 
336*12b4a537SBarry Smith   Level: developer
337d886c4f4SPeter Brune 
338*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDASetOffset()`, `DMDAVecGetArray()`
339d886c4f4SPeter Brune @*/
340d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetOffset(DM da, PetscInt *xo, PetscInt *yo, PetscInt *zo, PetscInt *Mo, PetscInt *No, PetscInt *Po)
341d71ae5a4SJacob Faibussowitsch {
342d886c4f4SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
343d886c4f4SPeter Brune 
344d886c4f4SPeter Brune   PetscFunctionBegin;
345a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
346d886c4f4SPeter Brune   if (xo) *xo = dd->xo;
347d886c4f4SPeter Brune   if (yo) *yo = dd->yo;
348d886c4f4SPeter Brune   if (zo) *zo = dd->zo;
34995c13181SPeter Brune   if (Mo) *Mo = dd->Mo;
35095c13181SPeter Brune   if (No) *No = dd->No;
35195c13181SPeter Brune   if (Po) *Po = dd->Po;
3523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
353d886c4f4SPeter Brune }
354d886c4f4SPeter Brune 
35540234c92SPeter Brune /*@
356*12b4a537SBarry Smith   DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain `DMDA`.
35740234c92SPeter Brune 
35820f4b53cSBarry Smith   Not Collective
35940234c92SPeter Brune 
36040234c92SPeter Brune   Input Parameter:
361dce8aebaSBarry Smith . da - The `DMDA`
36240234c92SPeter Brune 
36340234c92SPeter Brune   Output Parameters:
36440234c92SPeter Brune + xs - The start of the region in x
36540234c92SPeter Brune . ys - The start of the region in y
36640234c92SPeter Brune . zs - The start of the region in z
367a4e35b19SJacob Faibussowitsch . xm - The size of the region in x
368a4e35b19SJacob Faibussowitsch . ym - The size of the region in y
369a4e35b19SJacob Faibussowitsch - zm - The size of the region in z
37040234c92SPeter Brune 
37140234c92SPeter Brune   Level: intermediate
37240234c92SPeter Brune 
373*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()`
37440234c92SPeter Brune @*/
375d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm)
376d71ae5a4SJacob Faibussowitsch {
37740234c92SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
37840234c92SPeter Brune 
37940234c92SPeter Brune   PetscFunctionBegin;
380a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
38140234c92SPeter Brune   if (xs) *xs = dd->nonxs;
38240234c92SPeter Brune   if (ys) *ys = dd->nonys;
38340234c92SPeter Brune   if (zs) *zs = dd->nonzs;
38440234c92SPeter Brune   if (xm) *xm = dd->nonxm;
38540234c92SPeter Brune   if (ym) *ym = dd->nonym;
38640234c92SPeter Brune   if (zm) *zm = dd->nonzm;
3873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38840234c92SPeter Brune }
38940234c92SPeter Brune 
39040234c92SPeter Brune /*@
391*12b4a537SBarry Smith   DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain `DMDA`.
39240234c92SPeter Brune 
39320f4b53cSBarry Smith   Collective
39440234c92SPeter Brune 
395d8d19677SJose E. Roman   Input Parameters:
396dce8aebaSBarry Smith + da - The `DMDA`
39740234c92SPeter Brune . xs - The start of the region in x
39840234c92SPeter Brune . ys - The start of the region in y
39940234c92SPeter Brune . zs - The start of the region in z
400a4e35b19SJacob Faibussowitsch . xm - The size of the region in x
401a4e35b19SJacob Faibussowitsch . ym - The size of the region in y
402a4e35b19SJacob Faibussowitsch - zm - The size of the region in z
40340234c92SPeter Brune 
40440234c92SPeter Brune   Level: intermediate
40540234c92SPeter Brune 
406*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()`
40740234c92SPeter Brune @*/
408d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm)
409d71ae5a4SJacob Faibussowitsch {
41040234c92SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
41140234c92SPeter Brune 
41240234c92SPeter Brune   PetscFunctionBegin;
413a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
41440234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, xs, 2);
41540234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, ys, 3);
41640234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, zs, 4);
41740234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, xm, 5);
41840234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, ym, 6);
41940234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, zm, 7);
42040234c92SPeter Brune   dd->nonxs = xs;
42140234c92SPeter Brune   dd->nonys = ys;
42240234c92SPeter Brune   dd->nonzs = zs;
42340234c92SPeter Brune   dd->nonxm = xm;
42440234c92SPeter Brune   dd->nonym = ym;
42540234c92SPeter Brune   dd->nonzm = zm;
42640234c92SPeter Brune 
4273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42840234c92SPeter Brune }
42988661749SPeter Brune 
43047c6ae99SBarry Smith /*@
431aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
43247c6ae99SBarry Smith 
43320f4b53cSBarry Smith   Logically Collective
43447c6ae99SBarry Smith 
435d8d19677SJose E. Roman   Input Parameters:
436dce8aebaSBarry Smith + da    - The `DMDA`
437dce8aebaSBarry Smith - stype - The stencil type, use either `DMDA_STENCIL_BOX` or `DMDA_STENCIL_STAR`.
43847c6ae99SBarry Smith 
43947c6ae99SBarry Smith   Level: intermediate
44047c6ae99SBarry Smith 
441*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMDAStencilType`, `DMDA_STENCIL_BOX`, `DMDA_STENCIL_STAR.`
44247c6ae99SBarry Smith @*/
443d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype)
444d71ae5a4SJacob Faibussowitsch {
44547c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
44647c6ae99SBarry Smith 
44747c6ae99SBarry Smith   PetscFunctionBegin;
448a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
44947c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da, stype, 2);
4507a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
45147c6ae99SBarry Smith   dd->stencil_type = stype;
4523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45347c6ae99SBarry Smith }
45447c6ae99SBarry Smith 
455fb6725baSMatthew G. Knepley /*@
456fb6725baSMatthew G. Knepley   DMDAGetStencilType - Gets the type of the communication stencil
457fb6725baSMatthew G. Knepley 
45820f4b53cSBarry Smith   Not Collective
459fb6725baSMatthew G. Knepley 
460fb6725baSMatthew G. Knepley   Input Parameter:
461dce8aebaSBarry Smith . da - The `DMDA`
462fb6725baSMatthew G. Knepley 
463fb6725baSMatthew G. Knepley   Output Parameter:
464dce8aebaSBarry Smith . stype - The stencil type, use either `DMDA_STENCIL_BOX` or `DMDA_STENCIL_STAR`.
465fb6725baSMatthew G. Knepley 
466fb6725baSMatthew G. Knepley   Level: intermediate
467fb6725baSMatthew G. Knepley 
468*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMDAStencilType`, `DMDA_STENCIL_BOX`, `DMDA_STENCIL_STAR.`
469fb6725baSMatthew G. Knepley @*/
470d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype)
471d71ae5a4SJacob Faibussowitsch {
472fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *)da->data;
473fb6725baSMatthew G. Knepley 
474fb6725baSMatthew G. Knepley   PetscFunctionBegin;
475a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
4764f572ea9SToby Isaac   PetscAssertPointer(stype, 2);
477fb6725baSMatthew G. Knepley   *stype = dd->stencil_type;
4783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
479fb6725baSMatthew G. Knepley }
480fb6725baSMatthew G. Knepley 
48147c6ae99SBarry Smith /*@
482aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
48347c6ae99SBarry Smith 
48420f4b53cSBarry Smith   Logically Collective
48547c6ae99SBarry Smith 
486d8d19677SJose E. Roman   Input Parameters:
487dce8aebaSBarry Smith + da    - The `DMDA`
48847c6ae99SBarry Smith - width - The stencil width
48947c6ae99SBarry Smith 
49047c6ae99SBarry Smith   Level: intermediate
49147c6ae99SBarry Smith 
492*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMDAStencilType`, `DMDA_STENCIL_BOX`, `DMDA_STENCIL_STAR.`
49347c6ae99SBarry Smith @*/
494d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width)
495d71ae5a4SJacob Faibussowitsch {
49647c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
49747c6ae99SBarry Smith 
49847c6ae99SBarry Smith   PetscFunctionBegin;
499a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
50047c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, width, 2);
5017a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
50247c6ae99SBarry Smith   dd->s = width;
5033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50447c6ae99SBarry Smith }
50547c6ae99SBarry Smith 
506fb6725baSMatthew G. Knepley /*@
507fb6725baSMatthew G. Knepley   DMDAGetStencilWidth - Gets the width of the communication stencil
508fb6725baSMatthew G. Knepley 
50920f4b53cSBarry Smith   Not Collective
510fb6725baSMatthew G. Knepley 
511fb6725baSMatthew G. Knepley   Input Parameter:
512dce8aebaSBarry Smith . da - The `DMDA`
513fb6725baSMatthew G. Knepley 
514fb6725baSMatthew G. Knepley   Output Parameter:
515fb6725baSMatthew G. Knepley . width - The stencil width
516fb6725baSMatthew G. Knepley 
517fb6725baSMatthew G. Knepley   Level: intermediate
518fb6725baSMatthew G. Knepley 
519*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMDAStencilType`, `DMDA_STENCIL_BOX`, `DMDA_STENCIL_STAR.`
520fb6725baSMatthew G. Knepley @*/
521d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width)
522d71ae5a4SJacob Faibussowitsch {
523fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *)da->data;
524fb6725baSMatthew G. Knepley 
525fb6725baSMatthew G. Knepley   PetscFunctionBegin;
526a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
5274f572ea9SToby Isaac   PetscAssertPointer(width, 2);
528fb6725baSMatthew G. Knepley   *width = dd->s;
5293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
530fb6725baSMatthew G. Knepley }
531fb6725baSMatthew G. Knepley 
532d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da, PetscInt M, PetscInt m, const PetscInt lx[])
533d71ae5a4SJacob Faibussowitsch {
53447c6ae99SBarry Smith   PetscInt i, sum;
53547c6ae99SBarry Smith 
53647c6ae99SBarry Smith   PetscFunctionBegin;
5377a8be351SBarry Smith   PetscCheck(M >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Global dimension not set");
53847c6ae99SBarry Smith   for (i = sum = 0; i < m; i++) sum += lx[i];
53963a3b9bcSJacob 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);
5403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54147c6ae99SBarry Smith }
54247c6ae99SBarry Smith 
54347c6ae99SBarry Smith /*@
544aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
54547c6ae99SBarry Smith 
54620f4b53cSBarry Smith   Logically Collective
54747c6ae99SBarry Smith 
548d8d19677SJose E. Roman   Input Parameters:
549dce8aebaSBarry Smith + da - The `DMDA`
550*12b4a537SBarry 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
551*12b4a537SBarry 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
552*12b4a537SBarry 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.
55347c6ae99SBarry Smith 
55447c6ae99SBarry Smith   Level: intermediate
55547c6ae99SBarry Smith 
556a4e35b19SJacob Faibussowitsch   Note:
557a4e35b19SJacob Faibussowitsch   These numbers are NOT multiplied by the number of dof per node.
558e3f3e4b6SBarry Smith 
559*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`
56047c6ae99SBarry Smith @*/
561d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
562d71ae5a4SJacob Faibussowitsch {
56347c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
56447c6ae99SBarry Smith 
56547c6ae99SBarry Smith   PetscFunctionBegin;
566a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
5677a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
56847c6ae99SBarry Smith   if (lx) {
569*12b4a537SBarry Smith     PetscCheck(dd->m >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of processes");
5709566063dSJacob Faibussowitsch     PetscCall(DMDACheckOwnershipRanges_Private(da, dd->M, dd->m, lx));
57148a46eb9SPierre Jolivet     if (!dd->lx) PetscCall(PetscMalloc1(dd->m, &dd->lx));
5729566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd->lx, lx, dd->m));
57347c6ae99SBarry Smith   }
57447c6ae99SBarry Smith   if (ly) {
575*12b4a537SBarry Smith     PetscCheck(dd->n >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of processes");
5769566063dSJacob Faibussowitsch     PetscCall(DMDACheckOwnershipRanges_Private(da, dd->N, dd->n, ly));
57748a46eb9SPierre Jolivet     if (!dd->ly) PetscCall(PetscMalloc1(dd->n, &dd->ly));
5789566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd->ly, ly, dd->n));
57947c6ae99SBarry Smith   }
58047c6ae99SBarry Smith   if (lz) {
581*12b4a537SBarry Smith     PetscCheck(dd->p >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of processes");
5829566063dSJacob Faibussowitsch     PetscCall(DMDACheckOwnershipRanges_Private(da, dd->P, dd->p, lz));
58348a46eb9SPierre Jolivet     if (!dd->lz) PetscCall(PetscMalloc1(dd->p, &dd->lz));
5849566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd->lz, lz, dd->p));
58547c6ae99SBarry Smith   }
5863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58747c6ae99SBarry Smith }
58847c6ae99SBarry Smith 
58947c6ae99SBarry Smith /*@
590aa219208SBarry Smith   DMDASetInterpolationType - Sets the type of interpolation that will be
591dce8aebaSBarry Smith   returned by `DMCreateInterpolation()`
59247c6ae99SBarry Smith 
59320f4b53cSBarry Smith   Logically Collective
59447c6ae99SBarry Smith 
595d8d19677SJose E. Roman   Input Parameters:
59647c6ae99SBarry Smith + da    - initial distributed array
597dce8aebaSBarry Smith - ctype - `DMDA_Q1` and `DMDA_Q0` are currently the only supported forms
59847c6ae99SBarry Smith 
59947c6ae99SBarry Smith   Level: intermediate
60047c6ae99SBarry Smith 
601dce8aebaSBarry Smith   Note:
602dce8aebaSBarry Smith   You should call this on the coarser of the two `DMDA` you pass to `DMCreateInterpolation()`
60347c6ae99SBarry Smith 
604*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDestroy()`, `DMDAInterpolationType`,
605*12b4a537SBarry Smith           `DMDA_Q1`, `DMDA_Q0`
60647c6ae99SBarry Smith @*/
607d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetInterpolationType(DM da, DMDAInterpolationType ctype)
608d71ae5a4SJacob Faibussowitsch {
60947c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
61047c6ae99SBarry Smith 
61147c6ae99SBarry Smith   PetscFunctionBegin;
612a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
61347c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da, ctype, 2);
61447c6ae99SBarry Smith   dd->interptype = ctype;
6153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
61647c6ae99SBarry Smith }
61747c6ae99SBarry Smith 
6182dde6fd4SLisandro Dalcin /*@
6192dde6fd4SLisandro Dalcin   DMDAGetInterpolationType - Gets the type of interpolation that will be
620dce8aebaSBarry Smith   used by `DMCreateInterpolation()`
6212dde6fd4SLisandro Dalcin 
6222dde6fd4SLisandro Dalcin   Not Collective
6232dde6fd4SLisandro Dalcin 
6242dde6fd4SLisandro Dalcin   Input Parameter:
6252dde6fd4SLisandro Dalcin . da - distributed array
6262dde6fd4SLisandro Dalcin 
6272dde6fd4SLisandro Dalcin   Output Parameter:
628dce8aebaSBarry Smith . ctype - interpolation type (`DMDA_Q1` and `DMDA_Q0` are currently the only supported forms)
6292dde6fd4SLisandro Dalcin 
6302dde6fd4SLisandro Dalcin   Level: intermediate
6312dde6fd4SLisandro Dalcin 
632*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAInterpolationType`, `DMDASetInterpolationType()`, `DMCreateInterpolation()`,
633*12b4a537SBarry Smith           `DMDA_Q1`, `DMDA_Q0`
6342dde6fd4SLisandro Dalcin @*/
635d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetInterpolationType(DM da, DMDAInterpolationType *ctype)
636d71ae5a4SJacob Faibussowitsch {
6372dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA *)da->data;
6382dde6fd4SLisandro Dalcin 
6392dde6fd4SLisandro Dalcin   PetscFunctionBegin;
640a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
6414f572ea9SToby Isaac   PetscAssertPointer(ctype, 2);
6422dde6fd4SLisandro Dalcin   *ctype = dd->interptype;
6433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6442dde6fd4SLisandro Dalcin }
64547c6ae99SBarry Smith 
6466a119db4SBarry Smith /*@C
647aa219208SBarry Smith   DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
64847c6ae99SBarry Smith   processes neighbors.
64947c6ae99SBarry Smith 
65047c6ae99SBarry Smith   Not Collective
65147c6ae99SBarry Smith 
65247c6ae99SBarry Smith   Input Parameter:
653dce8aebaSBarry Smith . da - the `DMDA` object
65447c6ae99SBarry Smith 
6552fe279fdSBarry Smith   Output Parameter:
656*12b4a537SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. The process itself is in the list
65747c6ae99SBarry Smith 
65847c6ae99SBarry Smith   Level: intermediate
65947c6ae99SBarry Smith 
660dce8aebaSBarry Smith   Notes:
661dce8aebaSBarry Smith   In 2d the array is of length 9, in 3d of length 27
662dce8aebaSBarry Smith 
663dce8aebaSBarry Smith   Not supported in 1d
664dce8aebaSBarry Smith 
665dce8aebaSBarry Smith   Do not free the array, it is freed when the `DMDA` is destroyed.
666dce8aebaSBarry Smith 
667*12b4a537SBarry Smith   Fortran Note:
668*12b4a537SBarry Smith   Pass in an array of the appropriate length to contain the values
669dce8aebaSBarry Smith 
670*12b4a537SBarry Smith .seealso: [](sec_struct), `DMDA`, `DM`
67147c6ae99SBarry Smith @*/
672d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNeighbors(DM da, const PetscMPIInt *ranks[])
673d71ae5a4SJacob Faibussowitsch {
67447c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
6755fd66863SKarl Rupp 
67647c6ae99SBarry Smith   PetscFunctionBegin;
677a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
67847c6ae99SBarry Smith   *ranks = dd->neighbors;
6793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68047c6ae99SBarry Smith }
68147c6ae99SBarry Smith 
68247c6ae99SBarry Smith /*@C
683aa219208SBarry Smith   DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
68447c6ae99SBarry Smith 
68547c6ae99SBarry Smith   Not Collective
68647c6ae99SBarry Smith 
68747c6ae99SBarry Smith   Input Parameter:
688dce8aebaSBarry Smith . da - the `DMDA` object
68947c6ae99SBarry Smith 
690d8d19677SJose E. Roman   Output Parameters:
69147c6ae99SBarry Smith + lx - ownership along x direction (optional)
69247c6ae99SBarry Smith . ly - ownership along y direction (optional)
69347c6ae99SBarry Smith - lz - ownership along z direction (optional)
69447c6ae99SBarry Smith 
69547c6ae99SBarry Smith   Level: intermediate
69647c6ae99SBarry Smith 
697dce8aebaSBarry Smith   Note:
698dce8aebaSBarry Smith   These correspond to the optional final arguments passed to `DMDACreate()`, `DMDACreate2d()`, `DMDACreate3d()`
69947c6ae99SBarry Smith 
70047c6ae99SBarry Smith   In C you should not free these arrays, nor change the values in them. They will only have valid values while the
701dce8aebaSBarry Smith   `DMDA` they came from still exists (has not been destroyed).
70247c6ae99SBarry Smith 
703e3f3e4b6SBarry Smith   These numbers are NOT multiplied by the number of dof per node.
704e3f3e4b6SBarry Smith 
705*12b4a537SBarry Smith   Fortran Note:
706*12b4a537SBarry Smith   Pass in arrays `lx`, `ly`, and `lz` of the appropriate length to hold the values; the sixth, seventh and
707dce8aebaSBarry Smith   eighth arguments from `DMDAGetInfo()`
708dce8aebaSBarry Smith 
709*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetCorners()`, `DMDAGetGhostCorners()`, `DMDACreate()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `VecGetOwnershipRanges()`
71047c6ae99SBarry Smith @*/
711d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetOwnershipRanges(DM da, const PetscInt *lx[], const PetscInt *ly[], const PetscInt *lz[])
712d71ae5a4SJacob Faibussowitsch {
71347c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
71447c6ae99SBarry Smith 
71547c6ae99SBarry Smith   PetscFunctionBegin;
716a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
71747c6ae99SBarry Smith   if (lx) *lx = dd->lx;
71847c6ae99SBarry Smith   if (ly) *ly = dd->ly;
71947c6ae99SBarry Smith   if (lz) *lz = dd->lz;
7203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
72147c6ae99SBarry Smith }
72247c6ae99SBarry Smith 
72347c6ae99SBarry Smith /*@
724dce8aebaSBarry Smith   DMDASetRefinementFactor - Set the ratios that the `DMDA` grid is refined
72547c6ae99SBarry Smith 
72620f4b53cSBarry Smith   Logically Collective
72747c6ae99SBarry Smith 
72847c6ae99SBarry Smith   Input Parameters:
729dce8aebaSBarry Smith + da       - the `DMDA` object
73047c6ae99SBarry Smith . refine_x - ratio of fine grid to coarse in x direction (2 by default)
73147c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default)
73247c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default)
73347c6ae99SBarry Smith 
734dce8aebaSBarry Smith   Options Database Keys:
73548eeb7c8SBarry Smith + -da_refine_x refine_x - refinement ratio in x direction
73648eeb7c8SBarry Smith . -da_refine_y rafine_y - refinement ratio in y direction
73748eeb7c8SBarry Smith . -da_refine_z refine_z - refinement ratio in z direction
738*12b4a537SBarry Smith - -da_refine <n>        - refine the `DMDA` object n times when it is created.
73947c6ae99SBarry Smith 
74047c6ae99SBarry Smith   Level: intermediate
74147c6ae99SBarry Smith 
742dce8aebaSBarry Smith   Note:
743dce8aebaSBarry Smith   Pass `PETSC_IGNORE` to leave a value unchanged
74447c6ae99SBarry Smith 
745*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMRefine()`, `DMDAGetRefinementFactor()`
74647c6ae99SBarry Smith @*/
747d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y, PetscInt refine_z)
748d71ae5a4SJacob Faibussowitsch {
74947c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
75047c6ae99SBarry Smith 
75147c6ae99SBarry Smith   PetscFunctionBegin;
752a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
75347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, refine_x, 2);
75447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, refine_y, 3);
75547c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, refine_z, 4);
75647c6ae99SBarry Smith 
75747c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
75847c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
75947c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
7603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76147c6ae99SBarry Smith }
76247c6ae99SBarry Smith 
76347c6ae99SBarry Smith /*@C
764dce8aebaSBarry Smith   DMDAGetRefinementFactor - Gets the ratios that the `DMDA` grid is refined
76547c6ae99SBarry Smith 
76647c6ae99SBarry Smith   Not Collective
76747c6ae99SBarry Smith 
76847c6ae99SBarry Smith   Input Parameter:
769dce8aebaSBarry Smith . da - the `DMDA` object
77047c6ae99SBarry Smith 
77147c6ae99SBarry Smith   Output Parameters:
77247c6ae99SBarry Smith + refine_x - ratio of fine grid to coarse in x direction (2 by default)
77347c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default)
77447c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default)
77547c6ae99SBarry Smith 
77647c6ae99SBarry Smith   Level: intermediate
77747c6ae99SBarry Smith 
778dce8aebaSBarry Smith   Note:
77920f4b53cSBarry Smith   Pass `NULL` for values you do not need
78047c6ae99SBarry Smith 
781*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMRefine()`, `DMDASetRefinementFactor()`
78247c6ae99SBarry Smith @*/
783d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y, PetscInt *refine_z)
784d71ae5a4SJacob Faibussowitsch {
78547c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
78647c6ae99SBarry Smith 
78747c6ae99SBarry Smith   PetscFunctionBegin;
788a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
78947c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
79047c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
79147c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
7923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79347c6ae99SBarry Smith }
79447c6ae99SBarry Smith 
79547c6ae99SBarry Smith /*@C
796dce8aebaSBarry Smith   DMDASetGetMatrix - Sets the routine used by the `DMDA` to allocate a matrix.
79747c6ae99SBarry Smith 
79820f4b53cSBarry Smith   Logically Collective; No Fortran Support
79947c6ae99SBarry Smith 
80047c6ae99SBarry Smith   Input Parameters:
801dce8aebaSBarry Smith + da - the `DMDA` object
802*12b4a537SBarry Smith - f  - the function that allocates the matrix for that specific `DMDA`
803*12b4a537SBarry Smith 
804*12b4a537SBarry Smith   Calling sequence of `f`:
805*12b4a537SBarry Smith + da - the `DMDA` object
806*12b4a537SBarry Smith - A  - the created matrix
80747c6ae99SBarry Smith 
80847c6ae99SBarry Smith   Level: developer
80947c6ae99SBarry Smith 
810*12b4a537SBarry Smith   Notes:
811*12b4a537SBarry Smith   If the function is not provided a default function is used that uses the `DMDAStencilType`, `DMBoundaryType`, and value of `DMDASetStencilWidth()`
812*12b4a537SBarry Smith   to construct the matrix.
81347c6ae99SBarry Smith 
814*12b4a537SBarry Smith   See `DMDASetBlockFills()` that provides a simple way to provide the nonzero structure for
815*12b4a537SBarry Smith   the diagonal and off-diagonal blocks of the matrix without providing a custom function
816*12b4a537SBarry Smith 
817*12b4a537SBarry Smith   Developer Note:
818*12b4a537SBarry Smith   This should be called `DMDASetCreateMatrix()`
819*12b4a537SBarry Smith 
820*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateMatrix()`, `DMDASetBlockFills()`
82147c6ae99SBarry Smith @*/
822*12b4a537SBarry Smith PetscErrorCode DMDASetGetMatrix(DM da, PetscErrorCode (*f)(DM da, Mat *A))
823d71ae5a4SJacob Faibussowitsch {
82447c6ae99SBarry Smith   PetscFunctionBegin;
825a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
82625296bd5SBarry Smith   da->ops->creatematrix = f;
8273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
82847c6ae99SBarry Smith }
82947c6ae99SBarry Smith 
83038fb4e8eSJunchao Zhang /*@
831dce8aebaSBarry Smith   DMDAMapMatStencilToGlobal - Map a list of `MatStencil` on a grid to global indices.
83238fb4e8eSJunchao Zhang 
83338fb4e8eSJunchao Zhang   Not Collective
83438fb4e8eSJunchao Zhang 
83538fb4e8eSJunchao Zhang   Input Parameters:
836dce8aebaSBarry Smith + da   - the `DMDA` object
837*12b4a537SBarry Smith . m    - number of `MatStencil` to map
83838fb4e8eSJunchao Zhang - idxm - grid points (and component number when dof > 1)
83938fb4e8eSJunchao Zhang 
8401179163eSBarry Smith   Output Parameter:
8411179163eSBarry Smith . gidxm - global row indices
8421179163eSBarry Smith 
8431179163eSBarry Smith   Level: intermediate
84438fb4e8eSJunchao Zhang 
845*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `MatStencil`
84638fb4e8eSJunchao Zhang @*/
847d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAMapMatStencilToGlobal(DM da, PetscInt m, const MatStencil idxm[], PetscInt gidxm[])
848d71ae5a4SJacob Faibussowitsch {
84938fb4e8eSJunchao Zhang   const DM_DA           *dd  = (const DM_DA *)da->data;
85038fb4e8eSJunchao Zhang   const PetscInt        *dxm = (const PetscInt *)idxm;
85138fb4e8eSJunchao Zhang   PetscInt               i, j, sdim, tmp, dim;
85238fb4e8eSJunchao Zhang   PetscInt               dims[4], starts[4], dims2[3], starts2[3], dof = dd->w;
85338fb4e8eSJunchao Zhang   ISLocalToGlobalMapping ltog;
85438fb4e8eSJunchao Zhang 
85538fb4e8eSJunchao Zhang   PetscFunctionBegin;
8563ba16761SJacob Faibussowitsch   if (m <= 0) PetscFunctionReturn(PETSC_SUCCESS);
85738fb4e8eSJunchao Zhang 
8582f27f4d1SJunchao Zhang   /* Code adapted from DMDAGetGhostCorners() */
85938fb4e8eSJunchao Zhang   starts2[0] = dd->Xs / dof + dd->xo;
86038fb4e8eSJunchao Zhang   starts2[1] = dd->Ys + dd->yo;
86138fb4e8eSJunchao Zhang   starts2[2] = dd->Zs + dd->zo;
86238fb4e8eSJunchao Zhang   dims2[0]   = (dd->Xe - dd->Xs) / dof;
86338fb4e8eSJunchao Zhang   dims2[1]   = (dd->Ye - dd->Ys);
86438fb4e8eSJunchao Zhang   dims2[2]   = (dd->Ze - dd->Zs);
86538fb4e8eSJunchao Zhang 
8662f27f4d1SJunchao Zhang   /* As if we do MatSetStencil() to get dims[]/starts[] of mat->stencil */
8672f27f4d1SJunchao Zhang   dim  = da->dim;                   /* DA dim: 1 to 3 */
8682f27f4d1SJunchao Zhang   sdim = dim + (dof > 1 ? 1 : 0);   /* Dimensions in MatStencil's (k,j,i,c) view */
8692f27f4d1SJunchao Zhang   for (i = 0; i < dim; i++) {       /* Reverse the order and also skip the unused dimensions */
8702f27f4d1SJunchao 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 */
87138fb4e8eSJunchao Zhang     starts[i] = starts2[dim - i - 1];
87238fb4e8eSJunchao Zhang   }
8732f27f4d1SJunchao Zhang   starts[dim] = 0; /* Append the extra dim for dof (won't be used below if dof=1) */
87438fb4e8eSJunchao Zhang   dims[dim]   = dof;
87538fb4e8eSJunchao Zhang 
87638fb4e8eSJunchao Zhang   /* Map stencils to local indices (code adapted from MatSetValuesStencil()) */
87738fb4e8eSJunchao Zhang   for (i = 0; i < m; i++) {
8782f27f4d1SJunchao Zhang     dxm += 3 - dim; /* Input is {k,j,i,c}; move the pointer to the first used index, e.g., j in 2D */
8792f27f4d1SJunchao Zhang     tmp = 0;
8802f27f4d1SJunchao Zhang     for (j = 0; j < sdim; j++) {                                                      /* Iter over, ex. j,i or j,i,c in 2D */
8812f27f4d1SJunchao 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 */
8822f27f4d1SJunchao Zhang       else tmp = tmp * dims[j] + (dxm[j] - starts[j]);
88338fb4e8eSJunchao Zhang     }
88438fb4e8eSJunchao Zhang     gidxm[i] = tmp;
8852f27f4d1SJunchao Zhang     /* Move to the next MatStencil point */
8862f27f4d1SJunchao Zhang     if (dof > 1) dxm += sdim; /* c is already counted in sdim */
8872f27f4d1SJunchao Zhang     else dxm += sdim + 1;     /* skip the unused c */
88838fb4e8eSJunchao Zhang   }
88938fb4e8eSJunchao Zhang 
89038fb4e8eSJunchao Zhang   /* Map local indices to global indices */
89138fb4e8eSJunchao Zhang   PetscCall(DMGetLocalToGlobalMapping(da, &ltog));
89238fb4e8eSJunchao Zhang   PetscCall(ISLocalToGlobalMappingApply(ltog, m, gidxm, gidxm));
8933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
89438fb4e8eSJunchao Zhang }
89538fb4e8eSJunchao Zhang 
89647c6ae99SBarry Smith /*
89747c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
89847c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
89947c6ae99SBarry Smith 
90047c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
90147c6ae99SBarry Smith */
902d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDARefineOwnershipRanges(DM da, PetscBool periodic, PetscInt stencil_width, PetscInt ratio, PetscInt m, const PetscInt lc[], PetscInt lf[])
903d71ae5a4SJacob Faibussowitsch {
90447c6ae99SBarry Smith   PetscInt i, totalc = 0, remaining, startc = 0, startf = 0;
90547c6ae99SBarry Smith 
90647c6ae99SBarry Smith   PetscFunctionBegin;
90763a3b9bcSJacob Faibussowitsch   PetscCheck(ratio >= 1, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Requested refinement ratio %" PetscInt_FMT " must be at least 1", ratio);
90847c6ae99SBarry Smith   if (ratio == 1) {
9099566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(lf, lc, m));
9103ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
91147c6ae99SBarry Smith   }
91247c6ae99SBarry Smith   for (i = 0; i < m; i++) totalc += lc[i];
91347c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
91447c6ae99SBarry Smith   for (i = 0; i < m; i++) {
91547c6ae99SBarry Smith     PetscInt want = remaining / (m - i) + !!(remaining % (m - i));
91647c6ae99SBarry Smith     if (i == m - 1) lf[i] = want;
91747c6ae99SBarry Smith     else {
9187aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
9197aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
9207aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
9217aca7175SJed Brown       while ((startf + want) / ratio < nextc - stencil_width) want++;
9227aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
9237aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
9247aca7175SJed Brown       while ((startf + want - 1 + ratio - 1) / ratio > nextc - 1 + stencil_width) want--;
9257aca7175SJed Brown       /* Make sure all constraints are satisfied */
9269371c9d4SSatish Balay       if (want < 0 || want > remaining || ((startf + want) / ratio < nextc - stencil_width) || ((startf + want - 1 + ratio - 1) / ratio > nextc - 1 + stencil_width))
9279371c9d4SSatish Balay         SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Could not find a compatible refined ownership range");
92847c6ae99SBarry Smith     }
92947c6ae99SBarry Smith     lf[i] = want;
93047c6ae99SBarry Smith     startc += lc[i];
93147c6ae99SBarry Smith     startf += lf[i];
93247c6ae99SBarry Smith     remaining -= lf[i];
93347c6ae99SBarry Smith   }
9343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
93547c6ae99SBarry Smith }
93647c6ae99SBarry Smith 
9372be375d4SJed Brown /*
9382be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
9392be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
9402be375d4SJed Brown 
9412be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
9422be375d4SJed Brown */
943d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDACoarsenOwnershipRanges(DM da, PetscBool periodic, PetscInt stencil_width, PetscInt ratio, PetscInt m, const PetscInt lf[], PetscInt lc[])
944d71ae5a4SJacob Faibussowitsch {
9452be375d4SJed Brown   PetscInt i, totalf, remaining, startc, startf;
9462be375d4SJed Brown 
9472be375d4SJed Brown   PetscFunctionBegin;
94863a3b9bcSJacob Faibussowitsch   PetscCheck(ratio >= 1, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Requested refinement ratio %" PetscInt_FMT " must be at least 1", ratio);
9492be375d4SJed Brown   if (ratio == 1) {
9509566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(lc, lf, m));
9513ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
9522be375d4SJed Brown   }
9532be375d4SJed Brown   for (i = 0, totalf = 0; i < m; i++) totalf += lf[i];
9542be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
9552be375d4SJed Brown   for (i = 0, startc = 0, startf = 0; i < m; i++) {
9562be375d4SJed Brown     PetscInt want = remaining / (m - i) + !!(remaining % (m - i));
9572be375d4SJed Brown     if (i == m - 1) lc[i] = want;
9582be375d4SJed Brown     else {
9592be375d4SJed Brown       const PetscInt nextf = startf + lf[i];
9602be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
9612be375d4SJed Brown        * node is within one stencil width. */
9622be375d4SJed Brown       while (nextf / ratio < startc + want - stencil_width) want--;
9632be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
9642be375d4SJed Brown        * fine node is within one stencil width. */
9652be375d4SJed Brown       while ((nextf - 1 + ratio - 1) / ratio > startc + want - 1 + stencil_width) want++;
9669371c9d4SSatish Balay       if (want < 0 || want > remaining || (nextf / ratio < startc + want - stencil_width) || ((nextf - 1 + ratio - 1) / ratio > startc + want - 1 + stencil_width))
9679371c9d4SSatish Balay         SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Could not find a compatible coarsened ownership range");
9682be375d4SJed Brown     }
9692be375d4SJed Brown     lc[i] = want;
9702be375d4SJed Brown     startc += lc[i];
9712be375d4SJed Brown     startf += lf[i];
9722be375d4SJed Brown     remaining -= lc[i];
9732be375d4SJed Brown   }
9743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9752be375d4SJed Brown }
9762be375d4SJed Brown 
977d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_DA(DM da, MPI_Comm comm, DM *daref)
978d71ae5a4SJacob Faibussowitsch {
979c73cfb54SMatthew G. Knepley   PetscInt M, N, P, i, dim;
9806858538eSMatthew G. Knepley   Vec      coordsc, coordsf;
9819a42bb27SBarry Smith   DM       da2;
98247c6ae99SBarry Smith   DM_DA   *dd = (DM_DA *)da->data, *dd2;
98347c6ae99SBarry Smith 
98447c6ae99SBarry Smith   PetscFunctionBegin;
985a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
9864f572ea9SToby Isaac   PetscAssertPointer(daref, 3);
98747c6ae99SBarry Smith 
9889566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(da, &dim));
989bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
99047c6ae99SBarry Smith     M = dd->refine_x * dd->M;
99147c6ae99SBarry Smith   } else {
99247c6ae99SBarry Smith     M = 1 + dd->refine_x * (dd->M - 1);
99347c6ae99SBarry Smith   }
994bff4a2f0SMatthew G. Knepley   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
995c73cfb54SMatthew G. Knepley     if (dim > 1) {
99647c6ae99SBarry Smith       N = dd->refine_y * dd->N;
99747c6ae99SBarry Smith     } else {
9981860e6e9SBarry Smith       N = 1;
9991860e6e9SBarry Smith     }
10001860e6e9SBarry Smith   } else {
100147c6ae99SBarry Smith     N = 1 + dd->refine_y * (dd->N - 1);
100247c6ae99SBarry Smith   }
1003bff4a2f0SMatthew G. Knepley   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1004c73cfb54SMatthew G. Knepley     if (dim > 2) {
100547c6ae99SBarry Smith       P = dd->refine_z * dd->P;
100647c6ae99SBarry Smith     } else {
10071860e6e9SBarry Smith       P = 1;
10081860e6e9SBarry Smith     }
10091860e6e9SBarry Smith   } else {
101047c6ae99SBarry Smith     P = 1 + dd->refine_z * (dd->P - 1);
101147c6ae99SBarry Smith   }
10129566063dSJacob Faibussowitsch   PetscCall(DMDACreate(PetscObjectComm((PetscObject)da), &da2));
10139566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(da2, ((PetscObject)da)->prefix));
10149566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(da2, dim));
10159566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(da2, M, N, P));
10169566063dSJacob Faibussowitsch   PetscCall(DMDASetNumProcs(da2, dd->m, dd->n, dd->p));
10179566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(da2, dd->bx, dd->by, dd->bz));
10189566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(da2, dd->w));
10199566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilType(da2, dd->stencil_type));
10209566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(da2, dd->s));
1021c73cfb54SMatthew G. Knepley   if (dim == 3) {
102247c6ae99SBarry Smith     PetscInt *lx, *ly, *lz;
10239566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(dd->m, &lx, dd->n, &ly, dd->p, &lz));
10249566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx));
10259566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_y, dd->n, dd->ly, ly));
10269566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_z, dd->p, dd->lz, lz));
10279566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(da2, lx, ly, lz));
10289566063dSJacob Faibussowitsch     PetscCall(PetscFree3(lx, ly, lz));
1029c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
103047c6ae99SBarry Smith     PetscInt *lx, *ly;
10319566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(dd->m, &lx, dd->n, &ly));
10329566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx));
10339566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_y, dd->n, dd->ly, ly));
10349566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(da2, lx, ly, NULL));
10359566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lx, ly));
1036c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
103747c6ae99SBarry Smith     PetscInt *lx;
10389566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->m, &lx));
10399566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx));
10409566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(da2, lx, NULL, NULL));
10419566063dSJacob Faibussowitsch     PetscCall(PetscFree(lx));
104247c6ae99SBarry Smith   }
104347c6ae99SBarry Smith   dd2 = (DM_DA *)da2->data;
104447c6ae99SBarry Smith 
104547c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
104625296bd5SBarry Smith   da2->ops->creatematrix = da->ops->creatematrix;
104725296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
104847c6ae99SBarry Smith   da2->ops->getcoloring = da->ops->getcoloring;
104947c6ae99SBarry Smith   dd2->interptype       = dd->interptype;
105047c6ae99SBarry Smith 
105147c6ae99SBarry Smith   /* copy fill information if given */
105247c6ae99SBarry Smith   if (dd->dfill) {
10539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->dfill[dd->w] + dd->w + 1, &dd2->dfill));
10549566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->dfill, dd->dfill, dd->dfill[dd->w] + dd->w + 1));
105547c6ae99SBarry Smith   }
105647c6ae99SBarry Smith   if (dd->ofill) {
10579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->ofill[dd->w] + dd->w + 1, &dd2->ofill));
10589566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->ofill, dd->ofill, dd->ofill[dd->w] + dd->w + 1));
105947c6ae99SBarry Smith   }
106047c6ae99SBarry Smith   /* copy the refine information */
1061397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
1062397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
1063397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
106447c6ae99SBarry Smith 
1065897f7067SBarry Smith   if (dd->refine_z_hier) {
1066ad540459SPierre Jolivet     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];
1067ad540459SPierre Jolivet     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];
1068897f7067SBarry Smith     dd2->refine_z_hier_n = dd->refine_z_hier_n;
10699566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_z_hier_n, &dd2->refine_z_hier));
10709566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_z_hier, dd->refine_z_hier, dd2->refine_z_hier_n));
1071897f7067SBarry Smith   }
1072897f7067SBarry Smith   if (dd->refine_y_hier) {
1073ad540459SPierre Jolivet     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];
1074ad540459SPierre Jolivet     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];
1075897f7067SBarry Smith     dd2->refine_y_hier_n = dd->refine_y_hier_n;
10769566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_y_hier_n, &dd2->refine_y_hier));
10779566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_y_hier, dd->refine_y_hier, dd2->refine_y_hier_n));
1078897f7067SBarry Smith   }
1079897f7067SBarry Smith   if (dd->refine_x_hier) {
1080ad540459SPierre Jolivet     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];
1081ad540459SPierre Jolivet     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];
1082897f7067SBarry Smith     dd2->refine_x_hier_n = dd->refine_x_hier_n;
10839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_x_hier_n, &dd2->refine_x_hier));
10849566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_x_hier, dd->refine_x_hier, dd2->refine_x_hier_n));
1085897f7067SBarry Smith   }
1086897f7067SBarry Smith 
108747c6ae99SBarry Smith   /* copy vector type information */
10889566063dSJacob Faibussowitsch   PetscCall(DMSetVecType(da2, da->vectype));
1089ddcf8b74SDave May 
1090efd51863SBarry Smith   dd2->lf = dd->lf;
1091efd51863SBarry Smith   dd2->lj = dd->lj;
1092efd51863SBarry Smith 
10936e87535bSJed Brown   da2->leveldown = da->leveldown;
10946e87535bSJed Brown   da2->levelup   = da->levelup + 1;
10958865f1eaSKarl Rupp 
10969566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da2));
10976e87535bSJed Brown 
1098ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
10996858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(da, &coordsc));
11006858538eSMatthew G. Knepley   if (coordsc) {
1101ddcf8b74SDave May     DM  cdaf, cdac;
1102ddcf8b74SDave May     Mat II;
1103ddcf8b74SDave May 
11049566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(da, &cdac));
11059566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(da2, &cdaf));
1106b61d3410SDave May     /* force creation of the coordinate vector */
11079566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(da2, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
11089566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(da2, &coordsf));
11099566063dSJacob Faibussowitsch     PetscCall(DMCreateInterpolation(cdac, cdaf, &II, NULL));
11109566063dSJacob Faibussowitsch     PetscCall(MatInterpolate(II, coordsc, coordsf));
11119566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&II));
1112ddcf8b74SDave May   }
1113397b6216SJed Brown 
1114f3141302SJed Brown   for (i = 0; i < da->bs; i++) {
1115f3141302SJed Brown     const char *fieldname;
11169566063dSJacob Faibussowitsch     PetscCall(DMDAGetFieldName(da, i, &fieldname));
11179566063dSJacob Faibussowitsch     PetscCall(DMDASetFieldName(da2, i, fieldname));
1118f3141302SJed Brown   }
1119397b6216SJed Brown 
112047c6ae99SBarry Smith   *daref = da2;
11213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
112247c6ae99SBarry Smith }
112347c6ae99SBarry Smith 
1124d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsen_DA(DM dmf, MPI_Comm comm, DM *dmc)
1125d71ae5a4SJacob Faibussowitsch {
1126c73cfb54SMatthew G. Knepley   PetscInt M, N, P, i, dim;
11276858538eSMatthew G. Knepley   Vec      coordsc, coordsf;
1128a5bc1bf3SBarry Smith   DM       dmc2;
1129a5bc1bf3SBarry Smith   DM_DA   *dd = (DM_DA *)dmf->data, *dd2;
113047c6ae99SBarry Smith 
113147c6ae99SBarry Smith   PetscFunctionBegin;
1132a5bc1bf3SBarry Smith   PetscValidHeaderSpecificType(dmf, DM_CLASSID, 1, DMDA);
11334f572ea9SToby Isaac   PetscAssertPointer(dmc, 3);
113447c6ae99SBarry Smith 
11359566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dmf, &dim));
1136bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1137397b6216SJed Brown     M = dd->M / dd->coarsen_x;
113847c6ae99SBarry Smith   } else {
1139397b6216SJed Brown     M = 1 + (dd->M - 1) / dd->coarsen_x;
114047c6ae99SBarry Smith   }
1141bff4a2f0SMatthew G. Knepley   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1142c73cfb54SMatthew G. Knepley     if (dim > 1) {
1143397b6216SJed Brown       N = dd->N / dd->coarsen_y;
114447c6ae99SBarry Smith     } else {
11451860e6e9SBarry Smith       N = 1;
11461860e6e9SBarry Smith     }
11471860e6e9SBarry Smith   } else {
1148397b6216SJed Brown     N = 1 + (dd->N - 1) / dd->coarsen_y;
114947c6ae99SBarry Smith   }
1150bff4a2f0SMatthew G. Knepley   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1151c73cfb54SMatthew G. Knepley     if (dim > 2) {
1152397b6216SJed Brown       P = dd->P / dd->coarsen_z;
115347c6ae99SBarry Smith     } else {
11541860e6e9SBarry Smith       P = 1;
11551860e6e9SBarry Smith     }
11561860e6e9SBarry Smith   } else {
1157397b6216SJed Brown     P = 1 + (dd->P - 1) / dd->coarsen_z;
115847c6ae99SBarry Smith   }
11599566063dSJacob Faibussowitsch   PetscCall(DMDACreate(PetscObjectComm((PetscObject)dmf), &dmc2));
11609566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(dmc2, ((PetscObject)dmf)->prefix));
11619566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dmc2, dim));
11629566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(dmc2, M, N, P));
11639566063dSJacob Faibussowitsch   PetscCall(DMDASetNumProcs(dmc2, dd->m, dd->n, dd->p));
11649566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(dmc2, dd->bx, dd->by, dd->bz));
11659566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(dmc2, dd->w));
11669566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilType(dmc2, dd->stencil_type));
11679566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(dmc2, dd->s));
1168c73cfb54SMatthew G. Knepley   if (dim == 3) {
11692be375d4SJed Brown     PetscInt *lx, *ly, *lz;
11709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(dd->m, &lx, dd->n, &ly, dd->p, &lz));
11719566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx));
11729566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_y, dd->n, dd->ly, ly));
11739566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_z, dd->p, dd->lz, lz));
11749566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(dmc2, lx, ly, lz));
11759566063dSJacob Faibussowitsch     PetscCall(PetscFree3(lx, ly, lz));
1176c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
11772be375d4SJed Brown     PetscInt *lx, *ly;
11789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(dd->m, &lx, dd->n, &ly));
11799566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx));
11809566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_y, dd->n, dd->ly, ly));
11819566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(dmc2, lx, ly, NULL));
11829566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lx, ly));
1183c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
11842be375d4SJed Brown     PetscInt *lx;
11859566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->m, &lx));
11869566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx));
11879566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(dmc2, lx, NULL, NULL));
11889566063dSJacob Faibussowitsch     PetscCall(PetscFree(lx));
118947c6ae99SBarry Smith   }
1190a5bc1bf3SBarry Smith   dd2 = (DM_DA *)dmc2->data;
119147c6ae99SBarry Smith 
11924dcab191SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
1193a5bc1bf3SBarry Smith   /* dmc2->ops->createinterpolation = dmf->ops->createinterpolation; copying this one causes trouble for DMSetVI */
1194a5bc1bf3SBarry Smith   dmc2->ops->creatematrix = dmf->ops->creatematrix;
1195a5bc1bf3SBarry Smith   dmc2->ops->getcoloring  = dmf->ops->getcoloring;
119647c6ae99SBarry Smith   dd2->interptype         = dd->interptype;
119747c6ae99SBarry Smith 
119847c6ae99SBarry Smith   /* copy fill information if given */
119947c6ae99SBarry Smith   if (dd->dfill) {
12009566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->dfill[dd->w] + dd->w + 1, &dd2->dfill));
12019566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->dfill, dd->dfill, dd->dfill[dd->w] + dd->w + 1));
120247c6ae99SBarry Smith   }
120347c6ae99SBarry Smith   if (dd->ofill) {
12049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->ofill[dd->w] + dd->w + 1, &dd2->ofill));
12059566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->ofill, dd->ofill, dd->ofill[dd->w] + dd->w + 1));
120647c6ae99SBarry Smith   }
120747c6ae99SBarry Smith   /* copy the refine information */
1208397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
1209397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
1210397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
121147c6ae99SBarry Smith 
1212897f7067SBarry Smith   if (dd->refine_z_hier) {
1213ad540459SPierre Jolivet     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];
1214ad540459SPierre Jolivet     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];
1215897f7067SBarry Smith     dd2->refine_z_hier_n = dd->refine_z_hier_n;
12169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_z_hier_n, &dd2->refine_z_hier));
12179566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_z_hier, dd->refine_z_hier, dd2->refine_z_hier_n));
1218897f7067SBarry Smith   }
1219897f7067SBarry Smith   if (dd->refine_y_hier) {
1220ad540459SPierre Jolivet     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];
1221ad540459SPierre Jolivet     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];
1222897f7067SBarry Smith     dd2->refine_y_hier_n = dd->refine_y_hier_n;
12239566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_y_hier_n, &dd2->refine_y_hier));
12249566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_y_hier, dd->refine_y_hier, dd2->refine_y_hier_n));
1225897f7067SBarry Smith   }
1226897f7067SBarry Smith   if (dd->refine_x_hier) {
1227ad540459SPierre Jolivet     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];
1228ad540459SPierre Jolivet     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];
1229897f7067SBarry Smith     dd2->refine_x_hier_n = dd->refine_x_hier_n;
12309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_x_hier_n, &dd2->refine_x_hier));
12319566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_x_hier, dd->refine_x_hier, dd2->refine_x_hier_n));
1232897f7067SBarry Smith   }
1233897f7067SBarry Smith 
123447c6ae99SBarry Smith   /* copy vector type information */
12359566063dSJacob Faibussowitsch   PetscCall(DMSetVecType(dmc2, dmf->vectype));
123647c6ae99SBarry Smith 
1237644e2e5bSBarry Smith   dd2->lf = dd->lf;
1238644e2e5bSBarry Smith   dd2->lj = dd->lj;
1239644e2e5bSBarry Smith 
1240a5bc1bf3SBarry Smith   dmc2->leveldown = dmf->leveldown + 1;
1241a5bc1bf3SBarry Smith   dmc2->levelup   = dmf->levelup;
12428865f1eaSKarl Rupp 
12439566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmc2));
12446e87535bSJed Brown 
1245ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
12466858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dmf, &coordsf));
12476858538eSMatthew G. Knepley   if (coordsf) {
1248ddcf8b74SDave May     DM         cdaf, cdac;
12496dbf9973SLawrence Mitchell     Mat        inject;
12506dbf9973SLawrence Mitchell     VecScatter vscat;
1251ddcf8b74SDave May 
12529566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmf, &cdaf));
12539566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmc2, &cdac));
1254b61d3410SDave May     /* force creation of the coordinate vector */
12559566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(dmc2, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
12569566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmc2, &coordsc));
1257ddcf8b74SDave May 
12589566063dSJacob Faibussowitsch     PetscCall(DMCreateInjection(cdac, cdaf, &inject));
12599566063dSJacob Faibussowitsch     PetscCall(MatScatterGetVecScatter(inject, &vscat));
12609566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(vscat, coordsf, coordsc, INSERT_VALUES, SCATTER_FORWARD));
12619566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(vscat, coordsf, coordsc, INSERT_VALUES, SCATTER_FORWARD));
12629566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&inject));
1263ddcf8b74SDave May   }
1264f98405f7SJed Brown 
1265a5bc1bf3SBarry Smith   for (i = 0; i < dmf->bs; i++) {
1266f98405f7SJed Brown     const char *fieldname;
12679566063dSJacob Faibussowitsch     PetscCall(DMDAGetFieldName(dmf, i, &fieldname));
12689566063dSJacob Faibussowitsch     PetscCall(DMDASetFieldName(dmc2, i, fieldname));
1269f98405f7SJed Brown   }
1270f98405f7SJed Brown 
1271a5bc1bf3SBarry Smith   *dmc = dmc2;
12723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
127347c6ae99SBarry Smith }
127447c6ae99SBarry Smith 
1275d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefineHierarchy_DA(DM da, PetscInt nlevels, DM daf[])
1276d71ae5a4SJacob Faibussowitsch {
127747c6ae99SBarry Smith   PetscInt i, n, *refx, *refy, *refz;
127847c6ae99SBarry Smith 
127947c6ae99SBarry Smith   PetscFunctionBegin;
128047c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
12817a8be351SBarry Smith   PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative");
12823ba16761SJacob Faibussowitsch   if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS);
12834f572ea9SToby Isaac   PetscAssertPointer(daf, 3);
128447c6ae99SBarry Smith 
1285aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
12869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(nlevels, &refx, nlevels, &refy, nlevels, &refz));
128748a46eb9SPierre Jolivet   for (i = 0; i < nlevels; i++) PetscCall(DMDAGetRefinementFactor(da, &refx[i], &refy[i], &refz[i]));
128847c6ae99SBarry Smith   n = nlevels;
12899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_x", refx, &n, NULL));
129047c6ae99SBarry Smith   n = nlevels;
12919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_y", refy, &n, NULL));
129247c6ae99SBarry Smith   n = nlevels;
12939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_z", refz, &n, NULL));
129447c6ae99SBarry Smith 
12959566063dSJacob Faibussowitsch   PetscCall(DMDASetRefinementFactor(da, refx[0], refy[0], refz[0]));
12969566063dSJacob Faibussowitsch   PetscCall(DMRefine(da, PetscObjectComm((PetscObject)da), &daf[0]));
129747c6ae99SBarry Smith   for (i = 1; i < nlevels; i++) {
12989566063dSJacob Faibussowitsch     PetscCall(DMDASetRefinementFactor(daf[i - 1], refx[i], refy[i], refz[i]));
12999566063dSJacob Faibussowitsch     PetscCall(DMRefine(daf[i - 1], PetscObjectComm((PetscObject)da), &daf[i]));
130047c6ae99SBarry Smith   }
13019566063dSJacob Faibussowitsch   PetscCall(PetscFree3(refx, refy, refz));
13023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
130347c6ae99SBarry Smith }
130447c6ae99SBarry Smith 
1305d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsenHierarchy_DA(DM da, PetscInt nlevels, DM dac[])
1306d71ae5a4SJacob Faibussowitsch {
130747c6ae99SBarry Smith   PetscInt i;
130847c6ae99SBarry Smith 
130947c6ae99SBarry Smith   PetscFunctionBegin;
131047c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
13117a8be351SBarry Smith   PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative");
13123ba16761SJacob Faibussowitsch   if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS);
13134f572ea9SToby Isaac   PetscAssertPointer(dac, 3);
13149566063dSJacob Faibussowitsch   PetscCall(DMCoarsen(da, PetscObjectComm((PetscObject)da), &dac[0]));
131548a46eb9SPierre Jolivet   for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dac[i - 1], PetscObjectComm((PetscObject)da), &dac[i]));
13163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131747c6ae99SBarry Smith }
131862197512SBarry Smith 
131966976f2fSJacob Faibussowitsch static PetscErrorCode DMDASetGLLCoordinates_1d(DM dm, PetscInt n, PetscReal *nodes)
1320d71ae5a4SJacob Faibussowitsch {
13218272889dSSatish Balay   PetscInt     i, j, xs, xn, q;
132262197512SBarry Smith   PetscScalar *xx;
132362197512SBarry Smith   PetscReal    h;
132462197512SBarry Smith   Vec          x;
132562197512SBarry Smith   DM_DA       *da = (DM_DA *)dm->data;
132662197512SBarry Smith 
132762197512SBarry Smith   PetscFunctionBegin;
132862197512SBarry Smith   if (da->bx != DM_BOUNDARY_PERIODIC) {
13299566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(dm, NULL, &q, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
133062197512SBarry Smith     q = (q - 1) / (n - 1); /* number of spectral elements */
133162197512SBarry Smith     h = 2.0 / q;
13329566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xn, NULL, NULL));
133362197512SBarry Smith     xs = xs / (n - 1);
133462197512SBarry Smith     xn = xn / (n - 1);
13359566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(dm, -1., 1., 0., 0., 0., 0.));
13369566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dm, &x));
13379566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, x, &xx));
133862197512SBarry Smith 
133962197512SBarry Smith     /* loop over local spectral elements */
134062197512SBarry Smith     for (j = xs; j < xs + xn; j++) {
134162197512SBarry Smith       /*
134262197512SBarry Smith        Except for the first process, each process starts on the second GLL point of the first element on that process
134362197512SBarry Smith        */
1344ad540459SPierre Jolivet       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.;
134562197512SBarry Smith     }
13469566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, x, &xx));
134762197512SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Not yet implemented for periodic");
13483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
134962197512SBarry Smith }
135062197512SBarry Smith 
13511fd49c25SBarry Smith /*@
135262197512SBarry Smith 
135362197512SBarry 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
135462197512SBarry Smith 
135520f4b53cSBarry Smith   Collective
135662197512SBarry Smith 
135762197512SBarry Smith   Input Parameters:
1358dce8aebaSBarry Smith + da    - the `DMDA` object
13592fe279fdSBarry Smith . n     - the number of GLL nodes
13608272889dSSatish Balay - nodes - the GLL nodes
136162197512SBarry Smith 
1362edc382c3SSatish Balay   Level: advanced
1363edc382c3SSatish Balay 
1364dce8aebaSBarry Smith   Note:
1365dce8aebaSBarry Smith   The parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points
1366*12b4a537SBarry Smith   on each process much be divisible by the number of GLL elements needed per process. This depends on whether the `DMDA` is
1367dce8aebaSBarry Smith   periodic or not.
1368dce8aebaSBarry Smith 
1369*12b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `PetscDTGaussLobattoLegendreQuadrature()`, `DMGetCoordinates()`
137062197512SBarry Smith @*/
1371d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetGLLCoordinates(DM da, PetscInt n, PetscReal *nodes)
1372d71ae5a4SJacob Faibussowitsch {
137362197512SBarry Smith   PetscFunctionBegin;
137462197512SBarry Smith   if (da->dim == 1) {
13759566063dSJacob Faibussowitsch     PetscCall(DMDASetGLLCoordinates_1d(da, n, nodes));
137662197512SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Not yet implemented for 2 or 3d");
13773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
137862197512SBarry Smith }
13797c3cd84eSPatrick Sanan 
138066976f2fSJacob Faibussowitsch PetscErrorCode DMGetCompatibility_DA(DM da1, DM dm2, PetscBool *compatible, PetscBool *set)
1381d71ae5a4SJacob Faibussowitsch {
13827c3cd84eSPatrick Sanan   DM_DA    *dd1 = (DM_DA *)da1->data, *dd2;
13837c3cd84eSPatrick Sanan   DM        da2;
13847c3cd84eSPatrick Sanan   DMType    dmtype2;
13857c3cd84eSPatrick Sanan   PetscBool isda, compatibleLocal;
13867c3cd84eSPatrick Sanan   PetscInt  i;
13877c3cd84eSPatrick Sanan 
13887c3cd84eSPatrick Sanan   PetscFunctionBegin;
13897a8be351SBarry Smith   PetscCheck(da1->setupcalled, PetscObjectComm((PetscObject)da1), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() must be called on first DM before DMGetCompatibility()");
13909566063dSJacob Faibussowitsch   PetscCall(DMGetType(dm2, &dmtype2));
13919566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(dmtype2, DMDA, &isda));
13927c3cd84eSPatrick Sanan   if (isda) {
13937c3cd84eSPatrick Sanan     da2 = dm2;
13947c3cd84eSPatrick Sanan     dd2 = (DM_DA *)da2->data;
13957a8be351SBarry Smith     PetscCheck(da2->setupcalled, PetscObjectComm((PetscObject)da2), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() must be called on second DM before DMGetCompatibility()");
1396110623a0SKarl Rupp     compatibleLocal = (PetscBool)(da1->dim == da2->dim);
1397c790a739SKarl Rupp     if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->s == dd2->s)); /* Stencil width */
1398110623a0SKarl Rupp     /*                                                                           Global size              ranks               Boundary type */
1399c790a739SKarl Rupp     if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->M == dd2->M) && (dd1->m == dd2->m) && (dd1->bx == dd2->bx));
1400c790a739SKarl Rupp     if (compatibleLocal && da1->dim > 1) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->N == dd2->N) && (dd1->n == dd2->n) && (dd1->by == dd2->by));
1401c790a739SKarl Rupp     if (compatibleLocal && da1->dim > 2) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->P == dd2->P) && (dd1->p == dd2->p) && (dd1->bz == dd2->bz));
14027c3cd84eSPatrick Sanan     if (compatibleLocal) {
14039371c9d4SSatish Balay       for (i = 0; i < dd1->m; ++i) { compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lx[i] == dd2->lx[i])); /* Local size     */ }
14047c3cd84eSPatrick Sanan     }
14057c3cd84eSPatrick Sanan     if (compatibleLocal && da1->dim > 1) {
1406ad540459SPierre Jolivet       for (i = 0; i < dd1->n; ++i) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->ly[i] == dd2->ly[i]));
14077c3cd84eSPatrick Sanan     }
14087c3cd84eSPatrick Sanan     if (compatibleLocal && da1->dim > 2) {
1409ad540459SPierre Jolivet       for (i = 0; i < dd1->p; ++i) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lz[i] == dd2->lz[i]));
14107c3cd84eSPatrick Sanan     }
14117c3cd84eSPatrick Sanan     *compatible = compatibleLocal;
14127c3cd84eSPatrick Sanan     *set        = PETSC_TRUE;
14137c3cd84eSPatrick Sanan   } else {
14147c3cd84eSPatrick Sanan     /* Decline to determine compatibility with other DM types */
14157c3cd84eSPatrick Sanan     *set = PETSC_FALSE;
14167c3cd84eSPatrick Sanan   }
14173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14187c3cd84eSPatrick Sanan }
1419