xref: /petsc/src/dm/impls/da/da.c (revision c73cfb54f34510b9dcefc99e398efa3b88d5dde5)
14035e84dSBarry Smith #include <petsc-private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
247c6ae99SBarry Smith 
347c6ae99SBarry Smith #undef __FUNCT__
4aa219208SBarry Smith #define __FUNCT__ "DMDASetSizes"
547c6ae99SBarry Smith /*@
6aa219208SBarry Smith   DMDASetSizes - Sets the global sizes
747c6ae99SBarry Smith 
8aa219208SBarry Smith   Logically Collective on DMDA
947c6ae99SBarry Smith 
1047c6ae99SBarry Smith   Input Parameters:
11aa219208SBarry Smith + da - the DMDA
1247c6ae99SBarry Smith . M - the global X size (or PETSC_DECIDE)
1347c6ae99SBarry Smith . N - the global Y size (or PETSC_DECIDE)
1447c6ae99SBarry Smith - P - the global Z size (or PETSC_DECIDE)
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith   Level: intermediate
1747c6ae99SBarry Smith 
18aa219208SBarry Smith .seealso: DMDAGetSize(), PetscSplitOwnership()
1947c6ae99SBarry Smith @*/
207087cfbeSBarry Smith PetscErrorCode  DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P)
2147c6ae99SBarry Smith {
2247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   PetscFunctionBegin;
2547c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
2647c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,M,2);
2747c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,N,3);
2847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,P,4);
29ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
3047c6ae99SBarry Smith 
3147c6ae99SBarry Smith   dd->M = M;
3247c6ae99SBarry Smith   dd->N = N;
3347c6ae99SBarry Smith   dd->P = P;
3447c6ae99SBarry Smith   PetscFunctionReturn(0);
3547c6ae99SBarry Smith }
3647c6ae99SBarry Smith 
3747c6ae99SBarry Smith #undef __FUNCT__
38aa219208SBarry Smith #define __FUNCT__ "DMDASetNumProcs"
3947c6ae99SBarry Smith /*@
40aa219208SBarry Smith   DMDASetNumProcs - Sets the number of processes in each dimension
4147c6ae99SBarry Smith 
42aa219208SBarry Smith   Logically Collective on DMDA
4347c6ae99SBarry Smith 
4447c6ae99SBarry Smith   Input Parameters:
45aa219208SBarry Smith + da - the DMDA
4647c6ae99SBarry Smith . m - the number of X procs (or PETSC_DECIDE)
4747c6ae99SBarry Smith . n - the number of Y procs (or PETSC_DECIDE)
4847c6ae99SBarry Smith - p - the number of Z procs (or PETSC_DECIDE)
4947c6ae99SBarry Smith 
5047c6ae99SBarry Smith   Level: intermediate
5147c6ae99SBarry Smith 
52aa219208SBarry Smith .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership()
5347c6ae99SBarry Smith @*/
547087cfbeSBarry Smith PetscErrorCode  DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p)
5547c6ae99SBarry Smith {
5647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
57e3f3e4b6SBarry Smith   PetscErrorCode ierr;
5847c6ae99SBarry Smith 
5947c6ae99SBarry Smith   PetscFunctionBegin;
6047c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
6147c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,m,2);
6247c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,n,3);
6347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,p,4);
64ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
6547c6ae99SBarry Smith   dd->m = m;
6647c6ae99SBarry Smith   dd->n = n;
6747c6ae99SBarry Smith   dd->p = p;
68*c73cfb54SMatthew G. Knepley   if (da->dim == 2) {
69e3f3e4b6SBarry Smith     PetscMPIInt size;
70ce94432eSBarry Smith     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr);
71e3f3e4b6SBarry Smith     if ((dd->m > 0) && (dd->n < 0)) {
72e3f3e4b6SBarry Smith       dd->n = size/dd->m;
73ce94432eSBarry Smith       if (dd->n*dd->m != size) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"%D processes in X direction not divisible into comm size %d",m,size);
74e3f3e4b6SBarry Smith     }
75e3f3e4b6SBarry Smith     if ((dd->n > 0) && (dd->m < 0)) {
76e3f3e4b6SBarry Smith       dd->m = size/dd->n;
77ce94432eSBarry Smith       if (dd->n*dd->m != size) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"%D processes in Y direction not divisible into comm size %d",n,size);
78e3f3e4b6SBarry Smith     }
79e3f3e4b6SBarry Smith   }
8047c6ae99SBarry Smith   PetscFunctionReturn(0);
8147c6ae99SBarry Smith }
8247c6ae99SBarry Smith 
8347c6ae99SBarry Smith #undef __FUNCT__
843e24977eSEthan Coon #define __FUNCT__ "DMDASetBoundaryType"
8547c6ae99SBarry Smith /*@
861321219cSEthan Coon   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries.
8747c6ae99SBarry Smith 
8847c6ae99SBarry Smith   Not collective
8947c6ae99SBarry Smith 
9047c6ae99SBarry Smith   Input Parameter:
91aa219208SBarry Smith + da    - The DMDA
92bff4a2f0SMatthew G. Knepley - bx,by,bz - One of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC
9347c6ae99SBarry Smith 
9447c6ae99SBarry Smith   Level: intermediate
9547c6ae99SBarry Smith 
9647c6ae99SBarry Smith .keywords:  distributed array, periodicity
97bff4a2f0SMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA, DMBoundaryType
9847c6ae99SBarry Smith @*/
99bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDASetBoundaryType(DM da,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz)
10047c6ae99SBarry Smith {
10147c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
10247c6ae99SBarry Smith 
10347c6ae99SBarry Smith   PetscFunctionBegin;
10447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1055a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,bx,2);
1065a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,by,3);
1075a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,bz,4);
108ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1091321219cSEthan Coon   dd->bx = bx;
1101321219cSEthan Coon   dd->by = by;
1111321219cSEthan Coon   dd->bz = bz;
11247c6ae99SBarry Smith   PetscFunctionReturn(0);
11347c6ae99SBarry Smith }
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith #undef __FUNCT__
116aa219208SBarry Smith #define __FUNCT__ "DMDASetDof"
11747c6ae99SBarry Smith /*@
118aa219208SBarry Smith   DMDASetDof - Sets the number of degrees of freedom per vertex
11947c6ae99SBarry Smith 
12047c6ae99SBarry Smith   Not collective
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith   Input Parameter:
123aa219208SBarry Smith + da  - The DMDA
12447c6ae99SBarry Smith - dof - Number of degrees of freedom
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith   Level: intermediate
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith .keywords:  distributed array, degrees of freedom
12996e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
13047c6ae99SBarry Smith @*/
13154cfb0beSLisandro Dalcin PetscErrorCode  DMDASetDof(DM da, PetscInt dof)
13247c6ae99SBarry Smith {
13347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
13447c6ae99SBarry Smith 
13547c6ae99SBarry Smith   PetscFunctionBegin;
13647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
13754cfb0beSLisandro Dalcin   PetscValidLogicalCollectiveInt(da,dof,2);
138ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
13947c6ae99SBarry Smith   dd->w  = dof;
1401411c6eeSJed Brown   da->bs = dof;
14147c6ae99SBarry Smith   PetscFunctionReturn(0);
14247c6ae99SBarry Smith }
14347c6ae99SBarry Smith 
14447c6ae99SBarry Smith #undef __FUNCT__
1457ddda789SPeter Brune #define __FUNCT__ "DMDAGetOverlap"
1467ddda789SPeter Brune /*@
1477ddda789SPeter Brune   DMDAGetOverlap - Gets the size of the per-processor overlap.
1487ddda789SPeter Brune 
1497ddda789SPeter Brune   Not collective
1507ddda789SPeter Brune 
1517ddda789SPeter Brune   Input Parameters:
1527ddda789SPeter Brune . da  - The DMDA
1537ddda789SPeter Brune 
1547ddda789SPeter Brune   Output Parameters:
1557ddda789SPeter Brune + x   - Overlap in the x direction
1567ddda789SPeter Brune . y   - Overlap in the y direction
1577ddda789SPeter Brune - z   - Overlap in the z direction
1587ddda789SPeter Brune 
1597ddda789SPeter Brune   Level: intermediate
1607ddda789SPeter Brune 
1617ddda789SPeter Brune .keywords:  distributed array, overlap, domain decomposition
1627ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetOverlap(), DMDA
1637ddda789SPeter Brune @*/
1647ddda789SPeter Brune PetscErrorCode  DMDAGetOverlap(DM da,PetscInt *x,PetscInt *y,PetscInt *z)
1657ddda789SPeter Brune {
1667ddda789SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
1677ddda789SPeter Brune 
1687ddda789SPeter Brune   PetscFunctionBegin;
1697ddda789SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1707ddda789SPeter Brune   if (x) *x = dd->xol;
1717ddda789SPeter Brune   if (y) *y = dd->yol;
1727ddda789SPeter Brune   if (z) *z = dd->zol;
1737ddda789SPeter Brune   PetscFunctionReturn(0);
1747ddda789SPeter Brune }
1757ddda789SPeter Brune 
1767ddda789SPeter Brune #undef __FUNCT__
17788661749SPeter Brune #define __FUNCT__ "DMDASetOverlap"
17888661749SPeter Brune /*@
17988661749SPeter Brune   DMDASetOverlap - Sets the size of the per-processor overlap.
18088661749SPeter Brune 
18188661749SPeter Brune   Not collective
18288661749SPeter Brune 
1837ddda789SPeter Brune   Input Parameters:
18488661749SPeter Brune + da  - The DMDA
1857ddda789SPeter Brune . x   - Overlap in the x direction
1867ddda789SPeter Brune . y   - Overlap in the y direction
1877ddda789SPeter Brune - z   - Overlap in the z direction
18888661749SPeter Brune 
18988661749SPeter Brune   Level: intermediate
19088661749SPeter Brune 
1917ddda789SPeter Brune .keywords:  distributed array, overlap, domain decomposition
1927ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetOverlap(), DMDA
19388661749SPeter Brune @*/
1947ddda789SPeter Brune PetscErrorCode  DMDASetOverlap(DM da,PetscInt x,PetscInt y,PetscInt z)
19588661749SPeter Brune {
19688661749SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
19788661749SPeter Brune 
19888661749SPeter Brune   PetscFunctionBegin;
19988661749SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
2007ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da,x,2);
2017ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da,y,3);
2027ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da,z,4);
2037ddda789SPeter Brune   dd->xol = x;
2047ddda789SPeter Brune   dd->yol = y;
2057ddda789SPeter Brune   dd->zol = z;
20688661749SPeter Brune   PetscFunctionReturn(0);
20788661749SPeter Brune }
20888661749SPeter Brune 
2093e7870d2SPeter Brune 
2103e7870d2SPeter Brune #undef __FUNCT__
2113e7870d2SPeter Brune #define __FUNCT__ "DMDAGetNumLocalSubDomains"
2123e7870d2SPeter Brune /*@
2133e7870d2SPeter Brune   DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition.
2143e7870d2SPeter Brune 
2153e7870d2SPeter Brune   Not collective
2163e7870d2SPeter Brune 
2173e7870d2SPeter Brune   Input Parameters:
2183e7870d2SPeter Brune . da  - The DMDA
2193e7870d2SPeter Brune 
2203e7870d2SPeter Brune   Output Parameters:
2213e7870d2SPeter Brune + Nsub   - Number of local subdomains created upon decomposition
2223e7870d2SPeter Brune 
2233e7870d2SPeter Brune   Level: intermediate
2243e7870d2SPeter Brune 
2253e7870d2SPeter Brune .keywords:  distributed array, domain decomposition
2263e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetNumLocalSubDomains(), DMDA
2273e7870d2SPeter Brune @*/
2283e7870d2SPeter Brune PetscErrorCode  DMDAGetNumLocalSubDomains(DM da,PetscInt *Nsub)
2293e7870d2SPeter Brune {
2303e7870d2SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
2313e7870d2SPeter Brune 
2323e7870d2SPeter Brune   PetscFunctionBegin;
2333e7870d2SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
2343e7870d2SPeter Brune   if (Nsub) *Nsub = dd->Nsub;
2353e7870d2SPeter Brune   PetscFunctionReturn(0);
2363e7870d2SPeter Brune }
2373e7870d2SPeter Brune 
2383e7870d2SPeter Brune #undef __FUNCT__
2393e7870d2SPeter Brune #define __FUNCT__ "DMDASetNumLocalSubDomains"
2403e7870d2SPeter Brune /*@
2413e7870d2SPeter Brune   DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition.
2423e7870d2SPeter Brune 
2433e7870d2SPeter Brune   Not collective
2443e7870d2SPeter Brune 
2453e7870d2SPeter Brune   Input Parameters:
2463e7870d2SPeter Brune + da  - The DMDA
2473e7870d2SPeter Brune - Nsub - The number of local subdomains requested
2483e7870d2SPeter Brune 
2493e7870d2SPeter Brune   Level: intermediate
2503e7870d2SPeter Brune 
2513e7870d2SPeter Brune .keywords:  distributed array, domain decomposition
2523e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetNumLocalSubDomains(), DMDA
2533e7870d2SPeter Brune @*/
2543e7870d2SPeter Brune PetscErrorCode  DMDASetNumLocalSubDomains(DM da,PetscInt Nsub)
2553e7870d2SPeter Brune {
2563e7870d2SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
2573e7870d2SPeter Brune 
2583e7870d2SPeter Brune   PetscFunctionBegin;
2593e7870d2SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
2603e7870d2SPeter Brune   PetscValidLogicalCollectiveInt(da,Nsub,2);
2613e7870d2SPeter Brune   dd->Nsub = Nsub;
2623e7870d2SPeter Brune   PetscFunctionReturn(0);
2633e7870d2SPeter Brune }
2643e7870d2SPeter Brune 
265d886c4f4SPeter Brune #undef __FUNCT__
266d886c4f4SPeter Brune #define __FUNCT__ "DMDASetOffset"
267d886c4f4SPeter Brune /*@
268d886c4f4SPeter Brune   DMDASetOffset - Sets the index offset of the DA.
269d886c4f4SPeter Brune 
270d886c4f4SPeter Brune   Collective on DA
271d886c4f4SPeter Brune 
272d886c4f4SPeter Brune   Input Parameter:
273d886c4f4SPeter Brune + da  - The DMDA
274d886c4f4SPeter Brune . xo  - The offset in the x direction
275d886c4f4SPeter Brune . yo  - The offset in the y direction
276d886c4f4SPeter Brune - zo  - The offset in the z direction
277d886c4f4SPeter Brune 
278d886c4f4SPeter Brune   Level: intermediate
279d886c4f4SPeter Brune 
280d886c4f4SPeter Brune   Notes: This is used primarily to overlap a computation on a local DA with that on a global DA without
281d886c4f4SPeter Brune   changing boundary conditions or subdomain features that depend upon the global offsets.
282d886c4f4SPeter Brune 
283d886c4f4SPeter Brune .keywords:  distributed array, degrees of freedom
284d886c4f4SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray()
285d886c4f4SPeter Brune @*/
28695c13181SPeter Brune PetscErrorCode  DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po)
287d886c4f4SPeter Brune {
28895c13181SPeter Brune   PetscErrorCode ierr;
289d886c4f4SPeter Brune   DM_DA          *dd = (DM_DA*)da->data;
290d886c4f4SPeter Brune 
291d886c4f4SPeter Brune   PetscFunctionBegin;
292d886c4f4SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
293d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da,xo,2);
29495c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,yo,3);
29595c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,zo,4);
29695c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,Mo,5);
29795c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,No,6);
29895c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,Po,7);
299d886c4f4SPeter Brune   dd->xo = xo;
300d886c4f4SPeter Brune   dd->yo = yo;
301d886c4f4SPeter Brune   dd->zo = zo;
30295c13181SPeter Brune   dd->Mo = Mo;
30395c13181SPeter Brune   dd->No = No;
30495c13181SPeter Brune   dd->Po = Po;
30595c13181SPeter Brune 
30695c13181SPeter Brune   if (da->coordinateDM) {
30795c13181SPeter Brune     ierr = DMDASetOffset(da->coordinateDM,xo,yo,zo,Mo,No,Po);CHKERRQ(ierr);
30895c13181SPeter Brune   }
309d886c4f4SPeter Brune   PetscFunctionReturn(0);
310d886c4f4SPeter Brune }
311d886c4f4SPeter Brune 
312d886c4f4SPeter Brune #undef __FUNCT__
313d886c4f4SPeter Brune #define __FUNCT__ "DMDAGetOffset"
314d886c4f4SPeter Brune /*@
315d886c4f4SPeter Brune   DMDAGetOffset - Gets the index offset of the DA.
316d886c4f4SPeter Brune 
317d886c4f4SPeter Brune   Not collective
318d886c4f4SPeter Brune 
319d886c4f4SPeter Brune   Input Parameter:
320d886c4f4SPeter Brune . da  - The DMDA
321d886c4f4SPeter Brune 
322d886c4f4SPeter Brune   Output Parameters:
323d886c4f4SPeter Brune + xo  - The offset in the x direction
324d886c4f4SPeter Brune . yo  - The offset in the y direction
32595c13181SPeter Brune . zo  - The offset in the z direction
32695c13181SPeter Brune . Mo  - The global size in the x direction
32795c13181SPeter Brune . No  - The global size in the y direction
32895c13181SPeter Brune - Po  - The global size in the z direction
329d886c4f4SPeter Brune 
330d886c4f4SPeter Brune   Level: intermediate
331d886c4f4SPeter Brune 
332d886c4f4SPeter Brune .keywords:  distributed array, degrees of freedom
333d886c4f4SPeter Brune .seealso: DMDASetOffset(), DMDAVecGetArray()
334d886c4f4SPeter Brune @*/
33595c13181SPeter Brune PetscErrorCode  DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po)
336d886c4f4SPeter Brune {
337d886c4f4SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
338d886c4f4SPeter Brune 
339d886c4f4SPeter Brune   PetscFunctionBegin;
340d886c4f4SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
341d886c4f4SPeter Brune   if (xo) *xo = dd->xo;
342d886c4f4SPeter Brune   if (yo) *yo = dd->yo;
343d886c4f4SPeter Brune   if (zo) *zo = dd->zo;
34495c13181SPeter Brune   if (Mo) *Mo = dd->Mo;
34595c13181SPeter Brune   if (No) *No = dd->No;
34695c13181SPeter Brune   if (Po) *Po = dd->Po;
347d886c4f4SPeter Brune   PetscFunctionReturn(0);
348d886c4f4SPeter Brune }
349d886c4f4SPeter Brune 
35040234c92SPeter Brune #undef __FUNCT__
35140234c92SPeter Brune #define __FUNCT__ "DMDAGetNonOverlappingRegion"
35240234c92SPeter Brune /*@
35340234c92SPeter Brune   DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain DM.
35440234c92SPeter Brune 
35540234c92SPeter Brune   Not collective
35640234c92SPeter Brune 
35740234c92SPeter Brune   Input Parameter:
35840234c92SPeter Brune . da  - The DMDA
35940234c92SPeter Brune 
36040234c92SPeter Brune   Output Parameters:
36140234c92SPeter Brune + xs  - The start of the region in x
36240234c92SPeter Brune . ys  - The start of the region in y
36340234c92SPeter Brune . zs  - The start of the region in z
36440234c92SPeter Brune . xs  - The size of the region in x
36540234c92SPeter Brune . ys  - The size of the region in y
36640234c92SPeter Brune . zs  - The size of the region in z
36740234c92SPeter Brune 
36840234c92SPeter Brune   Level: intermediate
36940234c92SPeter Brune 
37040234c92SPeter Brune .keywords:  distributed array, degrees of freedom
37140234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray()
37240234c92SPeter Brune @*/
37340234c92SPeter Brune PetscErrorCode  DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm)
37440234c92SPeter Brune {
37540234c92SPeter Brune   DM_DA          *dd = (DM_DA*)da->data;
37640234c92SPeter Brune 
37740234c92SPeter Brune   PetscFunctionBegin;
37840234c92SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
37940234c92SPeter Brune   if (xs) *xs = dd->nonxs;
38040234c92SPeter Brune   if (ys) *ys = dd->nonys;
38140234c92SPeter Brune   if (zs) *zs = dd->nonzs;
38240234c92SPeter Brune   if (xm) *xm = dd->nonxm;
38340234c92SPeter Brune   if (ym) *ym = dd->nonym;
38440234c92SPeter Brune   if (zm) *zm = dd->nonzm;
38540234c92SPeter Brune   PetscFunctionReturn(0);
38640234c92SPeter Brune }
38740234c92SPeter Brune 
38840234c92SPeter Brune 
38940234c92SPeter Brune #undef __FUNCT__
39040234c92SPeter Brune #define __FUNCT__ "DMDASetNonOverlappingRegion"
39140234c92SPeter Brune /*@
39240234c92SPeter Brune   DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain DM.
39340234c92SPeter Brune 
39440234c92SPeter Brune   Collective on DA
39540234c92SPeter Brune 
39640234c92SPeter Brune   Input Parameter:
39740234c92SPeter Brune + da  - The DMDA
39840234c92SPeter Brune . xs  - The start of the region in x
39940234c92SPeter Brune . ys  - The start of the region in y
40040234c92SPeter Brune . zs  - The start of the region in z
40140234c92SPeter Brune . xs  - The size of the region in x
40240234c92SPeter Brune . ys  - The size of the region in y
40340234c92SPeter Brune . zs  - The size of the region in z
40440234c92SPeter Brune 
40540234c92SPeter Brune   Level: intermediate
40640234c92SPeter Brune 
40740234c92SPeter Brune .keywords:  distributed array, degrees of freedom
40840234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray()
40940234c92SPeter Brune @*/
41040234c92SPeter Brune PetscErrorCode  DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm)
41140234c92SPeter Brune {
41240234c92SPeter Brune   DM_DA          *dd = (DM_DA*)da->data;
41340234c92SPeter Brune 
41440234c92SPeter Brune   PetscFunctionBegin;
41540234c92SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
41640234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,xs,2);
41740234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,ys,3);
41840234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,zs,4);
41940234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,xm,5);
42040234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,ym,6);
42140234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,zm,7);
42240234c92SPeter Brune   dd->nonxs = xs;
42340234c92SPeter Brune   dd->nonys = ys;
42440234c92SPeter Brune   dd->nonzs = zs;
42540234c92SPeter Brune   dd->nonxm = xm;
42640234c92SPeter Brune   dd->nonym = ym;
42740234c92SPeter Brune   dd->nonzm = zm;
42840234c92SPeter Brune 
42940234c92SPeter Brune   PetscFunctionReturn(0);
43040234c92SPeter Brune }
43188661749SPeter Brune 
43288661749SPeter Brune #undef __FUNCT__
433aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType"
43447c6ae99SBarry Smith /*@
435aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
43647c6ae99SBarry Smith 
437aa219208SBarry Smith   Logically Collective on DMDA
43847c6ae99SBarry Smith 
43947c6ae99SBarry Smith   Input Parameter:
440aa219208SBarry Smith + da    - The DMDA
441aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
44247c6ae99SBarry Smith 
44347c6ae99SBarry Smith   Level: intermediate
44447c6ae99SBarry Smith 
44547c6ae99SBarry Smith .keywords:  distributed array, stencil
44696e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
44747c6ae99SBarry Smith @*/
4487087cfbeSBarry Smith PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
44947c6ae99SBarry Smith {
45047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
45147c6ae99SBarry Smith 
45247c6ae99SBarry Smith   PetscFunctionBegin;
45347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
45447c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,stype,2);
455ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
45647c6ae99SBarry Smith   dd->stencil_type = stype;
45747c6ae99SBarry Smith   PetscFunctionReturn(0);
45847c6ae99SBarry Smith }
45947c6ae99SBarry Smith 
46047c6ae99SBarry Smith #undef __FUNCT__
461aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth"
46247c6ae99SBarry Smith /*@
463aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
46447c6ae99SBarry Smith 
465aa219208SBarry Smith   Logically Collective on DMDA
46647c6ae99SBarry Smith 
46747c6ae99SBarry Smith   Input Parameter:
468aa219208SBarry Smith + da    - The DMDA
46947c6ae99SBarry Smith - width - The stencil width
47047c6ae99SBarry Smith 
47147c6ae99SBarry Smith   Level: intermediate
47247c6ae99SBarry Smith 
47347c6ae99SBarry Smith .keywords:  distributed array, stencil
47496e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
47547c6ae99SBarry Smith @*/
4767087cfbeSBarry Smith PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
47747c6ae99SBarry Smith {
47847c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
47947c6ae99SBarry Smith 
48047c6ae99SBarry Smith   PetscFunctionBegin;
48147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
48247c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,width,2);
483ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
48447c6ae99SBarry Smith   dd->s = width;
48547c6ae99SBarry Smith   PetscFunctionReturn(0);
48647c6ae99SBarry Smith }
48747c6ae99SBarry Smith 
48847c6ae99SBarry Smith #undef __FUNCT__
489aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
490aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
49147c6ae99SBarry Smith {
49247c6ae99SBarry Smith   PetscInt i,sum;
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith   PetscFunctionBegin;
495ce94432eSBarry Smith   if (M < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
49647c6ae99SBarry Smith   for (i=sum=0; i<m; i++) sum += lx[i];
497ce94432eSBarry Smith   if (sum != M) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
49847c6ae99SBarry Smith   PetscFunctionReturn(0);
49947c6ae99SBarry Smith }
50047c6ae99SBarry Smith 
50147c6ae99SBarry Smith #undef __FUNCT__
502aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges"
50347c6ae99SBarry Smith /*@
504aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
50547c6ae99SBarry Smith 
506aa219208SBarry Smith   Logically Collective on DMDA
50747c6ae99SBarry Smith 
50847c6ae99SBarry Smith   Input Parameter:
509aa219208SBarry Smith + da - The DMDA
5100298fd71SBarry 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
5110298fd71SBarry 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
5120298fd71SBarry 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.
51347c6ae99SBarry Smith 
51447c6ae99SBarry Smith   Level: intermediate
51547c6ae99SBarry Smith 
516e3f3e4b6SBarry Smith   Note: these numbers are NOT multiplied by the number of dof per node.
517e3f3e4b6SBarry Smith 
51847c6ae99SBarry Smith .keywords:  distributed array
51996e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
52047c6ae99SBarry Smith @*/
5217087cfbeSBarry Smith PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
52247c6ae99SBarry Smith {
52347c6ae99SBarry Smith   PetscErrorCode ierr;
52447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
52547c6ae99SBarry Smith 
52647c6ae99SBarry Smith   PetscFunctionBegin;
52747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
528ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
52947c6ae99SBarry Smith   if (lx) {
530ce94432eSBarry Smith     if (dd->m < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
531aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
53247c6ae99SBarry Smith     if (!dd->lx) {
533785e854fSJed Brown       ierr = PetscMalloc1(dd->m, &dd->lx);CHKERRQ(ierr);
53447c6ae99SBarry Smith     }
53547c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
53647c6ae99SBarry Smith   }
53747c6ae99SBarry Smith   if (ly) {
538ce94432eSBarry Smith     if (dd->n < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
539aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
54047c6ae99SBarry Smith     if (!dd->ly) {
541785e854fSJed Brown       ierr = PetscMalloc1(dd->n, &dd->ly);CHKERRQ(ierr);
54247c6ae99SBarry Smith     }
54347c6ae99SBarry Smith     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
54447c6ae99SBarry Smith   }
54547c6ae99SBarry Smith   if (lz) {
546ce94432eSBarry Smith     if (dd->p < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
547aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
54847c6ae99SBarry Smith     if (!dd->lz) {
549785e854fSJed Brown       ierr = PetscMalloc1(dd->p, &dd->lz);CHKERRQ(ierr);
55047c6ae99SBarry Smith     }
55147c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
55247c6ae99SBarry Smith   }
55347c6ae99SBarry Smith   PetscFunctionReturn(0);
55447c6ae99SBarry Smith }
55547c6ae99SBarry Smith 
55647c6ae99SBarry Smith #undef __FUNCT__
557aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType"
55847c6ae99SBarry Smith /*@
559aa219208SBarry Smith        DMDASetInterpolationType - Sets the type of interpolation that will be
560e727c939SJed Brown           returned by DMCreateInterpolation()
56147c6ae99SBarry Smith 
562aa219208SBarry Smith    Logically Collective on DMDA
56347c6ae99SBarry Smith 
56447c6ae99SBarry Smith    Input Parameter:
56547c6ae99SBarry Smith +  da - initial distributed array
566aa219208SBarry Smith .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
56747c6ae99SBarry Smith 
56847c6ae99SBarry Smith    Level: intermediate
56947c6ae99SBarry Smith 
570e727c939SJed Brown    Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation()
57147c6ae99SBarry Smith 
57247c6ae99SBarry Smith .keywords:  distributed array, interpolation
57347c6ae99SBarry Smith 
57496e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
57547c6ae99SBarry Smith @*/
5767087cfbeSBarry Smith PetscErrorCode  DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
57747c6ae99SBarry Smith {
57847c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
57947c6ae99SBarry Smith 
58047c6ae99SBarry Smith   PetscFunctionBegin;
58147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
58247c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,ctype,2);
58347c6ae99SBarry Smith   dd->interptype = ctype;
58447c6ae99SBarry Smith   PetscFunctionReturn(0);
58547c6ae99SBarry Smith }
58647c6ae99SBarry Smith 
5872dde6fd4SLisandro Dalcin #undef __FUNCT__
5882dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetInterpolationType"
5892dde6fd4SLisandro Dalcin /*@
5902dde6fd4SLisandro Dalcin        DMDAGetInterpolationType - Gets the type of interpolation that will be
591e727c939SJed Brown           used by DMCreateInterpolation()
5922dde6fd4SLisandro Dalcin 
5932dde6fd4SLisandro Dalcin    Not Collective
5942dde6fd4SLisandro Dalcin 
5952dde6fd4SLisandro Dalcin    Input Parameter:
5962dde6fd4SLisandro Dalcin .  da - distributed array
5972dde6fd4SLisandro Dalcin 
5982dde6fd4SLisandro Dalcin    Output Parameter:
5992dde6fd4SLisandro Dalcin .  ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms)
6002dde6fd4SLisandro Dalcin 
6012dde6fd4SLisandro Dalcin    Level: intermediate
6022dde6fd4SLisandro Dalcin 
6032dde6fd4SLisandro Dalcin .keywords:  distributed array, interpolation
6042dde6fd4SLisandro Dalcin 
605e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation()
6062dde6fd4SLisandro Dalcin @*/
6072dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype)
6082dde6fd4SLisandro Dalcin {
6092dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA*)da->data;
6102dde6fd4SLisandro Dalcin 
6112dde6fd4SLisandro Dalcin   PetscFunctionBegin;
6122dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
6132dde6fd4SLisandro Dalcin   PetscValidPointer(ctype,2);
6142dde6fd4SLisandro Dalcin   *ctype = dd->interptype;
6152dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
6162dde6fd4SLisandro Dalcin }
61747c6ae99SBarry Smith 
61847c6ae99SBarry Smith #undef __FUNCT__
619aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors"
62047c6ae99SBarry Smith /*@C
621aa219208SBarry Smith       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
62247c6ae99SBarry Smith         processes neighbors.
62347c6ae99SBarry Smith 
62447c6ae99SBarry Smith     Not Collective
62547c6ae99SBarry Smith 
62647c6ae99SBarry Smith    Input Parameter:
627aa219208SBarry Smith .     da - the DMDA object
62847c6ae99SBarry Smith 
62947c6ae99SBarry Smith    Output Parameters:
63047c6ae99SBarry Smith .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
63147c6ae99SBarry Smith               this process itself is in the list
63247c6ae99SBarry Smith 
63347c6ae99SBarry Smith    Notes: In 2d the array is of length 9, in 3d of length 27
63447c6ae99SBarry Smith           Not supported in 1d
635aa219208SBarry Smith           Do not free the array, it is freed when the DMDA is destroyed.
63647c6ae99SBarry Smith 
63747c6ae99SBarry Smith    Fortran Notes: In fortran you must pass in an array of the appropriate length.
63847c6ae99SBarry Smith 
63947c6ae99SBarry Smith    Level: intermediate
64047c6ae99SBarry Smith 
64147c6ae99SBarry Smith @*/
6427087cfbeSBarry Smith PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
64347c6ae99SBarry Smith {
64447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
6455fd66863SKarl Rupp 
64647c6ae99SBarry Smith   PetscFunctionBegin;
64747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
64847c6ae99SBarry Smith   *ranks = dd->neighbors;
64947c6ae99SBarry Smith   PetscFunctionReturn(0);
65047c6ae99SBarry Smith }
65147c6ae99SBarry Smith 
65247c6ae99SBarry Smith #undef __FUNCT__
653aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges"
65447c6ae99SBarry Smith /*@C
655aa219208SBarry Smith       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
65647c6ae99SBarry Smith 
65747c6ae99SBarry Smith     Not Collective
65847c6ae99SBarry Smith 
65947c6ae99SBarry Smith    Input Parameter:
660aa219208SBarry Smith .     da - the DMDA object
66147c6ae99SBarry Smith 
66247c6ae99SBarry Smith    Output Parameter:
66347c6ae99SBarry Smith +     lx - ownership along x direction (optional)
66447c6ae99SBarry Smith .     ly - ownership along y direction (optional)
66547c6ae99SBarry Smith -     lz - ownership along z direction (optional)
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith    Level: intermediate
66847c6ae99SBarry Smith 
669aa219208SBarry Smith     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
67047c6ae99SBarry Smith 
67147c6ae99SBarry Smith     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
672aa219208SBarry Smith     eighth arguments from DMDAGetInfo()
67347c6ae99SBarry Smith 
67447c6ae99SBarry Smith      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
675aa219208SBarry Smith     DMDA they came from still exists (has not been destroyed).
67647c6ae99SBarry Smith 
677e3f3e4b6SBarry Smith     These numbers are NOT multiplied by the number of dof per node.
678e3f3e4b6SBarry Smith 
679aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
68047c6ae99SBarry Smith @*/
6817087cfbeSBarry Smith PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
68247c6ae99SBarry Smith {
68347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
68447c6ae99SBarry Smith 
68547c6ae99SBarry Smith   PetscFunctionBegin;
68647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
68747c6ae99SBarry Smith   if (lx) *lx = dd->lx;
68847c6ae99SBarry Smith   if (ly) *ly = dd->ly;
68947c6ae99SBarry Smith   if (lz) *lz = dd->lz;
69047c6ae99SBarry Smith   PetscFunctionReturn(0);
69147c6ae99SBarry Smith }
69247c6ae99SBarry Smith 
69347c6ae99SBarry Smith #undef __FUNCT__
694aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor"
69547c6ae99SBarry Smith /*@
696aa219208SBarry Smith      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
69747c6ae99SBarry Smith 
698aa219208SBarry Smith     Logically Collective on DMDA
69947c6ae99SBarry Smith 
70047c6ae99SBarry Smith   Input Parameters:
701aa219208SBarry Smith +    da - the DMDA object
70247c6ae99SBarry Smith .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
70347c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
70447c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
70547c6ae99SBarry Smith 
70647c6ae99SBarry Smith   Options Database:
70747c6ae99SBarry Smith +  -da_refine_x - refinement ratio in x direction
70847c6ae99SBarry Smith .  -da_refine_y - refinement ratio in y direction
70947c6ae99SBarry Smith -  -da_refine_z - refinement ratio in z direction
71047c6ae99SBarry Smith 
71147c6ae99SBarry Smith   Level: intermediate
71247c6ae99SBarry Smith 
71347c6ae99SBarry Smith     Notes: Pass PETSC_IGNORE to leave a value unchanged
71447c6ae99SBarry Smith 
715aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor()
71647c6ae99SBarry Smith @*/
7177087cfbeSBarry Smith PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
71847c6ae99SBarry Smith {
71947c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
72047c6ae99SBarry Smith 
72147c6ae99SBarry Smith   PetscFunctionBegin;
72247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
72347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_x,2);
72447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_y,3);
72547c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_z,4);
72647c6ae99SBarry Smith 
72747c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
72847c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
72947c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
73047c6ae99SBarry Smith   PetscFunctionReturn(0);
73147c6ae99SBarry Smith }
73247c6ae99SBarry Smith 
73347c6ae99SBarry Smith #undef __FUNCT__
734aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor"
73547c6ae99SBarry Smith /*@C
736aa219208SBarry Smith      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
73747c6ae99SBarry Smith 
73847c6ae99SBarry Smith     Not Collective
73947c6ae99SBarry Smith 
74047c6ae99SBarry Smith   Input Parameter:
741aa219208SBarry Smith .    da - the DMDA object
74247c6ae99SBarry Smith 
74347c6ae99SBarry Smith   Output Parameters:
74447c6ae99SBarry Smith +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
74547c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
74647c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
74747c6ae99SBarry Smith 
74847c6ae99SBarry Smith   Level: intermediate
74947c6ae99SBarry Smith 
7500298fd71SBarry Smith     Notes: Pass NULL for values you do not need
75147c6ae99SBarry Smith 
752aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor()
75347c6ae99SBarry Smith @*/
7547087cfbeSBarry Smith PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
75547c6ae99SBarry Smith {
75647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
75747c6ae99SBarry Smith 
75847c6ae99SBarry Smith   PetscFunctionBegin;
75947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
76047c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
76147c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
76247c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
76347c6ae99SBarry Smith   PetscFunctionReturn(0);
76447c6ae99SBarry Smith }
76547c6ae99SBarry Smith 
76647c6ae99SBarry Smith #undef __FUNCT__
767aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix"
76847c6ae99SBarry Smith /*@C
769aa219208SBarry Smith      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
77047c6ae99SBarry Smith 
771aa219208SBarry Smith     Logically Collective on DMDA
77247c6ae99SBarry Smith 
77347c6ae99SBarry Smith   Input Parameters:
774aa219208SBarry Smith +    da - the DMDA object
775aa219208SBarry Smith -    f - the function that allocates the matrix for that specific DMDA
77647c6ae99SBarry Smith 
77747c6ae99SBarry Smith   Level: developer
77847c6ae99SBarry Smith 
779aa219208SBarry Smith    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
78047c6ae99SBarry Smith        the diagonal and off-diagonal blocks of the matrix
78147c6ae99SBarry Smith 
782950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills()
78347c6ae99SBarry Smith @*/
784b412c318SBarry Smith PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, Mat*))
78547c6ae99SBarry Smith {
78647c6ae99SBarry Smith   PetscFunctionBegin;
78747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
78825296bd5SBarry Smith   da->ops->creatematrix = f;
78947c6ae99SBarry Smith   PetscFunctionReturn(0);
79047c6ae99SBarry Smith }
79147c6ae99SBarry Smith 
79247c6ae99SBarry Smith #undef __FUNCT__
79394013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges"
79447c6ae99SBarry Smith /*
79547c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
79647c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
79747c6ae99SBarry Smith 
79847c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
79947c6ae99SBarry Smith */
80094013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
80147c6ae99SBarry Smith {
80247c6ae99SBarry Smith   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
80347c6ae99SBarry Smith   PetscErrorCode ierr;
80447c6ae99SBarry Smith 
80547c6ae99SBarry Smith   PetscFunctionBegin;
806ce94432eSBarry Smith   if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
80747c6ae99SBarry Smith   if (ratio == 1) {
80847c6ae99SBarry Smith     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
80947c6ae99SBarry Smith     PetscFunctionReturn(0);
81047c6ae99SBarry Smith   }
81147c6ae99SBarry Smith   for (i=0; i<m; i++) totalc += lc[i];
81247c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
81347c6ae99SBarry Smith   for (i=0; i<m; i++) {
81447c6ae99SBarry Smith     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
81547c6ae99SBarry Smith     if (i == m-1) lf[i] = want;
81647c6ae99SBarry Smith     else {
8177aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
8187aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
8197aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
8207aca7175SJed Brown       while ((startf+want)/ratio < nextc - stencil_width) want++;
8217aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
8227aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
8237aca7175SJed Brown       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
8247aca7175SJed Brown       /* Make sure all constraints are satisfied */
82530729d88SBarry Smith       if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width)
826ce94432eSBarry Smith           || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range");
82747c6ae99SBarry Smith     }
82847c6ae99SBarry Smith     lf[i]      = want;
82947c6ae99SBarry Smith     startc    += lc[i];
83047c6ae99SBarry Smith     startf    += lf[i];
83147c6ae99SBarry Smith     remaining -= lf[i];
83247c6ae99SBarry Smith   }
83347c6ae99SBarry Smith   PetscFunctionReturn(0);
83447c6ae99SBarry Smith }
83547c6ae99SBarry Smith 
83647c6ae99SBarry Smith #undef __FUNCT__
8372be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges"
8382be375d4SJed Brown /*
8392be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
8402be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
8412be375d4SJed Brown 
8422be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
8432be375d4SJed Brown */
8442be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
8452be375d4SJed Brown {
8462be375d4SJed Brown   PetscInt       i,totalf,remaining,startc,startf;
8472be375d4SJed Brown   PetscErrorCode ierr;
8482be375d4SJed Brown 
8492be375d4SJed Brown   PetscFunctionBegin;
850ce94432eSBarry Smith   if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
8512be375d4SJed Brown   if (ratio == 1) {
8522be375d4SJed Brown     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
8532be375d4SJed Brown     PetscFunctionReturn(0);
8542be375d4SJed Brown   }
8552be375d4SJed Brown   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
8562be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
8572be375d4SJed Brown   for (i=0,startc=0,startf=0; i<m; i++) {
8582be375d4SJed Brown     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
8592be375d4SJed Brown     if (i == m-1) lc[i] = want;
8602be375d4SJed Brown     else {
8612be375d4SJed Brown       const PetscInt nextf = startf+lf[i];
8622be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
8632be375d4SJed Brown        * node is within one stencil width. */
8642be375d4SJed Brown       while (nextf/ratio < startc+want-stencil_width) want--;
8652be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
8662be375d4SJed Brown        * fine node is within one stencil width. */
8672be375d4SJed Brown       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
8682be375d4SJed Brown       if (want < 0 || want > remaining
869ce94432eSBarry Smith           || (nextf/ratio < startc+want-stencil_width) || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range");
8702be375d4SJed Brown     }
8712be375d4SJed Brown     lc[i]      = want;
8722be375d4SJed Brown     startc    += lc[i];
8732be375d4SJed Brown     startf    += lf[i];
8742be375d4SJed Brown     remaining -= lc[i];
8752be375d4SJed Brown   }
8762be375d4SJed Brown   PetscFunctionReturn(0);
8772be375d4SJed Brown }
8782be375d4SJed Brown 
8792be375d4SJed Brown #undef __FUNCT__
88094013140SBarry Smith #define __FUNCT__ "DMRefine_DA"
8817087cfbeSBarry Smith PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
88247c6ae99SBarry Smith {
88347c6ae99SBarry Smith   PetscErrorCode ierr;
884*c73cfb54SMatthew G. Knepley   PetscInt       M,N,P,i,dim;
8859a42bb27SBarry Smith   DM             da2;
88647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
88747c6ae99SBarry Smith 
88847c6ae99SBarry Smith   PetscFunctionBegin;
88947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
89047c6ae99SBarry Smith   PetscValidPointer(daref,3);
89147c6ae99SBarry Smith 
892*c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(da, &dim);CHKERRQ(ierr);
893bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
89447c6ae99SBarry Smith     M = dd->refine_x*dd->M;
89547c6ae99SBarry Smith   } else {
89647c6ae99SBarry Smith     M = 1 + dd->refine_x*(dd->M - 1);
89747c6ae99SBarry Smith   }
898bff4a2f0SMatthew G. Knepley   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
899*c73cfb54SMatthew G. Knepley     if (dim > 1) {
90047c6ae99SBarry Smith       N = dd->refine_y*dd->N;
90147c6ae99SBarry Smith     } else {
9021860e6e9SBarry Smith       N = 1;
9031860e6e9SBarry Smith     }
9041860e6e9SBarry Smith   } else {
90547c6ae99SBarry Smith     N = 1 + dd->refine_y*(dd->N - 1);
90647c6ae99SBarry Smith   }
907bff4a2f0SMatthew G. Knepley   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
908*c73cfb54SMatthew G. Knepley     if (dim > 2) {
90947c6ae99SBarry Smith       P = dd->refine_z*dd->P;
91047c6ae99SBarry Smith     } else {
9111860e6e9SBarry Smith       P = 1;
9121860e6e9SBarry Smith     }
9131860e6e9SBarry Smith   } else {
91447c6ae99SBarry Smith     P = 1 + dd->refine_z*(dd->P - 1);
91547c6ae99SBarry Smith   }
916ce94432eSBarry Smith   ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr);
917397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
918*c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(da2,dim);CHKERRQ(ierr);
919397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
920397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
921397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
922397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
923397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
924397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
925*c73cfb54SMatthew G. Knepley   if (dim == 3) {
92647c6ae99SBarry Smith     PetscInt *lx,*ly,*lz;
927dcca6d9dSJed Brown     ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr);
928bff4a2f0SMatthew G. Knepley     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
929bff4a2f0SMatthew G. Knepley     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
930bff4a2f0SMatthew G. Knepley     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
931397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
93247c6ae99SBarry Smith     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
933*c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
93447c6ae99SBarry Smith     PetscInt *lx,*ly;
935dcca6d9dSJed Brown     ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr);
936bff4a2f0SMatthew G. Knepley     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
937bff4a2f0SMatthew G. Knepley     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
9380298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr);
93947c6ae99SBarry Smith     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
940*c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
94147c6ae99SBarry Smith     PetscInt *lx;
942785e854fSJed Brown     ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr);
943bff4a2f0SMatthew G. Knepley     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
9440298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr);
94547c6ae99SBarry Smith     ierr = PetscFree(lx);CHKERRQ(ierr);
94647c6ae99SBarry Smith   }
94747c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
94847c6ae99SBarry Smith 
94947c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
95025296bd5SBarry Smith   da2->ops->creatematrix = da->ops->creatematrix;
95125296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
95247c6ae99SBarry Smith   da2->ops->getcoloring = da->ops->getcoloring;
95347c6ae99SBarry Smith   dd2->interptype       = dd->interptype;
95447c6ae99SBarry Smith 
95547c6ae99SBarry Smith   /* copy fill information if given */
95647c6ae99SBarry Smith   if (dd->dfill) {
957785e854fSJed Brown     ierr = PetscMalloc1((dd->dfill[dd->w]+dd->w+1),&dd2->dfill);CHKERRQ(ierr);
95847c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
95947c6ae99SBarry Smith   }
96047c6ae99SBarry Smith   if (dd->ofill) {
961785e854fSJed Brown     ierr = PetscMalloc1((dd->ofill[dd->w]+dd->w+1),&dd2->ofill);CHKERRQ(ierr);
96247c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
96347c6ae99SBarry Smith   }
96447c6ae99SBarry Smith   /* copy the refine information */
965397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
966397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
967397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
96847c6ae99SBarry Smith 
96947c6ae99SBarry Smith   /* copy vector type information */
970c0dedaeaSBarry Smith   ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr);
971ddcf8b74SDave May 
972efd51863SBarry Smith   dd2->lf = dd->lf;
973efd51863SBarry Smith   dd2->lj = dd->lj;
974efd51863SBarry Smith 
9756e87535bSJed Brown   da2->leveldown = da->leveldown;
9766e87535bSJed Brown   da2->levelup   = da->levelup + 1;
9778865f1eaSKarl Rupp 
9786e87535bSJed Brown   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
9796e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
9806e87535bSJed Brown 
981ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
9826636e97aSMatthew G Knepley   if (da->coordinates) {
983ddcf8b74SDave May     DM  cdaf,cdac;
984ddcf8b74SDave May     Vec coordsc,coordsf;
985ddcf8b74SDave May     Mat II;
986ddcf8b74SDave May 
9876636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr);
9886636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr);
9896636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr);
990b61d3410SDave May     /* force creation of the coordinate vector */
991b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
9926636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
9930298fd71SBarry Smith     ierr = DMCreateInterpolation(cdac,cdaf,&II,NULL);CHKERRQ(ierr);
994ddcf8b74SDave May     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
995fcfd50ebSBarry Smith     ierr = MatDestroy(&II);CHKERRQ(ierr);
996ddcf8b74SDave May   }
997397b6216SJed Brown 
998f3141302SJed Brown   for (i=0; i<da->bs; i++) {
999f3141302SJed Brown     const char *fieldname;
1000f3141302SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
1001f3141302SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
1002f3141302SJed Brown   }
1003397b6216SJed Brown 
100447c6ae99SBarry Smith   *daref = da2;
100547c6ae99SBarry Smith   PetscFunctionReturn(0);
100647c6ae99SBarry Smith }
100747c6ae99SBarry Smith 
1008397b6216SJed Brown 
100947c6ae99SBarry Smith #undef __FUNCT__
101094013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA"
10117087cfbeSBarry Smith PetscErrorCode  DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
101247c6ae99SBarry Smith {
101347c6ae99SBarry Smith   PetscErrorCode ierr;
1014*c73cfb54SMatthew G. Knepley   PetscInt       M,N,P,i,dim;
10159a42bb27SBarry Smith   DM             da2;
101647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
101747c6ae99SBarry Smith 
101847c6ae99SBarry Smith   PetscFunctionBegin;
101947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
102047c6ae99SBarry Smith   PetscValidPointer(daref,3);
102147c6ae99SBarry Smith 
1022*c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(da, &dim);CHKERRQ(ierr);
1023bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1024397b6216SJed Brown     M = dd->M / dd->coarsen_x;
102547c6ae99SBarry Smith   } else {
1026397b6216SJed Brown     M = 1 + (dd->M - 1) / dd->coarsen_x;
102747c6ae99SBarry Smith   }
1028bff4a2f0SMatthew G. Knepley   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1029*c73cfb54SMatthew G. Knepley     if (dim > 1) {
1030397b6216SJed Brown       N = dd->N / dd->coarsen_y;
103147c6ae99SBarry Smith     } else {
10321860e6e9SBarry Smith       N = 1;
10331860e6e9SBarry Smith     }
10341860e6e9SBarry Smith   } else {
1035397b6216SJed Brown     N = 1 + (dd->N - 1) / dd->coarsen_y;
103647c6ae99SBarry Smith   }
1037bff4a2f0SMatthew G. Knepley   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1038*c73cfb54SMatthew G. Knepley     if (dim > 2) {
1039397b6216SJed Brown       P = dd->P / dd->coarsen_z;
104047c6ae99SBarry Smith     } else {
10411860e6e9SBarry Smith       P = 1;
10421860e6e9SBarry Smith     }
10431860e6e9SBarry Smith   } else {
1044397b6216SJed Brown     P = 1 + (dd->P - 1) / dd->coarsen_z;
104547c6ae99SBarry Smith   }
1046ce94432eSBarry Smith   ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr);
1047397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
1048*c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(da2,dim);CHKERRQ(ierr);
1049397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
1050397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
1051397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
1052397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
1053397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
1054397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
1055*c73cfb54SMatthew G. Knepley   if (dim == 3) {
10562be375d4SJed Brown     PetscInt *lx,*ly,*lz;
1057dcca6d9dSJed Brown     ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr);
1058bff4a2f0SMatthew G. Knepley     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
1059bff4a2f0SMatthew G. Knepley     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
1060bff4a2f0SMatthew G. Knepley     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
1061397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
10622be375d4SJed Brown     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
1063*c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
10642be375d4SJed Brown     PetscInt *lx,*ly;
1065dcca6d9dSJed Brown     ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr);
1066bff4a2f0SMatthew G. Knepley     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
1067bff4a2f0SMatthew G. Knepley     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
10680298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr);
10692be375d4SJed Brown     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
1070*c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
10712be375d4SJed Brown     PetscInt *lx;
1072785e854fSJed Brown     ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr);
1073bff4a2f0SMatthew G. Knepley     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
10740298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr);
10752be375d4SJed Brown     ierr = PetscFree(lx);CHKERRQ(ierr);
107647c6ae99SBarry Smith   }
107747c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
107847c6ae99SBarry Smith 
10794dcab191SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
108025296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */
108125296bd5SBarry Smith   da2->ops->creatematrix = da->ops->creatematrix;
108247c6ae99SBarry Smith   da2->ops->getcoloring  = da->ops->getcoloring;
108347c6ae99SBarry Smith   dd2->interptype        = dd->interptype;
108447c6ae99SBarry Smith 
108547c6ae99SBarry Smith   /* copy fill information if given */
108647c6ae99SBarry Smith   if (dd->dfill) {
1087785e854fSJed Brown     ierr = PetscMalloc1((dd->dfill[dd->w]+dd->w+1),&dd2->dfill);CHKERRQ(ierr);
108847c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
108947c6ae99SBarry Smith   }
109047c6ae99SBarry Smith   if (dd->ofill) {
1091785e854fSJed Brown     ierr = PetscMalloc1((dd->ofill[dd->w]+dd->w+1),&dd2->ofill);CHKERRQ(ierr);
109247c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
109347c6ae99SBarry Smith   }
109447c6ae99SBarry Smith   /* copy the refine information */
1095397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
1096397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
1097397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
109847c6ae99SBarry Smith 
109947c6ae99SBarry Smith   /* copy vector type information */
1100c0dedaeaSBarry Smith   ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr);
110147c6ae99SBarry Smith 
1102644e2e5bSBarry Smith   dd2->lf = dd->lf;
1103644e2e5bSBarry Smith   dd2->lj = dd->lj;
1104644e2e5bSBarry Smith 
11056e87535bSJed Brown   da2->leveldown = da->leveldown + 1;
11066e87535bSJed Brown   da2->levelup   = da->levelup;
11078865f1eaSKarl Rupp 
11086e87535bSJed Brown   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
11096e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
11106e87535bSJed Brown 
1111ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
11126636e97aSMatthew G Knepley   if (da->coordinates) {
1113ddcf8b74SDave May     DM         cdaf,cdac;
1114ddcf8b74SDave May     Vec        coordsc,coordsf;
1115ddcf8b74SDave May     VecScatter inject;
1116ddcf8b74SDave May 
11176636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr);
11186636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr);
11196636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr);
1120b61d3410SDave May     /* force creation of the coordinate vector */
1121b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
11226636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
1123ddcf8b74SDave May 
1124e727c939SJed Brown     ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
1125ddcf8b74SDave May     ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1126ddcf8b74SDave May     ierr = VecScatterEnd(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1127fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&inject);CHKERRQ(ierr);
1128ddcf8b74SDave May   }
1129f98405f7SJed Brown 
1130f98405f7SJed Brown   for (i=0; i<da->bs; i++) {
1131f98405f7SJed Brown     const char *fieldname;
1132f98405f7SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
1133f98405f7SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
1134f98405f7SJed Brown   }
1135f98405f7SJed Brown 
113647c6ae99SBarry Smith   *daref = da2;
113747c6ae99SBarry Smith   PetscFunctionReturn(0);
113847c6ae99SBarry Smith }
113947c6ae99SBarry Smith 
114047c6ae99SBarry Smith #undef __FUNCT__
114194013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA"
11427087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
114347c6ae99SBarry Smith {
114447c6ae99SBarry Smith   PetscErrorCode ierr;
114547c6ae99SBarry Smith   PetscInt       i,n,*refx,*refy,*refz;
114647c6ae99SBarry Smith 
114747c6ae99SBarry Smith   PetscFunctionBegin;
114847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1149ce94432eSBarry Smith   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
115047c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
115147c6ae99SBarry Smith   PetscValidPointer(daf,3);
115247c6ae99SBarry Smith 
1153aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
1154dcca6d9dSJed Brown   ierr = PetscMalloc3(nlevels,&refx,nlevels,&refy,nlevels,&refz);CHKERRQ(ierr);
115547c6ae99SBarry Smith   for (i=0; i<nlevels; i++) {
1156aa219208SBarry Smith     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
115747c6ae99SBarry Smith   }
115847c6ae99SBarry Smith   n    = nlevels;
11590298fd71SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr);
116047c6ae99SBarry Smith   n    = nlevels;
11610298fd71SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr);
116247c6ae99SBarry Smith   n    = nlevels;
11630298fd71SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr);
116447c6ae99SBarry Smith 
1165aa219208SBarry Smith   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
1166ce94432eSBarry Smith   ierr = DMRefine(da,PetscObjectComm((PetscObject)da),&daf[0]);CHKERRQ(ierr);
116747c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
1168aa219208SBarry Smith     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
1169ce94432eSBarry Smith     ierr = DMRefine(daf[i-1],PetscObjectComm((PetscObject)da),&daf[i]);CHKERRQ(ierr);
117047c6ae99SBarry Smith   }
117147c6ae99SBarry Smith   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
117247c6ae99SBarry Smith   PetscFunctionReturn(0);
117347c6ae99SBarry Smith }
117447c6ae99SBarry Smith 
117547c6ae99SBarry Smith #undef __FUNCT__
117694013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA"
11777087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
117847c6ae99SBarry Smith {
117947c6ae99SBarry Smith   PetscErrorCode ierr;
118047c6ae99SBarry Smith   PetscInt       i;
118147c6ae99SBarry Smith 
118247c6ae99SBarry Smith   PetscFunctionBegin;
118347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1184ce94432eSBarry Smith   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
118547c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
118647c6ae99SBarry Smith   PetscValidPointer(dac,3);
1187ce94432eSBarry Smith   ierr = DMCoarsen(da,PetscObjectComm((PetscObject)da),&dac[0]);CHKERRQ(ierr);
118847c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
1189ce94432eSBarry Smith     ierr = DMCoarsen(dac[i-1],PetscObjectComm((PetscObject)da),&dac[i]);CHKERRQ(ierr);
119047c6ae99SBarry Smith   }
119147c6ae99SBarry Smith   PetscFunctionReturn(0);
119247c6ae99SBarry Smith }
1193