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