xref: /petsc/src/dm/impls/da/da.c (revision d3be247eaf305b2db07a63a37f3a405e71b8058f)
1af0996ceSBarry 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;
68c73cfb54SMatthew G. Knepley   if (da->dim == 2) {
69*d3be247eSBarry Smith     PetscMPIInt size;
70*d3be247eSBarry 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 
12259f3ab6dSMatthew G. Knepley   Input Parameters:
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
129fb6725baSMatthew G. Knepley .seealso: DMDAGetDof(), 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__
145fb6725baSMatthew G. Knepley #define __FUNCT__ "DMDAGetDof"
146fb6725baSMatthew G. Knepley /*@
147fb6725baSMatthew G. Knepley   DMDAGetDof - Gets the number of degrees of freedom per vertex
148fb6725baSMatthew G. Knepley 
149fb6725baSMatthew G. Knepley   Not collective
150fb6725baSMatthew G. Knepley 
151fb6725baSMatthew G. Knepley   Input Parameter:
152fb6725baSMatthew G. Knepley . da  - The DMDA
153fb6725baSMatthew G. Knepley 
154fb6725baSMatthew G. Knepley   Output Parameter:
155fb6725baSMatthew G. Knepley . dof - Number of degrees of freedom
156fb6725baSMatthew G. Knepley 
157fb6725baSMatthew G. Knepley   Level: intermediate
158fb6725baSMatthew G. Knepley 
159fb6725baSMatthew G. Knepley .keywords:  distributed array, degrees of freedom
160fb6725baSMatthew G. Knepley .seealso: DMDASetDof(), DMDACreate(), DMDestroy(), DMDA
161fb6725baSMatthew G. Knepley @*/
162fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetDof(DM da, PetscInt *dof)
163fb6725baSMatthew G. Knepley {
164fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *) da->data;
165fb6725baSMatthew G. Knepley 
166fb6725baSMatthew G. Knepley   PetscFunctionBegin;
167fb6725baSMatthew G. Knepley   PetscValidHeaderSpecific(da,DM_CLASSID,1);
168fb6725baSMatthew G. Knepley   PetscValidPointer(dof,2);
169fb6725baSMatthew G. Knepley   *dof = dd->w;
170fb6725baSMatthew G. Knepley   PetscFunctionReturn(0);
171fb6725baSMatthew G. Knepley }
172fb6725baSMatthew G. Knepley 
173fb6725baSMatthew G. Knepley #undef __FUNCT__
1747ddda789SPeter Brune #define __FUNCT__ "DMDAGetOverlap"
1757ddda789SPeter Brune /*@
1767ddda789SPeter Brune   DMDAGetOverlap - Gets the size of the per-processor overlap.
1777ddda789SPeter Brune 
1787ddda789SPeter Brune   Not collective
1797ddda789SPeter Brune 
1807ddda789SPeter Brune   Input Parameters:
1817ddda789SPeter Brune . da  - The DMDA
1827ddda789SPeter Brune 
1837ddda789SPeter Brune   Output Parameters:
1847ddda789SPeter Brune + x   - Overlap in the x direction
1857ddda789SPeter Brune . y   - Overlap in the y direction
1867ddda789SPeter Brune - z   - Overlap in the z direction
1877ddda789SPeter Brune 
1887ddda789SPeter Brune   Level: intermediate
1897ddda789SPeter Brune 
1907ddda789SPeter Brune .keywords:  distributed array, overlap, domain decomposition
1917ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetOverlap(), DMDA
1927ddda789SPeter Brune @*/
1937ddda789SPeter Brune PetscErrorCode  DMDAGetOverlap(DM da,PetscInt *x,PetscInt *y,PetscInt *z)
1947ddda789SPeter Brune {
1957ddda789SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
1967ddda789SPeter Brune 
1977ddda789SPeter Brune   PetscFunctionBegin;
1987ddda789SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1997ddda789SPeter Brune   if (x) *x = dd->xol;
2007ddda789SPeter Brune   if (y) *y = dd->yol;
2017ddda789SPeter Brune   if (z) *z = dd->zol;
2027ddda789SPeter Brune   PetscFunctionReturn(0);
2037ddda789SPeter Brune }
2047ddda789SPeter Brune 
2057ddda789SPeter Brune #undef __FUNCT__
20688661749SPeter Brune #define __FUNCT__ "DMDASetOverlap"
20788661749SPeter Brune /*@
20888661749SPeter Brune   DMDASetOverlap - Sets the size of the per-processor overlap.
20988661749SPeter Brune 
21088661749SPeter Brune   Not collective
21188661749SPeter Brune 
2127ddda789SPeter Brune   Input Parameters:
21388661749SPeter Brune + da  - The DMDA
2147ddda789SPeter Brune . x   - Overlap in the x direction
2157ddda789SPeter Brune . y   - Overlap in the y direction
2167ddda789SPeter Brune - z   - Overlap in the z direction
21788661749SPeter Brune 
21888661749SPeter Brune   Level: intermediate
21988661749SPeter Brune 
2207ddda789SPeter Brune .keywords:  distributed array, overlap, domain decomposition
2217ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetOverlap(), DMDA
22288661749SPeter Brune @*/
2237ddda789SPeter Brune PetscErrorCode  DMDASetOverlap(DM da,PetscInt x,PetscInt y,PetscInt z)
22488661749SPeter Brune {
22588661749SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
22688661749SPeter Brune 
22788661749SPeter Brune   PetscFunctionBegin;
22888661749SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
2297ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da,x,2);
2307ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da,y,3);
2317ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da,z,4);
2327ddda789SPeter Brune   dd->xol = x;
2337ddda789SPeter Brune   dd->yol = y;
2347ddda789SPeter Brune   dd->zol = z;
23588661749SPeter Brune   PetscFunctionReturn(0);
23688661749SPeter Brune }
23788661749SPeter Brune 
2383e7870d2SPeter Brune 
2393e7870d2SPeter Brune #undef __FUNCT__
2403e7870d2SPeter Brune #define __FUNCT__ "DMDAGetNumLocalSubDomains"
2413e7870d2SPeter Brune /*@
2423e7870d2SPeter Brune   DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition.
2433e7870d2SPeter Brune 
2443e7870d2SPeter Brune   Not collective
2453e7870d2SPeter Brune 
2463e7870d2SPeter Brune   Input Parameters:
2473e7870d2SPeter Brune . da  - The DMDA
2483e7870d2SPeter Brune 
2493e7870d2SPeter Brune   Output Parameters:
2503e7870d2SPeter Brune + Nsub   - Number of local subdomains created upon decomposition
2513e7870d2SPeter Brune 
2523e7870d2SPeter Brune   Level: intermediate
2533e7870d2SPeter Brune 
2543e7870d2SPeter Brune .keywords:  distributed array, domain decomposition
2553e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetNumLocalSubDomains(), DMDA
2563e7870d2SPeter Brune @*/
2573e7870d2SPeter Brune PetscErrorCode  DMDAGetNumLocalSubDomains(DM da,PetscInt *Nsub)
2583e7870d2SPeter Brune {
2593e7870d2SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
2603e7870d2SPeter Brune 
2613e7870d2SPeter Brune   PetscFunctionBegin;
2623e7870d2SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
2633e7870d2SPeter Brune   if (Nsub) *Nsub = dd->Nsub;
2643e7870d2SPeter Brune   PetscFunctionReturn(0);
2653e7870d2SPeter Brune }
2663e7870d2SPeter Brune 
2673e7870d2SPeter Brune #undef __FUNCT__
2683e7870d2SPeter Brune #define __FUNCT__ "DMDASetNumLocalSubDomains"
2693e7870d2SPeter Brune /*@
2703e7870d2SPeter Brune   DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition.
2713e7870d2SPeter Brune 
2723e7870d2SPeter Brune   Not collective
2733e7870d2SPeter Brune 
2743e7870d2SPeter Brune   Input Parameters:
2753e7870d2SPeter Brune + da  - The DMDA
2763e7870d2SPeter Brune - Nsub - The number of local subdomains requested
2773e7870d2SPeter Brune 
2783e7870d2SPeter Brune   Level: intermediate
2793e7870d2SPeter Brune 
2803e7870d2SPeter Brune .keywords:  distributed array, domain decomposition
2813e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetNumLocalSubDomains(), DMDA
2823e7870d2SPeter Brune @*/
2833e7870d2SPeter Brune PetscErrorCode  DMDASetNumLocalSubDomains(DM da,PetscInt Nsub)
2843e7870d2SPeter Brune {
2853e7870d2SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
2863e7870d2SPeter Brune 
2873e7870d2SPeter Brune   PetscFunctionBegin;
2883e7870d2SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
2893e7870d2SPeter Brune   PetscValidLogicalCollectiveInt(da,Nsub,2);
2903e7870d2SPeter Brune   dd->Nsub = Nsub;
2913e7870d2SPeter Brune   PetscFunctionReturn(0);
2923e7870d2SPeter Brune }
2933e7870d2SPeter Brune 
294d886c4f4SPeter Brune #undef __FUNCT__
295d886c4f4SPeter Brune #define __FUNCT__ "DMDASetOffset"
296d886c4f4SPeter Brune /*@
297d886c4f4SPeter Brune   DMDASetOffset - Sets the index offset of the DA.
298d886c4f4SPeter Brune 
299d886c4f4SPeter Brune   Collective on DA
300d886c4f4SPeter Brune 
301d886c4f4SPeter Brune   Input Parameter:
302d886c4f4SPeter Brune + da  - The DMDA
303d886c4f4SPeter Brune . xo  - The offset in the x direction
304d886c4f4SPeter Brune . yo  - The offset in the y direction
305d886c4f4SPeter Brune - zo  - The offset in the z direction
306d886c4f4SPeter Brune 
307d886c4f4SPeter Brune   Level: intermediate
308d886c4f4SPeter Brune 
309d886c4f4SPeter Brune   Notes: This is used primarily to overlap a computation on a local DA with that on a global DA without
310d886c4f4SPeter Brune   changing boundary conditions or subdomain features that depend upon the global offsets.
311d886c4f4SPeter Brune 
312d886c4f4SPeter Brune .keywords:  distributed array, degrees of freedom
313d886c4f4SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray()
314d886c4f4SPeter Brune @*/
31595c13181SPeter Brune PetscErrorCode  DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po)
316d886c4f4SPeter Brune {
31795c13181SPeter Brune   PetscErrorCode ierr;
318d886c4f4SPeter Brune   DM_DA          *dd = (DM_DA*)da->data;
319d886c4f4SPeter Brune 
320d886c4f4SPeter Brune   PetscFunctionBegin;
321d886c4f4SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
322d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da,xo,2);
32395c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,yo,3);
32495c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,zo,4);
32595c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,Mo,5);
32695c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,No,6);
32795c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,Po,7);
328d886c4f4SPeter Brune   dd->xo = xo;
329d886c4f4SPeter Brune   dd->yo = yo;
330d886c4f4SPeter Brune   dd->zo = zo;
33195c13181SPeter Brune   dd->Mo = Mo;
33295c13181SPeter Brune   dd->No = No;
33395c13181SPeter Brune   dd->Po = Po;
33495c13181SPeter Brune 
33595c13181SPeter Brune   if (da->coordinateDM) {
33695c13181SPeter Brune     ierr = DMDASetOffset(da->coordinateDM,xo,yo,zo,Mo,No,Po);CHKERRQ(ierr);
33795c13181SPeter Brune   }
338d886c4f4SPeter Brune   PetscFunctionReturn(0);
339d886c4f4SPeter Brune }
340d886c4f4SPeter Brune 
341d886c4f4SPeter Brune #undef __FUNCT__
342d886c4f4SPeter Brune #define __FUNCT__ "DMDAGetOffset"
343d886c4f4SPeter Brune /*@
344d886c4f4SPeter Brune   DMDAGetOffset - Gets the index offset of the DA.
345d886c4f4SPeter Brune 
346d886c4f4SPeter Brune   Not collective
347d886c4f4SPeter Brune 
348d886c4f4SPeter Brune   Input Parameter:
349d886c4f4SPeter Brune . da  - The DMDA
350d886c4f4SPeter Brune 
351d886c4f4SPeter Brune   Output Parameters:
352d886c4f4SPeter Brune + xo  - The offset in the x direction
353d886c4f4SPeter Brune . yo  - The offset in the y direction
35495c13181SPeter Brune . zo  - The offset in the z direction
35595c13181SPeter Brune . Mo  - The global size in the x direction
35695c13181SPeter Brune . No  - The global size in the y direction
35795c13181SPeter Brune - Po  - The global size in the z direction
358d886c4f4SPeter Brune 
359d886c4f4SPeter Brune   Level: intermediate
360d886c4f4SPeter Brune 
361d886c4f4SPeter Brune .keywords:  distributed array, degrees of freedom
362d886c4f4SPeter Brune .seealso: DMDASetOffset(), DMDAVecGetArray()
363d886c4f4SPeter Brune @*/
36495c13181SPeter Brune PetscErrorCode  DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po)
365d886c4f4SPeter Brune {
366d886c4f4SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
367d886c4f4SPeter Brune 
368d886c4f4SPeter Brune   PetscFunctionBegin;
369d886c4f4SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
370d886c4f4SPeter Brune   if (xo) *xo = dd->xo;
371d886c4f4SPeter Brune   if (yo) *yo = dd->yo;
372d886c4f4SPeter Brune   if (zo) *zo = dd->zo;
37395c13181SPeter Brune   if (Mo) *Mo = dd->Mo;
37495c13181SPeter Brune   if (No) *No = dd->No;
37595c13181SPeter Brune   if (Po) *Po = dd->Po;
376d886c4f4SPeter Brune   PetscFunctionReturn(0);
377d886c4f4SPeter Brune }
378d886c4f4SPeter Brune 
37940234c92SPeter Brune #undef __FUNCT__
38040234c92SPeter Brune #define __FUNCT__ "DMDAGetNonOverlappingRegion"
38140234c92SPeter Brune /*@
38240234c92SPeter Brune   DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain DM.
38340234c92SPeter Brune 
38440234c92SPeter Brune   Not collective
38540234c92SPeter Brune 
38640234c92SPeter Brune   Input Parameter:
38740234c92SPeter Brune . da  - The DMDA
38840234c92SPeter Brune 
38940234c92SPeter Brune   Output Parameters:
39040234c92SPeter Brune + xs  - The start of the region in x
39140234c92SPeter Brune . ys  - The start of the region in y
39240234c92SPeter Brune . zs  - The start of the region in z
39340234c92SPeter Brune . xs  - The size of the region in x
39440234c92SPeter Brune . ys  - The size of the region in y
39540234c92SPeter Brune . zs  - The size of the region in z
39640234c92SPeter Brune 
39740234c92SPeter Brune   Level: intermediate
39840234c92SPeter Brune 
39940234c92SPeter Brune .keywords:  distributed array, degrees of freedom
40040234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray()
40140234c92SPeter Brune @*/
40240234c92SPeter Brune PetscErrorCode  DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm)
40340234c92SPeter Brune {
40440234c92SPeter Brune   DM_DA          *dd = (DM_DA*)da->data;
40540234c92SPeter Brune 
40640234c92SPeter Brune   PetscFunctionBegin;
40740234c92SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
40840234c92SPeter Brune   if (xs) *xs = dd->nonxs;
40940234c92SPeter Brune   if (ys) *ys = dd->nonys;
41040234c92SPeter Brune   if (zs) *zs = dd->nonzs;
41140234c92SPeter Brune   if (xm) *xm = dd->nonxm;
41240234c92SPeter Brune   if (ym) *ym = dd->nonym;
41340234c92SPeter Brune   if (zm) *zm = dd->nonzm;
41440234c92SPeter Brune   PetscFunctionReturn(0);
41540234c92SPeter Brune }
41640234c92SPeter Brune 
41740234c92SPeter Brune 
41840234c92SPeter Brune #undef __FUNCT__
41940234c92SPeter Brune #define __FUNCT__ "DMDASetNonOverlappingRegion"
42040234c92SPeter Brune /*@
42140234c92SPeter Brune   DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain DM.
42240234c92SPeter Brune 
42340234c92SPeter Brune   Collective on DA
42440234c92SPeter Brune 
42540234c92SPeter Brune   Input Parameter:
42640234c92SPeter Brune + da  - The DMDA
42740234c92SPeter Brune . xs  - The start of the region in x
42840234c92SPeter Brune . ys  - The start of the region in y
42940234c92SPeter Brune . zs  - The start of the region in z
43040234c92SPeter Brune . xs  - The size of the region in x
43140234c92SPeter Brune . ys  - The size of the region in y
43240234c92SPeter Brune . zs  - The size of the region in z
43340234c92SPeter Brune 
43440234c92SPeter Brune   Level: intermediate
43540234c92SPeter Brune 
43640234c92SPeter Brune .keywords:  distributed array, degrees of freedom
43740234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray()
43840234c92SPeter Brune @*/
43940234c92SPeter Brune PetscErrorCode  DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm)
44040234c92SPeter Brune {
44140234c92SPeter Brune   DM_DA          *dd = (DM_DA*)da->data;
44240234c92SPeter Brune 
44340234c92SPeter Brune   PetscFunctionBegin;
44440234c92SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
44540234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,xs,2);
44640234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,ys,3);
44740234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,zs,4);
44840234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,xm,5);
44940234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,ym,6);
45040234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,zm,7);
45140234c92SPeter Brune   dd->nonxs = xs;
45240234c92SPeter Brune   dd->nonys = ys;
45340234c92SPeter Brune   dd->nonzs = zs;
45440234c92SPeter Brune   dd->nonxm = xm;
45540234c92SPeter Brune   dd->nonym = ym;
45640234c92SPeter Brune   dd->nonzm = zm;
45740234c92SPeter Brune 
45840234c92SPeter Brune   PetscFunctionReturn(0);
45940234c92SPeter Brune }
46088661749SPeter Brune 
46188661749SPeter Brune #undef __FUNCT__
462aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType"
46347c6ae99SBarry Smith /*@
464aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
46547c6ae99SBarry Smith 
466aa219208SBarry Smith   Logically Collective on DMDA
46747c6ae99SBarry Smith 
46847c6ae99SBarry Smith   Input Parameter:
469aa219208SBarry Smith + da    - The DMDA
470aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
47147c6ae99SBarry Smith 
47247c6ae99SBarry Smith   Level: intermediate
47347c6ae99SBarry Smith 
47447c6ae99SBarry Smith .keywords:  distributed array, stencil
47596e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
47647c6ae99SBarry Smith @*/
4777087cfbeSBarry Smith PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
47847c6ae99SBarry Smith {
47947c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
48047c6ae99SBarry Smith 
48147c6ae99SBarry Smith   PetscFunctionBegin;
48247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
48347c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,stype,2);
484ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
48547c6ae99SBarry Smith   dd->stencil_type = stype;
48647c6ae99SBarry Smith   PetscFunctionReturn(0);
48747c6ae99SBarry Smith }
48847c6ae99SBarry Smith 
48947c6ae99SBarry Smith #undef __FUNCT__
490fb6725baSMatthew G. Knepley #define __FUNCT__ "DMDAGetStencilType"
491fb6725baSMatthew G. Knepley /*@
492fb6725baSMatthew G. Knepley   DMDAGetStencilType - Gets the type of the communication stencil
493fb6725baSMatthew G. Knepley 
494fb6725baSMatthew G. Knepley   Not collective
495fb6725baSMatthew G. Knepley 
496fb6725baSMatthew G. Knepley   Input Parameter:
497fb6725baSMatthew G. Knepley . da    - The DMDA
498fb6725baSMatthew G. Knepley 
499fb6725baSMatthew G. Knepley   Output Parameter:
500fb6725baSMatthew G. Knepley . stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
501fb6725baSMatthew G. Knepley 
502fb6725baSMatthew G. Knepley   Level: intermediate
503fb6725baSMatthew G. Knepley 
504fb6725baSMatthew G. Knepley .keywords:  distributed array, stencil
505fb6725baSMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA
506fb6725baSMatthew G. Knepley @*/
507fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype)
508fb6725baSMatthew G. Knepley {
509fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA*)da->data;
510fb6725baSMatthew G. Knepley 
511fb6725baSMatthew G. Knepley   PetscFunctionBegin;
512fb6725baSMatthew G. Knepley   PetscValidHeaderSpecific(da,DM_CLASSID,1);
513fb6725baSMatthew G. Knepley   PetscValidPointer(stype,2);
514fb6725baSMatthew G. Knepley   *stype = dd->stencil_type;
515fb6725baSMatthew G. Knepley   PetscFunctionReturn(0);
516fb6725baSMatthew G. Knepley }
517fb6725baSMatthew G. Knepley 
518fb6725baSMatthew G. Knepley #undef __FUNCT__
519aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth"
52047c6ae99SBarry Smith /*@
521aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
52247c6ae99SBarry Smith 
523aa219208SBarry Smith   Logically Collective on DMDA
52447c6ae99SBarry Smith 
52547c6ae99SBarry Smith   Input Parameter:
526aa219208SBarry Smith + da    - The DMDA
52747c6ae99SBarry Smith - width - The stencil width
52847c6ae99SBarry Smith 
52947c6ae99SBarry Smith   Level: intermediate
53047c6ae99SBarry Smith 
53147c6ae99SBarry Smith .keywords:  distributed array, stencil
53296e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
53347c6ae99SBarry Smith @*/
5347087cfbeSBarry Smith PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
53547c6ae99SBarry Smith {
53647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
53747c6ae99SBarry Smith 
53847c6ae99SBarry Smith   PetscFunctionBegin;
53947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
54047c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,width,2);
541ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
54247c6ae99SBarry Smith   dd->s = width;
54347c6ae99SBarry Smith   PetscFunctionReturn(0);
54447c6ae99SBarry Smith }
54547c6ae99SBarry Smith 
54647c6ae99SBarry Smith #undef __FUNCT__
547fb6725baSMatthew G. Knepley #define __FUNCT__ "DMDAGetStencilWidth"
548fb6725baSMatthew G. Knepley /*@
549fb6725baSMatthew G. Knepley   DMDAGetStencilWidth - Gets the width of the communication stencil
550fb6725baSMatthew G. Knepley 
551fb6725baSMatthew G. Knepley   Not collective
552fb6725baSMatthew G. Knepley 
553fb6725baSMatthew G. Knepley   Input Parameter:
554fb6725baSMatthew G. Knepley . da    - The DMDA
555fb6725baSMatthew G. Knepley 
556fb6725baSMatthew G. Knepley   Output Parameter:
557fb6725baSMatthew G. Knepley . width - The stencil width
558fb6725baSMatthew G. Knepley 
559fb6725baSMatthew G. Knepley   Level: intermediate
560fb6725baSMatthew G. Knepley 
561fb6725baSMatthew G. Knepley .keywords:  distributed array, stencil
562fb6725baSMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA
563fb6725baSMatthew G. Knepley @*/
564fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width)
565fb6725baSMatthew G. Knepley {
566fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *) da->data;
567fb6725baSMatthew G. Knepley 
568fb6725baSMatthew G. Knepley   PetscFunctionBegin;
569fb6725baSMatthew G. Knepley   PetscValidHeaderSpecific(da,DM_CLASSID,1);
570fb6725baSMatthew G. Knepley   PetscValidPointer(width,2);
571fb6725baSMatthew G. Knepley   *width = dd->s;
572fb6725baSMatthew G. Knepley   PetscFunctionReturn(0);
573fb6725baSMatthew G. Knepley }
574fb6725baSMatthew G. Knepley 
575fb6725baSMatthew G. Knepley #undef __FUNCT__
576aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
577aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
57847c6ae99SBarry Smith {
57947c6ae99SBarry Smith   PetscInt i,sum;
58047c6ae99SBarry Smith 
58147c6ae99SBarry Smith   PetscFunctionBegin;
582ce94432eSBarry Smith   if (M < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
58347c6ae99SBarry Smith   for (i=sum=0; i<m; i++) sum += lx[i];
584ce94432eSBarry Smith   if (sum != M) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
58547c6ae99SBarry Smith   PetscFunctionReturn(0);
58647c6ae99SBarry Smith }
58747c6ae99SBarry Smith 
58847c6ae99SBarry Smith #undef __FUNCT__
589aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges"
59047c6ae99SBarry Smith /*@
591aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
59247c6ae99SBarry Smith 
593aa219208SBarry Smith   Logically Collective on DMDA
59447c6ae99SBarry Smith 
59547c6ae99SBarry Smith   Input Parameter:
596aa219208SBarry Smith + da - The DMDA
5970298fd71SBarry 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
5980298fd71SBarry 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
5990298fd71SBarry 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.
60047c6ae99SBarry Smith 
60147c6ae99SBarry Smith   Level: intermediate
60247c6ae99SBarry Smith 
603e3f3e4b6SBarry Smith   Note: these numbers are NOT multiplied by the number of dof per node.
604e3f3e4b6SBarry Smith 
60547c6ae99SBarry Smith .keywords:  distributed array
60696e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
60747c6ae99SBarry Smith @*/
6087087cfbeSBarry Smith PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
60947c6ae99SBarry Smith {
61047c6ae99SBarry Smith   PetscErrorCode ierr;
61147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
61247c6ae99SBarry Smith 
61347c6ae99SBarry Smith   PetscFunctionBegin;
61447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
615ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
61647c6ae99SBarry Smith   if (lx) {
617ce94432eSBarry Smith     if (dd->m < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
618aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
61947c6ae99SBarry Smith     if (!dd->lx) {
620785e854fSJed Brown       ierr = PetscMalloc1(dd->m, &dd->lx);CHKERRQ(ierr);
62147c6ae99SBarry Smith     }
62247c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
62347c6ae99SBarry Smith   }
62447c6ae99SBarry Smith   if (ly) {
625ce94432eSBarry Smith     if (dd->n < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
626aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
62747c6ae99SBarry Smith     if (!dd->ly) {
628785e854fSJed Brown       ierr = PetscMalloc1(dd->n, &dd->ly);CHKERRQ(ierr);
62947c6ae99SBarry Smith     }
63047c6ae99SBarry Smith     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
63147c6ae99SBarry Smith   }
63247c6ae99SBarry Smith   if (lz) {
633ce94432eSBarry Smith     if (dd->p < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
634aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
63547c6ae99SBarry Smith     if (!dd->lz) {
636785e854fSJed Brown       ierr = PetscMalloc1(dd->p, &dd->lz);CHKERRQ(ierr);
63747c6ae99SBarry Smith     }
63847c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
63947c6ae99SBarry Smith   }
64047c6ae99SBarry Smith   PetscFunctionReturn(0);
64147c6ae99SBarry Smith }
64247c6ae99SBarry Smith 
64347c6ae99SBarry Smith #undef __FUNCT__
644aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType"
64547c6ae99SBarry Smith /*@
646aa219208SBarry Smith        DMDASetInterpolationType - Sets the type of interpolation that will be
647e727c939SJed Brown           returned by DMCreateInterpolation()
64847c6ae99SBarry Smith 
649aa219208SBarry Smith    Logically Collective on DMDA
65047c6ae99SBarry Smith 
65147c6ae99SBarry Smith    Input Parameter:
65247c6ae99SBarry Smith +  da - initial distributed array
653aa219208SBarry Smith .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
65447c6ae99SBarry Smith 
65547c6ae99SBarry Smith    Level: intermediate
65647c6ae99SBarry Smith 
657e727c939SJed Brown    Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation()
65847c6ae99SBarry Smith 
65947c6ae99SBarry Smith .keywords:  distributed array, interpolation
66047c6ae99SBarry Smith 
66196e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
66247c6ae99SBarry Smith @*/
6637087cfbeSBarry Smith PetscErrorCode  DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
66447c6ae99SBarry Smith {
66547c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith   PetscFunctionBegin;
66847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
66947c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,ctype,2);
67047c6ae99SBarry Smith   dd->interptype = ctype;
67147c6ae99SBarry Smith   PetscFunctionReturn(0);
67247c6ae99SBarry Smith }
67347c6ae99SBarry Smith 
6742dde6fd4SLisandro Dalcin #undef __FUNCT__
6752dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetInterpolationType"
6762dde6fd4SLisandro Dalcin /*@
6772dde6fd4SLisandro Dalcin        DMDAGetInterpolationType - Gets the type of interpolation that will be
678e727c939SJed Brown           used by DMCreateInterpolation()
6792dde6fd4SLisandro Dalcin 
6802dde6fd4SLisandro Dalcin    Not Collective
6812dde6fd4SLisandro Dalcin 
6822dde6fd4SLisandro Dalcin    Input Parameter:
6832dde6fd4SLisandro Dalcin .  da - distributed array
6842dde6fd4SLisandro Dalcin 
6852dde6fd4SLisandro Dalcin    Output Parameter:
6862dde6fd4SLisandro Dalcin .  ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms)
6872dde6fd4SLisandro Dalcin 
6882dde6fd4SLisandro Dalcin    Level: intermediate
6892dde6fd4SLisandro Dalcin 
6902dde6fd4SLisandro Dalcin .keywords:  distributed array, interpolation
6912dde6fd4SLisandro Dalcin 
692e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation()
6932dde6fd4SLisandro Dalcin @*/
6942dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype)
6952dde6fd4SLisandro Dalcin {
6962dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA*)da->data;
6972dde6fd4SLisandro Dalcin 
6982dde6fd4SLisandro Dalcin   PetscFunctionBegin;
6992dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
7002dde6fd4SLisandro Dalcin   PetscValidPointer(ctype,2);
7012dde6fd4SLisandro Dalcin   *ctype = dd->interptype;
7022dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
7032dde6fd4SLisandro Dalcin }
70447c6ae99SBarry Smith 
70547c6ae99SBarry Smith #undef __FUNCT__
706aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors"
70747c6ae99SBarry Smith /*@C
708aa219208SBarry Smith       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
70947c6ae99SBarry Smith         processes neighbors.
71047c6ae99SBarry Smith 
71147c6ae99SBarry Smith     Not Collective
71247c6ae99SBarry Smith 
71347c6ae99SBarry Smith    Input Parameter:
714aa219208SBarry Smith .     da - the DMDA object
71547c6ae99SBarry Smith 
71647c6ae99SBarry Smith    Output Parameters:
71747c6ae99SBarry Smith .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
71847c6ae99SBarry Smith               this process itself is in the list
71947c6ae99SBarry Smith 
72047c6ae99SBarry Smith    Notes: In 2d the array is of length 9, in 3d of length 27
72147c6ae99SBarry Smith           Not supported in 1d
722aa219208SBarry Smith           Do not free the array, it is freed when the DMDA is destroyed.
72347c6ae99SBarry Smith 
72447c6ae99SBarry Smith    Fortran Notes: In fortran you must pass in an array of the appropriate length.
72547c6ae99SBarry Smith 
72647c6ae99SBarry Smith    Level: intermediate
72747c6ae99SBarry Smith 
72847c6ae99SBarry Smith @*/
7297087cfbeSBarry Smith PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
73047c6ae99SBarry Smith {
73147c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
7325fd66863SKarl Rupp 
73347c6ae99SBarry Smith   PetscFunctionBegin;
73447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
73547c6ae99SBarry Smith   *ranks = dd->neighbors;
73647c6ae99SBarry Smith   PetscFunctionReturn(0);
73747c6ae99SBarry Smith }
73847c6ae99SBarry Smith 
73947c6ae99SBarry Smith #undef __FUNCT__
740aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges"
74147c6ae99SBarry Smith /*@C
742aa219208SBarry Smith       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
74347c6ae99SBarry Smith 
74447c6ae99SBarry Smith     Not Collective
74547c6ae99SBarry Smith 
74647c6ae99SBarry Smith    Input Parameter:
747aa219208SBarry Smith .     da - the DMDA object
74847c6ae99SBarry Smith 
74947c6ae99SBarry Smith    Output Parameter:
75047c6ae99SBarry Smith +     lx - ownership along x direction (optional)
75147c6ae99SBarry Smith .     ly - ownership along y direction (optional)
75247c6ae99SBarry Smith -     lz - ownership along z direction (optional)
75347c6ae99SBarry Smith 
75447c6ae99SBarry Smith    Level: intermediate
75547c6ae99SBarry Smith 
756aa219208SBarry Smith     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
75747c6ae99SBarry Smith 
75847c6ae99SBarry Smith     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
759aa219208SBarry Smith     eighth arguments from DMDAGetInfo()
76047c6ae99SBarry Smith 
76147c6ae99SBarry Smith      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
762aa219208SBarry Smith     DMDA they came from still exists (has not been destroyed).
76347c6ae99SBarry Smith 
764e3f3e4b6SBarry Smith     These numbers are NOT multiplied by the number of dof per node.
765e3f3e4b6SBarry Smith 
766aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
76747c6ae99SBarry Smith @*/
7687087cfbeSBarry Smith PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
76947c6ae99SBarry Smith {
77047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
77147c6ae99SBarry Smith 
77247c6ae99SBarry Smith   PetscFunctionBegin;
77347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
77447c6ae99SBarry Smith   if (lx) *lx = dd->lx;
77547c6ae99SBarry Smith   if (ly) *ly = dd->ly;
77647c6ae99SBarry Smith   if (lz) *lz = dd->lz;
77747c6ae99SBarry Smith   PetscFunctionReturn(0);
77847c6ae99SBarry Smith }
77947c6ae99SBarry Smith 
78047c6ae99SBarry Smith #undef __FUNCT__
781aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor"
78247c6ae99SBarry Smith /*@
783aa219208SBarry Smith      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
78447c6ae99SBarry Smith 
785aa219208SBarry Smith     Logically Collective on DMDA
78647c6ae99SBarry Smith 
78747c6ae99SBarry Smith   Input Parameters:
788aa219208SBarry Smith +    da - the DMDA object
78947c6ae99SBarry Smith .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
79047c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
79147c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
79247c6ae99SBarry Smith 
79347c6ae99SBarry Smith   Options Database:
79447c6ae99SBarry Smith +  -da_refine_x - refinement ratio in x direction
79547c6ae99SBarry Smith .  -da_refine_y - refinement ratio in y direction
79647c6ae99SBarry Smith -  -da_refine_z - refinement ratio in z direction
79747c6ae99SBarry Smith 
79847c6ae99SBarry Smith   Level: intermediate
79947c6ae99SBarry Smith 
80047c6ae99SBarry Smith     Notes: Pass PETSC_IGNORE to leave a value unchanged
80147c6ae99SBarry Smith 
802aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor()
80347c6ae99SBarry Smith @*/
8047087cfbeSBarry Smith PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
80547c6ae99SBarry Smith {
80647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
80747c6ae99SBarry Smith 
80847c6ae99SBarry Smith   PetscFunctionBegin;
80947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
81047c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_x,2);
81147c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_y,3);
81247c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_z,4);
81347c6ae99SBarry Smith 
81447c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
81547c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
81647c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
81747c6ae99SBarry Smith   PetscFunctionReturn(0);
81847c6ae99SBarry Smith }
81947c6ae99SBarry Smith 
82047c6ae99SBarry Smith #undef __FUNCT__
821aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor"
82247c6ae99SBarry Smith /*@C
823aa219208SBarry Smith      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
82447c6ae99SBarry Smith 
82547c6ae99SBarry Smith     Not Collective
82647c6ae99SBarry Smith 
82747c6ae99SBarry Smith   Input Parameter:
828aa219208SBarry Smith .    da - the DMDA object
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith   Output Parameters:
83147c6ae99SBarry Smith +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
83247c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
83347c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
83447c6ae99SBarry Smith 
83547c6ae99SBarry Smith   Level: intermediate
83647c6ae99SBarry Smith 
8370298fd71SBarry Smith     Notes: Pass NULL for values you do not need
83847c6ae99SBarry Smith 
839aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor()
84047c6ae99SBarry Smith @*/
8417087cfbeSBarry Smith PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
84247c6ae99SBarry Smith {
84347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
84447c6ae99SBarry Smith 
84547c6ae99SBarry Smith   PetscFunctionBegin;
84647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
84747c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
84847c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
84947c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
85047c6ae99SBarry Smith   PetscFunctionReturn(0);
85147c6ae99SBarry Smith }
85247c6ae99SBarry Smith 
85347c6ae99SBarry Smith #undef __FUNCT__
854aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix"
85547c6ae99SBarry Smith /*@C
856aa219208SBarry Smith      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
85747c6ae99SBarry Smith 
858aa219208SBarry Smith     Logically Collective on DMDA
85947c6ae99SBarry Smith 
86047c6ae99SBarry Smith   Input Parameters:
861aa219208SBarry Smith +    da - the DMDA object
862aa219208SBarry Smith -    f - the function that allocates the matrix for that specific DMDA
86347c6ae99SBarry Smith 
86447c6ae99SBarry Smith   Level: developer
86547c6ae99SBarry Smith 
866aa219208SBarry Smith    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
86747c6ae99SBarry Smith        the diagonal and off-diagonal blocks of the matrix
86847c6ae99SBarry Smith 
869950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills()
87047c6ae99SBarry Smith @*/
871b412c318SBarry Smith PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, Mat*))
87247c6ae99SBarry Smith {
87347c6ae99SBarry Smith   PetscFunctionBegin;
87447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
87525296bd5SBarry Smith   da->ops->creatematrix = f;
87647c6ae99SBarry Smith   PetscFunctionReturn(0);
87747c6ae99SBarry Smith }
87847c6ae99SBarry Smith 
87947c6ae99SBarry Smith #undef __FUNCT__
88094013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges"
88147c6ae99SBarry Smith /*
88247c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
88347c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
88447c6ae99SBarry Smith 
88547c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
88647c6ae99SBarry Smith */
88794013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
88847c6ae99SBarry Smith {
88947c6ae99SBarry Smith   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
89047c6ae99SBarry Smith   PetscErrorCode ierr;
89147c6ae99SBarry Smith 
89247c6ae99SBarry Smith   PetscFunctionBegin;
893ce94432eSBarry Smith   if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
89447c6ae99SBarry Smith   if (ratio == 1) {
89547c6ae99SBarry Smith     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
89647c6ae99SBarry Smith     PetscFunctionReturn(0);
89747c6ae99SBarry Smith   }
89847c6ae99SBarry Smith   for (i=0; i<m; i++) totalc += lc[i];
89947c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
90047c6ae99SBarry Smith   for (i=0; i<m; i++) {
90147c6ae99SBarry Smith     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
90247c6ae99SBarry Smith     if (i == m-1) lf[i] = want;
90347c6ae99SBarry Smith     else {
9047aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
9057aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
9067aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
9077aca7175SJed Brown       while ((startf+want)/ratio < nextc - stencil_width) want++;
9087aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
9097aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
9107aca7175SJed Brown       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
9117aca7175SJed Brown       /* Make sure all constraints are satisfied */
91230729d88SBarry Smith       if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width)
913ce94432eSBarry 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");
91447c6ae99SBarry Smith     }
91547c6ae99SBarry Smith     lf[i]      = want;
91647c6ae99SBarry Smith     startc    += lc[i];
91747c6ae99SBarry Smith     startf    += lf[i];
91847c6ae99SBarry Smith     remaining -= lf[i];
91947c6ae99SBarry Smith   }
92047c6ae99SBarry Smith   PetscFunctionReturn(0);
92147c6ae99SBarry Smith }
92247c6ae99SBarry Smith 
92347c6ae99SBarry Smith #undef __FUNCT__
9242be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges"
9252be375d4SJed Brown /*
9262be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
9272be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
9282be375d4SJed Brown 
9292be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
9302be375d4SJed Brown */
9312be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
9322be375d4SJed Brown {
9332be375d4SJed Brown   PetscInt       i,totalf,remaining,startc,startf;
9342be375d4SJed Brown   PetscErrorCode ierr;
9352be375d4SJed Brown 
9362be375d4SJed Brown   PetscFunctionBegin;
937ce94432eSBarry Smith   if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
9382be375d4SJed Brown   if (ratio == 1) {
9392be375d4SJed Brown     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
9402be375d4SJed Brown     PetscFunctionReturn(0);
9412be375d4SJed Brown   }
9422be375d4SJed Brown   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
9432be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
9442be375d4SJed Brown   for (i=0,startc=0,startf=0; i<m; i++) {
9452be375d4SJed Brown     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
9462be375d4SJed Brown     if (i == m-1) lc[i] = want;
9472be375d4SJed Brown     else {
9482be375d4SJed Brown       const PetscInt nextf = startf+lf[i];
9492be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
9502be375d4SJed Brown        * node is within one stencil width. */
9512be375d4SJed Brown       while (nextf/ratio < startc+want-stencil_width) want--;
9522be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
9532be375d4SJed Brown        * fine node is within one stencil width. */
9542be375d4SJed Brown       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
9552be375d4SJed Brown       if (want < 0 || want > remaining
956ce94432eSBarry 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");
9572be375d4SJed Brown     }
9582be375d4SJed Brown     lc[i]      = want;
9592be375d4SJed Brown     startc    += lc[i];
9602be375d4SJed Brown     startf    += lf[i];
9612be375d4SJed Brown     remaining -= lc[i];
9622be375d4SJed Brown   }
9632be375d4SJed Brown   PetscFunctionReturn(0);
9642be375d4SJed Brown }
9652be375d4SJed Brown 
9662be375d4SJed Brown #undef __FUNCT__
96794013140SBarry Smith #define __FUNCT__ "DMRefine_DA"
9687087cfbeSBarry Smith PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
96947c6ae99SBarry Smith {
97047c6ae99SBarry Smith   PetscErrorCode ierr;
971c73cfb54SMatthew G. Knepley   PetscInt       M,N,P,i,dim;
9729a42bb27SBarry Smith   DM             da2;
97347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
97447c6ae99SBarry Smith 
97547c6ae99SBarry Smith   PetscFunctionBegin;
97647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
97747c6ae99SBarry Smith   PetscValidPointer(daref,3);
97847c6ae99SBarry Smith 
979c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(da, &dim);CHKERRQ(ierr);
980bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
98147c6ae99SBarry Smith     M = dd->refine_x*dd->M;
98247c6ae99SBarry Smith   } else {
98347c6ae99SBarry Smith     M = 1 + dd->refine_x*(dd->M - 1);
98447c6ae99SBarry Smith   }
985bff4a2f0SMatthew G. Knepley   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
986c73cfb54SMatthew G. Knepley     if (dim > 1) {
98747c6ae99SBarry Smith       N = dd->refine_y*dd->N;
98847c6ae99SBarry Smith     } else {
9891860e6e9SBarry Smith       N = 1;
9901860e6e9SBarry Smith     }
9911860e6e9SBarry Smith   } else {
99247c6ae99SBarry Smith     N = 1 + dd->refine_y*(dd->N - 1);
99347c6ae99SBarry Smith   }
994bff4a2f0SMatthew G. Knepley   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
995c73cfb54SMatthew G. Knepley     if (dim > 2) {
99647c6ae99SBarry Smith       P = dd->refine_z*dd->P;
99747c6ae99SBarry Smith     } else {
9981860e6e9SBarry Smith       P = 1;
9991860e6e9SBarry Smith     }
10001860e6e9SBarry Smith   } else {
100147c6ae99SBarry Smith     P = 1 + dd->refine_z*(dd->P - 1);
100247c6ae99SBarry Smith   }
1003ce94432eSBarry Smith   ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr);
1004397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
1005c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(da2,dim);CHKERRQ(ierr);
1006397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
1007397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
1008397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
1009397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
1010397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
1011397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
1012c73cfb54SMatthew G. Knepley   if (dim == 3) {
101347c6ae99SBarry Smith     PetscInt *lx,*ly,*lz;
1014dcca6d9dSJed Brown     ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr);
1015bff4a2f0SMatthew 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);
1016bff4a2f0SMatthew 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);
1017bff4a2f0SMatthew 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);
1018397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
101947c6ae99SBarry Smith     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
1020c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
102147c6ae99SBarry Smith     PetscInt *lx,*ly;
1022dcca6d9dSJed Brown     ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr);
1023bff4a2f0SMatthew 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);
1024bff4a2f0SMatthew 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);
10250298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr);
102647c6ae99SBarry Smith     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
1027c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
102847c6ae99SBarry Smith     PetscInt *lx;
1029785e854fSJed Brown     ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr);
1030bff4a2f0SMatthew 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);
10310298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr);
103247c6ae99SBarry Smith     ierr = PetscFree(lx);CHKERRQ(ierr);
103347c6ae99SBarry Smith   }
103447c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
103547c6ae99SBarry Smith 
103647c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
103725296bd5SBarry Smith   da2->ops->creatematrix = da->ops->creatematrix;
103825296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
103947c6ae99SBarry Smith   da2->ops->getcoloring = da->ops->getcoloring;
104047c6ae99SBarry Smith   dd2->interptype       = dd->interptype;
104147c6ae99SBarry Smith 
104247c6ae99SBarry Smith   /* copy fill information if given */
104347c6ae99SBarry Smith   if (dd->dfill) {
1044854ce69bSBarry Smith     ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr);
104547c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
104647c6ae99SBarry Smith   }
104747c6ae99SBarry Smith   if (dd->ofill) {
1048854ce69bSBarry Smith     ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr);
104947c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
105047c6ae99SBarry Smith   }
105147c6ae99SBarry Smith   /* copy the refine information */
1052397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
1053397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
1054397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
105547c6ae99SBarry Smith 
1056897f7067SBarry Smith   if (dd->refine_z_hier) {
1057897f7067SBarry Smith     if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_z_hier_n) {
1058897f7067SBarry Smith       dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown + 1];
1059897f7067SBarry Smith     }
1060897f7067SBarry Smith     if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_z_hier_n) {
1061897f7067SBarry Smith       dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown];
1062897f7067SBarry Smith     }
1063897f7067SBarry Smith     dd2->refine_z_hier_n = dd->refine_z_hier_n;
1064897f7067SBarry Smith     ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr);
1065897f7067SBarry Smith     ierr = PetscMemcpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n*sizeof(PetscInt));CHKERRQ(ierr);
1066897f7067SBarry Smith   }
1067897f7067SBarry Smith   if (dd->refine_y_hier) {
1068897f7067SBarry Smith     if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_y_hier_n) {
1069897f7067SBarry Smith       dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown + 1];
1070897f7067SBarry Smith     }
1071897f7067SBarry Smith     if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_y_hier_n) {
1072897f7067SBarry Smith       dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown];
1073897f7067SBarry Smith     }
1074897f7067SBarry Smith     dd2->refine_y_hier_n = dd->refine_y_hier_n;
1075897f7067SBarry Smith     ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr);
1076897f7067SBarry Smith     ierr = PetscMemcpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n*sizeof(PetscInt));CHKERRQ(ierr);
1077897f7067SBarry Smith   }
1078897f7067SBarry Smith   if (dd->refine_x_hier) {
1079897f7067SBarry Smith     if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_x_hier_n) {
1080897f7067SBarry Smith       dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown + 1];
1081897f7067SBarry Smith     }
1082897f7067SBarry Smith     if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_x_hier_n) {
1083897f7067SBarry Smith       dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown];
1084897f7067SBarry Smith     }
1085897f7067SBarry Smith     dd2->refine_x_hier_n = dd->refine_x_hier_n;
1086897f7067SBarry Smith     ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr);
1087897f7067SBarry Smith     ierr = PetscMemcpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n*sizeof(PetscInt));CHKERRQ(ierr);
1088897f7067SBarry Smith   }
1089897f7067SBarry Smith 
1090897f7067SBarry Smith 
109147c6ae99SBarry Smith   /* copy vector type information */
1092c0dedaeaSBarry Smith   ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr);
1093ddcf8b74SDave May 
1094efd51863SBarry Smith   dd2->lf = dd->lf;
1095efd51863SBarry Smith   dd2->lj = dd->lj;
1096efd51863SBarry Smith 
10976e87535bSJed Brown   da2->leveldown = da->leveldown;
10986e87535bSJed Brown   da2->levelup   = da->levelup + 1;
10998865f1eaSKarl Rupp 
11006e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
11016e87535bSJed Brown 
1102ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
11036636e97aSMatthew G Knepley   if (da->coordinates) {
1104ddcf8b74SDave May     DM  cdaf,cdac;
1105ddcf8b74SDave May     Vec coordsc,coordsf;
1106ddcf8b74SDave May     Mat II;
1107ddcf8b74SDave May 
11086636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr);
11096636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr);
11106636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr);
1111b61d3410SDave May     /* force creation of the coordinate vector */
1112b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
11136636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
11140298fd71SBarry Smith     ierr = DMCreateInterpolation(cdac,cdaf,&II,NULL);CHKERRQ(ierr);
1115ddcf8b74SDave May     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
1116fcfd50ebSBarry Smith     ierr = MatDestroy(&II);CHKERRQ(ierr);
1117ddcf8b74SDave May   }
1118397b6216SJed Brown 
1119f3141302SJed Brown   for (i=0; i<da->bs; i++) {
1120f3141302SJed Brown     const char *fieldname;
1121f3141302SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
1122f3141302SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
1123f3141302SJed Brown   }
1124397b6216SJed Brown 
112547c6ae99SBarry Smith   *daref = da2;
112647c6ae99SBarry Smith   PetscFunctionReturn(0);
112747c6ae99SBarry Smith }
112847c6ae99SBarry Smith 
1129397b6216SJed Brown 
113047c6ae99SBarry Smith #undef __FUNCT__
113194013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA"
11327087cfbeSBarry Smith PetscErrorCode  DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
113347c6ae99SBarry Smith {
113447c6ae99SBarry Smith   PetscErrorCode ierr;
1135c73cfb54SMatthew G. Knepley   PetscInt       M,N,P,i,dim;
11369a42bb27SBarry Smith   DM             da2;
113747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
113847c6ae99SBarry Smith 
113947c6ae99SBarry Smith   PetscFunctionBegin;
114047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
114147c6ae99SBarry Smith   PetscValidPointer(daref,3);
114247c6ae99SBarry Smith 
1143c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(da, &dim);CHKERRQ(ierr);
1144bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1145397b6216SJed Brown     M = dd->M / dd->coarsen_x;
114647c6ae99SBarry Smith   } else {
1147397b6216SJed Brown     M = 1 + (dd->M - 1) / dd->coarsen_x;
114847c6ae99SBarry Smith   }
1149bff4a2f0SMatthew G. Knepley   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1150c73cfb54SMatthew G. Knepley     if (dim > 1) {
1151397b6216SJed Brown       N = dd->N / dd->coarsen_y;
115247c6ae99SBarry Smith     } else {
11531860e6e9SBarry Smith       N = 1;
11541860e6e9SBarry Smith     }
11551860e6e9SBarry Smith   } else {
1156397b6216SJed Brown     N = 1 + (dd->N - 1) / dd->coarsen_y;
115747c6ae99SBarry Smith   }
1158bff4a2f0SMatthew G. Knepley   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1159c73cfb54SMatthew G. Knepley     if (dim > 2) {
1160397b6216SJed Brown       P = dd->P / dd->coarsen_z;
116147c6ae99SBarry Smith     } else {
11621860e6e9SBarry Smith       P = 1;
11631860e6e9SBarry Smith     }
11641860e6e9SBarry Smith   } else {
1165397b6216SJed Brown     P = 1 + (dd->P - 1) / dd->coarsen_z;
116647c6ae99SBarry Smith   }
1167ce94432eSBarry Smith   ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr);
1168397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
1169c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(da2,dim);CHKERRQ(ierr);
1170397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
1171397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
1172397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
1173397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
1174397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
1175397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
1176c73cfb54SMatthew G. Knepley   if (dim == 3) {
11772be375d4SJed Brown     PetscInt *lx,*ly,*lz;
1178dcca6d9dSJed Brown     ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr);
1179bff4a2f0SMatthew 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);
1180bff4a2f0SMatthew 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);
1181bff4a2f0SMatthew 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);
1182397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
11832be375d4SJed Brown     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
1184c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
11852be375d4SJed Brown     PetscInt *lx,*ly;
1186dcca6d9dSJed Brown     ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr);
1187bff4a2f0SMatthew 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);
1188bff4a2f0SMatthew 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);
11890298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr);
11902be375d4SJed Brown     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
1191c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
11922be375d4SJed Brown     PetscInt *lx;
1193785e854fSJed Brown     ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr);
1194bff4a2f0SMatthew 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);
11950298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr);
11962be375d4SJed Brown     ierr = PetscFree(lx);CHKERRQ(ierr);
119747c6ae99SBarry Smith   }
119847c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
119947c6ae99SBarry Smith 
12004dcab191SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
120125296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */
120225296bd5SBarry Smith   da2->ops->creatematrix = da->ops->creatematrix;
120347c6ae99SBarry Smith   da2->ops->getcoloring  = da->ops->getcoloring;
120447c6ae99SBarry Smith   dd2->interptype        = dd->interptype;
120547c6ae99SBarry Smith 
120647c6ae99SBarry Smith   /* copy fill information if given */
120747c6ae99SBarry Smith   if (dd->dfill) {
1208854ce69bSBarry Smith     ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr);
120947c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
121047c6ae99SBarry Smith   }
121147c6ae99SBarry Smith   if (dd->ofill) {
1212854ce69bSBarry Smith     ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr);
121347c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
121447c6ae99SBarry Smith   }
121547c6ae99SBarry Smith   /* copy the refine information */
1216397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
1217397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
1218397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
121947c6ae99SBarry Smith 
1220897f7067SBarry Smith   if (dd->refine_z_hier) {
1221897f7067SBarry Smith     if (da->levelup - da->leveldown -1 > -1 && da->levelup - da->leveldown - 1< dd->refine_z_hier_n) {
1222897f7067SBarry Smith       dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown - 1];
1223897f7067SBarry Smith     }
1224897f7067SBarry Smith     if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_z_hier_n) {
1225897f7067SBarry Smith       dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown - 2];
1226897f7067SBarry Smith     }
1227897f7067SBarry Smith     dd2->refine_z_hier_n = dd->refine_z_hier_n;
1228897f7067SBarry Smith     ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr);
1229897f7067SBarry Smith     ierr = PetscMemcpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n*sizeof(PetscInt));CHKERRQ(ierr);
1230897f7067SBarry Smith   }
1231897f7067SBarry Smith   if (dd->refine_y_hier) {
1232897f7067SBarry Smith     if (da->levelup - da->leveldown - 1 > -1 && da->levelup - da->leveldown - 1< dd->refine_y_hier_n) {
1233897f7067SBarry Smith       dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown - 1];
1234897f7067SBarry Smith     }
1235897f7067SBarry Smith     if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_y_hier_n) {
1236897f7067SBarry Smith       dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown - 2];
1237897f7067SBarry Smith     }
1238897f7067SBarry Smith     dd2->refine_y_hier_n = dd->refine_y_hier_n;
1239897f7067SBarry Smith     ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr);
1240897f7067SBarry Smith     ierr = PetscMemcpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n*sizeof(PetscInt));CHKERRQ(ierr);
1241897f7067SBarry Smith   }
1242897f7067SBarry Smith   if (dd->refine_x_hier) {
1243897f7067SBarry Smith     if (da->levelup - da->leveldown - 1 > -1 && da->levelup - da->leveldown - 1 < dd->refine_x_hier_n) {
1244897f7067SBarry Smith       dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown - 1];
1245897f7067SBarry Smith     }
1246897f7067SBarry Smith     if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_x_hier_n) {
1247897f7067SBarry Smith       dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown - 2];
1248897f7067SBarry Smith     }
1249897f7067SBarry Smith     dd2->refine_x_hier_n = dd->refine_x_hier_n;
1250897f7067SBarry Smith     ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr);
1251897f7067SBarry Smith     ierr = PetscMemcpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n*sizeof(PetscInt));CHKERRQ(ierr);
1252897f7067SBarry Smith   }
1253897f7067SBarry Smith 
125447c6ae99SBarry Smith   /* copy vector type information */
1255c0dedaeaSBarry Smith   ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr);
125647c6ae99SBarry Smith 
1257644e2e5bSBarry Smith   dd2->lf = dd->lf;
1258644e2e5bSBarry Smith   dd2->lj = dd->lj;
1259644e2e5bSBarry Smith 
12606e87535bSJed Brown   da2->leveldown = da->leveldown + 1;
12616e87535bSJed Brown   da2->levelup   = da->levelup;
12628865f1eaSKarl Rupp 
12636e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
12646e87535bSJed Brown 
1265ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
12666636e97aSMatthew G Knepley   if (da->coordinates) {
1267ddcf8b74SDave May     DM         cdaf,cdac;
1268ddcf8b74SDave May     Vec        coordsc,coordsf;
12696dbf9973SLawrence Mitchell     Mat        inject;
12706dbf9973SLawrence Mitchell     VecScatter vscat;
1271ddcf8b74SDave May 
12726636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr);
12736636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr);
12746636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr);
1275b61d3410SDave May     /* force creation of the coordinate vector */
1276b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
12776636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
1278ddcf8b74SDave May 
1279e727c939SJed Brown     ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
12806dbf9973SLawrence Mitchell     ierr = MatScatterGetVecScatter(inject,&vscat);CHKERRQ(ierr);
12816dbf9973SLawrence Mitchell     ierr = VecScatterBegin(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12826dbf9973SLawrence Mitchell     ierr = VecScatterEnd(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12836dbf9973SLawrence Mitchell     ierr = MatDestroy(&inject);CHKERRQ(ierr);
1284ddcf8b74SDave May   }
1285f98405f7SJed Brown 
1286f98405f7SJed Brown   for (i=0; i<da->bs; i++) {
1287f98405f7SJed Brown     const char *fieldname;
1288f98405f7SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
1289f98405f7SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
1290f98405f7SJed Brown   }
1291f98405f7SJed Brown 
129247c6ae99SBarry Smith   *daref = da2;
129347c6ae99SBarry Smith   PetscFunctionReturn(0);
129447c6ae99SBarry Smith }
129547c6ae99SBarry Smith 
129647c6ae99SBarry Smith #undef __FUNCT__
129794013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA"
12987087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
129947c6ae99SBarry Smith {
130047c6ae99SBarry Smith   PetscErrorCode ierr;
130147c6ae99SBarry Smith   PetscInt       i,n,*refx,*refy,*refz;
130247c6ae99SBarry Smith 
130347c6ae99SBarry Smith   PetscFunctionBegin;
130447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1305ce94432eSBarry Smith   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
130647c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
130747c6ae99SBarry Smith   PetscValidPointer(daf,3);
130847c6ae99SBarry Smith 
1309aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
1310dcca6d9dSJed Brown   ierr = PetscMalloc3(nlevels,&refx,nlevels,&refy,nlevels,&refz);CHKERRQ(ierr);
131147c6ae99SBarry Smith   for (i=0; i<nlevels; i++) {
1312aa219208SBarry Smith     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
131347c6ae99SBarry Smith   }
131447c6ae99SBarry Smith   n    = nlevels;
1315c5929fdfSBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr);
131647c6ae99SBarry Smith   n    = nlevels;
1317c5929fdfSBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr);
131847c6ae99SBarry Smith   n    = nlevels;
1319c5929fdfSBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr);
132047c6ae99SBarry Smith 
1321aa219208SBarry Smith   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
1322ce94432eSBarry Smith   ierr = DMRefine(da,PetscObjectComm((PetscObject)da),&daf[0]);CHKERRQ(ierr);
132347c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
1324aa219208SBarry Smith     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
1325ce94432eSBarry Smith     ierr = DMRefine(daf[i-1],PetscObjectComm((PetscObject)da),&daf[i]);CHKERRQ(ierr);
132647c6ae99SBarry Smith   }
132747c6ae99SBarry Smith   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
132847c6ae99SBarry Smith   PetscFunctionReturn(0);
132947c6ae99SBarry Smith }
133047c6ae99SBarry Smith 
133147c6ae99SBarry Smith #undef __FUNCT__
133294013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA"
13337087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
133447c6ae99SBarry Smith {
133547c6ae99SBarry Smith   PetscErrorCode ierr;
133647c6ae99SBarry Smith   PetscInt       i;
133747c6ae99SBarry Smith 
133847c6ae99SBarry Smith   PetscFunctionBegin;
133947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1340ce94432eSBarry Smith   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
134147c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
134247c6ae99SBarry Smith   PetscValidPointer(dac,3);
1343ce94432eSBarry Smith   ierr = DMCoarsen(da,PetscObjectComm((PetscObject)da),&dac[0]);CHKERRQ(ierr);
134447c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
1345ce94432eSBarry Smith     ierr = DMCoarsen(dac[i-1],PetscObjectComm((PetscObject)da),&dac[i]);CHKERRQ(ierr);
134647c6ae99SBarry Smith   }
134747c6ae99SBarry Smith   PetscFunctionReturn(0);
134847c6ae99SBarry Smith }
1349