xref: /petsc/src/dm/impls/da/da.c (revision ca266f367873f97695dfce72da4e09458ac9ba8a)
1b45d2f2cSJed Brown #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
247c6ae99SBarry Smith 
347c6ae99SBarry Smith 
447c6ae99SBarry Smith #undef __FUNCT__
5aa219208SBarry Smith #define __FUNCT__ "DMDASetDim"
647c6ae99SBarry Smith /*@
7aa219208SBarry Smith   DMDASetDim - Sets the dimension
847c6ae99SBarry Smith 
9aa219208SBarry Smith   Collective on DMDA
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith   Input Parameters:
12aa219208SBarry Smith + da - the DMDA
1347c6ae99SBarry Smith - dim - the dimension (or PETSC_DECIDE)
1447c6ae99SBarry Smith 
1547c6ae99SBarry Smith   Level: intermediate
1647c6ae99SBarry Smith 
17aa219208SBarry Smith .seealso: DaGetDim(), DMDASetSizes()
1847c6ae99SBarry Smith @*/
197087cfbeSBarry Smith PetscErrorCode  DMDASetDim(DM da, PetscInt dim)
2047c6ae99SBarry Smith {
2147c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
2247c6ae99SBarry Smith 
2347c6ae99SBarry Smith   PetscFunctionBegin;
2447c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
25aa219208SBarry Smith   if (dd->dim > 0 && dim != dd->dim) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change DMDA dim from %D after it was set to %D",dd->dim,dim);
2647c6ae99SBarry Smith   dd->dim = dim;
2747c6ae99SBarry Smith   PetscFunctionReturn(0);
2847c6ae99SBarry Smith }
2947c6ae99SBarry Smith 
3047c6ae99SBarry Smith #undef __FUNCT__
31aa219208SBarry Smith #define __FUNCT__ "DMDASetSizes"
3247c6ae99SBarry Smith /*@
33aa219208SBarry Smith   DMDASetSizes - Sets the global sizes
3447c6ae99SBarry Smith 
35aa219208SBarry Smith   Logically Collective on DMDA
3647c6ae99SBarry Smith 
3747c6ae99SBarry Smith   Input Parameters:
38aa219208SBarry Smith + da - the DMDA
3947c6ae99SBarry Smith . M - the global X size (or PETSC_DECIDE)
4047c6ae99SBarry Smith . N - the global Y size (or PETSC_DECIDE)
4147c6ae99SBarry Smith - P - the global Z size (or PETSC_DECIDE)
4247c6ae99SBarry Smith 
4347c6ae99SBarry Smith   Level: intermediate
4447c6ae99SBarry Smith 
45aa219208SBarry Smith .seealso: DMDAGetSize(), PetscSplitOwnership()
4647c6ae99SBarry Smith @*/
477087cfbeSBarry Smith PetscErrorCode  DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P)
4847c6ae99SBarry Smith {
4947c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
5047c6ae99SBarry Smith 
5147c6ae99SBarry Smith   PetscFunctionBegin;
5247c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
5347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,M,2);
5447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,N,3);
5547c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,P,4);
56cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
5747c6ae99SBarry Smith 
5847c6ae99SBarry Smith   dd->M = M;
5947c6ae99SBarry Smith   dd->N = N;
6047c6ae99SBarry Smith   dd->P = P;
6147c6ae99SBarry Smith   PetscFunctionReturn(0);
6247c6ae99SBarry Smith }
6347c6ae99SBarry Smith 
6447c6ae99SBarry Smith #undef __FUNCT__
65aa219208SBarry Smith #define __FUNCT__ "DMDASetNumProcs"
6647c6ae99SBarry Smith /*@
67aa219208SBarry Smith   DMDASetNumProcs - Sets the number of processes in each dimension
6847c6ae99SBarry Smith 
69aa219208SBarry Smith   Logically Collective on DMDA
7047c6ae99SBarry Smith 
7147c6ae99SBarry Smith   Input Parameters:
72aa219208SBarry Smith + da - the DMDA
7347c6ae99SBarry Smith . m - the number of X procs (or PETSC_DECIDE)
7447c6ae99SBarry Smith . n - the number of Y procs (or PETSC_DECIDE)
7547c6ae99SBarry Smith - p - the number of Z procs (or PETSC_DECIDE)
7647c6ae99SBarry Smith 
7747c6ae99SBarry Smith   Level: intermediate
7847c6ae99SBarry Smith 
79aa219208SBarry Smith .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership()
8047c6ae99SBarry Smith @*/
817087cfbeSBarry Smith PetscErrorCode  DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p)
8247c6ae99SBarry Smith {
8347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
8447c6ae99SBarry Smith 
8547c6ae99SBarry Smith   PetscFunctionBegin;
8647c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
8747c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,m,2);
8847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,n,3);
8947c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,p,4);
90cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
9147c6ae99SBarry Smith   dd->m = m;
9247c6ae99SBarry Smith   dd->n = n;
9347c6ae99SBarry Smith   dd->p = p;
9447c6ae99SBarry Smith   PetscFunctionReturn(0);
9547c6ae99SBarry Smith }
9647c6ae99SBarry Smith 
9747c6ae99SBarry Smith #undef __FUNCT__
983e24977eSEthan Coon #define __FUNCT__ "DMDASetBoundaryType"
9947c6ae99SBarry Smith /*@
1001321219cSEthan Coon   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries.
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith   Not collective
10347c6ae99SBarry Smith 
10447c6ae99SBarry Smith   Input Parameter:
105aa219208SBarry Smith + da    - The DMDA
1061321219cSEthan Coon - bx,by,bz - One of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC
10747c6ae99SBarry Smith 
10847c6ae99SBarry Smith   Level: intermediate
10947c6ae99SBarry Smith 
11047c6ae99SBarry Smith .keywords:  distributed array, periodicity
11196e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA, DMDABoundaryType
11247c6ae99SBarry Smith @*/
1131321219cSEthan Coon PetscErrorCode  DMDASetBoundaryType(DM da,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz)
11447c6ae99SBarry Smith {
11547c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith   PetscFunctionBegin;
11847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1195a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,bx,2);
1205a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,by,3);
1215a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,bz,4);
122cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1231321219cSEthan Coon   dd->bx = bx;
1241321219cSEthan Coon   dd->by = by;
1251321219cSEthan Coon   dd->bz = bz;
12647c6ae99SBarry Smith   PetscFunctionReturn(0);
12747c6ae99SBarry Smith }
12847c6ae99SBarry Smith 
12947c6ae99SBarry Smith #undef __FUNCT__
130aa219208SBarry Smith #define __FUNCT__ "DMDASetDof"
13147c6ae99SBarry Smith /*@
132aa219208SBarry Smith   DMDASetDof - Sets the number of degrees of freedom per vertex
13347c6ae99SBarry Smith 
13447c6ae99SBarry Smith   Not collective
13547c6ae99SBarry Smith 
13647c6ae99SBarry Smith   Input Parameter:
137aa219208SBarry Smith + da  - The DMDA
13847c6ae99SBarry Smith - dof - Number of degrees of freedom
13947c6ae99SBarry Smith 
14047c6ae99SBarry Smith   Level: intermediate
14147c6ae99SBarry Smith 
14247c6ae99SBarry Smith .keywords:  distributed array, degrees of freedom
14396e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
14447c6ae99SBarry Smith @*/
14554cfb0beSLisandro Dalcin PetscErrorCode  DMDASetDof(DM da, PetscInt dof)
14647c6ae99SBarry Smith {
14747c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
14847c6ae99SBarry Smith 
14947c6ae99SBarry Smith   PetscFunctionBegin;
15047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
15154cfb0beSLisandro Dalcin   PetscValidLogicalCollectiveInt(da,dof,2);
152cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
15347c6ae99SBarry Smith   dd->w = dof;
1541411c6eeSJed Brown   da->bs = dof;
15547c6ae99SBarry Smith   PetscFunctionReturn(0);
15647c6ae99SBarry Smith }
15747c6ae99SBarry Smith 
15847c6ae99SBarry Smith #undef __FUNCT__
15988661749SPeter Brune #define __FUNCT__ "DMDASetOverlap"
16088661749SPeter Brune /*@
16188661749SPeter Brune   DMDASetOverlap - Sets the size of the per-processor overlap.
16288661749SPeter Brune 
16388661749SPeter Brune   Not collective
16488661749SPeter Brune 
16588661749SPeter Brune   Input Parameter:
16688661749SPeter Brune + da  - The DMDA
16788661749SPeter Brune - dof - Number of degrees of freedom
16888661749SPeter Brune 
16988661749SPeter Brune   Level: intermediate
17088661749SPeter Brune 
17188661749SPeter Brune .keywords:  distributed array, degrees of freedom
17288661749SPeter Brune .seealso: DMDACreate(), DMDestroy(), DMDA
17388661749SPeter Brune @*/
17488661749SPeter Brune PetscErrorCode  DMDASetOverlap(DM da, PetscInt overlap)
17588661749SPeter Brune {
17688661749SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
17788661749SPeter Brune 
17888661749SPeter Brune   PetscFunctionBegin;
17988661749SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
18088661749SPeter Brune   PetscValidLogicalCollectiveInt(da,overlap,2);
18188661749SPeter Brune   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
18288661749SPeter Brune   dd->overlap = overlap;
18388661749SPeter Brune   PetscFunctionReturn(0);
18488661749SPeter Brune }
18588661749SPeter Brune 
186d886c4f4SPeter Brune #undef __FUNCT__
187d886c4f4SPeter Brune #define __FUNCT__ "DMDASetOffset"
188d886c4f4SPeter Brune /*@
189d886c4f4SPeter Brune   DMDASetOffset - Sets the index offset of the DA.
190d886c4f4SPeter Brune 
191d886c4f4SPeter Brune   Collective on DA
192d886c4f4SPeter Brune 
193d886c4f4SPeter Brune   Input Parameter:
194d886c4f4SPeter Brune + da  - The DMDA
195d886c4f4SPeter Brune . xo  - The offset in the x direction
196d886c4f4SPeter Brune . yo  - The offset in the y direction
197d886c4f4SPeter Brune - zo  - The offset in the z direction
198d886c4f4SPeter Brune 
199d886c4f4SPeter Brune   Level: intermediate
200d886c4f4SPeter Brune 
201d886c4f4SPeter Brune   Notes: This is used primarily to overlap a computation on a local DA with that on a global DA without
202d886c4f4SPeter Brune   changing boundary conditions or subdomain features that depend upon the global offsets.
203d886c4f4SPeter Brune 
204d886c4f4SPeter Brune .keywords:  distributed array, degrees of freedom
205d886c4f4SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray()
206d886c4f4SPeter Brune @*/
207d886c4f4SPeter Brune PetscErrorCode  DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo)
208d886c4f4SPeter Brune {
209d886c4f4SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
210d886c4f4SPeter Brune 
211d886c4f4SPeter Brune   PetscFunctionBegin;
212d886c4f4SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
213d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da,xo,2);
214d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da,yo,2);
215d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da,zo,2);
216d886c4f4SPeter Brune   dd->xo = xo;
217d886c4f4SPeter Brune   dd->yo = yo;
218d886c4f4SPeter Brune   dd->zo = zo;
219d886c4f4SPeter Brune   PetscFunctionReturn(0);
220d886c4f4SPeter Brune }
221d886c4f4SPeter Brune 
222d886c4f4SPeter Brune #undef __FUNCT__
223d886c4f4SPeter Brune #define __FUNCT__ "DMDAGetOffset"
224d886c4f4SPeter Brune /*@
225d886c4f4SPeter Brune   DMDAGetOffset - Gets the index offset of the DA.
226d886c4f4SPeter Brune 
227d886c4f4SPeter Brune   Not collective
228d886c4f4SPeter Brune 
229d886c4f4SPeter Brune   Input Parameter:
230d886c4f4SPeter Brune . da  - The DMDA
231d886c4f4SPeter Brune 
232d886c4f4SPeter Brune   Output Parameters:
233d886c4f4SPeter Brune + xo  - The offset in the x direction
234d886c4f4SPeter Brune . yo  - The offset in the y direction
235d886c4f4SPeter Brune - zo  - The offset in the z direction
236d886c4f4SPeter Brune 
237d886c4f4SPeter Brune   Level: intermediate
238d886c4f4SPeter Brune 
239d886c4f4SPeter Brune .keywords:  distributed array, degrees of freedom
240d886c4f4SPeter Brune .seealso: DMDASetOffset(), DMDAVecGetArray()
241d886c4f4SPeter Brune @*/
242d886c4f4SPeter Brune PetscErrorCode  DMDAGetOffset(DM da, PetscInt *xo, PetscInt *yo, PetscInt *zo)
243d886c4f4SPeter Brune {
244d886c4f4SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
245d886c4f4SPeter Brune 
246d886c4f4SPeter Brune   PetscFunctionBegin;
247d886c4f4SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
248d886c4f4SPeter Brune   if (xo) *xo = dd->xo;
249d886c4f4SPeter Brune   if (yo) *yo = dd->yo;
250d886c4f4SPeter Brune   if (zo) *zo = dd->zo;
251d886c4f4SPeter Brune   PetscFunctionReturn(0);
252d886c4f4SPeter Brune }
253d886c4f4SPeter Brune 
25488661749SPeter Brune 
25588661749SPeter Brune #undef __FUNCT__
256aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType"
25747c6ae99SBarry Smith /*@
258aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
25947c6ae99SBarry Smith 
260aa219208SBarry Smith   Logically Collective on DMDA
26147c6ae99SBarry Smith 
26247c6ae99SBarry Smith   Input Parameter:
263aa219208SBarry Smith + da    - The DMDA
264aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
26547c6ae99SBarry Smith 
26647c6ae99SBarry Smith   Level: intermediate
26747c6ae99SBarry Smith 
26847c6ae99SBarry Smith .keywords:  distributed array, stencil
26996e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
27047c6ae99SBarry Smith @*/
2717087cfbeSBarry Smith PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
27247c6ae99SBarry Smith {
27347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
27447c6ae99SBarry Smith 
27547c6ae99SBarry Smith   PetscFunctionBegin;
27647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
27747c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,stype,2);
278cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
27947c6ae99SBarry Smith   dd->stencil_type = stype;
28047c6ae99SBarry Smith   PetscFunctionReturn(0);
28147c6ae99SBarry Smith }
28247c6ae99SBarry Smith 
28347c6ae99SBarry Smith #undef __FUNCT__
284aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth"
28547c6ae99SBarry Smith /*@
286aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
28747c6ae99SBarry Smith 
288aa219208SBarry Smith   Logically Collective on DMDA
28947c6ae99SBarry Smith 
29047c6ae99SBarry Smith   Input Parameter:
291aa219208SBarry Smith + da    - The DMDA
29247c6ae99SBarry Smith - width - The stencil width
29347c6ae99SBarry Smith 
29447c6ae99SBarry Smith   Level: intermediate
29547c6ae99SBarry Smith 
29647c6ae99SBarry Smith .keywords:  distributed array, stencil
29796e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
29847c6ae99SBarry Smith @*/
2997087cfbeSBarry Smith PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
30047c6ae99SBarry Smith {
30147c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
30247c6ae99SBarry Smith 
30347c6ae99SBarry Smith   PetscFunctionBegin;
30447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
30547c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,width,2);
306cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
30747c6ae99SBarry Smith   dd->s = width;
30847c6ae99SBarry Smith   PetscFunctionReturn(0);
30947c6ae99SBarry Smith }
31047c6ae99SBarry Smith 
31147c6ae99SBarry Smith #undef __FUNCT__
312aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
313aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
31447c6ae99SBarry Smith {
31547c6ae99SBarry Smith   PetscInt i,sum;
31647c6ae99SBarry Smith 
31747c6ae99SBarry Smith   PetscFunctionBegin;
31847c6ae99SBarry Smith   if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
31947c6ae99SBarry Smith   for (i=sum=0; i<m; i++) sum += lx[i];
32047c6ae99SBarry Smith   if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
32147c6ae99SBarry Smith   PetscFunctionReturn(0);
32247c6ae99SBarry Smith }
32347c6ae99SBarry Smith 
32447c6ae99SBarry Smith #undef __FUNCT__
325aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges"
32647c6ae99SBarry Smith /*@
327aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
32847c6ae99SBarry Smith 
329aa219208SBarry Smith   Logically Collective on DMDA
33047c6ae99SBarry Smith 
33147c6ae99SBarry Smith   Input Parameter:
332aa219208SBarry Smith + da - The DMDA
33347c6ae99SBarry Smith . lx - array containing number of nodes in the X direction on each process, or PETSC_NULL. If non-null, must be of length da->m
33447c6ae99SBarry Smith . ly - array containing number of nodes in the Y direction on each process, or PETSC_NULL. If non-null, must be of length da->n
33547c6ae99SBarry Smith - lz - array containing number of nodes in the Z direction on each process, or PETSC_NULL. If non-null, must be of length da->p.
33647c6ae99SBarry Smith 
33747c6ae99SBarry Smith   Level: intermediate
33847c6ae99SBarry Smith 
33947c6ae99SBarry Smith .keywords:  distributed array
34096e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
34147c6ae99SBarry Smith @*/
3427087cfbeSBarry Smith PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
34347c6ae99SBarry Smith {
34447c6ae99SBarry Smith   PetscErrorCode ierr;
34547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
34647c6ae99SBarry Smith 
34747c6ae99SBarry Smith   PetscFunctionBegin;
34847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
349cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
35047c6ae99SBarry Smith   if (lx) {
35147c6ae99SBarry Smith     if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
352aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
35347c6ae99SBarry Smith     if (!dd->lx) {
35447c6ae99SBarry Smith       ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
35547c6ae99SBarry Smith     }
35647c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
35747c6ae99SBarry Smith   }
35847c6ae99SBarry Smith   if (ly) {
35947c6ae99SBarry Smith     if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
360aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
36147c6ae99SBarry Smith     if (!dd->ly) {
36247c6ae99SBarry Smith       ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
36347c6ae99SBarry Smith     }
36447c6ae99SBarry Smith     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
36547c6ae99SBarry Smith   }
36647c6ae99SBarry Smith   if (lz) {
36747c6ae99SBarry Smith     if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
368aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
36947c6ae99SBarry Smith     if (!dd->lz) {
37047c6ae99SBarry Smith       ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
37147c6ae99SBarry Smith     }
37247c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
37347c6ae99SBarry Smith   }
37447c6ae99SBarry Smith   PetscFunctionReturn(0);
37547c6ae99SBarry Smith }
37647c6ae99SBarry Smith 
37747c6ae99SBarry Smith #undef __FUNCT__
378aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType"
37947c6ae99SBarry Smith /*@
380aa219208SBarry Smith        DMDASetInterpolationType - Sets the type of interpolation that will be
381e727c939SJed Brown           returned by DMCreateInterpolation()
38247c6ae99SBarry Smith 
383aa219208SBarry Smith    Logically Collective on DMDA
38447c6ae99SBarry Smith 
38547c6ae99SBarry Smith    Input Parameter:
38647c6ae99SBarry Smith +  da - initial distributed array
387aa219208SBarry Smith .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
38847c6ae99SBarry Smith 
38947c6ae99SBarry Smith    Level: intermediate
39047c6ae99SBarry Smith 
391e727c939SJed Brown    Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation()
39247c6ae99SBarry Smith 
39347c6ae99SBarry Smith .keywords:  distributed array, interpolation
39447c6ae99SBarry Smith 
39596e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
39647c6ae99SBarry Smith @*/
3977087cfbeSBarry Smith PetscErrorCode  DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
39847c6ae99SBarry Smith {
39947c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
40047c6ae99SBarry Smith 
40147c6ae99SBarry Smith   PetscFunctionBegin;
40247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
40347c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,ctype,2);
40447c6ae99SBarry Smith   dd->interptype = ctype;
40547c6ae99SBarry Smith   PetscFunctionReturn(0);
40647c6ae99SBarry Smith }
40747c6ae99SBarry Smith 
4082dde6fd4SLisandro Dalcin #undef __FUNCT__
4092dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetInterpolationType"
4102dde6fd4SLisandro Dalcin /*@
4112dde6fd4SLisandro Dalcin        DMDAGetInterpolationType - Gets the type of interpolation that will be
412e727c939SJed Brown           used by DMCreateInterpolation()
4132dde6fd4SLisandro Dalcin 
4142dde6fd4SLisandro Dalcin    Not Collective
4152dde6fd4SLisandro Dalcin 
4162dde6fd4SLisandro Dalcin    Input Parameter:
4172dde6fd4SLisandro Dalcin .  da - distributed array
4182dde6fd4SLisandro Dalcin 
4192dde6fd4SLisandro Dalcin    Output Parameter:
4202dde6fd4SLisandro Dalcin .  ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms)
4212dde6fd4SLisandro Dalcin 
4222dde6fd4SLisandro Dalcin    Level: intermediate
4232dde6fd4SLisandro Dalcin 
4242dde6fd4SLisandro Dalcin .keywords:  distributed array, interpolation
4252dde6fd4SLisandro Dalcin 
426e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation()
4272dde6fd4SLisandro Dalcin @*/
4282dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype)
4292dde6fd4SLisandro Dalcin {
4302dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA*)da->data;
4312dde6fd4SLisandro Dalcin 
4322dde6fd4SLisandro Dalcin   PetscFunctionBegin;
4332dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
4342dde6fd4SLisandro Dalcin   PetscValidPointer(ctype,2);
4352dde6fd4SLisandro Dalcin   *ctype = dd->interptype;
4362dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
4372dde6fd4SLisandro Dalcin }
43847c6ae99SBarry Smith 
43947c6ae99SBarry Smith #undef __FUNCT__
440aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors"
44147c6ae99SBarry Smith /*@C
442aa219208SBarry Smith       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
44347c6ae99SBarry Smith         processes neighbors.
44447c6ae99SBarry Smith 
44547c6ae99SBarry Smith     Not Collective
44647c6ae99SBarry Smith 
44747c6ae99SBarry Smith    Input Parameter:
448aa219208SBarry Smith .     da - the DMDA object
44947c6ae99SBarry Smith 
45047c6ae99SBarry Smith    Output Parameters:
45147c6ae99SBarry Smith .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
45247c6ae99SBarry Smith               this process itself is in the list
45347c6ae99SBarry Smith 
45447c6ae99SBarry Smith    Notes: In 2d the array is of length 9, in 3d of length 27
45547c6ae99SBarry Smith           Not supported in 1d
456aa219208SBarry Smith           Do not free the array, it is freed when the DMDA is destroyed.
45747c6ae99SBarry Smith 
45847c6ae99SBarry Smith    Fortran Notes: In fortran you must pass in an array of the appropriate length.
45947c6ae99SBarry Smith 
46047c6ae99SBarry Smith    Level: intermediate
46147c6ae99SBarry Smith 
46247c6ae99SBarry Smith @*/
4637087cfbeSBarry Smith PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
46447c6ae99SBarry Smith {
46547c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
46647c6ae99SBarry Smith   PetscFunctionBegin;
46747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
46847c6ae99SBarry Smith   *ranks = dd->neighbors;
46947c6ae99SBarry Smith   PetscFunctionReturn(0);
47047c6ae99SBarry Smith }
47147c6ae99SBarry Smith 
47247c6ae99SBarry Smith #undef __FUNCT__
473aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges"
47447c6ae99SBarry Smith /*@C
475aa219208SBarry Smith       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
47647c6ae99SBarry Smith 
47747c6ae99SBarry Smith     Not Collective
47847c6ae99SBarry Smith 
47947c6ae99SBarry Smith    Input Parameter:
480aa219208SBarry Smith .     da - the DMDA object
48147c6ae99SBarry Smith 
48247c6ae99SBarry Smith    Output Parameter:
48347c6ae99SBarry Smith +     lx - ownership along x direction (optional)
48447c6ae99SBarry Smith .     ly - ownership along y direction (optional)
48547c6ae99SBarry Smith -     lz - ownership along z direction (optional)
48647c6ae99SBarry Smith 
48747c6ae99SBarry Smith    Level: intermediate
48847c6ae99SBarry Smith 
489aa219208SBarry Smith     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
49047c6ae99SBarry Smith 
49147c6ae99SBarry Smith     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
492aa219208SBarry Smith     eighth arguments from DMDAGetInfo()
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
495aa219208SBarry Smith     DMDA they came from still exists (has not been destroyed).
49647c6ae99SBarry Smith 
497aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
49847c6ae99SBarry Smith @*/
4997087cfbeSBarry Smith PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
50047c6ae99SBarry Smith {
50147c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
50247c6ae99SBarry Smith 
50347c6ae99SBarry Smith   PetscFunctionBegin;
50447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
50547c6ae99SBarry Smith   if (lx) *lx = dd->lx;
50647c6ae99SBarry Smith   if (ly) *ly = dd->ly;
50747c6ae99SBarry Smith   if (lz) *lz = dd->lz;
50847c6ae99SBarry Smith   PetscFunctionReturn(0);
50947c6ae99SBarry Smith }
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith #undef __FUNCT__
512aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor"
51347c6ae99SBarry Smith /*@
514aa219208SBarry Smith      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
51547c6ae99SBarry Smith 
516aa219208SBarry Smith     Logically Collective on DMDA
51747c6ae99SBarry Smith 
51847c6ae99SBarry Smith   Input Parameters:
519aa219208SBarry Smith +    da - the DMDA object
52047c6ae99SBarry Smith .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
52147c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
52247c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
52347c6ae99SBarry Smith 
52447c6ae99SBarry Smith   Options Database:
52547c6ae99SBarry Smith +  -da_refine_x - refinement ratio in x direction
52647c6ae99SBarry Smith .  -da_refine_y - refinement ratio in y direction
52747c6ae99SBarry Smith -  -da_refine_z - refinement ratio in z direction
52847c6ae99SBarry Smith 
52947c6ae99SBarry Smith   Level: intermediate
53047c6ae99SBarry Smith 
53147c6ae99SBarry Smith     Notes: Pass PETSC_IGNORE to leave a value unchanged
53247c6ae99SBarry Smith 
533aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor()
53447c6ae99SBarry Smith @*/
5357087cfbeSBarry Smith PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
53647c6ae99SBarry Smith {
53747c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
53847c6ae99SBarry Smith 
53947c6ae99SBarry Smith   PetscFunctionBegin;
54047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
54147c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_x,2);
54247c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_y,3);
54347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_z,4);
54447c6ae99SBarry Smith 
54547c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
54647c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
54747c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
54847c6ae99SBarry Smith   PetscFunctionReturn(0);
54947c6ae99SBarry Smith }
55047c6ae99SBarry Smith 
55147c6ae99SBarry Smith #undef __FUNCT__
552aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor"
55347c6ae99SBarry Smith /*@C
554aa219208SBarry Smith      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
55547c6ae99SBarry Smith 
55647c6ae99SBarry Smith     Not Collective
55747c6ae99SBarry Smith 
55847c6ae99SBarry Smith   Input Parameter:
559aa219208SBarry Smith .    da - the DMDA object
56047c6ae99SBarry Smith 
56147c6ae99SBarry Smith   Output Parameters:
56247c6ae99SBarry Smith +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
56347c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
56447c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
56547c6ae99SBarry Smith 
56647c6ae99SBarry Smith   Level: intermediate
56747c6ae99SBarry Smith 
56847c6ae99SBarry Smith     Notes: Pass PETSC_NULL for values you do not need
56947c6ae99SBarry Smith 
570aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor()
57147c6ae99SBarry Smith @*/
5727087cfbeSBarry Smith PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
57347c6ae99SBarry Smith {
57447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
57547c6ae99SBarry Smith 
57647c6ae99SBarry Smith   PetscFunctionBegin;
57747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
57847c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
57947c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
58047c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
58147c6ae99SBarry Smith   PetscFunctionReturn(0);
58247c6ae99SBarry Smith }
58347c6ae99SBarry Smith 
58447c6ae99SBarry Smith #undef __FUNCT__
585aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix"
58647c6ae99SBarry Smith /*@C
587aa219208SBarry Smith      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
58847c6ae99SBarry Smith 
589aa219208SBarry Smith     Logically Collective on DMDA
59047c6ae99SBarry Smith 
59147c6ae99SBarry Smith   Input Parameters:
592aa219208SBarry Smith +    da - the DMDA object
593aa219208SBarry Smith -    f - the function that allocates the matrix for that specific DMDA
59447c6ae99SBarry Smith 
59547c6ae99SBarry Smith   Level: developer
59647c6ae99SBarry Smith 
597aa219208SBarry Smith    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
59847c6ae99SBarry Smith        the diagonal and off-diagonal blocks of the matrix
59947c6ae99SBarry Smith 
600950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills()
60147c6ae99SBarry Smith @*/
60219fd82e9SBarry Smith PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, MatType,Mat*))
60347c6ae99SBarry Smith {
60447c6ae99SBarry Smith   PetscFunctionBegin;
60547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
60625296bd5SBarry Smith   da->ops->creatematrix = f;
60747c6ae99SBarry Smith   PetscFunctionReturn(0);
60847c6ae99SBarry Smith }
60947c6ae99SBarry Smith 
61047c6ae99SBarry Smith #undef __FUNCT__
61194013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges"
61247c6ae99SBarry Smith /*
61347c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
61447c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
61547c6ae99SBarry Smith 
61647c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
61747c6ae99SBarry Smith */
61894013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
61947c6ae99SBarry Smith {
62047c6ae99SBarry Smith   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
62147c6ae99SBarry Smith   PetscErrorCode ierr;
62247c6ae99SBarry Smith 
62347c6ae99SBarry Smith   PetscFunctionBegin;
62447c6ae99SBarry Smith   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
62547c6ae99SBarry Smith   if (ratio == 1) {
62647c6ae99SBarry Smith     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
62747c6ae99SBarry Smith     PetscFunctionReturn(0);
62847c6ae99SBarry Smith   }
62947c6ae99SBarry Smith   for (i=0; i<m; i++) totalc += lc[i];
63047c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
63147c6ae99SBarry Smith   for (i=0; i<m; i++) {
63247c6ae99SBarry Smith     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
63347c6ae99SBarry Smith     if (i == m-1) lf[i] = want;
63447c6ae99SBarry Smith     else {
6357aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
6367aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
6377aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
6387aca7175SJed Brown       while ((startf+want)/ratio < nextc - stencil_width) want++;
6397aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
6407aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
6417aca7175SJed Brown       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
6427aca7175SJed Brown       /* Make sure all constraints are satisfied */
64330729d88SBarry Smith       if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width)
64430729d88SBarry Smith           || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range");
64547c6ae99SBarry Smith     }
64647c6ae99SBarry Smith     lf[i] = want;
64747c6ae99SBarry Smith     startc += lc[i];
64847c6ae99SBarry Smith     startf += lf[i];
64947c6ae99SBarry Smith     remaining -= lf[i];
65047c6ae99SBarry Smith   }
65147c6ae99SBarry Smith   PetscFunctionReturn(0);
65247c6ae99SBarry Smith }
65347c6ae99SBarry Smith 
65447c6ae99SBarry Smith #undef __FUNCT__
6552be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges"
6562be375d4SJed Brown /*
6572be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
6582be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
6592be375d4SJed Brown 
6602be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
6612be375d4SJed Brown */
6622be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
6632be375d4SJed Brown {
6642be375d4SJed Brown   PetscInt       i,totalf,remaining,startc,startf;
6652be375d4SJed Brown   PetscErrorCode ierr;
6662be375d4SJed Brown 
6672be375d4SJed Brown   PetscFunctionBegin;
6682be375d4SJed Brown   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
6692be375d4SJed Brown   if (ratio == 1) {
6702be375d4SJed Brown     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
6712be375d4SJed Brown     PetscFunctionReturn(0);
6722be375d4SJed Brown   }
6732be375d4SJed Brown   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
6742be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
6752be375d4SJed Brown   for (i=0,startc=0,startf=0; i<m; i++) {
6762be375d4SJed Brown     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
6772be375d4SJed Brown     if (i == m-1) lc[i] = want;
6782be375d4SJed Brown     else {
6792be375d4SJed Brown       const PetscInt nextf = startf+lf[i];
6802be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
6812be375d4SJed Brown        * node is within one stencil width. */
6822be375d4SJed Brown       while (nextf/ratio < startc+want-stencil_width) want--;
6832be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
6842be375d4SJed Brown        * fine node is within one stencil width. */
6852be375d4SJed Brown       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
6862be375d4SJed Brown       if (want < 0 || want > remaining
68730729d88SBarry Smith           || (nextf/ratio < startc+want-stencil_width) || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range");
6882be375d4SJed Brown     }
6892be375d4SJed Brown     lc[i] = want;
6902be375d4SJed Brown     startc += lc[i];
6912be375d4SJed Brown     startf += lf[i];
6922be375d4SJed Brown     remaining -= lc[i];
6932be375d4SJed Brown   }
6942be375d4SJed Brown   PetscFunctionReturn(0);
6952be375d4SJed Brown }
6962be375d4SJed Brown 
6972be375d4SJed Brown #undef __FUNCT__
69894013140SBarry Smith #define __FUNCT__ "DMRefine_DA"
6997087cfbeSBarry Smith PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
70047c6ae99SBarry Smith {
70147c6ae99SBarry Smith   PetscErrorCode ierr;
702f3141302SJed Brown   PetscInt       M,N,P,i;
7039a42bb27SBarry Smith   DM             da2;
70447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
70547c6ae99SBarry Smith 
70647c6ae99SBarry Smith   PetscFunctionBegin;
70747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
70847c6ae99SBarry Smith   PetscValidPointer(daref,3);
70947c6ae99SBarry Smith 
7101321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
71147c6ae99SBarry Smith     M = dd->refine_x*dd->M;
71247c6ae99SBarry Smith   } else {
71347c6ae99SBarry Smith     M = 1 + dd->refine_x*(dd->M - 1);
71447c6ae99SBarry Smith   }
7151321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
7161860e6e9SBarry Smith     if (dd->dim > 1) {
71747c6ae99SBarry Smith       N = dd->refine_y*dd->N;
71847c6ae99SBarry Smith     } else {
7191860e6e9SBarry Smith       N = 1;
7201860e6e9SBarry Smith     }
7211860e6e9SBarry Smith   } else {
72247c6ae99SBarry Smith     N = 1 + dd->refine_y*(dd->N - 1);
72347c6ae99SBarry Smith   }
7241321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
7251860e6e9SBarry Smith     if (dd->dim > 2) {
72647c6ae99SBarry Smith       P = dd->refine_z*dd->P;
72747c6ae99SBarry Smith     } else {
7281860e6e9SBarry Smith       P = 1;
7291860e6e9SBarry Smith     }
7301860e6e9SBarry Smith   } else {
73147c6ae99SBarry Smith     P = 1 + dd->refine_z*(dd->P - 1);
73247c6ae99SBarry Smith   }
733397b6216SJed Brown   ierr = DMDACreate(((PetscObject)da)->comm,&da2);CHKERRQ(ierr);
734397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
735397b6216SJed Brown   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
736397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
737397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
738397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
739397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
740397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
741397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
74247c6ae99SBarry Smith   if (dd->dim == 3) {
74347c6ae99SBarry Smith     PetscInt *lx,*ly,*lz;
74447c6ae99SBarry Smith     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
7451321219cSEthan Coon     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
7461321219cSEthan Coon     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
7471321219cSEthan Coon     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
748397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
74947c6ae99SBarry Smith     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
75047c6ae99SBarry Smith   } else if (dd->dim == 2) {
75147c6ae99SBarry Smith     PetscInt *lx,*ly;
75247c6ae99SBarry Smith     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
7531321219cSEthan Coon     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
7541321219cSEthan Coon     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
755397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,PETSC_NULL);CHKERRQ(ierr);
75647c6ae99SBarry Smith     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
75747c6ae99SBarry Smith   } else if (dd->dim == 1) {
75847c6ae99SBarry Smith     PetscInt *lx;
75947c6ae99SBarry Smith     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
7601321219cSEthan Coon     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
761397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
76247c6ae99SBarry Smith     ierr = PetscFree(lx);CHKERRQ(ierr);
76347c6ae99SBarry Smith   }
76447c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
76547c6ae99SBarry Smith 
76647c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
76725296bd5SBarry Smith   da2->ops->creatematrix        = da->ops->creatematrix;
76825296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
76947c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
77047c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
77147c6ae99SBarry Smith 
77247c6ae99SBarry Smith   /* copy fill information if given */
77347c6ae99SBarry Smith   if (dd->dfill) {
77447c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
77547c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
77647c6ae99SBarry Smith   }
77747c6ae99SBarry Smith   if (dd->ofill) {
77847c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
77947c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
78047c6ae99SBarry Smith   }
78147c6ae99SBarry Smith   /* copy the refine information */
782397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
783397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
784397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
78547c6ae99SBarry Smith 
78647c6ae99SBarry Smith   /* copy vector type information */
78747c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
78819fd82e9SBarry Smith   ierr = PetscStrallocpy(da->vectype,(char**)&da2->vectype);CHKERRQ(ierr);
789ddcf8b74SDave May 
790efd51863SBarry Smith   dd2->lf = dd->lf;
791efd51863SBarry Smith   dd2->lj = dd->lj;
792efd51863SBarry Smith 
7936e87535bSJed Brown   da2->leveldown = da->leveldown;
7946e87535bSJed Brown   da2->levelup   = da->levelup + 1;
7956e87535bSJed Brown   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
7966e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
797*ca266f36SBarry Smith   ierr = DMViewFromOptions(da2,"-dm_view");CHKERRQ(ierr);
7986e87535bSJed Brown 
799ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
8006636e97aSMatthew G Knepley   if (da->coordinates) {
801ddcf8b74SDave May     DM  cdaf,cdac;
802ddcf8b74SDave May     Vec coordsc,coordsf;
803ddcf8b74SDave May     Mat II;
804ddcf8b74SDave May 
8056636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr);
8066636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr);
8076636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr);
808b61d3410SDave May     /* force creation of the coordinate vector */
809b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
8106636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
811e727c939SJed Brown     ierr = DMCreateInterpolation(cdac,cdaf,&II,PETSC_NULL);CHKERRQ(ierr);
812ddcf8b74SDave May     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
813fcfd50ebSBarry Smith     ierr = MatDestroy(&II);CHKERRQ(ierr);
814ddcf8b74SDave May   }
815397b6216SJed Brown 
816f3141302SJed Brown   for (i=0; i<da->bs; i++) {
817f3141302SJed Brown     const char *fieldname;
818f3141302SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
819f3141302SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
820f3141302SJed Brown   }
821397b6216SJed Brown 
82247c6ae99SBarry Smith   *daref = da2;
82347c6ae99SBarry Smith   PetscFunctionReturn(0);
82447c6ae99SBarry Smith }
82547c6ae99SBarry Smith 
826397b6216SJed Brown 
82747c6ae99SBarry Smith #undef __FUNCT__
82894013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA"
8297087cfbeSBarry Smith PetscErrorCode  DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
83047c6ae99SBarry Smith {
83147c6ae99SBarry Smith   PetscErrorCode ierr;
832397b6216SJed Brown   PetscInt       M,N,P,i;
8339a42bb27SBarry Smith   DM             da2;
83447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
83547c6ae99SBarry Smith 
83647c6ae99SBarry Smith   PetscFunctionBegin;
83747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
83847c6ae99SBarry Smith   PetscValidPointer(daref,3);
83947c6ae99SBarry Smith 
8401321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
841397b6216SJed Brown     M = dd->M / dd->coarsen_x;
84247c6ae99SBarry Smith   } else {
843397b6216SJed Brown     M = 1 + (dd->M - 1) / dd->coarsen_x;
84447c6ae99SBarry Smith   }
8451321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
8461860e6e9SBarry Smith     if (dd->dim > 1) {
847397b6216SJed Brown       N = dd->N / dd->coarsen_y;
84847c6ae99SBarry Smith     } else {
8491860e6e9SBarry Smith       N = 1;
8501860e6e9SBarry Smith     }
8511860e6e9SBarry Smith   } else {
852397b6216SJed Brown     N = 1 + (dd->N - 1) / dd->coarsen_y;
85347c6ae99SBarry Smith   }
8541321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
8551860e6e9SBarry Smith     if (dd->dim > 2) {
856397b6216SJed Brown       P = dd->P / dd->coarsen_z;
85747c6ae99SBarry Smith     } else {
8581860e6e9SBarry Smith       P = 1;
8591860e6e9SBarry Smith     }
8601860e6e9SBarry Smith   } else {
861397b6216SJed Brown     P = 1 + (dd->P - 1) / dd->coarsen_z;
86247c6ae99SBarry Smith   }
863397b6216SJed Brown   ierr = DMDACreate(((PetscObject)da)->comm,&da2);CHKERRQ(ierr);
864397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
865397b6216SJed Brown   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
866397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
867397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
868397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
869397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
870397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
871397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
87247c6ae99SBarry Smith   if (dd->dim == 3) {
8732be375d4SJed Brown     PetscInt *lx,*ly,*lz;
8742be375d4SJed Brown     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
875397b6216SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
876397b6216SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
877397b6216SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
878397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
8792be375d4SJed Brown     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
88047c6ae99SBarry Smith   } else if (dd->dim == 2) {
8812be375d4SJed Brown     PetscInt *lx,*ly;
8822be375d4SJed Brown     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
883397b6216SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
884397b6216SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
885397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,PETSC_NULL);CHKERRQ(ierr);
8862be375d4SJed Brown     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
88747c6ae99SBarry Smith   } else if (dd->dim == 1) {
8882be375d4SJed Brown     PetscInt *lx;
8892be375d4SJed Brown     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
890397b6216SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
891397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
8922be375d4SJed Brown     ierr = PetscFree(lx);CHKERRQ(ierr);
89347c6ae99SBarry Smith   }
89447c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
89547c6ae99SBarry Smith 
8964dcab191SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
89725296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */
89825296bd5SBarry Smith   da2->ops->creatematrix        = da->ops->creatematrix;
89947c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
90047c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
90147c6ae99SBarry Smith 
90247c6ae99SBarry Smith   /* copy fill information if given */
90347c6ae99SBarry Smith   if (dd->dfill) {
90447c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
90547c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
90647c6ae99SBarry Smith   }
90747c6ae99SBarry Smith   if (dd->ofill) {
90847c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
90947c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
91047c6ae99SBarry Smith   }
91147c6ae99SBarry Smith   /* copy the refine information */
912397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
913397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
914397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
91547c6ae99SBarry Smith 
91647c6ae99SBarry Smith   /* copy vector type information */
91747c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
91819fd82e9SBarry Smith   ierr = PetscStrallocpy(da->vectype,(char**)&da2->vectype);CHKERRQ(ierr);
91947c6ae99SBarry Smith 
920644e2e5bSBarry Smith   dd2->lf = dd->lf;
921644e2e5bSBarry Smith   dd2->lj = dd->lj;
922644e2e5bSBarry Smith 
9236e87535bSJed Brown   da2->leveldown = da->leveldown + 1;
9246e87535bSJed Brown   da2->levelup   = da->levelup;
9256e87535bSJed Brown   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
9266e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
927*ca266f36SBarry Smith   ierr = DMViewFromOptions(da2,"-dm_view");CHKERRQ(ierr);
9286e87535bSJed Brown 
929ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
9306636e97aSMatthew G Knepley   if (da->coordinates) {
931ddcf8b74SDave May     DM         cdaf,cdac;
932ddcf8b74SDave May     Vec        coordsc,coordsf;
933ddcf8b74SDave May     VecScatter inject;
934ddcf8b74SDave May 
9356636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr);
9366636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr);
9376636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr);
938b61d3410SDave May     /* force creation of the coordinate vector */
939b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
9406636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
941ddcf8b74SDave May 
942e727c939SJed Brown     ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
943ddcf8b74SDave May     ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
944ddcf8b74SDave May     ierr = VecScatterEnd(inject  ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
945fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&inject);CHKERRQ(ierr);
946ddcf8b74SDave May   }
947f98405f7SJed Brown 
948f98405f7SJed Brown   for (i=0; i<da->bs; i++) {
949f98405f7SJed Brown     const char *fieldname;
950f98405f7SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
951f98405f7SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
952f98405f7SJed Brown   }
953f98405f7SJed Brown 
95447c6ae99SBarry Smith   *daref = da2;
95547c6ae99SBarry Smith   PetscFunctionReturn(0);
95647c6ae99SBarry Smith }
95747c6ae99SBarry Smith 
95847c6ae99SBarry Smith #undef __FUNCT__
95994013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA"
9607087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
96147c6ae99SBarry Smith {
96247c6ae99SBarry Smith   PetscErrorCode ierr;
96347c6ae99SBarry Smith   PetscInt       i,n,*refx,*refy,*refz;
96447c6ae99SBarry Smith 
96547c6ae99SBarry Smith   PetscFunctionBegin;
96647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
96747c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
96847c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
96947c6ae99SBarry Smith   PetscValidPointer(daf,3);
97047c6ae99SBarry Smith 
971aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
97247c6ae99SBarry Smith   ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr);
97347c6ae99SBarry Smith   for (i=0; i<nlevels; i++) {
974aa219208SBarry Smith     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
97547c6ae99SBarry Smith   }
97647c6ae99SBarry Smith   n = nlevels;
97747c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
97847c6ae99SBarry Smith   n = nlevels;
97947c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
98047c6ae99SBarry Smith   n = nlevels;
98147c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
98247c6ae99SBarry Smith 
983aa219208SBarry Smith   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
98494013140SBarry Smith   ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr);
98547c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
986aa219208SBarry Smith     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
98794013140SBarry Smith     ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr);
98847c6ae99SBarry Smith   }
98947c6ae99SBarry Smith   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
99047c6ae99SBarry Smith   PetscFunctionReturn(0);
99147c6ae99SBarry Smith }
99247c6ae99SBarry Smith 
99347c6ae99SBarry Smith #undef __FUNCT__
99494013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA"
9957087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
99647c6ae99SBarry Smith {
99747c6ae99SBarry Smith   PetscErrorCode ierr;
99847c6ae99SBarry Smith   PetscInt       i;
99947c6ae99SBarry Smith 
100047c6ae99SBarry Smith   PetscFunctionBegin;
100147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
100247c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
100347c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
100447c6ae99SBarry Smith   PetscValidPointer(dac,3);
100594013140SBarry Smith   ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr);
100647c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
100794013140SBarry Smith     ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr);
100847c6ae99SBarry Smith   }
100947c6ae99SBarry Smith   PetscFunctionReturn(0);
101047c6ae99SBarry Smith }
1011