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