xref: /petsc/src/dm/impls/da/da.c (revision e3f3e4b6e8c9f9971bdd3dd5f3154bef081df595)
1b45d2f2cSJed Brown #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
247c6ae99SBarry Smith 
347c6ae99SBarry Smith 
447c6ae99SBarry Smith #undef __FUNCT__
5aa219208SBarry Smith #define __FUNCT__ "DMDASetDim"
647c6ae99SBarry Smith /*@
7aa219208SBarry Smith   DMDASetDim - Sets the dimension
847c6ae99SBarry Smith 
9aa219208SBarry Smith   Collective on DMDA
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith   Input Parameters:
12aa219208SBarry Smith + da - the DMDA
1347c6ae99SBarry Smith - dim - the dimension (or PETSC_DECIDE)
1447c6ae99SBarry Smith 
1547c6ae99SBarry Smith   Level: intermediate
1647c6ae99SBarry Smith 
17aa219208SBarry Smith .seealso: DaGetDim(), DMDASetSizes()
1847c6ae99SBarry Smith @*/
197087cfbeSBarry Smith PetscErrorCode  DMDASetDim(DM da, PetscInt dim)
2047c6ae99SBarry Smith {
2147c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
2247c6ae99SBarry Smith 
2347c6ae99SBarry Smith   PetscFunctionBegin;
2447c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
25aa219208SBarry Smith   if (dd->dim > 0 && dim != dd->dim) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change DMDA dim from %D after it was set to %D",dd->dim,dim);
2647c6ae99SBarry Smith   dd->dim = dim;
2747c6ae99SBarry Smith   PetscFunctionReturn(0);
2847c6ae99SBarry Smith }
2947c6ae99SBarry Smith 
3047c6ae99SBarry Smith #undef __FUNCT__
31aa219208SBarry Smith #define __FUNCT__ "DMDASetSizes"
3247c6ae99SBarry Smith /*@
33aa219208SBarry Smith   DMDASetSizes - Sets the global sizes
3447c6ae99SBarry Smith 
35aa219208SBarry Smith   Logically Collective on DMDA
3647c6ae99SBarry Smith 
3747c6ae99SBarry Smith   Input Parameters:
38aa219208SBarry Smith + da - the DMDA
3947c6ae99SBarry Smith . M - the global X size (or PETSC_DECIDE)
4047c6ae99SBarry Smith . N - the global Y size (or PETSC_DECIDE)
4147c6ae99SBarry Smith - P - the global Z size (or PETSC_DECIDE)
4247c6ae99SBarry Smith 
4347c6ae99SBarry Smith   Level: intermediate
4447c6ae99SBarry Smith 
45aa219208SBarry Smith .seealso: DMDAGetSize(), PetscSplitOwnership()
4647c6ae99SBarry Smith @*/
477087cfbeSBarry Smith PetscErrorCode  DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P)
4847c6ae99SBarry Smith {
4947c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
5047c6ae99SBarry Smith 
5147c6ae99SBarry Smith   PetscFunctionBegin;
5247c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
5347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,M,2);
5447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,N,3);
5547c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,P,4);
56cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
5747c6ae99SBarry Smith 
5847c6ae99SBarry Smith   dd->M = M;
5947c6ae99SBarry Smith   dd->N = N;
6047c6ae99SBarry Smith   dd->P = P;
6147c6ae99SBarry Smith   PetscFunctionReturn(0);
6247c6ae99SBarry Smith }
6347c6ae99SBarry Smith 
6447c6ae99SBarry Smith #undef __FUNCT__
65aa219208SBarry Smith #define __FUNCT__ "DMDASetNumProcs"
6647c6ae99SBarry Smith /*@
67aa219208SBarry Smith   DMDASetNumProcs - Sets the number of processes in each dimension
6847c6ae99SBarry Smith 
69aa219208SBarry Smith   Logically Collective on DMDA
7047c6ae99SBarry Smith 
7147c6ae99SBarry Smith   Input Parameters:
72aa219208SBarry Smith + da - the DMDA
7347c6ae99SBarry Smith . m - the number of X procs (or PETSC_DECIDE)
7447c6ae99SBarry Smith . n - the number of Y procs (or PETSC_DECIDE)
7547c6ae99SBarry Smith - p - the number of Z procs (or PETSC_DECIDE)
7647c6ae99SBarry Smith 
7747c6ae99SBarry Smith   Level: intermediate
7847c6ae99SBarry Smith 
79aa219208SBarry Smith .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership()
8047c6ae99SBarry Smith @*/
817087cfbeSBarry Smith PetscErrorCode  DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p)
8247c6ae99SBarry Smith {
8347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
84*e3f3e4b6SBarry Smith   PetscErrorCode ierr;
8547c6ae99SBarry Smith 
8647c6ae99SBarry Smith   PetscFunctionBegin;
8747c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
8847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,m,2);
8947c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,n,3);
9047c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,p,4);
91cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
9247c6ae99SBarry Smith   dd->m = m;
9347c6ae99SBarry Smith   dd->n = n;
9447c6ae99SBarry Smith   dd->p = p;
95*e3f3e4b6SBarry Smith   if (dd->dim == 2) {
96*e3f3e4b6SBarry Smith     PetscMPIInt size;
97*e3f3e4b6SBarry Smith     ierr = MPI_Comm_size(((PetscObject)da)->comm,&size);CHKERRQ(ierr);
98*e3f3e4b6SBarry Smith     if ((dd->m > 0) && (dd->n < 0)) {
99*e3f3e4b6SBarry Smith       dd->n = size/dd->m;
100*e3f3e4b6SBarry Smith       if (dd->n*dd->m != size) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"%D processes in X direction not divisible into comm size %d",m,size);
101*e3f3e4b6SBarry Smith     }
102*e3f3e4b6SBarry Smith     if ((dd->n > 0) && (dd->m < 0)) {
103*e3f3e4b6SBarry Smith       dd->m = size/dd->n;
104*e3f3e4b6SBarry Smith       if (dd->n*dd->m != size) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"%D processes in Y direction not divisible into comm size %d",n,size);
105*e3f3e4b6SBarry Smith     }
106*e3f3e4b6SBarry Smith   }
10747c6ae99SBarry Smith   PetscFunctionReturn(0);
10847c6ae99SBarry Smith }
10947c6ae99SBarry Smith 
11047c6ae99SBarry Smith #undef __FUNCT__
1113e24977eSEthan Coon #define __FUNCT__ "DMDASetBoundaryType"
11247c6ae99SBarry Smith /*@
1131321219cSEthan Coon   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries.
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith   Not collective
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith   Input Parameter:
118aa219208SBarry Smith + da    - The DMDA
1191321219cSEthan Coon - bx,by,bz - One of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC
12047c6ae99SBarry Smith 
12147c6ae99SBarry Smith   Level: intermediate
12247c6ae99SBarry Smith 
12347c6ae99SBarry Smith .keywords:  distributed array, periodicity
12496e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA, DMDABoundaryType
12547c6ae99SBarry Smith @*/
1261321219cSEthan Coon PetscErrorCode  DMDASetBoundaryType(DM da,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz)
12747c6ae99SBarry Smith {
12847c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith   PetscFunctionBegin;
13147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1325a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,bx,2);
1335a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,by,3);
1345a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,bz,4);
135cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1361321219cSEthan Coon   dd->bx = bx;
1371321219cSEthan Coon   dd->by = by;
1381321219cSEthan Coon   dd->bz = bz;
13947c6ae99SBarry Smith   PetscFunctionReturn(0);
14047c6ae99SBarry Smith }
14147c6ae99SBarry Smith 
14247c6ae99SBarry Smith #undef __FUNCT__
143aa219208SBarry Smith #define __FUNCT__ "DMDASetDof"
14447c6ae99SBarry Smith /*@
145aa219208SBarry Smith   DMDASetDof - Sets the number of degrees of freedom per vertex
14647c6ae99SBarry Smith 
14747c6ae99SBarry Smith   Not collective
14847c6ae99SBarry Smith 
14947c6ae99SBarry Smith   Input Parameter:
150aa219208SBarry Smith + da  - The DMDA
15147c6ae99SBarry Smith - dof - Number of degrees of freedom
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith   Level: intermediate
15447c6ae99SBarry Smith 
15547c6ae99SBarry Smith .keywords:  distributed array, degrees of freedom
15696e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
15747c6ae99SBarry Smith @*/
15854cfb0beSLisandro Dalcin PetscErrorCode  DMDASetDof(DM da, PetscInt dof)
15947c6ae99SBarry Smith {
16047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
16147c6ae99SBarry Smith 
16247c6ae99SBarry Smith   PetscFunctionBegin;
16347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
16454cfb0beSLisandro Dalcin   PetscValidLogicalCollectiveInt(da,dof,2);
165cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
16647c6ae99SBarry Smith   dd->w = dof;
1671411c6eeSJed Brown   da->bs = dof;
16847c6ae99SBarry Smith   PetscFunctionReturn(0);
16947c6ae99SBarry Smith }
17047c6ae99SBarry Smith 
17147c6ae99SBarry Smith #undef __FUNCT__
17288661749SPeter Brune #define __FUNCT__ "DMDASetOverlap"
17388661749SPeter Brune /*@
17488661749SPeter Brune   DMDASetOverlap - Sets the size of the per-processor overlap.
17588661749SPeter Brune 
17688661749SPeter Brune   Not collective
17788661749SPeter Brune 
17888661749SPeter Brune   Input Parameter:
17988661749SPeter Brune + da  - The DMDA
18088661749SPeter Brune - dof - Number of degrees of freedom
18188661749SPeter Brune 
18288661749SPeter Brune   Level: intermediate
18388661749SPeter Brune 
18488661749SPeter Brune .keywords:  distributed array, degrees of freedom
18588661749SPeter Brune .seealso: DMDACreate(), DMDestroy(), DMDA
18688661749SPeter Brune @*/
18788661749SPeter Brune PetscErrorCode  DMDASetOverlap(DM da, PetscInt overlap)
18888661749SPeter Brune {
18988661749SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
19088661749SPeter Brune 
19188661749SPeter Brune   PetscFunctionBegin;
19288661749SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
19388661749SPeter Brune   PetscValidLogicalCollectiveInt(da,overlap,2);
19488661749SPeter Brune   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
19588661749SPeter Brune   dd->overlap = overlap;
19688661749SPeter Brune   PetscFunctionReturn(0);
19788661749SPeter Brune }
19888661749SPeter Brune 
199d886c4f4SPeter Brune #undef __FUNCT__
200d886c4f4SPeter Brune #define __FUNCT__ "DMDASetOffset"
201d886c4f4SPeter Brune /*@
202d886c4f4SPeter Brune   DMDASetOffset - Sets the index offset of the DA.
203d886c4f4SPeter Brune 
204d886c4f4SPeter Brune   Collective on DA
205d886c4f4SPeter Brune 
206d886c4f4SPeter Brune   Input Parameter:
207d886c4f4SPeter Brune + da  - The DMDA
208d886c4f4SPeter Brune . xo  - The offset in the x direction
209d886c4f4SPeter Brune . yo  - The offset in the y direction
210d886c4f4SPeter Brune - zo  - The offset in the z direction
211d886c4f4SPeter Brune 
212d886c4f4SPeter Brune   Level: intermediate
213d886c4f4SPeter Brune 
214d886c4f4SPeter Brune   Notes: This is used primarily to overlap a computation on a local DA with that on a global DA without
215d886c4f4SPeter Brune   changing boundary conditions or subdomain features that depend upon the global offsets.
216d886c4f4SPeter Brune 
217d886c4f4SPeter Brune .keywords:  distributed array, degrees of freedom
218d886c4f4SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray()
219d886c4f4SPeter Brune @*/
220d886c4f4SPeter Brune PetscErrorCode  DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo)
221d886c4f4SPeter Brune {
222d886c4f4SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
223d886c4f4SPeter Brune 
224d886c4f4SPeter Brune   PetscFunctionBegin;
225d886c4f4SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
226d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da,xo,2);
227d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da,yo,2);
228d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da,zo,2);
229d886c4f4SPeter Brune   dd->xo = xo;
230d886c4f4SPeter Brune   dd->yo = yo;
231d886c4f4SPeter Brune   dd->zo = zo;
232d886c4f4SPeter Brune   PetscFunctionReturn(0);
233d886c4f4SPeter Brune }
234d886c4f4SPeter Brune 
235d886c4f4SPeter Brune #undef __FUNCT__
236d886c4f4SPeter Brune #define __FUNCT__ "DMDAGetOffset"
237d886c4f4SPeter Brune /*@
238d886c4f4SPeter Brune   DMDAGetOffset - Gets the index offset of the DA.
239d886c4f4SPeter Brune 
240d886c4f4SPeter Brune   Not collective
241d886c4f4SPeter Brune 
242d886c4f4SPeter Brune   Input Parameter:
243d886c4f4SPeter Brune . da  - The DMDA
244d886c4f4SPeter Brune 
245d886c4f4SPeter Brune   Output Parameters:
246d886c4f4SPeter Brune + xo  - The offset in the x direction
247d886c4f4SPeter Brune . yo  - The offset in the y direction
248d886c4f4SPeter Brune - zo  - The offset in the z direction
249d886c4f4SPeter Brune 
250d886c4f4SPeter Brune   Level: intermediate
251d886c4f4SPeter Brune 
252d886c4f4SPeter Brune .keywords:  distributed array, degrees of freedom
253d886c4f4SPeter Brune .seealso: DMDASetOffset(), DMDAVecGetArray()
254d886c4f4SPeter Brune @*/
255d886c4f4SPeter Brune PetscErrorCode  DMDAGetOffset(DM da, PetscInt *xo, PetscInt *yo, PetscInt *zo)
256d886c4f4SPeter Brune {
257d886c4f4SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
258d886c4f4SPeter Brune 
259d886c4f4SPeter Brune   PetscFunctionBegin;
260d886c4f4SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
261d886c4f4SPeter Brune   if (xo) *xo = dd->xo;
262d886c4f4SPeter Brune   if (yo) *yo = dd->yo;
263d886c4f4SPeter Brune   if (zo) *zo = dd->zo;
264d886c4f4SPeter Brune   PetscFunctionReturn(0);
265d886c4f4SPeter Brune }
266d886c4f4SPeter Brune 
26788661749SPeter Brune 
26888661749SPeter Brune #undef __FUNCT__
269aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType"
27047c6ae99SBarry Smith /*@
271aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
27247c6ae99SBarry Smith 
273aa219208SBarry Smith   Logically Collective on DMDA
27447c6ae99SBarry Smith 
27547c6ae99SBarry Smith   Input Parameter:
276aa219208SBarry Smith + da    - The DMDA
277aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
27847c6ae99SBarry Smith 
27947c6ae99SBarry Smith   Level: intermediate
28047c6ae99SBarry Smith 
28147c6ae99SBarry Smith .keywords:  distributed array, stencil
28296e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
28347c6ae99SBarry Smith @*/
2847087cfbeSBarry Smith PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
28547c6ae99SBarry Smith {
28647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
28747c6ae99SBarry Smith 
28847c6ae99SBarry Smith   PetscFunctionBegin;
28947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
29047c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,stype,2);
291cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
29247c6ae99SBarry Smith   dd->stencil_type = stype;
29347c6ae99SBarry Smith   PetscFunctionReturn(0);
29447c6ae99SBarry Smith }
29547c6ae99SBarry Smith 
29647c6ae99SBarry Smith #undef __FUNCT__
297aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth"
29847c6ae99SBarry Smith /*@
299aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
30047c6ae99SBarry Smith 
301aa219208SBarry Smith   Logically Collective on DMDA
30247c6ae99SBarry Smith 
30347c6ae99SBarry Smith   Input Parameter:
304aa219208SBarry Smith + da    - The DMDA
30547c6ae99SBarry Smith - width - The stencil width
30647c6ae99SBarry Smith 
30747c6ae99SBarry Smith   Level: intermediate
30847c6ae99SBarry Smith 
30947c6ae99SBarry Smith .keywords:  distributed array, stencil
31096e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
31147c6ae99SBarry Smith @*/
3127087cfbeSBarry Smith PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
31347c6ae99SBarry Smith {
31447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
31547c6ae99SBarry Smith 
31647c6ae99SBarry Smith   PetscFunctionBegin;
31747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
31847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,width,2);
319cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
32047c6ae99SBarry Smith   dd->s = width;
32147c6ae99SBarry Smith   PetscFunctionReturn(0);
32247c6ae99SBarry Smith }
32347c6ae99SBarry Smith 
32447c6ae99SBarry Smith #undef __FUNCT__
325aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
326aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
32747c6ae99SBarry Smith {
32847c6ae99SBarry Smith   PetscInt i,sum;
32947c6ae99SBarry Smith 
33047c6ae99SBarry Smith   PetscFunctionBegin;
33147c6ae99SBarry Smith   if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
33247c6ae99SBarry Smith   for (i=sum=0; i<m; i++) sum += lx[i];
33347c6ae99SBarry Smith   if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
33447c6ae99SBarry Smith   PetscFunctionReturn(0);
33547c6ae99SBarry Smith }
33647c6ae99SBarry Smith 
33747c6ae99SBarry Smith #undef __FUNCT__
338aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges"
33947c6ae99SBarry Smith /*@
340aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
34147c6ae99SBarry Smith 
342aa219208SBarry Smith   Logically Collective on DMDA
34347c6ae99SBarry Smith 
34447c6ae99SBarry Smith   Input Parameter:
345aa219208SBarry Smith + da - The DMDA
34647c6ae99SBarry 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
34747c6ae99SBarry 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
34847c6ae99SBarry 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.
34947c6ae99SBarry Smith 
35047c6ae99SBarry Smith   Level: intermediate
35147c6ae99SBarry Smith 
352*e3f3e4b6SBarry Smith   Note: these numbers are NOT multiplied by the number of dof per node.
353*e3f3e4b6SBarry Smith 
35447c6ae99SBarry Smith .keywords:  distributed array
35596e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
35647c6ae99SBarry Smith @*/
3577087cfbeSBarry Smith PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
35847c6ae99SBarry Smith {
35947c6ae99SBarry Smith   PetscErrorCode ierr;
36047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
36147c6ae99SBarry Smith 
36247c6ae99SBarry Smith   PetscFunctionBegin;
36347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
364cb630486SJed Brown   if (da->setupcalled) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
36547c6ae99SBarry Smith   if (lx) {
36647c6ae99SBarry Smith     if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
367aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
36847c6ae99SBarry Smith     if (!dd->lx) {
36947c6ae99SBarry Smith       ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
37047c6ae99SBarry Smith     }
37147c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
37247c6ae99SBarry Smith   }
37347c6ae99SBarry Smith   if (ly) {
37447c6ae99SBarry Smith     if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
375aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
37647c6ae99SBarry Smith     if (!dd->ly) {
37747c6ae99SBarry Smith       ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
37847c6ae99SBarry Smith     }
37947c6ae99SBarry Smith     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
38047c6ae99SBarry Smith   }
38147c6ae99SBarry Smith   if (lz) {
38247c6ae99SBarry Smith     if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
383aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
38447c6ae99SBarry Smith     if (!dd->lz) {
38547c6ae99SBarry Smith       ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
38647c6ae99SBarry Smith     }
38747c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
38847c6ae99SBarry Smith   }
38947c6ae99SBarry Smith   PetscFunctionReturn(0);
39047c6ae99SBarry Smith }
39147c6ae99SBarry Smith 
39247c6ae99SBarry Smith #undef __FUNCT__
393aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType"
39447c6ae99SBarry Smith /*@
395aa219208SBarry Smith        DMDASetInterpolationType - Sets the type of interpolation that will be
396e727c939SJed Brown           returned by DMCreateInterpolation()
39747c6ae99SBarry Smith 
398aa219208SBarry Smith    Logically Collective on DMDA
39947c6ae99SBarry Smith 
40047c6ae99SBarry Smith    Input Parameter:
40147c6ae99SBarry Smith +  da - initial distributed array
402aa219208SBarry Smith .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
40347c6ae99SBarry Smith 
40447c6ae99SBarry Smith    Level: intermediate
40547c6ae99SBarry Smith 
406e727c939SJed Brown    Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation()
40747c6ae99SBarry Smith 
40847c6ae99SBarry Smith .keywords:  distributed array, interpolation
40947c6ae99SBarry Smith 
41096e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
41147c6ae99SBarry Smith @*/
4127087cfbeSBarry Smith PetscErrorCode  DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
41347c6ae99SBarry Smith {
41447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
41547c6ae99SBarry Smith 
41647c6ae99SBarry Smith   PetscFunctionBegin;
41747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
41847c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,ctype,2);
41947c6ae99SBarry Smith   dd->interptype = ctype;
42047c6ae99SBarry Smith   PetscFunctionReturn(0);
42147c6ae99SBarry Smith }
42247c6ae99SBarry Smith 
4232dde6fd4SLisandro Dalcin #undef __FUNCT__
4242dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetInterpolationType"
4252dde6fd4SLisandro Dalcin /*@
4262dde6fd4SLisandro Dalcin        DMDAGetInterpolationType - Gets the type of interpolation that will be
427e727c939SJed Brown           used by DMCreateInterpolation()
4282dde6fd4SLisandro Dalcin 
4292dde6fd4SLisandro Dalcin    Not Collective
4302dde6fd4SLisandro Dalcin 
4312dde6fd4SLisandro Dalcin    Input Parameter:
4322dde6fd4SLisandro Dalcin .  da - distributed array
4332dde6fd4SLisandro Dalcin 
4342dde6fd4SLisandro Dalcin    Output Parameter:
4352dde6fd4SLisandro Dalcin .  ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms)
4362dde6fd4SLisandro Dalcin 
4372dde6fd4SLisandro Dalcin    Level: intermediate
4382dde6fd4SLisandro Dalcin 
4392dde6fd4SLisandro Dalcin .keywords:  distributed array, interpolation
4402dde6fd4SLisandro Dalcin 
441e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation()
4422dde6fd4SLisandro Dalcin @*/
4432dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype)
4442dde6fd4SLisandro Dalcin {
4452dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA*)da->data;
4462dde6fd4SLisandro Dalcin 
4472dde6fd4SLisandro Dalcin   PetscFunctionBegin;
4482dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
4492dde6fd4SLisandro Dalcin   PetscValidPointer(ctype,2);
4502dde6fd4SLisandro Dalcin   *ctype = dd->interptype;
4512dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
4522dde6fd4SLisandro Dalcin }
45347c6ae99SBarry Smith 
45447c6ae99SBarry Smith #undef __FUNCT__
455aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors"
45647c6ae99SBarry Smith /*@C
457aa219208SBarry Smith       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
45847c6ae99SBarry Smith         processes neighbors.
45947c6ae99SBarry Smith 
46047c6ae99SBarry Smith     Not Collective
46147c6ae99SBarry Smith 
46247c6ae99SBarry Smith    Input Parameter:
463aa219208SBarry Smith .     da - the DMDA object
46447c6ae99SBarry Smith 
46547c6ae99SBarry Smith    Output Parameters:
46647c6ae99SBarry Smith .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
46747c6ae99SBarry Smith               this process itself is in the list
46847c6ae99SBarry Smith 
46947c6ae99SBarry Smith    Notes: In 2d the array is of length 9, in 3d of length 27
47047c6ae99SBarry Smith           Not supported in 1d
471aa219208SBarry Smith           Do not free the array, it is freed when the DMDA is destroyed.
47247c6ae99SBarry Smith 
47347c6ae99SBarry Smith    Fortran Notes: In fortran you must pass in an array of the appropriate length.
47447c6ae99SBarry Smith 
47547c6ae99SBarry Smith    Level: intermediate
47647c6ae99SBarry Smith 
47747c6ae99SBarry Smith @*/
4787087cfbeSBarry Smith PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
47947c6ae99SBarry Smith {
48047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
48147c6ae99SBarry Smith   PetscFunctionBegin;
48247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
48347c6ae99SBarry Smith   *ranks = dd->neighbors;
48447c6ae99SBarry Smith   PetscFunctionReturn(0);
48547c6ae99SBarry Smith }
48647c6ae99SBarry Smith 
48747c6ae99SBarry Smith #undef __FUNCT__
488aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges"
48947c6ae99SBarry Smith /*@C
490aa219208SBarry Smith       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
49147c6ae99SBarry Smith 
49247c6ae99SBarry Smith     Not Collective
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith    Input Parameter:
495aa219208SBarry Smith .     da - the DMDA object
49647c6ae99SBarry Smith 
49747c6ae99SBarry Smith    Output Parameter:
49847c6ae99SBarry Smith +     lx - ownership along x direction (optional)
49947c6ae99SBarry Smith .     ly - ownership along y direction (optional)
50047c6ae99SBarry Smith -     lz - ownership along z direction (optional)
50147c6ae99SBarry Smith 
50247c6ae99SBarry Smith    Level: intermediate
50347c6ae99SBarry Smith 
504aa219208SBarry Smith     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
507aa219208SBarry Smith     eighth arguments from DMDAGetInfo()
50847c6ae99SBarry Smith 
50947c6ae99SBarry Smith      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
510aa219208SBarry Smith     DMDA they came from still exists (has not been destroyed).
51147c6ae99SBarry Smith 
512*e3f3e4b6SBarry Smith     These numbers are NOT multiplied by the number of dof per node.
513*e3f3e4b6SBarry Smith 
514aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
51547c6ae99SBarry Smith @*/
5167087cfbeSBarry Smith PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
51747c6ae99SBarry Smith {
51847c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
51947c6ae99SBarry Smith 
52047c6ae99SBarry Smith   PetscFunctionBegin;
52147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
52247c6ae99SBarry Smith   if (lx) *lx = dd->lx;
52347c6ae99SBarry Smith   if (ly) *ly = dd->ly;
52447c6ae99SBarry Smith   if (lz) *lz = dd->lz;
52547c6ae99SBarry Smith   PetscFunctionReturn(0);
52647c6ae99SBarry Smith }
52747c6ae99SBarry Smith 
52847c6ae99SBarry Smith #undef __FUNCT__
529aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor"
53047c6ae99SBarry Smith /*@
531aa219208SBarry Smith      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
53247c6ae99SBarry Smith 
533aa219208SBarry Smith     Logically Collective on DMDA
53447c6ae99SBarry Smith 
53547c6ae99SBarry Smith   Input Parameters:
536aa219208SBarry Smith +    da - the DMDA object
53747c6ae99SBarry Smith .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
53847c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
53947c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
54047c6ae99SBarry Smith 
54147c6ae99SBarry Smith   Options Database:
54247c6ae99SBarry Smith +  -da_refine_x - refinement ratio in x direction
54347c6ae99SBarry Smith .  -da_refine_y - refinement ratio in y direction
54447c6ae99SBarry Smith -  -da_refine_z - refinement ratio in z direction
54547c6ae99SBarry Smith 
54647c6ae99SBarry Smith   Level: intermediate
54747c6ae99SBarry Smith 
54847c6ae99SBarry Smith     Notes: Pass PETSC_IGNORE to leave a value unchanged
54947c6ae99SBarry Smith 
550aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor()
55147c6ae99SBarry Smith @*/
5527087cfbeSBarry Smith PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
55347c6ae99SBarry Smith {
55447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
55547c6ae99SBarry Smith 
55647c6ae99SBarry Smith   PetscFunctionBegin;
55747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
55847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_x,2);
55947c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_y,3);
56047c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_z,4);
56147c6ae99SBarry Smith 
56247c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
56347c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
56447c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
56547c6ae99SBarry Smith   PetscFunctionReturn(0);
56647c6ae99SBarry Smith }
56747c6ae99SBarry Smith 
56847c6ae99SBarry Smith #undef __FUNCT__
569aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor"
57047c6ae99SBarry Smith /*@C
571aa219208SBarry Smith      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
57247c6ae99SBarry Smith 
57347c6ae99SBarry Smith     Not Collective
57447c6ae99SBarry Smith 
57547c6ae99SBarry Smith   Input Parameter:
576aa219208SBarry Smith .    da - the DMDA object
57747c6ae99SBarry Smith 
57847c6ae99SBarry Smith   Output Parameters:
57947c6ae99SBarry Smith +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
58047c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
58147c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
58247c6ae99SBarry Smith 
58347c6ae99SBarry Smith   Level: intermediate
58447c6ae99SBarry Smith 
58547c6ae99SBarry Smith     Notes: Pass PETSC_NULL for values you do not need
58647c6ae99SBarry Smith 
587aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor()
58847c6ae99SBarry Smith @*/
5897087cfbeSBarry Smith PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
59047c6ae99SBarry Smith {
59147c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
59247c6ae99SBarry Smith 
59347c6ae99SBarry Smith   PetscFunctionBegin;
59447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
59547c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
59647c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
59747c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
59847c6ae99SBarry Smith   PetscFunctionReturn(0);
59947c6ae99SBarry Smith }
60047c6ae99SBarry Smith 
60147c6ae99SBarry Smith #undef __FUNCT__
602aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix"
60347c6ae99SBarry Smith /*@C
604aa219208SBarry Smith      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
60547c6ae99SBarry Smith 
606aa219208SBarry Smith     Logically Collective on DMDA
60747c6ae99SBarry Smith 
60847c6ae99SBarry Smith   Input Parameters:
609aa219208SBarry Smith +    da - the DMDA object
610aa219208SBarry Smith -    f - the function that allocates the matrix for that specific DMDA
61147c6ae99SBarry Smith 
61247c6ae99SBarry Smith   Level: developer
61347c6ae99SBarry Smith 
614aa219208SBarry Smith    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
61547c6ae99SBarry Smith        the diagonal and off-diagonal blocks of the matrix
61647c6ae99SBarry Smith 
617950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills()
61847c6ae99SBarry Smith @*/
61919fd82e9SBarry Smith PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, MatType,Mat*))
62047c6ae99SBarry Smith {
62147c6ae99SBarry Smith   PetscFunctionBegin;
62247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
62325296bd5SBarry Smith   da->ops->creatematrix = f;
62447c6ae99SBarry Smith   PetscFunctionReturn(0);
62547c6ae99SBarry Smith }
62647c6ae99SBarry Smith 
62747c6ae99SBarry Smith #undef __FUNCT__
62894013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges"
62947c6ae99SBarry Smith /*
63047c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
63147c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
63247c6ae99SBarry Smith 
63347c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
63447c6ae99SBarry Smith */
63594013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
63647c6ae99SBarry Smith {
63747c6ae99SBarry Smith   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
63847c6ae99SBarry Smith   PetscErrorCode ierr;
63947c6ae99SBarry Smith 
64047c6ae99SBarry Smith   PetscFunctionBegin;
64147c6ae99SBarry Smith   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
64247c6ae99SBarry Smith   if (ratio == 1) {
64347c6ae99SBarry Smith     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
64447c6ae99SBarry Smith     PetscFunctionReturn(0);
64547c6ae99SBarry Smith   }
64647c6ae99SBarry Smith   for (i=0; i<m; i++) totalc += lc[i];
64747c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
64847c6ae99SBarry Smith   for (i=0; i<m; i++) {
64947c6ae99SBarry Smith     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
65047c6ae99SBarry Smith     if (i == m-1) lf[i] = want;
65147c6ae99SBarry Smith     else {
6527aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
6537aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
6547aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
6557aca7175SJed Brown       while ((startf+want)/ratio < nextc - stencil_width) want++;
6567aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
6577aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
6587aca7175SJed Brown       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
6597aca7175SJed Brown       /* Make sure all constraints are satisfied */
66030729d88SBarry Smith       if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width)
66130729d88SBarry Smith           || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range");
66247c6ae99SBarry Smith     }
66347c6ae99SBarry Smith     lf[i] = want;
66447c6ae99SBarry Smith     startc += lc[i];
66547c6ae99SBarry Smith     startf += lf[i];
66647c6ae99SBarry Smith     remaining -= lf[i];
66747c6ae99SBarry Smith   }
66847c6ae99SBarry Smith   PetscFunctionReturn(0);
66947c6ae99SBarry Smith }
67047c6ae99SBarry Smith 
67147c6ae99SBarry Smith #undef __FUNCT__
6722be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges"
6732be375d4SJed Brown /*
6742be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
6752be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
6762be375d4SJed Brown 
6772be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
6782be375d4SJed Brown */
6792be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
6802be375d4SJed Brown {
6812be375d4SJed Brown   PetscInt       i,totalf,remaining,startc,startf;
6822be375d4SJed Brown   PetscErrorCode ierr;
6832be375d4SJed Brown 
6842be375d4SJed Brown   PetscFunctionBegin;
6852be375d4SJed Brown   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
6862be375d4SJed Brown   if (ratio == 1) {
6872be375d4SJed Brown     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
6882be375d4SJed Brown     PetscFunctionReturn(0);
6892be375d4SJed Brown   }
6902be375d4SJed Brown   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
6912be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
6922be375d4SJed Brown   for (i=0,startc=0,startf=0; i<m; i++) {
6932be375d4SJed Brown     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
6942be375d4SJed Brown     if (i == m-1) lc[i] = want;
6952be375d4SJed Brown     else {
6962be375d4SJed Brown       const PetscInt nextf = startf+lf[i];
6972be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
6982be375d4SJed Brown        * node is within one stencil width. */
6992be375d4SJed Brown       while (nextf/ratio < startc+want-stencil_width) want--;
7002be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
7012be375d4SJed Brown        * fine node is within one stencil width. */
7022be375d4SJed Brown       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
7032be375d4SJed Brown       if (want < 0 || want > remaining
70430729d88SBarry Smith           || (nextf/ratio < startc+want-stencil_width) || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range");
7052be375d4SJed Brown     }
7062be375d4SJed Brown     lc[i] = want;
7072be375d4SJed Brown     startc += lc[i];
7082be375d4SJed Brown     startf += lf[i];
7092be375d4SJed Brown     remaining -= lc[i];
7102be375d4SJed Brown   }
7112be375d4SJed Brown   PetscFunctionReturn(0);
7122be375d4SJed Brown }
7132be375d4SJed Brown 
7142be375d4SJed Brown #undef __FUNCT__
71594013140SBarry Smith #define __FUNCT__ "DMRefine_DA"
7167087cfbeSBarry Smith PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
71747c6ae99SBarry Smith {
71847c6ae99SBarry Smith   PetscErrorCode ierr;
719f3141302SJed Brown   PetscInt       M,N,P,i;
7209a42bb27SBarry Smith   DM             da2;
72147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
72247c6ae99SBarry Smith 
72347c6ae99SBarry Smith   PetscFunctionBegin;
72447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
72547c6ae99SBarry Smith   PetscValidPointer(daref,3);
72647c6ae99SBarry Smith 
7271321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
72847c6ae99SBarry Smith     M = dd->refine_x*dd->M;
72947c6ae99SBarry Smith   } else {
73047c6ae99SBarry Smith     M = 1 + dd->refine_x*(dd->M - 1);
73147c6ae99SBarry Smith   }
7321321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
7331860e6e9SBarry Smith     if (dd->dim > 1) {
73447c6ae99SBarry Smith       N = dd->refine_y*dd->N;
73547c6ae99SBarry Smith     } else {
7361860e6e9SBarry Smith       N = 1;
7371860e6e9SBarry Smith     }
7381860e6e9SBarry Smith   } else {
73947c6ae99SBarry Smith     N = 1 + dd->refine_y*(dd->N - 1);
74047c6ae99SBarry Smith   }
7411321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
7421860e6e9SBarry Smith     if (dd->dim > 2) {
74347c6ae99SBarry Smith       P = dd->refine_z*dd->P;
74447c6ae99SBarry Smith     } else {
7451860e6e9SBarry Smith       P = 1;
7461860e6e9SBarry Smith     }
7471860e6e9SBarry Smith   } else {
74847c6ae99SBarry Smith     P = 1 + dd->refine_z*(dd->P - 1);
74947c6ae99SBarry Smith   }
750397b6216SJed Brown   ierr = DMDACreate(((PetscObject)da)->comm,&da2);CHKERRQ(ierr);
751397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
752397b6216SJed Brown   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
753397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
754397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
755397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
756397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
757397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
758397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
75947c6ae99SBarry Smith   if (dd->dim == 3) {
76047c6ae99SBarry Smith     PetscInt *lx,*ly,*lz;
76147c6ae99SBarry Smith     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
7621321219cSEthan 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);
7631321219cSEthan 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);
7641321219cSEthan 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);
765397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
76647c6ae99SBarry Smith     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
76747c6ae99SBarry Smith   } else if (dd->dim == 2) {
76847c6ae99SBarry Smith     PetscInt *lx,*ly;
76947c6ae99SBarry Smith     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
7701321219cSEthan 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);
7711321219cSEthan 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);
772397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,PETSC_NULL);CHKERRQ(ierr);
77347c6ae99SBarry Smith     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
77447c6ae99SBarry Smith   } else if (dd->dim == 1) {
77547c6ae99SBarry Smith     PetscInt *lx;
77647c6ae99SBarry Smith     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
7771321219cSEthan 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);
778397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
77947c6ae99SBarry Smith     ierr = PetscFree(lx);CHKERRQ(ierr);
78047c6ae99SBarry Smith   }
78147c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
78247c6ae99SBarry Smith 
78347c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
78425296bd5SBarry Smith   da2->ops->creatematrix        = da->ops->creatematrix;
78525296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
78647c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
78747c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
78847c6ae99SBarry Smith 
78947c6ae99SBarry Smith   /* copy fill information if given */
79047c6ae99SBarry Smith   if (dd->dfill) {
79147c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
79247c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
79347c6ae99SBarry Smith   }
79447c6ae99SBarry Smith   if (dd->ofill) {
79547c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
79647c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
79747c6ae99SBarry Smith   }
79847c6ae99SBarry Smith   /* copy the refine information */
799397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
800397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
801397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
80247c6ae99SBarry Smith 
80347c6ae99SBarry Smith   /* copy vector type information */
80447c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
80519fd82e9SBarry Smith   ierr = PetscStrallocpy(da->vectype,(char**)&da2->vectype);CHKERRQ(ierr);
806ddcf8b74SDave May 
807efd51863SBarry Smith   dd2->lf = dd->lf;
808efd51863SBarry Smith   dd2->lj = dd->lj;
809efd51863SBarry Smith 
8106e87535bSJed Brown   da2->leveldown = da->leveldown;
8116e87535bSJed Brown   da2->levelup   = da->levelup + 1;
8126e87535bSJed Brown   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
8136e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
814ca266f36SBarry Smith   ierr = DMViewFromOptions(da2,"-dm_view");CHKERRQ(ierr);
8156e87535bSJed Brown 
816ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
8176636e97aSMatthew G Knepley   if (da->coordinates) {
818ddcf8b74SDave May     DM  cdaf,cdac;
819ddcf8b74SDave May     Vec coordsc,coordsf;
820ddcf8b74SDave May     Mat II;
821ddcf8b74SDave May 
8226636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr);
8236636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr);
8246636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr);
825b61d3410SDave May     /* force creation of the coordinate vector */
826b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
8276636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
828e727c939SJed Brown     ierr = DMCreateInterpolation(cdac,cdaf,&II,PETSC_NULL);CHKERRQ(ierr);
829ddcf8b74SDave May     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
830fcfd50ebSBarry Smith     ierr = MatDestroy(&II);CHKERRQ(ierr);
831ddcf8b74SDave May   }
832397b6216SJed Brown 
833f3141302SJed Brown   for (i=0; i<da->bs; i++) {
834f3141302SJed Brown     const char *fieldname;
835f3141302SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
836f3141302SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
837f3141302SJed Brown   }
838397b6216SJed Brown 
83947c6ae99SBarry Smith   *daref = da2;
84047c6ae99SBarry Smith   PetscFunctionReturn(0);
84147c6ae99SBarry Smith }
84247c6ae99SBarry Smith 
843397b6216SJed Brown 
84447c6ae99SBarry Smith #undef __FUNCT__
84594013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA"
8467087cfbeSBarry Smith PetscErrorCode  DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
84747c6ae99SBarry Smith {
84847c6ae99SBarry Smith   PetscErrorCode ierr;
849397b6216SJed Brown   PetscInt       M,N,P,i;
8509a42bb27SBarry Smith   DM             da2;
85147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
85247c6ae99SBarry Smith 
85347c6ae99SBarry Smith   PetscFunctionBegin;
85447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
85547c6ae99SBarry Smith   PetscValidPointer(daref,3);
85647c6ae99SBarry Smith 
8571321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
858397b6216SJed Brown     M = dd->M / dd->coarsen_x;
85947c6ae99SBarry Smith   } else {
860397b6216SJed Brown     M = 1 + (dd->M - 1) / dd->coarsen_x;
86147c6ae99SBarry Smith   }
8621321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
8631860e6e9SBarry Smith     if (dd->dim > 1) {
864397b6216SJed Brown       N = dd->N / dd->coarsen_y;
86547c6ae99SBarry Smith     } else {
8661860e6e9SBarry Smith       N = 1;
8671860e6e9SBarry Smith     }
8681860e6e9SBarry Smith   } else {
869397b6216SJed Brown     N = 1 + (dd->N - 1) / dd->coarsen_y;
87047c6ae99SBarry Smith   }
8711321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
8721860e6e9SBarry Smith     if (dd->dim > 2) {
873397b6216SJed Brown       P = dd->P / dd->coarsen_z;
87447c6ae99SBarry Smith     } else {
8751860e6e9SBarry Smith       P = 1;
8761860e6e9SBarry Smith     }
8771860e6e9SBarry Smith   } else {
878397b6216SJed Brown     P = 1 + (dd->P - 1) / dd->coarsen_z;
87947c6ae99SBarry Smith   }
880397b6216SJed Brown   ierr = DMDACreate(((PetscObject)da)->comm,&da2);CHKERRQ(ierr);
881397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
882397b6216SJed Brown   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
883397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
884397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
885397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
886397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
887397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
888397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
88947c6ae99SBarry Smith   if (dd->dim == 3) {
8902be375d4SJed Brown     PetscInt *lx,*ly,*lz;
8912be375d4SJed Brown     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
892397b6216SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
893397b6216SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
894397b6216SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
895397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
8962be375d4SJed Brown     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
89747c6ae99SBarry Smith   } else if (dd->dim == 2) {
8982be375d4SJed Brown     PetscInt *lx,*ly;
8992be375d4SJed Brown     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
900397b6216SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
901397b6216SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
902397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,PETSC_NULL);CHKERRQ(ierr);
9032be375d4SJed Brown     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
90447c6ae99SBarry Smith   } else if (dd->dim == 1) {
9052be375d4SJed Brown     PetscInt *lx;
9062be375d4SJed Brown     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
907397b6216SJed Brown     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
908397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
9092be375d4SJed Brown     ierr = PetscFree(lx);CHKERRQ(ierr);
91047c6ae99SBarry Smith   }
91147c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
91247c6ae99SBarry Smith 
9134dcab191SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
91425296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */
91525296bd5SBarry Smith   da2->ops->creatematrix        = da->ops->creatematrix;
91647c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
91747c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
91847c6ae99SBarry Smith 
91947c6ae99SBarry Smith   /* copy fill information if given */
92047c6ae99SBarry Smith   if (dd->dfill) {
92147c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
92247c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
92347c6ae99SBarry Smith   }
92447c6ae99SBarry Smith   if (dd->ofill) {
92547c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
92647c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
92747c6ae99SBarry Smith   }
92847c6ae99SBarry Smith   /* copy the refine information */
929397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
930397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
931397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
93247c6ae99SBarry Smith 
93347c6ae99SBarry Smith   /* copy vector type information */
93447c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
93519fd82e9SBarry Smith   ierr = PetscStrallocpy(da->vectype,(char**)&da2->vectype);CHKERRQ(ierr);
93647c6ae99SBarry Smith 
937644e2e5bSBarry Smith   dd2->lf = dd->lf;
938644e2e5bSBarry Smith   dd2->lj = dd->lj;
939644e2e5bSBarry Smith 
9406e87535bSJed Brown   da2->leveldown = da->leveldown + 1;
9416e87535bSJed Brown   da2->levelup   = da->levelup;
9426e87535bSJed Brown   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
9436e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
944ca266f36SBarry Smith   ierr = DMViewFromOptions(da2,"-dm_view");CHKERRQ(ierr);
9456e87535bSJed Brown 
946ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
9476636e97aSMatthew G Knepley   if (da->coordinates) {
948ddcf8b74SDave May     DM         cdaf,cdac;
949ddcf8b74SDave May     Vec        coordsc,coordsf;
950ddcf8b74SDave May     VecScatter inject;
951ddcf8b74SDave May 
9526636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr);
9536636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr);
9546636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr);
955b61d3410SDave May     /* force creation of the coordinate vector */
956b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
9576636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
958ddcf8b74SDave May 
959e727c939SJed Brown     ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
960ddcf8b74SDave May     ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
961ddcf8b74SDave May     ierr = VecScatterEnd(inject  ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
962fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&inject);CHKERRQ(ierr);
963ddcf8b74SDave May   }
964f98405f7SJed Brown 
965f98405f7SJed Brown   for (i=0; i<da->bs; i++) {
966f98405f7SJed Brown     const char *fieldname;
967f98405f7SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
968f98405f7SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
969f98405f7SJed Brown   }
970f98405f7SJed Brown 
97147c6ae99SBarry Smith   *daref = da2;
97247c6ae99SBarry Smith   PetscFunctionReturn(0);
97347c6ae99SBarry Smith }
97447c6ae99SBarry Smith 
97547c6ae99SBarry Smith #undef __FUNCT__
97694013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA"
9777087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
97847c6ae99SBarry Smith {
97947c6ae99SBarry Smith   PetscErrorCode ierr;
98047c6ae99SBarry Smith   PetscInt       i,n,*refx,*refy,*refz;
98147c6ae99SBarry Smith 
98247c6ae99SBarry Smith   PetscFunctionBegin;
98347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
98447c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
98547c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
98647c6ae99SBarry Smith   PetscValidPointer(daf,3);
98747c6ae99SBarry Smith 
988aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
98947c6ae99SBarry Smith   ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr);
99047c6ae99SBarry Smith   for (i=0; i<nlevels; i++) {
991aa219208SBarry Smith     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
99247c6ae99SBarry Smith   }
99347c6ae99SBarry Smith   n = nlevels;
99447c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
99547c6ae99SBarry Smith   n = nlevels;
99647c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
99747c6ae99SBarry Smith   n = nlevels;
99847c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
99947c6ae99SBarry Smith 
1000aa219208SBarry Smith   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
100194013140SBarry Smith   ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr);
100247c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
1003aa219208SBarry Smith     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
100494013140SBarry Smith     ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr);
100547c6ae99SBarry Smith   }
100647c6ae99SBarry Smith   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
100747c6ae99SBarry Smith   PetscFunctionReturn(0);
100847c6ae99SBarry Smith }
100947c6ae99SBarry Smith 
101047c6ae99SBarry Smith #undef __FUNCT__
101194013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA"
10127087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
101347c6ae99SBarry Smith {
101447c6ae99SBarry Smith   PetscErrorCode ierr;
101547c6ae99SBarry Smith   PetscInt       i;
101647c6ae99SBarry Smith 
101747c6ae99SBarry Smith   PetscFunctionBegin;
101847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
101947c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
102047c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
102147c6ae99SBarry Smith   PetscValidPointer(dac,3);
102294013140SBarry Smith   ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr);
102347c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
102494013140SBarry Smith     ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr);
102547c6ae99SBarry Smith   }
102647c6ae99SBarry Smith   PetscFunctionReturn(0);
102747c6ae99SBarry Smith }
1028