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 PetscErrorCode DMCreateCoordinateDM_DA(DM dm, DM *cdm) 9 { 10 PetscErrorCode ierr; 11 PetscFunctionBegin; 12 ierr = DMDAGetReducedDMDA(dm,dm->dim,cdm);CHKERRQ(ierr); 13 PetscFunctionReturn(0); 14 } 15 16 /*@C 17 DMDASetFieldName - Sets the names of individual field components in multicomponent 18 vectors associated with a DMDA. 19 20 Not Collective 21 22 Input Parameters: 23 + da - the distributed array 24 . nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the 25 number of degrees of freedom per node within the DMDA 26 - names - the name of the field (component) 27 28 Notes: It must be called after having called DMSetUp(). 29 30 Level: intermediate 31 32 .keywords: distributed array, get, component name 33 34 .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMDASetFieldNames(), DMSetUp() 35 @*/ 36 PetscErrorCode DMDASetFieldName(DM da,PetscInt nf,const char name[]) 37 { 38 PetscErrorCode ierr; 39 DM_DA *dd = (DM_DA*)da->data; 40 41 PetscFunctionBegin; 42 PetscValidHeaderSpecific(da,DM_CLASSID,1); 43 if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf); 44 if (!dd->fieldname) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 45 ierr = PetscFree(dd->fieldname[nf]);CHKERRQ(ierr); 46 ierr = PetscStrallocpy(name,&dd->fieldname[nf]);CHKERRQ(ierr); 47 PetscFunctionReturn(0); 48 } 49 50 /*@C 51 DMDAGetFieldNames - Gets the name of each component in the vector associated with the DMDA 52 53 Collective on TS 54 55 Input Parameter: 56 . dm - the DMDA object 57 58 Output Parameter: 59 . names - the names of the components, final string is NULL, will have the same number of entries as the dof used in creating the DMDA 60 61 Level: intermediate 62 63 Not supported from Fortran, use DMDAGetFieldName() 64 65 .keywords: distributed array, get, component name 66 67 .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMDASetFieldName(), DMDASetFieldNames() 68 @*/ 69 PetscErrorCode DMDAGetFieldNames(DM da,const char * const **names) 70 { 71 DM_DA *dd = (DM_DA*)da->data; 72 73 PetscFunctionBegin; 74 *names = (const char * const *) dd->fieldname; 75 PetscFunctionReturn(0); 76 } 77 78 /*@C 79 DMDASetFieldNames - Sets the name of each component in the vector associated with the DMDA 80 81 Collective on TS 82 83 Input Parameters: 84 + dm - the DMDA object 85 - names - the names of the components, final string must be NULL, must have the same number of entries as the dof used in creating the DMDA 86 87 Notes: It must be called after having called DMSetUp(). 88 89 Level: intermediate 90 91 Not supported from Fortran, use DMDASetFieldName() 92 93 .keywords: distributed array, get, component name 94 95 .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMDASetFieldName(), DMSetUp() 96 @*/ 97 PetscErrorCode DMDASetFieldNames(DM da,const char * const *names) 98 { 99 PetscErrorCode ierr; 100 DM_DA *dd = (DM_DA*)da->data; 101 char **fieldname; 102 PetscInt nf = 0; 103 104 PetscFunctionBegin; 105 if (!dd->fieldname) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 106 while (names[nf++]) {}; 107 if (nf != dd->w+1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid number of fields %D",nf-1); 108 ierr = PetscStrArrayallocpy(names,&fieldname);CHKERRQ(ierr); 109 ierr = PetscStrArrayDestroy(&dd->fieldname);CHKERRQ(ierr); 110 dd->fieldname = fieldname; 111 PetscFunctionReturn(0); 112 } 113 114 /*@C 115 DMDAGetFieldName - Gets the names of individual field components in multicomponent 116 vectors associated with a DMDA. 117 118 Not Collective 119 120 Input Parameter: 121 + da - the distributed array 122 - nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the 123 number of degrees of freedom per node within the DMDA 124 125 Output Parameter: 126 . names - the name of the field (component) 127 128 Notes: It must be called after having called DMSetUp(). 129 130 Level: intermediate 131 132 .keywords: distributed array, get, component name 133 134 .seealso: DMDASetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMSetUp() 135 @*/ 136 PetscErrorCode DMDAGetFieldName(DM da,PetscInt nf,const char **name) 137 { 138 DM_DA *dd = (DM_DA*)da->data; 139 140 PetscFunctionBegin; 141 PetscValidHeaderSpecific(da,DM_CLASSID,1); 142 PetscValidPointer(name,3); 143 if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf); 144 if (!dd->fieldname) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 145 *name = dd->fieldname[nf]; 146 PetscFunctionReturn(0); 147 } 148 149 /*@C 150 DMDASetCoordinateName - Sets the name of the coordinate directions associated with a DMDA, for example "x" or "y" 151 152 Not Collective 153 154 Input Parameters: 155 + dm - the DM 156 . nf - coordinate number for the DMDA (0, 1, ... dim-1), 157 - name - the name of the coordinate 158 159 Notes: It must be called after having called DMSetUp(). 160 161 Level: intermediate 162 163 Not supported from Fortran 164 165 .keywords: distributed array, get, component name 166 167 .seealso: DMDAGetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName(), DMSetUp() 168 @*/ 169 PetscErrorCode DMDASetCoordinateName(DM dm,PetscInt nf,const char name[]) 170 { 171 PetscErrorCode ierr; 172 DM_DA *dd = (DM_DA*)dm->data; 173 174 PetscFunctionBegin; 175 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 176 if (nf < 0 || nf >= dm->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf); 177 if (!dd->coordinatename) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 178 ierr = PetscFree(dd->coordinatename[nf]);CHKERRQ(ierr); 179 ierr = PetscStrallocpy(name,&dd->coordinatename[nf]);CHKERRQ(ierr); 180 PetscFunctionReturn(0); 181 } 182 183 /*@C 184 DMDAGetCoordinateName - Gets the name of a coodinate direction associated with a DMDA. 185 186 Not Collective 187 188 Input Parameter: 189 + dm - the DM 190 - nf - number for the DMDA (0, 1, ... dim-1) 191 192 Output Parameter: 193 . names - the name of the coordinate direction 194 195 Notes: It must be called after having called DMSetUp(). 196 197 Level: intermediate 198 199 Not supported from Fortran 200 201 .keywords: distributed array, get, component name 202 203 .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName(), DMSetUp() 204 @*/ 205 PetscErrorCode DMDAGetCoordinateName(DM dm,PetscInt nf,const char **name) 206 { 207 DM_DA *dd = (DM_DA*)dm->data; 208 209 PetscFunctionBegin; 210 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 211 PetscValidPointer(name,3); 212 if (nf < 0 || nf >= dm->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf); 213 if (!dd->coordinatename) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 214 *name = dd->coordinatename[nf]; 215 PetscFunctionReturn(0); 216 } 217 218 /*@C 219 DMDAGetCorners - Returns the global (x,y,z) indices of the lower left 220 corner and size of the local region, excluding ghost points. 221 222 Not Collective 223 224 Input Parameter: 225 . da - the distributed array 226 227 Output Parameters: 228 + x,y,z - the corner indices (where y and z are optional; these are used 229 for 2D and 3D problems) 230 - m,n,p - widths in the corresponding directions (where n and p are optional; 231 these are used for 2D and 3D problems) 232 233 Note: 234 The corner information is independent of the number of degrees of 235 freedom per node set with the DMDACreateXX() routine. Thus the x, y, z, and 236 m, n, p can be thought of as coordinates on a logical grid, where each 237 grid point has (potentially) several degrees of freedom. 238 Any of y, z, n, and p can be passed in as NULL if not needed. 239 240 Level: beginner 241 242 .keywords: distributed array, get, corners, nodes, local indices 243 244 .seealso: DMDAGetGhostCorners(), DMDAGetOwnershipRanges() 245 @*/ 246 PetscErrorCode DMDAGetCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p) 247 { 248 PetscInt w; 249 DM_DA *dd = (DM_DA*)da->data; 250 251 PetscFunctionBegin; 252 PetscValidHeaderSpecific(da,DM_CLASSID,1); 253 /* since the xs, xe ... have all been multiplied by the number of degrees 254 of freedom per cell, w = dd->w, we divide that out before returning.*/ 255 w = dd->w; 256 if (x) *x = dd->xs/w + dd->xo; 257 /* the y and z have NOT been multiplied by w */ 258 if (y) *y = dd->ys + dd->yo; 259 if (z) *z = dd->zs + dd->zo; 260 if (m) *m = (dd->xe - dd->xs)/w; 261 if (n) *n = (dd->ye - dd->ys); 262 if (p) *p = (dd->ze - dd->zs); 263 PetscFunctionReturn(0); 264 } 265 266 /*@ 267 DMDAGetLocalBoundingBox - Returns the local bounding box for the DMDA. 268 269 Not Collective 270 271 Input Parameter: 272 . dm - the DM 273 274 Output Parameters: 275 + lmin - local minimum coordinates (length dim, optional) 276 - lmax - local maximim coordinates (length dim, optional) 277 278 Level: beginner 279 280 Not supported from Fortran 281 282 .keywords: distributed array, get, coordinates 283 284 .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetBoundingBox() 285 @*/ 286 PetscErrorCode DMDAGetLocalBoundingBox(DM dm,PetscReal lmin[],PetscReal lmax[]) 287 { 288 PetscErrorCode ierr; 289 Vec coords = NULL; 290 PetscInt dim,i,j; 291 const PetscScalar *local_coords; 292 PetscReal min[3]={PETSC_MAX_REAL,PETSC_MAX_REAL,PETSC_MAX_REAL},max[3]={PETSC_MIN_REAL,PETSC_MIN_REAL,PETSC_MIN_REAL}; 293 PetscInt N,Ni; 294 295 PetscFunctionBegin; 296 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 297 dim = dm->dim; 298 ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr); 299 if (coords) { 300 ierr = VecGetArrayRead(coords,&local_coords);CHKERRQ(ierr); 301 ierr = VecGetLocalSize(coords,&N);CHKERRQ(ierr); 302 Ni = N/dim; 303 for (i=0; i<Ni; i++) { 304 for (j=0; j<3; j++) { 305 min[j] = j < dim ? PetscMin(min[j],PetscRealPart(local_coords[i*dim+j])) : 0; 306 max[j] = j < dim ? PetscMax(max[j],PetscRealPart(local_coords[i*dim+j])) : 0; 307 } 308 } 309 ierr = VecRestoreArrayRead(coords,&local_coords);CHKERRQ(ierr); 310 } else { /* Just use grid indices */ 311 DMDALocalInfo info; 312 ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 313 min[0] = info.xs; 314 min[1] = info.ys; 315 min[2] = info.zs; 316 max[0] = info.xs + info.xm-1; 317 max[1] = info.ys + info.ym-1; 318 max[2] = info.zs + info.zm-1; 319 } 320 if (lmin) {ierr = PetscMemcpy(lmin,min,dim*sizeof(PetscReal));CHKERRQ(ierr);} 321 if (lmax) {ierr = PetscMemcpy(lmax,max,dim*sizeof(PetscReal));CHKERRQ(ierr);} 322 PetscFunctionReturn(0); 323 } 324 325 /*@ 326 DMDAGetBoundingBox - Returns the global bounding box for the DMDA. 327 328 Collective on DMDA 329 330 Input Parameter: 331 . dm - the DM 332 333 Output Parameters: 334 + gmin - global minimum coordinates (length dim, optional) 335 - gmax - global maximim coordinates (length dim, optional) 336 337 Level: beginner 338 339 .keywords: distributed array, get, coordinates 340 341 .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetLocalBoundingBox() 342 @*/ 343 PetscErrorCode DMDAGetBoundingBox(DM dm,PetscReal gmin[],PetscReal gmax[]) 344 { 345 PetscErrorCode ierr; 346 PetscMPIInt count; 347 PetscReal lmin[3],lmax[3]; 348 349 PetscFunctionBegin; 350 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 351 ierr = PetscMPIIntCast(dm->dim,&count);CHKERRQ(ierr); 352 ierr = DMDAGetLocalBoundingBox(dm,lmin,lmax);CHKERRQ(ierr); 353 if (gmin) {ierr = MPIU_Allreduce(lmin,gmin,count,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);} 354 if (gmax) {ierr = MPIU_Allreduce(lmax,gmax,count,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);} 355 PetscFunctionReturn(0); 356 } 357 358 /*@ 359 DMDAGetReducedDMDA - Gets the DMDA with the same layout but with fewer or more fields 360 361 Collective on DMDA 362 363 Input Parameters: 364 + da - the distributed array 365 - nfields - number of fields in new DMDA 366 367 Output Parameter: 368 . nda - the new DMDA 369 370 Level: intermediate 371 372 .keywords: distributed array, get, corners, nodes, local indices, coordinates 373 374 .seealso: DMDAGetGhostCorners(), DMSetCoordinates(), DMDASetUniformCoordinates(), DMGetCoordinates(), DMDAGetGhostedCoordinates() 375 @*/ 376 PetscErrorCode DMDAGetReducedDMDA(DM da,PetscInt nfields,DM *nda) 377 { 378 PetscErrorCode ierr; 379 DM_DA *dd = (DM_DA*)da->data; 380 PetscInt s,m,n,p,M,N,P,dim,Mo,No,Po; 381 const PetscInt *lx,*ly,*lz; 382 DMBoundaryType bx,by,bz; 383 DMDAStencilType stencil_type; 384 PetscInt ox,oy,oz; 385 PetscInt cl,rl; 386 387 PetscFunctionBegin; 388 dim = da->dim; 389 M = dd->M; 390 N = dd->N; 391 P = dd->P; 392 m = dd->m; 393 n = dd->n; 394 p = dd->p; 395 s = dd->s; 396 bx = dd->bx; 397 by = dd->by; 398 bz = dd->bz; 399 400 stencil_type = dd->stencil_type; 401 402 ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr); 403 if (dim == 1) { 404 ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),bx,M,nfields,s,dd->lx,nda);CHKERRQ(ierr); 405 } else if (dim == 2) { 406 ierr = DMDACreate2d(PetscObjectComm((PetscObject)da),bx,by,stencil_type,M,N,m,n,nfields,s,lx,ly,nda);CHKERRQ(ierr); 407 } else if (dim == 3) { 408 ierr = DMDACreate3d(PetscObjectComm((PetscObject)da),bx,by,bz,stencil_type,M,N,P,m,n,p,nfields,s,lx,ly,lz,nda);CHKERRQ(ierr); 409 } 410 ierr = DMSetUp(*nda);CHKERRQ(ierr); 411 if (da->coordinates) { 412 ierr = PetscObjectReference((PetscObject)da->coordinates);CHKERRQ(ierr); 413 (*nda)->coordinates = da->coordinates; 414 } 415 416 /* allow for getting a reduced DA corresponding to a domain decomposition */ 417 ierr = DMDAGetOffset(da,&ox,&oy,&oz,&Mo,&No,&Po);CHKERRQ(ierr); 418 ierr = DMDASetOffset(*nda,ox,oy,oz,Mo,No,Po);CHKERRQ(ierr); 419 420 /* allow for getting a reduced DA corresponding to a coarsened DA */ 421 ierr = DMGetCoarsenLevel(da,&cl);CHKERRQ(ierr); 422 ierr = DMGetRefineLevel(da,&rl);CHKERRQ(ierr); 423 424 (*nda)->levelup = rl; 425 (*nda)->leveldown = cl; 426 PetscFunctionReturn(0); 427 } 428 429 /*@C 430 DMDAGetCoordinateArray - Gets an array containing the coordinates of the DMDA 431 432 Not Collective 433 434 Input Parameter: 435 . dm - the DM 436 437 Output Parameter: 438 . xc - the coordinates 439 440 Level: intermediate 441 442 Not supported from Fortran 443 444 .keywords: distributed array, get, component name 445 446 .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName(), DMDARestoreCoordinateArray() 447 @*/ 448 PetscErrorCode DMDAGetCoordinateArray(DM dm,void *xc) 449 { 450 PetscErrorCode ierr; 451 DM cdm; 452 Vec x; 453 454 PetscFunctionBegin; 455 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 456 ierr = DMGetCoordinates(dm,&x);CHKERRQ(ierr); 457 ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr); 458 ierr = DMDAVecGetArray(cdm,x,xc);CHKERRQ(ierr); 459 PetscFunctionReturn(0); 460 } 461 462 /*@C 463 DMDARestoreCoordinateArray - Sets an array containing the coordinates of the DMDA 464 465 Not Collective 466 467 Input Parameter: 468 + dm - the DM 469 - xc - the coordinates 470 471 Level: intermediate 472 473 Not supported from Fortran 474 475 .keywords: distributed array, get, component name 476 477 .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName(), DMDAGetCoordinateArray() 478 @*/ 479 PetscErrorCode DMDARestoreCoordinateArray(DM dm,void *xc) 480 { 481 PetscErrorCode ierr; 482 DM cdm; 483 Vec x; 484 485 PetscFunctionBegin; 486 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 487 ierr = DMGetCoordinates(dm,&x);CHKERRQ(ierr); 488 ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr); 489 ierr = DMDAVecRestoreArray(cdm,x,xc);CHKERRQ(ierr); 490 PetscFunctionReturn(0); 491 } 492