xref: /petsc/src/dm/impls/da/da.c (revision ddcf8b7473737f7bd8efb4ed51fab5dc4c8fc0b1)
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 @*/
20aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT 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 @*/
48aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT 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 @*/
81aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT 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__
97aa219208SBarry Smith #define __FUNCT__ "DMDASetPeriodicity"
9847c6ae99SBarry Smith /*@
99aa219208SBarry Smith   DMDASetPeriodicity - Sets the type of periodicity
10047c6ae99SBarry Smith 
10147c6ae99SBarry Smith   Not collective
10247c6ae99SBarry Smith 
10347c6ae99SBarry Smith   Input Parameter:
104aa219208SBarry Smith + da    - The DMDA
105aa219208SBarry Smith - ptype - One of DMDA_NONPERIODIC, DMDA_XPERIODIC, DMDA_YPERIODIC, DMDA_ZPERIODIC, DMDA_XYPERIODIC, DMDA_XZPERIODIC, DMDA_YZPERIODIC, or DMDA_XYZPERIODIC
10647c6ae99SBarry Smith 
10747c6ae99SBarry Smith   Level: intermediate
10847c6ae99SBarry Smith 
10947c6ae99SBarry Smith .keywords:  distributed array, periodicity
110aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA, DMDAPeriodicType
11147c6ae99SBarry Smith @*/
112aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetPeriodicity(DM da, DMDAPeriodicType ptype)
11347c6ae99SBarry Smith {
11447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
11547c6ae99SBarry Smith 
11647c6ae99SBarry Smith   PetscFunctionBegin;
11747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
11847c6ae99SBarry Smith   dd->wrap = ptype;
11947c6ae99SBarry Smith   PetscFunctionReturn(0);
12047c6ae99SBarry Smith }
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith #undef __FUNCT__
123aa219208SBarry Smith #define __FUNCT__ "DMDASetDof"
12447c6ae99SBarry Smith /*@
125aa219208SBarry Smith   DMDASetDof - Sets the number of degrees of freedom per vertex
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith   Not collective
12847c6ae99SBarry Smith 
12947c6ae99SBarry Smith   Input Parameter:
130aa219208SBarry Smith + da  - The DMDA
13147c6ae99SBarry Smith - dof - Number of degrees of freedom
13247c6ae99SBarry Smith 
13347c6ae99SBarry Smith   Level: intermediate
13447c6ae99SBarry Smith 
13547c6ae99SBarry Smith .keywords:  distributed array, degrees of freedom
136aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
13747c6ae99SBarry Smith @*/
138aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetDof(DM da, int dof)
13947c6ae99SBarry Smith {
14047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
14147c6ae99SBarry Smith 
14247c6ae99SBarry Smith   PetscFunctionBegin;
14347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
14447c6ae99SBarry Smith   dd->w = dof;
1451411c6eeSJed Brown   da->bs = dof;
14647c6ae99SBarry Smith   PetscFunctionReturn(0);
14747c6ae99SBarry Smith }
14847c6ae99SBarry Smith 
14947c6ae99SBarry Smith #undef __FUNCT__
150aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType"
15147c6ae99SBarry Smith /*@
152aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
15347c6ae99SBarry Smith 
154aa219208SBarry Smith   Logically Collective on DMDA
15547c6ae99SBarry Smith 
15647c6ae99SBarry Smith   Input Parameter:
157aa219208SBarry Smith + da    - The DMDA
158aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
15947c6ae99SBarry Smith 
16047c6ae99SBarry Smith   Level: intermediate
16147c6ae99SBarry Smith 
16247c6ae99SBarry Smith .keywords:  distributed array, stencil
163aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
16447c6ae99SBarry Smith @*/
165aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetStencilType(DM da, DMDAStencilType stype)
16647c6ae99SBarry Smith {
16747c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
16847c6ae99SBarry Smith 
16947c6ae99SBarry Smith   PetscFunctionBegin;
17047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
17147c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,stype,2);
17247c6ae99SBarry Smith   dd->stencil_type = stype;
17347c6ae99SBarry Smith   PetscFunctionReturn(0);
17447c6ae99SBarry Smith }
17547c6ae99SBarry Smith 
17647c6ae99SBarry Smith #undef __FUNCT__
177aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth"
17847c6ae99SBarry Smith /*@
179aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
18047c6ae99SBarry Smith 
181aa219208SBarry Smith   Logically Collective on DMDA
18247c6ae99SBarry Smith 
18347c6ae99SBarry Smith   Input Parameter:
184aa219208SBarry Smith + da    - The DMDA
18547c6ae99SBarry Smith - width - The stencil width
18647c6ae99SBarry Smith 
18747c6ae99SBarry Smith   Level: intermediate
18847c6ae99SBarry Smith 
18947c6ae99SBarry Smith .keywords:  distributed array, stencil
190aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
19147c6ae99SBarry Smith @*/
192aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetStencilWidth(DM da, PetscInt width)
19347c6ae99SBarry Smith {
19447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
19547c6ae99SBarry Smith 
19647c6ae99SBarry Smith   PetscFunctionBegin;
19747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
19847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,width,2);
19947c6ae99SBarry Smith   dd->s = width;
20047c6ae99SBarry Smith   PetscFunctionReturn(0);
20147c6ae99SBarry Smith }
20247c6ae99SBarry Smith 
20347c6ae99SBarry Smith #undef __FUNCT__
204aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
205aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
20647c6ae99SBarry Smith {
20747c6ae99SBarry Smith   PetscInt i,sum;
20847c6ae99SBarry Smith 
20947c6ae99SBarry Smith   PetscFunctionBegin;
21047c6ae99SBarry Smith   if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
21147c6ae99SBarry Smith   for (i=sum=0; i<m; i++) sum += lx[i];
21247c6ae99SBarry Smith   if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
21347c6ae99SBarry Smith   PetscFunctionReturn(0);
21447c6ae99SBarry Smith }
21547c6ae99SBarry Smith 
21647c6ae99SBarry Smith #undef __FUNCT__
217aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges"
21847c6ae99SBarry Smith /*@
219aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
22047c6ae99SBarry Smith 
221aa219208SBarry Smith   Logically Collective on DMDA
22247c6ae99SBarry Smith 
22347c6ae99SBarry Smith   Input Parameter:
224aa219208SBarry Smith + da - The DMDA
22547c6ae99SBarry 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
22647c6ae99SBarry 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
22747c6ae99SBarry 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.
22847c6ae99SBarry Smith 
22947c6ae99SBarry Smith   Level: intermediate
23047c6ae99SBarry Smith 
23147c6ae99SBarry Smith .keywords:  distributed array
232aa219208SBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
23347c6ae99SBarry Smith @*/
234aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
23547c6ae99SBarry Smith {
23647c6ae99SBarry Smith   PetscErrorCode ierr;
23747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
23847c6ae99SBarry Smith 
23947c6ae99SBarry Smith   PetscFunctionBegin;
24047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
24147c6ae99SBarry Smith   if (lx) {
24247c6ae99SBarry Smith     if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
243aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
24447c6ae99SBarry Smith     if (!dd->lx) {
24547c6ae99SBarry Smith       ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
24647c6ae99SBarry Smith     }
24747c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
24847c6ae99SBarry Smith   }
24947c6ae99SBarry Smith   if (ly) {
25047c6ae99SBarry Smith     if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
251aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
25247c6ae99SBarry Smith     if (!dd->ly) {
25347c6ae99SBarry Smith       ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
25447c6ae99SBarry Smith     }
25547c6ae99SBarry Smith     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
25647c6ae99SBarry Smith   }
25747c6ae99SBarry Smith   if (lz) {
25847c6ae99SBarry Smith     if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
259aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
26047c6ae99SBarry Smith     if (!dd->lz) {
26147c6ae99SBarry Smith       ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
26247c6ae99SBarry Smith     }
26347c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
26447c6ae99SBarry Smith   }
26547c6ae99SBarry Smith   PetscFunctionReturn(0);
26647c6ae99SBarry Smith }
26747c6ae99SBarry Smith 
26847c6ae99SBarry Smith #undef __FUNCT__
269aa219208SBarry Smith #define __FUNCT__ "DMDACreateOwnershipRanges"
27047c6ae99SBarry Smith /*
271aa219208SBarry Smith  Ensure that da->lx, ly, and lz exist.  Collective on DMDA.
27247c6ae99SBarry Smith */
273aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDACreateOwnershipRanges(DM da)
27447c6ae99SBarry Smith {
27547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
27647c6ae99SBarry Smith   PetscErrorCode ierr;
27747c6ae99SBarry Smith   PetscInt       n;
27847c6ae99SBarry Smith   MPI_Comm       comm;
27947c6ae99SBarry Smith   PetscMPIInt    size;
28047c6ae99SBarry Smith 
28147c6ae99SBarry Smith   PetscFunctionBegin;
28247c6ae99SBarry Smith   if (!dd->lx) {
28347c6ae99SBarry Smith     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&dd->lx);CHKERRQ(ierr);
284aa219208SBarry Smith     ierr = DMDAGetProcessorSubset(da,DMDA_X,dd->xs,&comm);CHKERRQ(ierr);
28547c6ae99SBarry Smith     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
28647c6ae99SBarry Smith     n = (dd->xe-dd->xs)/dd->w;
28747c6ae99SBarry Smith     ierr = MPI_Allgather(&n,1,MPIU_INT,dd->lx,1,MPIU_INT,comm);CHKERRQ(ierr);
28847c6ae99SBarry Smith   }
28947c6ae99SBarry Smith   if (dd->dim > 1 && !dd->ly) {
29047c6ae99SBarry Smith     ierr = PetscMalloc(dd->n*sizeof(PetscInt),&dd->ly);CHKERRQ(ierr);
291aa219208SBarry Smith     ierr = DMDAGetProcessorSubset(da,DMDA_Y,dd->ys,&comm);CHKERRQ(ierr);
29247c6ae99SBarry Smith     n = dd->ye-dd->ys;
29347c6ae99SBarry Smith     ierr = MPI_Allgather(&n,1,MPIU_INT,dd->ly,1,MPIU_INT,comm);CHKERRQ(ierr);
29447c6ae99SBarry Smith   }
29547c6ae99SBarry Smith   if (dd->dim > 2 && !dd->lz) {
29647c6ae99SBarry Smith     ierr = PetscMalloc(dd->p*sizeof(PetscInt),&dd->lz);CHKERRQ(ierr);
297aa219208SBarry Smith     ierr = DMDAGetProcessorSubset(da,DMDA_Z,dd->zs,&comm);CHKERRQ(ierr);
29847c6ae99SBarry Smith     n = dd->ze-dd->zs;
29947c6ae99SBarry Smith     ierr = MPI_Allgather(&n,1,MPIU_INT,dd->lz,1,MPIU_INT,comm);CHKERRQ(ierr);
30047c6ae99SBarry Smith   }
30147c6ae99SBarry Smith   PetscFunctionReturn(0);
30247c6ae99SBarry Smith }
30347c6ae99SBarry Smith 
30447c6ae99SBarry Smith #undef __FUNCT__
305aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType"
30647c6ae99SBarry Smith /*@
307aa219208SBarry Smith        DMDASetInterpolationType - Sets the type of interpolation that will be
30894013140SBarry Smith           returned by DMGetInterpolation()
30947c6ae99SBarry Smith 
310aa219208SBarry Smith    Logically Collective on DMDA
31147c6ae99SBarry Smith 
31247c6ae99SBarry Smith    Input Parameter:
31347c6ae99SBarry Smith +  da - initial distributed array
314aa219208SBarry Smith .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
31547c6ae99SBarry Smith 
31647c6ae99SBarry Smith    Level: intermediate
31747c6ae99SBarry Smith 
318aa219208SBarry Smith    Notes: you should call this on the coarser of the two DMDAs you pass to DMGetInterpolation()
31947c6ae99SBarry Smith 
32047c6ae99SBarry Smith .keywords:  distributed array, interpolation
32147c6ae99SBarry Smith 
322aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
32347c6ae99SBarry Smith @*/
324aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
32547c6ae99SBarry Smith {
32647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
32747c6ae99SBarry Smith 
32847c6ae99SBarry Smith   PetscFunctionBegin;
32947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
33047c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,ctype,2);
33147c6ae99SBarry Smith   dd->interptype = ctype;
33247c6ae99SBarry Smith   PetscFunctionReturn(0);
33347c6ae99SBarry Smith }
33447c6ae99SBarry Smith 
33547c6ae99SBarry Smith 
33647c6ae99SBarry Smith #undef __FUNCT__
337aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors"
33847c6ae99SBarry Smith /*@C
339aa219208SBarry Smith       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
34047c6ae99SBarry Smith         processes neighbors.
34147c6ae99SBarry Smith 
34247c6ae99SBarry Smith     Not Collective
34347c6ae99SBarry Smith 
34447c6ae99SBarry Smith    Input Parameter:
345aa219208SBarry Smith .     da - the DMDA object
34647c6ae99SBarry Smith 
34747c6ae99SBarry Smith    Output Parameters:
34847c6ae99SBarry Smith .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
34947c6ae99SBarry Smith               this process itself is in the list
35047c6ae99SBarry Smith 
35147c6ae99SBarry Smith    Notes: In 2d the array is of length 9, in 3d of length 27
35247c6ae99SBarry Smith           Not supported in 1d
353aa219208SBarry Smith           Do not free the array, it is freed when the DMDA is destroyed.
35447c6ae99SBarry Smith 
35547c6ae99SBarry Smith    Fortran Notes: In fortran you must pass in an array of the appropriate length.
35647c6ae99SBarry Smith 
35747c6ae99SBarry Smith    Level: intermediate
35847c6ae99SBarry Smith 
35947c6ae99SBarry Smith @*/
360aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
36147c6ae99SBarry Smith {
36247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
36347c6ae99SBarry Smith   PetscFunctionBegin;
36447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
36547c6ae99SBarry Smith   *ranks = dd->neighbors;
36647c6ae99SBarry Smith   PetscFunctionReturn(0);
36747c6ae99SBarry Smith }
36847c6ae99SBarry Smith 
36947c6ae99SBarry Smith /*@C
370aa219208SBarry Smith       DMDASetElementType - Sets the element type to be returned by DMGetElements()
37147c6ae99SBarry Smith 
37247c6ae99SBarry Smith     Not Collective
37347c6ae99SBarry Smith 
37447c6ae99SBarry Smith    Input Parameter:
375aa219208SBarry Smith .     da - the DMDA object
37647c6ae99SBarry Smith 
37747c6ae99SBarry Smith    Output Parameters:
378aa219208SBarry Smith .     etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1
37947c6ae99SBarry Smith 
38047c6ae99SBarry Smith    Level: intermediate
38147c6ae99SBarry Smith 
382aa219208SBarry Smith .seealso: DMDAElementType, DMDAGetElementType(), DMGetElements(), DMRestoreElements()
38347c6ae99SBarry Smith @*/
38447c6ae99SBarry Smith #undef __FUNCT__
385aa219208SBarry Smith #define __FUNCT__ "DMDASetElementType"
386aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetElementType(DM da, DMDAElementType etype)
38747c6ae99SBarry Smith {
38847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
38947c6ae99SBarry Smith   PetscErrorCode ierr;
39047c6ae99SBarry Smith 
39147c6ae99SBarry Smith   PetscFunctionBegin;
39247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
39347c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,etype,2);
39447c6ae99SBarry Smith   if (dd->elementtype != etype) {
39547c6ae99SBarry Smith     ierr = PetscFree(dd->e);CHKERRQ(ierr);
39647c6ae99SBarry Smith     dd->elementtype = etype;
39747c6ae99SBarry Smith     dd->ne          = 0;
39847c6ae99SBarry Smith     dd->e           = PETSC_NULL;
39947c6ae99SBarry Smith   }
40047c6ae99SBarry Smith   PetscFunctionReturn(0);
40147c6ae99SBarry Smith }
40247c6ae99SBarry Smith 
40347c6ae99SBarry Smith /*@C
404aa219208SBarry Smith       DMDAGetElementType - Gets the element type to be returned by DMGetElements()
40547c6ae99SBarry Smith 
40647c6ae99SBarry Smith     Not Collective
40747c6ae99SBarry Smith 
40847c6ae99SBarry Smith    Input Parameter:
409aa219208SBarry Smith .     da - the DMDA object
41047c6ae99SBarry Smith 
41147c6ae99SBarry Smith    Output Parameters:
412aa219208SBarry Smith .     etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1
41347c6ae99SBarry Smith 
41447c6ae99SBarry Smith    Level: intermediate
41547c6ae99SBarry Smith 
416aa219208SBarry Smith .seealso: DMDAElementType, DMDASetElementType(), DMGetElements(), DMRestoreElements()
41747c6ae99SBarry Smith @*/
41847c6ae99SBarry Smith #undef __FUNCT__
419aa219208SBarry Smith #define __FUNCT__ "DMDAGetElementType"
420aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAGetElementType(DM da, DMDAElementType *etype)
42147c6ae99SBarry Smith {
42247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
42347c6ae99SBarry Smith   PetscFunctionBegin;
42447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
42547c6ae99SBarry Smith   PetscValidPointer(etype,2);
42647c6ae99SBarry Smith   *etype = dd->elementtype;
42747c6ae99SBarry Smith   PetscFunctionReturn(0);
42847c6ae99SBarry Smith }
42947c6ae99SBarry Smith 
43047c6ae99SBarry Smith #undef __FUNCT__
43147c6ae99SBarry Smith #define __FUNCT__ "DMGetElements"
43247c6ae99SBarry Smith /*@C
43347c6ae99SBarry Smith       DMGetElements - Gets an array containing the indices (in local coordinates)
43447c6ae99SBarry Smith                  of all the local elements
43547c6ae99SBarry Smith 
43647c6ae99SBarry Smith     Not Collective
43747c6ae99SBarry Smith 
43847c6ae99SBarry Smith    Input Parameter:
43947c6ae99SBarry Smith .     dm - the DM object
44047c6ae99SBarry Smith 
44147c6ae99SBarry Smith    Output Parameters:
442454e267fSLisandro Dalcin +     nel - number of local elements
443454e267fSLisandro Dalcin .     nen - number of element nodes
44447c6ae99SBarry Smith -     e - the indices of the elements vertices
44547c6ae99SBarry Smith 
44647c6ae99SBarry Smith    Level: intermediate
44747c6ae99SBarry Smith 
44847c6ae99SBarry Smith .seealso: DMElementType, DMSetElementType(), DMRestoreElements()
44947c6ae99SBarry Smith @*/
450454e267fSLisandro Dalcin PetscErrorCode PETSCDM_DLLEXPORT DMGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
45147c6ae99SBarry Smith {
45247c6ae99SBarry Smith   PetscErrorCode ierr;
45347c6ae99SBarry Smith   PetscFunctionBegin;
45447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
455454e267fSLisandro Dalcin   PetscValidIntPointer(nel,2);
456454e267fSLisandro Dalcin   PetscValidIntPointer(nen,3);
457454e267fSLisandro Dalcin   PetscValidPointer(e,4);
458454e267fSLisandro Dalcin   ierr = (dm->ops->getelements)(dm,nel,nen,e);CHKERRQ(ierr);
45947c6ae99SBarry Smith   PetscFunctionReturn(0);
46047c6ae99SBarry Smith }
46147c6ae99SBarry Smith 
46247c6ae99SBarry Smith #undef __FUNCT__
46347c6ae99SBarry Smith #define __FUNCT__ "DMRestoreElements"
46447c6ae99SBarry Smith /*@C
46547c6ae99SBarry Smith       DMRestoreElements - Returns an array containing the indices (in local coordinates)
46647c6ae99SBarry Smith                  of all the local elements obtained with DMGetElements()
46747c6ae99SBarry Smith 
46847c6ae99SBarry Smith     Not Collective
46947c6ae99SBarry Smith 
47047c6ae99SBarry Smith    Input Parameter:
47147c6ae99SBarry Smith +     dm - the DM object
472454e267fSLisandro Dalcin .     nel - number of local elements
473454e267fSLisandro Dalcin .     nen - number of element nodes
47447c6ae99SBarry Smith -     e - the indices of the elements vertices
47547c6ae99SBarry Smith 
47647c6ae99SBarry Smith    Level: intermediate
47747c6ae99SBarry Smith 
47847c6ae99SBarry Smith .seealso: DMElementType, DMSetElementType(), DMGetElements()
47947c6ae99SBarry Smith @*/
480454e267fSLisandro Dalcin PetscErrorCode PETSCDM_DLLEXPORT DMRestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
48147c6ae99SBarry Smith {
48247c6ae99SBarry Smith   PetscErrorCode ierr;
48347c6ae99SBarry Smith   PetscFunctionBegin;
48447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
485454e267fSLisandro Dalcin   PetscValidIntPointer(nel,2);
486454e267fSLisandro Dalcin   PetscValidIntPointer(nen,3);
487454e267fSLisandro Dalcin   PetscValidPointer(e,4);
48847c6ae99SBarry Smith   if (dm->ops->restoreelements) {
489454e267fSLisandro Dalcin     ierr = (dm->ops->restoreelements)(dm,nel,nen,e);CHKERRQ(ierr);
49047c6ae99SBarry Smith   }
49147c6ae99SBarry Smith   PetscFunctionReturn(0);
49247c6ae99SBarry Smith }
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith #undef __FUNCT__
495aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges"
49647c6ae99SBarry Smith /*@C
497aa219208SBarry Smith       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith     Not Collective
50047c6ae99SBarry Smith 
50147c6ae99SBarry Smith    Input Parameter:
502aa219208SBarry Smith .     da - the DMDA object
50347c6ae99SBarry Smith 
50447c6ae99SBarry Smith    Output Parameter:
50547c6ae99SBarry Smith +     lx - ownership along x direction (optional)
50647c6ae99SBarry Smith .     ly - ownership along y direction (optional)
50747c6ae99SBarry Smith -     lz - ownership along z direction (optional)
50847c6ae99SBarry Smith 
50947c6ae99SBarry Smith    Level: intermediate
51047c6ae99SBarry Smith 
511aa219208SBarry Smith     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
514aa219208SBarry Smith     eighth arguments from DMDAGetInfo()
51547c6ae99SBarry Smith 
51647c6ae99SBarry Smith      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
517aa219208SBarry Smith     DMDA they came from still exists (has not been destroyed).
51847c6ae99SBarry Smith 
519aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
52047c6ae99SBarry Smith @*/
521aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
52247c6ae99SBarry Smith {
52347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
52447c6ae99SBarry Smith 
52547c6ae99SBarry Smith   PetscFunctionBegin;
52647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
52747c6ae99SBarry Smith   if (lx) *lx = dd->lx;
52847c6ae99SBarry Smith   if (ly) *ly = dd->ly;
52947c6ae99SBarry Smith   if (lz) *lz = dd->lz;
53047c6ae99SBarry Smith   PetscFunctionReturn(0);
53147c6ae99SBarry Smith }
53247c6ae99SBarry Smith 
53347c6ae99SBarry Smith #undef __FUNCT__
534aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor"
53547c6ae99SBarry Smith /*@
536aa219208SBarry Smith      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
53747c6ae99SBarry Smith 
538aa219208SBarry Smith     Logically Collective on DMDA
53947c6ae99SBarry Smith 
54047c6ae99SBarry Smith   Input Parameters:
541aa219208SBarry Smith +    da - the DMDA object
54247c6ae99SBarry Smith .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
54347c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
54447c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
54547c6ae99SBarry Smith 
54647c6ae99SBarry Smith   Options Database:
54747c6ae99SBarry Smith +  -da_refine_x - refinement ratio in x direction
54847c6ae99SBarry Smith .  -da_refine_y - refinement ratio in y direction
54947c6ae99SBarry Smith -  -da_refine_z - refinement ratio in z direction
55047c6ae99SBarry Smith 
55147c6ae99SBarry Smith   Level: intermediate
55247c6ae99SBarry Smith 
55347c6ae99SBarry Smith     Notes: Pass PETSC_IGNORE to leave a value unchanged
55447c6ae99SBarry Smith 
555aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor()
55647c6ae99SBarry Smith @*/
557aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
55847c6ae99SBarry Smith {
55947c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
56047c6ae99SBarry Smith 
56147c6ae99SBarry Smith   PetscFunctionBegin;
56247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
56347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_x,2);
56447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_y,3);
56547c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_z,4);
56647c6ae99SBarry Smith 
56747c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
56847c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
56947c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
57047c6ae99SBarry Smith   PetscFunctionReturn(0);
57147c6ae99SBarry Smith }
57247c6ae99SBarry Smith 
57347c6ae99SBarry Smith #undef __FUNCT__
574aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor"
57547c6ae99SBarry Smith /*@C
576aa219208SBarry Smith      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
57747c6ae99SBarry Smith 
57847c6ae99SBarry Smith     Not Collective
57947c6ae99SBarry Smith 
58047c6ae99SBarry Smith   Input Parameter:
581aa219208SBarry Smith .    da - the DMDA object
58247c6ae99SBarry Smith 
58347c6ae99SBarry Smith   Output Parameters:
58447c6ae99SBarry Smith +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
58547c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
58647c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
58747c6ae99SBarry Smith 
58847c6ae99SBarry Smith   Level: intermediate
58947c6ae99SBarry Smith 
59047c6ae99SBarry Smith     Notes: Pass PETSC_NULL for values you do not need
59147c6ae99SBarry Smith 
592aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor()
59347c6ae99SBarry Smith @*/
594aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
59547c6ae99SBarry Smith {
59647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
59747c6ae99SBarry Smith 
59847c6ae99SBarry Smith   PetscFunctionBegin;
59947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
60047c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
60147c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
60247c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
60347c6ae99SBarry Smith   PetscFunctionReturn(0);
60447c6ae99SBarry Smith }
60547c6ae99SBarry Smith 
60647c6ae99SBarry Smith #undef __FUNCT__
607aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix"
60847c6ae99SBarry Smith /*@C
609aa219208SBarry Smith      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
61047c6ae99SBarry Smith 
611aa219208SBarry Smith     Logically Collective on DMDA
61247c6ae99SBarry Smith 
61347c6ae99SBarry Smith   Input Parameters:
614aa219208SBarry Smith +    da - the DMDA object
615aa219208SBarry Smith -    f - the function that allocates the matrix for that specific DMDA
61647c6ae99SBarry Smith 
61747c6ae99SBarry Smith   Level: developer
61847c6ae99SBarry Smith 
619aa219208SBarry Smith    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
62047c6ae99SBarry Smith        the diagonal and off-diagonal blocks of the matrix
62147c6ae99SBarry Smith 
622aa219208SBarry Smith .seealso: DMGetMatrix(), DMDASetBlockFills()
62347c6ae99SBarry Smith @*/
624aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, const MatType,Mat*))
62547c6ae99SBarry Smith {
62647c6ae99SBarry Smith   PetscFunctionBegin;
62747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
62847c6ae99SBarry Smith   da->ops->getmatrix = f;
62947c6ae99SBarry Smith   PetscFunctionReturn(0);
63047c6ae99SBarry Smith }
63147c6ae99SBarry Smith 
63247c6ae99SBarry Smith #undef __FUNCT__
63394013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges"
63447c6ae99SBarry Smith /*
63547c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
63647c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
63747c6ae99SBarry Smith 
63847c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
63947c6ae99SBarry Smith */
64094013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
64147c6ae99SBarry Smith {
64247c6ae99SBarry Smith   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
64347c6ae99SBarry Smith   PetscErrorCode ierr;
64447c6ae99SBarry Smith 
64547c6ae99SBarry Smith   PetscFunctionBegin;
64647c6ae99SBarry Smith   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
64747c6ae99SBarry Smith   if (ratio == 1) {
64847c6ae99SBarry Smith     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
64947c6ae99SBarry Smith     PetscFunctionReturn(0);
65047c6ae99SBarry Smith   }
65147c6ae99SBarry Smith   for (i=0; i<m; i++) totalc += lc[i];
65247c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
65347c6ae99SBarry Smith   for (i=0; i<m; i++) {
65447c6ae99SBarry Smith     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
65547c6ae99SBarry Smith     if (i == m-1) lf[i] = want;
65647c6ae99SBarry Smith     else {
6577aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
6587aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
6597aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
6607aca7175SJed Brown       while ((startf+want)/ratio < nextc - stencil_width) want++;
6617aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
6627aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
6637aca7175SJed Brown       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
6647aca7175SJed Brown       /* Make sure all constraints are satisfied */
6657aca7175SJed Brown       if (want < 0 || want > remaining
6667aca7175SJed Brown           || ((startf+want)/ratio < nextc - stencil_width)
6677aca7175SJed Brown           || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width))
6687aca7175SJed Brown         SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range");
66947c6ae99SBarry Smith     }
67047c6ae99SBarry Smith     lf[i] = want;
67147c6ae99SBarry Smith     startc += lc[i];
67247c6ae99SBarry Smith     startf += lf[i];
67347c6ae99SBarry Smith     remaining -= lf[i];
67447c6ae99SBarry Smith   }
67547c6ae99SBarry Smith   PetscFunctionReturn(0);
67647c6ae99SBarry Smith }
67747c6ae99SBarry Smith 
67847c6ae99SBarry Smith #undef __FUNCT__
6792be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges"
6802be375d4SJed Brown /*
6812be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
6822be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
6832be375d4SJed Brown 
6842be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
6852be375d4SJed Brown */
6862be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
6872be375d4SJed Brown {
6882be375d4SJed Brown   PetscInt       i,totalf,remaining,startc,startf;
6892be375d4SJed Brown   PetscErrorCode ierr;
6902be375d4SJed Brown 
6912be375d4SJed Brown   PetscFunctionBegin;
6922be375d4SJed Brown   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
6932be375d4SJed Brown   if (ratio == 1) {
6942be375d4SJed Brown     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
6952be375d4SJed Brown     PetscFunctionReturn(0);
6962be375d4SJed Brown   }
6972be375d4SJed Brown   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
6982be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
6992be375d4SJed Brown   for (i=0,startc=0,startf=0; i<m; i++) {
7002be375d4SJed Brown     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
7012be375d4SJed Brown     if (i == m-1) lc[i] = want;
7022be375d4SJed Brown     else {
7032be375d4SJed Brown       const PetscInt nextf = startf+lf[i];
7042be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
7052be375d4SJed Brown        * node is within one stencil width. */
7062be375d4SJed Brown       while (nextf/ratio < startc+want-stencil_width) want--;
7072be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
7082be375d4SJed Brown        * fine node is within one stencil width. */
7092be375d4SJed Brown       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
7102be375d4SJed Brown       if (want < 0 || want > remaining
7112be375d4SJed Brown           || (nextf/ratio < startc+want-stencil_width)
7122be375d4SJed Brown           || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width))
7132be375d4SJed Brown         SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range");
7142be375d4SJed Brown     }
7152be375d4SJed Brown     lc[i] = want;
7162be375d4SJed Brown     startc += lc[i];
7172be375d4SJed Brown     startf += lf[i];
7182be375d4SJed Brown     remaining -= lc[i];
7192be375d4SJed Brown   }
7202be375d4SJed Brown   PetscFunctionReturn(0);
7212be375d4SJed Brown }
7222be375d4SJed Brown 
7232be375d4SJed Brown #undef __FUNCT__
72494013140SBarry Smith #define __FUNCT__ "DMRefine_DA"
72594013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
72647c6ae99SBarry Smith {
72747c6ae99SBarry Smith   PetscErrorCode ierr;
72847c6ae99SBarry Smith   PetscInt       M,N,P;
7299a42bb27SBarry Smith   DM             da2;
73047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
73147c6ae99SBarry Smith 
73247c6ae99SBarry Smith   PetscFunctionBegin;
73347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
73447c6ae99SBarry Smith   PetscValidPointer(daref,3);
73547c6ae99SBarry Smith 
736aa219208SBarry Smith   if (DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){
73747c6ae99SBarry Smith     M = dd->refine_x*dd->M;
73847c6ae99SBarry Smith   } else {
73947c6ae99SBarry Smith     M = 1 + dd->refine_x*(dd->M - 1);
74047c6ae99SBarry Smith   }
741aa219208SBarry Smith   if (DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){
74247c6ae99SBarry Smith     N = dd->refine_y*dd->N;
74347c6ae99SBarry Smith   } else {
74447c6ae99SBarry Smith     N = 1 + dd->refine_y*(dd->N - 1);
74547c6ae99SBarry Smith   }
746aa219208SBarry Smith   if (DMDAZPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){
74747c6ae99SBarry Smith     P = dd->refine_z*dd->P;
74847c6ae99SBarry Smith   } else {
74947c6ae99SBarry Smith     P = 1 + dd->refine_z*(dd->P - 1);
75047c6ae99SBarry Smith   }
751aa219208SBarry Smith   ierr = DMDACreateOwnershipRanges(da);CHKERRQ(ierr);
75247c6ae99SBarry Smith   if (dd->dim == 3) {
75347c6ae99SBarry Smith     PetscInt *lx,*ly,*lz;
75447c6ae99SBarry Smith     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
755aa219208SBarry Smith     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
756aa219208SBarry Smith     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
757aa219208SBarry Smith     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAZPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
758aa219208SBarry Smith     ierr = DMDACreate3d(((PetscObject)da)->comm,dd->wrap,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr);
75947c6ae99SBarry Smith     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
76047c6ae99SBarry Smith   } else if (dd->dim == 2) {
76147c6ae99SBarry Smith     PetscInt *lx,*ly;
76247c6ae99SBarry Smith     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
763aa219208SBarry Smith     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
764aa219208SBarry Smith     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
765aa219208SBarry Smith     ierr = DMDACreate2d(((PetscObject)da)->comm,dd->wrap,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr);
76647c6ae99SBarry Smith     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
76747c6ae99SBarry Smith   } else if (dd->dim == 1) {
76847c6ae99SBarry Smith     PetscInt *lx;
76947c6ae99SBarry Smith     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
770aa219208SBarry Smith     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
771aa219208SBarry Smith     ierr = DMDACreate1d(((PetscObject)da)->comm,dd->wrap,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr);
77247c6ae99SBarry Smith     ierr = PetscFree(lx);CHKERRQ(ierr);
77347c6ae99SBarry Smith   }
77447c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
77547c6ae99SBarry Smith 
77647c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
77747c6ae99SBarry Smith   da2->ops->getmatrix        = da->ops->getmatrix;
77847c6ae99SBarry Smith   da2->ops->getinterpolation = da->ops->getinterpolation;
77947c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
78047c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
78147c6ae99SBarry Smith 
78247c6ae99SBarry Smith   /* copy fill information if given */
78347c6ae99SBarry Smith   if (dd->dfill) {
78447c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
78547c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
78647c6ae99SBarry Smith   }
78747c6ae99SBarry Smith   if (dd->ofill) {
78847c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
78947c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
79047c6ae99SBarry Smith   }
79147c6ae99SBarry Smith   /* copy the refine information */
79247c6ae99SBarry Smith   dd2->refine_x = dd->refine_x;
79347c6ae99SBarry Smith   dd2->refine_y = dd->refine_y;
79447c6ae99SBarry Smith   dd2->refine_z = dd->refine_z;
79547c6ae99SBarry Smith 
79647c6ae99SBarry Smith   /* copy vector type information */
79747c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
79847c6ae99SBarry Smith   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
799*ddcf8b74SDave May 
800*ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
801*ddcf8b74SDave May   if (dd->coordinates) {
802*ddcf8b74SDave May     DM  cdaf,cdac;
803*ddcf8b74SDave May     Vec coordsc,coordsf;
804*ddcf8b74SDave May     Mat II;
805*ddcf8b74SDave May 
806*ddcf8b74SDave May     ierr = DMDAGetCoordinateDA(da,&cdac);CHKERRQ(ierr);
807*ddcf8b74SDave May     ierr = DMDAGetCoordinateDA(da2,&cdaf);CHKERRQ(ierr);
808*ddcf8b74SDave May     ierr = DMDAGetCoordinates(da,&coordsc);CHKERRQ(ierr);
809*ddcf8b74SDave May     ierr = DMDAGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
810*ddcf8b74SDave May     ierr = DMGetInterpolation(cdac,cdaf,&II,PETSC_NULL);CHKERRQ(ierr);
811*ddcf8b74SDave May     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
812*ddcf8b74SDave May     ierr = MatDestroy(II);CHKERRQ(ierr);
813*ddcf8b74SDave May   }
81447c6ae99SBarry Smith   *daref = da2;
81547c6ae99SBarry Smith   PetscFunctionReturn(0);
81647c6ae99SBarry Smith }
81747c6ae99SBarry Smith 
81847c6ae99SBarry Smith #undef __FUNCT__
81994013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA"
82094013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
82147c6ae99SBarry Smith {
82247c6ae99SBarry Smith   PetscErrorCode ierr;
82347c6ae99SBarry Smith   PetscInt       M,N,P;
8249a42bb27SBarry Smith   DM             da2;
82547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
82647c6ae99SBarry Smith 
82747c6ae99SBarry Smith   PetscFunctionBegin;
82847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
82947c6ae99SBarry Smith   PetscValidPointer(daref,3);
83047c6ae99SBarry Smith 
831aa219208SBarry Smith   if (DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){
83247c6ae99SBarry Smith     M = dd->M / dd->refine_x;
83347c6ae99SBarry Smith   } else {
83447c6ae99SBarry Smith     M = 1 + (dd->M - 1) / dd->refine_x;
83547c6ae99SBarry Smith   }
836aa219208SBarry Smith   if (DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){
83747c6ae99SBarry Smith     N = dd->N / dd->refine_y;
83847c6ae99SBarry Smith   } else {
83947c6ae99SBarry Smith     N = 1 + (dd->N - 1) / dd->refine_y;
84047c6ae99SBarry Smith   }
841aa219208SBarry Smith   if (DMDAZPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){
84247c6ae99SBarry Smith     P = dd->P / dd->refine_z;
84347c6ae99SBarry Smith   } else {
84447c6ae99SBarry Smith     P = 1 + (dd->P - 1) / dd->refine_z;
84547c6ae99SBarry Smith   }
8462be375d4SJed Brown   ierr = DMDACreateOwnershipRanges(da);CHKERRQ(ierr);
84747c6ae99SBarry Smith   if (dd->dim == 3) {
8482be375d4SJed Brown     PetscInt *lx,*ly,*lz;
8492be375d4SJed Brown     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
8502be375d4SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
8512be375d4SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
8522be375d4SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAZPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
8532be375d4SJed Brown     ierr = DMDACreate3d(((PetscObject)da)->comm,dd->wrap,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr);
8542be375d4SJed Brown     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
85547c6ae99SBarry Smith   } else if (dd->dim == 2) {
8562be375d4SJed Brown     PetscInt *lx,*ly;
8572be375d4SJed Brown     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
8582be375d4SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
8592be375d4SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
8602be375d4SJed Brown     ierr = DMDACreate2d(((PetscObject)da)->comm,dd->wrap,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr);
8612be375d4SJed Brown     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
86247c6ae99SBarry Smith   } else if (dd->dim == 1) {
8632be375d4SJed Brown     PetscInt *lx;
8642be375d4SJed Brown     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
8652be375d4SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
8662be375d4SJed Brown     ierr = DMDACreate1d(((PetscObject)da)->comm,dd->wrap,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr);
8672be375d4SJed Brown     ierr = PetscFree(lx);CHKERRQ(ierr);
86847c6ae99SBarry Smith   }
86947c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
87047c6ae99SBarry Smith 
87147c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
87247c6ae99SBarry Smith   da2->ops->getmatrix        = da->ops->getmatrix;
87347c6ae99SBarry Smith   da2->ops->getinterpolation = da->ops->getinterpolation;
87447c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
87547c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
87647c6ae99SBarry Smith 
87747c6ae99SBarry Smith   /* copy fill information if given */
87847c6ae99SBarry Smith   if (dd->dfill) {
87947c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
88047c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
88147c6ae99SBarry Smith   }
88247c6ae99SBarry Smith   if (dd->ofill) {
88347c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
88447c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
88547c6ae99SBarry Smith   }
88647c6ae99SBarry Smith   /* copy the refine information */
88747c6ae99SBarry Smith   dd2->refine_x = dd->refine_x;
88847c6ae99SBarry Smith   dd2->refine_y = dd->refine_y;
88947c6ae99SBarry Smith   dd2->refine_z = dd->refine_z;
89047c6ae99SBarry Smith 
89147c6ae99SBarry Smith   /* copy vector type information */
89247c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
89347c6ae99SBarry Smith   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
89447c6ae99SBarry Smith 
895*ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
896*ddcf8b74SDave May   if (dd->coordinates) {
897*ddcf8b74SDave May     DM  cdaf,cdac;
898*ddcf8b74SDave May     Vec coordsc,coordsf;
899*ddcf8b74SDave May     VecScatter inject;
900*ddcf8b74SDave May 
901*ddcf8b74SDave May     ierr = DMDAGetCoordinateDA(da,&cdaf);CHKERRQ(ierr);
902*ddcf8b74SDave May     ierr = DMDAGetCoordinateDA(da2,&cdac);CHKERRQ(ierr);
903*ddcf8b74SDave May     ierr = DMDAGetCoordinates(da,&coordsf);CHKERRQ(ierr);
904*ddcf8b74SDave May     ierr = DMDAGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
905*ddcf8b74SDave May 
906*ddcf8b74SDave May     ierr = DMGetInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
907*ddcf8b74SDave May     ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
908*ddcf8b74SDave May     ierr = VecScatterEnd(inject  ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
909*ddcf8b74SDave May     ierr = VecScatterDestroy(inject);CHKERRQ(ierr);
910*ddcf8b74SDave May   }
91147c6ae99SBarry Smith   *daref = da2;
91247c6ae99SBarry Smith   PetscFunctionReturn(0);
91347c6ae99SBarry Smith }
91447c6ae99SBarry Smith 
91547c6ae99SBarry Smith #undef __FUNCT__
91694013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA"
91794013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
91847c6ae99SBarry Smith {
91947c6ae99SBarry Smith   PetscErrorCode ierr;
92047c6ae99SBarry Smith   PetscInt       i,n,*refx,*refy,*refz;
92147c6ae99SBarry Smith 
92247c6ae99SBarry Smith   PetscFunctionBegin;
92347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
92447c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
92547c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
92647c6ae99SBarry Smith   PetscValidPointer(daf,3);
92747c6ae99SBarry Smith 
928aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
92947c6ae99SBarry Smith   ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr);
93047c6ae99SBarry Smith   for (i=0; i<nlevels; i++) {
931aa219208SBarry Smith     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
93247c6ae99SBarry Smith   }
93347c6ae99SBarry Smith   n = nlevels;
93447c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
93547c6ae99SBarry Smith   n = nlevels;
93647c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
93747c6ae99SBarry Smith   n = nlevels;
93847c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
93947c6ae99SBarry Smith 
940aa219208SBarry Smith   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
94194013140SBarry Smith   ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr);
94247c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
943aa219208SBarry Smith     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
94494013140SBarry Smith     ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr);
94547c6ae99SBarry Smith   }
94647c6ae99SBarry Smith   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
94747c6ae99SBarry Smith   PetscFunctionReturn(0);
94847c6ae99SBarry Smith }
94947c6ae99SBarry Smith 
95047c6ae99SBarry Smith #undef __FUNCT__
95194013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA"
95294013140SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
95347c6ae99SBarry Smith {
95447c6ae99SBarry Smith   PetscErrorCode ierr;
95547c6ae99SBarry Smith   PetscInt       i;
95647c6ae99SBarry Smith 
95747c6ae99SBarry Smith   PetscFunctionBegin;
95847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
95947c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
96047c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
96147c6ae99SBarry Smith   PetscValidPointer(dac,3);
96294013140SBarry Smith   ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr);
96347c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
96494013140SBarry Smith     ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr);
96547c6ae99SBarry Smith   }
96647c6ae99SBarry Smith   PetscFunctionReturn(0);
96747c6ae99SBarry Smith }
968