xref: /petsc/src/dm/impls/da/da.c (revision 110623a03a130d18890d33f7df13a2b3e2bded77)
1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
247c6ae99SBarry Smith 
347c6ae99SBarry Smith /*@
4e43dc8daSBarry Smith   DMDASetSizes - Sets the number of grid points in the three dimensional directions
547c6ae99SBarry Smith 
6aa219208SBarry Smith   Logically Collective on DMDA
747c6ae99SBarry Smith 
847c6ae99SBarry Smith   Input Parameters:
9aa219208SBarry Smith + da - the DMDA
10e43dc8daSBarry Smith . M - the global X size
11e43dc8daSBarry Smith . N - the global Y size
12e43dc8daSBarry Smith - P - the global Z size
1347c6ae99SBarry Smith 
1447c6ae99SBarry Smith   Level: intermediate
1547c6ae99SBarry Smith 
16e43dc8daSBarry Smith   Developer Notes: Since the dimension may not yet have been set the code cannot error check for non-positive Y and Z number of grid points
17e43dc8daSBarry Smith 
18aa219208SBarry Smith .seealso: DMDAGetSize(), PetscSplitOwnership()
1947c6ae99SBarry Smith @*/
207087cfbeSBarry Smith PetscErrorCode  DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P)
2147c6ae99SBarry Smith {
2247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   PetscFunctionBegin;
25a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
2647c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,M,2);
2747c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,N,3);
2847c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,P,4);
29ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
30e43dc8daSBarry Smith   if (M < 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in X direction must be positive");
31e43dc8daSBarry Smith   if (N < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in Y direction must be positive");
32e43dc8daSBarry Smith   if (P < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in Z direction must be positive");
3347c6ae99SBarry Smith 
3447c6ae99SBarry Smith   dd->M = M;
3547c6ae99SBarry Smith   dd->N = N;
3647c6ae99SBarry Smith   dd->P = P;
3747c6ae99SBarry Smith   PetscFunctionReturn(0);
3847c6ae99SBarry Smith }
3947c6ae99SBarry Smith 
4047c6ae99SBarry Smith /*@
41aa219208SBarry Smith   DMDASetNumProcs - Sets the number of processes in each dimension
4247c6ae99SBarry Smith 
43aa219208SBarry Smith   Logically Collective on DMDA
4447c6ae99SBarry Smith 
4547c6ae99SBarry Smith   Input Parameters:
46aa219208SBarry Smith + da - the DMDA
4747c6ae99SBarry Smith . m - the number of X procs (or PETSC_DECIDE)
4847c6ae99SBarry Smith . n - the number of Y procs (or PETSC_DECIDE)
4947c6ae99SBarry Smith - p - the number of Z procs (or PETSC_DECIDE)
5047c6ae99SBarry Smith 
5147c6ae99SBarry Smith   Level: intermediate
5247c6ae99SBarry Smith 
53aa219208SBarry Smith .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership()
5447c6ae99SBarry Smith @*/
557087cfbeSBarry Smith PetscErrorCode  DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p)
5647c6ae99SBarry Smith {
5747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
58e3f3e4b6SBarry Smith   PetscErrorCode ierr;
5947c6ae99SBarry Smith 
6047c6ae99SBarry Smith   PetscFunctionBegin;
61a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
6247c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,m,2);
6347c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,n,3);
6447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,p,4);
65ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
6647c6ae99SBarry Smith   dd->m = m;
6747c6ae99SBarry Smith   dd->n = n;
6847c6ae99SBarry Smith   dd->p = p;
69c73cfb54SMatthew G. Knepley   if (da->dim == 2) {
70d3be247eSBarry Smith     PetscMPIInt size;
71d3be247eSBarry Smith     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr);
72e3f3e4b6SBarry Smith     if ((dd->m > 0) && (dd->n < 0)) {
73e3f3e4b6SBarry Smith       dd->n = size/dd->m;
74ce94432eSBarry 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);
75e3f3e4b6SBarry Smith     }
76e3f3e4b6SBarry Smith     if ((dd->n > 0) && (dd->m < 0)) {
77e3f3e4b6SBarry Smith       dd->m = size/dd->n;
78ce94432eSBarry 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);
79e3f3e4b6SBarry Smith     }
80e3f3e4b6SBarry Smith   }
8147c6ae99SBarry Smith   PetscFunctionReturn(0);
8247c6ae99SBarry Smith }
8347c6ae99SBarry Smith 
8447c6ae99SBarry Smith /*@
851321219cSEthan Coon   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries.
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith   Not collective
8847c6ae99SBarry Smith 
8947c6ae99SBarry Smith   Input Parameter:
90aa219208SBarry Smith + da    - The DMDA
91bff4a2f0SMatthew G. Knepley - bx,by,bz - One of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC
9247c6ae99SBarry Smith 
9347c6ae99SBarry Smith   Level: intermediate
9447c6ae99SBarry Smith 
9547c6ae99SBarry Smith .keywords:  distributed array, periodicity
96bff4a2f0SMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA, DMBoundaryType
9747c6ae99SBarry Smith @*/
98bff4a2f0SMatthew G. Knepley PetscErrorCode  DMDASetBoundaryType(DM da,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz)
9947c6ae99SBarry Smith {
10047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith   PetscFunctionBegin;
103a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
1045a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,bx,2);
1055a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,by,3);
1065a43f728SLisandro Dalcin   PetscValidLogicalCollectiveEnum(da,bz,4);
107ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
1081321219cSEthan Coon   dd->bx = bx;
1091321219cSEthan Coon   dd->by = by;
1101321219cSEthan Coon   dd->bz = bz;
11147c6ae99SBarry Smith   PetscFunctionReturn(0);
11247c6ae99SBarry Smith }
11347c6ae99SBarry Smith 
11447c6ae99SBarry Smith /*@
115aa219208SBarry Smith   DMDASetDof - Sets the number of degrees of freedom per vertex
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith   Not collective
11847c6ae99SBarry Smith 
11959f3ab6dSMatthew G. Knepley   Input Parameters:
120aa219208SBarry Smith + da  - The DMDA
12147c6ae99SBarry Smith - dof - Number of degrees of freedom
12247c6ae99SBarry Smith 
12347c6ae99SBarry Smith   Level: intermediate
12447c6ae99SBarry Smith 
12547c6ae99SBarry Smith .keywords:  distributed array, degrees of freedom
126fb6725baSMatthew G. Knepley .seealso: DMDAGetDof(), DMDACreate(), DMDestroy(), DMDA
12747c6ae99SBarry Smith @*/
12854cfb0beSLisandro Dalcin PetscErrorCode  DMDASetDof(DM da, PetscInt dof)
12947c6ae99SBarry Smith {
13047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith   PetscFunctionBegin;
133a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
13454cfb0beSLisandro Dalcin   PetscValidLogicalCollectiveInt(da,dof,2);
135ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
13647c6ae99SBarry Smith   dd->w  = dof;
1371411c6eeSJed Brown   da->bs = dof;
13847c6ae99SBarry Smith   PetscFunctionReturn(0);
13947c6ae99SBarry Smith }
14047c6ae99SBarry Smith 
141fb6725baSMatthew G. Knepley /*@
142fb6725baSMatthew G. Knepley   DMDAGetDof - Gets the number of degrees of freedom per vertex
143fb6725baSMatthew G. Knepley 
144fb6725baSMatthew G. Knepley   Not collective
145fb6725baSMatthew G. Knepley 
146fb6725baSMatthew G. Knepley   Input Parameter:
147fb6725baSMatthew G. Knepley . da  - The DMDA
148fb6725baSMatthew G. Knepley 
149fb6725baSMatthew G. Knepley   Output Parameter:
150fb6725baSMatthew G. Knepley . dof - Number of degrees of freedom
151fb6725baSMatthew G. Knepley 
152fb6725baSMatthew G. Knepley   Level: intermediate
153fb6725baSMatthew G. Knepley 
154fb6725baSMatthew G. Knepley .keywords:  distributed array, degrees of freedom
155fb6725baSMatthew G. Knepley .seealso: DMDASetDof(), DMDACreate(), DMDestroy(), DMDA
156fb6725baSMatthew G. Knepley @*/
157fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetDof(DM da, PetscInt *dof)
158fb6725baSMatthew G. Knepley {
159fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *) da->data;
160fb6725baSMatthew G. Knepley 
161fb6725baSMatthew G. Knepley   PetscFunctionBegin;
162a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
163fb6725baSMatthew G. Knepley   PetscValidPointer(dof,2);
164fb6725baSMatthew G. Knepley   *dof = dd->w;
165fb6725baSMatthew G. Knepley   PetscFunctionReturn(0);
166fb6725baSMatthew G. Knepley }
167fb6725baSMatthew G. Knepley 
1687ddda789SPeter Brune /*@
1697ddda789SPeter Brune   DMDAGetOverlap - Gets the size of the per-processor overlap.
1707ddda789SPeter Brune 
1717ddda789SPeter Brune   Not collective
1727ddda789SPeter Brune 
1737ddda789SPeter Brune   Input Parameters:
1747ddda789SPeter Brune . da  - The DMDA
1757ddda789SPeter Brune 
1767ddda789SPeter Brune   Output Parameters:
1777ddda789SPeter Brune + x   - Overlap in the x direction
1787ddda789SPeter Brune . y   - Overlap in the y direction
1797ddda789SPeter Brune - z   - Overlap in the z direction
1807ddda789SPeter Brune 
1817ddda789SPeter Brune   Level: intermediate
1827ddda789SPeter Brune 
1837ddda789SPeter Brune .keywords:  distributed array, overlap, domain decomposition
1847ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetOverlap(), DMDA
1857ddda789SPeter Brune @*/
1867ddda789SPeter Brune PetscErrorCode  DMDAGetOverlap(DM da,PetscInt *x,PetscInt *y,PetscInt *z)
1877ddda789SPeter Brune {
1887ddda789SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
1897ddda789SPeter Brune 
1907ddda789SPeter Brune   PetscFunctionBegin;
191a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
1927ddda789SPeter Brune   if (x) *x = dd->xol;
1937ddda789SPeter Brune   if (y) *y = dd->yol;
1947ddda789SPeter Brune   if (z) *z = dd->zol;
1957ddda789SPeter Brune   PetscFunctionReturn(0);
1967ddda789SPeter Brune }
1977ddda789SPeter Brune 
19888661749SPeter Brune /*@
19988661749SPeter Brune   DMDASetOverlap - Sets the size of the per-processor overlap.
20088661749SPeter Brune 
20188661749SPeter Brune   Not collective
20288661749SPeter Brune 
2037ddda789SPeter Brune   Input Parameters:
20488661749SPeter Brune + da  - The DMDA
2057ddda789SPeter Brune . x   - Overlap in the x direction
2067ddda789SPeter Brune . y   - Overlap in the y direction
2077ddda789SPeter Brune - z   - Overlap in the z direction
20888661749SPeter Brune 
20988661749SPeter Brune   Level: intermediate
21088661749SPeter Brune 
2117ddda789SPeter Brune .keywords:  distributed array, overlap, domain decomposition
2127ddda789SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetOverlap(), DMDA
21388661749SPeter Brune @*/
2147ddda789SPeter Brune PetscErrorCode  DMDASetOverlap(DM da,PetscInt x,PetscInt y,PetscInt z)
21588661749SPeter Brune {
21688661749SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
21788661749SPeter Brune 
21888661749SPeter Brune   PetscFunctionBegin;
219a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
2207ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da,x,2);
2217ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da,y,3);
2227ddda789SPeter Brune   PetscValidLogicalCollectiveInt(da,z,4);
2237ddda789SPeter Brune   dd->xol = x;
2247ddda789SPeter Brune   dd->yol = y;
2257ddda789SPeter Brune   dd->zol = z;
22688661749SPeter Brune   PetscFunctionReturn(0);
22788661749SPeter Brune }
22888661749SPeter Brune 
2293e7870d2SPeter Brune 
2303e7870d2SPeter Brune /*@
2313e7870d2SPeter Brune   DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition.
2323e7870d2SPeter Brune 
2333e7870d2SPeter Brune   Not collective
2343e7870d2SPeter Brune 
2353e7870d2SPeter Brune   Input Parameters:
2363e7870d2SPeter Brune . da  - The DMDA
2373e7870d2SPeter Brune 
2383e7870d2SPeter Brune   Output Parameters:
2393e7870d2SPeter Brune + Nsub   - Number of local subdomains created upon decomposition
2403e7870d2SPeter Brune 
2413e7870d2SPeter Brune   Level: intermediate
2423e7870d2SPeter Brune 
2433e7870d2SPeter Brune .keywords:  distributed array, domain decomposition
2443e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDASetNumLocalSubDomains(), DMDA
2453e7870d2SPeter Brune @*/
2463e7870d2SPeter Brune PetscErrorCode  DMDAGetNumLocalSubDomains(DM da,PetscInt *Nsub)
2473e7870d2SPeter Brune {
2483e7870d2SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
2493e7870d2SPeter Brune 
2503e7870d2SPeter Brune   PetscFunctionBegin;
251a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
2523e7870d2SPeter Brune   if (Nsub) *Nsub = dd->Nsub;
2533e7870d2SPeter Brune   PetscFunctionReturn(0);
2543e7870d2SPeter Brune }
2553e7870d2SPeter Brune 
2563e7870d2SPeter Brune /*@
2573e7870d2SPeter Brune   DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition.
2583e7870d2SPeter Brune 
2593e7870d2SPeter Brune   Not collective
2603e7870d2SPeter Brune 
2613e7870d2SPeter Brune   Input Parameters:
2623e7870d2SPeter Brune + da  - The DMDA
2633e7870d2SPeter Brune - Nsub - The number of local subdomains requested
2643e7870d2SPeter Brune 
2653e7870d2SPeter Brune   Level: intermediate
2663e7870d2SPeter Brune 
2673e7870d2SPeter Brune .keywords:  distributed array, domain decomposition
2683e7870d2SPeter Brune .seealso: DMDACreateDomainDecomposition(), DMDAGetNumLocalSubDomains(), DMDA
2693e7870d2SPeter Brune @*/
2703e7870d2SPeter Brune PetscErrorCode  DMDASetNumLocalSubDomains(DM da,PetscInt Nsub)
2713e7870d2SPeter Brune {
2723e7870d2SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
2733e7870d2SPeter Brune 
2743e7870d2SPeter Brune   PetscFunctionBegin;
275a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
2763e7870d2SPeter Brune   PetscValidLogicalCollectiveInt(da,Nsub,2);
2773e7870d2SPeter Brune   dd->Nsub = Nsub;
2783e7870d2SPeter Brune   PetscFunctionReturn(0);
2793e7870d2SPeter Brune }
2803e7870d2SPeter Brune 
281d886c4f4SPeter Brune /*@
282d886c4f4SPeter Brune   DMDASetOffset - Sets the index offset of the DA.
283d886c4f4SPeter Brune 
284d886c4f4SPeter Brune   Collective on DA
285d886c4f4SPeter Brune 
286d886c4f4SPeter Brune   Input Parameter:
287d886c4f4SPeter Brune + da  - The DMDA
288d886c4f4SPeter Brune . xo  - The offset in the x direction
289d886c4f4SPeter Brune . yo  - The offset in the y direction
290d886c4f4SPeter Brune - zo  - The offset in the z direction
291d886c4f4SPeter Brune 
292d886c4f4SPeter Brune   Level: intermediate
293d886c4f4SPeter Brune 
29495452b02SPatrick Sanan   Notes:
29595452b02SPatrick Sanan     This is used primarily to overlap a computation on a local DA with that on a global DA without
296d886c4f4SPeter Brune   changing boundary conditions or subdomain features that depend upon the global offsets.
297d886c4f4SPeter Brune 
298d886c4f4SPeter Brune .keywords:  distributed array, degrees of freedom
299d886c4f4SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray()
300d886c4f4SPeter Brune @*/
30195c13181SPeter Brune PetscErrorCode  DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po)
302d886c4f4SPeter Brune {
30395c13181SPeter Brune   PetscErrorCode ierr;
304d886c4f4SPeter Brune   DM_DA          *dd = (DM_DA*)da->data;
305d886c4f4SPeter Brune 
306d886c4f4SPeter Brune   PetscFunctionBegin;
307a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
308d886c4f4SPeter Brune   PetscValidLogicalCollectiveInt(da,xo,2);
30995c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,yo,3);
31095c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,zo,4);
31195c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,Mo,5);
31295c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,No,6);
31395c13181SPeter Brune   PetscValidLogicalCollectiveInt(da,Po,7);
314d886c4f4SPeter Brune   dd->xo = xo;
315d886c4f4SPeter Brune   dd->yo = yo;
316d886c4f4SPeter Brune   dd->zo = zo;
31795c13181SPeter Brune   dd->Mo = Mo;
31895c13181SPeter Brune   dd->No = No;
31995c13181SPeter Brune   dd->Po = Po;
32095c13181SPeter Brune 
32195c13181SPeter Brune   if (da->coordinateDM) {
32295c13181SPeter Brune     ierr = DMDASetOffset(da->coordinateDM,xo,yo,zo,Mo,No,Po);CHKERRQ(ierr);
32395c13181SPeter Brune   }
324d886c4f4SPeter Brune   PetscFunctionReturn(0);
325d886c4f4SPeter Brune }
326d886c4f4SPeter Brune 
327d886c4f4SPeter Brune /*@
328d886c4f4SPeter Brune   DMDAGetOffset - Gets the index offset of the DA.
329d886c4f4SPeter Brune 
330d886c4f4SPeter Brune   Not collective
331d886c4f4SPeter Brune 
332d886c4f4SPeter Brune   Input Parameter:
333d886c4f4SPeter Brune . da  - The DMDA
334d886c4f4SPeter Brune 
335d886c4f4SPeter Brune   Output Parameters:
336d886c4f4SPeter Brune + xo  - The offset in the x direction
337d886c4f4SPeter Brune . yo  - The offset in the y direction
33895c13181SPeter Brune . zo  - The offset in the z direction
33995c13181SPeter Brune . Mo  - The global size in the x direction
34095c13181SPeter Brune . No  - The global size in the y direction
34195c13181SPeter Brune - Po  - The global size in the z direction
342d886c4f4SPeter Brune 
343d886c4f4SPeter Brune   Level: intermediate
344d886c4f4SPeter Brune 
345d886c4f4SPeter Brune .keywords:  distributed array, degrees of freedom
346d886c4f4SPeter Brune .seealso: DMDASetOffset(), DMDAVecGetArray()
347d886c4f4SPeter Brune @*/
34895c13181SPeter Brune PetscErrorCode  DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po)
349d886c4f4SPeter Brune {
350d886c4f4SPeter Brune   DM_DA *dd = (DM_DA*)da->data;
351d886c4f4SPeter Brune 
352d886c4f4SPeter Brune   PetscFunctionBegin;
353a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
354d886c4f4SPeter Brune   if (xo) *xo = dd->xo;
355d886c4f4SPeter Brune   if (yo) *yo = dd->yo;
356d886c4f4SPeter Brune   if (zo) *zo = dd->zo;
35795c13181SPeter Brune   if (Mo) *Mo = dd->Mo;
35895c13181SPeter Brune   if (No) *No = dd->No;
35995c13181SPeter Brune   if (Po) *Po = dd->Po;
360d886c4f4SPeter Brune   PetscFunctionReturn(0);
361d886c4f4SPeter Brune }
362d886c4f4SPeter Brune 
36340234c92SPeter Brune /*@
36440234c92SPeter Brune   DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain DM.
36540234c92SPeter Brune 
36640234c92SPeter Brune   Not collective
36740234c92SPeter Brune 
36840234c92SPeter Brune   Input Parameter:
36940234c92SPeter Brune . da  - The DMDA
37040234c92SPeter Brune 
37140234c92SPeter Brune   Output Parameters:
37240234c92SPeter Brune + xs  - The start of the region in x
37340234c92SPeter Brune . ys  - The start of the region in y
37440234c92SPeter Brune . zs  - The start of the region in z
37540234c92SPeter Brune . xs  - The size of the region in x
37640234c92SPeter Brune . ys  - The size of the region in y
37740234c92SPeter Brune . zs  - The size of the region in z
37840234c92SPeter Brune 
37940234c92SPeter Brune   Level: intermediate
38040234c92SPeter Brune 
38140234c92SPeter Brune .keywords:  distributed array, degrees of freedom
38240234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray()
38340234c92SPeter Brune @*/
38440234c92SPeter Brune PetscErrorCode  DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm)
38540234c92SPeter Brune {
38640234c92SPeter Brune   DM_DA          *dd = (DM_DA*)da->data;
38740234c92SPeter Brune 
38840234c92SPeter Brune   PetscFunctionBegin;
389a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
39040234c92SPeter Brune   if (xs) *xs = dd->nonxs;
39140234c92SPeter Brune   if (ys) *ys = dd->nonys;
39240234c92SPeter Brune   if (zs) *zs = dd->nonzs;
39340234c92SPeter Brune   if (xm) *xm = dd->nonxm;
39440234c92SPeter Brune   if (ym) *ym = dd->nonym;
39540234c92SPeter Brune   if (zm) *zm = dd->nonzm;
39640234c92SPeter Brune   PetscFunctionReturn(0);
39740234c92SPeter Brune }
39840234c92SPeter Brune 
39940234c92SPeter Brune 
40040234c92SPeter Brune /*@
40140234c92SPeter Brune   DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain DM.
40240234c92SPeter Brune 
40340234c92SPeter Brune   Collective on DA
40440234c92SPeter Brune 
40540234c92SPeter Brune   Input Parameter:
40640234c92SPeter Brune + da  - The DMDA
40740234c92SPeter Brune . xs  - The start of the region in x
40840234c92SPeter Brune . ys  - The start of the region in y
40940234c92SPeter Brune . zs  - The start of the region in z
41040234c92SPeter Brune . xs  - The size of the region in x
41140234c92SPeter Brune . ys  - The size of the region in y
41240234c92SPeter Brune . zs  - The size of the region in z
41340234c92SPeter Brune 
41440234c92SPeter Brune   Level: intermediate
41540234c92SPeter Brune 
41640234c92SPeter Brune .keywords:  distributed array, degrees of freedom
41740234c92SPeter Brune .seealso: DMDAGetOffset(), DMDAVecGetArray()
41840234c92SPeter Brune @*/
41940234c92SPeter Brune PetscErrorCode  DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm)
42040234c92SPeter Brune {
42140234c92SPeter Brune   DM_DA          *dd = (DM_DA*)da->data;
42240234c92SPeter Brune 
42340234c92SPeter Brune   PetscFunctionBegin;
424a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
42540234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,xs,2);
42640234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,ys,3);
42740234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,zs,4);
42840234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,xm,5);
42940234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,ym,6);
43040234c92SPeter Brune   PetscValidLogicalCollectiveInt(da,zm,7);
43140234c92SPeter Brune   dd->nonxs = xs;
43240234c92SPeter Brune   dd->nonys = ys;
43340234c92SPeter Brune   dd->nonzs = zs;
43440234c92SPeter Brune   dd->nonxm = xm;
43540234c92SPeter Brune   dd->nonym = ym;
43640234c92SPeter Brune   dd->nonzm = zm;
43740234c92SPeter Brune 
43840234c92SPeter Brune   PetscFunctionReturn(0);
43940234c92SPeter Brune }
44088661749SPeter Brune 
44147c6ae99SBarry Smith /*@
442aa219208SBarry Smith   DMDASetStencilType - Sets the type of the communication stencil
44347c6ae99SBarry Smith 
444aa219208SBarry Smith   Logically Collective on DMDA
44547c6ae99SBarry Smith 
44647c6ae99SBarry Smith   Input Parameter:
447aa219208SBarry Smith + da    - The DMDA
448aa219208SBarry Smith - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
44947c6ae99SBarry Smith 
45047c6ae99SBarry Smith   Level: intermediate
45147c6ae99SBarry Smith 
45247c6ae99SBarry Smith .keywords:  distributed array, stencil
45396e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
45447c6ae99SBarry Smith @*/
4557087cfbeSBarry Smith PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
45647c6ae99SBarry Smith {
45747c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
45847c6ae99SBarry Smith 
45947c6ae99SBarry Smith   PetscFunctionBegin;
460a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
46147c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,stype,2);
462ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
46347c6ae99SBarry Smith   dd->stencil_type = stype;
46447c6ae99SBarry Smith   PetscFunctionReturn(0);
46547c6ae99SBarry Smith }
46647c6ae99SBarry Smith 
467fb6725baSMatthew G. Knepley /*@
468fb6725baSMatthew G. Knepley   DMDAGetStencilType - Gets the type of the communication stencil
469fb6725baSMatthew G. Knepley 
470fb6725baSMatthew G. Knepley   Not collective
471fb6725baSMatthew G. Knepley 
472fb6725baSMatthew G. Knepley   Input Parameter:
473fb6725baSMatthew G. Knepley . da    - The DMDA
474fb6725baSMatthew G. Knepley 
475fb6725baSMatthew G. Knepley   Output Parameter:
476fb6725baSMatthew G. Knepley . stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
477fb6725baSMatthew G. Knepley 
478fb6725baSMatthew G. Knepley   Level: intermediate
479fb6725baSMatthew G. Knepley 
480fb6725baSMatthew G. Knepley .keywords:  distributed array, stencil
481fb6725baSMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA
482fb6725baSMatthew G. Knepley @*/
483fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype)
484fb6725baSMatthew G. Knepley {
485fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA*)da->data;
486fb6725baSMatthew G. Knepley 
487fb6725baSMatthew G. Knepley   PetscFunctionBegin;
488a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
489fb6725baSMatthew G. Knepley   PetscValidPointer(stype,2);
490fb6725baSMatthew G. Knepley   *stype = dd->stencil_type;
491fb6725baSMatthew G. Knepley   PetscFunctionReturn(0);
492fb6725baSMatthew G. Knepley }
493fb6725baSMatthew G. Knepley 
49447c6ae99SBarry Smith /*@
495aa219208SBarry Smith   DMDASetStencilWidth - Sets the width of the communication stencil
49647c6ae99SBarry Smith 
497aa219208SBarry Smith   Logically Collective on DMDA
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith   Input Parameter:
500aa219208SBarry Smith + da    - The DMDA
50147c6ae99SBarry Smith - width - The stencil width
50247c6ae99SBarry Smith 
50347c6ae99SBarry Smith   Level: intermediate
50447c6ae99SBarry Smith 
50547c6ae99SBarry Smith .keywords:  distributed array, stencil
50696e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
50747c6ae99SBarry Smith @*/
5087087cfbeSBarry Smith PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
50947c6ae99SBarry Smith {
51047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
51147c6ae99SBarry Smith 
51247c6ae99SBarry Smith   PetscFunctionBegin;
513a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
51447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,width,2);
515ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
51647c6ae99SBarry Smith   dd->s = width;
51747c6ae99SBarry Smith   PetscFunctionReturn(0);
51847c6ae99SBarry Smith }
51947c6ae99SBarry Smith 
520fb6725baSMatthew G. Knepley /*@
521fb6725baSMatthew G. Knepley   DMDAGetStencilWidth - Gets the width of the communication stencil
522fb6725baSMatthew G. Knepley 
523fb6725baSMatthew G. Knepley   Not collective
524fb6725baSMatthew G. Knepley 
525fb6725baSMatthew G. Knepley   Input Parameter:
526fb6725baSMatthew G. Knepley . da    - The DMDA
527fb6725baSMatthew G. Knepley 
528fb6725baSMatthew G. Knepley   Output Parameter:
529fb6725baSMatthew G. Knepley . width - The stencil width
530fb6725baSMatthew G. Knepley 
531fb6725baSMatthew G. Knepley   Level: intermediate
532fb6725baSMatthew G. Knepley 
533fb6725baSMatthew G. Knepley .keywords:  distributed array, stencil
534fb6725baSMatthew G. Knepley .seealso: DMDACreate(), DMDestroy(), DMDA
535fb6725baSMatthew G. Knepley @*/
536fb6725baSMatthew G. Knepley PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width)
537fb6725baSMatthew G. Knepley {
538fb6725baSMatthew G. Knepley   DM_DA *dd = (DM_DA *) da->data;
539fb6725baSMatthew G. Knepley 
540fb6725baSMatthew G. Knepley   PetscFunctionBegin;
541a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
542fb6725baSMatthew G. Knepley   PetscValidPointer(width,2);
543fb6725baSMatthew G. Knepley   *width = dd->s;
544fb6725baSMatthew G. Knepley   PetscFunctionReturn(0);
545fb6725baSMatthew G. Knepley }
546fb6725baSMatthew G. Knepley 
547aa219208SBarry Smith static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
54847c6ae99SBarry Smith {
54947c6ae99SBarry Smith   PetscInt i,sum;
55047c6ae99SBarry Smith 
55147c6ae99SBarry Smith   PetscFunctionBegin;
552ce94432eSBarry Smith   if (M < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
55347c6ae99SBarry Smith   for (i=sum=0; i<m; i++) sum += lx[i];
554ce94432eSBarry Smith   if (sum != M) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
55547c6ae99SBarry Smith   PetscFunctionReturn(0);
55647c6ae99SBarry Smith }
55747c6ae99SBarry Smith 
55847c6ae99SBarry Smith /*@
559aa219208SBarry Smith   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
56047c6ae99SBarry Smith 
561aa219208SBarry Smith   Logically Collective on DMDA
56247c6ae99SBarry Smith 
56347c6ae99SBarry Smith   Input Parameter:
564aa219208SBarry Smith + da - The DMDA
5650298fd71SBarry 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
5660298fd71SBarry 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
5670298fd71SBarry 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.
56847c6ae99SBarry Smith 
56947c6ae99SBarry Smith   Level: intermediate
57047c6ae99SBarry Smith 
571e3f3e4b6SBarry Smith   Note: these numbers are NOT multiplied by the number of dof per node.
572e3f3e4b6SBarry Smith 
57347c6ae99SBarry Smith .keywords:  distributed array
57496e147daSBarry Smith .seealso: DMDACreate(), DMDestroy(), DMDA
57547c6ae99SBarry Smith @*/
5767087cfbeSBarry Smith PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
57747c6ae99SBarry Smith {
57847c6ae99SBarry Smith   PetscErrorCode ierr;
57947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
58047c6ae99SBarry Smith 
58147c6ae99SBarry Smith   PetscFunctionBegin;
582a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
583ce94432eSBarry Smith   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
58447c6ae99SBarry Smith   if (lx) {
585ce94432eSBarry Smith     if (dd->m < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
586aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
58747c6ae99SBarry Smith     if (!dd->lx) {
588785e854fSJed Brown       ierr = PetscMalloc1(dd->m, &dd->lx);CHKERRQ(ierr);
58947c6ae99SBarry Smith     }
59047c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
59147c6ae99SBarry Smith   }
59247c6ae99SBarry Smith   if (ly) {
593ce94432eSBarry Smith     if (dd->n < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
594aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
59547c6ae99SBarry Smith     if (!dd->ly) {
596785e854fSJed Brown       ierr = PetscMalloc1(dd->n, &dd->ly);CHKERRQ(ierr);
59747c6ae99SBarry Smith     }
59847c6ae99SBarry Smith     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
59947c6ae99SBarry Smith   }
60047c6ae99SBarry Smith   if (lz) {
601ce94432eSBarry Smith     if (dd->p < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
602aa219208SBarry Smith     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
60347c6ae99SBarry Smith     if (!dd->lz) {
604785e854fSJed Brown       ierr = PetscMalloc1(dd->p, &dd->lz);CHKERRQ(ierr);
60547c6ae99SBarry Smith     }
60647c6ae99SBarry Smith     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
60747c6ae99SBarry Smith   }
60847c6ae99SBarry Smith   PetscFunctionReturn(0);
60947c6ae99SBarry Smith }
61047c6ae99SBarry Smith 
61147c6ae99SBarry Smith /*@
612aa219208SBarry Smith        DMDASetInterpolationType - Sets the type of interpolation that will be
613e727c939SJed Brown           returned by DMCreateInterpolation()
61447c6ae99SBarry Smith 
615aa219208SBarry Smith    Logically Collective on DMDA
61647c6ae99SBarry Smith 
61747c6ae99SBarry Smith    Input Parameter:
61847c6ae99SBarry Smith +  da - initial distributed array
619aa219208SBarry Smith .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
62047c6ae99SBarry Smith 
62147c6ae99SBarry Smith    Level: intermediate
62247c6ae99SBarry Smith 
62395452b02SPatrick Sanan    Notes:
62495452b02SPatrick Sanan     you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation()
62547c6ae99SBarry Smith 
62647c6ae99SBarry Smith .keywords:  distributed array, interpolation
62747c6ae99SBarry Smith 
62896e147daSBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
62947c6ae99SBarry Smith @*/
6307087cfbeSBarry Smith PetscErrorCode  DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
63147c6ae99SBarry Smith {
63247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
63347c6ae99SBarry Smith 
63447c6ae99SBarry Smith   PetscFunctionBegin;
635a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
63647c6ae99SBarry Smith   PetscValidLogicalCollectiveEnum(da,ctype,2);
63747c6ae99SBarry Smith   dd->interptype = ctype;
63847c6ae99SBarry Smith   PetscFunctionReturn(0);
63947c6ae99SBarry Smith }
64047c6ae99SBarry Smith 
6412dde6fd4SLisandro Dalcin /*@
6422dde6fd4SLisandro Dalcin        DMDAGetInterpolationType - Gets the type of interpolation that will be
643e727c939SJed Brown           used by DMCreateInterpolation()
6442dde6fd4SLisandro Dalcin 
6452dde6fd4SLisandro Dalcin    Not Collective
6462dde6fd4SLisandro Dalcin 
6472dde6fd4SLisandro Dalcin    Input Parameter:
6482dde6fd4SLisandro Dalcin .  da - distributed array
6492dde6fd4SLisandro Dalcin 
6502dde6fd4SLisandro Dalcin    Output Parameter:
6512dde6fd4SLisandro Dalcin .  ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms)
6522dde6fd4SLisandro Dalcin 
6532dde6fd4SLisandro Dalcin    Level: intermediate
6542dde6fd4SLisandro Dalcin 
6552dde6fd4SLisandro Dalcin .keywords:  distributed array, interpolation
6562dde6fd4SLisandro Dalcin 
657e727c939SJed Brown .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation()
6582dde6fd4SLisandro Dalcin @*/
6592dde6fd4SLisandro Dalcin PetscErrorCode  DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype)
6602dde6fd4SLisandro Dalcin {
6612dde6fd4SLisandro Dalcin   DM_DA *dd = (DM_DA*)da->data;
6622dde6fd4SLisandro Dalcin 
6632dde6fd4SLisandro Dalcin   PetscFunctionBegin;
664a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
6652dde6fd4SLisandro Dalcin   PetscValidPointer(ctype,2);
6662dde6fd4SLisandro Dalcin   *ctype = dd->interptype;
6672dde6fd4SLisandro Dalcin   PetscFunctionReturn(0);
6682dde6fd4SLisandro Dalcin }
66947c6ae99SBarry Smith 
6706a119db4SBarry Smith /*@C
671aa219208SBarry Smith       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
67247c6ae99SBarry Smith         processes neighbors.
67347c6ae99SBarry Smith 
67447c6ae99SBarry Smith     Not Collective
67547c6ae99SBarry Smith 
67647c6ae99SBarry Smith    Input Parameter:
677aa219208SBarry Smith .     da - the DMDA object
67847c6ae99SBarry Smith 
67947c6ae99SBarry Smith    Output Parameters:
68047c6ae99SBarry Smith .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
68147c6ae99SBarry Smith               this process itself is in the list
68247c6ae99SBarry Smith 
68395452b02SPatrick Sanan    Notes:
68495452b02SPatrick Sanan     In 2d the array is of length 9, in 3d of length 27
68547c6ae99SBarry Smith           Not supported in 1d
686aa219208SBarry Smith           Do not free the array, it is freed when the DMDA is destroyed.
68747c6ae99SBarry Smith 
68895452b02SPatrick Sanan    Fortran Notes:
68995452b02SPatrick Sanan     In fortran you must pass in an array of the appropriate length.
69047c6ae99SBarry Smith 
69147c6ae99SBarry Smith    Level: intermediate
69247c6ae99SBarry Smith 
69347c6ae99SBarry Smith @*/
6947087cfbeSBarry Smith PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
69547c6ae99SBarry Smith {
69647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
6975fd66863SKarl Rupp 
69847c6ae99SBarry Smith   PetscFunctionBegin;
699a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
70047c6ae99SBarry Smith   *ranks = dd->neighbors;
70147c6ae99SBarry Smith   PetscFunctionReturn(0);
70247c6ae99SBarry Smith }
70347c6ae99SBarry Smith 
70447c6ae99SBarry Smith /*@C
705aa219208SBarry Smith       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
70647c6ae99SBarry Smith 
70747c6ae99SBarry Smith     Not Collective
70847c6ae99SBarry Smith 
70947c6ae99SBarry Smith    Input Parameter:
710aa219208SBarry Smith .     da - the DMDA object
71147c6ae99SBarry Smith 
71247c6ae99SBarry Smith    Output Parameter:
71347c6ae99SBarry Smith +     lx - ownership along x direction (optional)
71447c6ae99SBarry Smith .     ly - ownership along y direction (optional)
71547c6ae99SBarry Smith -     lz - ownership along z direction (optional)
71647c6ae99SBarry Smith 
71747c6ae99SBarry Smith    Level: intermediate
71847c6ae99SBarry Smith 
719aa219208SBarry Smith     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
72047c6ae99SBarry Smith 
72147c6ae99SBarry Smith     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
722aa219208SBarry Smith     eighth arguments from DMDAGetInfo()
72347c6ae99SBarry Smith 
72447c6ae99SBarry Smith      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
725aa219208SBarry Smith     DMDA they came from still exists (has not been destroyed).
72647c6ae99SBarry Smith 
727e3f3e4b6SBarry Smith     These numbers are NOT multiplied by the number of dof per node.
728e3f3e4b6SBarry Smith 
7291fd49c25SBarry Smith      Not available from Fortran
7301fd49c25SBarry Smith 
731aa219208SBarry Smith .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
73247c6ae99SBarry Smith @*/
7337087cfbeSBarry Smith PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
73447c6ae99SBarry Smith {
73547c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
73647c6ae99SBarry Smith 
73747c6ae99SBarry Smith   PetscFunctionBegin;
738a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
73947c6ae99SBarry Smith   if (lx) *lx = dd->lx;
74047c6ae99SBarry Smith   if (ly) *ly = dd->ly;
74147c6ae99SBarry Smith   if (lz) *lz = dd->lz;
74247c6ae99SBarry Smith   PetscFunctionReturn(0);
74347c6ae99SBarry Smith }
74447c6ae99SBarry Smith 
74547c6ae99SBarry Smith /*@
746aa219208SBarry Smith      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
74747c6ae99SBarry Smith 
748aa219208SBarry Smith     Logically Collective on DMDA
74947c6ae99SBarry Smith 
75047c6ae99SBarry Smith   Input Parameters:
751aa219208SBarry Smith +    da - the DMDA object
75247c6ae99SBarry Smith .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
75347c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
75447c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
75547c6ae99SBarry Smith 
75647c6ae99SBarry Smith   Options Database:
75747c6ae99SBarry Smith +  -da_refine_x - refinement ratio in x direction
75847c6ae99SBarry Smith .  -da_refine_y - refinement ratio in y direction
75947c6ae99SBarry Smith -  -da_refine_z - refinement ratio in z direction
76047c6ae99SBarry Smith 
76147c6ae99SBarry Smith   Level: intermediate
76247c6ae99SBarry Smith 
76395452b02SPatrick Sanan     Notes:
76495452b02SPatrick Sanan     Pass PETSC_IGNORE to leave a value unchanged
76547c6ae99SBarry Smith 
766aa219208SBarry Smith .seealso: DMRefine(), DMDAGetRefinementFactor()
76747c6ae99SBarry Smith @*/
7687087cfbeSBarry Smith PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
76947c6ae99SBarry Smith {
77047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
77147c6ae99SBarry Smith 
77247c6ae99SBarry Smith   PetscFunctionBegin;
773a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
77447c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_x,2);
77547c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_y,3);
77647c6ae99SBarry Smith   PetscValidLogicalCollectiveInt(da,refine_z,4);
77747c6ae99SBarry Smith 
77847c6ae99SBarry Smith   if (refine_x > 0) dd->refine_x = refine_x;
77947c6ae99SBarry Smith   if (refine_y > 0) dd->refine_y = refine_y;
78047c6ae99SBarry Smith   if (refine_z > 0) dd->refine_z = refine_z;
78147c6ae99SBarry Smith   PetscFunctionReturn(0);
78247c6ae99SBarry Smith }
78347c6ae99SBarry Smith 
78447c6ae99SBarry Smith /*@C
785aa219208SBarry Smith      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
78647c6ae99SBarry Smith 
78747c6ae99SBarry Smith     Not Collective
78847c6ae99SBarry Smith 
78947c6ae99SBarry Smith   Input Parameter:
790aa219208SBarry Smith .    da - the DMDA object
79147c6ae99SBarry Smith 
79247c6ae99SBarry Smith   Output Parameters:
79347c6ae99SBarry Smith +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
79447c6ae99SBarry Smith .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
79547c6ae99SBarry Smith -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
79647c6ae99SBarry Smith 
79747c6ae99SBarry Smith   Level: intermediate
79847c6ae99SBarry Smith 
79995452b02SPatrick Sanan     Notes:
80095452b02SPatrick Sanan     Pass NULL for values you do not need
80147c6ae99SBarry Smith 
802aa219208SBarry Smith .seealso: DMRefine(), DMDASetRefinementFactor()
80347c6ae99SBarry Smith @*/
8047087cfbeSBarry Smith PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
80547c6ae99SBarry Smith {
80647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
80747c6ae99SBarry Smith 
80847c6ae99SBarry Smith   PetscFunctionBegin;
809a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
81047c6ae99SBarry Smith   if (refine_x) *refine_x = dd->refine_x;
81147c6ae99SBarry Smith   if (refine_y) *refine_y = dd->refine_y;
81247c6ae99SBarry Smith   if (refine_z) *refine_z = dd->refine_z;
81347c6ae99SBarry Smith   PetscFunctionReturn(0);
81447c6ae99SBarry Smith }
81547c6ae99SBarry Smith 
81647c6ae99SBarry Smith /*@C
817aa219208SBarry Smith      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
81847c6ae99SBarry Smith 
819aa219208SBarry Smith     Logically Collective on DMDA
82047c6ae99SBarry Smith 
82147c6ae99SBarry Smith   Input Parameters:
822aa219208SBarry Smith +    da - the DMDA object
823aa219208SBarry Smith -    f - the function that allocates the matrix for that specific DMDA
82447c6ae99SBarry Smith 
82547c6ae99SBarry Smith   Level: developer
82647c6ae99SBarry Smith 
82795452b02SPatrick Sanan    Notes:
82895452b02SPatrick Sanan     See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
82947c6ae99SBarry Smith        the diagonal and off-diagonal blocks of the matrix
83047c6ae99SBarry Smith 
8311fd49c25SBarry Smith    Not supported from Fortran
8321fd49c25SBarry Smith 
833950540a4SJed Brown .seealso: DMCreateMatrix(), DMDASetBlockFills()
83447c6ae99SBarry Smith @*/
835b412c318SBarry Smith PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, Mat*))
83647c6ae99SBarry Smith {
83747c6ae99SBarry Smith   PetscFunctionBegin;
838a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
83925296bd5SBarry Smith   da->ops->creatematrix = f;
84047c6ae99SBarry Smith   PetscFunctionReturn(0);
84147c6ae99SBarry Smith }
84247c6ae99SBarry Smith 
84347c6ae99SBarry Smith /*
84447c6ae99SBarry Smith   Creates "balanced" ownership ranges after refinement, constrained by the need for the
84547c6ae99SBarry Smith   fine grid boundaries to fall within one stencil width of the coarse partition.
84647c6ae99SBarry Smith 
84747c6ae99SBarry Smith   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
84847c6ae99SBarry Smith */
84994013140SBarry Smith static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
85047c6ae99SBarry Smith {
85147c6ae99SBarry Smith   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
85247c6ae99SBarry Smith   PetscErrorCode ierr;
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith   PetscFunctionBegin;
855ce94432eSBarry Smith   if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
85647c6ae99SBarry Smith   if (ratio == 1) {
85747c6ae99SBarry Smith     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
85847c6ae99SBarry Smith     PetscFunctionReturn(0);
85947c6ae99SBarry Smith   }
86047c6ae99SBarry Smith   for (i=0; i<m; i++) totalc += lc[i];
86147c6ae99SBarry Smith   remaining = (!periodic) + ratio * (totalc - (!periodic));
86247c6ae99SBarry Smith   for (i=0; i<m; i++) {
86347c6ae99SBarry Smith     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
86447c6ae99SBarry Smith     if (i == m-1) lf[i] = want;
86547c6ae99SBarry Smith     else {
8667aca7175SJed Brown       const PetscInt nextc = startc + lc[i];
8677aca7175SJed Brown       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
8687aca7175SJed Brown        * coarse stencil width of the first coarse node in the next subdomain. */
8697aca7175SJed Brown       while ((startf+want)/ratio < nextc - stencil_width) want++;
8707aca7175SJed Brown       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
8717aca7175SJed Brown        * coarse stencil width of the last coarse node in the current subdomain. */
8727aca7175SJed Brown       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
8737aca7175SJed Brown       /* Make sure all constraints are satisfied */
87430729d88SBarry Smith       if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width)
875ce94432eSBarry 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");
87647c6ae99SBarry Smith     }
87747c6ae99SBarry Smith     lf[i]      = want;
87847c6ae99SBarry Smith     startc    += lc[i];
87947c6ae99SBarry Smith     startf    += lf[i];
88047c6ae99SBarry Smith     remaining -= lf[i];
88147c6ae99SBarry Smith   }
88247c6ae99SBarry Smith   PetscFunctionReturn(0);
88347c6ae99SBarry Smith }
88447c6ae99SBarry Smith 
8852be375d4SJed Brown /*
8862be375d4SJed Brown   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
8872be375d4SJed Brown   fine grid boundaries to fall within one stencil width of the coarse partition.
8882be375d4SJed Brown 
8892be375d4SJed Brown   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
8902be375d4SJed Brown */
8912be375d4SJed Brown static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
8922be375d4SJed Brown {
8932be375d4SJed Brown   PetscInt       i,totalf,remaining,startc,startf;
8942be375d4SJed Brown   PetscErrorCode ierr;
8952be375d4SJed Brown 
8962be375d4SJed Brown   PetscFunctionBegin;
897ce94432eSBarry Smith   if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
8982be375d4SJed Brown   if (ratio == 1) {
8992be375d4SJed Brown     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
9002be375d4SJed Brown     PetscFunctionReturn(0);
9012be375d4SJed Brown   }
9022be375d4SJed Brown   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
9032be375d4SJed Brown   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
9042be375d4SJed Brown   for (i=0,startc=0,startf=0; i<m; i++) {
9052be375d4SJed Brown     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
9062be375d4SJed Brown     if (i == m-1) lc[i] = want;
9072be375d4SJed Brown     else {
9082be375d4SJed Brown       const PetscInt nextf = startf+lf[i];
9092be375d4SJed Brown       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
9102be375d4SJed Brown        * node is within one stencil width. */
9112be375d4SJed Brown       while (nextf/ratio < startc+want-stencil_width) want--;
9122be375d4SJed Brown       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
9132be375d4SJed Brown        * fine node is within one stencil width. */
9142be375d4SJed Brown       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
9152be375d4SJed Brown       if (want < 0 || want > remaining
916ce94432eSBarry 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");
9172be375d4SJed Brown     }
9182be375d4SJed Brown     lc[i]      = want;
9192be375d4SJed Brown     startc    += lc[i];
9202be375d4SJed Brown     startf    += lf[i];
9212be375d4SJed Brown     remaining -= lc[i];
9222be375d4SJed Brown   }
9232be375d4SJed Brown   PetscFunctionReturn(0);
9242be375d4SJed Brown }
9252be375d4SJed Brown 
9267087cfbeSBarry Smith PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
92747c6ae99SBarry Smith {
92847c6ae99SBarry Smith   PetscErrorCode ierr;
929c73cfb54SMatthew G. Knepley   PetscInt       M,N,P,i,dim;
9309a42bb27SBarry Smith   DM             da2;
93147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
93247c6ae99SBarry Smith 
93347c6ae99SBarry Smith   PetscFunctionBegin;
934a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
93547c6ae99SBarry Smith   PetscValidPointer(daref,3);
93647c6ae99SBarry Smith 
937c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(da, &dim);CHKERRQ(ierr);
938bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
93947c6ae99SBarry Smith     M = dd->refine_x*dd->M;
94047c6ae99SBarry Smith   } else {
94147c6ae99SBarry Smith     M = 1 + dd->refine_x*(dd->M - 1);
94247c6ae99SBarry Smith   }
943bff4a2f0SMatthew G. Knepley   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
944c73cfb54SMatthew G. Knepley     if (dim > 1) {
94547c6ae99SBarry Smith       N = dd->refine_y*dd->N;
94647c6ae99SBarry Smith     } else {
9471860e6e9SBarry Smith       N = 1;
9481860e6e9SBarry Smith     }
9491860e6e9SBarry Smith   } else {
95047c6ae99SBarry Smith     N = 1 + dd->refine_y*(dd->N - 1);
95147c6ae99SBarry Smith   }
952bff4a2f0SMatthew G. Knepley   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
953c73cfb54SMatthew G. Knepley     if (dim > 2) {
95447c6ae99SBarry Smith       P = dd->refine_z*dd->P;
95547c6ae99SBarry Smith     } else {
9561860e6e9SBarry Smith       P = 1;
9571860e6e9SBarry Smith     }
9581860e6e9SBarry Smith   } else {
95947c6ae99SBarry Smith     P = 1 + dd->refine_z*(dd->P - 1);
96047c6ae99SBarry Smith   }
961ce94432eSBarry Smith   ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr);
962397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
963c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(da2,dim);CHKERRQ(ierr);
964397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
965397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
966397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
967397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
968397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
969397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
970c73cfb54SMatthew G. Knepley   if (dim == 3) {
97147c6ae99SBarry Smith     PetscInt *lx,*ly,*lz;
972dcca6d9dSJed Brown     ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr);
973bff4a2f0SMatthew G. Knepley     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
974bff4a2f0SMatthew G. Knepley     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
975bff4a2f0SMatthew G. Knepley     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
976397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
97747c6ae99SBarry Smith     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
978c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
97947c6ae99SBarry Smith     PetscInt *lx,*ly;
980dcca6d9dSJed Brown     ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr);
981bff4a2f0SMatthew G. Knepley     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
982bff4a2f0SMatthew G. Knepley     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
9830298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr);
98447c6ae99SBarry Smith     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
985c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
98647c6ae99SBarry Smith     PetscInt *lx;
987785e854fSJed Brown     ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr);
988bff4a2f0SMatthew G. Knepley     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
9890298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr);
99047c6ae99SBarry Smith     ierr = PetscFree(lx);CHKERRQ(ierr);
99147c6ae99SBarry Smith   }
99247c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
99347c6ae99SBarry Smith 
99447c6ae99SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
99525296bd5SBarry Smith   da2->ops->creatematrix = da->ops->creatematrix;
99625296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
99747c6ae99SBarry Smith   da2->ops->getcoloring = da->ops->getcoloring;
99847c6ae99SBarry Smith   dd2->interptype       = dd->interptype;
99947c6ae99SBarry Smith 
100047c6ae99SBarry Smith   /* copy fill information if given */
100147c6ae99SBarry Smith   if (dd->dfill) {
1002854ce69bSBarry Smith     ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr);
100347c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
100447c6ae99SBarry Smith   }
100547c6ae99SBarry Smith   if (dd->ofill) {
1006854ce69bSBarry Smith     ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr);
100747c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
100847c6ae99SBarry Smith   }
100947c6ae99SBarry Smith   /* copy the refine information */
1010397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
1011397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
1012397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
101347c6ae99SBarry Smith 
1014897f7067SBarry Smith   if (dd->refine_z_hier) {
1015897f7067SBarry Smith     if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_z_hier_n) {
1016897f7067SBarry Smith       dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown + 1];
1017897f7067SBarry Smith     }
1018897f7067SBarry Smith     if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_z_hier_n) {
1019897f7067SBarry Smith       dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown];
1020897f7067SBarry Smith     }
1021897f7067SBarry Smith     dd2->refine_z_hier_n = dd->refine_z_hier_n;
1022897f7067SBarry Smith     ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr);
1023897f7067SBarry Smith     ierr = PetscMemcpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n*sizeof(PetscInt));CHKERRQ(ierr);
1024897f7067SBarry Smith   }
1025897f7067SBarry Smith   if (dd->refine_y_hier) {
1026897f7067SBarry Smith     if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_y_hier_n) {
1027897f7067SBarry Smith       dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown + 1];
1028897f7067SBarry Smith     }
1029897f7067SBarry Smith     if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_y_hier_n) {
1030897f7067SBarry Smith       dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown];
1031897f7067SBarry Smith     }
1032897f7067SBarry Smith     dd2->refine_y_hier_n = dd->refine_y_hier_n;
1033897f7067SBarry Smith     ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr);
1034897f7067SBarry Smith     ierr = PetscMemcpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n*sizeof(PetscInt));CHKERRQ(ierr);
1035897f7067SBarry Smith   }
1036897f7067SBarry Smith   if (dd->refine_x_hier) {
1037897f7067SBarry Smith     if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_x_hier_n) {
1038897f7067SBarry Smith       dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown + 1];
1039897f7067SBarry Smith     }
1040897f7067SBarry Smith     if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_x_hier_n) {
1041897f7067SBarry Smith       dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown];
1042897f7067SBarry Smith     }
1043897f7067SBarry Smith     dd2->refine_x_hier_n = dd->refine_x_hier_n;
1044897f7067SBarry Smith     ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr);
1045897f7067SBarry Smith     ierr = PetscMemcpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n*sizeof(PetscInt));CHKERRQ(ierr);
1046897f7067SBarry Smith   }
1047897f7067SBarry Smith 
1048897f7067SBarry Smith 
104947c6ae99SBarry Smith   /* copy vector type information */
1050c0dedaeaSBarry Smith   ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr);
1051ddcf8b74SDave May 
1052efd51863SBarry Smith   dd2->lf = dd->lf;
1053efd51863SBarry Smith   dd2->lj = dd->lj;
1054efd51863SBarry Smith 
10556e87535bSJed Brown   da2->leveldown = da->leveldown;
10566e87535bSJed Brown   da2->levelup   = da->levelup + 1;
10578865f1eaSKarl Rupp 
10586e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
10596e87535bSJed Brown 
1060ddcf8b74SDave May   /* interpolate coordinates if they are set on the coarse grid */
10616636e97aSMatthew G Knepley   if (da->coordinates) {
1062ddcf8b74SDave May     DM  cdaf,cdac;
1063ddcf8b74SDave May     Vec coordsc,coordsf;
1064ddcf8b74SDave May     Mat II;
1065ddcf8b74SDave May 
10666636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr);
10676636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr);
10686636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr);
1069b61d3410SDave May     /* force creation of the coordinate vector */
1070b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
10716636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
10720298fd71SBarry Smith     ierr = DMCreateInterpolation(cdac,cdaf,&II,NULL);CHKERRQ(ierr);
1073ddcf8b74SDave May     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
1074fcfd50ebSBarry Smith     ierr = MatDestroy(&II);CHKERRQ(ierr);
1075ddcf8b74SDave May   }
1076397b6216SJed Brown 
1077f3141302SJed Brown   for (i=0; i<da->bs; i++) {
1078f3141302SJed Brown     const char *fieldname;
1079f3141302SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
1080f3141302SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
1081f3141302SJed Brown   }
1082397b6216SJed Brown 
108347c6ae99SBarry Smith   *daref = da2;
108447c6ae99SBarry Smith   PetscFunctionReturn(0);
108547c6ae99SBarry Smith }
108647c6ae99SBarry Smith 
1087397b6216SJed Brown 
10887087cfbeSBarry Smith PetscErrorCode  DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
108947c6ae99SBarry Smith {
109047c6ae99SBarry Smith   PetscErrorCode ierr;
1091c73cfb54SMatthew G. Knepley   PetscInt       M,N,P,i,dim;
10929a42bb27SBarry Smith   DM             da2;
109347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data,*dd2;
109447c6ae99SBarry Smith 
109547c6ae99SBarry Smith   PetscFunctionBegin;
1096a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
109747c6ae99SBarry Smith   PetscValidPointer(daref,3);
109847c6ae99SBarry Smith 
1099c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(da, &dim);CHKERRQ(ierr);
1100bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1101397b6216SJed Brown     M = dd->M / dd->coarsen_x;
110247c6ae99SBarry Smith   } else {
1103397b6216SJed Brown     M = 1 + (dd->M - 1) / dd->coarsen_x;
110447c6ae99SBarry Smith   }
1105bff4a2f0SMatthew G. Knepley   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1106c73cfb54SMatthew G. Knepley     if (dim > 1) {
1107397b6216SJed Brown       N = dd->N / dd->coarsen_y;
110847c6ae99SBarry Smith     } else {
11091860e6e9SBarry Smith       N = 1;
11101860e6e9SBarry Smith     }
11111860e6e9SBarry Smith   } else {
1112397b6216SJed Brown     N = 1 + (dd->N - 1) / dd->coarsen_y;
111347c6ae99SBarry Smith   }
1114bff4a2f0SMatthew G. Knepley   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1115c73cfb54SMatthew G. Knepley     if (dim > 2) {
1116397b6216SJed Brown       P = dd->P / dd->coarsen_z;
111747c6ae99SBarry Smith     } else {
11181860e6e9SBarry Smith       P = 1;
11191860e6e9SBarry Smith     }
11201860e6e9SBarry Smith   } else {
1121397b6216SJed Brown     P = 1 + (dd->P - 1) / dd->coarsen_z;
112247c6ae99SBarry Smith   }
1123ce94432eSBarry Smith   ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr);
1124397b6216SJed Brown   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
1125c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(da2,dim);CHKERRQ(ierr);
1126397b6216SJed Brown   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
1127397b6216SJed Brown   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
1128397b6216SJed Brown   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
1129397b6216SJed Brown   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
1130397b6216SJed Brown   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
1131397b6216SJed Brown   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
1132c73cfb54SMatthew G. Knepley   if (dim == 3) {
11332be375d4SJed Brown     PetscInt *lx,*ly,*lz;
1134dcca6d9dSJed Brown     ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr);
1135bff4a2f0SMatthew G. Knepley     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
1136bff4a2f0SMatthew G. Knepley     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
1137bff4a2f0SMatthew G. Knepley     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
1138397b6216SJed Brown     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
11392be375d4SJed Brown     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
1140c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
11412be375d4SJed Brown     PetscInt *lx,*ly;
1142dcca6d9dSJed Brown     ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr);
1143bff4a2f0SMatthew G. Knepley     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
1144bff4a2f0SMatthew G. Knepley     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
11450298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr);
11462be375d4SJed Brown     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
1147c73cfb54SMatthew G. Knepley   } else if (dim == 1) {
11482be375d4SJed Brown     PetscInt *lx;
1149785e854fSJed Brown     ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr);
1150bff4a2f0SMatthew G. Knepley     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
11510298fd71SBarry Smith     ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr);
11522be375d4SJed Brown     ierr = PetscFree(lx);CHKERRQ(ierr);
115347c6ae99SBarry Smith   }
115447c6ae99SBarry Smith   dd2 = (DM_DA*)da2->data;
115547c6ae99SBarry Smith 
11564dcab191SBarry Smith   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
115725296bd5SBarry Smith   /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */
115825296bd5SBarry Smith   da2->ops->creatematrix = da->ops->creatematrix;
115947c6ae99SBarry Smith   da2->ops->getcoloring  = da->ops->getcoloring;
116047c6ae99SBarry Smith   dd2->interptype        = dd->interptype;
116147c6ae99SBarry Smith 
116247c6ae99SBarry Smith   /* copy fill information if given */
116347c6ae99SBarry Smith   if (dd->dfill) {
1164854ce69bSBarry Smith     ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr);
116547c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
116647c6ae99SBarry Smith   }
116747c6ae99SBarry Smith   if (dd->ofill) {
1168854ce69bSBarry Smith     ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr);
116947c6ae99SBarry Smith     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
117047c6ae99SBarry Smith   }
117147c6ae99SBarry Smith   /* copy the refine information */
1172397b6216SJed Brown   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
1173397b6216SJed Brown   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
1174397b6216SJed Brown   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
117547c6ae99SBarry Smith 
1176897f7067SBarry Smith   if (dd->refine_z_hier) {
1177897f7067SBarry Smith     if (da->levelup - da->leveldown -1 > -1 && da->levelup - da->leveldown - 1< dd->refine_z_hier_n) {
1178897f7067SBarry Smith       dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown - 1];
1179897f7067SBarry Smith     }
1180897f7067SBarry Smith     if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_z_hier_n) {
1181897f7067SBarry Smith       dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown - 2];
1182897f7067SBarry Smith     }
1183897f7067SBarry Smith     dd2->refine_z_hier_n = dd->refine_z_hier_n;
1184897f7067SBarry Smith     ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr);
1185897f7067SBarry Smith     ierr = PetscMemcpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n*sizeof(PetscInt));CHKERRQ(ierr);
1186897f7067SBarry Smith   }
1187897f7067SBarry Smith   if (dd->refine_y_hier) {
1188897f7067SBarry Smith     if (da->levelup - da->leveldown - 1 > -1 && da->levelup - da->leveldown - 1< dd->refine_y_hier_n) {
1189897f7067SBarry Smith       dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown - 1];
1190897f7067SBarry Smith     }
1191897f7067SBarry Smith     if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_y_hier_n) {
1192897f7067SBarry Smith       dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown - 2];
1193897f7067SBarry Smith     }
1194897f7067SBarry Smith     dd2->refine_y_hier_n = dd->refine_y_hier_n;
1195897f7067SBarry Smith     ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr);
1196897f7067SBarry Smith     ierr = PetscMemcpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n*sizeof(PetscInt));CHKERRQ(ierr);
1197897f7067SBarry Smith   }
1198897f7067SBarry Smith   if (dd->refine_x_hier) {
1199897f7067SBarry Smith     if (da->levelup - da->leveldown - 1 > -1 && da->levelup - da->leveldown - 1 < dd->refine_x_hier_n) {
1200897f7067SBarry Smith       dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown - 1];
1201897f7067SBarry Smith     }
1202897f7067SBarry Smith     if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_x_hier_n) {
1203897f7067SBarry Smith       dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown - 2];
1204897f7067SBarry Smith     }
1205897f7067SBarry Smith     dd2->refine_x_hier_n = dd->refine_x_hier_n;
1206897f7067SBarry Smith     ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr);
1207897f7067SBarry Smith     ierr = PetscMemcpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n*sizeof(PetscInt));CHKERRQ(ierr);
1208897f7067SBarry Smith   }
1209897f7067SBarry Smith 
121047c6ae99SBarry Smith   /* copy vector type information */
1211c0dedaeaSBarry Smith   ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr);
121247c6ae99SBarry Smith 
1213644e2e5bSBarry Smith   dd2->lf = dd->lf;
1214644e2e5bSBarry Smith   dd2->lj = dd->lj;
1215644e2e5bSBarry Smith 
12166e87535bSJed Brown   da2->leveldown = da->leveldown + 1;
12176e87535bSJed Brown   da2->levelup   = da->levelup;
12188865f1eaSKarl Rupp 
12196e87535bSJed Brown   ierr = DMSetUp(da2);CHKERRQ(ierr);
12206e87535bSJed Brown 
1221ddcf8b74SDave May   /* inject coordinates if they are set on the fine grid */
12226636e97aSMatthew G Knepley   if (da->coordinates) {
1223ddcf8b74SDave May     DM         cdaf,cdac;
1224ddcf8b74SDave May     Vec        coordsc,coordsf;
12256dbf9973SLawrence Mitchell     Mat        inject;
12266dbf9973SLawrence Mitchell     VecScatter vscat;
1227ddcf8b74SDave May 
12286636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr);
12296636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr);
12306636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr);
1231b61d3410SDave May     /* force creation of the coordinate vector */
1232b61d3410SDave May     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
12336636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr);
1234ddcf8b74SDave May 
1235e727c939SJed Brown     ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
12366dbf9973SLawrence Mitchell     ierr = MatScatterGetVecScatter(inject,&vscat);CHKERRQ(ierr);
12376dbf9973SLawrence Mitchell     ierr = VecScatterBegin(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12386dbf9973SLawrence Mitchell     ierr = VecScatterEnd(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12396dbf9973SLawrence Mitchell     ierr = MatDestroy(&inject);CHKERRQ(ierr);
1240ddcf8b74SDave May   }
1241f98405f7SJed Brown 
1242f98405f7SJed Brown   for (i=0; i<da->bs; i++) {
1243f98405f7SJed Brown     const char *fieldname;
1244f98405f7SJed Brown     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
1245f98405f7SJed Brown     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
1246f98405f7SJed Brown   }
1247f98405f7SJed Brown 
124847c6ae99SBarry Smith   *daref = da2;
124947c6ae99SBarry Smith   PetscFunctionReturn(0);
125047c6ae99SBarry Smith }
125147c6ae99SBarry Smith 
12527087cfbeSBarry Smith PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
125347c6ae99SBarry Smith {
125447c6ae99SBarry Smith   PetscErrorCode ierr;
125547c6ae99SBarry Smith   PetscInt       i,n,*refx,*refy,*refz;
125647c6ae99SBarry Smith 
125747c6ae99SBarry Smith   PetscFunctionBegin;
125847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1259ce94432eSBarry Smith   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
126047c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
126147c6ae99SBarry Smith   PetscValidPointer(daf,3);
126247c6ae99SBarry Smith 
1263aa219208SBarry Smith   /* Get refinement factors, defaults taken from the coarse DMDA */
1264dcca6d9dSJed Brown   ierr = PetscMalloc3(nlevels,&refx,nlevels,&refy,nlevels,&refz);CHKERRQ(ierr);
126547c6ae99SBarry Smith   for (i=0; i<nlevels; i++) {
1266aa219208SBarry Smith     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
126747c6ae99SBarry Smith   }
126847c6ae99SBarry Smith   n    = nlevels;
1269c5929fdfSBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr);
127047c6ae99SBarry Smith   n    = nlevels;
1271c5929fdfSBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr);
127247c6ae99SBarry Smith   n    = nlevels;
1273c5929fdfSBarry Smith   ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr);
127447c6ae99SBarry Smith 
1275aa219208SBarry Smith   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
1276ce94432eSBarry Smith   ierr = DMRefine(da,PetscObjectComm((PetscObject)da),&daf[0]);CHKERRQ(ierr);
127747c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
1278aa219208SBarry Smith     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
1279ce94432eSBarry Smith     ierr = DMRefine(daf[i-1],PetscObjectComm((PetscObject)da),&daf[i]);CHKERRQ(ierr);
128047c6ae99SBarry Smith   }
128147c6ae99SBarry Smith   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
128247c6ae99SBarry Smith   PetscFunctionReturn(0);
128347c6ae99SBarry Smith }
128447c6ae99SBarry Smith 
12857087cfbeSBarry Smith PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
128647c6ae99SBarry Smith {
128747c6ae99SBarry Smith   PetscErrorCode ierr;
128847c6ae99SBarry Smith   PetscInt       i;
128947c6ae99SBarry Smith 
129047c6ae99SBarry Smith   PetscFunctionBegin;
129147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1292ce94432eSBarry Smith   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
129347c6ae99SBarry Smith   if (nlevels == 0) PetscFunctionReturn(0);
129447c6ae99SBarry Smith   PetscValidPointer(dac,3);
1295ce94432eSBarry Smith   ierr = DMCoarsen(da,PetscObjectComm((PetscObject)da),&dac[0]);CHKERRQ(ierr);
129647c6ae99SBarry Smith   for (i=1; i<nlevels; i++) {
1297ce94432eSBarry Smith     ierr = DMCoarsen(dac[i-1],PetscObjectComm((PetscObject)da),&dac[i]);CHKERRQ(ierr);
129847c6ae99SBarry Smith   }
129947c6ae99SBarry Smith   PetscFunctionReturn(0);
130047c6ae99SBarry Smith }
130162197512SBarry Smith 
130262197512SBarry Smith #include <petscgll.h>
130362197512SBarry Smith 
130462197512SBarry Smith PetscErrorCode DMDASetGLLCoordinates_1d(DM dm,PetscGLL *gll)
130562197512SBarry Smith {
130662197512SBarry Smith   PetscErrorCode ierr;
130762197512SBarry Smith   PetscInt       i,j,n = gll->n,xs,xn,q;
130862197512SBarry Smith   PetscScalar    *xx;
130962197512SBarry Smith   PetscReal      h;
131062197512SBarry Smith   Vec            x;
131162197512SBarry Smith   DM_DA          *da = (DM_DA*)dm->data;
131262197512SBarry Smith 
131362197512SBarry Smith   PetscFunctionBegin;
131462197512SBarry Smith   if (da->bx != DM_BOUNDARY_PERIODIC) {
131562197512SBarry Smith     ierr = DMDAGetInfo(dm,NULL,&q,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
131662197512SBarry Smith     q    = (q-1)/(n-1);  /* number of spectral elements */
131762197512SBarry Smith     h    = 2.0/q;
131862197512SBarry Smith     ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr);
131962197512SBarry Smith     xs   = xs/(n-1);
132062197512SBarry Smith     xn   = xn/(n-1);
132162197512SBarry Smith     ierr = DMDASetUniformCoordinates(dm,-1.,1.,0.,0.,0.,0.);CHKERRQ(ierr);
132262197512SBarry Smith     ierr = DMGetCoordinates(dm,&x);CHKERRQ(ierr);
132362197512SBarry Smith     ierr = DMDAVecGetArray(dm,x,&xx);CHKERRQ(ierr);
132462197512SBarry Smith 
132562197512SBarry Smith     /* loop over local spectral elements */
132662197512SBarry Smith     for (j=xs; j<xs+xn; j++) {
132762197512SBarry Smith       /*
132862197512SBarry Smith        Except for the first process, each process starts on the second GLL point of the first element on that process
132962197512SBarry Smith        */
133062197512SBarry Smith       for (i= (j == xs && xs > 0)? 1 : 0; i<n; i++) {
133162197512SBarry Smith         xx[j*(n-1) + i] = -1.0 + h*j + h*(gll->nodes[i]+1.0)/2.;
133262197512SBarry Smith       }
133362197512SBarry Smith     }
133462197512SBarry Smith     ierr = DMDAVecRestoreArray(dm,x,&xx);CHKERRQ(ierr);
133562197512SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not yet implemented for periodic");
133662197512SBarry Smith   PetscFunctionReturn(0);
133762197512SBarry Smith }
133862197512SBarry Smith 
13391fd49c25SBarry Smith /*@
134062197512SBarry Smith 
134162197512SBarry Smith      DMDASetGLLCoordinates - Sets the global coordinates from -1 to 1 to the GLL points of as many GLL elements that fit the number of grid points
134262197512SBarry Smith 
134362197512SBarry Smith    Collective on DM
134462197512SBarry Smith 
134562197512SBarry Smith    Input Parameters:
134662197512SBarry Smith +   da - the DMDA object
134762197512SBarry Smith -   gll - the GLL object
134862197512SBarry Smith 
134995452b02SPatrick Sanan    Notes:
135095452b02SPatrick Sanan     the parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points
135162197512SBarry Smith           on each process much be divisible by the number of GLL elements needed per process. This depends on whether the DM is
135262197512SBarry Smith           periodic or not.
135362197512SBarry Smith 
1354edc382c3SSatish Balay    Level: advanced
1355edc382c3SSatish Balay 
135662197512SBarry Smith .seealso:   DMDACreate(), PetscGLLCreate(), DMGetCoordinates()
135762197512SBarry Smith @*/
135862197512SBarry Smith PetscErrorCode DMDASetGLLCoordinates(DM da,PetscGLL *gll)
135962197512SBarry Smith {
136062197512SBarry Smith   PetscErrorCode ierr;
136162197512SBarry Smith 
136262197512SBarry Smith   PetscFunctionBegin;
136362197512SBarry Smith   if (da->dim == 1) {
136462197512SBarry Smith     ierr = DMDASetGLLCoordinates_1d(da,gll);CHKERRQ(ierr);
136562197512SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not yet implemented for 2 or 3d");
136662197512SBarry Smith   PetscFunctionReturn(0);
136762197512SBarry Smith }
13687c3cd84eSPatrick Sanan 
13697c3cd84eSPatrick Sanan PETSC_INTERN PetscErrorCode DMGetCompatibility_DA(DM da1,DM dm2,PetscBool *compatible,PetscBool *set)
13707c3cd84eSPatrick Sanan {
13717c3cd84eSPatrick Sanan   PetscErrorCode ierr;
13727c3cd84eSPatrick Sanan   DM_DA          *dd1 = (DM_DA*)da1->data,*dd2;
13737c3cd84eSPatrick Sanan   DM             da2;
13747c3cd84eSPatrick Sanan   DMType         dmtype2;
13757c3cd84eSPatrick Sanan   PetscBool      isda,compatibleLocal;
13767c3cd84eSPatrick Sanan   PetscInt       i;
13777c3cd84eSPatrick Sanan 
13787c3cd84eSPatrick Sanan   PetscFunctionBegin;
13797c3cd84eSPatrick Sanan   if (!da1->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da1),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() must be called on first DM before DMGetCompatibility()");
13807c3cd84eSPatrick Sanan   ierr = DMGetType(dm2,&dmtype2);CHKERRQ(ierr);
13817c3cd84eSPatrick Sanan   ierr = PetscStrcmp(dmtype2,DMDA,&isda);CHKERRQ(ierr);
13827c3cd84eSPatrick Sanan   if (isda) {
13837c3cd84eSPatrick Sanan     da2 = dm2;
13847c3cd84eSPatrick Sanan     dd2 = (DM_DA*)da2->data;
13857c3cd84eSPatrick Sanan     if (!da2->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da2),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() must be called on second DM before DMGetCompatibility()");
1386*110623a0SKarl Rupp     compatibleLocal = (PetscBool)(da1->dim == da2->dim);
1387*110623a0SKarl Rupp     if (compatibleLocal) compatibleLocal &= (PetscBool)(dd1->s == dd2->s); /* Stencil width */
1388*110623a0SKarl Rupp     /*                                                                      Global size              ranks               Boundary type */
1389*110623a0SKarl Rupp     if (compatibleLocal)                 compatibleLocal &= (PetscBool)((dd1->M  == dd2->M) && (dd1->m == dd2->m) && (dd1->bx == dd2->bx));
1390*110623a0SKarl Rupp     if (compatibleLocal && da1->dim > 1) compatibleLocal &= (PetscBool)((dd1->N  == dd2->N) && (dd1->n == dd2->n) && (dd1->by == dd2->by));
1391*110623a0SKarl Rupp     if (compatibleLocal && da1->dim > 2) compatibleLocal &= (PetscBool)((dd1->P  == dd2->P) && (dd1->p == dd2->p) && (dd1->bz == dd2->bz));
13927c3cd84eSPatrick Sanan     if (compatibleLocal) {
13937c3cd84eSPatrick Sanan       for (i=0; i<dd1->m; ++i) {
1394*110623a0SKarl Rupp         compatibleLocal &= (PetscBool)(dd1->lx[i] == dd2->lx[i]);           /* Local size     */
13957c3cd84eSPatrick Sanan       }
13967c3cd84eSPatrick Sanan     }
13977c3cd84eSPatrick Sanan     if (compatibleLocal && da1->dim > 1) {
13987c3cd84eSPatrick Sanan       for (i=0; i<dd1->n; ++i) {
1399*110623a0SKarl Rupp         compatibleLocal &= (PetscBool)(dd1->ly[i] == dd2->ly[i]);
14007c3cd84eSPatrick Sanan       }
14017c3cd84eSPatrick Sanan     }
14027c3cd84eSPatrick Sanan     if (compatibleLocal && da1->dim > 2) {
14037c3cd84eSPatrick Sanan       for (i=0; i<dd1->p; ++i) {
1404*110623a0SKarl Rupp         compatibleLocal &= (PetscBool)(dd1->lz[i] == dd2->lz[i]);
14057c3cd84eSPatrick Sanan       }
14067c3cd84eSPatrick Sanan     }
14077c3cd84eSPatrick Sanan     *compatible = compatibleLocal;
14087c3cd84eSPatrick Sanan     *set = PETSC_TRUE;
14097c3cd84eSPatrick Sanan   } else {
14107c3cd84eSPatrick Sanan     /* Decline to determine compatibility with other DM types */
14117c3cd84eSPatrick Sanan     *set = PETSC_FALSE;
14127c3cd84eSPatrick Sanan   }
14137c3cd84eSPatrick Sanan   PetscFunctionReturn(0);
14147c3cd84eSPatrick Sanan }
1415