xref: /petsc/src/dm/impls/da/dasub.c (revision 3ee9839e4674d2e0e9bcea975330f3458ebcb61e)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
347c6ae99SBarry Smith   Code for manipulating distributed regular arrays in parallel.
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
747c6ae99SBarry Smith 
869e49704SSatish Balay /*@
90a0169f3SBarry Smith    DMDAGetLogicalCoordinate - Returns a the i,j,k logical coordinate for the closest mesh point to a x,y,z point in the coordinates of the DMDA
100a0169f3SBarry Smith 
11d083f849SBarry Smith    Collective on da
120a0169f3SBarry Smith 
130a0169f3SBarry Smith    Input Parameters:
140a0169f3SBarry Smith +  da - the distributed array
15878cb397SSatish Balay -  x,y,z - the physical coordinates
160a0169f3SBarry Smith 
170a0169f3SBarry Smith    Output Parameters:
1833c53e3fSSatish Balay +   II, JJ, KK - the logical coordinate (-1 on processes that do not contain that point)
190a0169f3SBarry Smith -   X, Y, Z, - (optional) the coordinates of the located grid point
200a0169f3SBarry Smith 
210a0169f3SBarry Smith    Level: advanced
220a0169f3SBarry Smith 
230a0169f3SBarry Smith    Notes:
240a0169f3SBarry Smith    All processors that share the DMDA must call this with the same coordinate value
250a0169f3SBarry Smith 
260a0169f3SBarry Smith @*/
2733c53e3fSSatish Balay PetscErrorCode  DMDAGetLogicalCoordinate(DM da,PetscScalar x,PetscScalar y,PetscScalar z,PetscInt *II,PetscInt *JJ,PetscInt *KK,PetscScalar *X,PetscScalar *Y,PetscScalar *Z)
280a0169f3SBarry Smith {
290a0169f3SBarry Smith   PetscErrorCode ierr;
300a0169f3SBarry Smith   Vec            coors;
310a0169f3SBarry Smith   DM             dacoors;
320a0169f3SBarry Smith   DMDACoor2d     **c;
330a0169f3SBarry Smith   PetscInt       i,j,xs,xm,ys,ym;
34e97f6855SBarry Smith   PetscReal      d,D = PETSC_MAX_REAL,Dv;
35e97f6855SBarry Smith   PetscMPIInt    rank,root;
360a0169f3SBarry Smith 
370a0169f3SBarry Smith   PetscFunctionBegin;
38c73cfb54SMatthew G. Knepley   if (da->dim == 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot get point from 1d DMDA");
39c73cfb54SMatthew G. Knepley   if (da->dim == 3) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot get point from 3d DMDA");
400a0169f3SBarry Smith 
4133c53e3fSSatish Balay   *II = -1;
42e97f6855SBarry Smith   *JJ = -1;
430a0169f3SBarry Smith 
440a0169f3SBarry Smith   ierr = DMGetCoordinateDM(da,&dacoors);CHKERRQ(ierr);
450a0169f3SBarry Smith   ierr = DMDAGetCorners(dacoors,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
460a0169f3SBarry Smith   ierr = DMGetCoordinates(da,&coors);CHKERRQ(ierr);
475edff71fSBarry Smith   ierr = DMDAVecGetArrayRead(dacoors,coors,&c);CHKERRQ(ierr);
480a0169f3SBarry Smith   for (j=ys; j<ys+ym; j++) {
490a0169f3SBarry Smith     for (i=xs; i<xs+xm; i++) {
500a0169f3SBarry Smith       d = PetscSqrtReal(PetscRealPart( (c[j][i].x - x)*(c[j][i].x - x) + (c[j][i].y - y)*(c[j][i].y - y) ));
510a0169f3SBarry Smith       if (d < D) {
520a0169f3SBarry Smith         D   = d;
5333c53e3fSSatish Balay         *II = i;
5433c53e3fSSatish Balay         *JJ = j;
550a0169f3SBarry Smith       }
560a0169f3SBarry Smith     }
570a0169f3SBarry Smith   }
58b2566f29SBarry Smith   ierr = MPIU_Allreduce(&D,&Dv,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
59e97f6855SBarry Smith   if (D != Dv) {
60e97f6855SBarry Smith     *II  = -1;
61e97f6855SBarry Smith     *JJ  = -1;
62e97f6855SBarry Smith     rank = 0;
63e97f6855SBarry Smith   } else {
64e97f6855SBarry Smith     *X = c[*JJ][*II].x;
65e97f6855SBarry Smith     *Y = c[*JJ][*II].y;
66e97f6855SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr);
67e97f6855SBarry Smith     rank++;
68e97f6855SBarry Smith   }
69b2566f29SBarry Smith   ierr = MPIU_Allreduce(&rank,&root,1,MPI_INT,MPI_SUM,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
70e97f6855SBarry Smith   root--;
71ae95e784SBarry Smith   ierr = MPI_Bcast(X,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
72ae95e784SBarry Smith   ierr = MPI_Bcast(Y,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
735edff71fSBarry Smith   ierr = DMDAVecRestoreArrayRead(dacoors,coors,&c);CHKERRQ(ierr);
740a0169f3SBarry Smith   PetscFunctionReturn(0);
750a0169f3SBarry Smith }
760a0169f3SBarry Smith 
7769e49704SSatish Balay /*@
78db05f41bSBarry Smith    DMDAGetRay - Returns a vector on process zero that contains a row or column of the values in a DMDA vector
79db05f41bSBarry Smith 
80db05f41bSBarry Smith    Collective on DMDA
81db05f41bSBarry Smith 
82db05f41bSBarry Smith    Input Parameters:
83db05f41bSBarry Smith +  da - the distributed array
84db05f41bSBarry Smith .  vec - the vector
85*3ee9839eSMatthew G. Knepley .  dir - Cartesian direction, either DM_X, DM_Y, or DM_Z
86db05f41bSBarry Smith -  gp - global grid point number in this direction
87db05f41bSBarry Smith 
88db05f41bSBarry Smith    Output Parameters:
89db05f41bSBarry Smith +  newvec - the new vector that can hold the values (size zero on all processes except process 0)
90db05f41bSBarry Smith -  scatter - the VecScatter that will map from the original vector to the slice
91db05f41bSBarry Smith 
92db05f41bSBarry Smith    Level: advanced
93db05f41bSBarry Smith 
94db05f41bSBarry Smith    Notes:
95db05f41bSBarry Smith    All processors that share the DMDA must call this with the same gp value
96db05f41bSBarry Smith 
97db05f41bSBarry Smith @*/
98*3ee9839eSMatthew G. Knepley PetscErrorCode  DMDAGetRay(DM da,DMDirection dir,PetscInt gp,Vec *newvec,VecScatter *scatter)
99db05f41bSBarry Smith {
100db05f41bSBarry Smith   PetscMPIInt    rank;
101db05f41bSBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
102db05f41bSBarry Smith   PetscErrorCode ierr;
103db05f41bSBarry Smith   IS             is;
104db05f41bSBarry Smith   AO             ao;
105db05f41bSBarry Smith   Vec            vec;
106d1212d36SBarry Smith   PetscInt       *indices,i,j;
107db05f41bSBarry Smith 
108db05f41bSBarry Smith   PetscFunctionBegin;
109c73cfb54SMatthew G. Knepley   if (da->dim == 3) SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_SUP, "Cannot get slice from 3d DMDA");
110ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) da), &rank);CHKERRQ(ierr);
111c2001ae3SMatthew G. Knepley   ierr = DMDAGetAO(da, &ao);CHKERRQ(ierr);
112db05f41bSBarry Smith   if (!rank) {
113c73cfb54SMatthew G. Knepley     if (da->dim == 1) {
114*3ee9839eSMatthew G. Knepley       if (dir == DM_X) {
115785e854fSJed Brown         ierr = PetscMalloc1(dd->w, &indices);CHKERRQ(ierr);
116c2001ae3SMatthew G. Knepley         indices[0] = dd->w*gp;
117eca19243SMatthew G. Knepley         for (i = 1; i < dd->w; ++i) indices[i] = indices[i-1] + 1;
118c2001ae3SMatthew G. Knepley         ierr = AOApplicationToPetsc(ao, dd->w, indices);CHKERRQ(ierr);
119c2001ae3SMatthew G. Knepley         ierr = VecCreate(PETSC_COMM_SELF, newvec);CHKERRQ(ierr);
120c2001ae3SMatthew G. Knepley         ierr = VecSetBlockSize(*newvec, dd->w);CHKERRQ(ierr);
121c2001ae3SMatthew G. Knepley         ierr = VecSetSizes(*newvec, dd->w, PETSC_DETERMINE);CHKERRQ(ierr);
122c2001ae3SMatthew G. Knepley         ierr = VecSetType(*newvec, VECSEQ);CHKERRQ(ierr);
123c2001ae3SMatthew G. Knepley         ierr = ISCreateGeneral(PETSC_COMM_SELF, dd->w, indices, PETSC_OWN_POINTER, &is);CHKERRQ(ierr);
124*3ee9839eSMatthew G. Knepley       } else if (dir == DM_Y) SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_SUP, "Cannot get Y slice from 1d DMDA");
125*3ee9839eSMatthew G. Knepley       else SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_ARG_OUTOFRANGE, "Unknown DMDirection");
126c2001ae3SMatthew G. Knepley     } else {
127*3ee9839eSMatthew G. Knepley       if (dir == DM_Y) {
128785e854fSJed Brown         ierr       = PetscMalloc1(dd->w*dd->M,&indices);CHKERRQ(ierr);
129d1212d36SBarry Smith         indices[0] = gp*dd->M*dd->w;
1308865f1eaSKarl Rupp         for (i=1; i<dd->M*dd->w; i++) indices[i] = indices[i-1] + 1;
1318865f1eaSKarl Rupp 
132d1212d36SBarry Smith         ierr = AOApplicationToPetsc(ao,dd->M*dd->w,indices);CHKERRQ(ierr);
133db05f41bSBarry Smith         ierr = VecCreate(PETSC_COMM_SELF,newvec);CHKERRQ(ierr);
134db05f41bSBarry Smith         ierr = VecSetBlockSize(*newvec,dd->w);CHKERRQ(ierr);
135db05f41bSBarry Smith         ierr = VecSetSizes(*newvec,dd->M*dd->w,PETSC_DETERMINE);CHKERRQ(ierr);
136db05f41bSBarry Smith         ierr = VecSetType(*newvec,VECSEQ);CHKERRQ(ierr);
137d1212d36SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,dd->w*dd->M,indices,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
138*3ee9839eSMatthew G. Knepley       } else if (dir == DM_X) {
139785e854fSJed Brown         ierr       = PetscMalloc1(dd->w*dd->N,&indices);CHKERRQ(ierr);
140d1212d36SBarry Smith         indices[0] = dd->w*gp;
141d1212d36SBarry Smith         for (j=1; j<dd->w; j++) indices[j] = indices[j-1] + 1;
142d1212d36SBarry Smith         for (i=1; i<dd->N; i++) {
143d1212d36SBarry Smith           indices[i*dd->w] = indices[i*dd->w-1] + dd->w*dd->M - dd->w + 1;
144d1212d36SBarry Smith           for (j=1; j<dd->w; j++) indices[i*dd->w + j] = indices[i*dd->w + j - 1] + 1;
145d1212d36SBarry Smith         }
146d1212d36SBarry Smith         ierr = AOApplicationToPetsc(ao,dd->w*dd->N,indices);CHKERRQ(ierr);
147db05f41bSBarry Smith         ierr = VecCreate(PETSC_COMM_SELF,newvec);CHKERRQ(ierr);
148db05f41bSBarry Smith         ierr = VecSetBlockSize(*newvec,dd->w);CHKERRQ(ierr);
149db05f41bSBarry Smith         ierr = VecSetSizes(*newvec,dd->N*dd->w,PETSC_DETERMINE);CHKERRQ(ierr);
150db05f41bSBarry Smith         ierr = VecSetType(*newvec,VECSEQ);CHKERRQ(ierr);
151d1212d36SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,dd->w*dd->N,indices,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
152*3ee9839eSMatthew G. Knepley       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown DMDirection");
153c2001ae3SMatthew G. Knepley     }
154db05f41bSBarry Smith   } else {
155db05f41bSBarry Smith     ierr = VecCreateSeq(PETSC_COMM_SELF, 0, newvec);CHKERRQ(ierr);
156d1212d36SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, 0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
157db05f41bSBarry Smith   }
158db05f41bSBarry Smith   ierr = DMGetGlobalVector(da, &vec);CHKERRQ(ierr);
1599448b7f1SJunchao Zhang   ierr = VecScatterCreate(vec, is, *newvec, NULL, scatter);CHKERRQ(ierr);
160db05f41bSBarry Smith   ierr = DMRestoreGlobalVector(da, &vec);CHKERRQ(ierr);
161db05f41bSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
162db05f41bSBarry Smith   PetscFunctionReturn(0);
163db05f41bSBarry Smith }
164db05f41bSBarry Smith 
16547c6ae99SBarry Smith /*@C
166aa219208SBarry Smith    DMDAGetProcessorSubset - Returns a communicator consisting only of the
167aa219208SBarry Smith    processors in a DMDA that own a particular global x, y, or z grid point
16847c6ae99SBarry Smith    (corresponding to a logical plane in a 3D grid or a line in a 2D grid).
16947c6ae99SBarry Smith 
170d083f849SBarry Smith    Collective on da
17147c6ae99SBarry Smith 
17247c6ae99SBarry Smith    Input Parameters:
17347c6ae99SBarry Smith +  da - the distributed array
174*3ee9839eSMatthew G. Knepley .  dir - Cartesian direction, either DM_X, DM_Y, or DM_Z
17547c6ae99SBarry Smith -  gp - global grid point number in this direction
17647c6ae99SBarry Smith 
17747c6ae99SBarry Smith    Output Parameters:
17847c6ae99SBarry Smith .  comm - new communicator
17947c6ae99SBarry Smith 
18047c6ae99SBarry Smith    Level: advanced
18147c6ae99SBarry Smith 
18247c6ae99SBarry Smith    Notes:
183aa219208SBarry Smith    All processors that share the DMDA must call this with the same gp value
18447c6ae99SBarry Smith 
185bf3faf6aSSatish Balay    After use, comm should be freed with MPI_Comm_free()
186bf3faf6aSSatish Balay 
18747c6ae99SBarry Smith    This routine is particularly useful to compute boundary conditions
18847c6ae99SBarry Smith    or other application-specific calculations that require manipulating
18947c6ae99SBarry Smith    sets of data throughout a logical plane of grid points.
19047c6ae99SBarry Smith 
191f5f57ec0SBarry Smith    Not supported from Fortran
192f5f57ec0SBarry Smith 
19347c6ae99SBarry Smith @*/
194*3ee9839eSMatthew G. Knepley PetscErrorCode  DMDAGetProcessorSubset(DM da,DMDirection dir,PetscInt gp,MPI_Comm *comm)
19547c6ae99SBarry Smith {
19647c6ae99SBarry Smith   MPI_Group      group,subgroup;
19747c6ae99SBarry Smith   PetscErrorCode ierr;
19847c6ae99SBarry Smith   PetscInt       i,ict,flag,*owners,xs,xm,ys,ym,zs,zm;
1990298fd71SBarry Smith   PetscMPIInt    size,*ranks = NULL;
20047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith   PetscFunctionBegin;
203a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
20447c6ae99SBarry Smith   flag = 0;
205aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
206ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr);
207*3ee9839eSMatthew G. Knepley   if (dir == DM_Z) {
208*3ee9839eSMatthew G. Knepley     if (da->dim < 3) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"DM_Z invalid for DMDA dim < 3");
20947c6ae99SBarry Smith     if (gp < 0 || gp > dd->P) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
21047c6ae99SBarry Smith     if (gp >= zs && gp < zs+zm) flag = 1;
211*3ee9839eSMatthew G. Knepley   } else if (dir == DM_Y) {
212*3ee9839eSMatthew G. Knepley     if (da->dim == 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"DM_Y invalid for DMDA dim = 1");
21347c6ae99SBarry Smith     if (gp < 0 || gp > dd->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
21447c6ae99SBarry Smith     if (gp >= ys && gp < ys+ym) flag = 1;
215*3ee9839eSMatthew G. Knepley   } else if (dir == DM_X) {
21647c6ae99SBarry Smith     if (gp < 0 || gp > dd->M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
21747c6ae99SBarry Smith     if (gp >= xs && gp < xs+xm) flag = 1;
218ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction");
21947c6ae99SBarry Smith 
220dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&owners,size,&ranks);CHKERRQ(ierr);
221ce94432eSBarry Smith   ierr = MPI_Allgather(&flag,1,MPIU_INT,owners,1,MPIU_INT,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
22247c6ae99SBarry Smith   ict  = 0;
223c73cfb54SMatthew G. Knepley   ierr = PetscInfo2(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",da->dim,(int)dir);CHKERRQ(ierr);
22447c6ae99SBarry Smith   for (i=0; i<size; i++) {
22547c6ae99SBarry Smith     if (owners[i]) {
22647c6ae99SBarry Smith       ranks[ict] = i; ict++;
22747c6ae99SBarry Smith       ierr       = PetscInfo1(da,"%D ",i);CHKERRQ(ierr);
22847c6ae99SBarry Smith     }
22947c6ae99SBarry Smith   }
23047c6ae99SBarry Smith   ierr = PetscInfo(da,"\n");CHKERRQ(ierr);
231ce94432eSBarry Smith   ierr = MPI_Comm_group(PetscObjectComm((PetscObject)da),&group);CHKERRQ(ierr);
23247c6ae99SBarry Smith   ierr = MPI_Group_incl(group,ict,ranks,&subgroup);CHKERRQ(ierr);
233ce94432eSBarry Smith   ierr = MPI_Comm_create(PetscObjectComm((PetscObject)da),subgroup,comm);CHKERRQ(ierr);
23447c6ae99SBarry Smith   ierr = MPI_Group_free(&subgroup);CHKERRQ(ierr);
23547c6ae99SBarry Smith   ierr = MPI_Group_free(&group);CHKERRQ(ierr);
23647c6ae99SBarry Smith   ierr = PetscFree2(owners,ranks);CHKERRQ(ierr);
23747c6ae99SBarry Smith   PetscFunctionReturn(0);
23847c6ae99SBarry Smith }
23947c6ae99SBarry Smith 
24047c6ae99SBarry Smith /*@C
241aa219208SBarry Smith    DMDAGetProcessorSubsets - Returns communicators consisting only of the
242aa219208SBarry Smith    processors in a DMDA adjacent in a particular dimension,
24347c6ae99SBarry Smith    corresponding to a logical plane in a 3D grid or a line in a 2D grid.
24447c6ae99SBarry Smith 
245d083f849SBarry Smith    Collective on da
24647c6ae99SBarry Smith 
24747c6ae99SBarry Smith    Input Parameters:
24847c6ae99SBarry Smith +  da - the distributed array
249*3ee9839eSMatthew G. Knepley -  dir - Cartesian direction, either DM_X, DM_Y, or DM_Z
25047c6ae99SBarry Smith 
25147c6ae99SBarry Smith    Output Parameters:
25247c6ae99SBarry Smith .  subcomm - new communicator
25347c6ae99SBarry Smith 
25447c6ae99SBarry Smith    Level: advanced
25547c6ae99SBarry Smith 
25647c6ae99SBarry Smith    Notes:
25747c6ae99SBarry Smith    This routine is useful for distributing one-dimensional data in a tensor product grid.
25847c6ae99SBarry Smith 
259bf3faf6aSSatish Balay    After use, comm should be freed with MPI_Comm_free()
260bf3faf6aSSatish Balay 
261f5f57ec0SBarry Smith    Not supported from Fortran
262f5f57ec0SBarry Smith 
26347c6ae99SBarry Smith @*/
264*3ee9839eSMatthew G. Knepley PetscErrorCode  DMDAGetProcessorSubsets(DM da, DMDirection dir, MPI_Comm *subcomm)
26547c6ae99SBarry Smith {
26647c6ae99SBarry Smith   MPI_Comm       comm;
26747c6ae99SBarry Smith   MPI_Group      group, subgroup;
26847c6ae99SBarry Smith   PetscInt       subgroupSize = 0;
26947c6ae99SBarry Smith   PetscInt       *firstPoints;
2700298fd71SBarry Smith   PetscMPIInt    size, *subgroupRanks = NULL;
27147c6ae99SBarry Smith   PetscInt       xs, xm, ys, ym, zs, zm, firstPoint, p;
27247c6ae99SBarry Smith   PetscErrorCode ierr;
27347c6ae99SBarry Smith 
27447c6ae99SBarry Smith   PetscFunctionBegin;
275a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
27682f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
277aa219208SBarry Smith   ierr = DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm);CHKERRQ(ierr);
27847c6ae99SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
279*3ee9839eSMatthew G. Knepley   if (dir == DM_Z) {
280*3ee9839eSMatthew G. Knepley     if (da->dim < 3) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DM_Z invalid for DMDA dim < 3");
28147c6ae99SBarry Smith     firstPoint = zs;
282*3ee9839eSMatthew G. Knepley   } else if (dir == DM_Y) {
283*3ee9839eSMatthew G. Knepley     if (da->dim == 1) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DM_Y invalid for DMDA dim = 1");
28447c6ae99SBarry Smith     firstPoint = ys;
285*3ee9839eSMatthew G. Knepley   } else if (dir == DM_X) {
28647c6ae99SBarry Smith     firstPoint = xs;
28747c6ae99SBarry Smith   } else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction");
28847c6ae99SBarry Smith 
289dcca6d9dSJed Brown   ierr = PetscMalloc2(size, &firstPoints, size, &subgroupRanks);CHKERRQ(ierr);
29047c6ae99SBarry Smith   ierr = MPI_Allgather(&firstPoint, 1, MPIU_INT, firstPoints, 1, MPIU_INT, comm);CHKERRQ(ierr);
291c73cfb54SMatthew G. Knepley   ierr = PetscInfo2(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",da->dim,(int)dir);CHKERRQ(ierr);
29247c6ae99SBarry Smith   for (p = 0; p < size; ++p) {
29347c6ae99SBarry Smith     if (firstPoints[p] == firstPoint) {
29447c6ae99SBarry Smith       subgroupRanks[subgroupSize++] = p;
29547c6ae99SBarry Smith       ierr = PetscInfo1(da, "%D ", p);CHKERRQ(ierr);
29647c6ae99SBarry Smith     }
29747c6ae99SBarry Smith   }
29847c6ae99SBarry Smith   ierr = PetscInfo(da, "\n");CHKERRQ(ierr);
29947c6ae99SBarry Smith   ierr = MPI_Comm_group(comm, &group);CHKERRQ(ierr);
30047c6ae99SBarry Smith   ierr = MPI_Group_incl(group, subgroupSize, subgroupRanks, &subgroup);CHKERRQ(ierr);
30147c6ae99SBarry Smith   ierr = MPI_Comm_create(comm, subgroup, subcomm);CHKERRQ(ierr);
30247c6ae99SBarry Smith   ierr = MPI_Group_free(&subgroup);CHKERRQ(ierr);
30347c6ae99SBarry Smith   ierr = MPI_Group_free(&group);CHKERRQ(ierr);
30447c6ae99SBarry Smith   ierr = PetscFree2(firstPoints, subgroupRanks);CHKERRQ(ierr);
30547c6ae99SBarry Smith   PetscFunctionReturn(0);
30647c6ae99SBarry Smith }
307