xref: /petsc/src/dm/impls/da/dacorn.c (revision d886c4f4307575d4e45cf348c0d185d5b833ebc6)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
347c6ae99SBarry Smith   Code for manipulating distributed regular arrays in parallel.
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
6b45d2f2cSJed Brown #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
747c6ae99SBarry Smith 
847c6ae99SBarry Smith #undef __FUNCT__
96636e97aSMatthew G Knepley #define __FUNCT__ "DMCreateCoordinateDM_DA"
106636e97aSMatthew G Knepley PetscErrorCode DMCreateCoordinateDM_DA(DM dm, DM *cdm)
1147c6ae99SBarry Smith {
126636e97aSMatthew G Knepley   DM_DA         *da = (DM_DA *) dm->data;
1347c6ae99SBarry Smith   PetscMPIInt    size;
1447c6ae99SBarry Smith   PetscErrorCode ierr;
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith   PetscFunctionBegin;
176636e97aSMatthew G Knepley   ierr = MPI_Comm_size(((PetscObject) dm)->comm, &size);CHKERRQ(ierr);
186636e97aSMatthew G Knepley   if (da->dim == 1) {
1947c6ae99SBarry Smith     PetscInt         s,m,*lc,l;
201321219cSEthan Coon     DMDABoundaryType bx;
216636e97aSMatthew G Knepley 
226636e97aSMatthew G Knepley     ierr = DMDAGetInfo(dm,0,&m,0,0,0,0,0,0,&s,&bx,0,0,0);CHKERRQ(ierr);
236636e97aSMatthew G Knepley     ierr = DMDAGetCorners(dm,0,0,0,&l,0,0);CHKERRQ(ierr);
2447c6ae99SBarry Smith     ierr = PetscMalloc(size * sizeof(PetscInt), &lc);CHKERRQ(ierr);
256636e97aSMatthew G Knepley     ierr = MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
266636e97aSMatthew G Knepley     ierr = DMDACreate1d(((PetscObject)dm)->comm,bx,m,1,s,lc,cdm);CHKERRQ(ierr);
2747c6ae99SBarry Smith     ierr = PetscFree(lc);CHKERRQ(ierr);
286636e97aSMatthew G Knepley   } else if (da->dim == 2) {
2947c6ae99SBarry Smith     PetscInt         i,s,m,*lc,*ld,l,k,n,M,N;
301321219cSEthan Coon     DMDABoundaryType bx,by;
316636e97aSMatthew G Knepley 
326636e97aSMatthew G Knepley     ierr = DMDAGetInfo(dm,0,&m,&n,0,&M,&N,0,0,&s,&bx,&by,0,0);CHKERRQ(ierr);
336636e97aSMatthew G Knepley     ierr = DMDAGetCorners(dm,0,0,0,&l,&k,0);CHKERRQ(ierr);
3447c6ae99SBarry Smith     ierr = PetscMalloc2(size,PetscInt,&lc,size,PetscInt,&ld);CHKERRQ(ierr);
3547c6ae99SBarry Smith     /* only first M values in lc matter */
366636e97aSMatthew G Knepley     ierr = MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
3747c6ae99SBarry Smith     /* every Mth value in ld matters */
386636e97aSMatthew G Knepley       ierr = MPI_Allgather(&k,1,MPIU_INT,ld,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
396636e97aSMatthew G Knepley       for(i = 0; i < N; ++i) {
4047c6ae99SBarry Smith         ld[i] = ld[M*i];
4147c6ae99SBarry Smith       }
426636e97aSMatthew G Knepley       ierr = DMDACreate2d(((PetscObject)dm)->comm,bx,by,DMDA_STENCIL_BOX,m,n,M,N,2,s,lc,ld,cdm);CHKERRQ(ierr);
4347c6ae99SBarry Smith       ierr = PetscFree2(lc,ld);CHKERRQ(ierr);
446636e97aSMatthew G Knepley   } else if (da->dim == 3) {
4547c6ae99SBarry Smith     PetscInt         i,s,m,*lc,*ld,*le,l,k,q,n,M,N,P,p;
461321219cSEthan Coon     DMDABoundaryType bx,by,bz;
476636e97aSMatthew G Knepley 
486636e97aSMatthew G Knepley     ierr = DMDAGetInfo(dm,0,&m,&n,&p,&M,&N,&P,0,&s,&bx,&by,&bz,0);CHKERRQ(ierr);
496636e97aSMatthew G Knepley     ierr = DMDAGetCorners(dm,0,0,0,&l,&k,&q);CHKERRQ(ierr);
5047c6ae99SBarry Smith     ierr = PetscMalloc3(size,PetscInt,&lc,size,PetscInt,&ld,size,PetscInt,&le);CHKERRQ(ierr);
5147c6ae99SBarry Smith     /* only first M values in lc matter */
526636e97aSMatthew G Knepley     ierr = MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
5347c6ae99SBarry Smith     /* every Mth value in ld matters */
546636e97aSMatthew G Knepley     ierr = MPI_Allgather(&k,1,MPIU_INT,ld,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
556636e97aSMatthew G Knepley     for(i = 0; i < N; ++i) {
5647c6ae99SBarry Smith       ld[i] = ld[M*i];
5747c6ae99SBarry Smith     }
586636e97aSMatthew G Knepley     ierr = MPI_Allgather(&q,1,MPIU_INT,le,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
596636e97aSMatthew G Knepley     for(i = 0; i < P; ++i) {
6047c6ae99SBarry Smith       le[i] = le[M*N*i];
6147c6ae99SBarry Smith     }
626636e97aSMatthew G Knepley     ierr = DMDACreate3d(((PetscObject)dm)->comm,bx,by,bz,DMDA_STENCIL_BOX,m,n,p,M,N,P,3,s,lc,ld,le,cdm);CHKERRQ(ierr);
6347c6ae99SBarry Smith     ierr = PetscFree3(lc,ld,le);CHKERRQ(ierr);
6447c6ae99SBarry Smith   }
6547c6ae99SBarry Smith   PetscFunctionReturn(0);
6647c6ae99SBarry Smith }
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith #undef __FUNCT__
69aa219208SBarry Smith #define __FUNCT__ "DMDASetFieldName"
7047c6ae99SBarry Smith /*@C
71aa219208SBarry Smith    DMDASetFieldName - Sets the names of individual field components in multicomponent
72aa219208SBarry Smith    vectors associated with a DMDA.
7347c6ae99SBarry Smith 
7447c6ae99SBarry Smith    Not Collective
7547c6ae99SBarry Smith 
7647c6ae99SBarry Smith    Input Parameters:
7747c6ae99SBarry Smith +  da - the distributed array
78aa219208SBarry Smith .  nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the
79aa219208SBarry Smith         number of degrees of freedom per node within the DMDA
8047c6ae99SBarry Smith -  names - the name of the field (component)
8147c6ae99SBarry Smith 
8247c6ae99SBarry Smith   Level: intermediate
8347c6ae99SBarry Smith 
8447c6ae99SBarry Smith .keywords: distributed array, get, component name
8547c6ae99SBarry Smith 
86aa219208SBarry Smith .seealso: DMDAGetFieldName()
8747c6ae99SBarry Smith @*/
887087cfbeSBarry Smith PetscErrorCode  DMDASetFieldName(DM da,PetscInt nf,const char name[])
8947c6ae99SBarry Smith {
9047c6ae99SBarry Smith   PetscErrorCode ierr;
9147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
9247c6ae99SBarry Smith 
9347c6ae99SBarry Smith   PetscFunctionBegin;
9447c6ae99SBarry Smith    PetscValidHeaderSpecific(da,DM_CLASSID,1);
9547c6ae99SBarry Smith   if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf);
96c31cb41cSBarry Smith   ierr = PetscFree(dd->fieldname[nf]);CHKERRQ(ierr);
9747c6ae99SBarry Smith   ierr = PetscStrallocpy(name,&dd->fieldname[nf]);CHKERRQ(ierr);
9847c6ae99SBarry Smith   PetscFunctionReturn(0);
9947c6ae99SBarry Smith }
10047c6ae99SBarry Smith 
10147c6ae99SBarry Smith #undef __FUNCT__
102aa219208SBarry Smith #define __FUNCT__ "DMDAGetFieldName"
10347c6ae99SBarry Smith /*@C
104aa219208SBarry Smith    DMDAGetFieldName - Gets the names of individual field components in multicomponent
105aa219208SBarry Smith    vectors associated with a DMDA.
10647c6ae99SBarry Smith 
10747c6ae99SBarry Smith    Not Collective
10847c6ae99SBarry Smith 
10947c6ae99SBarry Smith    Input Parameter:
11047c6ae99SBarry Smith +  da - the distributed array
111aa219208SBarry Smith -  nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the
112aa219208SBarry Smith         number of degrees of freedom per node within the DMDA
11347c6ae99SBarry Smith 
11447c6ae99SBarry Smith    Output Parameter:
11547c6ae99SBarry Smith .  names - the name of the field (component)
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith   Level: intermediate
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith .keywords: distributed array, get, component name
12047c6ae99SBarry Smith 
121aa219208SBarry Smith .seealso: DMDASetFieldName()
12247c6ae99SBarry Smith @*/
1237087cfbeSBarry Smith PetscErrorCode  DMDAGetFieldName(DM da,PetscInt nf,const char **name)
12447c6ae99SBarry Smith {
12547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith   PetscFunctionBegin;
12847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
12947c6ae99SBarry Smith   PetscValidPointer(name,3);
13047c6ae99SBarry Smith   if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf);
13147c6ae99SBarry Smith   *name = dd->fieldname[nf];
13247c6ae99SBarry Smith   PetscFunctionReturn(0);
13347c6ae99SBarry Smith }
13447c6ae99SBarry Smith 
13547c6ae99SBarry Smith #undef __FUNCT__
136aa219208SBarry Smith #define __FUNCT__ "DMDAGetCorners"
13747c6ae99SBarry Smith /*@
138aa219208SBarry Smith    DMDAGetCorners - Returns the global (x,y,z) indices of the lower left
13947c6ae99SBarry Smith    corner of the local region, excluding ghost points.
14047c6ae99SBarry Smith 
14147c6ae99SBarry Smith    Not Collective
14247c6ae99SBarry Smith 
14347c6ae99SBarry Smith    Input Parameter:
14447c6ae99SBarry Smith .  da - the distributed array
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith    Output Parameters:
14747c6ae99SBarry Smith +  x,y,z - the corner indices (where y and z are optional; these are used
14847c6ae99SBarry Smith            for 2D and 3D problems)
14947c6ae99SBarry Smith -  m,n,p - widths in the corresponding directions (where n and p are optional;
15047c6ae99SBarry Smith            these are used for 2D and 3D problems)
15147c6ae99SBarry Smith 
15247c6ae99SBarry Smith    Note:
15347c6ae99SBarry Smith    The corner information is independent of the number of degrees of
154aa219208SBarry Smith    freedom per node set with the DMDACreateXX() routine. Thus the x, y, z, and
15547c6ae99SBarry Smith    m, n, p can be thought of as coordinates on a logical grid, where each
15647c6ae99SBarry Smith    grid point has (potentially) several degrees of freedom.
15747c6ae99SBarry Smith    Any of y, z, n, and p can be passed in as PETSC_NULL if not needed.
15847c6ae99SBarry Smith 
15947c6ae99SBarry Smith   Level: beginner
16047c6ae99SBarry Smith 
16147c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices
16247c6ae99SBarry Smith 
163aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetOwnershipRanges()
16447c6ae99SBarry Smith @*/
1657087cfbeSBarry Smith PetscErrorCode  DMDAGetCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p)
16647c6ae99SBarry Smith {
16747c6ae99SBarry Smith   PetscInt w;
16847c6ae99SBarry Smith   DM_DA    *dd = (DM_DA*)da->data;
16947c6ae99SBarry Smith 
17047c6ae99SBarry Smith   PetscFunctionBegin;
17147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
17247c6ae99SBarry Smith   /* since the xs, xe ... have all been multiplied by the number of degrees
17347c6ae99SBarry Smith      of freedom per cell, w = dd->w, we divide that out before returning.*/
17447c6ae99SBarry Smith   w = dd->w;
175*d886c4f4SPeter Brune   if (x) *x = dd->xs/w + dd->xo; if (m) *m = (dd->xe - dd->xs)/w;
17647c6ae99SBarry Smith   /* the y and z have NOT been multiplied by w */
177*d886c4f4SPeter Brune   if (y) *y = dd->ys + dd->yo;   if (n) *n = (dd->ye - dd->ys);
178*d886c4f4SPeter Brune   if (z) *z = dd->zs + dd->zo;   if (p) *p = (dd->ze - dd->zs);
17947c6ae99SBarry Smith   PetscFunctionReturn(0);
18047c6ae99SBarry Smith }
18147c6ae99SBarry Smith 
18247c6ae99SBarry Smith #undef __FUNCT__
183aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalBoundingBox"
18447c6ae99SBarry Smith /*@
185aa219208SBarry Smith    DMDAGetLocalBoundingBox - Returns the local bounding box for the DMDA.
18647c6ae99SBarry Smith 
18747c6ae99SBarry Smith    Not Collective
18847c6ae99SBarry Smith 
18947c6ae99SBarry Smith    Input Parameter:
19047c6ae99SBarry Smith .  da - the distributed array
19147c6ae99SBarry Smith 
19247c6ae99SBarry Smith    Output Parameters:
19347c6ae99SBarry Smith +  lmin - local minimum coordinates (length dim, optional)
19447c6ae99SBarry Smith -  lmax - local maximim coordinates (length dim, optional)
19547c6ae99SBarry Smith 
19647c6ae99SBarry Smith   Level: beginner
19747c6ae99SBarry Smith 
19847c6ae99SBarry Smith .keywords: distributed array, get, coordinates
19947c6ae99SBarry Smith 
2002150357eSBarry Smith .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetBoundingBox()
20147c6ae99SBarry Smith @*/
2027087cfbeSBarry Smith PetscErrorCode  DMDAGetLocalBoundingBox(DM da,PetscReal lmin[],PetscReal lmax[])
20347c6ae99SBarry Smith {
20447c6ae99SBarry Smith   PetscErrorCode    ierr;
20547c6ae99SBarry Smith   Vec               coords  = PETSC_NULL;
20647c6ae99SBarry Smith   PetscInt          dim,i,j;
20747c6ae99SBarry Smith   const PetscScalar *local_coords;
208ea345e14SBarry Smith   PetscReal         min[3]={PETSC_MAX_REAL,PETSC_MAX_REAL,PETSC_MAX_REAL},max[3]={PETSC_MIN_REAL,PETSC_MIN_REAL,PETSC_MIN_REAL};
20947c6ae99SBarry Smith   PetscInt          N,Ni;
21047c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
21147c6ae99SBarry Smith 
21247c6ae99SBarry Smith   PetscFunctionBegin;
21347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
21447c6ae99SBarry Smith   dim = dd->dim;
2156636e97aSMatthew G Knepley   ierr = DMGetCoordinates(da,&coords);CHKERRQ(ierr);
2167324c66bSJed Brown   if (coords) {
21747c6ae99SBarry Smith     ierr = VecGetArrayRead(coords,&local_coords);CHKERRQ(ierr);
21847c6ae99SBarry Smith     ierr = VecGetLocalSize(coords,&N);CHKERRQ(ierr);
21947c6ae99SBarry Smith     Ni = N/dim;
22047c6ae99SBarry Smith     for (i=0; i<Ni; i++) {
2217324c66bSJed Brown       for (j=0; j<3; j++) {
2227324c66bSJed Brown         min[j] = j < dim ? PetscMin(min[j],PetscRealPart(local_coords[i*dim+j])) : 0;
2237324c66bSJed Brown         max[j] = j < dim ? PetscMax(min[j],PetscRealPart(local_coords[i*dim+j])) : 0;
22447c6ae99SBarry Smith       }
22547c6ae99SBarry Smith     }
22647c6ae99SBarry Smith     ierr = VecRestoreArrayRead(coords,&local_coords);CHKERRQ(ierr);
2277324c66bSJed Brown   } else {                      /* Just use grid indices */
2287324c66bSJed Brown     DMDALocalInfo info;
2297324c66bSJed Brown     ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
2307324c66bSJed Brown     min[0] = info.xs;
2317324c66bSJed Brown     min[1] = info.ys;
2327324c66bSJed Brown     min[2] = info.zs;
2337324c66bSJed Brown     max[0] = info.xs + info.xm-1;
2347324c66bSJed Brown     max[1] = info.ys + info.ym-1;
2357324c66bSJed Brown     max[2] = info.zs + info.zm-1;
2367324c66bSJed Brown   }
23747c6ae99SBarry Smith   if (lmin) {ierr = PetscMemcpy(lmin,min,dim*sizeof(PetscReal));CHKERRQ(ierr);}
23847c6ae99SBarry Smith   if (lmax) {ierr = PetscMemcpy(lmax,max,dim*sizeof(PetscReal));CHKERRQ(ierr);}
23947c6ae99SBarry Smith   PetscFunctionReturn(0);
24047c6ae99SBarry Smith }
24147c6ae99SBarry Smith 
24247c6ae99SBarry Smith #undef __FUNCT__
243aa219208SBarry Smith #define __FUNCT__ "DMDAGetBoundingBox"
24447c6ae99SBarry Smith /*@
245aa219208SBarry Smith    DMDAGetBoundingBox - Returns the global bounding box for the DMDA.
24647c6ae99SBarry Smith 
247aa219208SBarry Smith    Collective on DMDA
24847c6ae99SBarry Smith 
24947c6ae99SBarry Smith    Input Parameter:
25047c6ae99SBarry Smith .  da - the distributed array
25147c6ae99SBarry Smith 
25247c6ae99SBarry Smith    Output Parameters:
25347c6ae99SBarry Smith +  gmin - global minimum coordinates (length dim, optional)
25447c6ae99SBarry Smith -  gmax - global maximim coordinates (length dim, optional)
25547c6ae99SBarry Smith 
25647c6ae99SBarry Smith   Level: beginner
25747c6ae99SBarry Smith 
25847c6ae99SBarry Smith .keywords: distributed array, get, coordinates
25947c6ae99SBarry Smith 
2602150357eSBarry Smith .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetLocalBoundingBox()
26147c6ae99SBarry Smith @*/
2627087cfbeSBarry Smith PetscErrorCode  DMDAGetBoundingBox(DM da,PetscReal gmin[],PetscReal gmax[])
26347c6ae99SBarry Smith {
26447c6ae99SBarry Smith   PetscErrorCode ierr;
26547c6ae99SBarry Smith   PetscMPIInt    count;
26647c6ae99SBarry Smith   PetscReal      lmin[3],lmax[3];
26747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
26847c6ae99SBarry Smith 
26947c6ae99SBarry Smith   PetscFunctionBegin;
27047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
27147c6ae99SBarry Smith   count = PetscMPIIntCast(dd->dim);
272aa219208SBarry Smith   ierr = DMDAGetLocalBoundingBox(da,lmin,lmax);CHKERRQ(ierr);
273d9822059SBarry Smith   if (gmin) {ierr = MPI_Allreduce(lmin,gmin,count,MPIU_REAL,MPIU_MIN,((PetscObject)da)->comm);CHKERRQ(ierr);}
274d9822059SBarry Smith   if (gmax) {ierr = MPI_Allreduce(lmax,gmax,count,MPIU_REAL,MPIU_MAX,((PetscObject)da)->comm);CHKERRQ(ierr);}
27547c6ae99SBarry Smith   PetscFunctionReturn(0);
27647c6ae99SBarry Smith }
277bc2bf880SBarry Smith 
278bc2bf880SBarry Smith #undef __FUNCT__
279bc2bf880SBarry Smith #define __FUNCT__ "DMDAGetReducedDA"
280bc2bf880SBarry Smith /*@
281bc2bf880SBarry Smith    DMDAGetReducedDA - Gets the DMDA with the same layout but with fewer or more fields
282bc2bf880SBarry Smith 
283bc2bf880SBarry Smith    Collective on DMDA
284bc2bf880SBarry Smith 
285bc2bf880SBarry Smith    Input Parameter:
286bc2bf880SBarry Smith +  da - the distributed array
287bc2bf880SBarry Smith .  nfields - number of fields in new DMDA
288bc2bf880SBarry Smith 
289bc2bf880SBarry Smith    Output Parameter:
290bc2bf880SBarry Smith .  nda - the new DMDA
291bc2bf880SBarry Smith 
292bc2bf880SBarry Smith   Level: intermediate
293bc2bf880SBarry Smith 
294bc2bf880SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
295bc2bf880SBarry Smith 
2962150357eSBarry Smith .seealso: DMDAGetGhostCorners(), DMSetCoordinates(), DMDASetUniformCoordinates(), DMGetCoordinates(), DMDAGetGhostedCoordinates()
297bc2bf880SBarry Smith @*/
298bc2bf880SBarry Smith PetscErrorCode  DMDAGetReducedDA(DM da,PetscInt nfields,DM *nda)
299bc2bf880SBarry Smith {
300bc2bf880SBarry Smith   PetscErrorCode   ierr;
301bc2bf880SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
302320964c4SBlaise Bourdin   PetscInt         s,m,n,p,M,N,P,dim;
303320964c4SBlaise Bourdin   const PetscInt   *lx,*ly,*lz;
304bc2bf880SBarry Smith   DMDABoundaryType bx,by,bz;
305320964c4SBlaise Bourdin   DMDAStencilType  stencil_type;
306320964c4SBlaise Bourdin 
307320964c4SBlaise Bourdin   PetscFunctionBegin;
308320964c4SBlaise Bourdin   ierr = DMDAGetInfo(da,&dim,&M,&N,&P,&m,&n,&p,0,&s,&bx,&by,&bz,&stencil_type);CHKERRQ(ierr);
309320964c4SBlaise Bourdin   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr);
310320964c4SBlaise Bourdin   if (dim == 1) {
311320964c4SBlaise Bourdin     ierr = DMDACreate1d(((PetscObject)da)->comm,bx,M,nfields,s,dd->lx,nda);CHKERRQ(ierr);
312320964c4SBlaise Bourdin   } else if (dim == 2) {
313320964c4SBlaise Bourdin     ierr = DMDACreate2d(((PetscObject)da)->comm,bx,by,stencil_type,M,N,m,n,nfields,s,lx,ly,nda);CHKERRQ(ierr);
314320964c4SBlaise Bourdin   } else if (dim == 3) {
315320964c4SBlaise Bourdin     ierr = DMDACreate3d(((PetscObject)da)->comm,bx,by,bz,stencil_type,M,N,P,m,n,p,nfields,s,lx,ly,lz,nda);CHKERRQ(ierr);
316bc2bf880SBarry Smith   }
3176636e97aSMatthew G Knepley   if (da->coordinates) {
3186636e97aSMatthew G Knepley     ierr        = PetscObjectReference((PetscObject)da->coordinates);CHKERRQ(ierr);
3196636e97aSMatthew G Knepley     (*nda)->coordinates = da->coordinates;
320bc2bf880SBarry Smith   }
321bc2bf880SBarry Smith   PetscFunctionReturn(0);
322bc2bf880SBarry Smith }
323bc2bf880SBarry Smith 
324