xref: /petsc/src/dm/impls/da/da.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
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 
1612b4a537SBarry 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 
1912b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `PetscSplitOwnership()`
2047c6ae99SBarry Smith @*/
DMDASetSizes(DM da,PetscInt M,PetscInt N,PetscInt P)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`
4812b4a537SBarry Smith . m  - the number of X processes (or `PETSC_DECIDE`)
4912b4a537SBarry Smith . n  - the number of Y processes (or `PETSC_DECIDE`)
5012b4a537SBarry Smith - p  - the number of Z processes (or `PETSC_DECIDE`)
5147c6ae99SBarry Smith 
5247c6ae99SBarry Smith   Level: intermediate
5347c6ae99SBarry Smith 
5412b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDASetSizes()`, `PetscSplitOwnership()`
5547c6ae99SBarry Smith @*/
DMDASetNumProcs(DM da,PetscInt m,PetscInt n,PetscInt p)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 /*@
85fe58071bSMatthew G. Knepley   DMDAGetBoundaryType - Gets the type of ghost nodes on domain boundaries.
86fe58071bSMatthew G. Knepley 
87fe58071bSMatthew G. Knepley   Not Collective
88fe58071bSMatthew G. Knepley 
89fe58071bSMatthew G. Knepley   Input Parameter:
90fe58071bSMatthew G. Knepley . da - The `DMDA`
91fe58071bSMatthew G. Knepley 
92fe58071bSMatthew G. Knepley   Output Parameters:
93fe58071bSMatthew G. Knepley + bx - x boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`
94fe58071bSMatthew G. Knepley . by - y boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`
95fe58071bSMatthew G. Knepley - bz - z boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`
96fe58071bSMatthew G. Knepley 
97fe58071bSMatthew G. Knepley   Level: intermediate
98fe58071bSMatthew G. Knepley 
99fe58071bSMatthew G. Knepley .seealso: [](sec_struct), `DMDASetBoundaryType()`, `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMBoundaryType`, `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`
100fe58071bSMatthew G. Knepley @*/
DMDAGetBoundaryType(DM da,PeOp DMBoundaryType * bx,PeOp DMBoundaryType * by,PeOp DMBoundaryType * bz)101ce78bad3SBarry Smith PetscErrorCode DMDAGetBoundaryType(DM da, PeOp DMBoundaryType *bx, PeOp DMBoundaryType *by, PeOp DMBoundaryType *bz)
102fe58071bSMatthew G. Knepley {
103fe58071bSMatthew G. Knepley   DM_DA *dd = (DM_DA *)da->data;
104fe58071bSMatthew G. Knepley 
105fe58071bSMatthew G. Knepley   PetscFunctionBegin;
106fe58071bSMatthew G. Knepley   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
107fe58071bSMatthew G. Knepley   if (bx) PetscAssertPointer(bx, 2);
108fe58071bSMatthew G. Knepley   if (by) PetscAssertPointer(by, 3);
109fe58071bSMatthew G. Knepley   if (bz) PetscAssertPointer(bz, 4);
110fe58071bSMatthew G. Knepley   if (bx) *bx = dd->bx;
111fe58071bSMatthew G. Knepley   if (by) *by = dd->by;
112fe58071bSMatthew G. Knepley   if (bz) *bz = dd->bz;
113fe58071bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
114fe58071bSMatthew G. Knepley }
115fe58071bSMatthew G. Knepley 
116fe58071bSMatthew G. Knepley /*@
1170b4b7b1cSBarry Smith   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries for a `DMDA` object.
11847c6ae99SBarry Smith 
11920f4b53cSBarry Smith   Not Collective
12047c6ae99SBarry Smith 
121d8d19677SJose E. Roman   Input Parameters:
122dce8aebaSBarry Smith + da - The `DMDA`
123a4e35b19SJacob Faibussowitsch . bx - x boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`
124a4e35b19SJacob Faibussowitsch . by - y boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`
125a4e35b19SJacob Faibussowitsch - bz - z boundary type, one of `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith   Level: intermediate
12847c6ae99SBarry Smith 
1290b4b7b1cSBarry Smith   Note:
1300b4b7b1cSBarry Smith   The default is `DM_BOUNDARY_NONE`
1310b4b7b1cSBarry Smith 
132fe58071bSMatthew G. Knepley .seealso: [](sec_struct), `DMDAGetBoundaryType()`, `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMBoundaryType`, `DM_BOUNDARY_NONE`, `DM_BOUNDARY_GHOSTED`, `DM_BOUNDARY_PERIODIC`
13347c6ae99SBarry Smith @*/
DMDASetBoundaryType(DM da,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz)134d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetBoundaryType(DM da, DMBoundaryType bx, DMBoundaryType by, DMBoundaryType bz)
135d71ae5a4SJacob Faibussowitsch {
13647c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
13747c6ae99SBarry Smith 
13847c6ae99SBarry Smith   PetscFunctionBegin;
139a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1405a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da, bx, 2);
1415a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da, by, 3);
1425a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da, bz, 4);
1437a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
1441321219cSEthan Coon   dd->bx = bx;
1451321219cSEthan Coon   dd->by = by;
1461321219cSEthan Coon   dd->bz = bz;
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14847c6ae99SBarry Smith }
14947c6ae99SBarry Smith 
15047c6ae99SBarry Smith /*@
151aa219208SBarry Smith   DMDASetDof - Sets the number of degrees of freedom per vertex
15247c6ae99SBarry Smith 
15320f4b53cSBarry Smith   Not Collective
15447c6ae99SBarry Smith 
15559f3ab6dSMatthew G. Knepley   Input Parameters:
156dce8aebaSBarry Smith + da  - The `DMDA`
15712b4a537SBarry Smith - dof - Number of degrees of freedom per vertex
15847c6ae99SBarry Smith 
15947c6ae99SBarry Smith   Level: intermediate
16047c6ae99SBarry Smith 
16112b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetDof()`, `DMDACreate()`, `DMDestroy()`
16247c6ae99SBarry Smith @*/
DMDASetDof(DM da,PetscInt dof)163d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetDof(DM da, PetscInt dof)
164d71ae5a4SJacob Faibussowitsch {
16547c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
16647c6ae99SBarry Smith 
16747c6ae99SBarry Smith   PetscFunctionBegin;
168a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
16954cfb0beSLisandro Dalcin   PetscValidLogicalCollectiveInt(da, dof, 2);
1707a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
17147c6ae99SBarry Smith   dd->w  = dof;
1721411c6eeSJed Brown   da->bs = dof;
1733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17447c6ae99SBarry Smith }
17547c6ae99SBarry Smith 
176fb6725baSMatthew G. Knepley /*@
177fb6725baSMatthew G. Knepley   DMDAGetDof - Gets the number of degrees of freedom per vertex
178fb6725baSMatthew G. Knepley 
17920f4b53cSBarry Smith   Not Collective
180fb6725baSMatthew G. Knepley 
181fb6725baSMatthew G. Knepley   Input Parameter:
182dce8aebaSBarry Smith . da - The `DMDA`
183fb6725baSMatthew G. Knepley 
184fb6725baSMatthew G. Knepley   Output Parameter:
18512b4a537SBarry Smith . dof - Number of degrees of freedom per vertex
186fb6725baSMatthew G. Knepley 
187fb6725baSMatthew G. Knepley   Level: intermediate
188fb6725baSMatthew G. Knepley 
18912b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDASetDof()`, `DMDACreate()`, `DMDestroy()`
190fb6725baSMatthew G. Knepley @*/
DMDAGetDof(DM da,PetscInt * dof)191d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetDof(DM da, PetscInt *dof)
192d71ae5a4SJacob Faibussowitsch {
193fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *)da->data;
194fb6725baSMatthew G. Knepley 
195fb6725baSMatthew G. Knepley   PetscFunctionBegin;
196a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1974f572ea9SToby Isaac   PetscAssertPointer(dof, 2);
198fb6725baSMatthew G. Knepley   *dof = dd->w;
1993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
200fb6725baSMatthew G. Knepley }
201fb6725baSMatthew G. Knepley 
2027ddda789SPeter Brune /*@
2037ddda789SPeter Brune   DMDAGetOverlap - Gets the size of the per-processor overlap.
2047ddda789SPeter Brune 
20520f4b53cSBarry Smith   Not Collective
2067ddda789SPeter Brune 
207f899ff85SJose E. Roman   Input Parameter:
208dce8aebaSBarry Smith . da - The `DMDA`
2097ddda789SPeter Brune 
2107ddda789SPeter Brune   Output Parameters:
2117ddda789SPeter Brune + x - Overlap in the x direction
2127ddda789SPeter Brune . y - Overlap in the y direction
2137ddda789SPeter Brune - z - Overlap in the z direction
2147ddda789SPeter Brune 
2157ddda789SPeter Brune   Level: intermediate
2167ddda789SPeter Brune 
21712b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDASetOverlap()`
2187ddda789SPeter Brune @*/
DMDAGetOverlap(DM da,PeOp PetscInt * x,PeOp PetscInt * y,PeOp PetscInt * z)219ce78bad3SBarry Smith PetscErrorCode DMDAGetOverlap(DM da, PeOp PetscInt *x, PeOp PetscInt *y, PeOp PetscInt *z)
220d71ae5a4SJacob Faibussowitsch {
2217ddda789SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
2227ddda789SPeter Brune 
2237ddda789SPeter Brune   PetscFunctionBegin;
224a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
2257ddda789SPeter Brune   if (x) *x = dd->xol;
2267ddda789SPeter Brune   if (y) *y = dd->yol;
2277ddda789SPeter Brune   if (z) *z = dd->zol;
2283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2297ddda789SPeter Brune }
2307ddda789SPeter Brune 
23188661749SPeter Brune /*@
23288661749SPeter Brune   DMDASetOverlap - Sets the size of the per-processor overlap.
23388661749SPeter Brune 
23420f4b53cSBarry Smith   Not Collective
23588661749SPeter Brune 
2367ddda789SPeter Brune   Input Parameters:
237dce8aebaSBarry Smith + da - The `DMDA`
2387ddda789SPeter Brune . x  - Overlap in the x direction
2397ddda789SPeter Brune . y  - Overlap in the y direction
2407ddda789SPeter Brune - z  - Overlap in the z direction
24188661749SPeter Brune 
24288661749SPeter Brune   Level: intermediate
24388661749SPeter Brune 
24412b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDAGetOverlap()`
24588661749SPeter Brune @*/
DMDASetOverlap(DM da,PetscInt x,PetscInt y,PetscInt z)246d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOverlap(DM da, PetscInt x, PetscInt y, PetscInt z)
247d71ae5a4SJacob Faibussowitsch {
24888661749SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
24988661749SPeter Brune 
25088661749SPeter Brune   PetscFunctionBegin;
251a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
2527ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da, x, 2);
2537ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da, y, 3);
2547ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da, z, 4);
2557ddda789SPeter Brune   dd->xol = x;
2567ddda789SPeter Brune   dd->yol = y;
2577ddda789SPeter Brune   dd->zol = z;
2583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25988661749SPeter Brune }
26088661749SPeter Brune 
2613e7870d2SPeter Brune /*@
26212b4a537SBarry Smith   DMDAGetNumLocalSubDomains - Gets the number of local subdomains that would be created upon decomposition.
2633e7870d2SPeter Brune 
26420f4b53cSBarry Smith   Not Collective
2653e7870d2SPeter Brune 
2662fe279fdSBarry Smith   Input Parameter:
267dce8aebaSBarry Smith . da - The `DMDA`
2683e7870d2SPeter Brune 
2692fe279fdSBarry Smith   Output Parameter:
270a2b725a8SWilliam Gropp . Nsub - Number of local subdomains created upon decomposition
2713e7870d2SPeter Brune 
2723e7870d2SPeter Brune   Level: intermediate
2733e7870d2SPeter Brune 
27412b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDASetNumLocalSubDomains()`
2753e7870d2SPeter Brune @*/
DMDAGetNumLocalSubDomains(DM da,PetscInt * Nsub)276d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumLocalSubDomains(DM da, PetscInt *Nsub)
277d71ae5a4SJacob Faibussowitsch {
2783e7870d2SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
2793e7870d2SPeter Brune 
2803e7870d2SPeter Brune   PetscFunctionBegin;
281a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
2823e7870d2SPeter Brune   if (Nsub) *Nsub = dd->Nsub;
2833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2843e7870d2SPeter Brune }
2853e7870d2SPeter Brune 
2863e7870d2SPeter Brune /*@
28712b4a537SBarry Smith   DMDASetNumLocalSubDomains - Sets the number of local subdomains to create when decomposing with `DMCreateDomainDecomposition()`
2883e7870d2SPeter Brune 
28920f4b53cSBarry Smith   Not Collective
2903e7870d2SPeter Brune 
2913e7870d2SPeter Brune   Input Parameters:
292dce8aebaSBarry Smith + da   - The `DMDA`
2933e7870d2SPeter Brune - Nsub - The number of local subdomains requested
2943e7870d2SPeter Brune 
2953e7870d2SPeter Brune   Level: intermediate
2963e7870d2SPeter Brune 
29712b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateDomainDecomposition()`, `DMDAGetNumLocalSubDomains()`
2983e7870d2SPeter Brune @*/
DMDASetNumLocalSubDomains(DM da,PetscInt Nsub)299d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetNumLocalSubDomains(DM da, PetscInt Nsub)
300d71ae5a4SJacob Faibussowitsch {
3013e7870d2SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
3023e7870d2SPeter Brune 
3033e7870d2SPeter Brune   PetscFunctionBegin;
304a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
3053e7870d2SPeter Brune   PetscValidLogicalCollectiveInt(da, Nsub, 2);
3063e7870d2SPeter Brune   dd->Nsub = Nsub;
3073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3083e7870d2SPeter Brune }
3093e7870d2SPeter Brune 
310d886c4f4SPeter Brune /*@
31112b4a537SBarry Smith   DMDASetOffset - Sets the index offset of the `DMDA`.
312d886c4f4SPeter Brune 
31320f4b53cSBarry Smith   Collective
314d886c4f4SPeter Brune 
315d8d19677SJose E. Roman   Input Parameters:
316dce8aebaSBarry Smith + da - The `DMDA`
317d886c4f4SPeter Brune . xo - The offset in the x direction
318d886c4f4SPeter Brune . yo - The offset in the y direction
319a4e35b19SJacob Faibussowitsch . zo - The offset in the z direction
320a4e35b19SJacob Faibussowitsch . Mo - The problem offset in the x direction
321a4e35b19SJacob Faibussowitsch . No - The problem offset in the y direction
322a4e35b19SJacob Faibussowitsch - Po - The problem offset in the z direction
323d886c4f4SPeter Brune 
32412b4a537SBarry Smith   Level: developer
325d886c4f4SPeter Brune 
326dce8aebaSBarry Smith   Note:
327dce8aebaSBarry Smith   This is used primarily to overlap a computation on a local `DMDA` with that on a global `DMDA` without
328d886c4f4SPeter Brune   changing boundary conditions or subdomain features that depend upon the global offsets.
329d886c4f4SPeter Brune 
33012b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()`
331d886c4f4SPeter Brune @*/
DMDASetOffset(DM da,PetscInt xo,PetscInt yo,PetscInt zo,PetscInt Mo,PetscInt No,PetscInt Po)332d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po)
333d71ae5a4SJacob Faibussowitsch {
334d886c4f4SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
335d886c4f4SPeter Brune 
336d886c4f4SPeter Brune   PetscFunctionBegin;
337a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
338d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da, xo, 2);
33995c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, yo, 3);
34095c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, zo, 4);
34195c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, Mo, 5);
34295c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, No, 6);
34395c13181SPeter Brune   PetscValidLogicalCollectiveInt(da, Po, 7);
344d886c4f4SPeter Brune   dd->xo = xo;
345d886c4f4SPeter Brune   dd->yo = yo;
346d886c4f4SPeter Brune   dd->zo = zo;
34795c13181SPeter Brune   dd->Mo = Mo;
34895c13181SPeter Brune   dd->No = No;
34995c13181SPeter Brune   dd->Po = Po;
35095c13181SPeter Brune 
3516858538eSMatthew G. Knepley   if (da->coordinates[0].dm) PetscCall(DMDASetOffset(da->coordinates[0].dm, xo, yo, zo, Mo, No, Po));
3523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
353d886c4f4SPeter Brune }
354d886c4f4SPeter Brune 
355d886c4f4SPeter Brune /*@
356dce8aebaSBarry Smith   DMDAGetOffset - Gets the index offset of the `DMDA`.
357d886c4f4SPeter Brune 
35820f4b53cSBarry Smith   Not Collective
359d886c4f4SPeter Brune 
360d886c4f4SPeter Brune   Input Parameter:
361dce8aebaSBarry Smith . da - The `DMDA`
362d886c4f4SPeter Brune 
363d886c4f4SPeter Brune   Output Parameters:
364d886c4f4SPeter Brune + xo - The offset in the x direction
365d886c4f4SPeter Brune . yo - The offset in the y direction
36695c13181SPeter Brune . zo - The offset in the z direction
36795c13181SPeter Brune . Mo - The global size in the x direction
36895c13181SPeter Brune . No - The global size in the y direction
36995c13181SPeter Brune - Po - The global size in the z direction
370d886c4f4SPeter Brune 
37112b4a537SBarry Smith   Level: developer
372d886c4f4SPeter Brune 
37312b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDASetOffset()`, `DMDAVecGetArray()`
374d886c4f4SPeter Brune @*/
DMDAGetOffset(DM da,PeOp PetscInt * xo,PeOp PetscInt * yo,PeOp PetscInt * zo,PeOp PetscInt * Mo,PeOp PetscInt * No,PeOp PetscInt * Po)375ce78bad3SBarry Smith PetscErrorCode DMDAGetOffset(DM da, PeOp PetscInt *xo, PeOp PetscInt *yo, PeOp PetscInt *zo, PeOp PetscInt *Mo, PeOp PetscInt *No, PeOp PetscInt *Po)
376d71ae5a4SJacob Faibussowitsch {
377d886c4f4SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
378d886c4f4SPeter Brune 
379d886c4f4SPeter Brune   PetscFunctionBegin;
380a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
381d886c4f4SPeter Brune   if (xo) *xo = dd->xo;
382d886c4f4SPeter Brune   if (yo) *yo = dd->yo;
383d886c4f4SPeter Brune   if (zo) *zo = dd->zo;
38495c13181SPeter Brune   if (Mo) *Mo = dd->Mo;
38595c13181SPeter Brune   if (No) *No = dd->No;
38695c13181SPeter Brune   if (Po) *Po = dd->Po;
3873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
388d886c4f4SPeter Brune }
389d886c4f4SPeter Brune 
39040234c92SPeter Brune /*@
39112b4a537SBarry Smith   DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain `DMDA`.
39240234c92SPeter Brune 
39320f4b53cSBarry Smith   Not Collective
39440234c92SPeter Brune 
39540234c92SPeter Brune   Input Parameter:
396dce8aebaSBarry Smith . da - The `DMDA`
39740234c92SPeter Brune 
39840234c92SPeter Brune   Output Parameters:
39940234c92SPeter Brune + xs - The start of the region in x
40040234c92SPeter Brune . ys - The start of the region in y
40140234c92SPeter Brune . zs - The start of the region in z
402a4e35b19SJacob Faibussowitsch . xm - The size of the region in x
403a4e35b19SJacob Faibussowitsch . ym - The size of the region in y
404a4e35b19SJacob Faibussowitsch - zm - The size of the region in z
40540234c92SPeter Brune 
40640234c92SPeter Brune   Level: intermediate
40740234c92SPeter Brune 
40812b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()`
40940234c92SPeter Brune @*/
DMDAGetNonOverlappingRegion(DM da,PeOp PetscInt * xs,PeOp PetscInt * ys,PeOp PetscInt * zs,PeOp PetscInt * xm,PeOp PetscInt * ym,PeOp PetscInt * zm)410ce78bad3SBarry Smith PetscErrorCode DMDAGetNonOverlappingRegion(DM da, PeOp PetscInt *xs, PeOp PetscInt *ys, PeOp PetscInt *zs, PeOp PetscInt *xm, PeOp PetscInt *ym, PeOp PetscInt *zm)
411d71ae5a4SJacob Faibussowitsch {
41240234c92SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
41340234c92SPeter Brune 
41440234c92SPeter Brune   PetscFunctionBegin;
415a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
41640234c92SPeter Brune   if (xs) *xs = dd->nonxs;
41740234c92SPeter Brune   if (ys) *ys = dd->nonys;
41840234c92SPeter Brune   if (zs) *zs = dd->nonzs;
41940234c92SPeter Brune   if (xm) *xm = dd->nonxm;
42040234c92SPeter Brune   if (ym) *ym = dd->nonym;
42140234c92SPeter Brune   if (zm) *zm = dd->nonzm;
4223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42340234c92SPeter Brune }
42440234c92SPeter Brune 
42540234c92SPeter Brune /*@
42612b4a537SBarry Smith   DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain `DMDA`.
42740234c92SPeter Brune 
42820f4b53cSBarry Smith   Collective
42940234c92SPeter Brune 
430d8d19677SJose E. Roman   Input Parameters:
431dce8aebaSBarry Smith + da - The `DMDA`
43240234c92SPeter Brune . xs - The start of the region in x
43340234c92SPeter Brune . ys - The start of the region in y
43440234c92SPeter Brune . zs - The start of the region in z
435a4e35b19SJacob Faibussowitsch . xm - The size of the region in x
436a4e35b19SJacob Faibussowitsch . ym - The size of the region in y
437a4e35b19SJacob Faibussowitsch - zm - The size of the region in z
43840234c92SPeter Brune 
43940234c92SPeter Brune   Level: intermediate
44040234c92SPeter Brune 
44112b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetOffset()`, `DMDAVecGetArray()`
44240234c92SPeter Brune @*/
DMDASetNonOverlappingRegion(DM da,PetscInt xs,PetscInt ys,PetscInt zs,PetscInt xm,PetscInt ym,PetscInt zm)443d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm)
444d71ae5a4SJacob Faibussowitsch {
44540234c92SPeter Brune   DM_DA *dd = (DM_DA *)da->data;
44640234c92SPeter Brune 
44740234c92SPeter Brune   PetscFunctionBegin;
448a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
44940234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, xs, 2);
45040234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, ys, 3);
45140234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, zs, 4);
45240234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, xm, 5);
45340234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, ym, 6);
45440234c92SPeter Brune   PetscValidLogicalCollectiveInt(da, zm, 7);
45540234c92SPeter Brune   dd->nonxs = xs;
45640234c92SPeter Brune   dd->nonys = ys;
45740234c92SPeter Brune   dd->nonzs = zs;
45840234c92SPeter Brune   dd->nonxm = xm;
45940234c92SPeter Brune   dd->nonym = ym;
46040234c92SPeter Brune   dd->nonzm = zm;
4613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46240234c92SPeter Brune }
46388661749SPeter Brune 
46447c6ae99SBarry Smith /*@
465aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
46647c6ae99SBarry Smith 
46720f4b53cSBarry Smith   Logically Collective
46847c6ae99SBarry Smith 
469d8d19677SJose E. Roman   Input Parameters:
470dce8aebaSBarry Smith + da    - The `DMDA`
471dce8aebaSBarry Smith - stype - The stencil type, use either `DMDA_STENCIL_BOX` or `DMDA_STENCIL_STAR`.
47247c6ae99SBarry Smith 
47347c6ae99SBarry Smith   Level: intermediate
47447c6ae99SBarry Smith 
47512b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMDAStencilType`, `DMDA_STENCIL_BOX`, `DMDA_STENCIL_STAR.`
47647c6ae99SBarry Smith @*/
DMDASetStencilType(DM da,DMDAStencilType stype)477d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype)
478d71ae5a4SJacob Faibussowitsch {
47947c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
48047c6ae99SBarry Smith 
48147c6ae99SBarry Smith   PetscFunctionBegin;
482a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
48347c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da, stype, 2);
4847a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
48547c6ae99SBarry Smith   dd->stencil_type = stype;
4863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48747c6ae99SBarry Smith }
48847c6ae99SBarry Smith 
489fb6725baSMatthew G. Knepley /*@
490fb6725baSMatthew G. Knepley   DMDAGetStencilType - Gets the type of the communication stencil
491fb6725baSMatthew G. Knepley 
49220f4b53cSBarry Smith   Not Collective
493fb6725baSMatthew G. Knepley 
494fb6725baSMatthew G. Knepley   Input Parameter:
495dce8aebaSBarry Smith . da - The `DMDA`
496fb6725baSMatthew G. Knepley 
497fb6725baSMatthew G. Knepley   Output Parameter:
498dce8aebaSBarry Smith . stype - The stencil type, use either `DMDA_STENCIL_BOX` or `DMDA_STENCIL_STAR`.
499fb6725baSMatthew G. Knepley 
500fb6725baSMatthew G. Knepley   Level: intermediate
501fb6725baSMatthew G. Knepley 
50212b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMDAStencilType`, `DMDA_STENCIL_BOX`, `DMDA_STENCIL_STAR.`
503fb6725baSMatthew G. Knepley @*/
DMDAGetStencilType(DM da,DMDAStencilType * stype)504d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype)
505d71ae5a4SJacob Faibussowitsch {
506fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *)da->data;
507fb6725baSMatthew G. Knepley 
508fb6725baSMatthew G. Knepley   PetscFunctionBegin;
509a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
5104f572ea9SToby Isaac   PetscAssertPointer(stype, 2);
511fb6725baSMatthew G. Knepley   *stype = dd->stencil_type;
5123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
513fb6725baSMatthew G. Knepley }
514fb6725baSMatthew G. Knepley 
51547c6ae99SBarry Smith /*@
516aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
51747c6ae99SBarry Smith 
51820f4b53cSBarry Smith   Logically Collective
51947c6ae99SBarry Smith 
520d8d19677SJose E. Roman   Input Parameters:
521dce8aebaSBarry Smith + da    - The `DMDA`
52247c6ae99SBarry Smith - width - The stencil width
52347c6ae99SBarry Smith 
52447c6ae99SBarry Smith   Level: intermediate
52547c6ae99SBarry Smith 
52612b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMDAStencilType`, `DMDA_STENCIL_BOX`, `DMDA_STENCIL_STAR.`
52747c6ae99SBarry Smith @*/
DMDASetStencilWidth(DM da,PetscInt width)528d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width)
529d71ae5a4SJacob Faibussowitsch {
53047c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
53147c6ae99SBarry Smith 
53247c6ae99SBarry Smith   PetscFunctionBegin;
533a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
53447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, width, 2);
5357a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
53647c6ae99SBarry Smith   dd->s = width;
5373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53847c6ae99SBarry Smith }
53947c6ae99SBarry Smith 
540fb6725baSMatthew G. Knepley /*@
541fb6725baSMatthew G. Knepley   DMDAGetStencilWidth - Gets the width of the communication stencil
542fb6725baSMatthew G. Knepley 
54320f4b53cSBarry Smith   Not Collective
544fb6725baSMatthew G. Knepley 
545fb6725baSMatthew G. Knepley   Input Parameter:
546dce8aebaSBarry Smith . da - The `DMDA`
547fb6725baSMatthew G. Knepley 
548fb6725baSMatthew G. Knepley   Output Parameter:
549fb6725baSMatthew G. Knepley . width - The stencil width
550fb6725baSMatthew G. Knepley 
551fb6725baSMatthew G. Knepley   Level: intermediate
552fb6725baSMatthew G. Knepley 
55312b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`, `DMDAStencilType`, `DMDA_STENCIL_BOX`, `DMDA_STENCIL_STAR.`
554fb6725baSMatthew G. Knepley @*/
DMDAGetStencilWidth(DM da,PetscInt * width)555d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width)
556d71ae5a4SJacob Faibussowitsch {
557fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *)da->data;
558fb6725baSMatthew G. Knepley 
559fb6725baSMatthew G. Knepley   PetscFunctionBegin;
560a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
5614f572ea9SToby Isaac   PetscAssertPointer(width, 2);
562fb6725baSMatthew G. Knepley   *width = dd->s;
5633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
564fb6725baSMatthew G. Knepley }
565fb6725baSMatthew G. Knepley 
DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])566d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da, PetscInt M, PetscInt m, const PetscInt lx[])
567d71ae5a4SJacob Faibussowitsch {
56847c6ae99SBarry Smith   PetscInt i, sum;
56947c6ae99SBarry Smith 
57047c6ae99SBarry Smith   PetscFunctionBegin;
5717a8be351SBarry Smith   PetscCheck(M >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Global dimension not set");
57247c6ae99SBarry Smith   for (i = sum = 0; i < m; i++) sum += lx[i];
57363a3b9bcSJacob 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);
5743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57547c6ae99SBarry Smith }
57647c6ae99SBarry Smith 
57747c6ae99SBarry Smith /*@
578aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
57947c6ae99SBarry Smith 
58020f4b53cSBarry Smith   Logically Collective
58147c6ae99SBarry Smith 
582d8d19677SJose E. Roman   Input Parameters:
583dce8aebaSBarry Smith + da - The `DMDA`
58412b4a537SBarry 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
58512b4a537SBarry 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
58612b4a537SBarry 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.
58747c6ae99SBarry Smith 
58847c6ae99SBarry Smith   Level: intermediate
58947c6ae99SBarry Smith 
590a4e35b19SJacob Faibussowitsch   Note:
591a4e35b19SJacob Faibussowitsch   These numbers are NOT multiplied by the number of dof per node.
592e3f3e4b6SBarry Smith 
59312b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `DMDestroy()`
59447c6ae99SBarry Smith @*/
DMDASetOwnershipRanges(DM da,const PetscInt lx[],const PetscInt ly[],const PetscInt lz[])595d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
596d71ae5a4SJacob Faibussowitsch {
59747c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
59847c6ae99SBarry Smith 
59947c6ae99SBarry Smith   PetscFunctionBegin;
600a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
6017a8be351SBarry Smith   PetscCheck(!da->setupcalled, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "This function must be called before DMSetUp()");
60247c6ae99SBarry Smith   if (lx) {
60312b4a537SBarry Smith     PetscCheck(dd->m >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of processes");
6049566063dSJacob Faibussowitsch     PetscCall(DMDACheckOwnershipRanges_Private(da, dd->M, dd->m, lx));
60548a46eb9SPierre Jolivet     if (!dd->lx) PetscCall(PetscMalloc1(dd->m, &dd->lx));
6069566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd->lx, lx, dd->m));
60747c6ae99SBarry Smith   }
60847c6ae99SBarry Smith   if (ly) {
60912b4a537SBarry Smith     PetscCheck(dd->n >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of processes");
6109566063dSJacob Faibussowitsch     PetscCall(DMDACheckOwnershipRanges_Private(da, dd->N, dd->n, ly));
61148a46eb9SPierre Jolivet     if (!dd->ly) PetscCall(PetscMalloc1(dd->n, &dd->ly));
6129566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd->ly, ly, dd->n));
61347c6ae99SBarry Smith   }
61447c6ae99SBarry Smith   if (lz) {
61512b4a537SBarry Smith     PetscCheck(dd->p >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set ownership ranges before setting number of processes");
6169566063dSJacob Faibussowitsch     PetscCall(DMDACheckOwnershipRanges_Private(da, dd->P, dd->p, lz));
61748a46eb9SPierre Jolivet     if (!dd->lz) PetscCall(PetscMalloc1(dd->p, &dd->lz));
6189566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd->lz, lz, dd->p));
61947c6ae99SBarry Smith   }
6203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
62147c6ae99SBarry Smith }
62247c6ae99SBarry Smith 
62347c6ae99SBarry Smith /*@
624aa219208SBarry Smith   DMDASetInterpolationType - Sets the type of interpolation that will be
625dce8aebaSBarry Smith   returned by `DMCreateInterpolation()`
62647c6ae99SBarry Smith 
62720f4b53cSBarry Smith   Logically Collective
62847c6ae99SBarry Smith 
629d8d19677SJose E. Roman   Input Parameters:
63047c6ae99SBarry Smith + da    - initial distributed array
631dce8aebaSBarry Smith - ctype - `DMDA_Q1` and `DMDA_Q0` are currently the only supported forms
63247c6ae99SBarry Smith 
63347c6ae99SBarry Smith   Level: intermediate
63447c6ae99SBarry Smith 
635dce8aebaSBarry Smith   Note:
636dce8aebaSBarry Smith   You should call this on the coarser of the two `DMDA` you pass to `DMCreateInterpolation()`
63747c6ae99SBarry Smith 
63812b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDestroy()`, `DMDAInterpolationType`,
63912b4a537SBarry Smith           `DMDA_Q1`, `DMDA_Q0`
64047c6ae99SBarry Smith @*/
DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)641d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetInterpolationType(DM da, DMDAInterpolationType ctype)
642d71ae5a4SJacob Faibussowitsch {
64347c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
64447c6ae99SBarry Smith 
64547c6ae99SBarry Smith   PetscFunctionBegin;
646a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
64747c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da, ctype, 2);
64847c6ae99SBarry Smith   dd->interptype = ctype;
6493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
65047c6ae99SBarry Smith }
65147c6ae99SBarry Smith 
6522dde6fd4SLisandro Dalcin /*@
6532dde6fd4SLisandro Dalcin   DMDAGetInterpolationType - Gets the type of interpolation that will be
654dce8aebaSBarry Smith   used by `DMCreateInterpolation()`
6552dde6fd4SLisandro Dalcin 
6562dde6fd4SLisandro Dalcin   Not Collective
6572dde6fd4SLisandro Dalcin 
6582dde6fd4SLisandro Dalcin   Input Parameter:
6592dde6fd4SLisandro Dalcin . da - distributed array
6602dde6fd4SLisandro Dalcin 
6612dde6fd4SLisandro Dalcin   Output Parameter:
662dce8aebaSBarry Smith . ctype - interpolation type (`DMDA_Q1` and `DMDA_Q0` are currently the only supported forms)
6632dde6fd4SLisandro Dalcin 
6642dde6fd4SLisandro Dalcin   Level: intermediate
6652dde6fd4SLisandro Dalcin 
66612b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAInterpolationType`, `DMDASetInterpolationType()`, `DMCreateInterpolation()`,
66712b4a537SBarry Smith           `DMDA_Q1`, `DMDA_Q0`
6682dde6fd4SLisandro Dalcin @*/
DMDAGetInterpolationType(DM da,DMDAInterpolationType * ctype)669d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetInterpolationType(DM da, DMDAInterpolationType *ctype)
670d71ae5a4SJacob Faibussowitsch {
6712dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA *)da->data;
6722dde6fd4SLisandro Dalcin 
6732dde6fd4SLisandro Dalcin   PetscFunctionBegin;
674a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
6754f572ea9SToby Isaac   PetscAssertPointer(ctype, 2);
6762dde6fd4SLisandro Dalcin   *ctype = dd->interptype;
6773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6782dde6fd4SLisandro Dalcin }
67947c6ae99SBarry Smith 
6806a119db4SBarry Smith /*@C
681aa219208SBarry Smith   DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
68247c6ae99SBarry Smith   processes neighbors.
68347c6ae99SBarry Smith 
68447c6ae99SBarry Smith   Not Collective
68547c6ae99SBarry Smith 
68647c6ae99SBarry Smith   Input Parameter:
687dce8aebaSBarry Smith . da - the `DMDA` object
68847c6ae99SBarry Smith 
6892fe279fdSBarry Smith   Output Parameter:
69012b4a537SBarry Smith . ranks - the neighbors ranks, stored with the x index increasing most rapidly. The process itself is in the list
69147c6ae99SBarry Smith 
69247c6ae99SBarry Smith   Level: intermediate
69347c6ae99SBarry Smith 
694dce8aebaSBarry Smith   Notes:
695f13dfd9eSBarry Smith   In 2d the `ranks` is of length 9, in 3d of length 27
696dce8aebaSBarry Smith 
697dce8aebaSBarry Smith   Not supported in 1d
698dce8aebaSBarry Smith 
699dce8aebaSBarry Smith   Do not free the array, it is freed when the `DMDA` is destroyed.
700dce8aebaSBarry Smith 
70112b4a537SBarry Smith   Fortran Note:
702ce78bad3SBarry Smith   Use `DMDARestoreNeighbors()` to return the array when no longer needed
703dce8aebaSBarry Smith 
70412b4a537SBarry Smith .seealso: [](sec_struct), `DMDA`, `DM`
70547c6ae99SBarry Smith @*/
DMDAGetNeighbors(DM da,const PetscMPIInt * ranks[])706d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNeighbors(DM da, const PetscMPIInt *ranks[])
707d71ae5a4SJacob Faibussowitsch {
70847c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
7095fd66863SKarl Rupp 
71047c6ae99SBarry Smith   PetscFunctionBegin;
711a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
71247c6ae99SBarry Smith   *ranks = dd->neighbors;
7133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71447c6ae99SBarry Smith }
71547c6ae99SBarry Smith 
71647c6ae99SBarry Smith /*@C
717ce78bad3SBarry Smith   DMDAGetOwnershipRanges - Gets the number of indices in the x, y and z direction that are owned by each process in that direction
71847c6ae99SBarry Smith 
71947c6ae99SBarry Smith   Not Collective
72047c6ae99SBarry Smith 
72147c6ae99SBarry Smith   Input Parameter:
722dce8aebaSBarry Smith . da - the `DMDA` object
72347c6ae99SBarry Smith 
724d8d19677SJose E. Roman   Output Parameters:
725ce78bad3SBarry Smith + lx - ownership along x direction (optional), its length is `m` the number of processes in the x-direction
726ce78bad3SBarry Smith . ly - ownership along y direction (optional), its length is `n` the number of processes in the y-direction
727ce78bad3SBarry Smith - lz - ownership along z direction (optional), its length is `p` the number of processes in the z-direction
72847c6ae99SBarry Smith 
72947c6ae99SBarry Smith   Level: intermediate
73047c6ae99SBarry Smith 
731ce78bad3SBarry Smith   Notes:
732dce8aebaSBarry Smith   These correspond to the optional final arguments passed to `DMDACreate()`, `DMDACreate2d()`, `DMDACreate3d()`
73347c6ae99SBarry Smith 
734ce78bad3SBarry Smith   You should not free these arrays, nor change the values in them. They will only have valid values while the
735dce8aebaSBarry Smith   `DMDA` they came from still exists (has not been destroyed).
73647c6ae99SBarry Smith 
737e3f3e4b6SBarry Smith   These numbers are NOT multiplied by the number of dof per node.
738e3f3e4b6SBarry Smith 
739ce78bad3SBarry Smith   The meaning of these is different than that returned by `VecGetOwnerShipRanges()`
740ce78bad3SBarry Smith 
741ce78bad3SBarry Smith   Fortran Notes:
742ce78bad3SBarry Smith   Pass `PETSC_NULL_INT_POINTER` for any array not needed.
743ce78bad3SBarry Smith 
744ce78bad3SBarry Smith   Use `DMDARestoreOwershipRange()` to return the arrays when no longer needed
745dce8aebaSBarry Smith 
74612b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetCorners()`, `DMDAGetGhostCorners()`, `DMDACreate()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `VecGetOwnershipRanges()`
74747c6ae99SBarry Smith @*/
DMDAGetOwnershipRanges(DM da,PeOp const PetscInt * lx[],PeOp const PetscInt * ly[],PeOp const PetscInt * lz[])748ce78bad3SBarry Smith PetscErrorCode DMDAGetOwnershipRanges(DM da, PeOp const PetscInt *lx[], PeOp const PetscInt *ly[], PeOp const PetscInt *lz[])
749d71ae5a4SJacob Faibussowitsch {
75047c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
75147c6ae99SBarry Smith 
75247c6ae99SBarry Smith   PetscFunctionBegin;
753a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
75447c6ae99SBarry Smith   if (lx) *lx = dd->lx;
75547c6ae99SBarry Smith   if (ly) *ly = dd->ly;
75647c6ae99SBarry Smith   if (lz) *lz = dd->lz;
7573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75847c6ae99SBarry Smith }
75947c6ae99SBarry Smith 
76047c6ae99SBarry Smith /*@
761dce8aebaSBarry Smith   DMDASetRefinementFactor - Set the ratios that the `DMDA` grid is refined
76247c6ae99SBarry Smith 
76320f4b53cSBarry Smith   Logically Collective
76447c6ae99SBarry Smith 
76547c6ae99SBarry Smith   Input Parameters:
766dce8aebaSBarry Smith + da       - the `DMDA` object
76747c6ae99SBarry Smith . refine_x - ratio of fine grid to coarse in x direction (2 by default)
76847c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default)
76947c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default)
77047c6ae99SBarry Smith 
771dce8aebaSBarry Smith   Options Database Keys:
77248eeb7c8SBarry Smith + -da_refine_x refine_x - refinement ratio in x direction
77348eeb7c8SBarry Smith . -da_refine_y rafine_y - refinement ratio in y direction
77448eeb7c8SBarry Smith . -da_refine_z refine_z - refinement ratio in z direction
77512b4a537SBarry Smith - -da_refine <n>        - refine the `DMDA` object n times when it is created.
77647c6ae99SBarry Smith 
77747c6ae99SBarry Smith   Level: intermediate
77847c6ae99SBarry Smith 
779dce8aebaSBarry Smith   Note:
780dce8aebaSBarry Smith   Pass `PETSC_IGNORE` to leave a value unchanged
78147c6ae99SBarry Smith 
78212b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMRefine()`, `DMDAGetRefinementFactor()`
78347c6ae99SBarry Smith @*/
DMDASetRefinementFactor(DM da,PetscInt refine_x,PetscInt refine_y,PetscInt refine_z)784d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y, PetscInt refine_z)
785d71ae5a4SJacob Faibussowitsch {
78647c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
78747c6ae99SBarry Smith 
78847c6ae99SBarry Smith   PetscFunctionBegin;
789a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
79047c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, refine_x, 2);
79147c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, refine_y, 3);
79247c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da, refine_z, 4);
79347c6ae99SBarry Smith 
79447c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
79547c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
79647c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
7973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79847c6ae99SBarry Smith }
79947c6ae99SBarry Smith 
800cc4c1da9SBarry Smith /*@
801dce8aebaSBarry Smith   DMDAGetRefinementFactor - Gets the ratios that the `DMDA` grid is refined
80247c6ae99SBarry Smith 
80347c6ae99SBarry Smith   Not Collective
80447c6ae99SBarry Smith 
80547c6ae99SBarry Smith   Input Parameter:
806dce8aebaSBarry Smith . da - the `DMDA` object
80747c6ae99SBarry Smith 
80847c6ae99SBarry Smith   Output Parameters:
80947c6ae99SBarry Smith + refine_x - ratio of fine grid to coarse in x direction (2 by default)
81047c6ae99SBarry Smith . refine_y - ratio of fine grid to coarse in y direction (2 by default)
81147c6ae99SBarry Smith - refine_z - ratio of fine grid to coarse in z direction (2 by default)
81247c6ae99SBarry Smith 
81347c6ae99SBarry Smith   Level: intermediate
81447c6ae99SBarry Smith 
815dce8aebaSBarry Smith   Note:
81620f4b53cSBarry Smith   Pass `NULL` for values you do not need
81747c6ae99SBarry Smith 
81812b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMRefine()`, `DMDASetRefinementFactor()`
81947c6ae99SBarry Smith @*/
DMDAGetRefinementFactor(DM da,PeOp PetscInt * refine_x,PeOp PetscInt * refine_y,PeOp PetscInt * refine_z)820ce78bad3SBarry Smith PetscErrorCode DMDAGetRefinementFactor(DM da, PeOp PetscInt *refine_x, PeOp PetscInt *refine_y, PeOp PetscInt *refine_z)
821d71ae5a4SJacob Faibussowitsch {
82247c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
82347c6ae99SBarry Smith 
82447c6ae99SBarry Smith   PetscFunctionBegin;
825a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
82647c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
82747c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
82847c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
8293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83047c6ae99SBarry Smith }
83147c6ae99SBarry Smith 
83247c6ae99SBarry Smith /*@C
833dce8aebaSBarry Smith   DMDASetGetMatrix - Sets the routine used by the `DMDA` to allocate a matrix.
83447c6ae99SBarry Smith 
83520f4b53cSBarry Smith   Logically Collective; No Fortran Support
83647c6ae99SBarry Smith 
83747c6ae99SBarry Smith   Input Parameters:
838dce8aebaSBarry Smith + da - the `DMDA` object
83912b4a537SBarry Smith - f  - the function that allocates the matrix for that specific `DMDA`
84012b4a537SBarry Smith 
84112b4a537SBarry Smith   Calling sequence of `f`:
84212b4a537SBarry Smith + da - the `DMDA` object
84312b4a537SBarry Smith - A  - the created matrix
84447c6ae99SBarry Smith 
84547c6ae99SBarry Smith   Level: developer
84647c6ae99SBarry Smith 
84712b4a537SBarry Smith   Notes:
84812b4a537SBarry Smith   If the function is not provided a default function is used that uses the `DMDAStencilType`, `DMBoundaryType`, and value of `DMDASetStencilWidth()`
84912b4a537SBarry Smith   to construct the matrix.
85047c6ae99SBarry Smith 
85112b4a537SBarry Smith   See `DMDASetBlockFills()` that provides a simple way to provide the nonzero structure for
85212b4a537SBarry Smith   the diagonal and off-diagonal blocks of the matrix without providing a custom function
85312b4a537SBarry Smith 
85412b4a537SBarry Smith   Developer Note:
85512b4a537SBarry Smith   This should be called `DMDASetCreateMatrix()`
85612b4a537SBarry Smith 
85712b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMCreateMatrix()`, `DMDASetBlockFills()`
85847c6ae99SBarry Smith @*/
DMDASetGetMatrix(DM da,PetscErrorCode (* f)(DM da,Mat * A))85912b4a537SBarry Smith PetscErrorCode DMDASetGetMatrix(DM da, PetscErrorCode (*f)(DM da, Mat *A))
860d71ae5a4SJacob Faibussowitsch {
86147c6ae99SBarry Smith   PetscFunctionBegin;
862a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
86325296bd5SBarry Smith   da->ops->creatematrix = f;
8643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86547c6ae99SBarry Smith }
86647c6ae99SBarry Smith 
86738fb4e8eSJunchao Zhang /*@
868dce8aebaSBarry Smith   DMDAMapMatStencilToGlobal - Map a list of `MatStencil` on a grid to global indices.
86938fb4e8eSJunchao Zhang 
87038fb4e8eSJunchao Zhang   Not Collective
87138fb4e8eSJunchao Zhang 
87238fb4e8eSJunchao Zhang   Input Parameters:
873dce8aebaSBarry Smith + da   - the `DMDA` object
87412b4a537SBarry Smith . m    - number of `MatStencil` to map
87538fb4e8eSJunchao Zhang - idxm - grid points (and component number when dof > 1)
87638fb4e8eSJunchao Zhang 
8771179163eSBarry Smith   Output Parameter:
8781179163eSBarry Smith . gidxm - global row indices
8791179163eSBarry Smith 
8801179163eSBarry Smith   Level: intermediate
88138fb4e8eSJunchao Zhang 
88212b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `MatStencil`
88338fb4e8eSJunchao Zhang @*/
DMDAMapMatStencilToGlobal(DM da,PetscInt m,const MatStencil idxm[],PetscInt gidxm[])884d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAMapMatStencilToGlobal(DM da, PetscInt m, const MatStencil idxm[], PetscInt gidxm[])
885d71ae5a4SJacob Faibussowitsch {
88638fb4e8eSJunchao Zhang   const DM_DA           *dd  = (const DM_DA *)da->data;
88738fb4e8eSJunchao Zhang   const PetscInt        *dxm = (const PetscInt *)idxm;
88838fb4e8eSJunchao Zhang   PetscInt               i, j, sdim, tmp, dim;
88938fb4e8eSJunchao Zhang   PetscInt               dims[4], starts[4], dims2[3], starts2[3], dof = dd->w;
89038fb4e8eSJunchao Zhang   ISLocalToGlobalMapping ltog;
89138fb4e8eSJunchao Zhang 
89238fb4e8eSJunchao Zhang   PetscFunctionBegin;
8933ba16761SJacob Faibussowitsch   if (m <= 0) PetscFunctionReturn(PETSC_SUCCESS);
89438fb4e8eSJunchao Zhang 
8952f27f4d1SJunchao Zhang   /* Code adapted from DMDAGetGhostCorners() */
89638fb4e8eSJunchao Zhang   starts2[0] = dd->Xs / dof + dd->xo;
89738fb4e8eSJunchao Zhang   starts2[1] = dd->Ys + dd->yo;
89838fb4e8eSJunchao Zhang   starts2[2] = dd->Zs + dd->zo;
89938fb4e8eSJunchao Zhang   dims2[0]   = (dd->Xe - dd->Xs) / dof;
90038fb4e8eSJunchao Zhang   dims2[1]   = (dd->Ye - dd->Ys);
90138fb4e8eSJunchao Zhang   dims2[2]   = (dd->Ze - dd->Zs);
90238fb4e8eSJunchao Zhang 
9032f27f4d1SJunchao Zhang   /* As if we do MatSetStencil() to get dims[]/starts[] of mat->stencil */
9042f27f4d1SJunchao Zhang   dim  = da->dim;                   /* DA dim: 1 to 3 */
9052f27f4d1SJunchao Zhang   sdim = dim + (dof > 1 ? 1 : 0);   /* Dimensions in MatStencil's (k,j,i,c) view */
9062f27f4d1SJunchao Zhang   for (i = 0; i < dim; i++) {       /* Reverse the order and also skip the unused dimensions */
9072f27f4d1SJunchao 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 */
90838fb4e8eSJunchao Zhang     starts[i] = starts2[dim - i - 1];
90938fb4e8eSJunchao Zhang   }
9102f27f4d1SJunchao Zhang   starts[dim] = 0; /* Append the extra dim for dof (won't be used below if dof=1) */
91138fb4e8eSJunchao Zhang   dims[dim]   = dof;
91238fb4e8eSJunchao Zhang 
91338fb4e8eSJunchao Zhang   /* Map stencils to local indices (code adapted from MatSetValuesStencil()) */
91438fb4e8eSJunchao Zhang   for (i = 0; i < m; i++) {
9152f27f4d1SJunchao Zhang     dxm += 3 - dim; /* Input is {k,j,i,c}; move the pointer to the first used index, e.g., j in 2D */
9162f27f4d1SJunchao Zhang     tmp = 0;
9172f27f4d1SJunchao Zhang     for (j = 0; j < sdim; j++) {                                                      /* Iter over, ex. j,i or j,i,c in 2D */
9182f27f4d1SJunchao 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 */
9192f27f4d1SJunchao Zhang       else tmp = tmp * dims[j] + (dxm[j] - starts[j]);
92038fb4e8eSJunchao Zhang     }
92138fb4e8eSJunchao Zhang     gidxm[i] = tmp;
9222f27f4d1SJunchao Zhang     /* Move to the next MatStencil point */
9232f27f4d1SJunchao Zhang     if (dof > 1) dxm += sdim; /* c is already counted in sdim */
9242f27f4d1SJunchao Zhang     else dxm += sdim + 1;     /* skip the unused c */
92538fb4e8eSJunchao Zhang   }
92638fb4e8eSJunchao Zhang 
92738fb4e8eSJunchao Zhang   /* Map local indices to global indices */
92838fb4e8eSJunchao Zhang   PetscCall(DMGetLocalToGlobalMapping(da, &ltog));
92938fb4e8eSJunchao Zhang   PetscCall(ISLocalToGlobalMappingApply(ltog, m, gidxm, gidxm));
9303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
93138fb4e8eSJunchao Zhang }
93238fb4e8eSJunchao Zhang 
93347c6ae99SBarry Smith /*
93447c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
93547c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
93647c6ae99SBarry Smith 
93747c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
93847c6ae99SBarry Smith */
DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])939d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDARefineOwnershipRanges(DM da, PetscBool periodic, PetscInt stencil_width, PetscInt ratio, PetscInt m, const PetscInt lc[], PetscInt lf[])
940d71ae5a4SJacob Faibussowitsch {
94147c6ae99SBarry Smith   PetscInt i, totalc = 0, remaining, startc = 0, startf = 0;
94247c6ae99SBarry Smith 
94347c6ae99SBarry Smith   PetscFunctionBegin;
94463a3b9bcSJacob Faibussowitsch   PetscCheck(ratio >= 1, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Requested refinement ratio %" PetscInt_FMT " must be at least 1", ratio);
94547c6ae99SBarry Smith   if (ratio == 1) {
9469566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(lf, lc, m));
9473ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
94847c6ae99SBarry Smith   }
94947c6ae99SBarry Smith   for (i = 0; i < m; i++) totalc += lc[i];
95047c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
95147c6ae99SBarry Smith   for (i = 0; i < m; i++) {
95247c6ae99SBarry Smith     PetscInt want = remaining / (m - i) + !!(remaining % (m - i));
95347c6ae99SBarry Smith     if (i == m - 1) lf[i] = want;
95447c6ae99SBarry Smith     else {
9557aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
9567aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
9577aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
9587aca7175SJed Brown       while ((startf + want) / ratio < nextc - stencil_width) want++;
9597aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
9607aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
9617aca7175SJed Brown       while ((startf + want - 1 + ratio - 1) / ratio > nextc - 1 + stencil_width) want--;
9627aca7175SJed Brown       /* Make sure all constraints are satisfied */
9639371c9d4SSatish Balay       if (want < 0 || want > remaining || ((startf + want) / ratio < nextc - stencil_width) || ((startf + want - 1 + ratio - 1) / ratio > nextc - 1 + stencil_width))
9649371c9d4SSatish Balay         SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Could not find a compatible refined ownership range");
96547c6ae99SBarry Smith     }
96647c6ae99SBarry Smith     lf[i] = want;
96747c6ae99SBarry Smith     startc += lc[i];
96847c6ae99SBarry Smith     startf += lf[i];
96947c6ae99SBarry Smith     remaining -= lf[i];
97047c6ae99SBarry Smith   }
9713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
97247c6ae99SBarry Smith }
97347c6ae99SBarry Smith 
9742be375d4SJed Brown /*
9752be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
9762be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
9772be375d4SJed Brown 
9782be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
9792be375d4SJed Brown */
DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])980d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDACoarsenOwnershipRanges(DM da, PetscBool periodic, PetscInt stencil_width, PetscInt ratio, PetscInt m, const PetscInt lf[], PetscInt lc[])
981d71ae5a4SJacob Faibussowitsch {
9822be375d4SJed Brown   PetscInt i, totalf, remaining, startc, startf;
9832be375d4SJed Brown 
9842be375d4SJed Brown   PetscFunctionBegin;
98563a3b9bcSJacob Faibussowitsch   PetscCheck(ratio >= 1, PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Requested refinement ratio %" PetscInt_FMT " must be at least 1", ratio);
9862be375d4SJed Brown   if (ratio == 1) {
9879566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(lc, lf, m));
9883ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
9892be375d4SJed Brown   }
9902be375d4SJed Brown   for (i = 0, totalf = 0; i < m; i++) totalf += lf[i];
9912be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
9922be375d4SJed Brown   for (i = 0, startc = 0, startf = 0; i < m; i++) {
9932be375d4SJed Brown     PetscInt want = remaining / (m - i) + !!(remaining % (m - i));
9942be375d4SJed Brown     if (i == m - 1) lc[i] = want;
9952be375d4SJed Brown     else {
9962be375d4SJed Brown       const PetscInt nextf = startf + lf[i];
9972be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
9982be375d4SJed Brown        * node is within one stencil width. */
9992be375d4SJed Brown       while (nextf / ratio < startc + want - stencil_width) want--;
10002be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
10012be375d4SJed Brown        * fine node is within one stencil width. */
10022be375d4SJed Brown       while ((nextf - 1 + ratio - 1) / ratio > startc + want - 1 + stencil_width) want++;
10039371c9d4SSatish Balay       if (want < 0 || want > remaining || (nextf / ratio < startc + want - stencil_width) || ((nextf - 1 + ratio - 1) / ratio > startc + want - 1 + stencil_width))
10049371c9d4SSatish Balay         SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_SIZ, "Could not find a compatible coarsened ownership range");
10052be375d4SJed Brown     }
10062be375d4SJed Brown     lc[i] = want;
10072be375d4SJed Brown     startc += lc[i];
10082be375d4SJed Brown     startf += lf[i];
10092be375d4SJed Brown     remaining -= lc[i];
10102be375d4SJed Brown   }
10113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10122be375d4SJed Brown }
10132be375d4SJed Brown 
DMRefine_DA(DM da,MPI_Comm comm,DM * daref)1014d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_DA(DM da, MPI_Comm comm, DM *daref)
1015d71ae5a4SJacob Faibussowitsch {
1016c73cfb54SMatthew G. Knepley   PetscInt M, N, P, i, dim;
10176858538eSMatthew G. Knepley   Vec      coordsc, coordsf;
10189a42bb27SBarry Smith   DM       da2;
101947c6ae99SBarry Smith   DM_DA   *dd = (DM_DA *)da->data, *dd2;
102047c6ae99SBarry Smith 
102147c6ae99SBarry Smith   PetscFunctionBegin;
1022a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
10234f572ea9SToby Isaac   PetscAssertPointer(daref, 3);
102447c6ae99SBarry Smith 
10259566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(da, &dim));
1026bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
102747c6ae99SBarry Smith     M = dd->refine_x * dd->M;
102847c6ae99SBarry Smith   } else {
102947c6ae99SBarry Smith     M = 1 + dd->refine_x * (dd->M - 1);
103047c6ae99SBarry Smith   }
1031bff4a2f0SMatthew G. Knepley   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1032c73cfb54SMatthew G. Knepley     if (dim > 1) {
103347c6ae99SBarry Smith       N = dd->refine_y * dd->N;
103447c6ae99SBarry Smith     } else {
10351860e6e9SBarry Smith       N = 1;
10361860e6e9SBarry Smith     }
10371860e6e9SBarry Smith   } else {
103847c6ae99SBarry Smith     N = 1 + dd->refine_y * (dd->N - 1);
103947c6ae99SBarry Smith   }
1040bff4a2f0SMatthew G. Knepley   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1041c73cfb54SMatthew G. Knepley     if (dim > 2) {
104247c6ae99SBarry Smith       P = dd->refine_z * dd->P;
104347c6ae99SBarry Smith     } else {
10441860e6e9SBarry Smith       P = 1;
10451860e6e9SBarry Smith     }
10461860e6e9SBarry Smith   } else {
104747c6ae99SBarry Smith     P = 1 + dd->refine_z * (dd->P - 1);
104847c6ae99SBarry Smith   }
10499566063dSJacob Faibussowitsch   PetscCall(DMDACreate(PetscObjectComm((PetscObject)da), &da2));
10509566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(da2, ((PetscObject)da)->prefix));
10519566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(da2, dim));
10529566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(da2, M, N, P));
10539566063dSJacob Faibussowitsch   PetscCall(DMDASetNumProcs(da2, dd->m, dd->n, dd->p));
10549566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(da2, dd->bx, dd->by, dd->bz));
10559566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(da2, dd->w));
10569566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilType(da2, dd->stencil_type));
10579566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(da2, dd->s));
1058c73cfb54SMatthew G. Knepley   if (dim == 3) {
105947c6ae99SBarry Smith     PetscInt *lx, *ly, *lz;
10609566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(dd->m, &lx, dd->n, &ly, dd->p, &lz));
10619566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx));
10629566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_y, dd->n, dd->ly, ly));
10639566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_z, dd->p, dd->lz, lz));
10649566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(da2, lx, ly, lz));
10659566063dSJacob Faibussowitsch     PetscCall(PetscFree3(lx, ly, lz));
1066c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
106747c6ae99SBarry Smith     PetscInt *lx, *ly;
10689566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(dd->m, &lx, dd->n, &ly));
10699566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx));
10709566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_y, dd->n, dd->ly, ly));
10719566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(da2, lx, ly, NULL));
10729566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lx, ly));
1073c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
107447c6ae99SBarry Smith     PetscInt *lx;
10759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->m, &lx));
10769566063dSJacob Faibussowitsch     PetscCall(DMDARefineOwnershipRanges(da, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->refine_x, dd->m, dd->lx, lx));
10779566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(da2, lx, NULL, NULL));
10789566063dSJacob Faibussowitsch     PetscCall(PetscFree(lx));
107947c6ae99SBarry Smith   }
108047c6ae99SBarry Smith   dd2 = (DM_DA *)da2->data;
108147c6ae99SBarry Smith 
108247c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
108325296bd5SBarry Smith   da2->ops->creatematrix = da->ops->creatematrix;
108425296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
108547c6ae99SBarry Smith   da2->ops->getcoloring = da->ops->getcoloring;
108647c6ae99SBarry Smith   dd2->interptype       = dd->interptype;
108747c6ae99SBarry Smith 
108847c6ae99SBarry Smith   /* copy fill information if given */
108947c6ae99SBarry Smith   if (dd->dfill) {
10909566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->dfill[dd->w] + dd->w + 1, &dd2->dfill));
10919566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->dfill, dd->dfill, dd->dfill[dd->w] + dd->w + 1));
109247c6ae99SBarry Smith   }
109347c6ae99SBarry Smith   if (dd->ofill) {
10949566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->ofill[dd->w] + dd->w + 1, &dd2->ofill));
10959566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->ofill, dd->ofill, dd->ofill[dd->w] + dd->w + 1));
109647c6ae99SBarry Smith   }
109747c6ae99SBarry Smith   /* copy the refine information */
1098397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
1099397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
1100397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
110147c6ae99SBarry Smith 
1102897f7067SBarry Smith   if (dd->refine_z_hier) {
1103ad540459SPierre 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];
1104ad540459SPierre 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];
1105897f7067SBarry Smith     dd2->refine_z_hier_n = dd->refine_z_hier_n;
11069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_z_hier_n, &dd2->refine_z_hier));
11079566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_z_hier, dd->refine_z_hier, dd2->refine_z_hier_n));
1108897f7067SBarry Smith   }
1109897f7067SBarry Smith   if (dd->refine_y_hier) {
1110ad540459SPierre 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];
1111ad540459SPierre 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];
1112897f7067SBarry Smith     dd2->refine_y_hier_n = dd->refine_y_hier_n;
11139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_y_hier_n, &dd2->refine_y_hier));
11149566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_y_hier, dd->refine_y_hier, dd2->refine_y_hier_n));
1115897f7067SBarry Smith   }
1116897f7067SBarry Smith   if (dd->refine_x_hier) {
1117ad540459SPierre 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];
1118ad540459SPierre 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];
1119897f7067SBarry Smith     dd2->refine_x_hier_n = dd->refine_x_hier_n;
11209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_x_hier_n, &dd2->refine_x_hier));
11219566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_x_hier, dd->refine_x_hier, dd2->refine_x_hier_n));
1122897f7067SBarry Smith   }
1123897f7067SBarry Smith 
112447c6ae99SBarry Smith   /* copy vector type information */
11259566063dSJacob Faibussowitsch   PetscCall(DMSetVecType(da2, da->vectype));
1126ddcf8b74SDave May 
1127efd51863SBarry Smith   dd2->lf = dd->lf;
1128efd51863SBarry Smith   dd2->lj = dd->lj;
1129efd51863SBarry Smith 
11306e87535bSJed Brown   da2->leveldown = da->leveldown;
11316e87535bSJed Brown   da2->levelup   = da->levelup + 1;
11328865f1eaSKarl Rupp 
11339566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da2));
11346e87535bSJed Brown 
1135ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
11366858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(da, &coordsc));
11376858538eSMatthew G. Knepley   if (coordsc) {
1138ddcf8b74SDave May     DM  cdaf, cdac;
1139ddcf8b74SDave May     Mat II;
1140ddcf8b74SDave May 
11419566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(da, &cdac));
11429566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(da2, &cdaf));
1143b61d3410SDave May     /* force creation of the coordinate vector */
11449566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(da2, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
11459566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(da2, &coordsf));
11469566063dSJacob Faibussowitsch     PetscCall(DMCreateInterpolation(cdac, cdaf, &II, NULL));
11479566063dSJacob Faibussowitsch     PetscCall(MatInterpolate(II, coordsc, coordsf));
11489566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&II));
1149ddcf8b74SDave May   }
1150397b6216SJed Brown 
1151f3141302SJed Brown   for (i = 0; i < da->bs; i++) {
1152f3141302SJed Brown     const char *fieldname;
11539566063dSJacob Faibussowitsch     PetscCall(DMDAGetFieldName(da, i, &fieldname));
11549566063dSJacob Faibussowitsch     PetscCall(DMDASetFieldName(da2, i, fieldname));
1155f3141302SJed Brown   }
1156397b6216SJed Brown 
115747c6ae99SBarry Smith   *daref = da2;
11583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115947c6ae99SBarry Smith }
116047c6ae99SBarry Smith 
DMCoarsen_DA(DM dmf,MPI_Comm comm,DM * dmc)1161d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsen_DA(DM dmf, MPI_Comm comm, DM *dmc)
1162d71ae5a4SJacob Faibussowitsch {
1163c73cfb54SMatthew G. Knepley   PetscInt M, N, P, i, dim;
11646858538eSMatthew G. Knepley   Vec      coordsc, coordsf;
1165a5bc1bf3SBarry Smith   DM       dmc2;
1166a5bc1bf3SBarry Smith   DM_DA   *dd = (DM_DA *)dmf->data, *dd2;
116747c6ae99SBarry Smith 
116847c6ae99SBarry Smith   PetscFunctionBegin;
1169a5bc1bf3SBarry Smith   PetscValidHeaderSpecificType(dmf, DM_CLASSID, 1, DMDA);
11704f572ea9SToby Isaac   PetscAssertPointer(dmc, 3);
117147c6ae99SBarry Smith 
11729566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dmf, &dim));
1173bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1174397b6216SJed Brown     M = dd->M / dd->coarsen_x;
117547c6ae99SBarry Smith   } else {
1176397b6216SJed Brown     M = 1 + (dd->M - 1) / dd->coarsen_x;
117747c6ae99SBarry Smith   }
1178bff4a2f0SMatthew G. Knepley   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1179c73cfb54SMatthew G. Knepley     if (dim > 1) {
1180397b6216SJed Brown       N = dd->N / dd->coarsen_y;
118147c6ae99SBarry Smith     } else {
11821860e6e9SBarry Smith       N = 1;
11831860e6e9SBarry Smith     }
11841860e6e9SBarry Smith   } else {
1185397b6216SJed Brown     N = 1 + (dd->N - 1) / dd->coarsen_y;
118647c6ae99SBarry Smith   }
1187bff4a2f0SMatthew G. Knepley   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1188c73cfb54SMatthew G. Knepley     if (dim > 2) {
1189397b6216SJed Brown       P = dd->P / dd->coarsen_z;
119047c6ae99SBarry Smith     } else {
11911860e6e9SBarry Smith       P = 1;
11921860e6e9SBarry Smith     }
11931860e6e9SBarry Smith   } else {
1194397b6216SJed Brown     P = 1 + (dd->P - 1) / dd->coarsen_z;
119547c6ae99SBarry Smith   }
11969566063dSJacob Faibussowitsch   PetscCall(DMDACreate(PetscObjectComm((PetscObject)dmf), &dmc2));
11979566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(dmc2, ((PetscObject)dmf)->prefix));
11989566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dmc2, dim));
11999566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(dmc2, M, N, P));
12009566063dSJacob Faibussowitsch   PetscCall(DMDASetNumProcs(dmc2, dd->m, dd->n, dd->p));
12019566063dSJacob Faibussowitsch   PetscCall(DMDASetBoundaryType(dmc2, dd->bx, dd->by, dd->bz));
12029566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(dmc2, dd->w));
12039566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilType(dmc2, dd->stencil_type));
12049566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(dmc2, dd->s));
1205c73cfb54SMatthew G. Knepley   if (dim == 3) {
12062be375d4SJed Brown     PetscInt *lx, *ly, *lz;
12079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(dd->m, &lx, dd->n, &ly, dd->p, &lz));
12089566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx));
12099566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_y, dd->n, dd->ly, ly));
12109566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_z, dd->p, dd->lz, lz));
12119566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(dmc2, lx, ly, lz));
12129566063dSJacob Faibussowitsch     PetscCall(PetscFree3(lx, ly, lz));
1213c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
12142be375d4SJed Brown     PetscInt *lx, *ly;
12159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(dd->m, &lx, dd->n, &ly));
12169566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx));
12179566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_y, dd->n, dd->ly, ly));
12189566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(dmc2, lx, ly, NULL));
12199566063dSJacob Faibussowitsch     PetscCall(PetscFree2(lx, ly));
1220c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
12212be375d4SJed Brown     PetscInt *lx;
12229566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->m, &lx));
12239566063dSJacob Faibussowitsch     PetscCall(DMDACoarsenOwnershipRanges(dmf, (PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0), dd->s, dd->coarsen_x, dd->m, dd->lx, lx));
12249566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(dmc2, lx, NULL, NULL));
12259566063dSJacob Faibussowitsch     PetscCall(PetscFree(lx));
122647c6ae99SBarry Smith   }
1227a5bc1bf3SBarry Smith   dd2 = (DM_DA *)dmc2->data;
122847c6ae99SBarry Smith 
12294dcab191SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
1230a5bc1bf3SBarry Smith   /* dmc2->ops->createinterpolation = dmf->ops->createinterpolation; copying this one causes trouble for DMSetVI */
1231a5bc1bf3SBarry Smith   dmc2->ops->creatematrix = dmf->ops->creatematrix;
1232a5bc1bf3SBarry Smith   dmc2->ops->getcoloring  = dmf->ops->getcoloring;
123347c6ae99SBarry Smith   dd2->interptype         = dd->interptype;
123447c6ae99SBarry Smith 
123547c6ae99SBarry Smith   /* copy fill information if given */
123647c6ae99SBarry Smith   if (dd->dfill) {
12379566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->dfill[dd->w] + dd->w + 1, &dd2->dfill));
12389566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->dfill, dd->dfill, dd->dfill[dd->w] + dd->w + 1));
123947c6ae99SBarry Smith   }
124047c6ae99SBarry Smith   if (dd->ofill) {
12419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->ofill[dd->w] + dd->w + 1, &dd2->ofill));
12429566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->ofill, dd->ofill, dd->ofill[dd->w] + dd->w + 1));
124347c6ae99SBarry Smith   }
124447c6ae99SBarry Smith   /* copy the refine information */
1245397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
1246397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
1247397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
124847c6ae99SBarry Smith 
1249897f7067SBarry Smith   if (dd->refine_z_hier) {
1250ad540459SPierre 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];
1251ad540459SPierre 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];
1252897f7067SBarry Smith     dd2->refine_z_hier_n = dd->refine_z_hier_n;
12539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_z_hier_n, &dd2->refine_z_hier));
12549566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_z_hier, dd->refine_z_hier, dd2->refine_z_hier_n));
1255897f7067SBarry Smith   }
1256897f7067SBarry Smith   if (dd->refine_y_hier) {
1257ad540459SPierre 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];
1258ad540459SPierre 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];
1259897f7067SBarry Smith     dd2->refine_y_hier_n = dd->refine_y_hier_n;
12609566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_y_hier_n, &dd2->refine_y_hier));
12619566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_y_hier, dd->refine_y_hier, dd2->refine_y_hier_n));
1262897f7067SBarry Smith   }
1263897f7067SBarry Smith   if (dd->refine_x_hier) {
1264ad540459SPierre 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];
1265ad540459SPierre 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];
1266897f7067SBarry Smith     dd2->refine_x_hier_n = dd->refine_x_hier_n;
12679566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd2->refine_x_hier_n, &dd2->refine_x_hier));
12689566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(dd2->refine_x_hier, dd->refine_x_hier, dd2->refine_x_hier_n));
1269897f7067SBarry Smith   }
1270897f7067SBarry Smith 
127147c6ae99SBarry Smith   /* copy vector type information */
12729566063dSJacob Faibussowitsch   PetscCall(DMSetVecType(dmc2, dmf->vectype));
127347c6ae99SBarry Smith 
1274644e2e5bSBarry Smith   dd2->lf = dd->lf;
1275644e2e5bSBarry Smith   dd2->lj = dd->lj;
1276644e2e5bSBarry Smith 
1277a5bc1bf3SBarry Smith   dmc2->leveldown = dmf->leveldown + 1;
1278a5bc1bf3SBarry Smith   dmc2->levelup   = dmf->levelup;
12798865f1eaSKarl Rupp 
12809566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmc2));
12816e87535bSJed Brown 
1282ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
12836858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dmf, &coordsf));
12846858538eSMatthew G. Knepley   if (coordsf) {
1285ddcf8b74SDave May     DM         cdaf, cdac;
12866dbf9973SLawrence Mitchell     Mat        inject;
12876dbf9973SLawrence Mitchell     VecScatter vscat;
1288ddcf8b74SDave May 
12899566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmf, &cdaf));
12909566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmc2, &cdac));
1291b61d3410SDave May     /* force creation of the coordinate vector */
12929566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(dmc2, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
12939566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmc2, &coordsc));
1294ddcf8b74SDave May 
12959566063dSJacob Faibussowitsch     PetscCall(DMCreateInjection(cdac, cdaf, &inject));
12969566063dSJacob Faibussowitsch     PetscCall(MatScatterGetVecScatter(inject, &vscat));
12979566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(vscat, coordsf, coordsc, INSERT_VALUES, SCATTER_FORWARD));
12989566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(vscat, coordsf, coordsc, INSERT_VALUES, SCATTER_FORWARD));
12999566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&inject));
1300ddcf8b74SDave May   }
1301f98405f7SJed Brown 
1302a5bc1bf3SBarry Smith   for (i = 0; i < dmf->bs; i++) {
1303f98405f7SJed Brown     const char *fieldname;
13049566063dSJacob Faibussowitsch     PetscCall(DMDAGetFieldName(dmf, i, &fieldname));
13059566063dSJacob Faibussowitsch     PetscCall(DMDASetFieldName(dmc2, i, fieldname));
1306f98405f7SJed Brown   }
1307f98405f7SJed Brown 
1308a5bc1bf3SBarry Smith   *dmc = dmc2;
13093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131047c6ae99SBarry Smith }
131147c6ae99SBarry Smith 
DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])1312d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefineHierarchy_DA(DM da, PetscInt nlevels, DM daf[])
1313d71ae5a4SJacob Faibussowitsch {
131447c6ae99SBarry Smith   PetscInt i, n, *refx, *refy, *refz;
131547c6ae99SBarry Smith 
131647c6ae99SBarry Smith   PetscFunctionBegin;
131747c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
13187a8be351SBarry Smith   PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative");
13193ba16761SJacob Faibussowitsch   if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS);
13204f572ea9SToby Isaac   PetscAssertPointer(daf, 3);
132147c6ae99SBarry Smith 
1322aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
13239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(nlevels, &refx, nlevels, &refy, nlevels, &refz));
132448a46eb9SPierre Jolivet   for (i = 0; i < nlevels; i++) PetscCall(DMDAGetRefinementFactor(da, &refx[i], &refy[i], &refz[i]));
132547c6ae99SBarry Smith   n = nlevels;
13269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_x", refx, &n, NULL));
132747c6ae99SBarry Smith   n = nlevels;
13289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_y", refy, &n, NULL));
132947c6ae99SBarry Smith   n = nlevels;
13309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_refine_hierarchy_z", refz, &n, NULL));
133147c6ae99SBarry Smith 
13329566063dSJacob Faibussowitsch   PetscCall(DMDASetRefinementFactor(da, refx[0], refy[0], refz[0]));
13339566063dSJacob Faibussowitsch   PetscCall(DMRefine(da, PetscObjectComm((PetscObject)da), &daf[0]));
133447c6ae99SBarry Smith   for (i = 1; i < nlevels; i++) {
13359566063dSJacob Faibussowitsch     PetscCall(DMDASetRefinementFactor(daf[i - 1], refx[i], refy[i], refz[i]));
13369566063dSJacob Faibussowitsch     PetscCall(DMRefine(daf[i - 1], PetscObjectComm((PetscObject)da), &daf[i]));
133747c6ae99SBarry Smith   }
13389566063dSJacob Faibussowitsch   PetscCall(PetscFree3(refx, refy, refz));
13393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
134047c6ae99SBarry Smith }
134147c6ae99SBarry Smith 
DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])1342d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsenHierarchy_DA(DM da, PetscInt nlevels, DM dac[])
1343d71ae5a4SJacob Faibussowitsch {
134447c6ae99SBarry Smith   PetscInt i;
134547c6ae99SBarry Smith 
134647c6ae99SBarry Smith   PetscFunctionBegin;
134747c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
13487a8be351SBarry Smith   PetscCheck(nlevels >= 0, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "nlevels cannot be negative");
13493ba16761SJacob Faibussowitsch   if (nlevels == 0) PetscFunctionReturn(PETSC_SUCCESS);
13504f572ea9SToby Isaac   PetscAssertPointer(dac, 3);
13519566063dSJacob Faibussowitsch   PetscCall(DMCoarsen(da, PetscObjectComm((PetscObject)da), &dac[0]));
135248a46eb9SPierre Jolivet   for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dac[i - 1], PetscObjectComm((PetscObject)da), &dac[i]));
13533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
135447c6ae99SBarry Smith }
135562197512SBarry Smith 
DMDASetGLLCoordinates_1d(DM dm,PetscInt n,PetscReal * nodes)135666976f2fSJacob Faibussowitsch static PetscErrorCode DMDASetGLLCoordinates_1d(DM dm, PetscInt n, PetscReal *nodes)
1357d71ae5a4SJacob Faibussowitsch {
13588272889dSSatish Balay   PetscInt     i, j, xs, xn, q;
135962197512SBarry Smith   PetscScalar *xx;
136062197512SBarry Smith   PetscReal    h;
136162197512SBarry Smith   Vec          x;
136262197512SBarry Smith   DM_DA       *da = (DM_DA *)dm->data;
136362197512SBarry Smith 
136462197512SBarry Smith   PetscFunctionBegin;
1365966bd95aSPierre Jolivet   PetscCheck(da->bx != DM_BOUNDARY_PERIODIC, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Not yet implemented for periodic");
13669566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, NULL, &q, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
136762197512SBarry Smith   q = (q - 1) / (n - 1); /* number of spectral elements */
136862197512SBarry Smith   h = 2.0 / q;
13699566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xn, NULL, NULL));
137062197512SBarry Smith   xs = xs / (n - 1);
137162197512SBarry Smith   xn = xn / (n - 1);
13729566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dm, -1., 1., 0., 0., 0., 0.));
13739566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(dm, &x));
13749566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, x, &xx));
137562197512SBarry Smith 
137662197512SBarry Smith   /* loop over local spectral elements */
137762197512SBarry Smith   for (j = xs; j < xs + xn; j++) {
137862197512SBarry Smith     /*
137962197512SBarry Smith      Except for the first process, each process starts on the second GLL point of the first element on that process
138062197512SBarry Smith      */
1381ad540459SPierre 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.;
138262197512SBarry Smith   }
13839566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, x, &xx));
13843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
138562197512SBarry Smith }
138662197512SBarry Smith 
13871fd49c25SBarry Smith /*@
138862197512SBarry Smith 
138962197512SBarry 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
139062197512SBarry Smith 
139120f4b53cSBarry Smith   Collective
139262197512SBarry Smith 
139362197512SBarry Smith   Input Parameters:
1394dce8aebaSBarry Smith + da    - the `DMDA` object
13952fe279fdSBarry Smith . n     - the number of GLL nodes
13968272889dSSatish Balay - nodes - the GLL nodes
139762197512SBarry Smith 
1398edc382c3SSatish Balay   Level: advanced
1399edc382c3SSatish Balay 
1400dce8aebaSBarry Smith   Note:
1401dce8aebaSBarry Smith   The parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points
140212b4a537SBarry Smith   on each process much be divisible by the number of GLL elements needed per process. This depends on whether the `DMDA` is
1403dce8aebaSBarry Smith   periodic or not.
1404dce8aebaSBarry Smith 
140512b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDACreate()`, `PetscDTGaussLobattoLegendreQuadrature()`, `DMGetCoordinates()`
140662197512SBarry Smith @*/
DMDASetGLLCoordinates(DM da,PetscInt n,PetscReal * nodes)1407d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetGLLCoordinates(DM da, PetscInt n, PetscReal *nodes)
1408d71ae5a4SJacob Faibussowitsch {
140962197512SBarry Smith   PetscFunctionBegin;
1410966bd95aSPierre Jolivet   PetscCheck(da->dim == 1, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Not yet implemented for 2 or 3d");
14119566063dSJacob Faibussowitsch   PetscCall(DMDASetGLLCoordinates_1d(da, n, nodes));
14123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
141362197512SBarry Smith }
14147c3cd84eSPatrick Sanan 
DMGetCompatibility_DA(DM da1,DM dm2,PetscBool * compatible,PetscBool * set)141566976f2fSJacob Faibussowitsch PetscErrorCode DMGetCompatibility_DA(DM da1, DM dm2, PetscBool *compatible, PetscBool *set)
1416d71ae5a4SJacob Faibussowitsch {
14177c3cd84eSPatrick Sanan   DM_DA    *dd1 = (DM_DA *)da1->data, *dd2;
14187c3cd84eSPatrick Sanan   DM        da2;
14197c3cd84eSPatrick Sanan   DMType    dmtype2;
14207c3cd84eSPatrick Sanan   PetscBool isda, compatibleLocal;
14217c3cd84eSPatrick Sanan   PetscInt  i;
14227c3cd84eSPatrick Sanan 
14237c3cd84eSPatrick Sanan   PetscFunctionBegin;
14247a8be351SBarry Smith   PetscCheck(da1->setupcalled, PetscObjectComm((PetscObject)da1), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() must be called on first DM before DMGetCompatibility()");
14259566063dSJacob Faibussowitsch   PetscCall(DMGetType(dm2, &dmtype2));
14269566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(dmtype2, DMDA, &isda));
14277c3cd84eSPatrick Sanan   if (isda) {
14287c3cd84eSPatrick Sanan     da2 = dm2;
14297c3cd84eSPatrick Sanan     dd2 = (DM_DA *)da2->data;
14307a8be351SBarry Smith     PetscCheck(da2->setupcalled, PetscObjectComm((PetscObject)da2), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() must be called on second DM before DMGetCompatibility()");
1431110623a0SKarl Rupp     compatibleLocal = (PetscBool)(da1->dim == da2->dim);
1432c790a739SKarl Rupp     if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->s == dd2->s)); /* Stencil width */
1433110623a0SKarl Rupp     /*                                                                           Global size              ranks               Boundary type */
1434c790a739SKarl Rupp     if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->M == dd2->M) && (dd1->m == dd2->m) && (dd1->bx == dd2->bx));
1435c790a739SKarl Rupp     if (compatibleLocal && da1->dim > 1) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->N == dd2->N) && (dd1->n == dd2->n) && (dd1->by == dd2->by));
1436c790a739SKarl Rupp     if (compatibleLocal && da1->dim > 2) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->P == dd2->P) && (dd1->p == dd2->p) && (dd1->bz == dd2->bz));
14377c3cd84eSPatrick Sanan     if (compatibleLocal) {
1438*ac530a7eSPierre Jolivet       for (i = 0; i < dd1->m; ++i) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lx[i] == dd2->lx[i])); /* Local size     */
14397c3cd84eSPatrick Sanan     }
14407c3cd84eSPatrick Sanan     if (compatibleLocal && da1->dim > 1) {
1441ad540459SPierre Jolivet       for (i = 0; i < dd1->n; ++i) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->ly[i] == dd2->ly[i]));
14427c3cd84eSPatrick Sanan     }
14437c3cd84eSPatrick Sanan     if (compatibleLocal && da1->dim > 2) {
1444ad540459SPierre Jolivet       for (i = 0; i < dd1->p; ++i) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lz[i] == dd2->lz[i]));
14457c3cd84eSPatrick Sanan     }
14467c3cd84eSPatrick Sanan     *compatible = compatibleLocal;
14477c3cd84eSPatrick Sanan     *set        = PETSC_TRUE;
14487c3cd84eSPatrick Sanan   } else {
14497c3cd84eSPatrick Sanan     /* Decline to determine compatibility with other DM types */
14507c3cd84eSPatrick Sanan     *set = PETSC_FALSE;
14517c3cd84eSPatrick Sanan   }
14523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14537c3cd84eSPatrick Sanan }
1454