1 2 /* 3 Code for manipulating distributed regular arrays in parallel. 4 */ 5 6 #include "private/daimpl.h" /*I "petscdm.h" I*/ 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "DMDAGetProcessorSubset" 10 /*@C 11 DMDAGetProcessorSubset - Returns a communicator consisting only of the 12 processors in a DMDA that own a particular global x, y, or z grid point 13 (corresponding to a logical plane in a 3D grid or a line in a 2D grid). 14 15 Collective on DMDA 16 17 Input Parameters: 18 + da - the distributed array 19 . dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z 20 - gp - global grid point number in this direction 21 22 Output Parameters: 23 . comm - new communicator 24 25 Level: advanced 26 27 Notes: 28 All processors that share the DMDA must call this with the same gp value 29 30 This routine is particularly useful to compute boundary conditions 31 or other application-specific calculations that require manipulating 32 sets of data throughout a logical plane of grid points. 33 34 .keywords: distributed array, get, processor subset 35 @*/ 36 PetscErrorCode DMDAGetProcessorSubset(DM da,DMDADirection dir,PetscInt gp,MPI_Comm *comm) 37 { 38 MPI_Group group,subgroup; 39 PetscErrorCode ierr; 40 PetscInt i,ict,flag,*owners,xs,xm,ys,ym,zs,zm; 41 PetscMPIInt size,*ranks = PETSC_NULL; 42 DM_DA *dd = (DM_DA*)da->data; 43 44 PetscFunctionBegin; 45 PetscValidHeaderSpecific(da,DM_CLASSID,1); 46 flag = 0; 47 ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 48 ierr = MPI_Comm_size(((PetscObject)da)->comm,&size);CHKERRQ(ierr); 49 if (dir == DMDA_Z) { 50 if (dd->dim < 3) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Z invalid for DMDA dim < 3"); 51 if (gp < 0 || gp > dd->P) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point"); 52 if (gp >= zs && gp < zs+zm) flag = 1; 53 } else if (dir == DMDA_Y) { 54 if (dd->dim == 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Y invalid for DMDA dim = 1"); 55 if (gp < 0 || gp > dd->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point"); 56 if (gp >= ys && gp < ys+ym) flag = 1; 57 } else if (dir == DMDA_X) { 58 if (gp < 0 || gp > dd->M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point"); 59 if (gp >= xs && gp < xs+xm) flag = 1; 60 } else SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction"); 61 62 ierr = PetscMalloc2(size,PetscInt,&owners,size,PetscMPIInt,&ranks);CHKERRQ(ierr); 63 ierr = MPI_Allgather(&flag,1,MPIU_INT,owners,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr); 64 ict = 0; 65 ierr = PetscInfo2(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",dd->dim,(int)dir);CHKERRQ(ierr); 66 for (i=0; i<size; i++) { 67 if (owners[i]) { 68 ranks[ict] = i; ict++; 69 ierr = PetscInfo1(da,"%D ",i);CHKERRQ(ierr); 70 } 71 } 72 ierr = PetscInfo(da,"\n");CHKERRQ(ierr); 73 ierr = MPI_Comm_group(((PetscObject)da)->comm,&group);CHKERRQ(ierr); 74 ierr = MPI_Group_incl(group,ict,ranks,&subgroup);CHKERRQ(ierr); 75 ierr = MPI_Comm_create(((PetscObject)da)->comm,subgroup,comm);CHKERRQ(ierr); 76 ierr = MPI_Group_free(&subgroup);CHKERRQ(ierr); 77 ierr = MPI_Group_free(&group);CHKERRQ(ierr); 78 ierr = PetscFree2(owners,ranks);CHKERRQ(ierr); 79 PetscFunctionReturn(0); 80 } 81 82 #undef __FUNCT__ 83 #define __FUNCT__ "DMDAGetProcessorSubsets" 84 /*@C 85 DMDAGetProcessorSubsets - Returns communicators consisting only of the 86 processors in a DMDA adjacent in a particular dimension, 87 corresponding to a logical plane in a 3D grid or a line in a 2D grid. 88 89 Collective on DMDA 90 91 Input Parameters: 92 + da - the distributed array 93 - dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z 94 95 Output Parameters: 96 . subcomm - new communicator 97 98 Level: advanced 99 100 Notes: 101 This routine is useful for distributing one-dimensional data in a tensor product grid. 102 103 .keywords: distributed array, get, processor subset 104 @*/ 105 PetscErrorCode DMDAGetProcessorSubsets(DM da, DMDADirection dir, MPI_Comm *subcomm) 106 { 107 MPI_Comm comm; 108 MPI_Group group, subgroup; 109 PetscInt subgroupSize = 0; 110 PetscInt *firstPoints; 111 PetscMPIInt size, *subgroupRanks = PETSC_NULL; 112 PetscInt xs, xm, ys, ym, zs, zm, firstPoint, p; 113 PetscErrorCode ierr; 114 DM_DA *dd = (DM_DA*)da->data; 115 116 PetscFunctionBegin; 117 PetscValidHeaderSpecific(da, DM_CLASSID, 1); 118 comm = ((PetscObject) da)->comm; 119 ierr = DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm);CHKERRQ(ierr); 120 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 121 if (dir == DMDA_Z) { 122 if (dd->dim < 3) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Z invalid for DMDA dim < 3"); 123 firstPoint = zs; 124 } else if (dir == DMDA_Y) { 125 if (dd->dim == 1) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Y invalid for DMDA dim = 1"); 126 firstPoint = ys; 127 } else if (dir == DMDA_X) { 128 firstPoint = xs; 129 } else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction"); 130 131 ierr = PetscMalloc2(size, PetscInt, &firstPoints, size, PetscMPIInt, &subgroupRanks);CHKERRQ(ierr); 132 ierr = MPI_Allgather(&firstPoint, 1, MPIU_INT, firstPoints, 1, MPIU_INT, comm);CHKERRQ(ierr); 133 ierr = PetscInfo2(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",dd->dim,(int)dir);CHKERRQ(ierr); 134 for(p = 0; p < size; ++p) { 135 if (firstPoints[p] == firstPoint) { 136 subgroupRanks[subgroupSize++] = p; 137 ierr = PetscInfo1(da, "%D ", p);CHKERRQ(ierr); 138 } 139 } 140 ierr = PetscInfo(da, "\n");CHKERRQ(ierr); 141 ierr = MPI_Comm_group(comm, &group);CHKERRQ(ierr); 142 ierr = MPI_Group_incl(group, subgroupSize, subgroupRanks, &subgroup);CHKERRQ(ierr); 143 ierr = MPI_Comm_create(comm, subgroup, subcomm);CHKERRQ(ierr); 144 ierr = MPI_Group_free(&subgroup);CHKERRQ(ierr); 145 ierr = MPI_Group_free(&group);CHKERRQ(ierr); 146 ierr = PetscFree2(firstPoints, subgroupRanks);CHKERRQ(ierr); 147 PetscFunctionReturn(0); 148 } 149