xref: /petsc/src/dm/impls/da/da.c (revision 1321219c1d5e9ec6e55fd54fae073108c73e7932)
147c6ae99SBarry Smith #define PETSCDM_DLL
2e1589f56SBarry Smith #include "private/daimpl.h"    /*I   "petscdm.h"   I*/
347c6ae99SBarry Smith 
447c6ae99SBarry Smith 
547c6ae99SBarry Smith #undef __FUNCT__
6aa219208SBarry Smith #define __FUNCT__ "DMDASetDim"
747c6ae99SBarry Smith /*@
8aa219208SBarry Smith   DMDASetDim - Sets the dimension
947c6ae99SBarry Smith 
10aa219208SBarry Smith   Collective on DMDA
1147c6ae99SBarry Smith 
1247c6ae99SBarry Smith   Input Parameters:
13aa219208SBarry Smith + da - the DMDA
1447c6ae99SBarry Smith - dim - the dimension (or PETSC_DECIDE)
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith   Level: intermediate
1747c6ae99SBarry Smith 
18aa219208SBarry Smith .seealso: DaGetDim(), DMDASetSizes()
1947c6ae99SBarry Smith @*/
207087cfbeSBarry Smith PetscErrorCode  DMDASetDim(DM da, PetscInt dim)
2147c6ae99SBarry Smith {
2247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   PetscFunctionBegin;
2547c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
26aa219208SBarry 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);
2747c6ae99SBarry Smith   dd->dim = dim;
2847c6ae99SBarry Smith   PetscFunctionReturn(0);
2947c6ae99SBarry Smith }
3047c6ae99SBarry Smith 
3147c6ae99SBarry Smith #undef __FUNCT__
32aa219208SBarry Smith #define __FUNCT__ "DMDASetSizes"
3347c6ae99SBarry Smith /*@
34aa219208SBarry Smith   DMDASetSizes - Sets the global sizes
3547c6ae99SBarry Smith 
36aa219208SBarry Smith   Logically Collective on DMDA
3747c6ae99SBarry Smith 
3847c6ae99SBarry Smith   Input Parameters:
39aa219208SBarry Smith + da - the DMDA
4047c6ae99SBarry Smith . M - the global X size (or PETSC_DECIDE)
4147c6ae99SBarry Smith . N - the global Y size (or PETSC_DECIDE)
4247c6ae99SBarry Smith - P - the global Z size (or PETSC_DECIDE)
4347c6ae99SBarry Smith 
4447c6ae99SBarry Smith   Level: intermediate
4547c6ae99SBarry Smith 
46aa219208SBarry Smith .seealso: DMDAGetSize(), PetscSplitOwnership()
4747c6ae99SBarry Smith @*/
487087cfbeSBarry Smith PetscErrorCode  DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P)
4947c6ae99SBarry Smith {
5047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
5147c6ae99SBarry Smith 
5247c6ae99SBarry Smith   PetscFunctionBegin;
5347c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
5447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,M,2);
5547c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,N,3);
5647c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,P,4);
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);
9047c6ae99SBarry Smith   dd->m = m;
9147c6ae99SBarry Smith   dd->n = n;
9247c6ae99SBarry Smith   dd->p = p;
9347c6ae99SBarry Smith   PetscFunctionReturn(0);
9447c6ae99SBarry Smith }
9547c6ae99SBarry Smith 
9647c6ae99SBarry Smith #undef __FUNCT__
973e24977eSEthan Coon #define __FUNCT__ "DMDASetBoundaryType"
9847c6ae99SBarry Smith /*@
99*1321219cSEthan Coon   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries.
10047c6ae99SBarry Smith 
10147c6ae99SBarry Smith   Not collective
10247c6ae99SBarry Smith 
10347c6ae99SBarry Smith   Input Parameter:
104aa219208SBarry Smith + da    - The DMDA
105*1321219cSEthan Coon - bx,by,bz - One of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC
10647c6ae99SBarry Smith 
10747c6ae99SBarry Smith   Level: intermediate
10847c6ae99SBarry Smith 
10947c6ae99SBarry Smith .keywords:  distributed array, periodicity
110db87c5ecSEthan Coon .seealso: DMDACreate(), DMDestroy(), DMDA, DMDABoundaryType
11147c6ae99SBarry Smith @*/
112*1321219cSEthan Coon PetscErrorCode  DMDASetBoundaryType(DM da,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz)
11347c6ae99SBarry Smith {
11447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
11547c6ae99SBarry Smith 
11647c6ae99SBarry Smith   PetscFunctionBegin;
11747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
118*1321219cSEthan Coon   PetscValidLogicalCollectiveInt(da,bx,2);
119*1321219cSEthan Coon   PetscValidLogicalCollectiveInt(da,by,3);
120*1321219cSEthan Coon   PetscValidLogicalCollectiveInt(da,bz,4);
121*1321219cSEthan Coon   dd->bx = bx;
122*1321219cSEthan Coon   dd->by = by;
123*1321219cSEthan Coon   dd->bz = bz;
12447c6ae99SBarry Smith   PetscFunctionReturn(0);
12547c6ae99SBarry Smith }
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith #undef __FUNCT__
128aa219208SBarry Smith #define __FUNCT__ "DMDASetDof"
12947c6ae99SBarry Smith /*@
130aa219208SBarry Smith   DMDASetDof - Sets the number of degrees of freedom per vertex
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith   Not collective
13347c6ae99SBarry Smith 
13447c6ae99SBarry Smith   Input Parameter:
135aa219208SBarry Smith + da  - The DMDA
13647c6ae99SBarry Smith - dof - Number of degrees of freedom
13747c6ae99SBarry Smith 
13847c6ae99SBarry Smith   Level: intermediate
13947c6ae99SBarry Smith 
14047c6ae99SBarry Smith .keywords:  distributed array, degrees of freedom
141aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
14247c6ae99SBarry Smith @*/
1437087cfbeSBarry Smith PetscErrorCode  DMDASetDof(DM da, int dof)
14447c6ae99SBarry Smith {
14547c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
14647c6ae99SBarry Smith 
14747c6ae99SBarry Smith   PetscFunctionBegin;
14847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
14947c6ae99SBarry Smith   dd->w = dof;
1501411c6eeSJed Brown   da->bs = dof;
15147c6ae99SBarry Smith   PetscFunctionReturn(0);
15247c6ae99SBarry Smith }
15347c6ae99SBarry Smith 
15447c6ae99SBarry Smith #undef __FUNCT__
155aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType"
15647c6ae99SBarry Smith /*@
157aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
15847c6ae99SBarry Smith 
159aa219208SBarry Smith   Logically Collective on DMDA
16047c6ae99SBarry Smith 
16147c6ae99SBarry Smith   Input Parameter:
162aa219208SBarry Smith + da    - The DMDA
163aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
16447c6ae99SBarry Smith 
16547c6ae99SBarry Smith   Level: intermediate
16647c6ae99SBarry Smith 
16747c6ae99SBarry Smith .keywords:  distributed array, stencil
168aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
16947c6ae99SBarry Smith @*/
1707087cfbeSBarry Smith PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
17147c6ae99SBarry Smith {
17247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
17347c6ae99SBarry Smith 
17447c6ae99SBarry Smith   PetscFunctionBegin;
17547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
17647c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,stype,2);
17747c6ae99SBarry Smith   dd->stencil_type = stype;
17847c6ae99SBarry Smith   PetscFunctionReturn(0);
17947c6ae99SBarry Smith }
18047c6ae99SBarry Smith 
18147c6ae99SBarry Smith #undef __FUNCT__
182aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth"
18347c6ae99SBarry Smith /*@
184aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
18547c6ae99SBarry Smith 
186aa219208SBarry Smith   Logically Collective on DMDA
18747c6ae99SBarry Smith 
18847c6ae99SBarry Smith   Input Parameter:
189aa219208SBarry Smith + da    - The DMDA
19047c6ae99SBarry Smith - width - The stencil width
19147c6ae99SBarry Smith 
19247c6ae99SBarry Smith   Level: intermediate
19347c6ae99SBarry Smith 
19447c6ae99SBarry Smith .keywords:  distributed array, stencil
195aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
19647c6ae99SBarry Smith @*/
1977087cfbeSBarry Smith PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
19847c6ae99SBarry Smith {
19947c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
20047c6ae99SBarry Smith 
20147c6ae99SBarry Smith   PetscFunctionBegin;
20247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
20347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,width,2);
20447c6ae99SBarry Smith   dd->s = width;
20547c6ae99SBarry Smith   PetscFunctionReturn(0);
20647c6ae99SBarry Smith }
20747c6ae99SBarry Smith 
20847c6ae99SBarry Smith #undef __FUNCT__
209aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
210aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
21147c6ae99SBarry Smith {
21247c6ae99SBarry Smith   PetscInt i,sum;
21347c6ae99SBarry Smith 
21447c6ae99SBarry Smith   PetscFunctionBegin;
21547c6ae99SBarry Smith   if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
21647c6ae99SBarry Smith   for (i=sum=0; i<m; i++) sum += lx[i];
21747c6ae99SBarry Smith   if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
21847c6ae99SBarry Smith   PetscFunctionReturn(0);
21947c6ae99SBarry Smith }
22047c6ae99SBarry Smith 
22147c6ae99SBarry Smith #undef __FUNCT__
222aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges"
22347c6ae99SBarry Smith /*@
224aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
22547c6ae99SBarry Smith 
226aa219208SBarry Smith   Logically Collective on DMDA
22747c6ae99SBarry Smith 
22847c6ae99SBarry Smith   Input Parameter:
229aa219208SBarry Smith + da - The DMDA
23047c6ae99SBarry 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
23147c6ae99SBarry 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
23247c6ae99SBarry 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.
23347c6ae99SBarry Smith 
23447c6ae99SBarry Smith   Level: intermediate
23547c6ae99SBarry Smith 
23647c6ae99SBarry Smith .keywords:  distributed array
237aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
23847c6ae99SBarry Smith @*/
2397087cfbeSBarry Smith PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
24047c6ae99SBarry Smith {
24147c6ae99SBarry Smith   PetscErrorCode ierr;
24247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
24347c6ae99SBarry Smith 
24447c6ae99SBarry Smith   PetscFunctionBegin;
24547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
24647c6ae99SBarry Smith   if (lx) {
24747c6ae99SBarry Smith     if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
248aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
24947c6ae99SBarry Smith     if (!dd->lx) {
25047c6ae99SBarry Smith       ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
25147c6ae99SBarry Smith     }
25247c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
25347c6ae99SBarry Smith   }
25447c6ae99SBarry Smith   if (ly) {
25547c6ae99SBarry Smith     if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
256aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
25747c6ae99SBarry Smith     if (!dd->ly) {
25847c6ae99SBarry Smith       ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
25947c6ae99SBarry Smith     }
26047c6ae99SBarry Smith     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
26147c6ae99SBarry Smith   }
26247c6ae99SBarry Smith   if (lz) {
26347c6ae99SBarry Smith     if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
264aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
26547c6ae99SBarry Smith     if (!dd->lz) {
26647c6ae99SBarry Smith       ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
26747c6ae99SBarry Smith     }
26847c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
26947c6ae99SBarry Smith   }
27047c6ae99SBarry Smith   PetscFunctionReturn(0);
27147c6ae99SBarry Smith }
27247c6ae99SBarry Smith 
27347c6ae99SBarry Smith #undef __FUNCT__
274aa219208SBarry Smith #define __FUNCT__ "DMDACreateOwnershipRanges"
27547c6ae99SBarry Smith /*
276aa219208SBarry Smith  Ensure that da->lx, ly, and lz exist.  Collective on DMDA.
27747c6ae99SBarry Smith */
2787087cfbeSBarry Smith PetscErrorCode  DMDACreateOwnershipRanges(DM da)
27947c6ae99SBarry Smith {
28047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
28147c6ae99SBarry Smith   PetscErrorCode ierr;
28247c6ae99SBarry Smith   PetscInt       n;
28347c6ae99SBarry Smith   MPI_Comm       comm;
28447c6ae99SBarry Smith   PetscMPIInt    size;
28547c6ae99SBarry Smith 
28647c6ae99SBarry Smith   PetscFunctionBegin;
28747c6ae99SBarry Smith   if (!dd->lx) {
28847c6ae99SBarry Smith     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&dd->lx);CHKERRQ(ierr);
289aa219208SBarry Smith     ierr = DMDAGetProcessorSubset(da,DMDA_X,dd->xs,&comm);CHKERRQ(ierr);
29047c6ae99SBarry Smith     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
29147c6ae99SBarry Smith     n = (dd->xe-dd->xs)/dd->w;
29247c6ae99SBarry Smith     ierr = MPI_Allgather(&n,1,MPIU_INT,dd->lx,1,MPIU_INT,comm);CHKERRQ(ierr);
29347c6ae99SBarry Smith   }
29447c6ae99SBarry Smith   if (dd->dim > 1 && !dd->ly) {
29547c6ae99SBarry Smith     ierr = PetscMalloc(dd->n*sizeof(PetscInt),&dd->ly);CHKERRQ(ierr);
296aa219208SBarry Smith     ierr = DMDAGetProcessorSubset(da,DMDA_Y,dd->ys,&comm);CHKERRQ(ierr);
29747c6ae99SBarry Smith     n = dd->ye-dd->ys;
29847c6ae99SBarry Smith     ierr = MPI_Allgather(&n,1,MPIU_INT,dd->ly,1,MPIU_INT,comm);CHKERRQ(ierr);
29947c6ae99SBarry Smith   }
30047c6ae99SBarry Smith   if (dd->dim > 2 && !dd->lz) {
30147c6ae99SBarry Smith     ierr = PetscMalloc(dd->p*sizeof(PetscInt),&dd->lz);CHKERRQ(ierr);
302aa219208SBarry Smith     ierr = DMDAGetProcessorSubset(da,DMDA_Z,dd->zs,&comm);CHKERRQ(ierr);
30347c6ae99SBarry Smith     n = dd->ze-dd->zs;
30447c6ae99SBarry Smith     ierr = MPI_Allgather(&n,1,MPIU_INT,dd->lz,1,MPIU_INT,comm);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
31394013140SBarry Smith           returned by DMGetInterpolation()
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 
323aa219208SBarry Smith    Notes: you should call this on the coarser of the two DMDAs you pass to DMGetInterpolation()
32447c6ae99SBarry Smith 
32547c6ae99SBarry Smith .keywords:  distributed array, interpolation
32647c6ae99SBarry Smith 
327aa219208SBarry 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 
34047c6ae99SBarry Smith 
34147c6ae99SBarry Smith #undef __FUNCT__
342aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors"
34347c6ae99SBarry Smith /*@C
344aa219208SBarry Smith       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
34547c6ae99SBarry Smith         processes neighbors.
34647c6ae99SBarry Smith 
34747c6ae99SBarry Smith     Not Collective
34847c6ae99SBarry Smith 
34947c6ae99SBarry Smith    Input Parameter:
350aa219208SBarry Smith .     da - the DMDA object
35147c6ae99SBarry Smith 
35247c6ae99SBarry Smith    Output Parameters:
35347c6ae99SBarry Smith .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
35447c6ae99SBarry Smith               this process itself is in the list
35547c6ae99SBarry Smith 
35647c6ae99SBarry Smith    Notes: In 2d the array is of length 9, in 3d of length 27
35747c6ae99SBarry Smith           Not supported in 1d
358aa219208SBarry Smith           Do not free the array, it is freed when the DMDA is destroyed.
35947c6ae99SBarry Smith 
36047c6ae99SBarry Smith    Fortran Notes: In fortran you must pass in an array of the appropriate length.
36147c6ae99SBarry Smith 
36247c6ae99SBarry Smith    Level: intermediate
36347c6ae99SBarry Smith 
36447c6ae99SBarry Smith @*/
3657087cfbeSBarry Smith PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
36647c6ae99SBarry Smith {
36747c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
36847c6ae99SBarry Smith   PetscFunctionBegin;
36947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
37047c6ae99SBarry Smith   *ranks = dd->neighbors;
37147c6ae99SBarry Smith   PetscFunctionReturn(0);
37247c6ae99SBarry Smith }
37347c6ae99SBarry Smith 
37447c6ae99SBarry Smith /*@C
375aa219208SBarry Smith       DMDASetElementType - Sets the element type to be returned by DMGetElements()
37647c6ae99SBarry Smith 
37747c6ae99SBarry Smith     Not Collective
37847c6ae99SBarry Smith 
37947c6ae99SBarry Smith    Input Parameter:
380aa219208SBarry Smith .     da - the DMDA object
38147c6ae99SBarry Smith 
38247c6ae99SBarry Smith    Output Parameters:
383aa219208SBarry Smith .     etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1
38447c6ae99SBarry Smith 
38547c6ae99SBarry Smith    Level: intermediate
38647c6ae99SBarry Smith 
387aa219208SBarry Smith .seealso: DMDAElementType, DMDAGetElementType(), DMGetElements(), DMRestoreElements()
38847c6ae99SBarry Smith @*/
38947c6ae99SBarry Smith #undef __FUNCT__
390aa219208SBarry Smith #define __FUNCT__ "DMDASetElementType"
3917087cfbeSBarry Smith PetscErrorCode  DMDASetElementType(DM da, DMDAElementType etype)
39247c6ae99SBarry Smith {
39347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
39447c6ae99SBarry Smith   PetscErrorCode ierr;
39547c6ae99SBarry Smith 
39647c6ae99SBarry Smith   PetscFunctionBegin;
39747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
39847c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,etype,2);
39947c6ae99SBarry Smith   if (dd->elementtype != etype) {
40047c6ae99SBarry Smith     ierr = PetscFree(dd->e);CHKERRQ(ierr);
40147c6ae99SBarry Smith     dd->elementtype = etype;
40247c6ae99SBarry Smith     dd->ne          = 0;
40347c6ae99SBarry Smith     dd->e           = PETSC_NULL;
40447c6ae99SBarry Smith   }
40547c6ae99SBarry Smith   PetscFunctionReturn(0);
40647c6ae99SBarry Smith }
40747c6ae99SBarry Smith 
40847c6ae99SBarry Smith /*@C
409aa219208SBarry Smith       DMDAGetElementType - Gets the element type to be returned by DMGetElements()
41047c6ae99SBarry Smith 
41147c6ae99SBarry Smith     Not Collective
41247c6ae99SBarry Smith 
41347c6ae99SBarry Smith    Input Parameter:
414aa219208SBarry Smith .     da - the DMDA object
41547c6ae99SBarry Smith 
41647c6ae99SBarry Smith    Output Parameters:
417aa219208SBarry Smith .     etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1
41847c6ae99SBarry Smith 
41947c6ae99SBarry Smith    Level: intermediate
42047c6ae99SBarry Smith 
421aa219208SBarry Smith .seealso: DMDAElementType, DMDASetElementType(), DMGetElements(), DMRestoreElements()
42247c6ae99SBarry Smith @*/
42347c6ae99SBarry Smith #undef __FUNCT__
424aa219208SBarry Smith #define __FUNCT__ "DMDAGetElementType"
4257087cfbeSBarry Smith PetscErrorCode  DMDAGetElementType(DM da, DMDAElementType *etype)
42647c6ae99SBarry Smith {
42747c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
42847c6ae99SBarry Smith   PetscFunctionBegin;
42947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
43047c6ae99SBarry Smith   PetscValidPointer(etype,2);
43147c6ae99SBarry Smith   *etype = dd->elementtype;
43247c6ae99SBarry Smith   PetscFunctionReturn(0);
43347c6ae99SBarry Smith }
43447c6ae99SBarry Smith 
43547c6ae99SBarry Smith #undef __FUNCT__
43647c6ae99SBarry Smith #define __FUNCT__ "DMGetElements"
43747c6ae99SBarry Smith /*@C
43847c6ae99SBarry Smith       DMGetElements - Gets an array containing the indices (in local coordinates)
43947c6ae99SBarry Smith                  of all the local elements
44047c6ae99SBarry Smith 
44147c6ae99SBarry Smith     Not Collective
44247c6ae99SBarry Smith 
44347c6ae99SBarry Smith    Input Parameter:
44447c6ae99SBarry Smith .     dm - the DM object
44547c6ae99SBarry Smith 
44647c6ae99SBarry Smith    Output Parameters:
447454e267fSLisandro Dalcin +     nel - number of local elements
448454e267fSLisandro Dalcin .     nen - number of element nodes
44947c6ae99SBarry Smith -     e - the indices of the elements vertices
45047c6ae99SBarry Smith 
45147c6ae99SBarry Smith    Level: intermediate
45247c6ae99SBarry Smith 
45347c6ae99SBarry Smith .seealso: DMElementType, DMSetElementType(), DMRestoreElements()
45447c6ae99SBarry Smith @*/
4557087cfbeSBarry Smith PetscErrorCode  DMGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
45647c6ae99SBarry Smith {
45747c6ae99SBarry Smith   PetscErrorCode ierr;
45847c6ae99SBarry Smith   PetscFunctionBegin;
45947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
460454e267fSLisandro Dalcin   PetscValidIntPointer(nel,2);
461454e267fSLisandro Dalcin   PetscValidIntPointer(nen,3);
462454e267fSLisandro Dalcin   PetscValidPointer(e,4);
463454e267fSLisandro Dalcin   ierr = (dm->ops->getelements)(dm,nel,nen,e);CHKERRQ(ierr);
46447c6ae99SBarry Smith   PetscFunctionReturn(0);
46547c6ae99SBarry Smith }
46647c6ae99SBarry Smith 
46747c6ae99SBarry Smith #undef __FUNCT__
46847c6ae99SBarry Smith #define __FUNCT__ "DMRestoreElements"
46947c6ae99SBarry Smith /*@C
47047c6ae99SBarry Smith       DMRestoreElements - Returns an array containing the indices (in local coordinates)
47147c6ae99SBarry Smith                  of all the local elements obtained with DMGetElements()
47247c6ae99SBarry Smith 
47347c6ae99SBarry Smith     Not Collective
47447c6ae99SBarry Smith 
47547c6ae99SBarry Smith    Input Parameter:
47647c6ae99SBarry Smith +     dm - the DM object
477454e267fSLisandro Dalcin .     nel - number of local elements
478454e267fSLisandro Dalcin .     nen - number of element nodes
47947c6ae99SBarry Smith -     e - the indices of the elements vertices
48047c6ae99SBarry Smith 
48147c6ae99SBarry Smith    Level: intermediate
48247c6ae99SBarry Smith 
48347c6ae99SBarry Smith .seealso: DMElementType, DMSetElementType(), DMGetElements()
48447c6ae99SBarry Smith @*/
4857087cfbeSBarry Smith PetscErrorCode  DMRestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
48647c6ae99SBarry Smith {
48747c6ae99SBarry Smith   PetscErrorCode ierr;
48847c6ae99SBarry Smith   PetscFunctionBegin;
48947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
490454e267fSLisandro Dalcin   PetscValidIntPointer(nel,2);
491454e267fSLisandro Dalcin   PetscValidIntPointer(nen,3);
492454e267fSLisandro Dalcin   PetscValidPointer(e,4);
49347c6ae99SBarry Smith   if (dm->ops->restoreelements) {
494454e267fSLisandro Dalcin     ierr = (dm->ops->restoreelements)(dm,nel,nen,e);CHKERRQ(ierr);
49547c6ae99SBarry Smith   }
49647c6ae99SBarry Smith   PetscFunctionReturn(0);
49747c6ae99SBarry Smith }
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith #undef __FUNCT__
500aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges"
50147c6ae99SBarry Smith /*@C
502aa219208SBarry Smith       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
50347c6ae99SBarry Smith 
50447c6ae99SBarry Smith     Not Collective
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith    Input Parameter:
507aa219208SBarry Smith .     da - the DMDA object
50847c6ae99SBarry Smith 
50947c6ae99SBarry Smith    Output Parameter:
51047c6ae99SBarry Smith +     lx - ownership along x direction (optional)
51147c6ae99SBarry Smith .     ly - ownership along y direction (optional)
51247c6ae99SBarry Smith -     lz - ownership along z direction (optional)
51347c6ae99SBarry Smith 
51447c6ae99SBarry Smith    Level: intermediate
51547c6ae99SBarry Smith 
516aa219208SBarry Smith     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
51747c6ae99SBarry Smith 
51847c6ae99SBarry Smith     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
519aa219208SBarry Smith     eighth arguments from DMDAGetInfo()
52047c6ae99SBarry Smith 
52147c6ae99SBarry Smith      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
522aa219208SBarry Smith     DMDA they came from still exists (has not been destroyed).
52347c6ae99SBarry Smith 
524aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
52547c6ae99SBarry Smith @*/
5267087cfbeSBarry Smith PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
52747c6ae99SBarry Smith {
52847c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
52947c6ae99SBarry Smith 
53047c6ae99SBarry Smith   PetscFunctionBegin;
53147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
53247c6ae99SBarry Smith   if (lx) *lx = dd->lx;
53347c6ae99SBarry Smith   if (ly) *ly = dd->ly;
53447c6ae99SBarry Smith   if (lz) *lz = dd->lz;
53547c6ae99SBarry Smith   PetscFunctionReturn(0);
53647c6ae99SBarry Smith }
53747c6ae99SBarry Smith 
53847c6ae99SBarry Smith #undef __FUNCT__
539aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor"
54047c6ae99SBarry Smith /*@
541aa219208SBarry Smith      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
54247c6ae99SBarry Smith 
543aa219208SBarry Smith     Logically Collective on DMDA
54447c6ae99SBarry Smith 
54547c6ae99SBarry Smith   Input Parameters:
546aa219208SBarry Smith +    da - the DMDA object
54747c6ae99SBarry Smith .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
54847c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
54947c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
55047c6ae99SBarry Smith 
55147c6ae99SBarry Smith   Options Database:
55247c6ae99SBarry Smith +  -da_refine_x - refinement ratio in x direction
55347c6ae99SBarry Smith .  -da_refine_y - refinement ratio in y direction
55447c6ae99SBarry Smith -  -da_refine_z - refinement ratio in z direction
55547c6ae99SBarry Smith 
55647c6ae99SBarry Smith   Level: intermediate
55747c6ae99SBarry Smith 
55847c6ae99SBarry Smith     Notes: Pass PETSC_IGNORE to leave a value unchanged
55947c6ae99SBarry Smith 
560aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor()
56147c6ae99SBarry Smith @*/
5627087cfbeSBarry Smith PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
56347c6ae99SBarry Smith {
56447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
56547c6ae99SBarry Smith 
56647c6ae99SBarry Smith   PetscFunctionBegin;
56747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
56847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_x,2);
56947c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_y,3);
57047c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_z,4);
57147c6ae99SBarry Smith 
57247c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
57347c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
57447c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
57547c6ae99SBarry Smith   PetscFunctionReturn(0);
57647c6ae99SBarry Smith }
57747c6ae99SBarry Smith 
57847c6ae99SBarry Smith #undef __FUNCT__
579aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor"
58047c6ae99SBarry Smith /*@C
581aa219208SBarry Smith      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
58247c6ae99SBarry Smith 
58347c6ae99SBarry Smith     Not Collective
58447c6ae99SBarry Smith 
58547c6ae99SBarry Smith   Input Parameter:
586aa219208SBarry Smith .    da - the DMDA object
58747c6ae99SBarry Smith 
58847c6ae99SBarry Smith   Output Parameters:
58947c6ae99SBarry Smith +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
59047c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
59147c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
59247c6ae99SBarry Smith 
59347c6ae99SBarry Smith   Level: intermediate
59447c6ae99SBarry Smith 
59547c6ae99SBarry Smith     Notes: Pass PETSC_NULL for values you do not need
59647c6ae99SBarry Smith 
597aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor()
59847c6ae99SBarry Smith @*/
5997087cfbeSBarry Smith PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
60047c6ae99SBarry Smith {
60147c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
60247c6ae99SBarry Smith 
60347c6ae99SBarry Smith   PetscFunctionBegin;
60447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
60547c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
60647c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
60747c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
60847c6ae99SBarry Smith   PetscFunctionReturn(0);
60947c6ae99SBarry Smith }
61047c6ae99SBarry Smith 
61147c6ae99SBarry Smith #undef __FUNCT__
612aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix"
61347c6ae99SBarry Smith /*@C
614aa219208SBarry Smith      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
61547c6ae99SBarry Smith 
616aa219208SBarry Smith     Logically Collective on DMDA
61747c6ae99SBarry Smith 
61847c6ae99SBarry Smith   Input Parameters:
619aa219208SBarry Smith +    da - the DMDA object
620aa219208SBarry Smith -    f - the function that allocates the matrix for that specific DMDA
62147c6ae99SBarry Smith 
62247c6ae99SBarry Smith   Level: developer
62347c6ae99SBarry Smith 
624aa219208SBarry Smith    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
62547c6ae99SBarry Smith        the diagonal and off-diagonal blocks of the matrix
62647c6ae99SBarry Smith 
627aa219208SBarry Smith .seealso: DMGetMatrix(), DMDASetBlockFills()
62847c6ae99SBarry Smith @*/
6297087cfbeSBarry Smith PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, const MatType,Mat*))
63047c6ae99SBarry Smith {
63147c6ae99SBarry Smith   PetscFunctionBegin;
63247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
63347c6ae99SBarry Smith   da->ops->getmatrix = f;
63447c6ae99SBarry Smith   PetscFunctionReturn(0);
63547c6ae99SBarry Smith }
63647c6ae99SBarry Smith 
63747c6ae99SBarry Smith #undef __FUNCT__
63894013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges"
63947c6ae99SBarry Smith /*
64047c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
64147c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
64247c6ae99SBarry Smith 
64347c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
64447c6ae99SBarry Smith */
64594013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
64647c6ae99SBarry Smith {
64747c6ae99SBarry Smith   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
64847c6ae99SBarry Smith   PetscErrorCode ierr;
64947c6ae99SBarry Smith 
65047c6ae99SBarry Smith   PetscFunctionBegin;
65147c6ae99SBarry Smith   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
65247c6ae99SBarry Smith   if (ratio == 1) {
65347c6ae99SBarry Smith     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
65447c6ae99SBarry Smith     PetscFunctionReturn(0);
65547c6ae99SBarry Smith   }
65647c6ae99SBarry Smith   for (i=0; i<m; i++) totalc += lc[i];
65747c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
65847c6ae99SBarry Smith   for (i=0; i<m; i++) {
65947c6ae99SBarry Smith     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
66047c6ae99SBarry Smith     if (i == m-1) lf[i] = want;
66147c6ae99SBarry Smith     else {
6627aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
6637aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
6647aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
6657aca7175SJed Brown       while ((startf+want)/ratio < nextc - stencil_width) want++;
6667aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
6677aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
6687aca7175SJed Brown       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
6697aca7175SJed Brown       /* Make sure all constraints are satisfied */
6707aca7175SJed Brown       if (want < 0 || want > remaining
6717aca7175SJed Brown           || ((startf+want)/ratio < nextc - stencil_width)
6727aca7175SJed Brown           || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width))
6737aca7175SJed Brown         SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range");
67447c6ae99SBarry Smith     }
67547c6ae99SBarry Smith     lf[i] = want;
67647c6ae99SBarry Smith     startc += lc[i];
67747c6ae99SBarry Smith     startf += lf[i];
67847c6ae99SBarry Smith     remaining -= lf[i];
67947c6ae99SBarry Smith   }
68047c6ae99SBarry Smith   PetscFunctionReturn(0);
68147c6ae99SBarry Smith }
68247c6ae99SBarry Smith 
68347c6ae99SBarry Smith #undef __FUNCT__
6842be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges"
6852be375d4SJed Brown /*
6862be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
6872be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
6882be375d4SJed Brown 
6892be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
6902be375d4SJed Brown */
6912be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
6922be375d4SJed Brown {
6932be375d4SJed Brown   PetscInt       i,totalf,remaining,startc,startf;
6942be375d4SJed Brown   PetscErrorCode ierr;
6952be375d4SJed Brown 
6962be375d4SJed Brown   PetscFunctionBegin;
6972be375d4SJed Brown   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
6982be375d4SJed Brown   if (ratio == 1) {
6992be375d4SJed Brown     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
7002be375d4SJed Brown     PetscFunctionReturn(0);
7012be375d4SJed Brown   }
7022be375d4SJed Brown   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
7032be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
7042be375d4SJed Brown   for (i=0,startc=0,startf=0; i<m; i++) {
7052be375d4SJed Brown     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
7062be375d4SJed Brown     if (i == m-1) lc[i] = want;
7072be375d4SJed Brown     else {
7082be375d4SJed Brown       const PetscInt nextf = startf+lf[i];
7092be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
7102be375d4SJed Brown        * node is within one stencil width. */
7112be375d4SJed Brown       while (nextf/ratio < startc+want-stencil_width) want--;
7122be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
7132be375d4SJed Brown        * fine node is within one stencil width. */
7142be375d4SJed Brown       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
7152be375d4SJed Brown       if (want < 0 || want > remaining
7162be375d4SJed Brown           || (nextf/ratio < startc+want-stencil_width)
7172be375d4SJed Brown           || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width))
7182be375d4SJed Brown         SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range");
7192be375d4SJed Brown     }
7202be375d4SJed Brown     lc[i] = want;
7212be375d4SJed Brown     startc += lc[i];
7222be375d4SJed Brown     startf += lf[i];
7232be375d4SJed Brown     remaining -= lc[i];
7242be375d4SJed Brown   }
7252be375d4SJed Brown   PetscFunctionReturn(0);
7262be375d4SJed Brown }
7272be375d4SJed Brown 
7282be375d4SJed Brown #undef __FUNCT__
72994013140SBarry Smith #define __FUNCT__ "DMRefine_DA"
7307087cfbeSBarry Smith PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
73147c6ae99SBarry Smith {
73247c6ae99SBarry Smith   PetscErrorCode ierr;
73347c6ae99SBarry Smith   PetscInt       M,N,P;
7349a42bb27SBarry Smith   DM             da2;
73547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
73647c6ae99SBarry Smith 
73747c6ae99SBarry Smith   PetscFunctionBegin;
73847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
73947c6ae99SBarry Smith   PetscValidPointer(daref,3);
74047c6ae99SBarry Smith 
741*1321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
74247c6ae99SBarry Smith     M = dd->refine_x*dd->M;
74347c6ae99SBarry Smith   } else {
74447c6ae99SBarry Smith     M = 1 + dd->refine_x*(dd->M - 1);
74547c6ae99SBarry Smith   }
746*1321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
74747c6ae99SBarry Smith     N = dd->refine_y*dd->N;
74847c6ae99SBarry Smith   } else {
74947c6ae99SBarry Smith     N = 1 + dd->refine_y*(dd->N - 1);
75047c6ae99SBarry Smith   }
751*1321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
75247c6ae99SBarry Smith     P = dd->refine_z*dd->P;
75347c6ae99SBarry Smith   } else {
75447c6ae99SBarry Smith     P = 1 + dd->refine_z*(dd->P - 1);
75547c6ae99SBarry Smith   }
756aa219208SBarry Smith   ierr = DMDACreateOwnershipRanges(da);CHKERRQ(ierr);
75747c6ae99SBarry Smith   if (dd->dim == 3) {
75847c6ae99SBarry Smith     PetscInt *lx,*ly,*lz;
75947c6ae99SBarry Smith     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
760*1321219cSEthan 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);
761*1321219cSEthan 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);
762*1321219cSEthan 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);
763*1321219cSEthan Coon     ierr = DMDACreate3d(((PetscObject)da)->comm,dd->bx,dd->by,dd->bz,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr);
76447c6ae99SBarry Smith     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
76547c6ae99SBarry Smith   } else if (dd->dim == 2) {
76647c6ae99SBarry Smith     PetscInt *lx,*ly;
76747c6ae99SBarry Smith     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
768*1321219cSEthan 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);
769*1321219cSEthan 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);
770*1321219cSEthan Coon     ierr = DMDACreate2d(((PetscObject)da)->comm,dd->bx,dd->by,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr);
77147c6ae99SBarry Smith     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
77247c6ae99SBarry Smith   } else if (dd->dim == 1) {
77347c6ae99SBarry Smith     PetscInt *lx;
77447c6ae99SBarry Smith     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
775*1321219cSEthan 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);
776*1321219cSEthan Coon     ierr = DMDACreate1d(((PetscObject)da)->comm,dd->bx,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr);
77747c6ae99SBarry Smith     ierr = PetscFree(lx);CHKERRQ(ierr);
77847c6ae99SBarry Smith   }
77947c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
78047c6ae99SBarry Smith 
78147c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
78247c6ae99SBarry Smith   da2->ops->getmatrix        = da->ops->getmatrix;
78347c6ae99SBarry Smith   da2->ops->getinterpolation = da->ops->getinterpolation;
78447c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
78547c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
78647c6ae99SBarry Smith 
78747c6ae99SBarry Smith   /* copy fill information if given */
78847c6ae99SBarry Smith   if (dd->dfill) {
78947c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
79047c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
79147c6ae99SBarry Smith   }
79247c6ae99SBarry Smith   if (dd->ofill) {
79347c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
79447c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
79547c6ae99SBarry Smith   }
79647c6ae99SBarry Smith   /* copy the refine information */
79747c6ae99SBarry Smith   dd2->refine_x = dd->refine_x;
79847c6ae99SBarry Smith   dd2->refine_y = dd->refine_y;
79947c6ae99SBarry Smith   dd2->refine_z = dd->refine_z;
80047c6ae99SBarry Smith 
80147c6ae99SBarry Smith   /* copy vector type information */
80247c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
80347c6ae99SBarry Smith   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
804ddcf8b74SDave May 
805ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
806ddcf8b74SDave May   if (dd->coordinates) {
807ddcf8b74SDave May     DM  cdaf,cdac;
808ddcf8b74SDave May     Vec coordsc,coordsf;
809ddcf8b74SDave May     Mat II;
810ddcf8b74SDave May 
811ddcf8b74SDave May     ierr = DMDAGetCoordinateDA(da,&cdac);CHKERRQ(ierr);
812ddcf8b74SDave May     ierr = DMDAGetCoordinates(da,&coordsc);CHKERRQ(ierr);
813b61d3410SDave May     ierr = DMDAGetCoordinateDA(da2,&cdaf);CHKERRQ(ierr);
814b61d3410SDave May     /* force creation of the coordinate vector */
815b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
816ddcf8b74SDave May     ierr = DMDAGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
817ddcf8b74SDave May     ierr = DMGetInterpolation(cdac,cdaf,&II,PETSC_NULL);CHKERRQ(ierr);
818ddcf8b74SDave May     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
819ddcf8b74SDave May     ierr = MatDestroy(II);CHKERRQ(ierr);
820ddcf8b74SDave May   }
82147c6ae99SBarry Smith   *daref = da2;
82247c6ae99SBarry Smith   PetscFunctionReturn(0);
82347c6ae99SBarry Smith }
82447c6ae99SBarry Smith 
82547c6ae99SBarry Smith #undef __FUNCT__
82694013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA"
8277087cfbeSBarry Smith PetscErrorCode  DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
82847c6ae99SBarry Smith {
82947c6ae99SBarry Smith   PetscErrorCode ierr;
83047c6ae99SBarry Smith   PetscInt       M,N,P;
8319a42bb27SBarry Smith   DM             da2;
83247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
83347c6ae99SBarry Smith 
83447c6ae99SBarry Smith   PetscFunctionBegin;
83547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
83647c6ae99SBarry Smith   PetscValidPointer(daref,3);
83747c6ae99SBarry Smith 
838*1321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
83947c6ae99SBarry Smith     M = dd->M / dd->refine_x;
84047c6ae99SBarry Smith   } else {
84147c6ae99SBarry Smith     M = 1 + (dd->M - 1) / dd->refine_x;
84247c6ae99SBarry Smith   }
843*1321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
84447c6ae99SBarry Smith     N = dd->N / dd->refine_y;
84547c6ae99SBarry Smith   } else {
84647c6ae99SBarry Smith     N = 1 + (dd->N - 1) / dd->refine_y;
84747c6ae99SBarry Smith   }
848*1321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
84947c6ae99SBarry Smith     P = dd->P / dd->refine_z;
85047c6ae99SBarry Smith   } else {
85147c6ae99SBarry Smith     P = 1 + (dd->P - 1) / dd->refine_z;
85247c6ae99SBarry Smith   }
8532be375d4SJed Brown   ierr = DMDACreateOwnershipRanges(da);CHKERRQ(ierr);
85447c6ae99SBarry Smith   if (dd->dim == 3) {
8552be375d4SJed Brown     PetscInt *lx,*ly,*lz;
8562be375d4SJed Brown     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
857*1321219cSEthan Coon     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
858*1321219cSEthan Coon     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
859*1321219cSEthan Coon     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
860*1321219cSEthan Coon     ierr = DMDACreate3d(((PetscObject)da)->comm,dd->bx,dd->by,dd->bz,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr);
8612be375d4SJed Brown     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
86247c6ae99SBarry Smith   } else if (dd->dim == 2) {
8632be375d4SJed Brown     PetscInt *lx,*ly;
8642be375d4SJed Brown     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
865*1321219cSEthan Coon     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
866*1321219cSEthan Coon     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
867*1321219cSEthan Coon     ierr = DMDACreate2d(((PetscObject)da)->comm,dd->bx,dd->by,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr);
8682be375d4SJed Brown     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
86947c6ae99SBarry Smith   } else if (dd->dim == 1) {
8702be375d4SJed Brown     PetscInt *lx;
8712be375d4SJed Brown     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
872*1321219cSEthan Coon     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
873*1321219cSEthan Coon     ierr = DMDACreate1d(((PetscObject)da)->comm,dd->bx,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr);
8742be375d4SJed Brown     ierr = PetscFree(lx);CHKERRQ(ierr);
87547c6ae99SBarry Smith   }
87647c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
87747c6ae99SBarry Smith 
87847c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
87947c6ae99SBarry Smith   da2->ops->getmatrix        = da->ops->getmatrix;
88047c6ae99SBarry Smith   da2->ops->getinterpolation = da->ops->getinterpolation;
88147c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
88247c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
88347c6ae99SBarry Smith 
88447c6ae99SBarry Smith   /* copy fill information if given */
88547c6ae99SBarry Smith   if (dd->dfill) {
88647c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
88747c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
88847c6ae99SBarry Smith   }
88947c6ae99SBarry Smith   if (dd->ofill) {
89047c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
89147c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
89247c6ae99SBarry Smith   }
89347c6ae99SBarry Smith   /* copy the refine information */
89447c6ae99SBarry Smith   dd2->refine_x = dd->refine_x;
89547c6ae99SBarry Smith   dd2->refine_y = dd->refine_y;
89647c6ae99SBarry Smith   dd2->refine_z = dd->refine_z;
89747c6ae99SBarry Smith 
89847c6ae99SBarry Smith   /* copy vector type information */
89947c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
90047c6ae99SBarry Smith   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
90147c6ae99SBarry Smith 
902ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
903ddcf8b74SDave May   if (dd->coordinates) {
904ddcf8b74SDave May     DM  cdaf,cdac;
905ddcf8b74SDave May     Vec coordsc,coordsf;
906ddcf8b74SDave May     VecScatter inject;
907ddcf8b74SDave May 
908ddcf8b74SDave May     ierr = DMDAGetCoordinateDA(da,&cdaf);CHKERRQ(ierr);
909ddcf8b74SDave May     ierr = DMDAGetCoordinates(da,&coordsf);CHKERRQ(ierr);
910b61d3410SDave May     ierr = DMDAGetCoordinateDA(da2,&cdac);CHKERRQ(ierr);
911b61d3410SDave May     /* force creation of the coordinate vector */
912b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
913ddcf8b74SDave May     ierr = DMDAGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
914ddcf8b74SDave May 
915ddcf8b74SDave May     ierr = DMGetInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
916ddcf8b74SDave May     ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
917ddcf8b74SDave May     ierr = VecScatterEnd(inject  ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
918ddcf8b74SDave May     ierr = VecScatterDestroy(inject);CHKERRQ(ierr);
919ddcf8b74SDave May   }
92047c6ae99SBarry Smith   *daref = da2;
92147c6ae99SBarry Smith   PetscFunctionReturn(0);
92247c6ae99SBarry Smith }
92347c6ae99SBarry Smith 
92447c6ae99SBarry Smith #undef __FUNCT__
92594013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA"
9267087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
92747c6ae99SBarry Smith {
92847c6ae99SBarry Smith   PetscErrorCode ierr;
92947c6ae99SBarry Smith   PetscInt       i,n,*refx,*refy,*refz;
93047c6ae99SBarry Smith 
93147c6ae99SBarry Smith   PetscFunctionBegin;
93247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
93347c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
93447c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
93547c6ae99SBarry Smith   PetscValidPointer(daf,3);
93647c6ae99SBarry Smith 
937aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
93847c6ae99SBarry Smith   ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr);
93947c6ae99SBarry Smith   for (i=0; i<nlevels; i++) {
940aa219208SBarry Smith     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
94147c6ae99SBarry Smith   }
94247c6ae99SBarry Smith   n = nlevels;
94347c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
94447c6ae99SBarry Smith   n = nlevels;
94547c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
94647c6ae99SBarry Smith   n = nlevels;
94747c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
94847c6ae99SBarry Smith 
949aa219208SBarry Smith   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
95094013140SBarry Smith   ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr);
95147c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
952aa219208SBarry Smith     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
95394013140SBarry Smith     ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr);
95447c6ae99SBarry Smith   }
95547c6ae99SBarry Smith   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
95647c6ae99SBarry Smith   PetscFunctionReturn(0);
95747c6ae99SBarry Smith }
95847c6ae99SBarry Smith 
95947c6ae99SBarry Smith #undef __FUNCT__
96094013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA"
9617087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
96247c6ae99SBarry Smith {
96347c6ae99SBarry Smith   PetscErrorCode ierr;
96447c6ae99SBarry Smith   PetscInt       i;
96547c6ae99SBarry Smith 
96647c6ae99SBarry Smith   PetscFunctionBegin;
96747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
96847c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
96947c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
97047c6ae99SBarry Smith   PetscValidPointer(dac,3);
97194013140SBarry Smith   ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr);
97247c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
97394013140SBarry Smith     ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr);
97447c6ae99SBarry Smith   }
97547c6ae99SBarry Smith   PetscFunctionReturn(0);
97647c6ae99SBarry Smith }
977