xref: /petsc/src/dm/impls/da/da.c (revision f31413023db3672df454ca19e1331d88063c4fa6)
1c6db04a5SJed Brown #include <private/daimpl.h>    /*I   "petscdmda.h"   I*/
247c6ae99SBarry Smith 
347c6ae99SBarry Smith 
447c6ae99SBarry Smith #undef __FUNCT__
5aa219208SBarry Smith #define __FUNCT__ "DMDASetDim"
647c6ae99SBarry Smith /*@
7aa219208SBarry Smith   DMDASetDim - Sets the dimension
847c6ae99SBarry Smith 
9aa219208SBarry Smith   Collective on DMDA
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith   Input Parameters:
12aa219208SBarry Smith + da - the DMDA
1347c6ae99SBarry Smith - dim - the dimension (or PETSC_DECIDE)
1447c6ae99SBarry Smith 
1547c6ae99SBarry Smith   Level: intermediate
1647c6ae99SBarry Smith 
17aa219208SBarry Smith .seealso: DaGetDim(), DMDASetSizes()
1847c6ae99SBarry Smith @*/
197087cfbeSBarry Smith PetscErrorCode  DMDASetDim(DM da, PetscInt dim)
2047c6ae99SBarry Smith {
2147c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
2247c6ae99SBarry Smith 
2347c6ae99SBarry Smith   PetscFunctionBegin;
2447c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
25aa219208SBarry Smith   if (dd->dim > 0 && dim != dd->dim) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change DMDA dim from %D after it was set to %D",dd->dim,dim);
2647c6ae99SBarry Smith   dd->dim = dim;
2747c6ae99SBarry Smith   PetscFunctionReturn(0);
2847c6ae99SBarry Smith }
2947c6ae99SBarry Smith 
3047c6ae99SBarry Smith #undef __FUNCT__
31aa219208SBarry Smith #define __FUNCT__ "DMDASetSizes"
3247c6ae99SBarry Smith /*@
33aa219208SBarry Smith   DMDASetSizes - Sets the global sizes
3447c6ae99SBarry Smith 
35aa219208SBarry Smith   Logically Collective on DMDA
3647c6ae99SBarry Smith 
3747c6ae99SBarry Smith   Input Parameters:
38aa219208SBarry Smith + da - the DMDA
3947c6ae99SBarry Smith . M - the global X size (or PETSC_DECIDE)
4047c6ae99SBarry Smith . N - the global Y size (or PETSC_DECIDE)
4147c6ae99SBarry Smith - P - the global Z size (or PETSC_DECIDE)
4247c6ae99SBarry Smith 
4347c6ae99SBarry Smith   Level: intermediate
4447c6ae99SBarry Smith 
45aa219208SBarry Smith .seealso: DMDAGetSize(), PetscSplitOwnership()
4647c6ae99SBarry Smith @*/
477087cfbeSBarry Smith PetscErrorCode  DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P)
4847c6ae99SBarry Smith {
4947c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
5047c6ae99SBarry Smith 
5147c6ae99SBarry Smith   PetscFunctionBegin;
5247c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
5347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,M,2);
5447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,N,3);
5547c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,P,4);
5647c6ae99SBarry Smith 
5747c6ae99SBarry Smith   dd->M = M;
5847c6ae99SBarry Smith   dd->N = N;
5947c6ae99SBarry Smith   dd->P = P;
6047c6ae99SBarry Smith   PetscFunctionReturn(0);
6147c6ae99SBarry Smith }
6247c6ae99SBarry Smith 
6347c6ae99SBarry Smith #undef __FUNCT__
64aa219208SBarry Smith #define __FUNCT__ "DMDASetNumProcs"
6547c6ae99SBarry Smith /*@
66aa219208SBarry Smith   DMDASetNumProcs - Sets the number of processes in each dimension
6747c6ae99SBarry Smith 
68aa219208SBarry Smith   Logically Collective on DMDA
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith   Input Parameters:
71aa219208SBarry Smith + da - the DMDA
7247c6ae99SBarry Smith . m - the number of X procs (or PETSC_DECIDE)
7347c6ae99SBarry Smith . n - the number of Y procs (or PETSC_DECIDE)
7447c6ae99SBarry Smith - p - the number of Z procs (or PETSC_DECIDE)
7547c6ae99SBarry Smith 
7647c6ae99SBarry Smith   Level: intermediate
7747c6ae99SBarry Smith 
78aa219208SBarry Smith .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership()
7947c6ae99SBarry Smith @*/
807087cfbeSBarry Smith PetscErrorCode  DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p)
8147c6ae99SBarry Smith {
8247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
8347c6ae99SBarry Smith 
8447c6ae99SBarry Smith   PetscFunctionBegin;
8547c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
8647c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,m,2);
8747c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,n,3);
8847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,p,4);
8947c6ae99SBarry Smith   dd->m = m;
9047c6ae99SBarry Smith   dd->n = n;
9147c6ae99SBarry Smith   dd->p = p;
9247c6ae99SBarry Smith   PetscFunctionReturn(0);
9347c6ae99SBarry Smith }
9447c6ae99SBarry Smith 
9547c6ae99SBarry Smith #undef __FUNCT__
963e24977eSEthan Coon #define __FUNCT__ "DMDASetBoundaryType"
9747c6ae99SBarry Smith /*@
981321219cSEthan Coon   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries.
9947c6ae99SBarry Smith 
10047c6ae99SBarry Smith   Not collective
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith   Input Parameter:
103aa219208SBarry Smith + da    - The DMDA
1041321219cSEthan Coon - bx,by,bz - One of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC
10547c6ae99SBarry Smith 
10647c6ae99SBarry Smith   Level: intermediate
10747c6ae99SBarry Smith 
10847c6ae99SBarry Smith .keywords:  distributed array, periodicity
10996e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA, DMDABoundaryType
11047c6ae99SBarry Smith @*/
1111321219cSEthan Coon PetscErrorCode  DMDASetBoundaryType(DM da,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz)
11247c6ae99SBarry Smith {
11347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith   PetscFunctionBegin;
11647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1175a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,bx,2);
1185a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,by,3);
1195a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,bz,4);
1201321219cSEthan Coon   dd->bx = bx;
1211321219cSEthan Coon   dd->by = by;
1221321219cSEthan Coon   dd->bz = bz;
12347c6ae99SBarry Smith   PetscFunctionReturn(0);
12447c6ae99SBarry Smith }
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith #undef __FUNCT__
127aa219208SBarry Smith #define __FUNCT__ "DMDASetDof"
12847c6ae99SBarry Smith /*@
129aa219208SBarry Smith   DMDASetDof - Sets the number of degrees of freedom per vertex
13047c6ae99SBarry Smith 
13147c6ae99SBarry Smith   Not collective
13247c6ae99SBarry Smith 
13347c6ae99SBarry Smith   Input Parameter:
134aa219208SBarry Smith + da  - The DMDA
13547c6ae99SBarry Smith - dof - Number of degrees of freedom
13647c6ae99SBarry Smith 
13747c6ae99SBarry Smith   Level: intermediate
13847c6ae99SBarry Smith 
13947c6ae99SBarry Smith .keywords:  distributed array, degrees of freedom
14096e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
14147c6ae99SBarry Smith @*/
14254cfb0beSLisandro Dalcin PetscErrorCode  DMDASetDof(DM da, PetscInt dof)
14347c6ae99SBarry Smith {
14447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith   PetscFunctionBegin;
14747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
14854cfb0beSLisandro Dalcin   PetscValidLogicalCollectiveInt(da,dof,2);
14947c6ae99SBarry Smith   dd->w = dof;
1501411c6eeSJed Brown   da->bs = dof;
15147c6ae99SBarry Smith   PetscFunctionReturn(0);
15247c6ae99SBarry Smith }
15347c6ae99SBarry Smith 
15447c6ae99SBarry Smith #undef __FUNCT__
155aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType"
15647c6ae99SBarry Smith /*@
157aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
15847c6ae99SBarry Smith 
159aa219208SBarry Smith   Logically Collective on DMDA
16047c6ae99SBarry Smith 
16147c6ae99SBarry Smith   Input Parameter:
162aa219208SBarry Smith + da    - The DMDA
163aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
16447c6ae99SBarry Smith 
16547c6ae99SBarry Smith   Level: intermediate
16647c6ae99SBarry Smith 
16747c6ae99SBarry Smith .keywords:  distributed array, stencil
16896e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
16947c6ae99SBarry Smith @*/
1707087cfbeSBarry Smith PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
17147c6ae99SBarry Smith {
17247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
17347c6ae99SBarry Smith 
17447c6ae99SBarry Smith   PetscFunctionBegin;
17547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
17647c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,stype,2);
17747c6ae99SBarry Smith   dd->stencil_type = stype;
17847c6ae99SBarry Smith   PetscFunctionReturn(0);
17947c6ae99SBarry Smith }
18047c6ae99SBarry Smith 
18147c6ae99SBarry Smith #undef __FUNCT__
182aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth"
18347c6ae99SBarry Smith /*@
184aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
18547c6ae99SBarry Smith 
186aa219208SBarry Smith   Logically Collective on DMDA
18747c6ae99SBarry Smith 
18847c6ae99SBarry Smith   Input Parameter:
189aa219208SBarry Smith + da    - The DMDA
19047c6ae99SBarry Smith - width - The stencil width
19147c6ae99SBarry Smith 
19247c6ae99SBarry Smith   Level: intermediate
19347c6ae99SBarry Smith 
19447c6ae99SBarry Smith .keywords:  distributed array, stencil
19596e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
19647c6ae99SBarry Smith @*/
1977087cfbeSBarry Smith PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
19847c6ae99SBarry Smith {
19947c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
20047c6ae99SBarry Smith 
20147c6ae99SBarry Smith   PetscFunctionBegin;
20247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
20347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,width,2);
20447c6ae99SBarry Smith   dd->s = width;
20547c6ae99SBarry Smith   PetscFunctionReturn(0);
20647c6ae99SBarry Smith }
20747c6ae99SBarry Smith 
20847c6ae99SBarry Smith #undef __FUNCT__
209aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
210aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
21147c6ae99SBarry Smith {
21247c6ae99SBarry Smith   PetscInt i,sum;
21347c6ae99SBarry Smith 
21447c6ae99SBarry Smith   PetscFunctionBegin;
21547c6ae99SBarry Smith   if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
21647c6ae99SBarry Smith   for (i=sum=0; i<m; i++) sum += lx[i];
21747c6ae99SBarry Smith   if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
21847c6ae99SBarry Smith   PetscFunctionReturn(0);
21947c6ae99SBarry Smith }
22047c6ae99SBarry Smith 
22147c6ae99SBarry Smith #undef __FUNCT__
222aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges"
22347c6ae99SBarry Smith /*@
224aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
22547c6ae99SBarry Smith 
226aa219208SBarry Smith   Logically Collective on DMDA
22747c6ae99SBarry Smith 
22847c6ae99SBarry Smith   Input Parameter:
229aa219208SBarry Smith + da - The DMDA
23047c6ae99SBarry Smith . lx - array containing number of nodes in the X direction on each process, or PETSC_NULL. If non-null, must be of length da->m
23147c6ae99SBarry Smith . ly - array containing number of nodes in the Y direction on each process, or PETSC_NULL. If non-null, must be of length da->n
23247c6ae99SBarry Smith - lz - array containing number of nodes in the Z direction on each process, or PETSC_NULL. If non-null, must be of length da->p.
23347c6ae99SBarry Smith 
23447c6ae99SBarry Smith   Level: intermediate
23547c6ae99SBarry Smith 
23647c6ae99SBarry Smith .keywords:  distributed array
23796e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
23847c6ae99SBarry Smith @*/
2397087cfbeSBarry Smith PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
24047c6ae99SBarry Smith {
24147c6ae99SBarry Smith   PetscErrorCode ierr;
24247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
24347c6ae99SBarry Smith 
24447c6ae99SBarry Smith   PetscFunctionBegin;
24547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
24647c6ae99SBarry Smith   if (lx) {
24747c6ae99SBarry Smith     if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
248aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
24947c6ae99SBarry Smith     if (!dd->lx) {
25047c6ae99SBarry Smith       ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
25147c6ae99SBarry Smith     }
25247c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
25347c6ae99SBarry Smith   }
25447c6ae99SBarry Smith   if (ly) {
25547c6ae99SBarry Smith     if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
256aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
25747c6ae99SBarry Smith     if (!dd->ly) {
25847c6ae99SBarry Smith       ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
25947c6ae99SBarry Smith     }
26047c6ae99SBarry Smith     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
26147c6ae99SBarry Smith   }
26247c6ae99SBarry Smith   if (lz) {
26347c6ae99SBarry Smith     if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
264aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
26547c6ae99SBarry Smith     if (!dd->lz) {
26647c6ae99SBarry Smith       ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
26747c6ae99SBarry Smith     }
26847c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
26947c6ae99SBarry Smith   }
27047c6ae99SBarry Smith   PetscFunctionReturn(0);
27147c6ae99SBarry Smith }
27247c6ae99SBarry Smith 
27347c6ae99SBarry Smith #undef __FUNCT__
274aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType"
27547c6ae99SBarry Smith /*@
276aa219208SBarry Smith        DMDASetInterpolationType - Sets the type of interpolation that will be
27794013140SBarry Smith           returned by DMGetInterpolation()
27847c6ae99SBarry Smith 
279aa219208SBarry Smith    Logically Collective on DMDA
28047c6ae99SBarry Smith 
28147c6ae99SBarry Smith    Input Parameter:
28247c6ae99SBarry Smith +  da - initial distributed array
283aa219208SBarry Smith .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
28447c6ae99SBarry Smith 
28547c6ae99SBarry Smith    Level: intermediate
28647c6ae99SBarry Smith 
287aa219208SBarry Smith    Notes: you should call this on the coarser of the two DMDAs you pass to DMGetInterpolation()
28847c6ae99SBarry Smith 
28947c6ae99SBarry Smith .keywords:  distributed array, interpolation
29047c6ae99SBarry Smith 
29196e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
29247c6ae99SBarry Smith @*/
2937087cfbeSBarry Smith PetscErrorCode  DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
29447c6ae99SBarry Smith {
29547c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
29647c6ae99SBarry Smith 
29747c6ae99SBarry Smith   PetscFunctionBegin;
29847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
29947c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,ctype,2);
30047c6ae99SBarry Smith   dd->interptype = ctype;
30147c6ae99SBarry Smith   PetscFunctionReturn(0);
30247c6ae99SBarry Smith }
30347c6ae99SBarry Smith 
3042dde6fd4SLisandro Dalcin #undef __FUNCT__
3052dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetInterpolationType"
3062dde6fd4SLisandro Dalcin /*@
3072dde6fd4SLisandro Dalcin        DMDAGetInterpolationType - Gets the type of interpolation that will be
3082dde6fd4SLisandro Dalcin           used by DMGetInterpolation()
3092dde6fd4SLisandro Dalcin 
3102dde6fd4SLisandro Dalcin    Not Collective
3112dde6fd4SLisandro Dalcin 
3122dde6fd4SLisandro Dalcin    Input Parameter:
3132dde6fd4SLisandro Dalcin .  da - distributed array
3142dde6fd4SLisandro Dalcin 
3152dde6fd4SLisandro Dalcin    Output Parameter:
3162dde6fd4SLisandro Dalcin .  ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms)
3172dde6fd4SLisandro Dalcin 
3182dde6fd4SLisandro Dalcin    Level: intermediate
3192dde6fd4SLisandro Dalcin 
3202dde6fd4SLisandro Dalcin .keywords:  distributed array, interpolation
3212dde6fd4SLisandro Dalcin 
3222dde6fd4SLisandro Dalcin .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMGetInterpolation()
3232dde6fd4SLisandro Dalcin @*/
3242dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype)
3252dde6fd4SLisandro Dalcin {
3262dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA*)da->data;
3272dde6fd4SLisandro Dalcin 
3282dde6fd4SLisandro Dalcin   PetscFunctionBegin;
3292dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
3302dde6fd4SLisandro Dalcin   PetscValidPointer(ctype,2);
3312dde6fd4SLisandro Dalcin   *ctype = dd->interptype;
3322dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
3332dde6fd4SLisandro Dalcin }
33447c6ae99SBarry Smith 
33547c6ae99SBarry Smith #undef __FUNCT__
336aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors"
33747c6ae99SBarry Smith /*@C
338aa219208SBarry Smith       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
33947c6ae99SBarry Smith         processes neighbors.
34047c6ae99SBarry Smith 
34147c6ae99SBarry Smith     Not Collective
34247c6ae99SBarry Smith 
34347c6ae99SBarry Smith    Input Parameter:
344aa219208SBarry Smith .     da - the DMDA object
34547c6ae99SBarry Smith 
34647c6ae99SBarry Smith    Output Parameters:
34747c6ae99SBarry Smith .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
34847c6ae99SBarry Smith               this process itself is in the list
34947c6ae99SBarry Smith 
35047c6ae99SBarry Smith    Notes: In 2d the array is of length 9, in 3d of length 27
35147c6ae99SBarry Smith           Not supported in 1d
352aa219208SBarry Smith           Do not free the array, it is freed when the DMDA is destroyed.
35347c6ae99SBarry Smith 
35447c6ae99SBarry Smith    Fortran Notes: In fortran you must pass in an array of the appropriate length.
35547c6ae99SBarry Smith 
35647c6ae99SBarry Smith    Level: intermediate
35747c6ae99SBarry Smith 
35847c6ae99SBarry Smith @*/
3597087cfbeSBarry Smith PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
36047c6ae99SBarry Smith {
36147c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
36247c6ae99SBarry Smith   PetscFunctionBegin;
36347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
36447c6ae99SBarry Smith   *ranks = dd->neighbors;
36547c6ae99SBarry Smith   PetscFunctionReturn(0);
36647c6ae99SBarry Smith }
36747c6ae99SBarry Smith 
36847c6ae99SBarry Smith #undef __FUNCT__
369aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges"
37047c6ae99SBarry Smith /*@C
371aa219208SBarry Smith       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
37247c6ae99SBarry Smith 
37347c6ae99SBarry Smith     Not Collective
37447c6ae99SBarry Smith 
37547c6ae99SBarry Smith    Input Parameter:
376aa219208SBarry Smith .     da - the DMDA object
37747c6ae99SBarry Smith 
37847c6ae99SBarry Smith    Output Parameter:
37947c6ae99SBarry Smith +     lx - ownership along x direction (optional)
38047c6ae99SBarry Smith .     ly - ownership along y direction (optional)
38147c6ae99SBarry Smith -     lz - ownership along z direction (optional)
38247c6ae99SBarry Smith 
38347c6ae99SBarry Smith    Level: intermediate
38447c6ae99SBarry Smith 
385aa219208SBarry Smith     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
38647c6ae99SBarry Smith 
38747c6ae99SBarry Smith     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
388aa219208SBarry Smith     eighth arguments from DMDAGetInfo()
38947c6ae99SBarry Smith 
39047c6ae99SBarry Smith      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
391aa219208SBarry Smith     DMDA they came from still exists (has not been destroyed).
39247c6ae99SBarry Smith 
393aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
39447c6ae99SBarry Smith @*/
3957087cfbeSBarry Smith PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
39647c6ae99SBarry Smith {
39747c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
39847c6ae99SBarry Smith 
39947c6ae99SBarry Smith   PetscFunctionBegin;
40047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
40147c6ae99SBarry Smith   if (lx) *lx = dd->lx;
40247c6ae99SBarry Smith   if (ly) *ly = dd->ly;
40347c6ae99SBarry Smith   if (lz) *lz = dd->lz;
40447c6ae99SBarry Smith   PetscFunctionReturn(0);
40547c6ae99SBarry Smith }
40647c6ae99SBarry Smith 
40747c6ae99SBarry Smith #undef __FUNCT__
408aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor"
40947c6ae99SBarry Smith /*@
410aa219208SBarry Smith      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
41147c6ae99SBarry Smith 
412aa219208SBarry Smith     Logically Collective on DMDA
41347c6ae99SBarry Smith 
41447c6ae99SBarry Smith   Input Parameters:
415aa219208SBarry Smith +    da - the DMDA object
41647c6ae99SBarry Smith .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
41747c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
41847c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
41947c6ae99SBarry Smith 
42047c6ae99SBarry Smith   Options Database:
42147c6ae99SBarry Smith +  -da_refine_x - refinement ratio in x direction
42247c6ae99SBarry Smith .  -da_refine_y - refinement ratio in y direction
42347c6ae99SBarry Smith -  -da_refine_z - refinement ratio in z direction
42447c6ae99SBarry Smith 
42547c6ae99SBarry Smith   Level: intermediate
42647c6ae99SBarry Smith 
42747c6ae99SBarry Smith     Notes: Pass PETSC_IGNORE to leave a value unchanged
42847c6ae99SBarry Smith 
429aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor()
43047c6ae99SBarry Smith @*/
4317087cfbeSBarry Smith PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
43247c6ae99SBarry Smith {
43347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
43447c6ae99SBarry Smith 
43547c6ae99SBarry Smith   PetscFunctionBegin;
43647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
43747c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_x,2);
43847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_y,3);
43947c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_z,4);
44047c6ae99SBarry Smith 
44147c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
44247c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
44347c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
44447c6ae99SBarry Smith   PetscFunctionReturn(0);
44547c6ae99SBarry Smith }
44647c6ae99SBarry Smith 
44747c6ae99SBarry Smith #undef __FUNCT__
448aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor"
44947c6ae99SBarry Smith /*@C
450aa219208SBarry Smith      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
45147c6ae99SBarry Smith 
45247c6ae99SBarry Smith     Not Collective
45347c6ae99SBarry Smith 
45447c6ae99SBarry Smith   Input Parameter:
455aa219208SBarry Smith .    da - the DMDA object
45647c6ae99SBarry Smith 
45747c6ae99SBarry Smith   Output Parameters:
45847c6ae99SBarry Smith +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
45947c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
46047c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
46147c6ae99SBarry Smith 
46247c6ae99SBarry Smith   Level: intermediate
46347c6ae99SBarry Smith 
46447c6ae99SBarry Smith     Notes: Pass PETSC_NULL for values you do not need
46547c6ae99SBarry Smith 
466aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor()
46747c6ae99SBarry Smith @*/
4687087cfbeSBarry Smith PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
46947c6ae99SBarry Smith {
47047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
47147c6ae99SBarry Smith 
47247c6ae99SBarry Smith   PetscFunctionBegin;
47347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
47447c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
47547c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
47647c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
47747c6ae99SBarry Smith   PetscFunctionReturn(0);
47847c6ae99SBarry Smith }
47947c6ae99SBarry Smith 
48047c6ae99SBarry Smith #undef __FUNCT__
481aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix"
48247c6ae99SBarry Smith /*@C
483aa219208SBarry Smith      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
48447c6ae99SBarry Smith 
485aa219208SBarry Smith     Logically Collective on DMDA
48647c6ae99SBarry Smith 
48747c6ae99SBarry Smith   Input Parameters:
488aa219208SBarry Smith +    da - the DMDA object
489aa219208SBarry Smith -    f - the function that allocates the matrix for that specific DMDA
49047c6ae99SBarry Smith 
49147c6ae99SBarry Smith   Level: developer
49247c6ae99SBarry Smith 
493aa219208SBarry Smith    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
49447c6ae99SBarry Smith        the diagonal and off-diagonal blocks of the matrix
49547c6ae99SBarry Smith 
496aa219208SBarry Smith .seealso: DMGetMatrix(), DMDASetBlockFills()
49747c6ae99SBarry Smith @*/
4987087cfbeSBarry Smith PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, const MatType,Mat*))
49947c6ae99SBarry Smith {
50047c6ae99SBarry Smith   PetscFunctionBegin;
50147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
50247c6ae99SBarry Smith   da->ops->getmatrix = f;
50347c6ae99SBarry Smith   PetscFunctionReturn(0);
50447c6ae99SBarry Smith }
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith #undef __FUNCT__
50794013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges"
50847c6ae99SBarry Smith /*
50947c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
51047c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
51147c6ae99SBarry Smith 
51247c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
51347c6ae99SBarry Smith */
51494013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
51547c6ae99SBarry Smith {
51647c6ae99SBarry Smith   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
51747c6ae99SBarry Smith   PetscErrorCode ierr;
51847c6ae99SBarry Smith 
51947c6ae99SBarry Smith   PetscFunctionBegin;
52047c6ae99SBarry Smith   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
52147c6ae99SBarry Smith   if (ratio == 1) {
52247c6ae99SBarry Smith     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
52347c6ae99SBarry Smith     PetscFunctionReturn(0);
52447c6ae99SBarry Smith   }
52547c6ae99SBarry Smith   for (i=0; i<m; i++) totalc += lc[i];
52647c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
52747c6ae99SBarry Smith   for (i=0; i<m; i++) {
52847c6ae99SBarry Smith     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
52947c6ae99SBarry Smith     if (i == m-1) lf[i] = want;
53047c6ae99SBarry Smith     else {
5317aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
5327aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
5337aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
5347aca7175SJed Brown       while ((startf+want)/ratio < nextc - stencil_width) want++;
5357aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
5367aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
5377aca7175SJed Brown       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
5387aca7175SJed Brown       /* Make sure all constraints are satisfied */
5397aca7175SJed Brown       if (want < 0 || want > remaining
5407aca7175SJed Brown           || ((startf+want)/ratio < nextc - stencil_width)
5417aca7175SJed Brown           || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width))
5427aca7175SJed Brown         SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range");
54347c6ae99SBarry Smith     }
54447c6ae99SBarry Smith     lf[i] = want;
54547c6ae99SBarry Smith     startc += lc[i];
54647c6ae99SBarry Smith     startf += lf[i];
54747c6ae99SBarry Smith     remaining -= lf[i];
54847c6ae99SBarry Smith   }
54947c6ae99SBarry Smith   PetscFunctionReturn(0);
55047c6ae99SBarry Smith }
55147c6ae99SBarry Smith 
55247c6ae99SBarry Smith #undef __FUNCT__
5532be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges"
5542be375d4SJed Brown /*
5552be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
5562be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
5572be375d4SJed Brown 
5582be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
5592be375d4SJed Brown */
5602be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
5612be375d4SJed Brown {
5622be375d4SJed Brown   PetscInt       i,totalf,remaining,startc,startf;
5632be375d4SJed Brown   PetscErrorCode ierr;
5642be375d4SJed Brown 
5652be375d4SJed Brown   PetscFunctionBegin;
5662be375d4SJed Brown   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
5672be375d4SJed Brown   if (ratio == 1) {
5682be375d4SJed Brown     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
5692be375d4SJed Brown     PetscFunctionReturn(0);
5702be375d4SJed Brown   }
5712be375d4SJed Brown   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
5722be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
5732be375d4SJed Brown   for (i=0,startc=0,startf=0; i<m; i++) {
5742be375d4SJed Brown     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
5752be375d4SJed Brown     if (i == m-1) lc[i] = want;
5762be375d4SJed Brown     else {
5772be375d4SJed Brown       const PetscInt nextf = startf+lf[i];
5782be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
5792be375d4SJed Brown        * node is within one stencil width. */
5802be375d4SJed Brown       while (nextf/ratio < startc+want-stencil_width) want--;
5812be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
5822be375d4SJed Brown        * fine node is within one stencil width. */
5832be375d4SJed Brown       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
5842be375d4SJed Brown       if (want < 0 || want > remaining
5852be375d4SJed Brown           || (nextf/ratio < startc+want-stencil_width)
5862be375d4SJed Brown           || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width))
5872be375d4SJed Brown         SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range");
5882be375d4SJed Brown     }
5892be375d4SJed Brown     lc[i] = want;
5902be375d4SJed Brown     startc += lc[i];
5912be375d4SJed Brown     startf += lf[i];
5922be375d4SJed Brown     remaining -= lc[i];
5932be375d4SJed Brown   }
5942be375d4SJed Brown   PetscFunctionReturn(0);
5952be375d4SJed Brown }
5962be375d4SJed Brown 
5972be375d4SJed Brown #undef __FUNCT__
59894013140SBarry Smith #define __FUNCT__ "DMRefine_DA"
5997087cfbeSBarry Smith PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
60047c6ae99SBarry Smith {
60147c6ae99SBarry Smith   PetscErrorCode ierr;
602*f3141302SJed Brown   PetscInt       M,N,P,i;
6039a42bb27SBarry Smith   DM             da2;
60447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
60547c6ae99SBarry Smith 
60647c6ae99SBarry Smith   PetscFunctionBegin;
60747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
60847c6ae99SBarry Smith   PetscValidPointer(daref,3);
60947c6ae99SBarry Smith 
6101321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
61147c6ae99SBarry Smith     M = dd->refine_x*dd->M;
61247c6ae99SBarry Smith   } else {
61347c6ae99SBarry Smith     M = 1 + dd->refine_x*(dd->M - 1);
61447c6ae99SBarry Smith   }
6151321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
61647c6ae99SBarry Smith     N = dd->refine_y*dd->N;
61747c6ae99SBarry Smith   } else {
61847c6ae99SBarry Smith     N = 1 + dd->refine_y*(dd->N - 1);
61947c6ae99SBarry Smith   }
6201321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
62147c6ae99SBarry Smith     P = dd->refine_z*dd->P;
62247c6ae99SBarry Smith   } else {
62347c6ae99SBarry Smith     P = 1 + dd->refine_z*(dd->P - 1);
62447c6ae99SBarry Smith   }
62547c6ae99SBarry Smith   if (dd->dim == 3) {
62647c6ae99SBarry Smith     PetscInt *lx,*ly,*lz;
62747c6ae99SBarry Smith     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
6281321219cSEthan 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);
6291321219cSEthan 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);
6301321219cSEthan 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);
6311321219cSEthan Coon     ierr = DMDACreate3d(((PetscObject)da)->comm,dd->bx,dd->by,dd->bz,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr);
63247c6ae99SBarry Smith     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
63347c6ae99SBarry Smith   } else if (dd->dim == 2) {
63447c6ae99SBarry Smith     PetscInt *lx,*ly;
63547c6ae99SBarry Smith     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
6361321219cSEthan 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);
6371321219cSEthan 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);
6381321219cSEthan Coon     ierr = DMDACreate2d(((PetscObject)da)->comm,dd->bx,dd->by,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr);
63947c6ae99SBarry Smith     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
64047c6ae99SBarry Smith   } else if (dd->dim == 1) {
64147c6ae99SBarry Smith     PetscInt *lx;
64247c6ae99SBarry Smith     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
6431321219cSEthan 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);
6441321219cSEthan Coon     ierr = DMDACreate1d(((PetscObject)da)->comm,dd->bx,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr);
64547c6ae99SBarry Smith     ierr = PetscFree(lx);CHKERRQ(ierr);
64647c6ae99SBarry Smith   }
64747c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
64847c6ae99SBarry Smith 
64947c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
65047c6ae99SBarry Smith   da2->ops->getmatrix        = da->ops->getmatrix;
6514c661befSBarry Smith   /* da2->ops->getinterpolation = da->ops->getinterpolation; this causes problem with SNESVI */
65247c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
65347c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
65447c6ae99SBarry Smith 
65547c6ae99SBarry Smith   /* copy fill information if given */
65647c6ae99SBarry Smith   if (dd->dfill) {
65747c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
65847c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
65947c6ae99SBarry Smith   }
66047c6ae99SBarry Smith   if (dd->ofill) {
66147c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
66247c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
66347c6ae99SBarry Smith   }
66447c6ae99SBarry Smith   /* copy the refine information */
66547c6ae99SBarry Smith   dd2->refine_x = dd->refine_x;
66647c6ae99SBarry Smith   dd2->refine_y = dd->refine_y;
66747c6ae99SBarry Smith   dd2->refine_z = dd->refine_z;
66847c6ae99SBarry Smith 
66947c6ae99SBarry Smith   /* copy vector type information */
67047c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
67147c6ae99SBarry Smith   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
672ddcf8b74SDave May 
673efd51863SBarry Smith   dd2->lf = dd->lf;
674efd51863SBarry Smith   dd2->lj = dd->lj;
675efd51863SBarry Smith 
676ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
677ddcf8b74SDave May   if (dd->coordinates) {
678ddcf8b74SDave May     DM  cdaf,cdac;
679ddcf8b74SDave May     Vec coordsc,coordsf;
680ddcf8b74SDave May     Mat II;
681ddcf8b74SDave May 
682ddcf8b74SDave May     ierr = DMDAGetCoordinateDA(da,&cdac);CHKERRQ(ierr);
683ddcf8b74SDave May     ierr = DMDAGetCoordinates(da,&coordsc);CHKERRQ(ierr);
684b61d3410SDave May     ierr = DMDAGetCoordinateDA(da2,&cdaf);CHKERRQ(ierr);
685b61d3410SDave May     /* force creation of the coordinate vector */
686b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
687ddcf8b74SDave May     ierr = DMDAGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
688ddcf8b74SDave May     ierr = DMGetInterpolation(cdac,cdaf,&II,PETSC_NULL);CHKERRQ(ierr);
689ddcf8b74SDave May     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
690fcfd50ebSBarry Smith     ierr = MatDestroy(&II);CHKERRQ(ierr);
691ddcf8b74SDave May   }
692*f3141302SJed Brown   for (i=0; i<da->bs; i++) {
693*f3141302SJed Brown     const char *fieldname;
694*f3141302SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
695*f3141302SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
696*f3141302SJed Brown   }
69747c6ae99SBarry Smith   *daref = da2;
69847c6ae99SBarry Smith   PetscFunctionReturn(0);
69947c6ae99SBarry Smith }
70047c6ae99SBarry Smith 
70147c6ae99SBarry Smith #undef __FUNCT__
70294013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA"
7037087cfbeSBarry Smith PetscErrorCode  DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
70447c6ae99SBarry Smith {
70547c6ae99SBarry Smith   PetscErrorCode ierr;
70647c6ae99SBarry Smith   PetscInt       M,N,P;
7079a42bb27SBarry Smith   DM             da2;
70847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
70947c6ae99SBarry Smith 
71047c6ae99SBarry Smith   PetscFunctionBegin;
71147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
71247c6ae99SBarry Smith   PetscValidPointer(daref,3);
71347c6ae99SBarry Smith 
7141321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
71547c6ae99SBarry Smith     M = dd->M / dd->refine_x;
71647c6ae99SBarry Smith   } else {
71747c6ae99SBarry Smith     M = 1 + (dd->M - 1) / dd->refine_x;
71847c6ae99SBarry Smith   }
7191321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
72047c6ae99SBarry Smith     N = dd->N / dd->refine_y;
72147c6ae99SBarry Smith   } else {
72247c6ae99SBarry Smith     N = 1 + (dd->N - 1) / dd->refine_y;
72347c6ae99SBarry Smith   }
7241321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
72547c6ae99SBarry Smith     P = dd->P / dd->refine_z;
72647c6ae99SBarry Smith   } else {
72747c6ae99SBarry Smith     P = 1 + (dd->P - 1) / dd->refine_z;
72847c6ae99SBarry Smith   }
72947c6ae99SBarry Smith   if (dd->dim == 3) {
7302be375d4SJed Brown     PetscInt *lx,*ly,*lz;
7312be375d4SJed Brown     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
7321321219cSEthan Coon     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
7331321219cSEthan Coon     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
7341321219cSEthan Coon     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
7351321219cSEthan Coon     ierr = DMDACreate3d(((PetscObject)da)->comm,dd->bx,dd->by,dd->bz,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr);
7362be375d4SJed Brown     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
73747c6ae99SBarry Smith   } else if (dd->dim == 2) {
7382be375d4SJed Brown     PetscInt *lx,*ly;
7392be375d4SJed Brown     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
7401321219cSEthan Coon     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
7411321219cSEthan Coon     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
7421321219cSEthan Coon     ierr = DMDACreate2d(((PetscObject)da)->comm,dd->bx,dd->by,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr);
7432be375d4SJed Brown     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
74447c6ae99SBarry Smith   } else if (dd->dim == 1) {
7452be375d4SJed Brown     PetscInt *lx;
7462be375d4SJed Brown     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
7471321219cSEthan Coon     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
7481321219cSEthan Coon     ierr = DMDACreate1d(((PetscObject)da)->comm,dd->bx,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr);
7492be375d4SJed Brown     ierr = PetscFree(lx);CHKERRQ(ierr);
75047c6ae99SBarry Smith   }
75147c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
75247c6ae99SBarry Smith 
7534dcab191SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
7544dcab191SBarry Smith   /* da2->ops->getinterpolation = da->ops->getinterpolation; copying this one causes trouble for DMSetVI */
75547c6ae99SBarry Smith   da2->ops->getmatrix        = da->ops->getmatrix;
75647c6ae99SBarry Smith   da2->ops->getcoloring      = da->ops->getcoloring;
75747c6ae99SBarry Smith   dd2->interptype            = dd->interptype;
75847c6ae99SBarry Smith 
75947c6ae99SBarry Smith   /* copy fill information if given */
76047c6ae99SBarry Smith   if (dd->dfill) {
76147c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
76247c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
76347c6ae99SBarry Smith   }
76447c6ae99SBarry Smith   if (dd->ofill) {
76547c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
76647c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
76747c6ae99SBarry Smith   }
76847c6ae99SBarry Smith   /* copy the refine information */
76947c6ae99SBarry Smith   dd2->refine_x = dd->refine_x;
77047c6ae99SBarry Smith   dd2->refine_y = dd->refine_y;
77147c6ae99SBarry Smith   dd2->refine_z = dd->refine_z;
77247c6ae99SBarry Smith 
77347c6ae99SBarry Smith   /* copy vector type information */
77447c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
77547c6ae99SBarry Smith   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
77647c6ae99SBarry Smith 
777644e2e5bSBarry Smith   dd2->lf = dd->lf;
778644e2e5bSBarry Smith   dd2->lj = dd->lj;
779644e2e5bSBarry Smith 
780ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
781ddcf8b74SDave May   if (dd->coordinates) {
782ddcf8b74SDave May     DM         cdaf,cdac;
783ddcf8b74SDave May     Vec        coordsc,coordsf;
784ddcf8b74SDave May     VecScatter inject;
785ddcf8b74SDave May 
786ddcf8b74SDave May     ierr = DMDAGetCoordinateDA(da,&cdaf);CHKERRQ(ierr);
787ddcf8b74SDave May     ierr = DMDAGetCoordinates(da,&coordsf);CHKERRQ(ierr);
788b61d3410SDave May     ierr = DMDAGetCoordinateDA(da2,&cdac);CHKERRQ(ierr);
789b61d3410SDave May     /* force creation of the coordinate vector */
790b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
791ddcf8b74SDave May     ierr = DMDAGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
792ddcf8b74SDave May 
793ddcf8b74SDave May     ierr = DMGetInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
794ddcf8b74SDave May     ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
795ddcf8b74SDave May     ierr = VecScatterEnd(inject  ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
796fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&inject);CHKERRQ(ierr);
797ddcf8b74SDave May   }
79847c6ae99SBarry Smith   *daref = da2;
79947c6ae99SBarry Smith   PetscFunctionReturn(0);
80047c6ae99SBarry Smith }
80147c6ae99SBarry Smith 
80247c6ae99SBarry Smith #undef __FUNCT__
80394013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA"
8047087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
80547c6ae99SBarry Smith {
80647c6ae99SBarry Smith   PetscErrorCode ierr;
80747c6ae99SBarry Smith   PetscInt       i,n,*refx,*refy,*refz;
80847c6ae99SBarry Smith 
80947c6ae99SBarry Smith   PetscFunctionBegin;
81047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
81147c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
81247c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
81347c6ae99SBarry Smith   PetscValidPointer(daf,3);
81447c6ae99SBarry Smith 
815aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
81647c6ae99SBarry Smith   ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr);
81747c6ae99SBarry Smith   for (i=0; i<nlevels; i++) {
818aa219208SBarry Smith     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
81947c6ae99SBarry Smith   }
82047c6ae99SBarry Smith   n = nlevels;
82147c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
82247c6ae99SBarry Smith   n = nlevels;
82347c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
82447c6ae99SBarry Smith   n = nlevels;
82547c6ae99SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
82647c6ae99SBarry Smith 
827aa219208SBarry Smith   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
82894013140SBarry Smith   ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr);
82947c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
830aa219208SBarry Smith     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
83194013140SBarry Smith     ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr);
83247c6ae99SBarry Smith   }
83347c6ae99SBarry Smith   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
83447c6ae99SBarry Smith   PetscFunctionReturn(0);
83547c6ae99SBarry Smith }
83647c6ae99SBarry Smith 
83747c6ae99SBarry Smith #undef __FUNCT__
83894013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA"
8397087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
84047c6ae99SBarry Smith {
84147c6ae99SBarry Smith   PetscErrorCode ierr;
84247c6ae99SBarry Smith   PetscInt       i;
84347c6ae99SBarry Smith 
84447c6ae99SBarry Smith   PetscFunctionBegin;
84547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
84647c6ae99SBarry Smith   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
84747c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
84847c6ae99SBarry Smith   PetscValidPointer(dac,3);
84994013140SBarry Smith   ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr);
85047c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
85194013140SBarry Smith     ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr);
85247c6ae99SBarry Smith   }
85347c6ae99SBarry Smith   PetscFunctionReturn(0);
85447c6ae99SBarry Smith }
855