xref: /petsc/src/dm/impls/da/dasub.c (revision 69e497046ff11ae3505710514a3758ee19418d4b)
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 
8*69e49704SSatish 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 
110a0169f3SBarry Smith    Collective on DMDA
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 .keywords: distributed array, get, processor subset
270a0169f3SBarry Smith @*/
2833c53e3fSSatish Balay PetscErrorCode  DMDAGetLogicalCoordinate(DM da,PetscScalar x,PetscScalar y,PetscScalar z,PetscInt *II,PetscInt *JJ,PetscInt *KK,PetscScalar *X,PetscScalar *Y,PetscScalar *Z)
290a0169f3SBarry Smith {
300a0169f3SBarry Smith   PetscErrorCode ierr;
310a0169f3SBarry Smith   Vec            coors;
320a0169f3SBarry Smith   DM             dacoors;
330a0169f3SBarry Smith   DMDACoor2d     **c;
340a0169f3SBarry Smith   PetscInt       i,j,xs,xm,ys,ym;
35e97f6855SBarry Smith   PetscReal      d,D = PETSC_MAX_REAL,Dv;
36e97f6855SBarry Smith   PetscMPIInt    rank,root;
370a0169f3SBarry Smith 
380a0169f3SBarry Smith   PetscFunctionBegin;
39c73cfb54SMatthew G. Knepley   if (da->dim == 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot get point from 1d DMDA");
40c73cfb54SMatthew G. Knepley   if (da->dim == 3) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot get point from 3d DMDA");
410a0169f3SBarry Smith 
4233c53e3fSSatish Balay   *II = -1;
43e97f6855SBarry Smith   *JJ = -1;
440a0169f3SBarry Smith 
450a0169f3SBarry Smith   ierr = DMGetCoordinateDM(da,&dacoors);CHKERRQ(ierr);
460a0169f3SBarry Smith   ierr = DMDAGetCorners(dacoors,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
470a0169f3SBarry Smith   ierr = DMGetCoordinates(da,&coors);CHKERRQ(ierr);
485edff71fSBarry Smith   ierr = DMDAVecGetArrayRead(dacoors,coors,&c);CHKERRQ(ierr);
490a0169f3SBarry Smith   for (j=ys; j<ys+ym; j++) {
500a0169f3SBarry Smith     for (i=xs; i<xs+xm; i++) {
510a0169f3SBarry Smith       d = PetscSqrtReal(PetscRealPart( (c[j][i].x - x)*(c[j][i].x - x) + (c[j][i].y - y)*(c[j][i].y - y) ));
520a0169f3SBarry Smith       if (d < D) {
530a0169f3SBarry Smith         D   = d;
5433c53e3fSSatish Balay         *II = i;
5533c53e3fSSatish Balay         *JJ = j;
560a0169f3SBarry Smith       }
570a0169f3SBarry Smith     }
580a0169f3SBarry Smith   }
59b2566f29SBarry Smith   ierr = MPIU_Allreduce(&D,&Dv,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
60e97f6855SBarry Smith   if (D != Dv) {
61e97f6855SBarry Smith     *II  = -1;
62e97f6855SBarry Smith     *JJ  = -1;
63e97f6855SBarry Smith     rank = 0;
64e97f6855SBarry Smith   } else {
65e97f6855SBarry Smith     *X = c[*JJ][*II].x;
66e97f6855SBarry Smith     *Y = c[*JJ][*II].y;
67e97f6855SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr);
68e97f6855SBarry Smith     rank++;
69e97f6855SBarry Smith   }
70b2566f29SBarry Smith   ierr = MPIU_Allreduce(&rank,&root,1,MPI_INT,MPI_SUM,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
71e97f6855SBarry Smith   root--;
72ae95e784SBarry Smith   ierr = MPI_Bcast(X,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
73ae95e784SBarry Smith   ierr = MPI_Bcast(Y,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
745edff71fSBarry Smith   ierr = DMDAVecRestoreArrayRead(dacoors,coors,&c);CHKERRQ(ierr);
750a0169f3SBarry Smith   PetscFunctionReturn(0);
760a0169f3SBarry Smith }
770a0169f3SBarry Smith 
78*69e49704SSatish Balay /*@
79db05f41bSBarry Smith    DMDAGetRay - Returns a vector on process zero that contains a row or column of the values in a DMDA vector
80db05f41bSBarry Smith 
81db05f41bSBarry Smith    Collective on DMDA
82db05f41bSBarry Smith 
83db05f41bSBarry Smith    Input Parameters:
84db05f41bSBarry Smith +  da - the distributed array
85db05f41bSBarry Smith .  vec - the vector
86db05f41bSBarry Smith .  dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z
87db05f41bSBarry Smith -  gp - global grid point number in this direction
88db05f41bSBarry Smith 
89db05f41bSBarry Smith    Output Parameters:
90db05f41bSBarry Smith +  newvec - the new vector that can hold the values (size zero on all processes except process 0)
91db05f41bSBarry Smith -  scatter - the VecScatter that will map from the original vector to the slice
92db05f41bSBarry Smith 
93db05f41bSBarry Smith    Level: advanced
94db05f41bSBarry Smith 
95db05f41bSBarry Smith    Notes:
96db05f41bSBarry Smith    All processors that share the DMDA must call this with the same gp value
97db05f41bSBarry Smith 
98db05f41bSBarry Smith .keywords: distributed array, get, processor subset
99db05f41bSBarry Smith @*/
100db05f41bSBarry Smith PetscErrorCode  DMDAGetRay(DM da,DMDADirection dir,PetscInt gp,Vec *newvec,VecScatter *scatter)
101db05f41bSBarry Smith {
102db05f41bSBarry Smith   PetscMPIInt    rank;
103db05f41bSBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
104db05f41bSBarry Smith   PetscErrorCode ierr;
105db05f41bSBarry Smith   IS             is;
106db05f41bSBarry Smith   AO             ao;
107db05f41bSBarry Smith   Vec            vec;
108d1212d36SBarry Smith   PetscInt       *indices,i,j;
109db05f41bSBarry Smith 
110db05f41bSBarry Smith   PetscFunctionBegin;
111c73cfb54SMatthew G. Knepley   if (da->dim == 3) SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_SUP, "Cannot get slice from 3d DMDA");
112ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) da), &rank);CHKERRQ(ierr);
113c2001ae3SMatthew G. Knepley   ierr = DMDAGetAO(da, &ao);CHKERRQ(ierr);
114db05f41bSBarry Smith   if (!rank) {
115c73cfb54SMatthew G. Knepley     if (da->dim == 1) {
116c2001ae3SMatthew G. Knepley       if (dir == DMDA_X) {
117785e854fSJed Brown         ierr = PetscMalloc1(dd->w, &indices);CHKERRQ(ierr);
118c2001ae3SMatthew G. Knepley         indices[0] = dd->w*gp;
119eca19243SMatthew G. Knepley         for (i = 1; i < dd->w; ++i) indices[i] = indices[i-1] + 1;
120c2001ae3SMatthew G. Knepley         ierr = AOApplicationToPetsc(ao, dd->w, indices);CHKERRQ(ierr);
121c2001ae3SMatthew G. Knepley         ierr = VecCreate(PETSC_COMM_SELF, newvec);CHKERRQ(ierr);
122c2001ae3SMatthew G. Knepley         ierr = VecSetBlockSize(*newvec, dd->w);CHKERRQ(ierr);
123c2001ae3SMatthew G. Knepley         ierr = VecSetSizes(*newvec, dd->w, PETSC_DETERMINE);CHKERRQ(ierr);
124c2001ae3SMatthew G. Knepley         ierr = VecSetType(*newvec, VECSEQ);CHKERRQ(ierr);
125c2001ae3SMatthew G. Knepley         ierr = ISCreateGeneral(PETSC_COMM_SELF, dd->w, indices, PETSC_OWN_POINTER, &is);CHKERRQ(ierr);
126c2001ae3SMatthew G. Knepley       } else if (dir == DMDA_Y) SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_SUP, "Cannot get Y slice from 1d DMDA");
127c2001ae3SMatthew G. Knepley       else SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_ARG_OUTOFRANGE, "Unknown DMDADirection");
128c2001ae3SMatthew G. Knepley     } else {
129db05f41bSBarry Smith       if (dir == DMDA_Y) {
130785e854fSJed Brown         ierr       = PetscMalloc1(dd->w*dd->M,&indices);CHKERRQ(ierr);
131d1212d36SBarry Smith         indices[0] = gp*dd->M*dd->w;
1328865f1eaSKarl Rupp         for (i=1; i<dd->M*dd->w; i++) indices[i] = indices[i-1] + 1;
1338865f1eaSKarl Rupp 
134d1212d36SBarry Smith         ierr = AOApplicationToPetsc(ao,dd->M*dd->w,indices);CHKERRQ(ierr);
135db05f41bSBarry Smith         ierr = VecCreate(PETSC_COMM_SELF,newvec);CHKERRQ(ierr);
136db05f41bSBarry Smith         ierr = VecSetBlockSize(*newvec,dd->w);CHKERRQ(ierr);
137db05f41bSBarry Smith         ierr = VecSetSizes(*newvec,dd->M*dd->w,PETSC_DETERMINE);CHKERRQ(ierr);
138db05f41bSBarry Smith         ierr = VecSetType(*newvec,VECSEQ);CHKERRQ(ierr);
139d1212d36SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,dd->w*dd->M,indices,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
140db05f41bSBarry Smith       } else if (dir == DMDA_X) {
141785e854fSJed Brown         ierr       = PetscMalloc1(dd->w*dd->N,&indices);CHKERRQ(ierr);
142d1212d36SBarry Smith         indices[0] = dd->w*gp;
143d1212d36SBarry Smith         for (j=1; j<dd->w; j++) indices[j] = indices[j-1] + 1;
144d1212d36SBarry Smith         for (i=1; i<dd->N; i++) {
145d1212d36SBarry Smith           indices[i*dd->w] = indices[i*dd->w-1] + dd->w*dd->M - dd->w + 1;
146d1212d36SBarry Smith           for (j=1; j<dd->w; j++) indices[i*dd->w + j] = indices[i*dd->w + j - 1] + 1;
147d1212d36SBarry Smith         }
148d1212d36SBarry Smith         ierr = AOApplicationToPetsc(ao,dd->w*dd->N,indices);CHKERRQ(ierr);
149db05f41bSBarry Smith         ierr = VecCreate(PETSC_COMM_SELF,newvec);CHKERRQ(ierr);
150db05f41bSBarry Smith         ierr = VecSetBlockSize(*newvec,dd->w);CHKERRQ(ierr);
151db05f41bSBarry Smith         ierr = VecSetSizes(*newvec,dd->N*dd->w,PETSC_DETERMINE);CHKERRQ(ierr);
152db05f41bSBarry Smith         ierr = VecSetType(*newvec,VECSEQ);CHKERRQ(ierr);
153d1212d36SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,dd->w*dd->N,indices,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
154db05f41bSBarry Smith       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown DMDADirection");
155c2001ae3SMatthew G. Knepley     }
156db05f41bSBarry Smith   } else {
157db05f41bSBarry Smith     ierr = VecCreateSeq(PETSC_COMM_SELF, 0, newvec);CHKERRQ(ierr);
158d1212d36SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, 0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
159db05f41bSBarry Smith   }
160db05f41bSBarry Smith   ierr = DMGetGlobalVector(da, &vec);CHKERRQ(ierr);
1610298fd71SBarry Smith   ierr = VecScatterCreate(vec, is, *newvec, NULL, scatter);CHKERRQ(ierr);
162db05f41bSBarry Smith   ierr = DMRestoreGlobalVector(da, &vec);CHKERRQ(ierr);
163db05f41bSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
164db05f41bSBarry Smith   PetscFunctionReturn(0);
165db05f41bSBarry Smith }
166db05f41bSBarry Smith 
16747c6ae99SBarry Smith /*@C
168aa219208SBarry Smith    DMDAGetProcessorSubset - Returns a communicator consisting only of the
169aa219208SBarry Smith    processors in a DMDA that own a particular global x, y, or z grid point
17047c6ae99SBarry Smith    (corresponding to a logical plane in a 3D grid or a line in a 2D grid).
17147c6ae99SBarry Smith 
172aa219208SBarry Smith    Collective on DMDA
17347c6ae99SBarry Smith 
17447c6ae99SBarry Smith    Input Parameters:
17547c6ae99SBarry Smith +  da - the distributed array
176aa219208SBarry Smith .  dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z
17747c6ae99SBarry Smith -  gp - global grid point number in this direction
17847c6ae99SBarry Smith 
17947c6ae99SBarry Smith    Output Parameters:
18047c6ae99SBarry Smith .  comm - new communicator
18147c6ae99SBarry Smith 
18247c6ae99SBarry Smith    Level: advanced
18347c6ae99SBarry Smith 
18447c6ae99SBarry Smith    Notes:
185aa219208SBarry Smith    All processors that share the DMDA must call this with the same gp value
18647c6ae99SBarry Smith 
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 
19147c6ae99SBarry Smith .keywords: distributed array, get, processor subset
19247c6ae99SBarry Smith @*/
1937087cfbeSBarry Smith PetscErrorCode  DMDAGetProcessorSubset(DM da,DMDADirection dir,PetscInt gp,MPI_Comm *comm)
19447c6ae99SBarry Smith {
19547c6ae99SBarry Smith   MPI_Group      group,subgroup;
19647c6ae99SBarry Smith   PetscErrorCode ierr;
19747c6ae99SBarry Smith   PetscInt       i,ict,flag,*owners,xs,xm,ys,ym,zs,zm;
1980298fd71SBarry Smith   PetscMPIInt    size,*ranks = NULL;
19947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
20047c6ae99SBarry Smith 
20147c6ae99SBarry Smith   PetscFunctionBegin;
20247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
20347c6ae99SBarry Smith   flag = 0;
204aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
205ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr);
206aa219208SBarry Smith   if (dir == DMDA_Z) {
207c73cfb54SMatthew G. Knepley     if (da->dim < 3) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Z invalid for DMDA dim < 3");
20847c6ae99SBarry Smith     if (gp < 0 || gp > dd->P) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
20947c6ae99SBarry Smith     if (gp >= zs && gp < zs+zm) flag = 1;
210aa219208SBarry Smith   } else if (dir == DMDA_Y) {
211c73cfb54SMatthew G. Knepley     if (da->dim == 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Y invalid for DMDA dim = 1");
21247c6ae99SBarry Smith     if (gp < 0 || gp > dd->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
21347c6ae99SBarry Smith     if (gp >= ys && gp < ys+ym) flag = 1;
214aa219208SBarry Smith   } else if (dir == DMDA_X) {
21547c6ae99SBarry Smith     if (gp < 0 || gp > dd->M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
21647c6ae99SBarry Smith     if (gp >= xs && gp < xs+xm) flag = 1;
217ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction");
21847c6ae99SBarry Smith 
219dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&owners,size,&ranks);CHKERRQ(ierr);
220ce94432eSBarry Smith   ierr = MPI_Allgather(&flag,1,MPIU_INT,owners,1,MPIU_INT,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
22147c6ae99SBarry Smith   ict  = 0;
222c73cfb54SMatthew G. Knepley   ierr = PetscInfo2(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",da->dim,(int)dir);CHKERRQ(ierr);
22347c6ae99SBarry Smith   for (i=0; i<size; i++) {
22447c6ae99SBarry Smith     if (owners[i]) {
22547c6ae99SBarry Smith       ranks[ict] = i; ict++;
22647c6ae99SBarry Smith       ierr       = PetscInfo1(da,"%D ",i);CHKERRQ(ierr);
22747c6ae99SBarry Smith     }
22847c6ae99SBarry Smith   }
22947c6ae99SBarry Smith   ierr = PetscInfo(da,"\n");CHKERRQ(ierr);
230ce94432eSBarry Smith   ierr = MPI_Comm_group(PetscObjectComm((PetscObject)da),&group);CHKERRQ(ierr);
23147c6ae99SBarry Smith   ierr = MPI_Group_incl(group,ict,ranks,&subgroup);CHKERRQ(ierr);
232ce94432eSBarry Smith   ierr = MPI_Comm_create(PetscObjectComm((PetscObject)da),subgroup,comm);CHKERRQ(ierr);
23347c6ae99SBarry Smith   ierr = MPI_Group_free(&subgroup);CHKERRQ(ierr);
23447c6ae99SBarry Smith   ierr = MPI_Group_free(&group);CHKERRQ(ierr);
23547c6ae99SBarry Smith   ierr = PetscFree2(owners,ranks);CHKERRQ(ierr);
23647c6ae99SBarry Smith   PetscFunctionReturn(0);
23747c6ae99SBarry Smith }
23847c6ae99SBarry Smith 
23947c6ae99SBarry Smith /*@C
240aa219208SBarry Smith    DMDAGetProcessorSubsets - Returns communicators consisting only of the
241aa219208SBarry Smith    processors in a DMDA adjacent in a particular dimension,
24247c6ae99SBarry Smith    corresponding to a logical plane in a 3D grid or a line in a 2D grid.
24347c6ae99SBarry Smith 
244aa219208SBarry Smith    Collective on DMDA
24547c6ae99SBarry Smith 
24647c6ae99SBarry Smith    Input Parameters:
24747c6ae99SBarry Smith +  da - the distributed array
248aa219208SBarry Smith -  dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z
24947c6ae99SBarry Smith 
25047c6ae99SBarry Smith    Output Parameters:
25147c6ae99SBarry Smith .  subcomm - new communicator
25247c6ae99SBarry Smith 
25347c6ae99SBarry Smith    Level: advanced
25447c6ae99SBarry Smith 
25547c6ae99SBarry Smith    Notes:
25647c6ae99SBarry Smith    This routine is useful for distributing one-dimensional data in a tensor product grid.
25747c6ae99SBarry Smith 
25847c6ae99SBarry Smith .keywords: distributed array, get, processor subset
25947c6ae99SBarry Smith @*/
2607087cfbeSBarry Smith PetscErrorCode  DMDAGetProcessorSubsets(DM da, DMDADirection dir, MPI_Comm *subcomm)
26147c6ae99SBarry Smith {
26247c6ae99SBarry Smith   MPI_Comm       comm;
26347c6ae99SBarry Smith   MPI_Group      group, subgroup;
26447c6ae99SBarry Smith   PetscInt       subgroupSize = 0;
26547c6ae99SBarry Smith   PetscInt       *firstPoints;
2660298fd71SBarry Smith   PetscMPIInt    size, *subgroupRanks = NULL;
26747c6ae99SBarry Smith   PetscInt       xs, xm, ys, ym, zs, zm, firstPoint, p;
26847c6ae99SBarry Smith   PetscErrorCode ierr;
26947c6ae99SBarry Smith 
27047c6ae99SBarry Smith   PetscFunctionBegin;
27147c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
27282f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
273aa219208SBarry Smith   ierr = DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm);CHKERRQ(ierr);
27447c6ae99SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
275aa219208SBarry Smith   if (dir == DMDA_Z) {
276c73cfb54SMatthew G. Knepley     if (da->dim < 3) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Z invalid for DMDA dim < 3");
27747c6ae99SBarry Smith     firstPoint = zs;
278aa219208SBarry Smith   } else if (dir == DMDA_Y) {
279c73cfb54SMatthew G. Knepley     if (da->dim == 1) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Y invalid for DMDA dim = 1");
28047c6ae99SBarry Smith     firstPoint = ys;
281aa219208SBarry Smith   } else if (dir == DMDA_X) {
28247c6ae99SBarry Smith     firstPoint = xs;
28347c6ae99SBarry Smith   } else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction");
28447c6ae99SBarry Smith 
285dcca6d9dSJed Brown   ierr = PetscMalloc2(size, &firstPoints, size, &subgroupRanks);CHKERRQ(ierr);
28647c6ae99SBarry Smith   ierr = MPI_Allgather(&firstPoint, 1, MPIU_INT, firstPoints, 1, MPIU_INT, comm);CHKERRQ(ierr);
287c73cfb54SMatthew G. Knepley   ierr = PetscInfo2(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",da->dim,(int)dir);CHKERRQ(ierr);
28847c6ae99SBarry Smith   for (p = 0; p < size; ++p) {
28947c6ae99SBarry Smith     if (firstPoints[p] == firstPoint) {
29047c6ae99SBarry Smith       subgroupRanks[subgroupSize++] = p;
29147c6ae99SBarry Smith       ierr = PetscInfo1(da, "%D ", p);CHKERRQ(ierr);
29247c6ae99SBarry Smith     }
29347c6ae99SBarry Smith   }
29447c6ae99SBarry Smith   ierr = PetscInfo(da, "\n");CHKERRQ(ierr);
29547c6ae99SBarry Smith   ierr = MPI_Comm_group(comm, &group);CHKERRQ(ierr);
29647c6ae99SBarry Smith   ierr = MPI_Group_incl(group, subgroupSize, subgroupRanks, &subgroup);CHKERRQ(ierr);
29747c6ae99SBarry Smith   ierr = MPI_Comm_create(comm, subgroup, subcomm);CHKERRQ(ierr);
29847c6ae99SBarry Smith   ierr = MPI_Group_free(&subgroup);CHKERRQ(ierr);
29947c6ae99SBarry Smith   ierr = MPI_Group_free(&group);CHKERRQ(ierr);
30047c6ae99SBarry Smith   ierr = PetscFree2(firstPoints, subgroupRanks);CHKERRQ(ierr);
30147c6ae99SBarry Smith   PetscFunctionReturn(0);
30247c6ae99SBarry Smith }
303