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 Vec coors; 360a0169f3SBarry Smith DM dacoors; 370a0169f3SBarry Smith DMDACoor2d **c; 380a0169f3SBarry Smith PetscInt i,j,xs,xm,ys,ym; 39e97f6855SBarry Smith PetscReal d,D = PETSC_MAX_REAL,Dv; 40e97f6855SBarry Smith PetscMPIInt rank,root; 410a0169f3SBarry Smith 420a0169f3SBarry Smith PetscFunctionBegin; 4308401ef6SPierre Jolivet PetscCheck(da->dim != 1,PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot get point from 1d DMDA"); 4408401ef6SPierre Jolivet PetscCheck(da->dim != 3,PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot get point from 3d DMDA"); 450a0169f3SBarry Smith 4633c53e3fSSatish Balay *II = -1; 47e97f6855SBarry Smith *JJ = -1; 480a0169f3SBarry Smith 499566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da,&dacoors)); 509566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dacoors,&xs,&ys,NULL,&xm,&ym,NULL)); 519566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da,&coors)); 529566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(dacoors,coors,&c)); 530a0169f3SBarry Smith for (j=ys; j<ys+ym; j++) { 540a0169f3SBarry Smith for (i=xs; i<xs+xm; i++) { 550a0169f3SBarry Smith d = PetscSqrtReal(PetscRealPart((c[j][i].x - x)*(c[j][i].x - x) + (c[j][i].y - y)*(c[j][i].y - y))); 560a0169f3SBarry Smith if (d < D) { 570a0169f3SBarry Smith D = d; 5833c53e3fSSatish Balay *II = i; 5933c53e3fSSatish Balay *JJ = j; 600a0169f3SBarry Smith } 610a0169f3SBarry Smith } 620a0169f3SBarry Smith } 631c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&D,&Dv,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)da))); 64e97f6855SBarry Smith if (D != Dv) { 65e97f6855SBarry Smith *II = -1; 66e97f6855SBarry Smith *JJ = -1; 67e97f6855SBarry Smith rank = 0; 68e97f6855SBarry Smith } else { 69e97f6855SBarry Smith *X = c[*JJ][*II].x; 70e97f6855SBarry Smith *Y = c[*JJ][*II].y; 719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank)); 72e97f6855SBarry Smith rank++; 73e97f6855SBarry Smith } 741c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&rank,&root,1,MPI_INT,MPI_SUM,PetscObjectComm((PetscObject)da))); 75e97f6855SBarry Smith root--; 769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(X,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)da))); 779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(Y,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)da))); 789566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(dacoors,coors,&c)); 790a0169f3SBarry Smith PetscFunctionReturn(0); 800a0169f3SBarry Smith } 810a0169f3SBarry Smith 8269e49704SSatish Balay /*@ 83db05f41bSBarry Smith DMDAGetRay - Returns a vector on process zero that contains a row or column of the values in a DMDA vector 84db05f41bSBarry Smith 85db05f41bSBarry Smith Collective on DMDA 86db05f41bSBarry Smith 87db05f41bSBarry Smith Input Parameters: 88db05f41bSBarry Smith + da - the distributed array 893ee9839eSMatthew G. Knepley . dir - Cartesian direction, either DM_X, DM_Y, or DM_Z 90db05f41bSBarry Smith - gp - global grid point number in this direction 91db05f41bSBarry Smith 92db05f41bSBarry Smith Output Parameters: 93db05f41bSBarry Smith + newvec - the new vector that can hold the values (size zero on all processes except process 0) 94db05f41bSBarry Smith - scatter - the VecScatter that will map from the original vector to the slice 95db05f41bSBarry Smith 96db05f41bSBarry Smith Level: advanced 97db05f41bSBarry Smith 98db05f41bSBarry Smith Notes: 99db05f41bSBarry Smith All processors that share the DMDA must call this with the same gp value 100db05f41bSBarry Smith 101db05f41bSBarry Smith @*/ 1023ee9839eSMatthew G. Knepley PetscErrorCode DMDAGetRay(DM da,DMDirection dir,PetscInt gp,Vec *newvec,VecScatter *scatter) 103db05f41bSBarry Smith { 104db05f41bSBarry Smith PetscMPIInt rank; 105db05f41bSBarry Smith DM_DA *dd = (DM_DA*)da->data; 106db05f41bSBarry Smith IS is; 107db05f41bSBarry Smith AO ao; 108db05f41bSBarry Smith Vec vec; 109d1212d36SBarry Smith PetscInt *indices,i,j; 110db05f41bSBarry Smith 111db05f41bSBarry Smith PetscFunctionBegin; 11208401ef6SPierre Jolivet PetscCheck(da->dim != 3,PetscObjectComm((PetscObject) da), PETSC_ERR_SUP, "Cannot get slice from 3d DMDA"); 1139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) da), &rank)); 1149566063dSJacob Faibussowitsch PetscCall(DMDAGetAO(da, &ao)); 115dd400576SPatrick Sanan if (rank == 0) { 116c73cfb54SMatthew G. Knepley if (da->dim == 1) { 1173ee9839eSMatthew G. Knepley if (dir == DM_X) { 1189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->w, &indices)); 119c2001ae3SMatthew G. Knepley indices[0] = dd->w*gp; 120eca19243SMatthew G. Knepley for (i = 1; i < dd->w; ++i) indices[i] = indices[i-1] + 1; 1219566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, dd->w, indices)); 1229566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, newvec)); 1239566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*newvec, dd->w)); 1249566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*newvec, dd->w, PETSC_DETERMINE)); 1259566063dSJacob Faibussowitsch PetscCall(VecSetType(*newvec, VECSEQ)); 1269566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, dd->w, indices, PETSC_OWN_POINTER, &is)); 127*f7d195e4SLawrence Mitchell } else { 128*f7d195e4SLawrence Mitchell PetscCheck(dir != DM_Y,PetscObjectComm((PetscObject) da), PETSC_ERR_SUP, "Cannot get Y slice from 1d DMDA"); 129*f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_ARG_OUTOFRANGE, "Unknown DMDirection"); 130*f7d195e4SLawrence Mitchell } 131c2001ae3SMatthew G. Knepley } else { 1323ee9839eSMatthew G. Knepley if (dir == DM_Y) { 1339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->w*dd->M,&indices)); 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 1379566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao,dd->M*dd->w,indices)); 1389566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,newvec)); 1399566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*newvec,dd->w)); 1409566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*newvec,dd->M*dd->w,PETSC_DETERMINE)); 1419566063dSJacob Faibussowitsch PetscCall(VecSetType(*newvec,VECSEQ)); 1429566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,dd->w*dd->M,indices,PETSC_OWN_POINTER,&is)); 1433ee9839eSMatthew G. Knepley } else if (dir == DM_X) { 1449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->w*dd->N,&indices)); 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 } 1519566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao,dd->w*dd->N,indices)); 1529566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,newvec)); 1539566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*newvec,dd->w)); 1549566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*newvec,dd->N*dd->w,PETSC_DETERMINE)); 1559566063dSJacob Faibussowitsch PetscCall(VecSetType(*newvec,VECSEQ)); 1569566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,dd->w*dd->N,indices,PETSC_OWN_POINTER,&is)); 1573ee9839eSMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown DMDirection"); 158c2001ae3SMatthew G. Knepley } 159db05f41bSBarry Smith } else { 1609566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, newvec)); 1619566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is)); 162db05f41bSBarry Smith } 1639566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &vec)); 1649566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vec, is, *newvec, NULL, scatter)); 1659566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &vec)); 1669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 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 PetscInt i,ict,flag,*owners,xs,xm,ys,ym,zs,zm; 2030298fd71SBarry Smith PetscMPIInt size,*ranks = NULL; 20447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 20547c6ae99SBarry Smith 20647c6ae99SBarry Smith PetscFunctionBegin; 207a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 20847c6ae99SBarry Smith flag = 0; 2099566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 2109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da),&size)); 2113ee9839eSMatthew G. Knepley if (dir == DM_Z) { 21208401ef6SPierre Jolivet PetscCheck(da->dim >= 3,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"DM_Z invalid for DMDA dim < 3"); 2131dca8a05SBarry Smith PetscCheck(gp >= 0 && gp <= dd->P,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point"); 21447c6ae99SBarry Smith if (gp >= zs && gp < zs+zm) flag = 1; 2153ee9839eSMatthew G. Knepley } else if (dir == DM_Y) { 21608401ef6SPierre Jolivet PetscCheck(da->dim != 1,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"DM_Y invalid for DMDA dim = 1"); 2171dca8a05SBarry Smith PetscCheck(gp >= 0 && gp <= dd->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point"); 21847c6ae99SBarry Smith if (gp >= ys && gp < ys+ym) flag = 1; 2193ee9839eSMatthew G. Knepley } else if (dir == DM_X) { 2201dca8a05SBarry Smith PetscCheck(gp >= 0 && gp <= dd->M,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point"); 22147c6ae99SBarry Smith if (gp >= xs && gp < xs+xm) flag = 1; 222ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction"); 22347c6ae99SBarry Smith 2249566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size,&owners,size,&ranks)); 2259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&flag,1,MPIU_INT,owners,1,MPIU_INT,PetscObjectComm((PetscObject)da))); 22647c6ae99SBarry Smith ict = 0; 22763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(da,"DMDAGetProcessorSubset: dim=%" PetscInt_FMT ", direction=%d, procs: ",da->dim,(int)dir)); 22847c6ae99SBarry Smith for (i=0; i<size; i++) { 22947c6ae99SBarry Smith if (owners[i]) { 23047c6ae99SBarry Smith ranks[ict] = i; ict++; 23163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(da,"%" PetscInt_FMT " ",i)); 23247c6ae99SBarry Smith } 23347c6ae99SBarry Smith } 2349566063dSJacob Faibussowitsch PetscCall(PetscInfo(da,"\n")); 2359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)da),&group)); 2369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_incl(group,ict,ranks,&subgroup)); 2379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create(PetscObjectComm((PetscObject)da),subgroup,comm)); 2389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&subgroup)); 2399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&group)); 2409566063dSJacob Faibussowitsch PetscCall(PetscFree2(owners,ranks)); 24147c6ae99SBarry Smith PetscFunctionReturn(0); 24247c6ae99SBarry Smith } 24347c6ae99SBarry Smith 24447c6ae99SBarry Smith /*@C 245aa219208SBarry Smith DMDAGetProcessorSubsets - Returns communicators consisting only of the 246aa219208SBarry Smith processors in a DMDA adjacent in a particular dimension, 24747c6ae99SBarry Smith corresponding to a logical plane in a 3D grid or a line in a 2D grid. 24847c6ae99SBarry Smith 249d083f849SBarry Smith Collective on da 25047c6ae99SBarry Smith 25147c6ae99SBarry Smith Input Parameters: 25247c6ae99SBarry Smith + da - the distributed array 2533ee9839eSMatthew G. Knepley - dir - Cartesian direction, either DM_X, DM_Y, or DM_Z 25447c6ae99SBarry Smith 25597bb3fdcSJose E. Roman Output Parameter: 25647c6ae99SBarry Smith . subcomm - new communicator 25747c6ae99SBarry Smith 25847c6ae99SBarry Smith Level: advanced 25947c6ae99SBarry Smith 26047c6ae99SBarry Smith Notes: 26147c6ae99SBarry Smith This routine is useful for distributing one-dimensional data in a tensor product grid. 26247c6ae99SBarry Smith 263bf3faf6aSSatish Balay After use, comm should be freed with MPI_Comm_free() 264bf3faf6aSSatish Balay 265f5f57ec0SBarry Smith Not supported from Fortran 266f5f57ec0SBarry Smith 26747c6ae99SBarry Smith @*/ 2683ee9839eSMatthew G. Knepley PetscErrorCode DMDAGetProcessorSubsets(DM da, DMDirection dir, MPI_Comm *subcomm) 26947c6ae99SBarry Smith { 27047c6ae99SBarry Smith MPI_Comm comm; 27147c6ae99SBarry Smith MPI_Group group, subgroup; 27247c6ae99SBarry Smith PetscInt subgroupSize = 0; 27347c6ae99SBarry Smith PetscInt *firstPoints; 2740298fd71SBarry Smith PetscMPIInt size, *subgroupRanks = NULL; 27547c6ae99SBarry Smith PetscInt xs, xm, ys, ym, zs, zm, firstPoint, p; 27647c6ae99SBarry Smith 27747c6ae99SBarry Smith PetscFunctionBegin; 278a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 2799566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da,&comm)); 2809566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 2819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2823ee9839eSMatthew G. Knepley if (dir == DM_Z) { 28308401ef6SPierre Jolivet PetscCheck(da->dim >= 3,comm,PETSC_ERR_ARG_OUTOFRANGE,"DM_Z invalid for DMDA dim < 3"); 28447c6ae99SBarry Smith firstPoint = zs; 2853ee9839eSMatthew G. Knepley } else if (dir == DM_Y) { 28608401ef6SPierre Jolivet PetscCheck(da->dim != 1,comm,PETSC_ERR_ARG_OUTOFRANGE,"DM_Y invalid for DMDA dim = 1"); 28747c6ae99SBarry Smith firstPoint = ys; 2883ee9839eSMatthew G. Knepley } else if (dir == DM_X) { 28947c6ae99SBarry Smith firstPoint = xs; 29047c6ae99SBarry Smith } else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction"); 29147c6ae99SBarry Smith 2929566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &firstPoints, size, &subgroupRanks)); 2939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&firstPoint, 1, MPIU_INT, firstPoints, 1, MPIU_INT, comm)); 29463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(da,"DMDAGetProcessorSubset: dim=%" PetscInt_FMT ", direction=%d, procs: ",da->dim,(int)dir)); 29547c6ae99SBarry Smith for (p = 0; p < size; ++p) { 29647c6ae99SBarry Smith if (firstPoints[p] == firstPoint) { 29747c6ae99SBarry Smith subgroupRanks[subgroupSize++] = p; 29863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(da, "%" PetscInt_FMT " ", p)); 29947c6ae99SBarry Smith } 30047c6ae99SBarry Smith } 3019566063dSJacob Faibussowitsch PetscCall(PetscInfo(da, "\n")); 3029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(comm, &group)); 3039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_incl(group, subgroupSize, subgroupRanks, &subgroup)); 3049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create(comm, subgroup, subcomm)); 3059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&subgroup)); 3069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&group)); 3079566063dSJacob Faibussowitsch PetscCall(PetscFree2(firstPoints, subgroupRanks)); 30847c6ae99SBarry Smith PetscFunctionReturn(0); 30947c6ae99SBarry Smith } 310