xref: /petsc/src/dm/impls/da/dacorn.c (revision 9a42bb27a39f0cdf3306a1e22d33cd9809484eaa)
147c6ae99SBarry Smith #define PETSCDM_DLL
247c6ae99SBarry Smith 
347c6ae99SBarry Smith /*
447c6ae99SBarry Smith   Code for manipulating distributed regular arrays in parallel.
547c6ae99SBarry Smith */
647c6ae99SBarry Smith 
747c6ae99SBarry Smith #include "private/daimpl.h"    /*I   "petscda.h"   I*/
847c6ae99SBarry Smith 
947c6ae99SBarry Smith #undef __FUNCT__
1047c6ae99SBarry Smith #define __FUNCT__ "DASetCoordinates"
1147c6ae99SBarry Smith /*@
1247c6ae99SBarry Smith    DASetCoordinates - Sets into the DA a vector that indicates the
1347c6ae99SBarry Smith       coordinates of the local nodes (NOT including ghost nodes).
1447c6ae99SBarry Smith 
1547c6ae99SBarry Smith    Collective on DA
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith    Input Parameter:
1847c6ae99SBarry Smith +  da - the distributed array
1947c6ae99SBarry Smith -  c - coordinate vector
2047c6ae99SBarry Smith 
2147c6ae99SBarry Smith    Note:
2247c6ae99SBarry Smith     The coordinates should NOT include those for all ghost points
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   Level: intermediate
2547c6ae99SBarry Smith 
2647c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith .seealso: DAGetGhostCorners(), DAGetCoordinates(), DASetUniformCoordinates(). DAGetGhostedCoordinates(), DAGetCoordinateDA()
2947c6ae99SBarry Smith @*/
30*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetCoordinates(DM da,Vec c)
3147c6ae99SBarry Smith {
3247c6ae99SBarry Smith   PetscErrorCode ierr;
3347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
3447c6ae99SBarry Smith 
3547c6ae99SBarry Smith   PetscFunctionBegin;
3647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
3747c6ae99SBarry Smith   PetscValidHeaderSpecific(c,VEC_CLASSID,2);
3847c6ae99SBarry Smith   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
3947c6ae99SBarry Smith   if (dd->coordinates) {ierr = VecDestroy(dd->coordinates);CHKERRQ(ierr);}
4047c6ae99SBarry Smith   dd->coordinates = c;
4147c6ae99SBarry Smith   ierr = VecSetBlockSize(c,dd->dim);CHKERRQ(ierr);
4247c6ae99SBarry Smith   if (dd->ghosted_coordinates) { /* The ghosted coordinates are no longer valid */
4347c6ae99SBarry Smith     ierr = VecDestroy(dd->ghosted_coordinates);CHKERRQ(ierr);
4447c6ae99SBarry Smith     dd->ghosted_coordinates = PETSC_NULL;
4547c6ae99SBarry Smith   }
4647c6ae99SBarry Smith   PetscFunctionReturn(0);
4747c6ae99SBarry Smith }
4847c6ae99SBarry Smith 
4947c6ae99SBarry Smith #undef __FUNCT__
5047c6ae99SBarry Smith #define __FUNCT__ "DAGetCoordinates"
5147c6ae99SBarry Smith /*@
5247c6ae99SBarry Smith    DAGetCoordinates - Gets the node coordinates associated with a DA.
5347c6ae99SBarry Smith 
5447c6ae99SBarry Smith    Not Collective
5547c6ae99SBarry Smith 
5647c6ae99SBarry Smith    Input Parameter:
5747c6ae99SBarry Smith .  da - the distributed array
5847c6ae99SBarry Smith 
5947c6ae99SBarry Smith    Output Parameter:
6047c6ae99SBarry Smith .  c - coordinate vector
6147c6ae99SBarry Smith 
6247c6ae99SBarry Smith    Note:
6347c6ae99SBarry Smith     Each process has only the coordinates for its local nodes (does NOT have the
6447c6ae99SBarry Smith   coordinates for the ghost nodes).
6547c6ae99SBarry Smith 
6647c6ae99SBarry Smith     For two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
6747c6ae99SBarry Smith     and (x_0,y_0,z_0,x_1,y_1,z_1...)
6847c6ae99SBarry Smith 
6947c6ae99SBarry Smith   Level: intermediate
7047c6ae99SBarry Smith 
7147c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
7247c6ae99SBarry Smith 
7347c6ae99SBarry Smith .seealso: DAGetGhostCorners(), DASetCoordinates(), DASetUniformCoordinates(), DAGetGhostedCoordinates(), DAGetCoordinateDA()
7447c6ae99SBarry Smith @*/
75*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetCoordinates(DM da,Vec *c)
7647c6ae99SBarry Smith {
7747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
7847c6ae99SBarry Smith   PetscFunctionBegin;
7947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
8047c6ae99SBarry Smith   PetscValidPointer(c,2);
8147c6ae99SBarry Smith   *c = dd->coordinates;
8247c6ae99SBarry Smith   PetscFunctionReturn(0);
8347c6ae99SBarry Smith }
8447c6ae99SBarry Smith 
8547c6ae99SBarry Smith #undef __FUNCT__
8647c6ae99SBarry Smith #define __FUNCT__ "DAGetCoordinateDA"
8747c6ae99SBarry Smith /*@
8847c6ae99SBarry Smith    DAGetCoordinateDA - Gets the DA that scatters between global and local DA coordinates
8947c6ae99SBarry Smith 
9047c6ae99SBarry Smith    Collective on DA
9147c6ae99SBarry Smith 
9247c6ae99SBarry Smith    Input Parameter:
9347c6ae99SBarry Smith .  da - the distributed array
9447c6ae99SBarry Smith 
9547c6ae99SBarry Smith    Output Parameter:
9647c6ae99SBarry Smith .  dac - coordinate DA
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith   Level: intermediate
9947c6ae99SBarry Smith 
10047c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith .seealso: DAGetGhostCorners(), DASetCoordinates(), DASetUniformCoordinates(), DAGetCoordinates(), DAGetGhostedCoordinates()
10347c6ae99SBarry Smith @*/
104*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetCoordinateDA(DM da,DM *cda)
10547c6ae99SBarry Smith {
10647c6ae99SBarry Smith   PetscMPIInt    size;
10747c6ae99SBarry Smith   PetscErrorCode ierr;
10847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
10947c6ae99SBarry Smith 
11047c6ae99SBarry Smith   PetscFunctionBegin;
11147c6ae99SBarry Smith   if (!dd->da_coordinates) {
11247c6ae99SBarry Smith     ierr = MPI_Comm_size(((PetscObject)da)->comm,&size);CHKERRQ(ierr);
11347c6ae99SBarry Smith     if (dd->dim == 1) {
11447c6ae99SBarry Smith       PetscInt            s,m,*lc,l;
11547c6ae99SBarry Smith       DAPeriodicType pt;
11647c6ae99SBarry Smith       ierr = DAGetInfo(da,0,&m,0,0,0,0,0,0,&s,&pt,0);CHKERRQ(ierr);
11747c6ae99SBarry Smith       ierr = DAGetCorners(da,0,0,0,&l,0,0);CHKERRQ(ierr);
11847c6ae99SBarry Smith       ierr = PetscMalloc(size*sizeof(PetscInt),&lc);CHKERRQ(ierr);
11947c6ae99SBarry Smith       ierr = MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
12047c6ae99SBarry Smith       ierr = DACreate1d(((PetscObject)da)->comm,pt,m,1,s,lc,&dd->da_coordinates);CHKERRQ(ierr);
12147c6ae99SBarry Smith       ierr = PetscFree(lc);CHKERRQ(ierr);
12247c6ae99SBarry Smith     } else if (dd->dim == 2) {
12347c6ae99SBarry Smith       PetscInt            i,s,m,*lc,*ld,l,k,n,M,N;
12447c6ae99SBarry Smith       DAPeriodicType pt;
12547c6ae99SBarry Smith       ierr = DAGetInfo(da,0,&m,&n,0,&M,&N,0,0,&s,&pt,0);CHKERRQ(ierr);
12647c6ae99SBarry Smith       ierr = DAGetCorners(da,0,0,0,&l,&k,0);CHKERRQ(ierr);
12747c6ae99SBarry Smith       ierr = PetscMalloc2(size,PetscInt,&lc,size,PetscInt,&ld);CHKERRQ(ierr);
12847c6ae99SBarry Smith       /* only first M values in lc matter */
12947c6ae99SBarry Smith       ierr = MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
13047c6ae99SBarry Smith       /* every Mth value in ld matters */
13147c6ae99SBarry Smith       ierr = MPI_Allgather(&k,1,MPIU_INT,ld,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
13247c6ae99SBarry Smith       for ( i=0; i<N; i++) {
13347c6ae99SBarry Smith         ld[i] = ld[M*i];
13447c6ae99SBarry Smith       }
13547c6ae99SBarry Smith       ierr = DACreate2d(((PetscObject)da)->comm,pt,DA_STENCIL_BOX,m,n,M,N,2,s,lc,ld,&dd->da_coordinates);CHKERRQ(ierr);
13647c6ae99SBarry Smith       ierr = PetscFree2(lc,ld);CHKERRQ(ierr);
13747c6ae99SBarry Smith     } else if (dd->dim == 3) {
13847c6ae99SBarry Smith       PetscInt            i,s,m,*lc,*ld,*le,l,k,q,n,M,N,P,p;
13947c6ae99SBarry Smith       DAPeriodicType pt;
14047c6ae99SBarry Smith       ierr = DAGetInfo(da,0,&m,&n,&p,&M,&N,&P,0,&s,&pt,0);CHKERRQ(ierr);
14147c6ae99SBarry Smith       ierr = DAGetCorners(da,0,0,0,&l,&k,&q);CHKERRQ(ierr);
14247c6ae99SBarry Smith       ierr = PetscMalloc3(size,PetscInt,&lc,size,PetscInt,&ld,size,PetscInt,&le);CHKERRQ(ierr);
14347c6ae99SBarry Smith       /* only first M values in lc matter */
14447c6ae99SBarry Smith       ierr = MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
14547c6ae99SBarry Smith       /* every Mth value in ld matters */
14647c6ae99SBarry Smith       ierr = MPI_Allgather(&k,1,MPIU_INT,ld,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
14747c6ae99SBarry Smith       for ( i=0; i<N; i++) {
14847c6ae99SBarry Smith         ld[i] = ld[M*i];
14947c6ae99SBarry Smith       }
15047c6ae99SBarry Smith       ierr = MPI_Allgather(&q,1,MPIU_INT,le,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr);
15147c6ae99SBarry Smith       for ( i=0; i<P; i++) {
15247c6ae99SBarry Smith         le[i] = le[M*N*i];
15347c6ae99SBarry Smith       }
15447c6ae99SBarry Smith       ierr = DACreate3d(((PetscObject)da)->comm,pt,DA_STENCIL_BOX,m,n,p,M,N,P,3,s,lc,ld,le,&dd->da_coordinates);CHKERRQ(ierr);
15547c6ae99SBarry Smith       ierr = PetscFree3(lc,ld,le);CHKERRQ(ierr);
15647c6ae99SBarry Smith     }
15747c6ae99SBarry Smith   }
15847c6ae99SBarry Smith   *cda = dd->da_coordinates;
15947c6ae99SBarry Smith   PetscFunctionReturn(0);
16047c6ae99SBarry Smith }
16147c6ae99SBarry Smith 
16247c6ae99SBarry Smith 
16347c6ae99SBarry Smith #undef __FUNCT__
16447c6ae99SBarry Smith #define __FUNCT__ "DAGetGhostedCoordinates"
16547c6ae99SBarry Smith /*@
16647c6ae99SBarry Smith    DAGetGhostedCoordinates - Gets the node coordinates associated with a DA.
16747c6ae99SBarry Smith 
16847c6ae99SBarry Smith    Collective on DA
16947c6ae99SBarry Smith 
17047c6ae99SBarry Smith    Input Parameter:
17147c6ae99SBarry Smith .  da - the distributed array
17247c6ae99SBarry Smith 
17347c6ae99SBarry Smith    Output Parameter:
17447c6ae99SBarry Smith .  c - coordinate vector
17547c6ae99SBarry Smith 
17647c6ae99SBarry Smith    Note:
17747c6ae99SBarry Smith     Each process has only the coordinates for its local AND ghost nodes
17847c6ae99SBarry Smith 
17947c6ae99SBarry Smith     For two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
18047c6ae99SBarry Smith     and (x_0,y_0,z_0,x_1,y_1,z_1...)
18147c6ae99SBarry Smith 
18247c6ae99SBarry Smith   Level: intermediate
18347c6ae99SBarry Smith 
18447c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
18547c6ae99SBarry Smith 
18647c6ae99SBarry Smith .seealso: DAGetGhostCorners(), DASetCoordinates(), DASetUniformCoordinates(), DAGetCoordinates(), DAGetCoordinateDA()
18747c6ae99SBarry Smith @*/
188*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetGhostedCoordinates(DM da,Vec *c)
18947c6ae99SBarry Smith {
19047c6ae99SBarry Smith   PetscErrorCode ierr;
19147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
19247c6ae99SBarry Smith 
19347c6ae99SBarry Smith   PetscFunctionBegin;
19447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
19547c6ae99SBarry Smith   PetscValidPointer(c,2);
19647c6ae99SBarry Smith   if (!dd->coordinates) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call DASetCoordinates() before this call");
19747c6ae99SBarry Smith   if (!dd->ghosted_coordinates) {
198*9a42bb27SBarry Smith     DM dac;
19947c6ae99SBarry Smith     ierr = DAGetCoordinateDA(da,&dac);CHKERRQ(ierr);
20047c6ae99SBarry Smith     ierr = DACreateLocalVector(dac,&dd->ghosted_coordinates);CHKERRQ(ierr);
201*9a42bb27SBarry Smith     ierr = DMGlobalToLocalBegin(dac,dd->coordinates,INSERT_VALUES,dd->ghosted_coordinates);CHKERRQ(ierr);
202*9a42bb27SBarry Smith     ierr = DMGlobalToLocalEnd(dac,dd->coordinates,INSERT_VALUES,dd->ghosted_coordinates);CHKERRQ(ierr);
20347c6ae99SBarry Smith   }
20447c6ae99SBarry Smith   *c = dd->ghosted_coordinates;
20547c6ae99SBarry Smith   PetscFunctionReturn(0);
20647c6ae99SBarry Smith }
20747c6ae99SBarry Smith 
20847c6ae99SBarry Smith #undef __FUNCT__
20947c6ae99SBarry Smith #define __FUNCT__ "DASetFieldName"
21047c6ae99SBarry Smith /*@C
21147c6ae99SBarry Smith    DASetFieldName - Sets the names of individual field components in multicomponent
21247c6ae99SBarry Smith    vectors associated with a DA.
21347c6ae99SBarry Smith 
21447c6ae99SBarry Smith    Not Collective
21547c6ae99SBarry Smith 
21647c6ae99SBarry Smith    Input Parameters:
21747c6ae99SBarry Smith +  da - the distributed array
21847c6ae99SBarry Smith .  nf - field number for the DA (0, 1, ... dof-1), where dof indicates the
21947c6ae99SBarry Smith         number of degrees of freedom per node within the DA
22047c6ae99SBarry Smith -  names - the name of the field (component)
22147c6ae99SBarry Smith 
22247c6ae99SBarry Smith   Level: intermediate
22347c6ae99SBarry Smith 
22447c6ae99SBarry Smith .keywords: distributed array, get, component name
22547c6ae99SBarry Smith 
22647c6ae99SBarry Smith .seealso: DAGetFieldName()
22747c6ae99SBarry Smith @*/
228*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetFieldName(DM da,PetscInt nf,const char name[])
22947c6ae99SBarry Smith {
23047c6ae99SBarry Smith   PetscErrorCode ierr;
23147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
23247c6ae99SBarry Smith 
23347c6ae99SBarry Smith   PetscFunctionBegin;
23447c6ae99SBarry Smith    PetscValidHeaderSpecific(da,DM_CLASSID,1);
23547c6ae99SBarry Smith   if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf);
23647c6ae99SBarry Smith   if (dd->fieldname[nf]) {ierr = PetscFree(dd->fieldname[nf]);CHKERRQ(ierr);}
23747c6ae99SBarry Smith    ierr = PetscStrallocpy(name,&dd->fieldname[nf]);CHKERRQ(ierr);
23847c6ae99SBarry Smith   PetscFunctionReturn(0);
23947c6ae99SBarry Smith }
24047c6ae99SBarry Smith 
24147c6ae99SBarry Smith #undef __FUNCT__
24247c6ae99SBarry Smith #define __FUNCT__ "DAGetFieldName"
24347c6ae99SBarry Smith /*@C
24447c6ae99SBarry Smith    DAGetFieldName - Gets the names of individual field components in multicomponent
24547c6ae99SBarry Smith    vectors associated with a DA.
24647c6ae99SBarry Smith 
24747c6ae99SBarry Smith    Not Collective
24847c6ae99SBarry Smith 
24947c6ae99SBarry Smith    Input Parameter:
25047c6ae99SBarry Smith +  da - the distributed array
25147c6ae99SBarry Smith -  nf - field number for the DA (0, 1, ... dof-1), where dof indicates the
25247c6ae99SBarry Smith         number of degrees of freedom per node within the DA
25347c6ae99SBarry Smith 
25447c6ae99SBarry Smith    Output Parameter:
25547c6ae99SBarry Smith .  names - the name of the field (component)
25647c6ae99SBarry Smith 
25747c6ae99SBarry Smith   Level: intermediate
25847c6ae99SBarry Smith 
25947c6ae99SBarry Smith .keywords: distributed array, get, component name
26047c6ae99SBarry Smith 
26147c6ae99SBarry Smith .seealso: DASetFieldName()
26247c6ae99SBarry Smith @*/
263*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetFieldName(DM da,PetscInt nf,const char **name)
26447c6ae99SBarry Smith {
26547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
26647c6ae99SBarry Smith 
26747c6ae99SBarry Smith   PetscFunctionBegin;
26847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
26947c6ae99SBarry Smith   PetscValidPointer(name,3);
27047c6ae99SBarry Smith   if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf);
27147c6ae99SBarry Smith   *name = dd->fieldname[nf];
27247c6ae99SBarry Smith   PetscFunctionReturn(0);
27347c6ae99SBarry Smith }
27447c6ae99SBarry Smith 
27547c6ae99SBarry Smith #undef __FUNCT__
27647c6ae99SBarry Smith #define __FUNCT__ "DAGetCorners"
27747c6ae99SBarry Smith /*@
27847c6ae99SBarry Smith    DAGetCorners - Returns the global (x,y,z) indices of the lower left
27947c6ae99SBarry Smith    corner of the local region, excluding ghost points.
28047c6ae99SBarry Smith 
28147c6ae99SBarry Smith    Not Collective
28247c6ae99SBarry Smith 
28347c6ae99SBarry Smith    Input Parameter:
28447c6ae99SBarry Smith .  da - the distributed array
28547c6ae99SBarry Smith 
28647c6ae99SBarry Smith    Output Parameters:
28747c6ae99SBarry Smith +  x,y,z - the corner indices (where y and z are optional; these are used
28847c6ae99SBarry Smith            for 2D and 3D problems)
28947c6ae99SBarry Smith -  m,n,p - widths in the corresponding directions (where n and p are optional;
29047c6ae99SBarry Smith            these are used for 2D and 3D problems)
29147c6ae99SBarry Smith 
29247c6ae99SBarry Smith    Note:
29347c6ae99SBarry Smith    The corner information is independent of the number of degrees of
29447c6ae99SBarry Smith    freedom per node set with the DACreateXX() routine. Thus the x, y, z, and
29547c6ae99SBarry Smith    m, n, p can be thought of as coordinates on a logical grid, where each
29647c6ae99SBarry Smith    grid point has (potentially) several degrees of freedom.
29747c6ae99SBarry Smith    Any of y, z, n, and p can be passed in as PETSC_NULL if not needed.
29847c6ae99SBarry Smith 
29947c6ae99SBarry Smith   Level: beginner
30047c6ae99SBarry Smith 
30147c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices
30247c6ae99SBarry Smith 
30347c6ae99SBarry Smith .seealso: DAGetGhostCorners(), DAGetOwnershipRanges()
30447c6ae99SBarry Smith @*/
305*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p)
30647c6ae99SBarry Smith {
30747c6ae99SBarry Smith   PetscInt w;
30847c6ae99SBarry Smith   DM_DA    *dd = (DM_DA*)da->data;
30947c6ae99SBarry Smith 
31047c6ae99SBarry Smith   PetscFunctionBegin;
31147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
31247c6ae99SBarry Smith   /* since the xs, xe ... have all been multiplied by the number of degrees
31347c6ae99SBarry Smith      of freedom per cell, w = dd->w, we divide that out before returning.*/
31447c6ae99SBarry Smith   w = dd->w;
31547c6ae99SBarry Smith   if (x) *x = dd->xs/w; if(m) *m = (dd->xe - dd->xs)/w;
31647c6ae99SBarry Smith   /* the y and z have NOT been multiplied by w */
31747c6ae99SBarry Smith   if (y) *y = dd->ys;   if (n) *n = (dd->ye - dd->ys);
31847c6ae99SBarry Smith   if (z) *z = dd->zs;   if (p) *p = (dd->ze - dd->zs);
31947c6ae99SBarry Smith   PetscFunctionReturn(0);
32047c6ae99SBarry Smith }
32147c6ae99SBarry Smith 
32247c6ae99SBarry Smith #undef __FUNCT__
32347c6ae99SBarry Smith #define __FUNCT__ "DAGetLocalBoundingBox"
32447c6ae99SBarry Smith /*@
32547c6ae99SBarry Smith    DAGetLocalBoundingBox - Returns the local bounding box for the DA.
32647c6ae99SBarry Smith 
32747c6ae99SBarry Smith    Not Collective
32847c6ae99SBarry Smith 
32947c6ae99SBarry Smith    Input Parameter:
33047c6ae99SBarry Smith .  da - the distributed array
33147c6ae99SBarry Smith 
33247c6ae99SBarry Smith    Output Parameters:
33347c6ae99SBarry Smith +  lmin - local minimum coordinates (length dim, optional)
33447c6ae99SBarry Smith -  lmax - local maximim coordinates (length dim, optional)
33547c6ae99SBarry Smith 
33647c6ae99SBarry Smith   Level: beginner
33747c6ae99SBarry Smith 
33847c6ae99SBarry Smith .keywords: distributed array, get, coordinates
33947c6ae99SBarry Smith 
34047c6ae99SBarry Smith .seealso: DAGetCoordinateDA(), DAGetCoordinates(), DAGetBoundingBox()
34147c6ae99SBarry Smith @*/
342*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetLocalBoundingBox(DM da,PetscReal lmin[],PetscReal lmax[])
34347c6ae99SBarry Smith {
34447c6ae99SBarry Smith   PetscErrorCode    ierr;
34547c6ae99SBarry Smith   Vec               coords  = PETSC_NULL;
34647c6ae99SBarry Smith   PetscInt          dim,i,j;
34747c6ae99SBarry Smith   const PetscScalar *local_coords;
34847c6ae99SBarry Smith   PetscReal         min[3]={PETSC_MAX,PETSC_MAX,PETSC_MAX},max[3]={PETSC_MIN,PETSC_MIN,PETSC_MIN};
34947c6ae99SBarry Smith   PetscInt          N,Ni;
35047c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
35147c6ae99SBarry Smith 
35247c6ae99SBarry Smith   PetscFunctionBegin;
35347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
35447c6ae99SBarry Smith   dim = dd->dim;
35547c6ae99SBarry Smith   ierr = DAGetCoordinates(da,&coords);CHKERRQ(ierr);
35647c6ae99SBarry Smith   ierr = VecGetArrayRead(coords,&local_coords);CHKERRQ(ierr);
35747c6ae99SBarry Smith   ierr = VecGetLocalSize(coords,&N);CHKERRQ(ierr);
35847c6ae99SBarry Smith   Ni = N/dim;
35947c6ae99SBarry Smith   for (i=0; i<Ni; i++) {
36047c6ae99SBarry Smith     for (j=0; j<dim; j++) {
36147c6ae99SBarry Smith       min[j] = PetscMin(min[j],PetscRealPart(local_coords[i*dim+j]));CHKERRQ(ierr);
36247c6ae99SBarry Smith       max[j] = PetscMax(min[j],PetscRealPart(local_coords[i*dim+j]));CHKERRQ(ierr);
36347c6ae99SBarry Smith     }
36447c6ae99SBarry Smith   }
36547c6ae99SBarry Smith   ierr = VecRestoreArrayRead(coords,&local_coords);CHKERRQ(ierr);
36647c6ae99SBarry Smith   if (lmin) {ierr = PetscMemcpy(lmin,min,dim*sizeof(PetscReal));CHKERRQ(ierr);}
36747c6ae99SBarry Smith   if (lmax) {ierr = PetscMemcpy(lmax,max,dim*sizeof(PetscReal));CHKERRQ(ierr);}
36847c6ae99SBarry Smith   PetscFunctionReturn(0);
36947c6ae99SBarry Smith }
37047c6ae99SBarry Smith 
37147c6ae99SBarry Smith #undef __FUNCT__
37247c6ae99SBarry Smith #define __FUNCT__ "DAGetBoundingBox"
37347c6ae99SBarry Smith /*@
37447c6ae99SBarry Smith    DAGetBoundingBox - Returns the global bounding box for the DA.
37547c6ae99SBarry Smith 
37647c6ae99SBarry Smith    Collective on DA
37747c6ae99SBarry Smith 
37847c6ae99SBarry Smith    Input Parameter:
37947c6ae99SBarry Smith .  da - the distributed array
38047c6ae99SBarry Smith 
38147c6ae99SBarry Smith    Output Parameters:
38247c6ae99SBarry Smith +  gmin - global minimum coordinates (length dim, optional)
38347c6ae99SBarry Smith -  gmax - global maximim coordinates (length dim, optional)
38447c6ae99SBarry Smith 
38547c6ae99SBarry Smith   Level: beginner
38647c6ae99SBarry Smith 
38747c6ae99SBarry Smith .keywords: distributed array, get, coordinates
38847c6ae99SBarry Smith 
38947c6ae99SBarry Smith .seealso: DAGetCoordinateDA(), DAGetCoordinates(), DAGetLocalBoundingBox()
39047c6ae99SBarry Smith @*/
391*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetBoundingBox(DM da,PetscReal gmin[],PetscReal gmax[])
39247c6ae99SBarry Smith {
39347c6ae99SBarry Smith   PetscErrorCode ierr;
39447c6ae99SBarry Smith   PetscMPIInt    count;
39547c6ae99SBarry Smith   PetscReal      lmin[3],lmax[3];
39647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
39747c6ae99SBarry Smith 
39847c6ae99SBarry Smith   PetscFunctionBegin;
39947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
40047c6ae99SBarry Smith   count = PetscMPIIntCast(dd->dim);
40147c6ae99SBarry Smith   ierr = DAGetLocalBoundingBox(da,lmin,lmax);CHKERRQ(ierr);
40247c6ae99SBarry Smith   if (gmin) {ierr = MPI_Allreduce(lmin,gmin,count,MPIU_REAL,MPI_MIN,((PetscObject)da)->comm);CHKERRQ(ierr);}
40347c6ae99SBarry Smith   if (gmax) {ierr = MPI_Allreduce(lmax,gmax,count,MPIU_REAL,MPI_MAX,((PetscObject)da)->comm);CHKERRQ(ierr);}
40447c6ae99SBarry Smith   PetscFunctionReturn(0);
40547c6ae99SBarry Smith }
406