xref: /petsc/src/dm/impls/da/da.c (revision 3e7870d263ac1815e57c04a2ae8ebb11bfa82906)
14035e84dSBarry Smith #include <petsc-private/dmdaimpl.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);
25ce94432eSBarry Smith   if (dd->dim > 0 && dim != dd->dim) SETERRQ2(PetscObjectComm((PetscObject)da),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);
56ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),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;
84e3f3e4b6SBarry 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);
91ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),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;
95e3f3e4b6SBarry Smith   if (dd->dim == 2) {
96e3f3e4b6SBarry Smith     PetscMPIInt size;
97ce94432eSBarry Smith     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr);
98e3f3e4b6SBarry Smith     if ((dd->m > 0) && (dd->n < 0)) {
99e3f3e4b6SBarry Smith       dd->n = size/dd->m;
100ce94432eSBarry Smith       if (dd->n*dd->m != size) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"%D processes in X direction not divisible into comm size %d",m,size);
101e3f3e4b6SBarry Smith     }
102e3f3e4b6SBarry Smith     if ((dd->n > 0) && (dd->m < 0)) {
103e3f3e4b6SBarry Smith       dd->m = size/dd->n;
104ce94432eSBarry Smith       if (dd->n*dd->m != size) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"%D processes in Y direction not divisible into comm size %d",n,size);
105e3f3e4b6SBarry Smith     }
106e3f3e4b6SBarry 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);
135ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),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);
165ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),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__
1727ddda789SPeter Brune #define __FUNCT__ "DMDAGetOverlap"
1737ddda789SPeter Brune /*@
1747ddda789SPeter Brune   DMDAGetOverlap - Gets the size of the per-processor overlap.
1757ddda789SPeter Brune 
1767ddda789SPeter Brune   Not collective
1777ddda789SPeter Brune 
1787ddda789SPeter Brune   Input Parameters:
1797ddda789SPeter Brune . da  - The DMDA
1807ddda789SPeter Brune 
1817ddda789SPeter Brune   Output Parameters:
1827ddda789SPeter Brune + x   - Overlap in the x direction
1837ddda789SPeter Brune . y   - Overlap in the y direction
1847ddda789SPeter Brune - z   - Overlap in the z direction
1857ddda789SPeter Brune 
1867ddda789SPeter Brune   Level: intermediate
1877ddda789SPeter Brune 
1887ddda789SPeter Brune .keywords:  distributed array, overlap, domain decomposition
1897ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetOverlap(), DMDA
1907ddda789SPeter Brune @*/
1917ddda789SPeter Brune PetscErrorCode  DMDAGetOverlap(DM da,PetscInt *x,PetscInt *y,PetscInt *z)
1927ddda789SPeter Brune {
1937ddda789SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
1947ddda789SPeter Brune 
1957ddda789SPeter Brune   PetscFunctionBegin;
1967ddda789SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1977ddda789SPeter Brune   if (x) *x = dd->xol;
1987ddda789SPeter Brune   if (y) *y = dd->yol;
1997ddda789SPeter Brune   if (z) *z = dd->zol;
2007ddda789SPeter Brune   PetscFunctionReturn(0);
2017ddda789SPeter Brune }
2027ddda789SPeter Brune 
2037ddda789SPeter Brune #undef __FUNCT__
20488661749SPeter Brune #define __FUNCT__ "DMDASetOverlap"
20588661749SPeter Brune /*@
20688661749SPeter Brune   DMDASetOverlap - Sets the size of the per-processor overlap.
20788661749SPeter Brune 
20888661749SPeter Brune   Not collective
20988661749SPeter Brune 
2107ddda789SPeter Brune   Input Parameters:
21188661749SPeter Brune + da  - The DMDA
2127ddda789SPeter Brune . x   - Overlap in the x direction
2137ddda789SPeter Brune . y   - Overlap in the y direction
2147ddda789SPeter Brune - z   - Overlap in the z direction
21588661749SPeter Brune 
21688661749SPeter Brune   Level: intermediate
21788661749SPeter Brune 
2187ddda789SPeter Brune .keywords:  distributed array, overlap, domain decomposition
2197ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetOverlap(), DMDA
22088661749SPeter Brune @*/
2217ddda789SPeter Brune PetscErrorCode  DMDASetOverlap(DM da,PetscInt x,PetscInt y,PetscInt z)
22288661749SPeter Brune {
22388661749SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
22488661749SPeter Brune 
22588661749SPeter Brune   PetscFunctionBegin;
22688661749SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
2277ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da,x,2);
2287ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da,y,3);
2297ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da,z,4);
2307ddda789SPeter Brune   dd->xol = x;
2317ddda789SPeter Brune   dd->yol = y;
2327ddda789SPeter Brune   dd->zol = z;
23388661749SPeter Brune   PetscFunctionReturn(0);
23488661749SPeter Brune }
23588661749SPeter Brune 
236*3e7870d2SPeter Brune 
237*3e7870d2SPeter Brune #undef __FUNCT__
238*3e7870d2SPeter Brune #define __FUNCT__ "DMDAGetNumLocalSubDomains"
239*3e7870d2SPeter Brune /*@
240*3e7870d2SPeter Brune   DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition.
241*3e7870d2SPeter Brune 
242*3e7870d2SPeter Brune   Not collective
243*3e7870d2SPeter Brune 
244*3e7870d2SPeter Brune   Input Parameters:
245*3e7870d2SPeter Brune . da  - The DMDA
246*3e7870d2SPeter Brune 
247*3e7870d2SPeter Brune   Output Parameters:
248*3e7870d2SPeter Brune + Nsub   - Number of local subdomains created upon decomposition
249*3e7870d2SPeter Brune 
250*3e7870d2SPeter Brune   Level: intermediate
251*3e7870d2SPeter Brune 
252*3e7870d2SPeter Brune .keywords:  distributed array, domain decomposition
253*3e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetNumLocalSubDomains(), DMDA
254*3e7870d2SPeter Brune @*/
255*3e7870d2SPeter Brune PetscErrorCode  DMDAGetNumLocalSubDomains(DM da,PetscInt *Nsub)
256*3e7870d2SPeter Brune {
257*3e7870d2SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
258*3e7870d2SPeter Brune 
259*3e7870d2SPeter Brune   PetscFunctionBegin;
260*3e7870d2SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
261*3e7870d2SPeter Brune   if (Nsub) *Nsub = dd->Nsub;
262*3e7870d2SPeter Brune   PetscFunctionReturn(0);
263*3e7870d2SPeter Brune }
264*3e7870d2SPeter Brune 
265*3e7870d2SPeter Brune #undef __FUNCT__
266*3e7870d2SPeter Brune #define __FUNCT__ "DMDASetNumLocalSubDomains"
267*3e7870d2SPeter Brune /*@
268*3e7870d2SPeter Brune   DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition.
269*3e7870d2SPeter Brune 
270*3e7870d2SPeter Brune   Not collective
271*3e7870d2SPeter Brune 
272*3e7870d2SPeter Brune   Input Parameters:
273*3e7870d2SPeter Brune + da  - The DMDA
274*3e7870d2SPeter Brune - Nsub - The number of local subdomains requested
275*3e7870d2SPeter Brune 
276*3e7870d2SPeter Brune   Level: intermediate
277*3e7870d2SPeter Brune 
278*3e7870d2SPeter Brune .keywords:  distributed array, domain decomposition
279*3e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetNumLocalSubDomains(), DMDA
280*3e7870d2SPeter Brune @*/
281*3e7870d2SPeter Brune PetscErrorCode  DMDASetNumLocalSubDomains(DM da,PetscInt Nsub)
282*3e7870d2SPeter Brune {
283*3e7870d2SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
284*3e7870d2SPeter Brune 
285*3e7870d2SPeter Brune   PetscFunctionBegin;
286*3e7870d2SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
287*3e7870d2SPeter Brune   PetscValidLogicalCollectiveInt(da,Nsub,2);
288*3e7870d2SPeter Brune   dd->Nsub = Nsub;
289*3e7870d2SPeter Brune   PetscFunctionReturn(0);
290*3e7870d2SPeter Brune }
291*3e7870d2SPeter Brune 
292d886c4f4SPeter Brune #undef __FUNCT__
293d886c4f4SPeter Brune #define __FUNCT__ "DMDASetOffset"
294d886c4f4SPeter Brune /*@
295d886c4f4SPeter Brune   DMDASetOffset - Sets the index offset of the DA.
296d886c4f4SPeter Brune 
297d886c4f4SPeter Brune   Collective on DA
298d886c4f4SPeter Brune 
299d886c4f4SPeter Brune   Input Parameter:
300d886c4f4SPeter Brune + da  - The DMDA
301d886c4f4SPeter Brune . xo  - The offset in the x direction
302d886c4f4SPeter Brune . yo  - The offset in the y direction
303d886c4f4SPeter Brune - zo  - The offset in the z direction
304d886c4f4SPeter Brune 
305d886c4f4SPeter Brune   Level: intermediate
306d886c4f4SPeter Brune 
307d886c4f4SPeter Brune   Notes: This is used primarily to overlap a computation on a local DA with that on a global DA without
308d886c4f4SPeter Brune   changing boundary conditions or subdomain features that depend upon the global offsets.
309d886c4f4SPeter Brune 
310d886c4f4SPeter Brune .keywords:  distributed array, degrees of freedom
311d886c4f4SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray()
312d886c4f4SPeter Brune @*/
31395c13181SPeter Brune PetscErrorCode  DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po)
314d886c4f4SPeter Brune {
31595c13181SPeter Brune   PetscErrorCode ierr;
316d886c4f4SPeter Brune   DM_DA          *dd = (DM_DA*)da->data;
317d886c4f4SPeter Brune 
318d886c4f4SPeter Brune   PetscFunctionBegin;
319d886c4f4SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
320d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da,xo,2);
32195c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,yo,3);
32295c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,zo,4);
32395c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,Mo,5);
32495c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,No,6);
32595c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,Po,7);
326d886c4f4SPeter Brune   dd->xo = xo;
327d886c4f4SPeter Brune   dd->yo = yo;
328d886c4f4SPeter Brune   dd->zo = zo;
32995c13181SPeter Brune   dd->Mo = Mo;
33095c13181SPeter Brune   dd->No = No;
33195c13181SPeter Brune   dd->Po = Po;
33295c13181SPeter Brune 
33395c13181SPeter Brune   if (da->coordinateDM) {
33495c13181SPeter Brune     ierr = DMDASetOffset(da->coordinateDM,xo,yo,zo,Mo,No,Po);CHKERRQ(ierr);
33595c13181SPeter Brune   }
336d886c4f4SPeter Brune   PetscFunctionReturn(0);
337d886c4f4SPeter Brune }
338d886c4f4SPeter Brune 
339d886c4f4SPeter Brune #undef __FUNCT__
340d886c4f4SPeter Brune #define __FUNCT__ "DMDAGetOffset"
341d886c4f4SPeter Brune /*@
342d886c4f4SPeter Brune   DMDAGetOffset - Gets the index offset of the DA.
343d886c4f4SPeter Brune 
344d886c4f4SPeter Brune   Not collective
345d886c4f4SPeter Brune 
346d886c4f4SPeter Brune   Input Parameter:
347d886c4f4SPeter Brune . da  - The DMDA
348d886c4f4SPeter Brune 
349d886c4f4SPeter Brune   Output Parameters:
350d886c4f4SPeter Brune + xo  - The offset in the x direction
351d886c4f4SPeter Brune . yo  - The offset in the y direction
35295c13181SPeter Brune . zo  - The offset in the z direction
35395c13181SPeter Brune . Mo  - The global size in the x direction
35495c13181SPeter Brune . No  - The global size in the y direction
35595c13181SPeter Brune - Po  - The global size in the z direction
356d886c4f4SPeter Brune 
357d886c4f4SPeter Brune   Level: intermediate
358d886c4f4SPeter Brune 
359d886c4f4SPeter Brune .keywords:  distributed array, degrees of freedom
360d886c4f4SPeter Brune .seealso: DMDASetOffset(), DMDAVecGetArray()
361d886c4f4SPeter Brune @*/
36295c13181SPeter Brune PetscErrorCode  DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po)
363d886c4f4SPeter Brune {
364d886c4f4SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
365d886c4f4SPeter Brune 
366d886c4f4SPeter Brune   PetscFunctionBegin;
367d886c4f4SPeter Brune   PetscValidHeaderSpecific(da,DM_CLASSID,1);
368d886c4f4SPeter Brune   if (xo) *xo = dd->xo;
369d886c4f4SPeter Brune   if (yo) *yo = dd->yo;
370d886c4f4SPeter Brune   if (zo) *zo = dd->zo;
37195c13181SPeter Brune   if (Mo) *Mo = dd->Mo;
37295c13181SPeter Brune   if (No) *No = dd->No;
37395c13181SPeter Brune   if (Po) *Po = dd->Po;
374d886c4f4SPeter Brune   PetscFunctionReturn(0);
375d886c4f4SPeter Brune }
376d886c4f4SPeter Brune 
37788661749SPeter Brune 
37888661749SPeter Brune #undef __FUNCT__
379aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilType"
38047c6ae99SBarry Smith /*@
381aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
38247c6ae99SBarry Smith 
383aa219208SBarry Smith   Logically Collective on DMDA
38447c6ae99SBarry Smith 
38547c6ae99SBarry Smith   Input Parameter:
386aa219208SBarry Smith + da    - The DMDA
387aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
38847c6ae99SBarry Smith 
38947c6ae99SBarry Smith   Level: intermediate
39047c6ae99SBarry Smith 
39147c6ae99SBarry Smith .keywords:  distributed array, stencil
39296e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
39347c6ae99SBarry Smith @*/
3947087cfbeSBarry Smith PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
39547c6ae99SBarry Smith {
39647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
39747c6ae99SBarry Smith 
39847c6ae99SBarry Smith   PetscFunctionBegin;
39947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
40047c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,stype,2);
401ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
40247c6ae99SBarry Smith   dd->stencil_type = stype;
40347c6ae99SBarry Smith   PetscFunctionReturn(0);
40447c6ae99SBarry Smith }
40547c6ae99SBarry Smith 
40647c6ae99SBarry Smith #undef __FUNCT__
407aa219208SBarry Smith #define __FUNCT__ "DMDASetStencilWidth"
40847c6ae99SBarry Smith /*@
409aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
41047c6ae99SBarry Smith 
411aa219208SBarry Smith   Logically Collective on DMDA
41247c6ae99SBarry Smith 
41347c6ae99SBarry Smith   Input Parameter:
414aa219208SBarry Smith + da    - The DMDA
41547c6ae99SBarry Smith - width - The stencil width
41647c6ae99SBarry Smith 
41747c6ae99SBarry Smith   Level: intermediate
41847c6ae99SBarry Smith 
41947c6ae99SBarry Smith .keywords:  distributed array, stencil
42096e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
42147c6ae99SBarry Smith @*/
4227087cfbeSBarry Smith PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
42347c6ae99SBarry Smith {
42447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
42547c6ae99SBarry Smith 
42647c6ae99SBarry Smith   PetscFunctionBegin;
42747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
42847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,width,2);
429ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
43047c6ae99SBarry Smith   dd->s = width;
43147c6ae99SBarry Smith   PetscFunctionReturn(0);
43247c6ae99SBarry Smith }
43347c6ae99SBarry Smith 
43447c6ae99SBarry Smith #undef __FUNCT__
435aa219208SBarry Smith #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
436aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
43747c6ae99SBarry Smith {
43847c6ae99SBarry Smith   PetscInt i,sum;
43947c6ae99SBarry Smith 
44047c6ae99SBarry Smith   PetscFunctionBegin;
441ce94432eSBarry Smith   if (M < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
44247c6ae99SBarry Smith   for (i=sum=0; i<m; i++) sum += lx[i];
443ce94432eSBarry Smith   if (sum != M) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
44447c6ae99SBarry Smith   PetscFunctionReturn(0);
44547c6ae99SBarry Smith }
44647c6ae99SBarry Smith 
44747c6ae99SBarry Smith #undef __FUNCT__
448aa219208SBarry Smith #define __FUNCT__ "DMDASetOwnershipRanges"
44947c6ae99SBarry Smith /*@
450aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
45147c6ae99SBarry Smith 
452aa219208SBarry Smith   Logically Collective on DMDA
45347c6ae99SBarry Smith 
45447c6ae99SBarry Smith   Input Parameter:
455aa219208SBarry Smith + da - The DMDA
4560298fd71SBarry Smith . lx - array containing number of nodes in the X direction on each process, or NULL. If non-null, must be of length da->m
4570298fd71SBarry Smith . ly - array containing number of nodes in the Y direction on each process, or NULL. If non-null, must be of length da->n
4580298fd71SBarry Smith - lz - array containing number of nodes in the Z direction on each process, or NULL. If non-null, must be of length da->p.
45947c6ae99SBarry Smith 
46047c6ae99SBarry Smith   Level: intermediate
46147c6ae99SBarry Smith 
462e3f3e4b6SBarry Smith   Note: these numbers are NOT multiplied by the number of dof per node.
463e3f3e4b6SBarry Smith 
46447c6ae99SBarry Smith .keywords:  distributed array
46596e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
46647c6ae99SBarry Smith @*/
4677087cfbeSBarry Smith PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
46847c6ae99SBarry Smith {
46947c6ae99SBarry Smith   PetscErrorCode ierr;
47047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
47147c6ae99SBarry Smith 
47247c6ae99SBarry Smith   PetscFunctionBegin;
47347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
474ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
47547c6ae99SBarry Smith   if (lx) {
476ce94432eSBarry Smith     if (dd->m < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
477aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
47847c6ae99SBarry Smith     if (!dd->lx) {
47947c6ae99SBarry Smith       ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
48047c6ae99SBarry Smith     }
48147c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
48247c6ae99SBarry Smith   }
48347c6ae99SBarry Smith   if (ly) {
484ce94432eSBarry Smith     if (dd->n < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
485aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
48647c6ae99SBarry Smith     if (!dd->ly) {
48747c6ae99SBarry Smith       ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
48847c6ae99SBarry Smith     }
48947c6ae99SBarry Smith     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
49047c6ae99SBarry Smith   }
49147c6ae99SBarry Smith   if (lz) {
492ce94432eSBarry Smith     if (dd->p < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
493aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
49447c6ae99SBarry Smith     if (!dd->lz) {
49547c6ae99SBarry Smith       ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
49647c6ae99SBarry Smith     }
49747c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
49847c6ae99SBarry Smith   }
49947c6ae99SBarry Smith   PetscFunctionReturn(0);
50047c6ae99SBarry Smith }
50147c6ae99SBarry Smith 
50247c6ae99SBarry Smith #undef __FUNCT__
503aa219208SBarry Smith #define __FUNCT__ "DMDASetInterpolationType"
50447c6ae99SBarry Smith /*@
505aa219208SBarry Smith        DMDASetInterpolationType - Sets the type of interpolation that will be
506e727c939SJed Brown           returned by DMCreateInterpolation()
50747c6ae99SBarry Smith 
508aa219208SBarry Smith    Logically Collective on DMDA
50947c6ae99SBarry Smith 
51047c6ae99SBarry Smith    Input Parameter:
51147c6ae99SBarry Smith +  da - initial distributed array
512aa219208SBarry Smith .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
51347c6ae99SBarry Smith 
51447c6ae99SBarry Smith    Level: intermediate
51547c6ae99SBarry Smith 
516e727c939SJed Brown    Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation()
51747c6ae99SBarry Smith 
51847c6ae99SBarry Smith .keywords:  distributed array, interpolation
51947c6ae99SBarry Smith 
52096e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
52147c6ae99SBarry Smith @*/
5227087cfbeSBarry Smith PetscErrorCode  DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
52347c6ae99SBarry Smith {
52447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
52547c6ae99SBarry Smith 
52647c6ae99SBarry Smith   PetscFunctionBegin;
52747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
52847c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,ctype,2);
52947c6ae99SBarry Smith   dd->interptype = ctype;
53047c6ae99SBarry Smith   PetscFunctionReturn(0);
53147c6ae99SBarry Smith }
53247c6ae99SBarry Smith 
5332dde6fd4SLisandro Dalcin #undef __FUNCT__
5342dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetInterpolationType"
5352dde6fd4SLisandro Dalcin /*@
5362dde6fd4SLisandro Dalcin        DMDAGetInterpolationType - Gets the type of interpolation that will be
537e727c939SJed Brown           used by DMCreateInterpolation()
5382dde6fd4SLisandro Dalcin 
5392dde6fd4SLisandro Dalcin    Not Collective
5402dde6fd4SLisandro Dalcin 
5412dde6fd4SLisandro Dalcin    Input Parameter:
5422dde6fd4SLisandro Dalcin .  da - distributed array
5432dde6fd4SLisandro Dalcin 
5442dde6fd4SLisandro Dalcin    Output Parameter:
5452dde6fd4SLisandro Dalcin .  ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms)
5462dde6fd4SLisandro Dalcin 
5472dde6fd4SLisandro Dalcin    Level: intermediate
5482dde6fd4SLisandro Dalcin 
5492dde6fd4SLisandro Dalcin .keywords:  distributed array, interpolation
5502dde6fd4SLisandro Dalcin 
551e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation()
5522dde6fd4SLisandro Dalcin @*/
5532dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype)
5542dde6fd4SLisandro Dalcin {
5552dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA*)da->data;
5562dde6fd4SLisandro Dalcin 
5572dde6fd4SLisandro Dalcin   PetscFunctionBegin;
5582dde6fd4SLisandro Dalcin   PetscValidHeaderSpecific(da,DM_CLASSID,1);
5592dde6fd4SLisandro Dalcin   PetscValidPointer(ctype,2);
5602dde6fd4SLisandro Dalcin   *ctype = dd->interptype;
5612dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
5622dde6fd4SLisandro Dalcin }
56347c6ae99SBarry Smith 
56447c6ae99SBarry Smith #undef __FUNCT__
565aa219208SBarry Smith #define __FUNCT__ "DMDAGetNeighbors"
56647c6ae99SBarry Smith /*@C
567aa219208SBarry Smith       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
56847c6ae99SBarry Smith         processes neighbors.
56947c6ae99SBarry Smith 
57047c6ae99SBarry Smith     Not Collective
57147c6ae99SBarry Smith 
57247c6ae99SBarry Smith    Input Parameter:
573aa219208SBarry Smith .     da - the DMDA object
57447c6ae99SBarry Smith 
57547c6ae99SBarry Smith    Output Parameters:
57647c6ae99SBarry Smith .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
57747c6ae99SBarry Smith               this process itself is in the list
57847c6ae99SBarry Smith 
57947c6ae99SBarry Smith    Notes: In 2d the array is of length 9, in 3d of length 27
58047c6ae99SBarry Smith           Not supported in 1d
581aa219208SBarry Smith           Do not free the array, it is freed when the DMDA is destroyed.
58247c6ae99SBarry Smith 
58347c6ae99SBarry Smith    Fortran Notes: In fortran you must pass in an array of the appropriate length.
58447c6ae99SBarry Smith 
58547c6ae99SBarry Smith    Level: intermediate
58647c6ae99SBarry Smith 
58747c6ae99SBarry Smith @*/
5887087cfbeSBarry Smith PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
58947c6ae99SBarry Smith {
59047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
5915fd66863SKarl Rupp 
59247c6ae99SBarry Smith   PetscFunctionBegin;
59347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
59447c6ae99SBarry Smith   *ranks = dd->neighbors;
59547c6ae99SBarry Smith   PetscFunctionReturn(0);
59647c6ae99SBarry Smith }
59747c6ae99SBarry Smith 
59847c6ae99SBarry Smith #undef __FUNCT__
599aa219208SBarry Smith #define __FUNCT__ "DMDAGetOwnershipRanges"
60047c6ae99SBarry Smith /*@C
601aa219208SBarry Smith       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
60247c6ae99SBarry Smith 
60347c6ae99SBarry Smith     Not Collective
60447c6ae99SBarry Smith 
60547c6ae99SBarry Smith    Input Parameter:
606aa219208SBarry Smith .     da - the DMDA object
60747c6ae99SBarry Smith 
60847c6ae99SBarry Smith    Output Parameter:
60947c6ae99SBarry Smith +     lx - ownership along x direction (optional)
61047c6ae99SBarry Smith .     ly - ownership along y direction (optional)
61147c6ae99SBarry Smith -     lz - ownership along z direction (optional)
61247c6ae99SBarry Smith 
61347c6ae99SBarry Smith    Level: intermediate
61447c6ae99SBarry Smith 
615aa219208SBarry Smith     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
61647c6ae99SBarry Smith 
61747c6ae99SBarry Smith     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
618aa219208SBarry Smith     eighth arguments from DMDAGetInfo()
61947c6ae99SBarry Smith 
62047c6ae99SBarry Smith      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
621aa219208SBarry Smith     DMDA they came from still exists (has not been destroyed).
62247c6ae99SBarry Smith 
623e3f3e4b6SBarry Smith     These numbers are NOT multiplied by the number of dof per node.
624e3f3e4b6SBarry Smith 
625aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
62647c6ae99SBarry Smith @*/
6277087cfbeSBarry Smith PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
62847c6ae99SBarry Smith {
62947c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
63047c6ae99SBarry Smith 
63147c6ae99SBarry Smith   PetscFunctionBegin;
63247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
63347c6ae99SBarry Smith   if (lx) *lx = dd->lx;
63447c6ae99SBarry Smith   if (ly) *ly = dd->ly;
63547c6ae99SBarry Smith   if (lz) *lz = dd->lz;
63647c6ae99SBarry Smith   PetscFunctionReturn(0);
63747c6ae99SBarry Smith }
63847c6ae99SBarry Smith 
63947c6ae99SBarry Smith #undef __FUNCT__
640aa219208SBarry Smith #define __FUNCT__ "DMDASetRefinementFactor"
64147c6ae99SBarry Smith /*@
642aa219208SBarry Smith      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
64347c6ae99SBarry Smith 
644aa219208SBarry Smith     Logically Collective on DMDA
64547c6ae99SBarry Smith 
64647c6ae99SBarry Smith   Input Parameters:
647aa219208SBarry Smith +    da - the DMDA object
64847c6ae99SBarry Smith .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
64947c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
65047c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
65147c6ae99SBarry Smith 
65247c6ae99SBarry Smith   Options Database:
65347c6ae99SBarry Smith +  -da_refine_x - refinement ratio in x direction
65447c6ae99SBarry Smith .  -da_refine_y - refinement ratio in y direction
65547c6ae99SBarry Smith -  -da_refine_z - refinement ratio in z direction
65647c6ae99SBarry Smith 
65747c6ae99SBarry Smith   Level: intermediate
65847c6ae99SBarry Smith 
65947c6ae99SBarry Smith     Notes: Pass PETSC_IGNORE to leave a value unchanged
66047c6ae99SBarry Smith 
661aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor()
66247c6ae99SBarry Smith @*/
6637087cfbeSBarry Smith PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
66447c6ae99SBarry Smith {
66547c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith   PetscFunctionBegin;
66847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
66947c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_x,2);
67047c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_y,3);
67147c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_z,4);
67247c6ae99SBarry Smith 
67347c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
67447c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
67547c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
67647c6ae99SBarry Smith   PetscFunctionReturn(0);
67747c6ae99SBarry Smith }
67847c6ae99SBarry Smith 
67947c6ae99SBarry Smith #undef __FUNCT__
680aa219208SBarry Smith #define __FUNCT__ "DMDAGetRefinementFactor"
68147c6ae99SBarry Smith /*@C
682aa219208SBarry Smith      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
68347c6ae99SBarry Smith 
68447c6ae99SBarry Smith     Not Collective
68547c6ae99SBarry Smith 
68647c6ae99SBarry Smith   Input Parameter:
687aa219208SBarry Smith .    da - the DMDA object
68847c6ae99SBarry Smith 
68947c6ae99SBarry Smith   Output Parameters:
69047c6ae99SBarry Smith +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
69147c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
69247c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
69347c6ae99SBarry Smith 
69447c6ae99SBarry Smith   Level: intermediate
69547c6ae99SBarry Smith 
6960298fd71SBarry Smith     Notes: Pass NULL for values you do not need
69747c6ae99SBarry Smith 
698aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor()
69947c6ae99SBarry Smith @*/
7007087cfbeSBarry Smith PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
70147c6ae99SBarry Smith {
70247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
70347c6ae99SBarry Smith 
70447c6ae99SBarry Smith   PetscFunctionBegin;
70547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
70647c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
70747c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
70847c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
70947c6ae99SBarry Smith   PetscFunctionReturn(0);
71047c6ae99SBarry Smith }
71147c6ae99SBarry Smith 
71247c6ae99SBarry Smith #undef __FUNCT__
713aa219208SBarry Smith #define __FUNCT__ "DMDASetGetMatrix"
71447c6ae99SBarry Smith /*@C
715aa219208SBarry Smith      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
71647c6ae99SBarry Smith 
717aa219208SBarry Smith     Logically Collective on DMDA
71847c6ae99SBarry Smith 
71947c6ae99SBarry Smith   Input Parameters:
720aa219208SBarry Smith +    da - the DMDA object
721aa219208SBarry Smith -    f - the function that allocates the matrix for that specific DMDA
72247c6ae99SBarry Smith 
72347c6ae99SBarry Smith   Level: developer
72447c6ae99SBarry Smith 
725aa219208SBarry Smith    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
72647c6ae99SBarry Smith        the diagonal and off-diagonal blocks of the matrix
72747c6ae99SBarry Smith 
728950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills()
72947c6ae99SBarry Smith @*/
73019fd82e9SBarry Smith PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, MatType,Mat*))
73147c6ae99SBarry Smith {
73247c6ae99SBarry Smith   PetscFunctionBegin;
73347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
73425296bd5SBarry Smith   da->ops->creatematrix = f;
73547c6ae99SBarry Smith   PetscFunctionReturn(0);
73647c6ae99SBarry Smith }
73747c6ae99SBarry Smith 
73847c6ae99SBarry Smith #undef __FUNCT__
73994013140SBarry Smith #define __FUNCT__ "DMDARefineOwnershipRanges"
74047c6ae99SBarry Smith /*
74147c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
74247c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
74347c6ae99SBarry Smith 
74447c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
74547c6ae99SBarry Smith */
74694013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
74747c6ae99SBarry Smith {
74847c6ae99SBarry Smith   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
74947c6ae99SBarry Smith   PetscErrorCode ierr;
75047c6ae99SBarry Smith 
75147c6ae99SBarry Smith   PetscFunctionBegin;
752ce94432eSBarry Smith   if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
75347c6ae99SBarry Smith   if (ratio == 1) {
75447c6ae99SBarry Smith     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
75547c6ae99SBarry Smith     PetscFunctionReturn(0);
75647c6ae99SBarry Smith   }
75747c6ae99SBarry Smith   for (i=0; i<m; i++) totalc += lc[i];
75847c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
75947c6ae99SBarry Smith   for (i=0; i<m; i++) {
76047c6ae99SBarry Smith     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
76147c6ae99SBarry Smith     if (i == m-1) lf[i] = want;
76247c6ae99SBarry Smith     else {
7637aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
7647aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
7657aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
7667aca7175SJed Brown       while ((startf+want)/ratio < nextc - stencil_width) want++;
7677aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
7687aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
7697aca7175SJed Brown       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
7707aca7175SJed Brown       /* Make sure all constraints are satisfied */
77130729d88SBarry Smith       if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width)
772ce94432eSBarry Smith           || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range");
77347c6ae99SBarry Smith     }
77447c6ae99SBarry Smith     lf[i]      = want;
77547c6ae99SBarry Smith     startc    += lc[i];
77647c6ae99SBarry Smith     startf    += lf[i];
77747c6ae99SBarry Smith     remaining -= lf[i];
77847c6ae99SBarry Smith   }
77947c6ae99SBarry Smith   PetscFunctionReturn(0);
78047c6ae99SBarry Smith }
78147c6ae99SBarry Smith 
78247c6ae99SBarry Smith #undef __FUNCT__
7832be375d4SJed Brown #define __FUNCT__ "DMDACoarsenOwnershipRanges"
7842be375d4SJed Brown /*
7852be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
7862be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
7872be375d4SJed Brown 
7882be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
7892be375d4SJed Brown */
7902be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
7912be375d4SJed Brown {
7922be375d4SJed Brown   PetscInt       i,totalf,remaining,startc,startf;
7932be375d4SJed Brown   PetscErrorCode ierr;
7942be375d4SJed Brown 
7952be375d4SJed Brown   PetscFunctionBegin;
796ce94432eSBarry Smith   if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
7972be375d4SJed Brown   if (ratio == 1) {
7982be375d4SJed Brown     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
7992be375d4SJed Brown     PetscFunctionReturn(0);
8002be375d4SJed Brown   }
8012be375d4SJed Brown   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
8022be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
8032be375d4SJed Brown   for (i=0,startc=0,startf=0; i<m; i++) {
8042be375d4SJed Brown     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
8052be375d4SJed Brown     if (i == m-1) lc[i] = want;
8062be375d4SJed Brown     else {
8072be375d4SJed Brown       const PetscInt nextf = startf+lf[i];
8082be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
8092be375d4SJed Brown        * node is within one stencil width. */
8102be375d4SJed Brown       while (nextf/ratio < startc+want-stencil_width) want--;
8112be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
8122be375d4SJed Brown        * fine node is within one stencil width. */
8132be375d4SJed Brown       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
8142be375d4SJed Brown       if (want < 0 || want > remaining
815ce94432eSBarry Smith           || (nextf/ratio < startc+want-stencil_width) || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range");
8162be375d4SJed Brown     }
8172be375d4SJed Brown     lc[i]      = want;
8182be375d4SJed Brown     startc    += lc[i];
8192be375d4SJed Brown     startf    += lf[i];
8202be375d4SJed Brown     remaining -= lc[i];
8212be375d4SJed Brown   }
8222be375d4SJed Brown   PetscFunctionReturn(0);
8232be375d4SJed Brown }
8242be375d4SJed Brown 
8252be375d4SJed Brown #undef __FUNCT__
82694013140SBarry Smith #define __FUNCT__ "DMRefine_DA"
8277087cfbeSBarry Smith PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
82847c6ae99SBarry Smith {
82947c6ae99SBarry Smith   PetscErrorCode ierr;
830f3141302SJed Brown   PetscInt       M,N,P,i;
8319a42bb27SBarry Smith   DM             da2;
83247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
83347c6ae99SBarry Smith 
83447c6ae99SBarry Smith   PetscFunctionBegin;
83547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
83647c6ae99SBarry Smith   PetscValidPointer(daref,3);
83747c6ae99SBarry Smith 
8381321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
83947c6ae99SBarry Smith     M = dd->refine_x*dd->M;
84047c6ae99SBarry Smith   } else {
84147c6ae99SBarry Smith     M = 1 + dd->refine_x*(dd->M - 1);
84247c6ae99SBarry Smith   }
8431321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
8441860e6e9SBarry Smith     if (dd->dim > 1) {
84547c6ae99SBarry Smith       N = dd->refine_y*dd->N;
84647c6ae99SBarry Smith     } else {
8471860e6e9SBarry Smith       N = 1;
8481860e6e9SBarry Smith     }
8491860e6e9SBarry Smith   } else {
85047c6ae99SBarry Smith     N = 1 + dd->refine_y*(dd->N - 1);
85147c6ae99SBarry Smith   }
8521321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
8531860e6e9SBarry Smith     if (dd->dim > 2) {
85447c6ae99SBarry Smith       P = dd->refine_z*dd->P;
85547c6ae99SBarry Smith     } else {
8561860e6e9SBarry Smith       P = 1;
8571860e6e9SBarry Smith     }
8581860e6e9SBarry Smith   } else {
85947c6ae99SBarry Smith     P = 1 + dd->refine_z*(dd->P - 1);
86047c6ae99SBarry Smith   }
861ce94432eSBarry Smith   ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr);
862397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
863397b6216SJed Brown   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
864397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
865397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
866397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
867397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
868397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
869397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
87047c6ae99SBarry Smith   if (dd->dim == 3) {
87147c6ae99SBarry Smith     PetscInt *lx,*ly,*lz;
87247c6ae99SBarry Smith     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
8731321219cSEthan 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);
8741321219cSEthan 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);
8751321219cSEthan 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);
876397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
87747c6ae99SBarry Smith     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
87847c6ae99SBarry Smith   } else if (dd->dim == 2) {
87947c6ae99SBarry Smith     PetscInt *lx,*ly;
88047c6ae99SBarry Smith     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
8811321219cSEthan 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);
8821321219cSEthan 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);
8830298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr);
88447c6ae99SBarry Smith     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
88547c6ae99SBarry Smith   } else if (dd->dim == 1) {
88647c6ae99SBarry Smith     PetscInt *lx;
88747c6ae99SBarry Smith     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
8881321219cSEthan 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);
8890298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr);
89047c6ae99SBarry Smith     ierr = PetscFree(lx);CHKERRQ(ierr);
89147c6ae99SBarry Smith   }
89247c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
89347c6ae99SBarry Smith 
89447c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
89525296bd5SBarry Smith   da2->ops->creatematrix = da->ops->creatematrix;
89625296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
89747c6ae99SBarry Smith   da2->ops->getcoloring = da->ops->getcoloring;
89847c6ae99SBarry Smith   dd2->interptype       = dd->interptype;
89947c6ae99SBarry Smith 
90047c6ae99SBarry Smith   /* copy fill information if given */
90147c6ae99SBarry Smith   if (dd->dfill) {
90247c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
90347c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
90447c6ae99SBarry Smith   }
90547c6ae99SBarry Smith   if (dd->ofill) {
90647c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
90747c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
90847c6ae99SBarry Smith   }
90947c6ae99SBarry Smith   /* copy the refine information */
910397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
911397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
912397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
91347c6ae99SBarry Smith 
91447c6ae99SBarry Smith   /* copy vector type information */
91547c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
91619fd82e9SBarry Smith   ierr = PetscStrallocpy(da->vectype,(char**)&da2->vectype);CHKERRQ(ierr);
917ddcf8b74SDave May 
918efd51863SBarry Smith   dd2->lf = dd->lf;
919efd51863SBarry Smith   dd2->lj = dd->lj;
920efd51863SBarry Smith 
9216e87535bSJed Brown   da2->leveldown = da->leveldown;
9226e87535bSJed Brown   da2->levelup   = da->levelup + 1;
9238865f1eaSKarl Rupp 
9246e87535bSJed Brown   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
9256e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
926ca266f36SBarry Smith   ierr = DMViewFromOptions(da2,"-dm_view");CHKERRQ(ierr);
9276e87535bSJed Brown 
928ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
9296636e97aSMatthew G Knepley   if (da->coordinates) {
930ddcf8b74SDave May     DM  cdaf,cdac;
931ddcf8b74SDave May     Vec coordsc,coordsf;
932ddcf8b74SDave May     Mat II;
933ddcf8b74SDave May 
9346636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr);
9356636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr);
9366636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr);
937b61d3410SDave May     /* force creation of the coordinate vector */
938b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
9396636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
9400298fd71SBarry Smith     ierr = DMCreateInterpolation(cdac,cdaf,&II,NULL);CHKERRQ(ierr);
941ddcf8b74SDave May     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
942fcfd50ebSBarry Smith     ierr = MatDestroy(&II);CHKERRQ(ierr);
943ddcf8b74SDave May   }
944397b6216SJed Brown 
945f3141302SJed Brown   for (i=0; i<da->bs; i++) {
946f3141302SJed Brown     const char *fieldname;
947f3141302SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
948f3141302SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
949f3141302SJed Brown   }
950397b6216SJed Brown 
95147c6ae99SBarry Smith   *daref = da2;
95247c6ae99SBarry Smith   PetscFunctionReturn(0);
95347c6ae99SBarry Smith }
95447c6ae99SBarry Smith 
955397b6216SJed Brown 
95647c6ae99SBarry Smith #undef __FUNCT__
95794013140SBarry Smith #define __FUNCT__ "DMCoarsen_DA"
9587087cfbeSBarry Smith PetscErrorCode  DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
95947c6ae99SBarry Smith {
96047c6ae99SBarry Smith   PetscErrorCode ierr;
961397b6216SJed Brown   PetscInt       M,N,P,i;
9629a42bb27SBarry Smith   DM             da2;
96347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
96447c6ae99SBarry Smith 
96547c6ae99SBarry Smith   PetscFunctionBegin;
96647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
96747c6ae99SBarry Smith   PetscValidPointer(daref,3);
96847c6ae99SBarry Smith 
9691321219cSEthan Coon   if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
970397b6216SJed Brown     M = dd->M / dd->coarsen_x;
97147c6ae99SBarry Smith   } else {
972397b6216SJed Brown     M = 1 + (dd->M - 1) / dd->coarsen_x;
97347c6ae99SBarry Smith   }
9741321219cSEthan Coon   if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
9751860e6e9SBarry Smith     if (dd->dim > 1) {
976397b6216SJed Brown       N = dd->N / dd->coarsen_y;
97747c6ae99SBarry Smith     } else {
9781860e6e9SBarry Smith       N = 1;
9791860e6e9SBarry Smith     }
9801860e6e9SBarry Smith   } else {
981397b6216SJed Brown     N = 1 + (dd->N - 1) / dd->coarsen_y;
98247c6ae99SBarry Smith   }
9831321219cSEthan Coon   if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
9841860e6e9SBarry Smith     if (dd->dim > 2) {
985397b6216SJed Brown       P = dd->P / dd->coarsen_z;
98647c6ae99SBarry Smith     } else {
9871860e6e9SBarry Smith       P = 1;
9881860e6e9SBarry Smith     }
9891860e6e9SBarry Smith   } else {
990397b6216SJed Brown     P = 1 + (dd->P - 1) / dd->coarsen_z;
99147c6ae99SBarry Smith   }
992ce94432eSBarry Smith   ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr);
993397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
994397b6216SJed Brown   ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr);
995397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
996397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
997397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
998397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
999397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
1000397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
100147c6ae99SBarry Smith   if (dd->dim == 3) {
10022be375d4SJed Brown     PetscInt *lx,*ly,*lz;
10032be375d4SJed Brown     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
1004397b6216SJed 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);
1005397b6216SJed 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);
1006397b6216SJed 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);
1007397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
10082be375d4SJed Brown     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
100947c6ae99SBarry Smith   } else if (dd->dim == 2) {
10102be375d4SJed Brown     PetscInt *lx,*ly;
10112be375d4SJed Brown     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
1012397b6216SJed 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);
1013397b6216SJed 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);
10140298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr);
10152be375d4SJed Brown     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
101647c6ae99SBarry Smith   } else if (dd->dim == 1) {
10172be375d4SJed Brown     PetscInt *lx;
10182be375d4SJed Brown     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
1019397b6216SJed 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);
10200298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr);
10212be375d4SJed Brown     ierr = PetscFree(lx);CHKERRQ(ierr);
102247c6ae99SBarry Smith   }
102347c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
102447c6ae99SBarry Smith 
10254dcab191SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
102625296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */
102725296bd5SBarry Smith   da2->ops->creatematrix = da->ops->creatematrix;
102847c6ae99SBarry Smith   da2->ops->getcoloring  = da->ops->getcoloring;
102947c6ae99SBarry Smith   dd2->interptype        = dd->interptype;
103047c6ae99SBarry Smith 
103147c6ae99SBarry Smith   /* copy fill information if given */
103247c6ae99SBarry Smith   if (dd->dfill) {
103347c6ae99SBarry Smith     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
103447c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
103547c6ae99SBarry Smith   }
103647c6ae99SBarry Smith   if (dd->ofill) {
103747c6ae99SBarry Smith     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
103847c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
103947c6ae99SBarry Smith   }
104047c6ae99SBarry Smith   /* copy the refine information */
1041397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
1042397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
1043397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
104447c6ae99SBarry Smith 
104547c6ae99SBarry Smith   /* copy vector type information */
104647c6ae99SBarry Smith   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
104719fd82e9SBarry Smith   ierr = PetscStrallocpy(da->vectype,(char**)&da2->vectype);CHKERRQ(ierr);
104847c6ae99SBarry Smith 
1049644e2e5bSBarry Smith   dd2->lf = dd->lf;
1050644e2e5bSBarry Smith   dd2->lj = dd->lj;
1051644e2e5bSBarry Smith 
10526e87535bSJed Brown   da2->leveldown = da->leveldown + 1;
10536e87535bSJed Brown   da2->levelup   = da->levelup;
10548865f1eaSKarl Rupp 
10556e87535bSJed Brown   ierr = DMSetFromOptions(da2);CHKERRQ(ierr);
10566e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
1057ca266f36SBarry Smith   ierr = DMViewFromOptions(da2,"-dm_view");CHKERRQ(ierr);
10586e87535bSJed Brown 
1059ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
10606636e97aSMatthew G Knepley   if (da->coordinates) {
1061ddcf8b74SDave May     DM         cdaf,cdac;
1062ddcf8b74SDave May     Vec        coordsc,coordsf;
1063ddcf8b74SDave May     VecScatter inject;
1064ddcf8b74SDave May 
10656636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr);
10666636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr);
10676636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr);
1068b61d3410SDave May     /* force creation of the coordinate vector */
1069b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
10706636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
1071ddcf8b74SDave May 
1072e727c939SJed Brown     ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
1073ddcf8b74SDave May     ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1074ddcf8b74SDave May     ierr = VecScatterEnd(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1075fcfd50ebSBarry Smith     ierr = VecScatterDestroy(&inject);CHKERRQ(ierr);
1076ddcf8b74SDave May   }
1077f98405f7SJed Brown 
1078f98405f7SJed Brown   for (i=0; i<da->bs; i++) {
1079f98405f7SJed Brown     const char *fieldname;
1080f98405f7SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
1081f98405f7SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
1082f98405f7SJed Brown   }
1083f98405f7SJed Brown 
108447c6ae99SBarry Smith   *daref = da2;
108547c6ae99SBarry Smith   PetscFunctionReturn(0);
108647c6ae99SBarry Smith }
108747c6ae99SBarry Smith 
108847c6ae99SBarry Smith #undef __FUNCT__
108994013140SBarry Smith #define __FUNCT__ "DMRefineHierarchy_DA"
10907087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
109147c6ae99SBarry Smith {
109247c6ae99SBarry Smith   PetscErrorCode ierr;
109347c6ae99SBarry Smith   PetscInt       i,n,*refx,*refy,*refz;
109447c6ae99SBarry Smith 
109547c6ae99SBarry Smith   PetscFunctionBegin;
109647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1097ce94432eSBarry Smith   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
109847c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
109947c6ae99SBarry Smith   PetscValidPointer(daf,3);
110047c6ae99SBarry Smith 
1101aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
110247c6ae99SBarry Smith   ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr);
110347c6ae99SBarry Smith   for (i=0; i<nlevels; i++) {
1104aa219208SBarry Smith     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
110547c6ae99SBarry Smith   }
110647c6ae99SBarry Smith   n    = nlevels;
11070298fd71SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr);
110847c6ae99SBarry Smith   n    = nlevels;
11090298fd71SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr);
111047c6ae99SBarry Smith   n    = nlevels;
11110298fd71SBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr);
111247c6ae99SBarry Smith 
1113aa219208SBarry Smith   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
1114ce94432eSBarry Smith   ierr = DMRefine(da,PetscObjectComm((PetscObject)da),&daf[0]);CHKERRQ(ierr);
111547c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
1116aa219208SBarry Smith     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
1117ce94432eSBarry Smith     ierr = DMRefine(daf[i-1],PetscObjectComm((PetscObject)da),&daf[i]);CHKERRQ(ierr);
111847c6ae99SBarry Smith   }
111947c6ae99SBarry Smith   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
112047c6ae99SBarry Smith   PetscFunctionReturn(0);
112147c6ae99SBarry Smith }
112247c6ae99SBarry Smith 
112347c6ae99SBarry Smith #undef __FUNCT__
112494013140SBarry Smith #define __FUNCT__ "DMCoarsenHierarchy_DA"
11257087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
112647c6ae99SBarry Smith {
112747c6ae99SBarry Smith   PetscErrorCode ierr;
112847c6ae99SBarry Smith   PetscInt       i;
112947c6ae99SBarry Smith 
113047c6ae99SBarry Smith   PetscFunctionBegin;
113147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1132ce94432eSBarry Smith   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
113347c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
113447c6ae99SBarry Smith   PetscValidPointer(dac,3);
1135ce94432eSBarry Smith   ierr = DMCoarsen(da,PetscObjectComm((PetscObject)da),&dac[0]);CHKERRQ(ierr);
113647c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
1137ce94432eSBarry Smith     ierr = DMCoarsen(dac[i-1],PetscObjectComm((PetscObject)da),&dac[i]);CHKERRQ(ierr);
113847c6ae99SBarry Smith   }
113947c6ae99SBarry Smith   PetscFunctionReturn(0);
114047c6ae99SBarry Smith }
1141