xref: /petsc/src/dm/impls/da/dasub.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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
156b867d5aSJose E. Roman .  x  - the first physical coordinate
166b867d5aSJose E. Roman .  y  - the second physical coordinate
176b867d5aSJose E. Roman -  z  - the third physical coordinate
180a0169f3SBarry Smith 
190a0169f3SBarry Smith    Output Parameters:
206b867d5aSJose E. Roman +  II - the first logical coordinate (-1 on processes that do not contain that point)
216b867d5aSJose E. Roman .  JJ - the second logical coordinate (-1 on processes that do not contain that point)
226b867d5aSJose E. Roman .  KK - the third logical coordinate (-1 on processes that do not contain that point)
236b867d5aSJose E. Roman .  X  - (optional) the first coordinate of the located grid point
246b867d5aSJose E. Roman .  Y  - (optional) the second coordinate of the located grid point
256b867d5aSJose E. Roman -  Z  - (optional) the third coordinate of the located grid point
260a0169f3SBarry Smith 
270a0169f3SBarry Smith    Level: advanced
280a0169f3SBarry Smith 
290a0169f3SBarry Smith    Notes:
300a0169f3SBarry Smith    All processors that share the DMDA must call this with the same coordinate value
310a0169f3SBarry Smith 
320a0169f3SBarry Smith @*/
3333c53e3fSSatish Balay PetscErrorCode  DMDAGetLogicalCoordinate(DM da,PetscScalar x,PetscScalar y,PetscScalar z,PetscInt *II,PetscInt *JJ,PetscInt *KK,PetscScalar *X,PetscScalar *Y,PetscScalar *Z)
340a0169f3SBarry Smith {
350a0169f3SBarry Smith   PetscErrorCode ierr;
360a0169f3SBarry Smith   Vec            coors;
370a0169f3SBarry Smith   DM             dacoors;
380a0169f3SBarry Smith   DMDACoor2d     **c;
390a0169f3SBarry Smith   PetscInt       i,j,xs,xm,ys,ym;
40e97f6855SBarry Smith   PetscReal      d,D = PETSC_MAX_REAL,Dv;
41e97f6855SBarry Smith   PetscMPIInt    rank,root;
420a0169f3SBarry Smith 
430a0169f3SBarry Smith   PetscFunctionBegin;
44*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(da->dim == 1,PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot get point from 1d DMDA");
45*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(da->dim == 3,PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot get point from 3d DMDA");
460a0169f3SBarry Smith 
4733c53e3fSSatish Balay   *II = -1;
48e97f6855SBarry Smith   *JJ = -1;
490a0169f3SBarry Smith 
500a0169f3SBarry Smith   ierr = DMGetCoordinateDM(da,&dacoors);CHKERRQ(ierr);
510a0169f3SBarry Smith   ierr = DMDAGetCorners(dacoors,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
520a0169f3SBarry Smith   ierr = DMGetCoordinates(da,&coors);CHKERRQ(ierr);
535edff71fSBarry Smith   ierr = DMDAVecGetArrayRead(dacoors,coors,&c);CHKERRQ(ierr);
540a0169f3SBarry Smith   for (j=ys; j<ys+ym; j++) {
550a0169f3SBarry Smith     for (i=xs; i<xs+xm; i++) {
560a0169f3SBarry Smith       d = PetscSqrtReal(PetscRealPart((c[j][i].x - x)*(c[j][i].x - x) + (c[j][i].y - y)*(c[j][i].y - y)));
570a0169f3SBarry Smith       if (d < D) {
580a0169f3SBarry Smith         D   = d;
5933c53e3fSSatish Balay         *II = i;
6033c53e3fSSatish Balay         *JJ = j;
610a0169f3SBarry Smith       }
620a0169f3SBarry Smith     }
630a0169f3SBarry Smith   }
64820f2d46SBarry Smith   ierr = MPIU_Allreduce(&D,&Dv,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr);
65e97f6855SBarry Smith   if (D != Dv) {
66e97f6855SBarry Smith     *II  = -1;
67e97f6855SBarry Smith     *JJ  = -1;
68e97f6855SBarry Smith     rank = 0;
69e97f6855SBarry Smith   } else {
70e97f6855SBarry Smith     *X = c[*JJ][*II].x;
71e97f6855SBarry Smith     *Y = c[*JJ][*II].y;
72ffc4695bSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRMPI(ierr);
73e97f6855SBarry Smith     rank++;
74e97f6855SBarry Smith   }
75820f2d46SBarry Smith   ierr = MPIU_Allreduce(&rank,&root,1,MPI_INT,MPI_SUM,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr);
76e97f6855SBarry Smith   root--;
77ffc4695bSBarry Smith   ierr = MPI_Bcast(X,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr);
78ffc4695bSBarry Smith   ierr = MPI_Bcast(Y,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr);
795edff71fSBarry Smith   ierr = DMDAVecRestoreArrayRead(dacoors,coors,&c);CHKERRQ(ierr);
800a0169f3SBarry Smith   PetscFunctionReturn(0);
810a0169f3SBarry Smith }
820a0169f3SBarry Smith 
8369e49704SSatish Balay /*@
84db05f41bSBarry Smith    DMDAGetRay - Returns a vector on process zero that contains a row or column of the values in a DMDA vector
85db05f41bSBarry Smith 
86db05f41bSBarry Smith    Collective on DMDA
87db05f41bSBarry Smith 
88db05f41bSBarry Smith    Input Parameters:
89db05f41bSBarry Smith +  da - the distributed array
903ee9839eSMatthew G. Knepley .  dir - Cartesian direction, either DM_X, DM_Y, or DM_Z
91db05f41bSBarry Smith -  gp - global grid point number in this direction
92db05f41bSBarry Smith 
93db05f41bSBarry Smith    Output Parameters:
94db05f41bSBarry Smith +  newvec - the new vector that can hold the values (size zero on all processes except process 0)
95db05f41bSBarry Smith -  scatter - the VecScatter that will map from the original vector to the slice
96db05f41bSBarry Smith 
97db05f41bSBarry Smith    Level: advanced
98db05f41bSBarry Smith 
99db05f41bSBarry Smith    Notes:
100db05f41bSBarry Smith    All processors that share the DMDA must call this with the same gp value
101db05f41bSBarry Smith 
102db05f41bSBarry Smith @*/
1033ee9839eSMatthew G. Knepley PetscErrorCode  DMDAGetRay(DM da,DMDirection dir,PetscInt gp,Vec *newvec,VecScatter *scatter)
104db05f41bSBarry Smith {
105db05f41bSBarry Smith   PetscMPIInt    rank;
106db05f41bSBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
107db05f41bSBarry Smith   PetscErrorCode ierr;
108db05f41bSBarry Smith   IS             is;
109db05f41bSBarry Smith   AO             ao;
110db05f41bSBarry Smith   Vec            vec;
111d1212d36SBarry Smith   PetscInt       *indices,i,j;
112db05f41bSBarry Smith 
113db05f41bSBarry Smith   PetscFunctionBegin;
114*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(da->dim == 3,PetscObjectComm((PetscObject) da), PETSC_ERR_SUP, "Cannot get slice from 3d DMDA");
115ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) da), &rank);CHKERRMPI(ierr);
116c2001ae3SMatthew G. Knepley   ierr = DMDAGetAO(da, &ao);CHKERRQ(ierr);
117dd400576SPatrick Sanan   if (rank == 0) {
118c73cfb54SMatthew G. Knepley     if (da->dim == 1) {
1193ee9839eSMatthew G. Knepley       if (dir == DM_X) {
120785e854fSJed Brown         ierr = PetscMalloc1(dd->w, &indices);CHKERRQ(ierr);
121c2001ae3SMatthew G. Knepley         indices[0] = dd->w*gp;
122eca19243SMatthew G. Knepley         for (i = 1; i < dd->w; ++i) indices[i] = indices[i-1] + 1;
123c2001ae3SMatthew G. Knepley         ierr = AOApplicationToPetsc(ao, dd->w, indices);CHKERRQ(ierr);
124c2001ae3SMatthew G. Knepley         ierr = VecCreate(PETSC_COMM_SELF, newvec);CHKERRQ(ierr);
125c2001ae3SMatthew G. Knepley         ierr = VecSetBlockSize(*newvec, dd->w);CHKERRQ(ierr);
126c2001ae3SMatthew G. Knepley         ierr = VecSetSizes(*newvec, dd->w, PETSC_DETERMINE);CHKERRQ(ierr);
127c2001ae3SMatthew G. Knepley         ierr = VecSetType(*newvec, VECSEQ);CHKERRQ(ierr);
128c2001ae3SMatthew G. Knepley         ierr = ISCreateGeneral(PETSC_COMM_SELF, dd->w, indices, PETSC_OWN_POINTER, &is);CHKERRQ(ierr);
129*2c71b3e2SJacob Faibussowitsch       } else PetscCheckFalse(dir == DM_Y,PetscObjectComm((PetscObject) da), PETSC_ERR_SUP, "Cannot get Y slice from 1d DMDA");
1303ee9839eSMatthew G. Knepley       else SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_ARG_OUTOFRANGE, "Unknown DMDirection");
131c2001ae3SMatthew G. Knepley     } else {
1323ee9839eSMatthew G. Knepley       if (dir == DM_Y) {
133785e854fSJed Brown         ierr       = PetscMalloc1(dd->w*dd->M,&indices);CHKERRQ(ierr);
134d1212d36SBarry Smith         indices[0] = gp*dd->M*dd->w;
1358865f1eaSKarl Rupp         for (i=1; i<dd->M*dd->w; i++) indices[i] = indices[i-1] + 1;
1368865f1eaSKarl Rupp 
137d1212d36SBarry Smith         ierr = AOApplicationToPetsc(ao,dd->M*dd->w,indices);CHKERRQ(ierr);
138db05f41bSBarry Smith         ierr = VecCreate(PETSC_COMM_SELF,newvec);CHKERRQ(ierr);
139db05f41bSBarry Smith         ierr = VecSetBlockSize(*newvec,dd->w);CHKERRQ(ierr);
140db05f41bSBarry Smith         ierr = VecSetSizes(*newvec,dd->M*dd->w,PETSC_DETERMINE);CHKERRQ(ierr);
141db05f41bSBarry Smith         ierr = VecSetType(*newvec,VECSEQ);CHKERRQ(ierr);
142d1212d36SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,dd->w*dd->M,indices,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
1433ee9839eSMatthew G. Knepley       } else if (dir == DM_X) {
144785e854fSJed Brown         ierr       = PetscMalloc1(dd->w*dd->N,&indices);CHKERRQ(ierr);
145d1212d36SBarry Smith         indices[0] = dd->w*gp;
146d1212d36SBarry Smith         for (j=1; j<dd->w; j++) indices[j] = indices[j-1] + 1;
147d1212d36SBarry Smith         for (i=1; i<dd->N; i++) {
148d1212d36SBarry Smith           indices[i*dd->w] = indices[i*dd->w-1] + dd->w*dd->M - dd->w + 1;
149d1212d36SBarry Smith           for (j=1; j<dd->w; j++) indices[i*dd->w + j] = indices[i*dd->w + j - 1] + 1;
150d1212d36SBarry Smith         }
151d1212d36SBarry Smith         ierr = AOApplicationToPetsc(ao,dd->w*dd->N,indices);CHKERRQ(ierr);
152db05f41bSBarry Smith         ierr = VecCreate(PETSC_COMM_SELF,newvec);CHKERRQ(ierr);
153db05f41bSBarry Smith         ierr = VecSetBlockSize(*newvec,dd->w);CHKERRQ(ierr);
154db05f41bSBarry Smith         ierr = VecSetSizes(*newvec,dd->N*dd->w,PETSC_DETERMINE);CHKERRQ(ierr);
155db05f41bSBarry Smith         ierr = VecSetType(*newvec,VECSEQ);CHKERRQ(ierr);
156d1212d36SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,dd->w*dd->N,indices,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
1573ee9839eSMatthew G. Knepley       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown DMDirection");
158c2001ae3SMatthew G. Knepley     }
159db05f41bSBarry Smith   } else {
160db05f41bSBarry Smith     ierr = VecCreateSeq(PETSC_COMM_SELF, 0, newvec);CHKERRQ(ierr);
161ea78f98cSLisandro Dalcin     ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is);CHKERRQ(ierr);
162db05f41bSBarry Smith   }
163db05f41bSBarry Smith   ierr = DMGetGlobalVector(da, &vec);CHKERRQ(ierr);
1649448b7f1SJunchao Zhang   ierr = VecScatterCreate(vec, is, *newvec, NULL, scatter);CHKERRQ(ierr);
165db05f41bSBarry Smith   ierr = DMRestoreGlobalVector(da, &vec);CHKERRQ(ierr);
166db05f41bSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
167db05f41bSBarry Smith   PetscFunctionReturn(0);
168db05f41bSBarry Smith }
169db05f41bSBarry Smith 
17047c6ae99SBarry Smith /*@C
171aa219208SBarry Smith    DMDAGetProcessorSubset - Returns a communicator consisting only of the
172aa219208SBarry Smith    processors in a DMDA that own a particular global x, y, or z grid point
17347c6ae99SBarry Smith    (corresponding to a logical plane in a 3D grid or a line in a 2D grid).
17447c6ae99SBarry Smith 
175d083f849SBarry Smith    Collective on da
17647c6ae99SBarry Smith 
17747c6ae99SBarry Smith    Input Parameters:
17847c6ae99SBarry Smith +  da - the distributed array
1793ee9839eSMatthew G. Knepley .  dir - Cartesian direction, either DM_X, DM_Y, or DM_Z
18047c6ae99SBarry Smith -  gp - global grid point number in this direction
18147c6ae99SBarry Smith 
18297bb3fdcSJose E. Roman    Output Parameter:
18347c6ae99SBarry Smith .  comm - new communicator
18447c6ae99SBarry Smith 
18547c6ae99SBarry Smith    Level: advanced
18647c6ae99SBarry Smith 
18747c6ae99SBarry Smith    Notes:
188aa219208SBarry Smith    All processors that share the DMDA must call this with the same gp value
18947c6ae99SBarry Smith 
190bf3faf6aSSatish Balay    After use, comm should be freed with MPI_Comm_free()
191bf3faf6aSSatish Balay 
19247c6ae99SBarry Smith    This routine is particularly useful to compute boundary conditions
19347c6ae99SBarry Smith    or other application-specific calculations that require manipulating
19447c6ae99SBarry Smith    sets of data throughout a logical plane of grid points.
19547c6ae99SBarry Smith 
196f5f57ec0SBarry Smith    Not supported from Fortran
197f5f57ec0SBarry Smith 
19847c6ae99SBarry Smith @*/
1993ee9839eSMatthew G. Knepley PetscErrorCode  DMDAGetProcessorSubset(DM da,DMDirection dir,PetscInt gp,MPI_Comm *comm)
20047c6ae99SBarry Smith {
20147c6ae99SBarry Smith   MPI_Group      group,subgroup;
20247c6ae99SBarry Smith   PetscErrorCode ierr;
20347c6ae99SBarry Smith   PetscInt       i,ict,flag,*owners,xs,xm,ys,ym,zs,zm;
2040298fd71SBarry Smith   PetscMPIInt    size,*ranks = NULL;
20547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
20647c6ae99SBarry Smith 
20747c6ae99SBarry Smith   PetscFunctionBegin;
208a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
20947c6ae99SBarry Smith   flag = 0;
210aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
211ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRMPI(ierr);
2123ee9839eSMatthew G. Knepley   if (dir == DM_Z) {
213*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(da->dim < 3,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"DM_Z invalid for DMDA dim < 3");
214*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(gp < 0 || gp > dd->P,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
21547c6ae99SBarry Smith     if (gp >= zs && gp < zs+zm) flag = 1;
2163ee9839eSMatthew G. Knepley   } else if (dir == DM_Y) {
217*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(da->dim == 1,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"DM_Y invalid for DMDA dim = 1");
218*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(gp < 0 || gp > dd->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
21947c6ae99SBarry Smith     if (gp >= ys && gp < ys+ym) flag = 1;
2203ee9839eSMatthew G. Knepley   } else if (dir == DM_X) {
221*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(gp < 0 || gp > dd->M,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point");
22247c6ae99SBarry Smith     if (gp >= xs && gp < xs+xm) flag = 1;
223ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction");
22447c6ae99SBarry Smith 
225dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&owners,size,&ranks);CHKERRQ(ierr);
226ffc4695bSBarry Smith   ierr = MPI_Allgather(&flag,1,MPIU_INT,owners,1,MPIU_INT,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr);
22747c6ae99SBarry Smith   ict  = 0;
2287d3de750SJacob Faibussowitsch   ierr = PetscInfo(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",da->dim,(int)dir);CHKERRQ(ierr);
22947c6ae99SBarry Smith   for (i=0; i<size; i++) {
23047c6ae99SBarry Smith     if (owners[i]) {
23147c6ae99SBarry Smith       ranks[ict] = i; ict++;
2327d3de750SJacob Faibussowitsch       ierr       = PetscInfo(da,"%D ",i);CHKERRQ(ierr);
23347c6ae99SBarry Smith     }
23447c6ae99SBarry Smith   }
23547c6ae99SBarry Smith   ierr = PetscInfo(da,"\n");CHKERRQ(ierr);
236ffc4695bSBarry Smith   ierr = MPI_Comm_group(PetscObjectComm((PetscObject)da),&group);CHKERRMPI(ierr);
237ffc4695bSBarry Smith   ierr = MPI_Group_incl(group,ict,ranks,&subgroup);CHKERRMPI(ierr);
238ffc4695bSBarry Smith   ierr = MPI_Comm_create(PetscObjectComm((PetscObject)da),subgroup,comm);CHKERRMPI(ierr);
239ffc4695bSBarry Smith   ierr = MPI_Group_free(&subgroup);CHKERRMPI(ierr);
240ffc4695bSBarry Smith   ierr = MPI_Group_free(&group);CHKERRMPI(ierr);
24147c6ae99SBarry Smith   ierr = PetscFree2(owners,ranks);CHKERRQ(ierr);
24247c6ae99SBarry Smith   PetscFunctionReturn(0);
24347c6ae99SBarry Smith }
24447c6ae99SBarry Smith 
24547c6ae99SBarry Smith /*@C
246aa219208SBarry Smith    DMDAGetProcessorSubsets - Returns communicators consisting only of the
247aa219208SBarry Smith    processors in a DMDA adjacent in a particular dimension,
24847c6ae99SBarry Smith    corresponding to a logical plane in a 3D grid or a line in a 2D grid.
24947c6ae99SBarry Smith 
250d083f849SBarry Smith    Collective on da
25147c6ae99SBarry Smith 
25247c6ae99SBarry Smith    Input Parameters:
25347c6ae99SBarry Smith +  da - the distributed array
2543ee9839eSMatthew G. Knepley -  dir - Cartesian direction, either DM_X, DM_Y, or DM_Z
25547c6ae99SBarry Smith 
25697bb3fdcSJose E. Roman    Output Parameter:
25747c6ae99SBarry Smith .  subcomm - new communicator
25847c6ae99SBarry Smith 
25947c6ae99SBarry Smith    Level: advanced
26047c6ae99SBarry Smith 
26147c6ae99SBarry Smith    Notes:
26247c6ae99SBarry Smith    This routine is useful for distributing one-dimensional data in a tensor product grid.
26347c6ae99SBarry Smith 
264bf3faf6aSSatish Balay    After use, comm should be freed with MPI_Comm_free()
265bf3faf6aSSatish Balay 
266f5f57ec0SBarry Smith    Not supported from Fortran
267f5f57ec0SBarry Smith 
26847c6ae99SBarry Smith @*/
2693ee9839eSMatthew G. Knepley PetscErrorCode  DMDAGetProcessorSubsets(DM da, DMDirection dir, MPI_Comm *subcomm)
27047c6ae99SBarry Smith {
27147c6ae99SBarry Smith   MPI_Comm       comm;
27247c6ae99SBarry Smith   MPI_Group      group, subgroup;
27347c6ae99SBarry Smith   PetscInt       subgroupSize = 0;
27447c6ae99SBarry Smith   PetscInt       *firstPoints;
2750298fd71SBarry Smith   PetscMPIInt    size, *subgroupRanks = NULL;
27647c6ae99SBarry Smith   PetscInt       xs, xm, ys, ym, zs, zm, firstPoint, p;
27747c6ae99SBarry Smith   PetscErrorCode ierr;
27847c6ae99SBarry Smith 
27947c6ae99SBarry Smith   PetscFunctionBegin;
280a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
28182f516ccSBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
282aa219208SBarry Smith   ierr = DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm);CHKERRQ(ierr);
283ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
2843ee9839eSMatthew G. Knepley   if (dir == DM_Z) {
285*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(da->dim < 3,comm,PETSC_ERR_ARG_OUTOFRANGE,"DM_Z invalid for DMDA dim < 3");
28647c6ae99SBarry Smith     firstPoint = zs;
2873ee9839eSMatthew G. Knepley   } else if (dir == DM_Y) {
288*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(da->dim == 1,comm,PETSC_ERR_ARG_OUTOFRANGE,"DM_Y invalid for DMDA dim = 1");
28947c6ae99SBarry Smith     firstPoint = ys;
2903ee9839eSMatthew G. Knepley   } else if (dir == DM_X) {
29147c6ae99SBarry Smith     firstPoint = xs;
29247c6ae99SBarry Smith   } else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction");
29347c6ae99SBarry Smith 
294dcca6d9dSJed Brown   ierr = PetscMalloc2(size, &firstPoints, size, &subgroupRanks);CHKERRQ(ierr);
295ffc4695bSBarry Smith   ierr = MPI_Allgather(&firstPoint, 1, MPIU_INT, firstPoints, 1, MPIU_INT, comm);CHKERRMPI(ierr);
2967d3de750SJacob Faibussowitsch   ierr = PetscInfo(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",da->dim,(int)dir);CHKERRQ(ierr);
29747c6ae99SBarry Smith   for (p = 0; p < size; ++p) {
29847c6ae99SBarry Smith     if (firstPoints[p] == firstPoint) {
29947c6ae99SBarry Smith       subgroupRanks[subgroupSize++] = p;
3007d3de750SJacob Faibussowitsch       ierr = PetscInfo(da, "%D ", p);CHKERRQ(ierr);
30147c6ae99SBarry Smith     }
30247c6ae99SBarry Smith   }
30347c6ae99SBarry Smith   ierr = PetscInfo(da, "\n");CHKERRQ(ierr);
304ffc4695bSBarry Smith   ierr = MPI_Comm_group(comm, &group);CHKERRMPI(ierr);
305ffc4695bSBarry Smith   ierr = MPI_Group_incl(group, subgroupSize, subgroupRanks, &subgroup);CHKERRMPI(ierr);
306ffc4695bSBarry Smith   ierr = MPI_Comm_create(comm, subgroup, subcomm);CHKERRMPI(ierr);
307ffc4695bSBarry Smith   ierr = MPI_Group_free(&subgroup);CHKERRMPI(ierr);
308ffc4695bSBarry Smith   ierr = MPI_Group_free(&group);CHKERRMPI(ierr);
30947c6ae99SBarry Smith   ierr = PetscFree2(firstPoints, subgroupRanks);CHKERRQ(ierr);
31047c6ae99SBarry Smith   PetscFunctionReturn(0);
31147c6ae99SBarry Smith }
312