1 2 /* 3 Code for manipulating distributed regular arrays in parallel. 4 */ 5 6 #include <petsc-private/dmdaimpl.h> /*I "petscdmda.h" I*/ 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "DMDAGetLogicalCoordinate" 10 /*@C 11 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 12 13 Collective on DMDA 14 15 Input Parameters: 16 + da - the distributed array 17 - x,y,z - the physical coordinates 18 19 Output Parameters: 20 + II, JJ, KK - the logical coordinate (-1 on processes that do not contain that point) 21 - X, Y, Z, - (optional) the coordinates of the located grid point 22 23 Level: advanced 24 25 Notes: 26 All processors that share the DMDA must call this with the same coordinate value 27 28 .keywords: distributed array, get, processor subset 29 @*/ 30 PetscErrorCode DMDAGetLogicalCoordinate(DM da,PetscScalar x,PetscScalar y,PetscScalar z,PetscInt *II,PetscInt *JJ,PetscInt *KK,PetscScalar *X,PetscScalar *Y,PetscScalar *Z) 31 { 32 DM_DA *dd = (DM_DA*)da->data; 33 PetscErrorCode ierr; 34 Vec coors; 35 DM dacoors; 36 DMDACoor2d **c; 37 PetscInt i,j,xs,xm,ys,ym; 38 PetscReal d,D = PETSC_MAX_REAL,Dv; 39 PetscMPIInt rank,root; 40 41 PetscFunctionBegin; 42 if (dd->dim == 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot get point from 1d DMDA"); 43 if (dd->dim == 3) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot get point from 3d DMDA"); 44 45 *II = -1; 46 *JJ = -1; 47 48 ierr = DMGetCoordinateDM(da,&dacoors);CHKERRQ(ierr); 49 ierr = DMDAGetCorners(dacoors,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 50 ierr = DMGetCoordinates(da,&coors);CHKERRQ(ierr); 51 ierr = DMDAVecGetArray(dacoors,coors,&c);CHKERRQ(ierr); 52 for (j=ys; j<ys+ym; j++) { 53 for (i=xs; i<xs+xm; i++) { 54 d = PetscSqrtReal(PetscRealPart( (c[j][i].x - x)*(c[j][i].x - x) + (c[j][i].y - y)*(c[j][i].y - y) )); 55 if (d < D) { 56 D = d; 57 *II = i; 58 *JJ = j; 59 } 60 } 61 } 62 ierr = MPI_Allreduce(&D,&Dv,1,MPIU_REAL,MPI_MIN,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 63 if (D != Dv) { 64 *II = -1; 65 *JJ = -1; 66 rank = 0; 67 } else { 68 *X = c[*JJ][*II].x; 69 *Y = c[*JJ][*II].y; 70 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr); 71 rank++; 72 } 73 ierr = MPI_Allreduce(&rank,&root,1,MPI_INT,MPI_SUM,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 74 root--; 75 ierr = MPI_Bcast(X,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 76 ierr = MPI_Bcast(Y,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 77 ierr = DMDAVecRestoreArray(dacoors,coors,&c);CHKERRQ(ierr); 78 PetscFunctionReturn(0); 79 } 80 81 #undef __FUNCT__ 82 #define __FUNCT__ "DMDAGetRay" 83 /*@C 84 DMDAGetRay - Returns a vector on process zero that contains a row or column of the values in a DMDA vector 85 86 Collective on DMDA 87 88 Input Parameters: 89 + da - the distributed array 90 . vec - the vector 91 . dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z 92 - gp - global grid point number in this direction 93 94 Output Parameters: 95 + newvec - the new vector that can hold the values (size zero on all processes except process 0) 96 - scatter - the VecScatter that will map from the original vector to the slice 97 98 Level: advanced 99 100 Notes: 101 All processors that share the DMDA must call this with the same gp value 102 103 .keywords: distributed array, get, processor subset 104 @*/ 105 PetscErrorCode DMDAGetRay(DM da,DMDADirection dir,PetscInt gp,Vec *newvec,VecScatter *scatter) 106 { 107 PetscMPIInt rank; 108 DM_DA *dd = (DM_DA*)da->data; 109 PetscErrorCode ierr; 110 IS is; 111 AO ao; 112 Vec vec; 113 PetscInt *indices,i,j; 114 115 PetscFunctionBegin; 116 if (dd->dim == 3) SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_SUP, "Cannot get slice from 3d DMDA"); 117 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) da), &rank);CHKERRQ(ierr); 118 ierr = DMDAGetAO(da, &ao);CHKERRQ(ierr); 119 if (!rank) { 120 if (dd->dim == 1) { 121 if (dir == DMDA_X) { 122 ierr = PetscMalloc1(dd->w, &indices);CHKERRQ(ierr); 123 indices[0] = dd->w*gp; 124 for (i = 1; i < dd->w; ++i) indices[i] = indices[i-1] + 1; 125 ierr = AOApplicationToPetsc(ao, dd->w, indices);CHKERRQ(ierr); 126 ierr = VecCreate(PETSC_COMM_SELF, newvec);CHKERRQ(ierr); 127 ierr = VecSetBlockSize(*newvec, dd->w);CHKERRQ(ierr); 128 ierr = VecSetSizes(*newvec, dd->w, PETSC_DETERMINE);CHKERRQ(ierr); 129 ierr = VecSetType(*newvec, VECSEQ);CHKERRQ(ierr); 130 ierr = ISCreateGeneral(PETSC_COMM_SELF, dd->w, indices, PETSC_OWN_POINTER, &is);CHKERRQ(ierr); 131 } else if (dir == DMDA_Y) SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_SUP, "Cannot get Y slice from 1d DMDA"); 132 else SETERRQ(PetscObjectComm((PetscObject) da), PETSC_ERR_ARG_OUTOFRANGE, "Unknown DMDADirection"); 133 } else { 134 if (dir == DMDA_Y) { 135 ierr = PetscMalloc1(dd->w*dd->M,&indices);CHKERRQ(ierr); 136 indices[0] = gp*dd->M*dd->w; 137 for (i=1; i<dd->M*dd->w; i++) indices[i] = indices[i-1] + 1; 138 139 ierr = AOApplicationToPetsc(ao,dd->M*dd->w,indices);CHKERRQ(ierr); 140 ierr = VecCreate(PETSC_COMM_SELF,newvec);CHKERRQ(ierr); 141 ierr = VecSetBlockSize(*newvec,dd->w);CHKERRQ(ierr); 142 ierr = VecSetSizes(*newvec,dd->M*dd->w,PETSC_DETERMINE);CHKERRQ(ierr); 143 ierr = VecSetType(*newvec,VECSEQ);CHKERRQ(ierr); 144 ierr = ISCreateGeneral(PETSC_COMM_SELF,dd->w*dd->M,indices,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 145 } else if (dir == DMDA_X) { 146 ierr = PetscMalloc1(dd->w*dd->N,&indices);CHKERRQ(ierr); 147 indices[0] = dd->w*gp; 148 for (j=1; j<dd->w; j++) indices[j] = indices[j-1] + 1; 149 for (i=1; i<dd->N; i++) { 150 indices[i*dd->w] = indices[i*dd->w-1] + dd->w*dd->M - dd->w + 1; 151 for (j=1; j<dd->w; j++) indices[i*dd->w + j] = indices[i*dd->w + j - 1] + 1; 152 } 153 ierr = AOApplicationToPetsc(ao,dd->w*dd->N,indices);CHKERRQ(ierr); 154 ierr = VecCreate(PETSC_COMM_SELF,newvec);CHKERRQ(ierr); 155 ierr = VecSetBlockSize(*newvec,dd->w);CHKERRQ(ierr); 156 ierr = VecSetSizes(*newvec,dd->N*dd->w,PETSC_DETERMINE);CHKERRQ(ierr); 157 ierr = VecSetType(*newvec,VECSEQ);CHKERRQ(ierr); 158 ierr = ISCreateGeneral(PETSC_COMM_SELF,dd->w*dd->N,indices,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 159 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown DMDADirection"); 160 } 161 } else { 162 ierr = VecCreateSeq(PETSC_COMM_SELF, 0, newvec);CHKERRQ(ierr); 163 ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, 0, PETSC_COPY_VALUES, &is);CHKERRQ(ierr); 164 } 165 ierr = DMGetGlobalVector(da, &vec);CHKERRQ(ierr); 166 ierr = VecScatterCreate(vec, is, *newvec, NULL, scatter);CHKERRQ(ierr); 167 ierr = DMRestoreGlobalVector(da, &vec);CHKERRQ(ierr); 168 ierr = ISDestroy(&is);CHKERRQ(ierr); 169 PetscFunctionReturn(0); 170 } 171 172 #undef __FUNCT__ 173 #define __FUNCT__ "DMDAGetProcessorSubset" 174 /*@C 175 DMDAGetProcessorSubset - Returns a communicator consisting only of the 176 processors in a DMDA that own a particular global x, y, or z grid point 177 (corresponding to a logical plane in a 3D grid or a line in a 2D grid). 178 179 Collective on DMDA 180 181 Input Parameters: 182 + da - the distributed array 183 . dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z 184 - gp - global grid point number in this direction 185 186 Output Parameters: 187 . comm - new communicator 188 189 Level: advanced 190 191 Notes: 192 All processors that share the DMDA must call this with the same gp value 193 194 This routine is particularly useful to compute boundary conditions 195 or other application-specific calculations that require manipulating 196 sets of data throughout a logical plane of grid points. 197 198 .keywords: distributed array, get, processor subset 199 @*/ 200 PetscErrorCode DMDAGetProcessorSubset(DM da,DMDADirection dir,PetscInt gp,MPI_Comm *comm) 201 { 202 MPI_Group group,subgroup; 203 PetscErrorCode ierr; 204 PetscInt i,ict,flag,*owners,xs,xm,ys,ym,zs,zm; 205 PetscMPIInt size,*ranks = NULL; 206 DM_DA *dd = (DM_DA*)da->data; 207 208 PetscFunctionBegin; 209 PetscValidHeaderSpecific(da,DM_CLASSID,1); 210 flag = 0; 211 ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 212 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr); 213 if (dir == DMDA_Z) { 214 if (dd->dim < 3) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Z invalid for DMDA dim < 3"); 215 if (gp < 0 || gp > dd->P) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point"); 216 if (gp >= zs && gp < zs+zm) flag = 1; 217 } else if (dir == DMDA_Y) { 218 if (dd->dim == 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Y invalid for DMDA dim = 1"); 219 if (gp < 0 || gp > dd->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point"); 220 if (gp >= ys && gp < ys+ym) flag = 1; 221 } else if (dir == DMDA_X) { 222 if (gp < 0 || gp > dd->M) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid grid point"); 223 if (gp >= xs && gp < xs+xm) flag = 1; 224 } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction"); 225 226 ierr = PetscMalloc2(size,&owners,size,&ranks);CHKERRQ(ierr); 227 ierr = MPI_Allgather(&flag,1,MPIU_INT,owners,1,MPIU_INT,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); 228 ict = 0; 229 ierr = PetscInfo2(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",dd->dim,(int)dir);CHKERRQ(ierr); 230 for (i=0; i<size; i++) { 231 if (owners[i]) { 232 ranks[ict] = i; ict++; 233 ierr = PetscInfo1(da,"%D ",i);CHKERRQ(ierr); 234 } 235 } 236 ierr = PetscInfo(da,"\n");CHKERRQ(ierr); 237 ierr = MPI_Comm_group(PetscObjectComm((PetscObject)da),&group);CHKERRQ(ierr); 238 ierr = MPI_Group_incl(group,ict,ranks,&subgroup);CHKERRQ(ierr); 239 ierr = MPI_Comm_create(PetscObjectComm((PetscObject)da),subgroup,comm);CHKERRQ(ierr); 240 ierr = MPI_Group_free(&subgroup);CHKERRQ(ierr); 241 ierr = MPI_Group_free(&group);CHKERRQ(ierr); 242 ierr = PetscFree2(owners,ranks);CHKERRQ(ierr); 243 PetscFunctionReturn(0); 244 } 245 246 #undef __FUNCT__ 247 #define __FUNCT__ "DMDAGetProcessorSubsets" 248 /*@C 249 DMDAGetProcessorSubsets - Returns communicators consisting only of the 250 processors in a DMDA adjacent in a particular dimension, 251 corresponding to a logical plane in a 3D grid or a line in a 2D grid. 252 253 Collective on DMDA 254 255 Input Parameters: 256 + da - the distributed array 257 - dir - Cartesian direction, either DMDA_X, DMDA_Y, or DMDA_Z 258 259 Output Parameters: 260 . subcomm - new communicator 261 262 Level: advanced 263 264 Notes: 265 This routine is useful for distributing one-dimensional data in a tensor product grid. 266 267 .keywords: distributed array, get, processor subset 268 @*/ 269 PetscErrorCode DMDAGetProcessorSubsets(DM da, DMDADirection dir, MPI_Comm *subcomm) 270 { 271 MPI_Comm comm; 272 MPI_Group group, subgroup; 273 PetscInt subgroupSize = 0; 274 PetscInt *firstPoints; 275 PetscMPIInt size, *subgroupRanks = NULL; 276 PetscInt xs, xm, ys, ym, zs, zm, firstPoint, p; 277 PetscErrorCode ierr; 278 DM_DA *dd = (DM_DA*)da->data; 279 280 PetscFunctionBegin; 281 PetscValidHeaderSpecific(da, DM_CLASSID, 1); 282 ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 283 ierr = DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm);CHKERRQ(ierr); 284 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 285 if (dir == DMDA_Z) { 286 if (dd->dim < 3) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Z invalid for DMDA dim < 3"); 287 firstPoint = zs; 288 } else if (dir == DMDA_Y) { 289 if (dd->dim == 1) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"DMDA_Y invalid for DMDA dim = 1"); 290 firstPoint = ys; 291 } else if (dir == DMDA_X) { 292 firstPoint = xs; 293 } else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid direction"); 294 295 ierr = PetscMalloc2(size, &firstPoints, size, &subgroupRanks);CHKERRQ(ierr); 296 ierr = MPI_Allgather(&firstPoint, 1, MPIU_INT, firstPoints, 1, MPIU_INT, comm);CHKERRQ(ierr); 297 ierr = PetscInfo2(da,"DMDAGetProcessorSubset: dim=%D, direction=%d, procs: ",dd->dim,(int)dir);CHKERRQ(ierr); 298 for (p = 0; p < size; ++p) { 299 if (firstPoints[p] == firstPoint) { 300 subgroupRanks[subgroupSize++] = p; 301 ierr = PetscInfo1(da, "%D ", p);CHKERRQ(ierr); 302 } 303 } 304 ierr = PetscInfo(da, "\n");CHKERRQ(ierr); 305 ierr = MPI_Comm_group(comm, &group);CHKERRQ(ierr); 306 ierr = MPI_Group_incl(group, subgroupSize, subgroupRanks, &subgroup);CHKERRQ(ierr); 307 ierr = MPI_Comm_create(comm, subgroup, subcomm);CHKERRQ(ierr); 308 ierr = MPI_Group_free(&subgroup);CHKERRQ(ierr); 309 ierr = MPI_Group_free(&group);CHKERRQ(ierr); 310 ierr = PetscFree2(firstPoints, subgroupRanks);CHKERRQ(ierr); 311 PetscFunctionReturn(0); 312 } 313