xref: /petsc/src/dm/impls/da/da.c (revision 6636e97a9dbd9dfa567780c5ec4fb4018c0e427e)
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 
18688661749SPeter Brune 
18788661749SPeter Brune #undef __FUNCT__
188aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType"
18947c6ae99SBarry Smith /*@
190aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
19147c6ae99SBarry Smith 
192aa219208SBarry Smith   Logically Collective on DMDA
19347c6ae99SBarry Smith 
19447c6ae99SBarry Smith   Input Parameter:
195aa219208SBarry Smith + da    - The DMDA
196aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
19747c6ae99SBarry Smith 
19847c6ae99SBarry Smith   Level: intermediate
19947c6ae99SBarry Smith 
20047c6ae99SBarry Smith .keywords:  distributed array, stencil
20196e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
20247c6ae99SBarry Smith @*/
2037087cfbeSBarry Smith PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
20447c6ae99SBarry Smith {
20547c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
20647c6ae99SBarry Smith 
20747c6ae99SBarry Smith   PetscFunctionBegin;
20847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
20947c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,stype,2);
210cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
21147c6ae99SBarry Smith   dd->stencil_type = stype;
21247c6ae99SBarry Smith   PetscFunctionReturn(0);
21347c6ae99SBarry Smith }
21447c6ae99SBarry Smith 
21547c6ae99SBarry Smith #undef __FUNCT__
216aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth"
21747c6ae99SBarry Smith /*@
218aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
21947c6ae99SBarry Smith 
220aa219208SBarry Smith   Logically Collective on DMDA
22147c6ae99SBarry Smith 
22247c6ae99SBarry Smith   Input Parameter:
223aa219208SBarry Smith + da    - The DMDA
22447c6ae99SBarry Smith - width - The stencil width
22547c6ae99SBarry Smith 
22647c6ae99SBarry Smith   Level: intermediate
22747c6ae99SBarry Smith 
22847c6ae99SBarry Smith .keywords:  distributed array, stencil
22996e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
23047c6ae99SBarry Smith @*/
2317087cfbeSBarry Smith PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
23247c6ae99SBarry Smith {
23347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
23447c6ae99SBarry Smith 
23547c6ae99SBarry Smith   PetscFunctionBegin;
23647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
23747c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,width,2);
238cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
23947c6ae99SBarry Smith   dd->s = width;
24047c6ae99SBarry Smith   PetscFunctionReturn(0);
24147c6ae99SBarry Smith }
24247c6ae99SBarry Smith 
24347c6ae99SBarry Smith #undef __FUNCT__
244aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
245aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
24647c6ae99SBarry Smith {
24747c6ae99SBarry Smith   PetscInt i,sum;
24847c6ae99SBarry Smith 
24947c6ae99SBarry Smith   PetscFunctionBegin;
25047c6ae99SBarry Smith   if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
25147c6ae99SBarry Smith   for (i=sum=0; i<m; i++) sum += lx[i];
25247c6ae99SBarry Smith   if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
25347c6ae99SBarry Smith   PetscFunctionReturn(0);
25447c6ae99SBarry Smith }
25547c6ae99SBarry Smith 
25647c6ae99SBarry Smith #undef __FUNCT__
257aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges"
25847c6ae99SBarry Smith /*@
259aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
26047c6ae99SBarry Smith 
261aa219208SBarry Smith   Logically Collective on DMDA
26247c6ae99SBarry Smith 
26347c6ae99SBarry Smith   Input Parameter:
264aa219208SBarry Smith + da - The DMDA
26547c6ae99SBarry 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
26647c6ae99SBarry 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
26747c6ae99SBarry 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.
26847c6ae99SBarry Smith 
26947c6ae99SBarry Smith   Level: intermediate
27047c6ae99SBarry Smith 
27147c6ae99SBarry Smith .keywords:  distributed array
27296e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
27347c6ae99SBarry Smith @*/
2747087cfbeSBarry Smith PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
27547c6ae99SBarry Smith {
27647c6ae99SBarry Smith   PetscErrorCode ierr;
27747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
27847c6ae99SBarry Smith 
27947c6ae99SBarry Smith   PetscFunctionBegin;
28047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
281cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
28247c6ae99SBarry Smith   if (lx) {
28347c6ae99SBarry Smith     if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
284aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
28547c6ae99SBarry Smith     if (!dd->lx) {
28647c6ae99SBarry Smith       ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
28747c6ae99SBarry Smith     }
28847c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
28947c6ae99SBarry Smith   }
29047c6ae99SBarry Smith   if (ly) {
29147c6ae99SBarry Smith     if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
292aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
29347c6ae99SBarry Smith     if (!dd->ly) {
29447c6ae99SBarry Smith       ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
29547c6ae99SBarry Smith     }
29647c6ae99SBarry Smith     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
29747c6ae99SBarry Smith   }
29847c6ae99SBarry Smith   if (lz) {
29947c6ae99SBarry Smith     if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
300aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
30147c6ae99SBarry Smith     if (!dd->lz) {
30247c6ae99SBarry Smith       ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
30347c6ae99SBarry Smith     }
30447c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
30547c6ae99SBarry Smith   }
30647c6ae99SBarry Smith   PetscFunctionReturn(0);
30747c6ae99SBarry Smith }
30847c6ae99SBarry Smith 
30947c6ae99SBarry Smith #undef __FUNCT__
310aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType"
31147c6ae99SBarry Smith /*@
312aa219208SBarry Smith        DMDASetInterpolationType - Sets the type of interpolation that will be
313e727c939SJed Brown           returned by DMCreateInterpolation()
31447c6ae99SBarry Smith 
315aa219208SBarry Smith    Logically Collective on DMDA
31647c6ae99SBarry Smith 
31747c6ae99SBarry Smith    Input Parameter:
31847c6ae99SBarry Smith +  da - initial distributed array
319aa219208SBarry Smith .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
32047c6ae99SBarry Smith 
32147c6ae99SBarry Smith    Level: intermediate
32247c6ae99SBarry Smith 
323e727c939SJed Brown    Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation()
32447c6ae99SBarry Smith 
32547c6ae99SBarry Smith .keywords:  distributed array, interpolation
32647c6ae99SBarry Smith 
32796e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
32847c6ae99SBarry Smith @*/
3297087cfbeSBarry Smith PetscErrorCode  DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
33047c6ae99SBarry Smith {
33147c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
33247c6ae99SBarry Smith 
33347c6ae99SBarry Smith   PetscFunctionBegin;
33447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
33547c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,ctype,2);
33647c6ae99SBarry Smith   dd->interptype = ctype;
33747c6ae99SBarry Smith   PetscFunctionReturn(0);
33847c6ae99SBarry Smith }
33947c6ae99SBarry Smith 
3402dde6fd4SLisandro Dalcin #undef __FUNCT__
3412dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetInterpolationType"
3422dde6fd4SLisandro Dalcin /*@
3432dde6fd4SLisandro Dalcin        DMDAGetInterpolationType - Gets the type of interpolation that will be
344e727c939SJed Brown           used by DMCreateInterpolation()
3452dde6fd4SLisandro Dalcin 
3462dde6fd4SLisandro Dalcin    Not Collective
3472dde6fd4SLisandro Dalcin 
3482dde6fd4SLisandro Dalcin    Input Parameter:
3492dde6fd4SLisandro Dalcin .  da - distributed array
3502dde6fd4SLisandro Dalcin 
3512dde6fd4SLisandro Dalcin    Output Parameter:
3522dde6fd4SLisandro Dalcin .  ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms)
3532dde6fd4SLisandro Dalcin 
3542dde6fd4SLisandro Dalcin    Level: intermediate
3552dde6fd4SLisandro Dalcin 
3562dde6fd4SLisandro Dalcin .keywords:  distributed array, interpolation
3572dde6fd4SLisandro Dalcin 
358e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation()
3592dde6fd4SLisandro Dalcin @*/
3602dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype)
3612dde6fd4SLisandro Dalcin {
3622dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA*)da->data;
3632dde6fd4SLisandro Dalcin 
3642dde6fd4SLisandro Dalcin   PetscFunctionBegin;
3652dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
3662dde6fd4SLisandro Dalcin   PetscValidPointer(ctype,2);
3672dde6fd4SLisandro Dalcin   *ctype = dd->interptype;
3682dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
3692dde6fd4SLisandro Dalcin }
37047c6ae99SBarry Smith 
37147c6ae99SBarry Smith #undef __FUNCT__
372aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors"
37347c6ae99SBarry Smith /*@C
374aa219208SBarry Smith       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
37547c6ae99SBarry Smith         processes neighbors.
37647c6ae99SBarry Smith 
37747c6ae99SBarry Smith     Not Collective
37847c6ae99SBarry Smith 
37947c6ae99SBarry Smith    Input Parameter:
380aa219208SBarry Smith .     da - the DMDA object
38147c6ae99SBarry Smith 
38247c6ae99SBarry Smith    Output Parameters:
38347c6ae99SBarry Smith .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
38447c6ae99SBarry Smith               this process itself is in the list
38547c6ae99SBarry Smith 
38647c6ae99SBarry Smith    Notes: In 2d the array is of length 9, in 3d of length 27
38747c6ae99SBarry Smith           Not supported in 1d
388aa219208SBarry Smith           Do not free the array, it is freed when the DMDA is destroyed.
38947c6ae99SBarry Smith 
39047c6ae99SBarry Smith    Fortran Notes: In fortran you must pass in an array of the appropriate length.
39147c6ae99SBarry Smith 
39247c6ae99SBarry Smith    Level: intermediate
39347c6ae99SBarry Smith 
39447c6ae99SBarry Smith @*/
3957087cfbeSBarry Smith PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
39647c6ae99SBarry Smith {
39747c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
39847c6ae99SBarry Smith   PetscFunctionBegin;
39947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
40047c6ae99SBarry Smith   *ranks = dd->neighbors;
40147c6ae99SBarry Smith   PetscFunctionReturn(0);
40247c6ae99SBarry Smith }
40347c6ae99SBarry Smith 
40447c6ae99SBarry Smith #undef __FUNCT__
405aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges"
40647c6ae99SBarry Smith /*@C
407aa219208SBarry Smith       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
40847c6ae99SBarry Smith 
40947c6ae99SBarry Smith     Not Collective
41047c6ae99SBarry Smith 
41147c6ae99SBarry Smith    Input Parameter:
412aa219208SBarry Smith .     da - the DMDA object
41347c6ae99SBarry Smith 
41447c6ae99SBarry Smith    Output Parameter:
41547c6ae99SBarry Smith +     lx - ownership along x direction (optional)
41647c6ae99SBarry Smith .     ly - ownership along y direction (optional)
41747c6ae99SBarry Smith -     lz - ownership along z direction (optional)
41847c6ae99SBarry Smith 
41947c6ae99SBarry Smith    Level: intermediate
42047c6ae99SBarry Smith 
421aa219208SBarry Smith     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
42247c6ae99SBarry Smith 
42347c6ae99SBarry Smith     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
424aa219208SBarry Smith     eighth arguments from DMDAGetInfo()
42547c6ae99SBarry Smith 
42647c6ae99SBarry Smith      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
427aa219208SBarry Smith     DMDA they came from still exists (has not been destroyed).
42847c6ae99SBarry Smith 
429aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
43047c6ae99SBarry Smith @*/
4317087cfbeSBarry Smith PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
43247c6ae99SBarry Smith {
43347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
43447c6ae99SBarry Smith 
43547c6ae99SBarry Smith   PetscFunctionBegin;
43647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
43747c6ae99SBarry Smith   if (lx) *lx = dd->lx;
43847c6ae99SBarry Smith   if (ly) *ly = dd->ly;
43947c6ae99SBarry Smith   if (lz) *lz = dd->lz;
44047c6ae99SBarry Smith   PetscFunctionReturn(0);
44147c6ae99SBarry Smith }
44247c6ae99SBarry Smith 
44347c6ae99SBarry Smith #undef __FUNCT__
444aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor"
44547c6ae99SBarry Smith /*@
446aa219208SBarry Smith      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
44747c6ae99SBarry Smith 
448aa219208SBarry Smith     Logically Collective on DMDA
44947c6ae99SBarry Smith 
45047c6ae99SBarry Smith   Input Parameters:
451aa219208SBarry Smith +    da - the DMDA object
45247c6ae99SBarry Smith .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
45347c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
45447c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
45547c6ae99SBarry Smith 
45647c6ae99SBarry Smith   Options Database:
45747c6ae99SBarry Smith +  -da_refine_x - refinement ratio in x direction
45847c6ae99SBarry Smith .  -da_refine_y - refinement ratio in y direction
45947c6ae99SBarry Smith -  -da_refine_z - refinement ratio in z direction
46047c6ae99SBarry Smith 
46147c6ae99SBarry Smith   Level: intermediate
46247c6ae99SBarry Smith 
46347c6ae99SBarry Smith     Notes: Pass PETSC_IGNORE to leave a value unchanged
46447c6ae99SBarry Smith 
465aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor()
46647c6ae99SBarry Smith @*/
4677087cfbeSBarry Smith PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
46847c6ae99SBarry Smith {
46947c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
47047c6ae99SBarry Smith 
47147c6ae99SBarry Smith   PetscFunctionBegin;
47247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
47347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_x,2);
47447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_y,3);
47547c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_z,4);
47647c6ae99SBarry Smith 
47747c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
47847c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
47947c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
48047c6ae99SBarry Smith   PetscFunctionReturn(0);
48147c6ae99SBarry Smith }
48247c6ae99SBarry Smith 
48347c6ae99SBarry Smith #undef __FUNCT__
484aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor"
48547c6ae99SBarry Smith /*@C
486aa219208SBarry Smith      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
48747c6ae99SBarry Smith 
48847c6ae99SBarry Smith     Not Collective
48947c6ae99SBarry Smith 
49047c6ae99SBarry Smith   Input Parameter:
491aa219208SBarry Smith .    da - the DMDA object
49247c6ae99SBarry Smith 
49347c6ae99SBarry Smith   Output Parameters:
49447c6ae99SBarry Smith +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
49547c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
49647c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
49747c6ae99SBarry Smith 
49847c6ae99SBarry Smith   Level: intermediate
49947c6ae99SBarry Smith 
50047c6ae99SBarry Smith     Notes: Pass PETSC_NULL for values you do not need
50147c6ae99SBarry Smith 
502aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor()
50347c6ae99SBarry Smith @*/
5047087cfbeSBarry Smith PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
50547c6ae99SBarry Smith {
50647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
50747c6ae99SBarry Smith 
50847c6ae99SBarry Smith   PetscFunctionBegin;
50947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
51047c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
51147c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
51247c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
51347c6ae99SBarry Smith   PetscFunctionReturn(0);
51447c6ae99SBarry Smith }
51547c6ae99SBarry Smith 
51647c6ae99SBarry Smith #undef __FUNCT__
517aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix"
51847c6ae99SBarry Smith /*@C
519aa219208SBarry Smith      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
52047c6ae99SBarry Smith 
521aa219208SBarry Smith     Logically Collective on DMDA
52247c6ae99SBarry Smith 
52347c6ae99SBarry Smith   Input Parameters:
524aa219208SBarry Smith +    da - the DMDA object
525aa219208SBarry Smith -    f - the function that allocates the matrix for that specific DMDA
52647c6ae99SBarry Smith 
52747c6ae99SBarry Smith   Level: developer
52847c6ae99SBarry Smith 
529aa219208SBarry Smith    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
53047c6ae99SBarry Smith        the diagonal and off-diagonal blocks of the matrix
53147c6ae99SBarry Smith 
532950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills()
53347c6ae99SBarry Smith @*/
5347087cfbeSBarry Smith PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, const MatType,Mat*))
53547c6ae99SBarry Smith {
53647c6ae99SBarry Smith   PetscFunctionBegin;
53747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
53825296bd5SBarry Smith   da->ops->creatematrix = f;
53947c6ae99SBarry Smith   PetscFunctionReturn(0);
54047c6ae99SBarry Smith }
54147c6ae99SBarry Smith 
54247c6ae99SBarry Smith #undef __FUNCT__
54394013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges"
54447c6ae99SBarry Smith /*
54547c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
54647c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
54747c6ae99SBarry Smith 
54847c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
54947c6ae99SBarry Smith */
55094013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
55147c6ae99SBarry Smith {
55247c6ae99SBarry Smith   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
55347c6ae99SBarry Smith   PetscErrorCode ierr;
55447c6ae99SBarry Smith 
55547c6ae99SBarry Smith   PetscFunctionBegin;
55647c6ae99SBarry Smith   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
55747c6ae99SBarry Smith   if (ratio == 1) {
55847c6ae99SBarry Smith     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
55947c6ae99SBarry Smith     PetscFunctionReturn(0);
56047c6ae99SBarry Smith   }
56147c6ae99SBarry Smith   for (i=0; i<m; i++) totalc += lc[i];
56247c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
56347c6ae99SBarry Smith   for (i=0; i<m; i++) {
56447c6ae99SBarry Smith     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
56547c6ae99SBarry Smith     if (i == m-1) lf[i] = want;
56647c6ae99SBarry Smith     else {
5677aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
5687aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
5697aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
5707aca7175SJed Brown       while ((startf+want)/ratio < nextc - stencil_width) want++;
5717aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
5727aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
5737aca7175SJed Brown       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
5747aca7175SJed Brown       /* Make sure all constraints are satisfied */
57530729d88SBarry Smith       if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width)
57630729d88SBarry 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");
57747c6ae99SBarry Smith     }
57847c6ae99SBarry Smith     lf[i] = want;
57947c6ae99SBarry Smith     startc += lc[i];
58047c6ae99SBarry Smith     startf += lf[i];
58147c6ae99SBarry Smith     remaining -= lf[i];
58247c6ae99SBarry Smith   }
58347c6ae99SBarry Smith   PetscFunctionReturn(0);
58447c6ae99SBarry Smith }
58547c6ae99SBarry Smith 
58647c6ae99SBarry Smith #undef __FUNCT__
5872be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges"
5882be375d4SJed Brown /*
5892be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
5902be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
5912be375d4SJed Brown 
5922be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
5932be375d4SJed Brown */
5942be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
5952be375d4SJed Brown {
5962be375d4SJed Brown   PetscInt       i,totalf,remaining,startc,startf;
5972be375d4SJed Brown   PetscErrorCode ierr;
5982be375d4SJed Brown 
5992be375d4SJed Brown   PetscFunctionBegin;
6002be375d4SJed Brown   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
6012be375d4SJed Brown   if (ratio == 1) {
6022be375d4SJed Brown     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
6032be375d4SJed Brown     PetscFunctionReturn(0);
6042be375d4SJed Brown   }
6052be375d4SJed Brown   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
6062be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
6072be375d4SJed Brown   for (i=0,startc=0,startf=0; i<m; i++) {
6082be375d4SJed Brown     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
6092be375d4SJed Brown     if (i == m-1) lc[i] = want;
6102be375d4SJed Brown     else {
6112be375d4SJed Brown       const PetscInt nextf = startf+lf[i];
6122be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
6132be375d4SJed Brown        * node is within one stencil width. */
6142be375d4SJed Brown       while (nextf/ratio < startc+want-stencil_width) want--;
6152be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
6162be375d4SJed Brown        * fine node is within one stencil width. */
6172be375d4SJed Brown       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
6182be375d4SJed Brown       if (want < 0 || want > remaining
61930729d88SBarry 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");
6202be375d4SJed Brown     }
6212be375d4SJed Brown     lc[i] = want;
6222be375d4SJed Brown     startc += lc[i];
6232be375d4SJed Brown     startf += lf[i];
6242be375d4SJed Brown     remaining -= lc[i];
6252be375d4SJed Brown   }
6262be375d4SJed Brown   PetscFunctionReturn(0);
6272be375d4SJed Brown }
6282be375d4SJed Brown 
6292be375d4SJed Brown #undef __FUNCT__
63094013140SBarry Smith #define __FUNCT__ "DMRefine_DA"
6317087cfbeSBarry Smith PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
63247c6ae99SBarry Smith {
63347c6ae99SBarry Smith   PetscErrorCode ierr;
634f3141302SJed Brown   PetscInt       M,N,P,i;
6359a42bb27SBarry Smith   DM             da2;
63647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
63747c6ae99SBarry Smith 
63847c6ae99SBarry Smith   PetscFunctionBegin;
63947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
64047c6ae99SBarry Smith   PetscValidPointer(daref,3);
64147c6ae99SBarry Smith 
6421321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
64347c6ae99SBarry Smith     M = dd->refine_x*dd->M;
64447c6ae99SBarry Smith   } else {
64547c6ae99SBarry Smith     M = 1 + dd->refine_x*(dd->M - 1);
64647c6ae99SBarry Smith   }
6471321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
6481860e6e9SBarry Smith     if (dd->dim > 1) {
64947c6ae99SBarry Smith       N = dd->refine_y*dd->N;
65047c6ae99SBarry Smith     } else {
6511860e6e9SBarry Smith       N = 1;
6521860e6e9SBarry Smith     }
6531860e6e9SBarry Smith   } else {
65447c6ae99SBarry Smith     N = 1 + dd->refine_y*(dd->N - 1);
65547c6ae99SBarry Smith   }
6561321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
6571860e6e9SBarry Smith     if (dd->dim > 2) {
65847c6ae99SBarry Smith       P = dd->refine_z*dd->P;
65947c6ae99SBarry Smith     } else {
6601860e6e9SBarry Smith       P = 1;
6611860e6e9SBarry Smith     }
6621860e6e9SBarry Smith   } else {
66347c6ae99SBarry Smith     P = 1 + dd->refine_z*(dd->P - 1);
66447c6ae99SBarry Smith   }
665397b6216SJed Brown   ierr = DMDACreate(((PetscObject)da)->comm,&da2);CHKERRQ(ierr);
666397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
667397b6216SJed Brown   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
668397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
669397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
670397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
671397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
672397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
673397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
67447c6ae99SBarry Smith   if (dd->dim == 3) {
67547c6ae99SBarry Smith     PetscInt *lx,*ly,*lz;
67647c6ae99SBarry Smith     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
6771321219cSEthan 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);
6781321219cSEthan 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);
6791321219cSEthan 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);
680397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
68147c6ae99SBarry Smith     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
68247c6ae99SBarry Smith   } else if (dd->dim == 2) {
68347c6ae99SBarry Smith     PetscInt *lx,*ly;
68447c6ae99SBarry Smith     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
6851321219cSEthan 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);
6861321219cSEthan 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);
687397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,PETSC_NULL);CHKERRQ(ierr);
68847c6ae99SBarry Smith     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
68947c6ae99SBarry Smith   } else if (dd->dim == 1) {
69047c6ae99SBarry Smith     PetscInt *lx;
69147c6ae99SBarry Smith     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
6921321219cSEthan 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);
693397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
69447c6ae99SBarry Smith     ierr = PetscFree(lx);CHKERRQ(ierr);
69547c6ae99SBarry Smith   }
69647c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
69747c6ae99SBarry Smith 
69847c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
69925296bd5SBarry Smith   da2->ops->creatematrix        = da->ops->creatematrix;
70025296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
70147c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
70247c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
70347c6ae99SBarry Smith 
70447c6ae99SBarry Smith   /* copy fill information if given */
70547c6ae99SBarry Smith   if (dd->dfill) {
70647c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
70747c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
70847c6ae99SBarry Smith   }
70947c6ae99SBarry Smith   if (dd->ofill) {
71047c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
71147c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
71247c6ae99SBarry Smith   }
71347c6ae99SBarry Smith   /* copy the refine information */
714397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
715397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
716397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
71747c6ae99SBarry Smith 
71847c6ae99SBarry Smith   /* copy vector type information */
71947c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
72047c6ae99SBarry Smith   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
721ddcf8b74SDave May 
722efd51863SBarry Smith   dd2->lf = dd->lf;
723efd51863SBarry Smith   dd2->lj = dd->lj;
724efd51863SBarry Smith 
7256e87535bSJed Brown   da2->leveldown = da->leveldown;
7266e87535bSJed Brown   da2->levelup   = da->levelup + 1;
7276e87535bSJed Brown   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
7286e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
7296e87535bSJed Brown   ierr = DMView_DA_Private(da2);CHKERRQ(ierr);
7306e87535bSJed Brown 
731ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
732*6636e97aSMatthew G Knepley   if (da->coordinates) {
733ddcf8b74SDave May     DM  cdaf,cdac;
734ddcf8b74SDave May     Vec coordsc,coordsf;
735ddcf8b74SDave May     Mat II;
736ddcf8b74SDave May 
737*6636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr);
738*6636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr);
739*6636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr);
740b61d3410SDave May     /* force creation of the coordinate vector */
741b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
742*6636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
743e727c939SJed Brown     ierr = DMCreateInterpolation(cdac,cdaf,&II,PETSC_NULL);CHKERRQ(ierr);
744ddcf8b74SDave May     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
745fcfd50ebSBarry Smith     ierr = MatDestroy(&II);CHKERRQ(ierr);
746ddcf8b74SDave May   }
747397b6216SJed Brown 
748f3141302SJed Brown   for (i=0; i<da->bs; i++) {
749f3141302SJed Brown     const char *fieldname;
750f3141302SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
751f3141302SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
752f3141302SJed Brown   }
753397b6216SJed Brown 
75447c6ae99SBarry Smith   *daref = da2;
75547c6ae99SBarry Smith   PetscFunctionReturn(0);
75647c6ae99SBarry Smith }
75747c6ae99SBarry Smith 
758397b6216SJed Brown 
75947c6ae99SBarry Smith #undef __FUNCT__
76094013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA"
7617087cfbeSBarry Smith PetscErrorCode  DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
76247c6ae99SBarry Smith {
76347c6ae99SBarry Smith   PetscErrorCode ierr;
764397b6216SJed Brown   PetscInt       M,N,P,i;
7659a42bb27SBarry Smith   DM             da2;
76647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
76747c6ae99SBarry Smith 
76847c6ae99SBarry Smith   PetscFunctionBegin;
76947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
77047c6ae99SBarry Smith   PetscValidPointer(daref,3);
77147c6ae99SBarry Smith 
7721321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
773397b6216SJed Brown     M = dd->M / dd->coarsen_x;
77447c6ae99SBarry Smith   } else {
775397b6216SJed Brown     M = 1 + (dd->M - 1) / dd->coarsen_x;
77647c6ae99SBarry Smith   }
7771321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
7781860e6e9SBarry Smith     if (dd->dim > 1) {
779397b6216SJed Brown       N = dd->N / dd->coarsen_y;
78047c6ae99SBarry Smith     } else {
7811860e6e9SBarry Smith       N = 1;
7821860e6e9SBarry Smith     }
7831860e6e9SBarry Smith   } else {
784397b6216SJed Brown     N = 1 + (dd->N - 1) / dd->coarsen_y;
78547c6ae99SBarry Smith   }
7861321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
7871860e6e9SBarry Smith     if (dd->dim > 2) {
788397b6216SJed Brown       P = dd->P / dd->coarsen_z;
78947c6ae99SBarry Smith     } else {
7901860e6e9SBarry Smith       P = 1;
7911860e6e9SBarry Smith     }
7921860e6e9SBarry Smith   } else {
793397b6216SJed Brown     P = 1 + (dd->P - 1) / dd->coarsen_z;
79447c6ae99SBarry Smith   }
795397b6216SJed Brown   ierr = DMDACreate(((PetscObject)da)->comm,&da2);CHKERRQ(ierr);
796397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
797397b6216SJed Brown   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
798397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
799397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
800397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
801397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
802397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
803397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
80447c6ae99SBarry Smith   if (dd->dim == 3) {
8052be375d4SJed Brown     PetscInt *lx,*ly,*lz;
8062be375d4SJed Brown     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
807397b6216SJed 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);
808397b6216SJed 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);
809397b6216SJed 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);
810397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
8112be375d4SJed Brown     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
81247c6ae99SBarry Smith   } else if (dd->dim == 2) {
8132be375d4SJed Brown     PetscInt *lx,*ly;
8142be375d4SJed Brown     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
815397b6216SJed 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);
816397b6216SJed 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);
817397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,PETSC_NULL);CHKERRQ(ierr);
8182be375d4SJed Brown     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
81947c6ae99SBarry Smith   } else if (dd->dim == 1) {
8202be375d4SJed Brown     PetscInt *lx;
8212be375d4SJed Brown     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
822397b6216SJed 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);
823397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
8242be375d4SJed Brown     ierr = PetscFree(lx);CHKERRQ(ierr);
82547c6ae99SBarry Smith   }
82647c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
82747c6ae99SBarry Smith 
8284dcab191SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
82925296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */
83025296bd5SBarry Smith   da2->ops->creatematrix        = da->ops->creatematrix;
83147c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
83247c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
83347c6ae99SBarry Smith 
83447c6ae99SBarry Smith   /* copy fill information if given */
83547c6ae99SBarry Smith   if (dd->dfill) {
83647c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
83747c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
83847c6ae99SBarry Smith   }
83947c6ae99SBarry Smith   if (dd->ofill) {
84047c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
84147c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
84247c6ae99SBarry Smith   }
84347c6ae99SBarry Smith   /* copy the refine information */
844397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
845397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
846397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
84747c6ae99SBarry Smith 
84847c6ae99SBarry Smith   /* copy vector type information */
84947c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
85047c6ae99SBarry Smith   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
85147c6ae99SBarry Smith 
852644e2e5bSBarry Smith   dd2->lf = dd->lf;
853644e2e5bSBarry Smith   dd2->lj = dd->lj;
854644e2e5bSBarry Smith 
8556e87535bSJed Brown   da2->leveldown = da->leveldown + 1;
8566e87535bSJed Brown   da2->levelup   = da->levelup;
8576e87535bSJed Brown   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
8586e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
8596e87535bSJed Brown   ierr = DMView_DA_Private(da2);CHKERRQ(ierr);
8606e87535bSJed Brown 
861ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
862*6636e97aSMatthew G Knepley   if (da->coordinates) {
863ddcf8b74SDave May     DM         cdaf,cdac;
864ddcf8b74SDave May     Vec        coordsc,coordsf;
865ddcf8b74SDave May     VecScatter inject;
866ddcf8b74SDave May 
867*6636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr);
868*6636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr);
869*6636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr);
870b61d3410SDave May     /* force creation of the coordinate vector */
871b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
872*6636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
873ddcf8b74SDave May 
874e727c939SJed Brown     ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
875ddcf8b74SDave May     ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
876ddcf8b74SDave May     ierr = VecScatterEnd(inject  ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
877fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&inject);CHKERRQ(ierr);
878ddcf8b74SDave May   }
879f98405f7SJed Brown 
880f98405f7SJed Brown   for (i=0; i<da->bs; i++) {
881f98405f7SJed Brown     const char *fieldname;
882f98405f7SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
883f98405f7SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
884f98405f7SJed Brown   }
885f98405f7SJed Brown 
88647c6ae99SBarry Smith   *daref = da2;
88747c6ae99SBarry Smith   PetscFunctionReturn(0);
88847c6ae99SBarry Smith }
88947c6ae99SBarry Smith 
89047c6ae99SBarry Smith #undef __FUNCT__
89194013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA"
8927087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
89347c6ae99SBarry Smith {
89447c6ae99SBarry Smith   PetscErrorCode ierr;
89547c6ae99SBarry Smith   PetscInt       i,n,*refx,*refy,*refz;
89647c6ae99SBarry Smith 
89747c6ae99SBarry Smith   PetscFunctionBegin;
89847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
89947c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
90047c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
90147c6ae99SBarry Smith   PetscValidPointer(daf,3);
90247c6ae99SBarry Smith 
903aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
90447c6ae99SBarry Smith   ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr);
90547c6ae99SBarry Smith   for (i=0; i<nlevels; i++) {
906aa219208SBarry Smith     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
90747c6ae99SBarry Smith   }
90847c6ae99SBarry Smith   n = nlevels;
90947c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
91047c6ae99SBarry Smith   n = nlevels;
91147c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
91247c6ae99SBarry Smith   n = nlevels;
91347c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
91447c6ae99SBarry Smith 
915aa219208SBarry Smith   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
91694013140SBarry Smith   ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr);
91747c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
918aa219208SBarry Smith     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
91994013140SBarry Smith     ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr);
92047c6ae99SBarry Smith   }
92147c6ae99SBarry Smith   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
92247c6ae99SBarry Smith   PetscFunctionReturn(0);
92347c6ae99SBarry Smith }
92447c6ae99SBarry Smith 
92547c6ae99SBarry Smith #undef __FUNCT__
92694013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA"
9277087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
92847c6ae99SBarry Smith {
92947c6ae99SBarry Smith   PetscErrorCode ierr;
93047c6ae99SBarry Smith   PetscInt       i;
93147c6ae99SBarry Smith 
93247c6ae99SBarry Smith   PetscFunctionBegin;
93347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
93447c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
93547c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
93647c6ae99SBarry Smith   PetscValidPointer(dac,3);
93794013140SBarry Smith   ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr);
93847c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
93994013140SBarry Smith     ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr);
94047c6ae99SBarry Smith   }
94147c6ae99SBarry Smith   PetscFunctionReturn(0);
94247c6ae99SBarry Smith }
943