xref: /petsc/src/dm/impls/da/dacorn.c (revision 8865f1ea4ea560cd84ab8db62e98b7095cdff96f)
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);
39*8865f1eaSKarl Rupp       for (i = 0; i < N; ++i) ld[i] = ld[M*i];
40cc76e2cfSBarry Smith       if (bx == DMDA_BOUNDARY_MIRROR || by == DMDA_BOUNDARY_MIRROR) {
41cc76e2cfSBarry Smith         ierr = DMDACreate2d(((PetscObject)dm)->comm,bx,by,DMDA_STENCIL_STAR,m,n,M,N,2,s,lc,ld,cdm);CHKERRQ(ierr);
42cc76e2cfSBarry Smith       } else {
436636e97aSMatthew G Knepley         ierr = DMDACreate2d(((PetscObject)dm)->comm,bx,by,DMDA_STENCIL_BOX,m,n,M,N,2,s,lc,ld,cdm);CHKERRQ(ierr);
44cc76e2cfSBarry Smith       }
4547c6ae99SBarry Smith       ierr = PetscFree2(lc,ld);CHKERRQ(ierr);
466636e97aSMatthew G Knepley   } else if (da->dim == 3) {
4747c6ae99SBarry Smith     PetscInt         i,s,m,*lc,*ld,*le,l,k,q,n,M,N,P,p;
481321219cSEthan Coon     DMDABoundaryType bx,by,bz;
496636e97aSMatthew G Knepley 
506636e97aSMatthew G Knepley     ierr = DMDAGetInfo(dm,0,&m,&n,&p,&M,&N,&P,0,&s,&bx,&by,&bz,0);CHKERRQ(ierr);
516636e97aSMatthew G Knepley     ierr = DMDAGetCorners(dm,0,0,0,&l,&k,&q);CHKERRQ(ierr);
5247c6ae99SBarry Smith     ierr = PetscMalloc3(size,PetscInt,&lc,size,PetscInt,&ld,size,PetscInt,&le);CHKERRQ(ierr);
5347c6ae99SBarry Smith     /* only first M values in lc matter */
546636e97aSMatthew G Knepley     ierr = MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
5547c6ae99SBarry Smith     /* every Mth value in ld matters */
566636e97aSMatthew G Knepley     ierr = MPI_Allgather(&k,1,MPIU_INT,ld,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
57*8865f1eaSKarl Rupp     for (i = 0; i < N; ++i) ld[i] = ld[M*i];
586636e97aSMatthew G Knepley     ierr = MPI_Allgather(&q,1,MPIU_INT,le,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
59*8865f1eaSKarl Rupp     for (i = 0; i < P; ++i) le[i] = le[M*N*i];
606636e97aSMatthew 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);
6147c6ae99SBarry Smith     ierr = PetscFree3(lc,ld,le);CHKERRQ(ierr);
6247c6ae99SBarry Smith   }
6347c6ae99SBarry Smith   PetscFunctionReturn(0);
6447c6ae99SBarry Smith }
6547c6ae99SBarry Smith 
6647c6ae99SBarry Smith #undef __FUNCT__
67aa219208SBarry Smith #define __FUNCT__ "DMDASetFieldName"
6847c6ae99SBarry Smith /*@C
69aa219208SBarry Smith    DMDASetFieldName - Sets the names of individual field components in multicomponent
70aa219208SBarry Smith    vectors associated with a DMDA.
7147c6ae99SBarry Smith 
7247c6ae99SBarry Smith    Not Collective
7347c6ae99SBarry Smith 
7447c6ae99SBarry Smith    Input Parameters:
7547c6ae99SBarry Smith +  da - the distributed array
76aa219208SBarry Smith .  nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the
77aa219208SBarry Smith         number of degrees of freedom per node within the DMDA
7847c6ae99SBarry Smith -  names - the name of the field (component)
7947c6ae99SBarry Smith 
8047c6ae99SBarry Smith   Level: intermediate
8147c6ae99SBarry Smith 
8247c6ae99SBarry Smith .keywords: distributed array, get, component name
8347c6ae99SBarry Smith 
84109c9344SBarry Smith .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName()
8547c6ae99SBarry Smith @*/
867087cfbeSBarry Smith PetscErrorCode  DMDASetFieldName(DM da,PetscInt nf,const char name[])
8747c6ae99SBarry Smith {
8847c6ae99SBarry Smith   PetscErrorCode ierr;
8947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
9047c6ae99SBarry Smith 
9147c6ae99SBarry Smith   PetscFunctionBegin;
9247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
9347c6ae99SBarry Smith   if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf);
94c31cb41cSBarry Smith   ierr = PetscFree(dd->fieldname[nf]);CHKERRQ(ierr);
9547c6ae99SBarry Smith   ierr = PetscStrallocpy(name,&dd->fieldname[nf]);CHKERRQ(ierr);
9647c6ae99SBarry Smith   PetscFunctionReturn(0);
9747c6ae99SBarry Smith }
9847c6ae99SBarry Smith 
9947c6ae99SBarry Smith #undef __FUNCT__
100aa219208SBarry Smith #define __FUNCT__ "DMDAGetFieldName"
10147c6ae99SBarry Smith /*@C
102aa219208SBarry Smith    DMDAGetFieldName - Gets the names of individual field components in multicomponent
103aa219208SBarry Smith    vectors associated with a DMDA.
10447c6ae99SBarry Smith 
10547c6ae99SBarry Smith    Not Collective
10647c6ae99SBarry Smith 
10747c6ae99SBarry Smith    Input Parameter:
10847c6ae99SBarry Smith +  da - the distributed array
109aa219208SBarry Smith -  nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the
110aa219208SBarry Smith         number of degrees of freedom per node within the DMDA
11147c6ae99SBarry Smith 
11247c6ae99SBarry Smith    Output Parameter:
11347c6ae99SBarry Smith .  names - the name of the field (component)
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith   Level: intermediate
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith .keywords: distributed array, get, component name
11847c6ae99SBarry Smith 
119109c9344SBarry Smith .seealso: DMDASetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName()
12047c6ae99SBarry Smith @*/
1217087cfbeSBarry Smith PetscErrorCode  DMDAGetFieldName(DM da,PetscInt nf,const char **name)
12247c6ae99SBarry Smith {
12347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
12447c6ae99SBarry Smith 
12547c6ae99SBarry Smith   PetscFunctionBegin;
12647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
12747c6ae99SBarry Smith   PetscValidPointer(name,3);
12847c6ae99SBarry Smith   if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf);
12947c6ae99SBarry Smith   *name = dd->fieldname[nf];
13047c6ae99SBarry Smith   PetscFunctionReturn(0);
13147c6ae99SBarry Smith }
13247c6ae99SBarry Smith 
13347c6ae99SBarry Smith #undef __FUNCT__
134109c9344SBarry Smith #define __FUNCT__ "DMDASetCoordinateName"
135109c9344SBarry Smith /*@C
136109c9344SBarry Smith    DMDASetCoordinateName - Sets the name of the coordinate directions associated with a DMDA, for example "x" or "y"
137109c9344SBarry Smith 
138109c9344SBarry Smith    Not Collective
139109c9344SBarry Smith 
140109c9344SBarry Smith    Input Parameters:
141109c9344SBarry Smith +  da - the distributed array
142109c9344SBarry Smith .  nf - coordinate number for the DMDA (0, 1, ... dim-1),
143109c9344SBarry Smith -  name - the name of the coordinate
144109c9344SBarry Smith 
145109c9344SBarry Smith   Level: intermediate
146109c9344SBarry Smith 
147109c9344SBarry Smith .keywords: distributed array, get, component name
148109c9344SBarry Smith 
149109c9344SBarry Smith .seealso: DMDAGetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName()
150109c9344SBarry Smith @*/
151109c9344SBarry Smith PetscErrorCode  DMDASetCoordinateName(DM da,PetscInt nf,const char name[])
152109c9344SBarry Smith {
153109c9344SBarry Smith   PetscErrorCode ierr;
154109c9344SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
155109c9344SBarry Smith 
156109c9344SBarry Smith   PetscFunctionBegin;
157109c9344SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
158109c9344SBarry Smith   if (nf < 0 || nf >= dd->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf);
159109c9344SBarry Smith   ierr = PetscFree(dd->coordinatename[nf]);CHKERRQ(ierr);
160109c9344SBarry Smith   ierr = PetscStrallocpy(name,&dd->coordinatename[nf]);CHKERRQ(ierr);
161109c9344SBarry Smith   PetscFunctionReturn(0);
162109c9344SBarry Smith }
163109c9344SBarry Smith 
164109c9344SBarry Smith #undef __FUNCT__
165109c9344SBarry Smith #define __FUNCT__ "DMDAGetCoordinateName"
166109c9344SBarry Smith /*@C
167109c9344SBarry Smith    DMDAGetCoordinateName - Gets the name of a coodinate direction associated with a DMDA.
168109c9344SBarry Smith 
169109c9344SBarry Smith    Not Collective
170109c9344SBarry Smith 
171109c9344SBarry Smith    Input Parameter:
172109c9344SBarry Smith +  da - the distributed array
173109c9344SBarry Smith -  nf -  number for the DMDA (0, 1, ... dim-1)
174109c9344SBarry Smith 
175109c9344SBarry Smith    Output Parameter:
176109c9344SBarry Smith .  names - the name of the coordinate direction
177109c9344SBarry Smith 
178109c9344SBarry Smith   Level: intermediate
179109c9344SBarry Smith 
180109c9344SBarry Smith .keywords: distributed array, get, component name
181109c9344SBarry Smith 
182109c9344SBarry Smith .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName()
183109c9344SBarry Smith @*/
184109c9344SBarry Smith PetscErrorCode  DMDAGetCoordinateName(DM da,PetscInt nf,const char **name)
185109c9344SBarry Smith {
186109c9344SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
187109c9344SBarry Smith 
188109c9344SBarry Smith   PetscFunctionBegin;
189109c9344SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
190109c9344SBarry Smith   PetscValidPointer(name,3);
191109c9344SBarry Smith   if (nf < 0 || nf >= dd->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf);
192109c9344SBarry Smith   *name = dd->coordinatename[nf];
193109c9344SBarry Smith   PetscFunctionReturn(0);
194109c9344SBarry Smith }
195109c9344SBarry Smith 
196109c9344SBarry Smith #undef __FUNCT__
197aa219208SBarry Smith #define __FUNCT__ "DMDAGetCorners"
19847c6ae99SBarry Smith /*@
199aa219208SBarry Smith    DMDAGetCorners - Returns the global (x,y,z) indices of the lower left
20047c6ae99SBarry Smith    corner of the local region, excluding ghost points.
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith    Not Collective
20347c6ae99SBarry Smith 
20447c6ae99SBarry Smith    Input Parameter:
20547c6ae99SBarry Smith .  da - the distributed array
20647c6ae99SBarry Smith 
20747c6ae99SBarry Smith    Output Parameters:
20847c6ae99SBarry Smith +  x,y,z - the corner indices (where y and z are optional; these are used
20947c6ae99SBarry Smith            for 2D and 3D problems)
21047c6ae99SBarry Smith -  m,n,p - widths in the corresponding directions (where n and p are optional;
21147c6ae99SBarry Smith            these are used for 2D and 3D problems)
21247c6ae99SBarry Smith 
21347c6ae99SBarry Smith    Note:
21447c6ae99SBarry Smith    The corner information is independent of the number of degrees of
215aa219208SBarry Smith    freedom per node set with the DMDACreateXX() routine. Thus the x, y, z, and
21647c6ae99SBarry Smith    m, n, p can be thought of as coordinates on a logical grid, where each
21747c6ae99SBarry Smith    grid point has (potentially) several degrees of freedom.
21847c6ae99SBarry Smith    Any of y, z, n, and p can be passed in as PETSC_NULL if not needed.
21947c6ae99SBarry Smith 
22047c6ae99SBarry Smith   Level: beginner
22147c6ae99SBarry Smith 
22247c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices
22347c6ae99SBarry Smith 
224aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetOwnershipRanges()
22547c6ae99SBarry Smith @*/
2267087cfbeSBarry Smith PetscErrorCode  DMDAGetCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p)
22747c6ae99SBarry Smith {
22847c6ae99SBarry Smith   PetscInt w;
22947c6ae99SBarry Smith   DM_DA    *dd = (DM_DA*)da->data;
23047c6ae99SBarry Smith 
23147c6ae99SBarry Smith   PetscFunctionBegin;
23247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
23347c6ae99SBarry Smith   /* since the xs, xe ... have all been multiplied by the number of degrees
23447c6ae99SBarry Smith      of freedom per cell, w = dd->w, we divide that out before returning.*/
23547c6ae99SBarry Smith   w = dd->w;
236d886c4f4SPeter Brune   if (x) *x = dd->xs/w + dd->xo; if (m) *m = (dd->xe - dd->xs)/w;
23747c6ae99SBarry Smith   /* the y and z have NOT been multiplied by w */
238d886c4f4SPeter Brune   if (y) *y = dd->ys + dd->yo; if (n) *n = (dd->ye - dd->ys);
239d886c4f4SPeter Brune   if (z) *z = dd->zs + dd->zo; if (p) *p = (dd->ze - dd->zs);
24047c6ae99SBarry Smith   PetscFunctionReturn(0);
24147c6ae99SBarry Smith }
24247c6ae99SBarry Smith 
24347c6ae99SBarry Smith #undef __FUNCT__
244aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalBoundingBox"
24547c6ae99SBarry Smith /*@
246aa219208SBarry Smith    DMDAGetLocalBoundingBox - Returns the local bounding box for the DMDA.
24747c6ae99SBarry Smith 
24847c6ae99SBarry Smith    Not Collective
24947c6ae99SBarry Smith 
25047c6ae99SBarry Smith    Input Parameter:
25147c6ae99SBarry Smith .  da - the distributed array
25247c6ae99SBarry Smith 
25347c6ae99SBarry Smith    Output Parameters:
25447c6ae99SBarry Smith +  lmin - local minimum coordinates (length dim, optional)
25547c6ae99SBarry Smith -  lmax - local maximim coordinates (length dim, optional)
25647c6ae99SBarry Smith 
25747c6ae99SBarry Smith   Level: beginner
25847c6ae99SBarry Smith 
25947c6ae99SBarry Smith .keywords: distributed array, get, coordinates
26047c6ae99SBarry Smith 
2612150357eSBarry Smith .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetBoundingBox()
26247c6ae99SBarry Smith @*/
2637087cfbeSBarry Smith PetscErrorCode  DMDAGetLocalBoundingBox(DM da,PetscReal lmin[],PetscReal lmax[])
26447c6ae99SBarry Smith {
26547c6ae99SBarry Smith   PetscErrorCode    ierr;
26647c6ae99SBarry Smith   Vec               coords = PETSC_NULL;
26747c6ae99SBarry Smith   PetscInt          dim,i,j;
26847c6ae99SBarry Smith   const PetscScalar *local_coords;
269ea345e14SBarry Smith   PetscReal         min[3]={PETSC_MAX_REAL,PETSC_MAX_REAL,PETSC_MAX_REAL},max[3]={PETSC_MIN_REAL,PETSC_MIN_REAL,PETSC_MIN_REAL};
27047c6ae99SBarry Smith   PetscInt          N,Ni;
27147c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
27247c6ae99SBarry Smith 
27347c6ae99SBarry Smith   PetscFunctionBegin;
27447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
27547c6ae99SBarry Smith   dim  = dd->dim;
2766636e97aSMatthew G Knepley   ierr = DMGetCoordinates(da,&coords);CHKERRQ(ierr);
2777324c66bSJed Brown   if (coords) {
27847c6ae99SBarry Smith     ierr = VecGetArrayRead(coords,&local_coords);CHKERRQ(ierr);
27947c6ae99SBarry Smith     ierr = VecGetLocalSize(coords,&N);CHKERRQ(ierr);
28047c6ae99SBarry Smith     Ni   = N/dim;
28147c6ae99SBarry Smith     for (i=0; i<Ni; i++) {
2827324c66bSJed Brown       for (j=0; j<3; j++) {
2837324c66bSJed Brown         min[j] = j < dim ? PetscMin(min[j],PetscRealPart(local_coords[i*dim+j])) : 0;
2842197688cSJed Brown         max[j] = j < dim ? PetscMax(max[j],PetscRealPart(local_coords[i*dim+j])) : 0;
28547c6ae99SBarry Smith       }
28647c6ae99SBarry Smith     }
28747c6ae99SBarry Smith     ierr = VecRestoreArrayRead(coords,&local_coords);CHKERRQ(ierr);
2887324c66bSJed Brown   } else {                      /* Just use grid indices */
2897324c66bSJed Brown     DMDALocalInfo info;
2907324c66bSJed Brown     ierr   = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
2917324c66bSJed Brown     min[0] = info.xs;
2927324c66bSJed Brown     min[1] = info.ys;
2937324c66bSJed Brown     min[2] = info.zs;
2947324c66bSJed Brown     max[0] = info.xs + info.xm-1;
2957324c66bSJed Brown     max[1] = info.ys + info.ym-1;
2967324c66bSJed Brown     max[2] = info.zs + info.zm-1;
2977324c66bSJed Brown   }
29847c6ae99SBarry Smith   if (lmin) {ierr = PetscMemcpy(lmin,min,dim*sizeof(PetscReal));CHKERRQ(ierr);}
29947c6ae99SBarry Smith   if (lmax) {ierr = PetscMemcpy(lmax,max,dim*sizeof(PetscReal));CHKERRQ(ierr);}
30047c6ae99SBarry Smith   PetscFunctionReturn(0);
30147c6ae99SBarry Smith }
30247c6ae99SBarry Smith 
30347c6ae99SBarry Smith #undef __FUNCT__
304aa219208SBarry Smith #define __FUNCT__ "DMDAGetBoundingBox"
30547c6ae99SBarry Smith /*@
306aa219208SBarry Smith    DMDAGetBoundingBox - Returns the global bounding box for the DMDA.
30747c6ae99SBarry Smith 
308aa219208SBarry Smith    Collective on DMDA
30947c6ae99SBarry Smith 
31047c6ae99SBarry Smith    Input Parameter:
31147c6ae99SBarry Smith .  da - the distributed array
31247c6ae99SBarry Smith 
31347c6ae99SBarry Smith    Output Parameters:
31447c6ae99SBarry Smith +  gmin - global minimum coordinates (length dim, optional)
31547c6ae99SBarry Smith -  gmax - global maximim coordinates (length dim, optional)
31647c6ae99SBarry Smith 
31747c6ae99SBarry Smith   Level: beginner
31847c6ae99SBarry Smith 
31947c6ae99SBarry Smith .keywords: distributed array, get, coordinates
32047c6ae99SBarry Smith 
3212150357eSBarry Smith .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetLocalBoundingBox()
32247c6ae99SBarry Smith @*/
3237087cfbeSBarry Smith PetscErrorCode  DMDAGetBoundingBox(DM da,PetscReal gmin[],PetscReal gmax[])
32447c6ae99SBarry Smith {
32547c6ae99SBarry Smith   PetscErrorCode ierr;
32647c6ae99SBarry Smith   PetscMPIInt    count;
32747c6ae99SBarry Smith   PetscReal      lmin[3],lmax[3];
32847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
32947c6ae99SBarry Smith 
33047c6ae99SBarry Smith   PetscFunctionBegin;
33147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
3324dc2109aSBarry Smith   ierr = PetscMPIIntCast(dd->dim,&count);CHKERRQ(ierr);
333aa219208SBarry Smith   ierr = DMDAGetLocalBoundingBox(da,lmin,lmax);CHKERRQ(ierr);
334d9822059SBarry Smith   if (gmin) {ierr = MPI_Allreduce(lmin,gmin,count,MPIU_REAL,MPIU_MIN,((PetscObject)da)->comm);CHKERRQ(ierr);}
335d9822059SBarry Smith   if (gmax) {ierr = MPI_Allreduce(lmax,gmax,count,MPIU_REAL,MPIU_MAX,((PetscObject)da)->comm);CHKERRQ(ierr);}
33647c6ae99SBarry Smith   PetscFunctionReturn(0);
33747c6ae99SBarry Smith }
338bc2bf880SBarry Smith 
339bc2bf880SBarry Smith #undef __FUNCT__
340db05f41bSBarry Smith #define __FUNCT__ "DMDAGetReducedDMDA"
341bc2bf880SBarry Smith /*@
342db05f41bSBarry Smith    DMDAGetReducedDMDA - Gets the DMDA with the same layout but with fewer or more fields
343bc2bf880SBarry Smith 
344bc2bf880SBarry Smith    Collective on DMDA
345bc2bf880SBarry Smith 
346bc2bf880SBarry Smith    Input Parameter:
347bc2bf880SBarry Smith +  da - the distributed array
348bc2bf880SBarry Smith .  nfields - number of fields in new DMDA
349bc2bf880SBarry Smith 
350bc2bf880SBarry Smith    Output Parameter:
351bc2bf880SBarry Smith .  nda - the new DMDA
352bc2bf880SBarry Smith 
353bc2bf880SBarry Smith   Level: intermediate
354bc2bf880SBarry Smith 
355bc2bf880SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
356bc2bf880SBarry Smith 
3572150357eSBarry Smith .seealso: DMDAGetGhostCorners(), DMSetCoordinates(), DMDASetUniformCoordinates(), DMGetCoordinates(), DMDAGetGhostedCoordinates()
358bc2bf880SBarry Smith @*/
359db05f41bSBarry Smith PetscErrorCode  DMDAGetReducedDMDA(DM da,PetscInt nfields,DM *nda)
360bc2bf880SBarry Smith {
361bc2bf880SBarry Smith   PetscErrorCode   ierr;
362bc2bf880SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
36395c13181SPeter Brune   PetscInt         s,m,n,p,M,N,P,dim,Mo,No,Po;
364320964c4SBlaise Bourdin   const PetscInt   *lx,*ly,*lz;
365bc2bf880SBarry Smith   DMDABoundaryType bx,by,bz;
366320964c4SBlaise Bourdin   DMDAStencilType  stencil_type;
36795c13181SPeter Brune   PetscInt         ox,oy,oz;
36895c13181SPeter Brune   PetscInt         cl,rl;
369320964c4SBlaise Bourdin 
370320964c4SBlaise Bourdin   PetscFunctionBegin;
37195c13181SPeter Brune   dim = dd->dim;
37295c13181SPeter Brune   M   = dd->M;
37395c13181SPeter Brune   N   = dd->N;
37495c13181SPeter Brune   P   = dd->P;
37595c13181SPeter Brune   m   = dd->m;
37695c13181SPeter Brune   n   = dd->n;
37795c13181SPeter Brune   p   = dd->p;
37895c13181SPeter Brune   s   = dd->s;
37995c13181SPeter Brune   bx  = dd->bx;
38095c13181SPeter Brune   by  = dd->by;
38195c13181SPeter Brune   bz  = dd->bz;
382*8865f1eaSKarl Rupp 
38395c13181SPeter Brune   stencil_type = dd->stencil_type;
384*8865f1eaSKarl Rupp 
385320964c4SBlaise Bourdin   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr);
386320964c4SBlaise Bourdin   if (dim == 1) {
387320964c4SBlaise Bourdin     ierr = DMDACreate1d(((PetscObject)da)->comm,bx,M,nfields,s,dd->lx,nda);CHKERRQ(ierr);
388320964c4SBlaise Bourdin   } else if (dim == 2) {
389320964c4SBlaise Bourdin     ierr = DMDACreate2d(((PetscObject)da)->comm,bx,by,stencil_type,M,N,m,n,nfields,s,lx,ly,nda);CHKERRQ(ierr);
390320964c4SBlaise Bourdin   } else if (dim == 3) {
391320964c4SBlaise 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);
392bc2bf880SBarry Smith   }
3936636e97aSMatthew G Knepley   if (da->coordinates) {
3946636e97aSMatthew G Knepley     ierr = PetscObjectReference((PetscObject)da->coordinates);CHKERRQ(ierr);
395*8865f1eaSKarl Rupp 
3966636e97aSMatthew G Knepley     (*nda)->coordinates = da->coordinates;
397bc2bf880SBarry Smith   }
39895c13181SPeter Brune 
39995c13181SPeter Brune   /* allow for getting a reduced DA corresponding to a domain decomposition */
40095c13181SPeter Brune   ierr = DMDAGetOffset(da,&ox,&oy,&oz,&Mo,&No,&Po);CHKERRQ(ierr);
40195c13181SPeter Brune   ierr = DMDASetOffset(*nda,ox,oy,oz,Mo,No,Po);CHKERRQ(ierr);
40295c13181SPeter Brune 
40395c13181SPeter Brune   /* allow for getting a reduced DA corresponding to a coarsened DA */
40495c13181SPeter Brune   ierr = DMGetCoarsenLevel(da,&cl);CHKERRQ(ierr);
40595c13181SPeter Brune   ierr = DMGetRefineLevel(da,&rl);CHKERRQ(ierr);
406*8865f1eaSKarl Rupp 
40795c13181SPeter Brune   (*nda)->levelup   = rl;
40895c13181SPeter Brune   (*nda)->leveldown = cl;
409bc2bf880SBarry Smith   PetscFunctionReturn(0);
410bc2bf880SBarry Smith }
411bc2bf880SBarry Smith 
412