xref: /petsc/src/dm/impls/da/da.c (revision c6db04a5321582041def2b1e244c75985478b3ef)
1*c6db04a5SJed Brown #include <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);
5647c6ae99SBarry Smith 
5747c6ae99SBarry Smith   dd->M = M;
5847c6ae99SBarry Smith   dd->N = N;
5947c6ae99SBarry Smith   dd->P = P;
6047c6ae99SBarry Smith   PetscFunctionReturn(0);
6147c6ae99SBarry Smith }
6247c6ae99SBarry Smith 
6347c6ae99SBarry Smith #undef __FUNCT__
64aa219208SBarry Smith #define __FUNCT__ "DMDASetNumProcs"
6547c6ae99SBarry Smith /*@
66aa219208SBarry Smith   DMDASetNumProcs - Sets the number of processes in each dimension
6747c6ae99SBarry Smith 
68aa219208SBarry Smith   Logically Collective on DMDA
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith   Input Parameters:
71aa219208SBarry Smith + da - the DMDA
7247c6ae99SBarry Smith . m - the number of X procs (or PETSC_DECIDE)
7347c6ae99SBarry Smith . n - the number of Y procs (or PETSC_DECIDE)
7447c6ae99SBarry Smith - p - the number of Z procs (or PETSC_DECIDE)
7547c6ae99SBarry Smith 
7647c6ae99SBarry Smith   Level: intermediate
7747c6ae99SBarry Smith 
78aa219208SBarry Smith .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership()
7947c6ae99SBarry Smith @*/
807087cfbeSBarry Smith PetscErrorCode  DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p)
8147c6ae99SBarry Smith {
8247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
8347c6ae99SBarry Smith 
8447c6ae99SBarry Smith   PetscFunctionBegin;
8547c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
8647c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,m,2);
8747c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,n,3);
8847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,p,4);
8947c6ae99SBarry Smith   dd->m = m;
9047c6ae99SBarry Smith   dd->n = n;
9147c6ae99SBarry Smith   dd->p = p;
9247c6ae99SBarry Smith   PetscFunctionReturn(0);
9347c6ae99SBarry Smith }
9447c6ae99SBarry Smith 
9547c6ae99SBarry Smith #undef __FUNCT__
963e24977eSEthan Coon #define __FUNCT__ "DMDASetBoundaryType"
9747c6ae99SBarry Smith /*@
981321219cSEthan Coon   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries.
9947c6ae99SBarry Smith 
10047c6ae99SBarry Smith   Not collective
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith   Input Parameter:
103aa219208SBarry Smith + da    - The DMDA
1041321219cSEthan Coon - bx,by,bz - One of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC
10547c6ae99SBarry Smith 
10647c6ae99SBarry Smith   Level: intermediate
10747c6ae99SBarry Smith 
10847c6ae99SBarry Smith .keywords:  distributed array, periodicity
109db87c5ecSEthan Coon .seealso: DMDACreate(), DMDestroy(), DMDA, DMDABoundaryType
11047c6ae99SBarry Smith @*/
1111321219cSEthan Coon PetscErrorCode  DMDASetBoundaryType(DM da,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz)
11247c6ae99SBarry Smith {
11347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith   PetscFunctionBegin;
11647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1171321219cSEthan Coon   PetscValidLogicalCollectiveInt(da,bx,2);
1181321219cSEthan Coon   PetscValidLogicalCollectiveInt(da,by,3);
1191321219cSEthan Coon   PetscValidLogicalCollectiveInt(da,bz,4);
1201321219cSEthan Coon   dd->bx = bx;
1211321219cSEthan Coon   dd->by = by;
1221321219cSEthan Coon   dd->bz = bz;
12347c6ae99SBarry Smith   PetscFunctionReturn(0);
12447c6ae99SBarry Smith }
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith #undef __FUNCT__
127aa219208SBarry Smith #define __FUNCT__ "DMDASetDof"
12847c6ae99SBarry Smith /*@
129aa219208SBarry Smith   DMDASetDof - Sets the number of degrees of freedom per vertex
13047c6ae99SBarry Smith 
13147c6ae99SBarry Smith   Not collective
13247c6ae99SBarry Smith 
13347c6ae99SBarry Smith   Input Parameter:
134aa219208SBarry Smith + da  - The DMDA
13547c6ae99SBarry Smith - dof - Number of degrees of freedom
13647c6ae99SBarry Smith 
13747c6ae99SBarry Smith   Level: intermediate
13847c6ae99SBarry Smith 
13947c6ae99SBarry Smith .keywords:  distributed array, degrees of freedom
140aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
14147c6ae99SBarry Smith @*/
1427087cfbeSBarry Smith PetscErrorCode  DMDASetDof(DM da, int dof)
14347c6ae99SBarry Smith {
14447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith   PetscFunctionBegin;
14747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
14847c6ae99SBarry Smith   dd->w = dof;
1491411c6eeSJed Brown   da->bs = dof;
15047c6ae99SBarry Smith   PetscFunctionReturn(0);
15147c6ae99SBarry Smith }
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith #undef __FUNCT__
154aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType"
15547c6ae99SBarry Smith /*@
156aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
15747c6ae99SBarry Smith 
158aa219208SBarry Smith   Logically Collective on DMDA
15947c6ae99SBarry Smith 
16047c6ae99SBarry Smith   Input Parameter:
161aa219208SBarry Smith + da    - The DMDA
162aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
16347c6ae99SBarry Smith 
16447c6ae99SBarry Smith   Level: intermediate
16547c6ae99SBarry Smith 
16647c6ae99SBarry Smith .keywords:  distributed array, stencil
167aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
16847c6ae99SBarry Smith @*/
1697087cfbeSBarry Smith PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
17047c6ae99SBarry Smith {
17147c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
17247c6ae99SBarry Smith 
17347c6ae99SBarry Smith   PetscFunctionBegin;
17447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
17547c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,stype,2);
17647c6ae99SBarry Smith   dd->stencil_type = stype;
17747c6ae99SBarry Smith   PetscFunctionReturn(0);
17847c6ae99SBarry Smith }
17947c6ae99SBarry Smith 
18047c6ae99SBarry Smith #undef __FUNCT__
181aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth"
18247c6ae99SBarry Smith /*@
183aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
18447c6ae99SBarry Smith 
185aa219208SBarry Smith   Logically Collective on DMDA
18647c6ae99SBarry Smith 
18747c6ae99SBarry Smith   Input Parameter:
188aa219208SBarry Smith + da    - The DMDA
18947c6ae99SBarry Smith - width - The stencil width
19047c6ae99SBarry Smith 
19147c6ae99SBarry Smith   Level: intermediate
19247c6ae99SBarry Smith 
19347c6ae99SBarry Smith .keywords:  distributed array, stencil
194aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
19547c6ae99SBarry Smith @*/
1967087cfbeSBarry Smith PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
19747c6ae99SBarry Smith {
19847c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
19947c6ae99SBarry Smith 
20047c6ae99SBarry Smith   PetscFunctionBegin;
20147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
20247c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,width,2);
20347c6ae99SBarry Smith   dd->s = width;
20447c6ae99SBarry Smith   PetscFunctionReturn(0);
20547c6ae99SBarry Smith }
20647c6ae99SBarry Smith 
20747c6ae99SBarry Smith #undef __FUNCT__
208aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
209aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
21047c6ae99SBarry Smith {
21147c6ae99SBarry Smith   PetscInt i,sum;
21247c6ae99SBarry Smith 
21347c6ae99SBarry Smith   PetscFunctionBegin;
21447c6ae99SBarry Smith   if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
21547c6ae99SBarry Smith   for (i=sum=0; i<m; i++) sum += lx[i];
21647c6ae99SBarry Smith   if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
21747c6ae99SBarry Smith   PetscFunctionReturn(0);
21847c6ae99SBarry Smith }
21947c6ae99SBarry Smith 
22047c6ae99SBarry Smith #undef __FUNCT__
221aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges"
22247c6ae99SBarry Smith /*@
223aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
22447c6ae99SBarry Smith 
225aa219208SBarry Smith   Logically Collective on DMDA
22647c6ae99SBarry Smith 
22747c6ae99SBarry Smith   Input Parameter:
228aa219208SBarry Smith + da - The DMDA
22947c6ae99SBarry 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
23047c6ae99SBarry 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
23147c6ae99SBarry 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.
23247c6ae99SBarry Smith 
23347c6ae99SBarry Smith   Level: intermediate
23447c6ae99SBarry Smith 
23547c6ae99SBarry Smith .keywords:  distributed array
236aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
23747c6ae99SBarry Smith @*/
2387087cfbeSBarry Smith PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
23947c6ae99SBarry Smith {
24047c6ae99SBarry Smith   PetscErrorCode ierr;
24147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
24247c6ae99SBarry Smith 
24347c6ae99SBarry Smith   PetscFunctionBegin;
24447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
24547c6ae99SBarry Smith   if (lx) {
24647c6ae99SBarry Smith     if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
247aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
24847c6ae99SBarry Smith     if (!dd->lx) {
24947c6ae99SBarry Smith       ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
25047c6ae99SBarry Smith     }
25147c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
25247c6ae99SBarry Smith   }
25347c6ae99SBarry Smith   if (ly) {
25447c6ae99SBarry Smith     if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
255aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
25647c6ae99SBarry Smith     if (!dd->ly) {
25747c6ae99SBarry Smith       ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
25847c6ae99SBarry Smith     }
25947c6ae99SBarry Smith     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
26047c6ae99SBarry Smith   }
26147c6ae99SBarry Smith   if (lz) {
26247c6ae99SBarry Smith     if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
263aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
26447c6ae99SBarry Smith     if (!dd->lz) {
26547c6ae99SBarry Smith       ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
26647c6ae99SBarry Smith     }
26747c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
26847c6ae99SBarry Smith   }
26947c6ae99SBarry Smith   PetscFunctionReturn(0);
27047c6ae99SBarry Smith }
27147c6ae99SBarry Smith 
27247c6ae99SBarry Smith #undef __FUNCT__
273aa219208SBarry Smith #define __FUNCT__ "DMDACreateOwnershipRanges"
27447c6ae99SBarry Smith /*
275aa219208SBarry Smith  Ensure that da->lx, ly, and lz exist.  Collective on DMDA.
27647c6ae99SBarry Smith */
2777087cfbeSBarry Smith PetscErrorCode  DMDACreateOwnershipRanges(DM da)
27847c6ae99SBarry Smith {
27947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
28047c6ae99SBarry Smith   PetscErrorCode ierr;
28147c6ae99SBarry Smith   PetscInt       n;
28247c6ae99SBarry Smith   MPI_Comm       comm;
28347c6ae99SBarry Smith   PetscMPIInt    size;
28447c6ae99SBarry Smith 
28547c6ae99SBarry Smith   PetscFunctionBegin;
28647c6ae99SBarry Smith   if (!dd->lx) {
28747c6ae99SBarry Smith     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&dd->lx);CHKERRQ(ierr);
288aa219208SBarry Smith     ierr = DMDAGetProcessorSubset(da,DMDA_X,dd->xs,&comm);CHKERRQ(ierr);
28947c6ae99SBarry Smith     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
29047c6ae99SBarry Smith     n = (dd->xe-dd->xs)/dd->w;
29147c6ae99SBarry Smith     ierr = MPI_Allgather(&n,1,MPIU_INT,dd->lx,1,MPIU_INT,comm);CHKERRQ(ierr);
29247c6ae99SBarry Smith   }
29347c6ae99SBarry Smith   if (dd->dim > 1 && !dd->ly) {
29447c6ae99SBarry Smith     ierr = PetscMalloc(dd->n*sizeof(PetscInt),&dd->ly);CHKERRQ(ierr);
295aa219208SBarry Smith     ierr = DMDAGetProcessorSubset(da,DMDA_Y,dd->ys,&comm);CHKERRQ(ierr);
29647c6ae99SBarry Smith     n = dd->ye-dd->ys;
29747c6ae99SBarry Smith     ierr = MPI_Allgather(&n,1,MPIU_INT,dd->ly,1,MPIU_INT,comm);CHKERRQ(ierr);
29847c6ae99SBarry Smith   }
29947c6ae99SBarry Smith   if (dd->dim > 2 && !dd->lz) {
30047c6ae99SBarry Smith     ierr = PetscMalloc(dd->p*sizeof(PetscInt),&dd->lz);CHKERRQ(ierr);
301aa219208SBarry Smith     ierr = DMDAGetProcessorSubset(da,DMDA_Z,dd->zs,&comm);CHKERRQ(ierr);
30247c6ae99SBarry Smith     n = dd->ze-dd->zs;
30347c6ae99SBarry Smith     ierr = MPI_Allgather(&n,1,MPIU_INT,dd->lz,1,MPIU_INT,comm);CHKERRQ(ierr);
30447c6ae99SBarry Smith   }
30547c6ae99SBarry Smith   PetscFunctionReturn(0);
30647c6ae99SBarry Smith }
30747c6ae99SBarry Smith 
30847c6ae99SBarry Smith #undef __FUNCT__
309aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType"
31047c6ae99SBarry Smith /*@
311aa219208SBarry Smith        DMDASetInterpolationType - Sets the type of interpolation that will be
31294013140SBarry Smith           returned by DMGetInterpolation()
31347c6ae99SBarry Smith 
314aa219208SBarry Smith    Logically Collective on DMDA
31547c6ae99SBarry Smith 
31647c6ae99SBarry Smith    Input Parameter:
31747c6ae99SBarry Smith +  da - initial distributed array
318aa219208SBarry Smith .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
31947c6ae99SBarry Smith 
32047c6ae99SBarry Smith    Level: intermediate
32147c6ae99SBarry Smith 
322aa219208SBarry Smith    Notes: you should call this on the coarser of the two DMDAs you pass to DMGetInterpolation()
32347c6ae99SBarry Smith 
32447c6ae99SBarry Smith .keywords:  distributed array, interpolation
32547c6ae99SBarry Smith 
326aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
32747c6ae99SBarry Smith @*/
3287087cfbeSBarry Smith PetscErrorCode  DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
32947c6ae99SBarry Smith {
33047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
33147c6ae99SBarry Smith 
33247c6ae99SBarry Smith   PetscFunctionBegin;
33347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
33447c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,ctype,2);
33547c6ae99SBarry Smith   dd->interptype = ctype;
33647c6ae99SBarry Smith   PetscFunctionReturn(0);
33747c6ae99SBarry Smith }
33847c6ae99SBarry Smith 
33947c6ae99SBarry Smith 
34047c6ae99SBarry Smith #undef __FUNCT__
341aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors"
34247c6ae99SBarry Smith /*@C
343aa219208SBarry Smith       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
34447c6ae99SBarry Smith         processes neighbors.
34547c6ae99SBarry Smith 
34647c6ae99SBarry Smith     Not Collective
34747c6ae99SBarry Smith 
34847c6ae99SBarry Smith    Input Parameter:
349aa219208SBarry Smith .     da - the DMDA object
35047c6ae99SBarry Smith 
35147c6ae99SBarry Smith    Output Parameters:
35247c6ae99SBarry Smith .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
35347c6ae99SBarry Smith               this process itself is in the list
35447c6ae99SBarry Smith 
35547c6ae99SBarry Smith    Notes: In 2d the array is of length 9, in 3d of length 27
35647c6ae99SBarry Smith           Not supported in 1d
357aa219208SBarry Smith           Do not free the array, it is freed when the DMDA is destroyed.
35847c6ae99SBarry Smith 
35947c6ae99SBarry Smith    Fortran Notes: In fortran you must pass in an array of the appropriate length.
36047c6ae99SBarry Smith 
36147c6ae99SBarry Smith    Level: intermediate
36247c6ae99SBarry Smith 
36347c6ae99SBarry Smith @*/
3647087cfbeSBarry Smith PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
36547c6ae99SBarry Smith {
36647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
36747c6ae99SBarry Smith   PetscFunctionBegin;
36847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
36947c6ae99SBarry Smith   *ranks = dd->neighbors;
37047c6ae99SBarry Smith   PetscFunctionReturn(0);
37147c6ae99SBarry Smith }
37247c6ae99SBarry Smith 
37347c6ae99SBarry Smith /*@C
374aa219208SBarry Smith       DMDASetElementType - Sets the element type to be returned by DMGetElements()
37547c6ae99SBarry Smith 
37647c6ae99SBarry Smith     Not Collective
37747c6ae99SBarry Smith 
37847c6ae99SBarry Smith    Input Parameter:
379aa219208SBarry Smith .     da - the DMDA object
38047c6ae99SBarry Smith 
38147c6ae99SBarry Smith    Output Parameters:
382aa219208SBarry Smith .     etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1
38347c6ae99SBarry Smith 
38447c6ae99SBarry Smith    Level: intermediate
38547c6ae99SBarry Smith 
386aa219208SBarry Smith .seealso: DMDAElementType, DMDAGetElementType(), DMGetElements(), DMRestoreElements()
38747c6ae99SBarry Smith @*/
38847c6ae99SBarry Smith #undef __FUNCT__
389aa219208SBarry Smith #define __FUNCT__ "DMDASetElementType"
3907087cfbeSBarry Smith PetscErrorCode  DMDASetElementType(DM da, DMDAElementType etype)
39147c6ae99SBarry Smith {
39247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
39347c6ae99SBarry Smith   PetscErrorCode ierr;
39447c6ae99SBarry Smith 
39547c6ae99SBarry Smith   PetscFunctionBegin;
39647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
39747c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,etype,2);
39847c6ae99SBarry Smith   if (dd->elementtype != etype) {
39947c6ae99SBarry Smith     ierr = PetscFree(dd->e);CHKERRQ(ierr);
40047c6ae99SBarry Smith     dd->elementtype = etype;
40147c6ae99SBarry Smith     dd->ne          = 0;
40247c6ae99SBarry Smith     dd->e           = PETSC_NULL;
40347c6ae99SBarry Smith   }
40447c6ae99SBarry Smith   PetscFunctionReturn(0);
40547c6ae99SBarry Smith }
40647c6ae99SBarry Smith 
40747c6ae99SBarry Smith /*@C
408aa219208SBarry Smith       DMDAGetElementType - Gets the element type to be returned by DMGetElements()
40947c6ae99SBarry Smith 
41047c6ae99SBarry Smith     Not Collective
41147c6ae99SBarry Smith 
41247c6ae99SBarry Smith    Input Parameter:
413aa219208SBarry Smith .     da - the DMDA object
41447c6ae99SBarry Smith 
41547c6ae99SBarry Smith    Output Parameters:
416aa219208SBarry Smith .     etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1
41747c6ae99SBarry Smith 
41847c6ae99SBarry Smith    Level: intermediate
41947c6ae99SBarry Smith 
420aa219208SBarry Smith .seealso: DMDAElementType, DMDASetElementType(), DMGetElements(), DMRestoreElements()
42147c6ae99SBarry Smith @*/
42247c6ae99SBarry Smith #undef __FUNCT__
423aa219208SBarry Smith #define __FUNCT__ "DMDAGetElementType"
4247087cfbeSBarry Smith PetscErrorCode  DMDAGetElementType(DM da, DMDAElementType *etype)
42547c6ae99SBarry Smith {
42647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
42747c6ae99SBarry Smith   PetscFunctionBegin;
42847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
42947c6ae99SBarry Smith   PetscValidPointer(etype,2);
43047c6ae99SBarry Smith   *etype = dd->elementtype;
43147c6ae99SBarry Smith   PetscFunctionReturn(0);
43247c6ae99SBarry Smith }
43347c6ae99SBarry Smith 
43447c6ae99SBarry Smith #undef __FUNCT__
43547c6ae99SBarry Smith #define __FUNCT__ "DMGetElements"
43647c6ae99SBarry Smith /*@C
43747c6ae99SBarry Smith       DMGetElements - Gets an array containing the indices (in local coordinates)
43847c6ae99SBarry Smith                  of all the local elements
43947c6ae99SBarry Smith 
44047c6ae99SBarry Smith     Not Collective
44147c6ae99SBarry Smith 
44247c6ae99SBarry Smith    Input Parameter:
44347c6ae99SBarry Smith .     dm - the DM object
44447c6ae99SBarry Smith 
44547c6ae99SBarry Smith    Output Parameters:
446454e267fSLisandro Dalcin +     nel - number of local elements
447454e267fSLisandro Dalcin .     nen - number of element nodes
44847c6ae99SBarry Smith -     e - the indices of the elements vertices
44947c6ae99SBarry Smith 
45047c6ae99SBarry Smith    Level: intermediate
45147c6ae99SBarry Smith 
45247c6ae99SBarry Smith .seealso: DMElementType, DMSetElementType(), DMRestoreElements()
45347c6ae99SBarry Smith @*/
4547087cfbeSBarry Smith PetscErrorCode  DMGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
45547c6ae99SBarry Smith {
45647c6ae99SBarry Smith   PetscErrorCode ierr;
45747c6ae99SBarry Smith   PetscFunctionBegin;
45847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
459454e267fSLisandro Dalcin   PetscValidIntPointer(nel,2);
460454e267fSLisandro Dalcin   PetscValidIntPointer(nen,3);
461454e267fSLisandro Dalcin   PetscValidPointer(e,4);
462454e267fSLisandro Dalcin   ierr = (dm->ops->getelements)(dm,nel,nen,e);CHKERRQ(ierr);
46347c6ae99SBarry Smith   PetscFunctionReturn(0);
46447c6ae99SBarry Smith }
46547c6ae99SBarry Smith 
46647c6ae99SBarry Smith #undef __FUNCT__
46747c6ae99SBarry Smith #define __FUNCT__ "DMRestoreElements"
46847c6ae99SBarry Smith /*@C
46947c6ae99SBarry Smith       DMRestoreElements - Returns an array containing the indices (in local coordinates)
47047c6ae99SBarry Smith                  of all the local elements obtained with DMGetElements()
47147c6ae99SBarry Smith 
47247c6ae99SBarry Smith     Not Collective
47347c6ae99SBarry Smith 
47447c6ae99SBarry Smith    Input Parameter:
47547c6ae99SBarry Smith +     dm - the DM object
476454e267fSLisandro Dalcin .     nel - number of local elements
477454e267fSLisandro Dalcin .     nen - number of element nodes
47847c6ae99SBarry Smith -     e - the indices of the elements vertices
47947c6ae99SBarry Smith 
48047c6ae99SBarry Smith    Level: intermediate
48147c6ae99SBarry Smith 
48247c6ae99SBarry Smith .seealso: DMElementType, DMSetElementType(), DMGetElements()
48347c6ae99SBarry Smith @*/
4847087cfbeSBarry Smith PetscErrorCode  DMRestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
48547c6ae99SBarry Smith {
48647c6ae99SBarry Smith   PetscErrorCode ierr;
48747c6ae99SBarry Smith   PetscFunctionBegin;
48847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
489454e267fSLisandro Dalcin   PetscValidIntPointer(nel,2);
490454e267fSLisandro Dalcin   PetscValidIntPointer(nen,3);
491454e267fSLisandro Dalcin   PetscValidPointer(e,4);
49247c6ae99SBarry Smith   if (dm->ops->restoreelements) {
493454e267fSLisandro Dalcin     ierr = (dm->ops->restoreelements)(dm,nel,nen,e);CHKERRQ(ierr);
49447c6ae99SBarry Smith   }
49547c6ae99SBarry Smith   PetscFunctionReturn(0);
49647c6ae99SBarry Smith }
49747c6ae99SBarry Smith 
49847c6ae99SBarry Smith #undef __FUNCT__
499aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges"
50047c6ae99SBarry Smith /*@C
501aa219208SBarry Smith       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
50247c6ae99SBarry Smith 
50347c6ae99SBarry Smith     Not Collective
50447c6ae99SBarry Smith 
50547c6ae99SBarry Smith    Input Parameter:
506aa219208SBarry Smith .     da - the DMDA object
50747c6ae99SBarry Smith 
50847c6ae99SBarry Smith    Output Parameter:
50947c6ae99SBarry Smith +     lx - ownership along x direction (optional)
51047c6ae99SBarry Smith .     ly - ownership along y direction (optional)
51147c6ae99SBarry Smith -     lz - ownership along z direction (optional)
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith    Level: intermediate
51447c6ae99SBarry Smith 
515aa219208SBarry Smith     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
51647c6ae99SBarry Smith 
51747c6ae99SBarry Smith     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
518aa219208SBarry Smith     eighth arguments from DMDAGetInfo()
51947c6ae99SBarry Smith 
52047c6ae99SBarry Smith      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
521aa219208SBarry Smith     DMDA they came from still exists (has not been destroyed).
52247c6ae99SBarry Smith 
523aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
52447c6ae99SBarry Smith @*/
5257087cfbeSBarry Smith PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
52647c6ae99SBarry Smith {
52747c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
52847c6ae99SBarry Smith 
52947c6ae99SBarry Smith   PetscFunctionBegin;
53047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
53147c6ae99SBarry Smith   if (lx) *lx = dd->lx;
53247c6ae99SBarry Smith   if (ly) *ly = dd->ly;
53347c6ae99SBarry Smith   if (lz) *lz = dd->lz;
53447c6ae99SBarry Smith   PetscFunctionReturn(0);
53547c6ae99SBarry Smith }
53647c6ae99SBarry Smith 
53747c6ae99SBarry Smith #undef __FUNCT__
538aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor"
53947c6ae99SBarry Smith /*@
540aa219208SBarry Smith      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
54147c6ae99SBarry Smith 
542aa219208SBarry Smith     Logically Collective on DMDA
54347c6ae99SBarry Smith 
54447c6ae99SBarry Smith   Input Parameters:
545aa219208SBarry Smith +    da - the DMDA object
54647c6ae99SBarry Smith .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
54747c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
54847c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
54947c6ae99SBarry Smith 
55047c6ae99SBarry Smith   Options Database:
55147c6ae99SBarry Smith +  -da_refine_x - refinement ratio in x direction
55247c6ae99SBarry Smith .  -da_refine_y - refinement ratio in y direction
55347c6ae99SBarry Smith -  -da_refine_z - refinement ratio in z direction
55447c6ae99SBarry Smith 
55547c6ae99SBarry Smith   Level: intermediate
55647c6ae99SBarry Smith 
55747c6ae99SBarry Smith     Notes: Pass PETSC_IGNORE to leave a value unchanged
55847c6ae99SBarry Smith 
559aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor()
56047c6ae99SBarry Smith @*/
5617087cfbeSBarry Smith PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
56247c6ae99SBarry Smith {
56347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
56447c6ae99SBarry Smith 
56547c6ae99SBarry Smith   PetscFunctionBegin;
56647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
56747c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_x,2);
56847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_y,3);
56947c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_z,4);
57047c6ae99SBarry Smith 
57147c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
57247c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
57347c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
57447c6ae99SBarry Smith   PetscFunctionReturn(0);
57547c6ae99SBarry Smith }
57647c6ae99SBarry Smith 
57747c6ae99SBarry Smith #undef __FUNCT__
578aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor"
57947c6ae99SBarry Smith /*@C
580aa219208SBarry Smith      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
58147c6ae99SBarry Smith 
58247c6ae99SBarry Smith     Not Collective
58347c6ae99SBarry Smith 
58447c6ae99SBarry Smith   Input Parameter:
585aa219208SBarry Smith .    da - the DMDA object
58647c6ae99SBarry Smith 
58747c6ae99SBarry Smith   Output Parameters:
58847c6ae99SBarry Smith +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
58947c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
59047c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
59147c6ae99SBarry Smith 
59247c6ae99SBarry Smith   Level: intermediate
59347c6ae99SBarry Smith 
59447c6ae99SBarry Smith     Notes: Pass PETSC_NULL for values you do not need
59547c6ae99SBarry Smith 
596aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor()
59747c6ae99SBarry Smith @*/
5987087cfbeSBarry Smith PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
59947c6ae99SBarry Smith {
60047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
60147c6ae99SBarry Smith 
60247c6ae99SBarry Smith   PetscFunctionBegin;
60347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
60447c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
60547c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
60647c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
60747c6ae99SBarry Smith   PetscFunctionReturn(0);
60847c6ae99SBarry Smith }
60947c6ae99SBarry Smith 
61047c6ae99SBarry Smith #undef __FUNCT__
611aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix"
61247c6ae99SBarry Smith /*@C
613aa219208SBarry Smith      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
61447c6ae99SBarry Smith 
615aa219208SBarry Smith     Logically Collective on DMDA
61647c6ae99SBarry Smith 
61747c6ae99SBarry Smith   Input Parameters:
618aa219208SBarry Smith +    da - the DMDA object
619aa219208SBarry Smith -    f - the function that allocates the matrix for that specific DMDA
62047c6ae99SBarry Smith 
62147c6ae99SBarry Smith   Level: developer
62247c6ae99SBarry Smith 
623aa219208SBarry Smith    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
62447c6ae99SBarry Smith        the diagonal and off-diagonal blocks of the matrix
62547c6ae99SBarry Smith 
626aa219208SBarry Smith .seealso: DMGetMatrix(), DMDASetBlockFills()
62747c6ae99SBarry Smith @*/
6287087cfbeSBarry Smith PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, const MatType,Mat*))
62947c6ae99SBarry Smith {
63047c6ae99SBarry Smith   PetscFunctionBegin;
63147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
63247c6ae99SBarry Smith   da->ops->getmatrix = f;
63347c6ae99SBarry Smith   PetscFunctionReturn(0);
63447c6ae99SBarry Smith }
63547c6ae99SBarry Smith 
63647c6ae99SBarry Smith #undef __FUNCT__
63794013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges"
63847c6ae99SBarry Smith /*
63947c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
64047c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
64147c6ae99SBarry Smith 
64247c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
64347c6ae99SBarry Smith */
64494013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
64547c6ae99SBarry Smith {
64647c6ae99SBarry Smith   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
64747c6ae99SBarry Smith   PetscErrorCode ierr;
64847c6ae99SBarry Smith 
64947c6ae99SBarry Smith   PetscFunctionBegin;
65047c6ae99SBarry Smith   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
65147c6ae99SBarry Smith   if (ratio == 1) {
65247c6ae99SBarry Smith     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
65347c6ae99SBarry Smith     PetscFunctionReturn(0);
65447c6ae99SBarry Smith   }
65547c6ae99SBarry Smith   for (i=0; i<m; i++) totalc += lc[i];
65647c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
65747c6ae99SBarry Smith   for (i=0; i<m; i++) {
65847c6ae99SBarry Smith     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
65947c6ae99SBarry Smith     if (i == m-1) lf[i] = want;
66047c6ae99SBarry Smith     else {
6617aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
6627aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
6637aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
6647aca7175SJed Brown       while ((startf+want)/ratio < nextc - stencil_width) want++;
6657aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
6667aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
6677aca7175SJed Brown       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
6687aca7175SJed Brown       /* Make sure all constraints are satisfied */
6697aca7175SJed Brown       if (want < 0 || want > remaining
6707aca7175SJed Brown           || ((startf+want)/ratio < nextc - stencil_width)
6717aca7175SJed Brown           || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width))
6727aca7175SJed Brown         SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range");
67347c6ae99SBarry Smith     }
67447c6ae99SBarry Smith     lf[i] = want;
67547c6ae99SBarry Smith     startc += lc[i];
67647c6ae99SBarry Smith     startf += lf[i];
67747c6ae99SBarry Smith     remaining -= lf[i];
67847c6ae99SBarry Smith   }
67947c6ae99SBarry Smith   PetscFunctionReturn(0);
68047c6ae99SBarry Smith }
68147c6ae99SBarry Smith 
68247c6ae99SBarry Smith #undef __FUNCT__
6832be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges"
6842be375d4SJed Brown /*
6852be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
6862be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
6872be375d4SJed Brown 
6882be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
6892be375d4SJed Brown */
6902be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
6912be375d4SJed Brown {
6922be375d4SJed Brown   PetscInt       i,totalf,remaining,startc,startf;
6932be375d4SJed Brown   PetscErrorCode ierr;
6942be375d4SJed Brown 
6952be375d4SJed Brown   PetscFunctionBegin;
6962be375d4SJed Brown   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
6972be375d4SJed Brown   if (ratio == 1) {
6982be375d4SJed Brown     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
6992be375d4SJed Brown     PetscFunctionReturn(0);
7002be375d4SJed Brown   }
7012be375d4SJed Brown   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
7022be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
7032be375d4SJed Brown   for (i=0,startc=0,startf=0; i<m; i++) {
7042be375d4SJed Brown     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
7052be375d4SJed Brown     if (i == m-1) lc[i] = want;
7062be375d4SJed Brown     else {
7072be375d4SJed Brown       const PetscInt nextf = startf+lf[i];
7082be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
7092be375d4SJed Brown        * node is within one stencil width. */
7102be375d4SJed Brown       while (nextf/ratio < startc+want-stencil_width) want--;
7112be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
7122be375d4SJed Brown        * fine node is within one stencil width. */
7132be375d4SJed Brown       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
7142be375d4SJed Brown       if (want < 0 || want > remaining
7152be375d4SJed Brown           || (nextf/ratio < startc+want-stencil_width)
7162be375d4SJed Brown           || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width))
7172be375d4SJed Brown         SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range");
7182be375d4SJed Brown     }
7192be375d4SJed Brown     lc[i] = want;
7202be375d4SJed Brown     startc += lc[i];
7212be375d4SJed Brown     startf += lf[i];
7222be375d4SJed Brown     remaining -= lc[i];
7232be375d4SJed Brown   }
7242be375d4SJed Brown   PetscFunctionReturn(0);
7252be375d4SJed Brown }
7262be375d4SJed Brown 
7272be375d4SJed Brown #undef __FUNCT__
72894013140SBarry Smith #define __FUNCT__ "DMRefine_DA"
7297087cfbeSBarry Smith PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
73047c6ae99SBarry Smith {
73147c6ae99SBarry Smith   PetscErrorCode ierr;
73247c6ae99SBarry Smith   PetscInt       M,N,P;
7339a42bb27SBarry Smith   DM             da2;
73447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
73547c6ae99SBarry Smith 
73647c6ae99SBarry Smith   PetscFunctionBegin;
73747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
73847c6ae99SBarry Smith   PetscValidPointer(daref,3);
73947c6ae99SBarry Smith 
7401321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
74147c6ae99SBarry Smith     M = dd->refine_x*dd->M;
74247c6ae99SBarry Smith   } else {
74347c6ae99SBarry Smith     M = 1 + dd->refine_x*(dd->M - 1);
74447c6ae99SBarry Smith   }
7451321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
74647c6ae99SBarry Smith     N = dd->refine_y*dd->N;
74747c6ae99SBarry Smith   } else {
74847c6ae99SBarry Smith     N = 1 + dd->refine_y*(dd->N - 1);
74947c6ae99SBarry Smith   }
7501321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
75147c6ae99SBarry Smith     P = dd->refine_z*dd->P;
75247c6ae99SBarry Smith   } else {
75347c6ae99SBarry Smith     P = 1 + dd->refine_z*(dd->P - 1);
75447c6ae99SBarry Smith   }
755aa219208SBarry Smith   ierr = DMDACreateOwnershipRanges(da);CHKERRQ(ierr);
75647c6ae99SBarry Smith   if (dd->dim == 3) {
75747c6ae99SBarry Smith     PetscInt *lx,*ly,*lz;
75847c6ae99SBarry Smith     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
7591321219cSEthan 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);
7601321219cSEthan 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);
7611321219cSEthan 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);
7621321219cSEthan 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);
76347c6ae99SBarry Smith     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
76447c6ae99SBarry Smith   } else if (dd->dim == 2) {
76547c6ae99SBarry Smith     PetscInt *lx,*ly;
76647c6ae99SBarry Smith     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
7671321219cSEthan 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);
7681321219cSEthan 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);
7691321219cSEthan 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);
77047c6ae99SBarry Smith     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
77147c6ae99SBarry Smith   } else if (dd->dim == 1) {
77247c6ae99SBarry Smith     PetscInt *lx;
77347c6ae99SBarry Smith     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
7741321219cSEthan 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);
7751321219cSEthan Coon     ierr = DMDACreate1d(((PetscObject)da)->comm,dd->bx,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr);
77647c6ae99SBarry Smith     ierr = PetscFree(lx);CHKERRQ(ierr);
77747c6ae99SBarry Smith   }
77847c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
77947c6ae99SBarry Smith 
78047c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
78147c6ae99SBarry Smith   da2->ops->getmatrix        = da->ops->getmatrix;
78247c6ae99SBarry Smith   da2->ops->getinterpolation = da->ops->getinterpolation;
78347c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
78447c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
78547c6ae99SBarry Smith 
78647c6ae99SBarry Smith   /* copy fill information if given */
78747c6ae99SBarry Smith   if (dd->dfill) {
78847c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
78947c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
79047c6ae99SBarry Smith   }
79147c6ae99SBarry Smith   if (dd->ofill) {
79247c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
79347c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
79447c6ae99SBarry Smith   }
79547c6ae99SBarry Smith   /* copy the refine information */
79647c6ae99SBarry Smith   dd2->refine_x = dd->refine_x;
79747c6ae99SBarry Smith   dd2->refine_y = dd->refine_y;
79847c6ae99SBarry Smith   dd2->refine_z = dd->refine_z;
79947c6ae99SBarry Smith 
80047c6ae99SBarry Smith   /* copy vector type information */
80147c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
80247c6ae99SBarry Smith   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
803ddcf8b74SDave May 
804ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
805ddcf8b74SDave May   if (dd->coordinates) {
806ddcf8b74SDave May     DM  cdaf,cdac;
807ddcf8b74SDave May     Vec coordsc,coordsf;
808ddcf8b74SDave May     Mat II;
809ddcf8b74SDave May 
810ddcf8b74SDave May     ierr = DMDAGetCoordinateDA(da,&cdac);CHKERRQ(ierr);
811ddcf8b74SDave May     ierr = DMDAGetCoordinates(da,&coordsc);CHKERRQ(ierr);
812b61d3410SDave May     ierr = DMDAGetCoordinateDA(da2,&cdaf);CHKERRQ(ierr);
813b61d3410SDave May     /* force creation of the coordinate vector */
814b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
815ddcf8b74SDave May     ierr = DMDAGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
816ddcf8b74SDave May     ierr = DMGetInterpolation(cdac,cdaf,&II,PETSC_NULL);CHKERRQ(ierr);
817ddcf8b74SDave May     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
818ddcf8b74SDave May     ierr = MatDestroy(II);CHKERRQ(ierr);
819ddcf8b74SDave May   }
82047c6ae99SBarry Smith   *daref = da2;
82147c6ae99SBarry Smith   PetscFunctionReturn(0);
82247c6ae99SBarry Smith }
82347c6ae99SBarry Smith 
82447c6ae99SBarry Smith #undef __FUNCT__
82594013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA"
8267087cfbeSBarry Smith PetscErrorCode  DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
82747c6ae99SBarry Smith {
82847c6ae99SBarry Smith   PetscErrorCode ierr;
82947c6ae99SBarry Smith   PetscInt       M,N,P;
8309a42bb27SBarry Smith   DM             da2;
83147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
83247c6ae99SBarry Smith 
83347c6ae99SBarry Smith   PetscFunctionBegin;
83447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
83547c6ae99SBarry Smith   PetscValidPointer(daref,3);
83647c6ae99SBarry Smith 
8371321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
83847c6ae99SBarry Smith     M = dd->M / dd->refine_x;
83947c6ae99SBarry Smith   } else {
84047c6ae99SBarry Smith     M = 1 + (dd->M - 1) / dd->refine_x;
84147c6ae99SBarry Smith   }
8421321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
84347c6ae99SBarry Smith     N = dd->N / dd->refine_y;
84447c6ae99SBarry Smith   } else {
84547c6ae99SBarry Smith     N = 1 + (dd->N - 1) / dd->refine_y;
84647c6ae99SBarry Smith   }
8471321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
84847c6ae99SBarry Smith     P = dd->P / dd->refine_z;
84947c6ae99SBarry Smith   } else {
85047c6ae99SBarry Smith     P = 1 + (dd->P - 1) / dd->refine_z;
85147c6ae99SBarry Smith   }
8522be375d4SJed Brown   ierr = DMDACreateOwnershipRanges(da);CHKERRQ(ierr);
85347c6ae99SBarry Smith   if (dd->dim == 3) {
8542be375d4SJed Brown     PetscInt *lx,*ly,*lz;
8552be375d4SJed Brown     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
8561321219cSEthan 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);
8571321219cSEthan 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);
8581321219cSEthan 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);
8591321219cSEthan 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);
8602be375d4SJed Brown     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
86147c6ae99SBarry Smith   } else if (dd->dim == 2) {
8622be375d4SJed Brown     PetscInt *lx,*ly;
8632be375d4SJed Brown     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
8641321219cSEthan 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);
8651321219cSEthan 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);
8661321219cSEthan 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);
8672be375d4SJed Brown     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
86847c6ae99SBarry Smith   } else if (dd->dim == 1) {
8692be375d4SJed Brown     PetscInt *lx;
8702be375d4SJed Brown     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
8711321219cSEthan 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);
8721321219cSEthan Coon     ierr = DMDACreate1d(((PetscObject)da)->comm,dd->bx,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr);
8732be375d4SJed Brown     ierr = PetscFree(lx);CHKERRQ(ierr);
87447c6ae99SBarry Smith   }
87547c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
87647c6ae99SBarry Smith 
87747c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
87847c6ae99SBarry Smith   da2->ops->getmatrix        = da->ops->getmatrix;
87947c6ae99SBarry Smith   da2->ops->getinterpolation = da->ops->getinterpolation;
88047c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
88147c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
88247c6ae99SBarry Smith 
88347c6ae99SBarry Smith   /* copy fill information if given */
88447c6ae99SBarry Smith   if (dd->dfill) {
88547c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
88647c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
88747c6ae99SBarry Smith   }
88847c6ae99SBarry Smith   if (dd->ofill) {
88947c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
89047c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
89147c6ae99SBarry Smith   }
89247c6ae99SBarry Smith   /* copy the refine information */
89347c6ae99SBarry Smith   dd2->refine_x = dd->refine_x;
89447c6ae99SBarry Smith   dd2->refine_y = dd->refine_y;
89547c6ae99SBarry Smith   dd2->refine_z = dd->refine_z;
89647c6ae99SBarry Smith 
89747c6ae99SBarry Smith   /* copy vector type information */
89847c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
89947c6ae99SBarry Smith   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
90047c6ae99SBarry Smith 
901ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
902ddcf8b74SDave May   if (dd->coordinates) {
903ddcf8b74SDave May     DM  cdaf,cdac;
904ddcf8b74SDave May     Vec coordsc,coordsf;
905ddcf8b74SDave May     VecScatter inject;
906ddcf8b74SDave May 
907ddcf8b74SDave May     ierr = DMDAGetCoordinateDA(da,&cdaf);CHKERRQ(ierr);
908ddcf8b74SDave May     ierr = DMDAGetCoordinates(da,&coordsf);CHKERRQ(ierr);
909b61d3410SDave May     ierr = DMDAGetCoordinateDA(da2,&cdac);CHKERRQ(ierr);
910b61d3410SDave May     /* force creation of the coordinate vector */
911b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
912ddcf8b74SDave May     ierr = DMDAGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
913ddcf8b74SDave May 
914ddcf8b74SDave May     ierr = DMGetInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
915ddcf8b74SDave May     ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
916ddcf8b74SDave May     ierr = VecScatterEnd(inject  ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
917ddcf8b74SDave May     ierr = VecScatterDestroy(inject);CHKERRQ(ierr);
918ddcf8b74SDave May   }
91947c6ae99SBarry Smith   *daref = da2;
92047c6ae99SBarry Smith   PetscFunctionReturn(0);
92147c6ae99SBarry Smith }
92247c6ae99SBarry Smith 
92347c6ae99SBarry Smith #undef __FUNCT__
92494013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA"
9257087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
92647c6ae99SBarry Smith {
92747c6ae99SBarry Smith   PetscErrorCode ierr;
92847c6ae99SBarry Smith   PetscInt       i,n,*refx,*refy,*refz;
92947c6ae99SBarry Smith 
93047c6ae99SBarry Smith   PetscFunctionBegin;
93147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
93247c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
93347c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
93447c6ae99SBarry Smith   PetscValidPointer(daf,3);
93547c6ae99SBarry Smith 
936aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
93747c6ae99SBarry Smith   ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr);
93847c6ae99SBarry Smith   for (i=0; i<nlevels; i++) {
939aa219208SBarry Smith     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
94047c6ae99SBarry Smith   }
94147c6ae99SBarry Smith   n = nlevels;
94247c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
94347c6ae99SBarry Smith   n = nlevels;
94447c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
94547c6ae99SBarry Smith   n = nlevels;
94647c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
94747c6ae99SBarry Smith 
948aa219208SBarry Smith   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
94994013140SBarry Smith   ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr);
95047c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
951aa219208SBarry Smith     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
95294013140SBarry Smith     ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr);
95347c6ae99SBarry Smith   }
95447c6ae99SBarry Smith   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
95547c6ae99SBarry Smith   PetscFunctionReturn(0);
95647c6ae99SBarry Smith }
95747c6ae99SBarry Smith 
95847c6ae99SBarry Smith #undef __FUNCT__
95994013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA"
9607087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
96147c6ae99SBarry Smith {
96247c6ae99SBarry Smith   PetscErrorCode ierr;
96347c6ae99SBarry Smith   PetscInt       i;
96447c6ae99SBarry Smith 
96547c6ae99SBarry Smith   PetscFunctionBegin;
96647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
96747c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
96847c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
96947c6ae99SBarry Smith   PetscValidPointer(dac,3);
97094013140SBarry Smith   ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr);
97147c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
97294013140SBarry Smith     ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr);
97347c6ae99SBarry Smith   }
97447c6ae99SBarry Smith   PetscFunctionReturn(0);
97547c6ae99SBarry Smith }
976