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__ "DMCreateCoordinateDM_DA" 10 PetscErrorCode DMCreateCoordinateDM_DA(DM dm, DM *cdm) 11 { 12 PetscErrorCode ierr; 13 PetscFunctionBegin; 14 ierr = DMDAGetReducedDMDA(dm,dm->dim,cdm);CHKERRQ(ierr); 15 PetscFunctionReturn(0); 16 } 17 18 #undef __FUNCT__ 19 #define __FUNCT__ "DMDASetFieldName" 20 /*@C 21 DMDASetFieldName - Sets the names of individual field components in multicomponent 22 vectors associated with a DMDA. 23 24 Not Collective 25 26 Input Parameters: 27 + da - the distributed array 28 . nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the 29 number of degrees of freedom per node within the DMDA 30 - names - the name of the field (component) 31 32 Level: intermediate 33 34 .keywords: distributed array, get, component name 35 36 .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName() 37 @*/ 38 PetscErrorCode DMDASetFieldName(DM da,PetscInt nf,const char name[]) 39 { 40 PetscErrorCode ierr; 41 DM_DA *dd = (DM_DA*)da->data; 42 43 PetscFunctionBegin; 44 PetscValidHeaderSpecific(da,DM_CLASSID,1); 45 if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf); 46 ierr = PetscFree(dd->fieldname[nf]);CHKERRQ(ierr); 47 ierr = PetscStrallocpy(name,&dd->fieldname[nf]);CHKERRQ(ierr); 48 PetscFunctionReturn(0); 49 } 50 51 #undef __FUNCT__ 52 #define __FUNCT__ "DMDAGetFieldName" 53 /*@C 54 DMDAGetFieldName - Gets the names of individual field components in multicomponent 55 vectors associated with a DMDA. 56 57 Not Collective 58 59 Input Parameter: 60 + da - the distributed array 61 - nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the 62 number of degrees of freedom per node within the DMDA 63 64 Output Parameter: 65 . names - the name of the field (component) 66 67 Level: intermediate 68 69 .keywords: distributed array, get, component name 70 71 .seealso: DMDASetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName() 72 @*/ 73 PetscErrorCode DMDAGetFieldName(DM da,PetscInt nf,const char **name) 74 { 75 DM_DA *dd = (DM_DA*)da->data; 76 77 PetscFunctionBegin; 78 PetscValidHeaderSpecific(da,DM_CLASSID,1); 79 PetscValidPointer(name,3); 80 if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf); 81 *name = dd->fieldname[nf]; 82 PetscFunctionReturn(0); 83 } 84 85 #undef __FUNCT__ 86 #define __FUNCT__ "DMDASetCoordinateName" 87 /*@C 88 DMDASetCoordinateName - Sets the name of the coordinate directions associated with a DMDA, for example "x" or "y" 89 90 Not Collective 91 92 Input Parameters: 93 + dm - the DM 94 . nf - coordinate number for the DMDA (0, 1, ... dim-1), 95 - name - the name of the coordinate 96 97 Level: intermediate 98 99 .keywords: distributed array, get, component name 100 101 .seealso: DMDAGetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName() 102 @*/ 103 PetscErrorCode DMDASetCoordinateName(DM dm,PetscInt nf,const char name[]) 104 { 105 PetscErrorCode ierr; 106 DM_DA *dd = (DM_DA*)dm->data; 107 108 PetscFunctionBegin; 109 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 110 if (nf < 0 || nf >= dm->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf); 111 ierr = PetscFree(dd->coordinatename[nf]);CHKERRQ(ierr); 112 ierr = PetscStrallocpy(name,&dd->coordinatename[nf]);CHKERRQ(ierr); 113 PetscFunctionReturn(0); 114 } 115 116 #undef __FUNCT__ 117 #define __FUNCT__ "DMDAGetCoordinateName" 118 /*@C 119 DMDAGetCoordinateName - Gets the name of a coodinate direction associated with a DMDA. 120 121 Not Collective 122 123 Input Parameter: 124 + dm - the DM 125 - nf - number for the DMDA (0, 1, ... dim-1) 126 127 Output Parameter: 128 . names - the name of the coordinate direction 129 130 Level: intermediate 131 132 .keywords: distributed array, get, component name 133 134 .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName() 135 @*/ 136 PetscErrorCode DMDAGetCoordinateName(DM dm,PetscInt nf,const char **name) 137 { 138 DM_DA *dd = (DM_DA*)dm->data; 139 140 PetscFunctionBegin; 141 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 142 PetscValidPointer(name,3); 143 if (nf < 0 || nf >= dm->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf); 144 *name = dd->coordinatename[nf]; 145 PetscFunctionReturn(0); 146 } 147 148 #undef __FUNCT__ 149 #define __FUNCT__ "DMDAGetCorners" 150 /*@ 151 DMDAGetCorners - Returns the global (x,y,z) indices of the lower left 152 corner of the local region, excluding ghost points. 153 154 Not Collective 155 156 Input Parameter: 157 . da - the distributed array 158 159 Output Parameters: 160 + x,y,z - the corner indices (where y and z are optional; these are used 161 for 2D and 3D problems) 162 - m,n,p - widths in the corresponding directions (where n and p are optional; 163 these are used for 2D and 3D problems) 164 165 Note: 166 The corner information is independent of the number of degrees of 167 freedom per node set with the DMDACreateXX() routine. Thus the x, y, z, and 168 m, n, p can be thought of as coordinates on a logical grid, where each 169 grid point has (potentially) several degrees of freedom. 170 Any of y, z, n, and p can be passed in as NULL if not needed. 171 172 Level: beginner 173 174 .keywords: distributed array, get, corners, nodes, local indices 175 176 .seealso: DMDAGetGhostCorners(), DMDAGetOwnershipRanges() 177 @*/ 178 PetscErrorCode DMDAGetCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p) 179 { 180 PetscInt w; 181 DM_DA *dd = (DM_DA*)da->data; 182 183 PetscFunctionBegin; 184 PetscValidHeaderSpecific(da,DM_CLASSID,1); 185 /* since the xs, xe ... have all been multiplied by the number of degrees 186 of freedom per cell, w = dd->w, we divide that out before returning.*/ 187 w = dd->w; 188 if (x) *x = dd->xs/w + dd->xo; if (m) *m = (dd->xe - dd->xs)/w; 189 /* the y and z have NOT been multiplied by w */ 190 if (y) *y = dd->ys + dd->yo; if (n) *n = (dd->ye - dd->ys); 191 if (z) *z = dd->zs + dd->zo; if (p) *p = (dd->ze - dd->zs); 192 PetscFunctionReturn(0); 193 } 194 195 #undef __FUNCT__ 196 #define __FUNCT__ "DMDAGetLocalBoundingBox" 197 /*@ 198 DMDAGetLocalBoundingBox - Returns the local bounding box for the DMDA. 199 200 Not Collective 201 202 Input Parameter: 203 . dm - the DM 204 205 Output Parameters: 206 + lmin - local minimum coordinates (length dim, optional) 207 - lmax - local maximim coordinates (length dim, optional) 208 209 Level: beginner 210 211 .keywords: distributed array, get, coordinates 212 213 .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetBoundingBox() 214 @*/ 215 PetscErrorCode DMDAGetLocalBoundingBox(DM dm,PetscReal lmin[],PetscReal lmax[]) 216 { 217 PetscErrorCode ierr; 218 Vec coords = NULL; 219 PetscInt dim,i,j; 220 const PetscScalar *local_coords; 221 PetscReal min[3]={PETSC_MAX_REAL,PETSC_MAX_REAL,PETSC_MAX_REAL},max[3]={PETSC_MIN_REAL,PETSC_MIN_REAL,PETSC_MIN_REAL}; 222 PetscInt N,Ni; 223 224 PetscFunctionBegin; 225 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 226 dim = dm->dim; 227 ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr); 228 if (coords) { 229 ierr = VecGetArrayRead(coords,&local_coords);CHKERRQ(ierr); 230 ierr = VecGetLocalSize(coords,&N);CHKERRQ(ierr); 231 Ni = N/dim; 232 for (i=0; i<Ni; i++) { 233 for (j=0; j<3; j++) { 234 min[j] = j < dim ? PetscMin(min[j],PetscRealPart(local_coords[i*dim+j])) : 0; 235 max[j] = j < dim ? PetscMax(max[j],PetscRealPart(local_coords[i*dim+j])) : 0; 236 } 237 } 238 ierr = VecRestoreArrayRead(coords,&local_coords);CHKERRQ(ierr); 239 } else { /* Just use grid indices */ 240 DMDALocalInfo info; 241 ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 242 min[0] = info.xs; 243 min[1] = info.ys; 244 min[2] = info.zs; 245 max[0] = info.xs + info.xm-1; 246 max[1] = info.ys + info.ym-1; 247 max[2] = info.zs + info.zm-1; 248 } 249 if (lmin) {ierr = PetscMemcpy(lmin,min,dim*sizeof(PetscReal));CHKERRQ(ierr);} 250 if (lmax) {ierr = PetscMemcpy(lmax,max,dim*sizeof(PetscReal));CHKERRQ(ierr);} 251 PetscFunctionReturn(0); 252 } 253 254 #undef __FUNCT__ 255 #define __FUNCT__ "DMDAGetBoundingBox" 256 /*@ 257 DMDAGetBoundingBox - Returns the global bounding box for the DMDA. 258 259 Collective on DMDA 260 261 Input Parameter: 262 . dm - the DM 263 264 Output Parameters: 265 + gmin - global minimum coordinates (length dim, optional) 266 - gmax - global maximim coordinates (length dim, optional) 267 268 Level: beginner 269 270 .keywords: distributed array, get, coordinates 271 272 .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetLocalBoundingBox() 273 @*/ 274 PetscErrorCode DMDAGetBoundingBox(DM dm,PetscReal gmin[],PetscReal gmax[]) 275 { 276 PetscErrorCode ierr; 277 PetscMPIInt count; 278 PetscReal lmin[3],lmax[3]; 279 280 PetscFunctionBegin; 281 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 282 ierr = PetscMPIIntCast(dm->dim,&count);CHKERRQ(ierr); 283 ierr = DMDAGetLocalBoundingBox(dm,lmin,lmax);CHKERRQ(ierr); 284 if (gmin) {ierr = MPI_Allreduce(lmin,gmin,count,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);} 285 if (gmax) {ierr = MPI_Allreduce(lmax,gmax,count,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);} 286 PetscFunctionReturn(0); 287 } 288 289 #undef __FUNCT__ 290 #define __FUNCT__ "DMDAGetReducedDMDA" 291 /*@ 292 DMDAGetReducedDMDA - Gets the DMDA with the same layout but with fewer or more fields 293 294 Collective on DMDA 295 296 Input Parameter: 297 + da - the distributed array 298 . nfields - number of fields in new DMDA 299 300 Output Parameter: 301 . nda - the new DMDA 302 303 Level: intermediate 304 305 .keywords: distributed array, get, corners, nodes, local indices, coordinates 306 307 .seealso: DMDAGetGhostCorners(), DMSetCoordinates(), DMDASetUniformCoordinates(), DMGetCoordinates(), DMDAGetGhostedCoordinates() 308 @*/ 309 PetscErrorCode DMDAGetReducedDMDA(DM da,PetscInt nfields,DM *nda) 310 { 311 PetscErrorCode ierr; 312 DM_DA *dd = (DM_DA*)da->data; 313 PetscInt s,m,n,p,M,N,P,dim,Mo,No,Po; 314 const PetscInt *lx,*ly,*lz; 315 DMBoundaryType bx,by,bz; 316 DMDAStencilType stencil_type; 317 PetscInt ox,oy,oz; 318 PetscInt cl,rl; 319 320 PetscFunctionBegin; 321 dim = da->dim; 322 M = dd->M; 323 N = dd->N; 324 P = dd->P; 325 m = dd->m; 326 n = dd->n; 327 p = dd->p; 328 s = dd->s; 329 bx = dd->bx; 330 by = dd->by; 331 bz = dd->bz; 332 333 stencil_type = dd->stencil_type; 334 335 ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr); 336 if (dim == 1) { 337 ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),bx,M,nfields,s,dd->lx,nda);CHKERRQ(ierr); 338 } else if (dim == 2) { 339 ierr = DMDACreate2d(PetscObjectComm((PetscObject)da),bx,by,stencil_type,M,N,m,n,nfields,s,lx,ly,nda);CHKERRQ(ierr); 340 } else if (dim == 3) { 341 ierr = DMDACreate3d(PetscObjectComm((PetscObject)da),bx,by,bz,stencil_type,M,N,P,m,n,p,nfields,s,lx,ly,lz,nda);CHKERRQ(ierr); 342 } 343 if (da->coordinates) { 344 ierr = PetscObjectReference((PetscObject)da->coordinates);CHKERRQ(ierr); 345 346 (*nda)->coordinates = da->coordinates; 347 } 348 349 /* allow for getting a reduced DA corresponding to a domain decomposition */ 350 ierr = DMDAGetOffset(da,&ox,&oy,&oz,&Mo,&No,&Po);CHKERRQ(ierr); 351 ierr = DMDASetOffset(*nda,ox,oy,oz,Mo,No,Po);CHKERRQ(ierr); 352 353 /* allow for getting a reduced DA corresponding to a coarsened DA */ 354 ierr = DMGetCoarsenLevel(da,&cl);CHKERRQ(ierr); 355 ierr = DMGetRefineLevel(da,&rl);CHKERRQ(ierr); 356 357 (*nda)->levelup = rl; 358 (*nda)->leveldown = cl; 359 PetscFunctionReturn(0); 360 } 361 362