147c6ae99SBarry Smith /* 247c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel. 347c6ae99SBarry Smith */ 447c6ae99SBarry Smith 5af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 647c6ae99SBarry Smith 769e49704SSatish Balay /*@ 8*dce8aebaSBarry 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` 90a0169f3SBarry Smith 10d083f849SBarry Smith Collective on da 110a0169f3SBarry Smith 120a0169f3SBarry Smith Input Parameters: 130a0169f3SBarry Smith + da - the distributed array 146b867d5aSJose E. Roman . x - the first physical coordinate 156b867d5aSJose E. Roman . y - the second physical coordinate 166b867d5aSJose E. Roman - z - the third physical coordinate 170a0169f3SBarry Smith 180a0169f3SBarry Smith Output Parameters: 196b867d5aSJose E. Roman + II - the first logical coordinate (-1 on processes that do not contain that point) 206b867d5aSJose E. Roman . JJ - the second logical coordinate (-1 on processes that do not contain that point) 216b867d5aSJose E. Roman . KK - the third logical coordinate (-1 on processes that do not contain that point) 226b867d5aSJose E. Roman . X - (optional) the first coordinate of the located grid point 236b867d5aSJose E. Roman . Y - (optional) the second coordinate of the located grid point 246b867d5aSJose E. Roman - Z - (optional) the third coordinate of the located grid point 250a0169f3SBarry Smith 260a0169f3SBarry Smith Level: advanced 270a0169f3SBarry Smith 28*dce8aebaSBarry Smith Note: 29*dce8aebaSBarry Smith All processors that share the `DMDA` must call this with the same coordinate value 300a0169f3SBarry Smith 31*dce8aebaSBarry Smith .seealso: `DM`, `DMDA` 320a0169f3SBarry Smith @*/ 33d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetLogicalCoordinate(DM da, PetscScalar x, PetscScalar y, PetscScalar z, PetscInt *II, PetscInt *JJ, PetscInt *KK, PetscScalar *X, PetscScalar *Y, PetscScalar *Z) 34d71ae5a4SJacob Faibussowitsch { 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 /*@ 83*dce8aebaSBarry Smith DMDAGetRay - Returns a vector on process zero that contains a row or column of the values in a `DMDA` vector 84db05f41bSBarry Smith 85*dce8aebaSBarry Smith Collective on da 86db05f41bSBarry Smith 87db05f41bSBarry Smith Input Parameters: 88db05f41bSBarry Smith + da - the distributed array 89*dce8aebaSBarry Smith . 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) 94*dce8aebaSBarry Smith - scatter - the `VecScatter` that will map from the original vector to the slice 95db05f41bSBarry Smith 96db05f41bSBarry Smith Level: advanced 97db05f41bSBarry Smith 98*dce8aebaSBarry Smith Note: 99*dce8aebaSBarry Smith All processors that share the `DMDA` must call this with the same gp value 100db05f41bSBarry Smith 101*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDirection`, `Vec`, `VecScatter` 102db05f41bSBarry Smith @*/ 103d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetRay(DM da, DMDirection dir, PetscInt gp, Vec *newvec, VecScatter *scatter) 104d71ae5a4SJacob Faibussowitsch { 105db05f41bSBarry Smith PetscMPIInt rank; 106db05f41bSBarry Smith DM_DA *dd = (DM_DA *)da->data; 107db05f41bSBarry Smith IS is; 108db05f41bSBarry Smith AO ao; 109db05f41bSBarry Smith Vec vec; 110d1212d36SBarry Smith PetscInt *indices, i, j; 111db05f41bSBarry Smith 112db05f41bSBarry Smith PetscFunctionBegin; 11308401ef6SPierre Jolivet PetscCheck(da->dim != 3, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot get slice from 3d DMDA"); 1149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)da), &rank)); 1159566063dSJacob Faibussowitsch PetscCall(DMDAGetAO(da, &ao)); 116dd400576SPatrick Sanan if (rank == 0) { 117c73cfb54SMatthew G. Knepley if (da->dim == 1) { 1183ee9839eSMatthew G. Knepley if (dir == DM_X) { 1199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->w, &indices)); 120c2001ae3SMatthew G. Knepley indices[0] = dd->w * gp; 121eca19243SMatthew G. Knepley for (i = 1; i < dd->w; ++i) indices[i] = indices[i - 1] + 1; 1229566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, dd->w, indices)); 1239566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, newvec)); 1249566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*newvec, dd->w)); 1259566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*newvec, dd->w, PETSC_DETERMINE)); 1269566063dSJacob Faibussowitsch PetscCall(VecSetType(*newvec, VECSEQ)); 1279566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, dd->w, indices, PETSC_OWN_POINTER, &is)); 128f7d195e4SLawrence Mitchell } else { 129f7d195e4SLawrence Mitchell PetscCheck(dir != DM_Y, PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot get Y slice from 1d DMDA"); 130f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Unknown DMDirection"); 131f7d195e4SLawrence Mitchell } 132c2001ae3SMatthew G. Knepley } else { 1333ee9839eSMatthew G. Knepley if (dir == DM_Y) { 1349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->w * dd->M, &indices)); 135d1212d36SBarry Smith indices[0] = gp * dd->M * dd->w; 1368865f1eaSKarl Rupp for (i = 1; i < dd->M * dd->w; i++) indices[i] = indices[i - 1] + 1; 1378865f1eaSKarl Rupp 1389566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, dd->M * dd->w, indices)); 1399566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, newvec)); 1409566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*newvec, dd->w)); 1419566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*newvec, dd->M * dd->w, PETSC_DETERMINE)); 1429566063dSJacob Faibussowitsch PetscCall(VecSetType(*newvec, VECSEQ)); 1439566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, dd->w * dd->M, indices, PETSC_OWN_POINTER, &is)); 1443ee9839eSMatthew G. Knepley } else if (dir == DM_X) { 1459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->w * dd->N, &indices)); 146d1212d36SBarry Smith indices[0] = dd->w * gp; 147d1212d36SBarry Smith for (j = 1; j < dd->w; j++) indices[j] = indices[j - 1] + 1; 148d1212d36SBarry Smith for (i = 1; i < dd->N; i++) { 149d1212d36SBarry Smith indices[i * dd->w] = indices[i * dd->w - 1] + dd->w * dd->M - dd->w + 1; 150d1212d36SBarry Smith for (j = 1; j < dd->w; j++) indices[i * dd->w + j] = indices[i * dd->w + j - 1] + 1; 151d1212d36SBarry Smith } 1529566063dSJacob Faibussowitsch PetscCall(AOApplicationToPetsc(ao, dd->w * dd->N, indices)); 1539566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, newvec)); 1549566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*newvec, dd->w)); 1559566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*newvec, dd->N * dd->w, PETSC_DETERMINE)); 1569566063dSJacob Faibussowitsch PetscCall(VecSetType(*newvec, VECSEQ)); 1579566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, dd->w * dd->N, indices, PETSC_OWN_POINTER, &is)); 1583ee9839eSMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown DMDirection"); 159c2001ae3SMatthew G. Knepley } 160db05f41bSBarry Smith } else { 1619566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, newvec)); 1629566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_COPY_VALUES, &is)); 163db05f41bSBarry Smith } 1649566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &vec)); 1659566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vec, is, *newvec, NULL, scatter)); 1669566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &vec)); 1679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 168db05f41bSBarry Smith PetscFunctionReturn(0); 169db05f41bSBarry Smith } 170db05f41bSBarry Smith 17147c6ae99SBarry Smith /*@C 172aa219208SBarry Smith DMDAGetProcessorSubset - Returns a communicator consisting only of the 173*dce8aebaSBarry Smith processors in a `DMDA` that own a particular global x, y, or z grid point 17447c6ae99SBarry Smith (corresponding to a logical plane in a 3D grid or a line in a 2D grid). 17547c6ae99SBarry Smith 176d083f849SBarry Smith Collective on da 17747c6ae99SBarry Smith 17847c6ae99SBarry Smith Input Parameters: 17947c6ae99SBarry Smith + da - the distributed array 180*dce8aebaSBarry Smith . dir - Cartesian direction, either `DM_X`, `DM_Y`, or `DM_Z` 18147c6ae99SBarry Smith - gp - global grid point number in this direction 18247c6ae99SBarry Smith 18397bb3fdcSJose E. Roman Output Parameter: 18447c6ae99SBarry Smith . comm - new communicator 18547c6ae99SBarry Smith 18647c6ae99SBarry Smith Level: advanced 18747c6ae99SBarry Smith 18847c6ae99SBarry Smith Notes: 189*dce8aebaSBarry Smith All processors that share the `DMDA` must call this with the same gp value 19047c6ae99SBarry Smith 191*dce8aebaSBarry Smith After use, comm should be freed with `MPI_Comm_free()` 192bf3faf6aSSatish Balay 19347c6ae99SBarry Smith This routine is particularly useful to compute boundary conditions 19447c6ae99SBarry Smith or other application-specific calculations that require manipulating 19547c6ae99SBarry Smith sets of data throughout a logical plane of grid points. 19647c6ae99SBarry Smith 197*dce8aebaSBarry Smith Fortran Note: 198f5f57ec0SBarry Smith Not supported from Fortran 199f5f57ec0SBarry Smith 200*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDirection` 20147c6ae99SBarry Smith @*/ 202d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetProcessorSubset(DM da, DMDirection dir, PetscInt gp, MPI_Comm *comm) 203d71ae5a4SJacob Faibussowitsch { 20447c6ae99SBarry Smith MPI_Group group, subgroup; 20547c6ae99SBarry Smith PetscInt i, ict, flag, *owners, xs, xm, ys, ym, zs, zm; 2060298fd71SBarry Smith PetscMPIInt size, *ranks = NULL; 20747c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 20847c6ae99SBarry Smith 20947c6ae99SBarry Smith PetscFunctionBegin; 210a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 21147c6ae99SBarry Smith flag = 0; 2129566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 2139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)da), &size)); 2143ee9839eSMatthew G. Knepley if (dir == DM_Z) { 21508401ef6SPierre Jolivet PetscCheck(da->dim >= 3, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "DM_Z invalid for DMDA dim < 3"); 2161dca8a05SBarry Smith PetscCheck(gp >= 0 && gp <= dd->P, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid grid point"); 21747c6ae99SBarry Smith if (gp >= zs && gp < zs + zm) flag = 1; 2183ee9839eSMatthew G. Knepley } else if (dir == DM_Y) { 21908401ef6SPierre Jolivet PetscCheck(da->dim != 1, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "DM_Y invalid for DMDA dim = 1"); 2201dca8a05SBarry Smith PetscCheck(gp >= 0 && gp <= dd->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid grid point"); 22147c6ae99SBarry Smith if (gp >= ys && gp < ys + ym) flag = 1; 2223ee9839eSMatthew G. Knepley } else if (dir == DM_X) { 2231dca8a05SBarry Smith PetscCheck(gp >= 0 && gp <= dd->M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid grid point"); 22447c6ae99SBarry Smith if (gp >= xs && gp < xs + xm) flag = 1; 225ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Invalid direction"); 22647c6ae99SBarry Smith 2279566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &owners, size, &ranks)); 2289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&flag, 1, MPIU_INT, owners, 1, MPIU_INT, PetscObjectComm((PetscObject)da))); 22947c6ae99SBarry Smith ict = 0; 23063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(da, "DMDAGetProcessorSubset: dim=%" PetscInt_FMT ", direction=%d, procs: ", da->dim, (int)dir)); 23147c6ae99SBarry Smith for (i = 0; i < size; i++) { 23247c6ae99SBarry Smith if (owners[i]) { 2339371c9d4SSatish Balay ranks[ict] = i; 2349371c9d4SSatish Balay ict++; 23563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(da, "%" PetscInt_FMT " ", i)); 23647c6ae99SBarry Smith } 23747c6ae99SBarry Smith } 2389566063dSJacob Faibussowitsch PetscCall(PetscInfo(da, "\n")); 2399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)da), &group)); 2409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_incl(group, ict, ranks, &subgroup)); 2419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create(PetscObjectComm((PetscObject)da), subgroup, comm)); 2429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&subgroup)); 2439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&group)); 2449566063dSJacob Faibussowitsch PetscCall(PetscFree2(owners, ranks)); 24547c6ae99SBarry Smith PetscFunctionReturn(0); 24647c6ae99SBarry Smith } 24747c6ae99SBarry Smith 24847c6ae99SBarry Smith /*@C 249aa219208SBarry Smith DMDAGetProcessorSubsets - Returns communicators consisting only of the 250*dce8aebaSBarry Smith processors in a `DMDA` adjacent in a particular dimension, 25147c6ae99SBarry Smith corresponding to a logical plane in a 3D grid or a line in a 2D grid. 25247c6ae99SBarry Smith 253d083f849SBarry Smith Collective on da 25447c6ae99SBarry Smith 25547c6ae99SBarry Smith Input Parameters: 25647c6ae99SBarry Smith + da - the distributed array 257*dce8aebaSBarry Smith - dir - Cartesian direction, either `DM_X`, `DM_Y`, or `DM_Z` 25847c6ae99SBarry Smith 25997bb3fdcSJose E. Roman Output Parameter: 26047c6ae99SBarry Smith . subcomm - new communicator 26147c6ae99SBarry Smith 26247c6ae99SBarry Smith Level: advanced 26347c6ae99SBarry Smith 26447c6ae99SBarry Smith Notes: 26547c6ae99SBarry Smith This routine is useful for distributing one-dimensional data in a tensor product grid. 26647c6ae99SBarry Smith 267*dce8aebaSBarry Smith After use, comm should be freed with` MPI_Comm_free()` 268bf3faf6aSSatish Balay 269*dce8aebaSBarry Smith Fortran Note: 270f5f57ec0SBarry Smith Not supported from Fortran 271f5f57ec0SBarry Smith 272*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDirection` 27347c6ae99SBarry Smith @*/ 274d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetProcessorSubsets(DM da, DMDirection dir, MPI_Comm *subcomm) 275d71ae5a4SJacob Faibussowitsch { 27647c6ae99SBarry Smith MPI_Comm comm; 27747c6ae99SBarry Smith MPI_Group group, subgroup; 27847c6ae99SBarry Smith PetscInt subgroupSize = 0; 27947c6ae99SBarry Smith PetscInt *firstPoints; 2800298fd71SBarry Smith PetscMPIInt size, *subgroupRanks = NULL; 28147c6ae99SBarry Smith PetscInt xs, xm, ys, ym, zs, zm, firstPoint, p; 28247c6ae99SBarry Smith 28347c6ae99SBarry Smith PetscFunctionBegin; 284a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 2859566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 2869566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 2879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2883ee9839eSMatthew G. Knepley if (dir == DM_Z) { 28908401ef6SPierre Jolivet PetscCheck(da->dim >= 3, comm, PETSC_ERR_ARG_OUTOFRANGE, "DM_Z invalid for DMDA dim < 3"); 29047c6ae99SBarry Smith firstPoint = zs; 2913ee9839eSMatthew G. Knepley } else if (dir == DM_Y) { 29208401ef6SPierre Jolivet PetscCheck(da->dim != 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "DM_Y invalid for DMDA dim = 1"); 29347c6ae99SBarry Smith firstPoint = ys; 2943ee9839eSMatthew G. Knepley } else if (dir == DM_X) { 29547c6ae99SBarry Smith firstPoint = xs; 29647c6ae99SBarry Smith } else SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid direction"); 29747c6ae99SBarry Smith 2989566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &firstPoints, size, &subgroupRanks)); 2999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&firstPoint, 1, MPIU_INT, firstPoints, 1, MPIU_INT, comm)); 30063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(da, "DMDAGetProcessorSubset: dim=%" PetscInt_FMT ", direction=%d, procs: ", da->dim, (int)dir)); 30147c6ae99SBarry Smith for (p = 0; p < size; ++p) { 30247c6ae99SBarry Smith if (firstPoints[p] == firstPoint) { 30347c6ae99SBarry Smith subgroupRanks[subgroupSize++] = p; 30463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(da, "%" PetscInt_FMT " ", p)); 30547c6ae99SBarry Smith } 30647c6ae99SBarry Smith } 3079566063dSJacob Faibussowitsch PetscCall(PetscInfo(da, "\n")); 3089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(comm, &group)); 3099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_incl(group, subgroupSize, subgroupRanks, &subgroup)); 3109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create(comm, subgroup, subcomm)); 3119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&subgroup)); 3129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&group)); 3139566063dSJacob Faibussowitsch PetscCall(PetscFree2(firstPoints, subgroupRanks)); 31447c6ae99SBarry Smith PetscFunctionReturn(0); 31547c6ae99SBarry Smith } 316