1 #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 2 3 /*@ 4 DMDASetSizes - Sets the number of grid points in the three dimensional directions 5 6 Logically Collective on DMDA 7 8 Input Parameters: 9 + da - the DMDA 10 . M - the global X size 11 . N - the global Y size 12 - P - the global Z size 13 14 Level: intermediate 15 16 Developer Notes: 17 Since the dimension may not yet have been set the code cannot error check for non-positive Y and Z number of grid points 18 19 .seealso: PetscSplitOwnership() 20 @*/ 21 PetscErrorCode DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P) 22 { 23 DM_DA *dd = (DM_DA*)da->data; 24 25 PetscFunctionBegin; 26 PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 27 PetscValidLogicalCollectiveInt(da,M,2); 28 PetscValidLogicalCollectiveInt(da,N,3); 29 PetscValidLogicalCollectiveInt(da,P,4); 30 if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 31 if (M < 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in X direction must be positive"); 32 if (N < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in Y direction must be positive"); 33 if (P < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in Z direction must be positive"); 34 35 dd->M = M; 36 dd->N = N; 37 dd->P = P; 38 PetscFunctionReturn(0); 39 } 40 41 /*@ 42 DMDASetNumProcs - Sets the number of processes in each dimension 43 44 Logically Collective on DMDA 45 46 Input Parameters: 47 + da - the DMDA 48 . m - the number of X procs (or PETSC_DECIDE) 49 . n - the number of Y procs (or PETSC_DECIDE) 50 - p - the number of Z procs (or PETSC_DECIDE) 51 52 Level: intermediate 53 54 .seealso: DMDASetSizes(), PetscSplitOwnership() 55 @*/ 56 PetscErrorCode DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p) 57 { 58 DM_DA *dd = (DM_DA*)da->data; 59 PetscErrorCode ierr; 60 61 PetscFunctionBegin; 62 PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 63 PetscValidLogicalCollectiveInt(da,m,2); 64 PetscValidLogicalCollectiveInt(da,n,3); 65 PetscValidLogicalCollectiveInt(da,p,4); 66 if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 67 dd->m = m; 68 dd->n = n; 69 dd->p = p; 70 if (da->dim == 2) { 71 PetscMPIInt size; 72 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr); 73 if ((dd->m > 0) && (dd->n < 0)) { 74 dd->n = size/dd->m; 75 if (dd->n*dd->m != size) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"%D processes in X direction not divisible into comm size %d",m,size); 76 } 77 if ((dd->n > 0) && (dd->m < 0)) { 78 dd->m = size/dd->n; 79 if (dd->n*dd->m != size) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"%D processes in Y direction not divisible into comm size %d",n,size); 80 } 81 } 82 PetscFunctionReturn(0); 83 } 84 85 /*@ 86 DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries. 87 88 Not collective 89 90 Input Parameter: 91 + da - The DMDA 92 - bx,by,bz - One of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC 93 94 Level: intermediate 95 96 .keywords: distributed array, periodicity 97 .seealso: DMDACreate(), DMDestroy(), DMDA, DMBoundaryType 98 @*/ 99 PetscErrorCode DMDASetBoundaryType(DM da,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz) 100 { 101 DM_DA *dd = (DM_DA*)da->data; 102 103 PetscFunctionBegin; 104 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 105 PetscValidLogicalCollectiveEnum(da,bx,2); 106 PetscValidLogicalCollectiveEnum(da,by,3); 107 PetscValidLogicalCollectiveEnum(da,bz,4); 108 if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 109 dd->bx = bx; 110 dd->by = by; 111 dd->bz = bz; 112 PetscFunctionReturn(0); 113 } 114 115 /*@ 116 DMDASetDof - Sets the number of degrees of freedom per vertex 117 118 Not collective 119 120 Input Parameters: 121 + da - The DMDA 122 - dof - Number of degrees of freedom 123 124 Level: intermediate 125 126 .keywords: distributed array, degrees of freedom 127 .seealso: DMDAGetDof(), DMDACreate(), DMDestroy(), DMDA 128 @*/ 129 PetscErrorCode DMDASetDof(DM da, PetscInt dof) 130 { 131 DM_DA *dd = (DM_DA*)da->data; 132 133 PetscFunctionBegin; 134 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 135 PetscValidLogicalCollectiveInt(da,dof,2); 136 if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 137 dd->w = dof; 138 da->bs = dof; 139 PetscFunctionReturn(0); 140 } 141 142 /*@ 143 DMDAGetDof - Gets the number of degrees of freedom per vertex 144 145 Not collective 146 147 Input Parameter: 148 . da - The DMDA 149 150 Output Parameter: 151 . dof - Number of degrees of freedom 152 153 Level: intermediate 154 155 .keywords: distributed array, degrees of freedom 156 .seealso: DMDASetDof(), DMDACreate(), DMDestroy(), DMDA 157 @*/ 158 PetscErrorCode DMDAGetDof(DM da, PetscInt *dof) 159 { 160 DM_DA *dd = (DM_DA *) da->data; 161 162 PetscFunctionBegin; 163 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 164 PetscValidPointer(dof,2); 165 *dof = dd->w; 166 PetscFunctionReturn(0); 167 } 168 169 /*@ 170 DMDAGetOverlap - Gets the size of the per-processor overlap. 171 172 Not collective 173 174 Input Parameters: 175 . da - The DMDA 176 177 Output Parameters: 178 + x - Overlap in the x direction 179 . y - Overlap in the y direction 180 - z - Overlap in the z direction 181 182 Level: intermediate 183 184 .keywords: distributed array, overlap, domain decomposition 185 .seealso: DMDACreateDomainDecomposition(), DMDASetOverlap(), DMDA 186 @*/ 187 PetscErrorCode DMDAGetOverlap(DM da,PetscInt *x,PetscInt *y,PetscInt *z) 188 { 189 DM_DA *dd = (DM_DA*)da->data; 190 191 PetscFunctionBegin; 192 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 193 if (x) *x = dd->xol; 194 if (y) *y = dd->yol; 195 if (z) *z = dd->zol; 196 PetscFunctionReturn(0); 197 } 198 199 /*@ 200 DMDASetOverlap - Sets the size of the per-processor overlap. 201 202 Not collective 203 204 Input Parameters: 205 + da - The DMDA 206 . x - Overlap in the x direction 207 . y - Overlap in the y direction 208 - z - Overlap in the z direction 209 210 Level: intermediate 211 212 .keywords: distributed array, overlap, domain decomposition 213 .seealso: DMDACreateDomainDecomposition(), DMDAGetOverlap(), DMDA 214 @*/ 215 PetscErrorCode DMDASetOverlap(DM da,PetscInt x,PetscInt y,PetscInt z) 216 { 217 DM_DA *dd = (DM_DA*)da->data; 218 219 PetscFunctionBegin; 220 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 221 PetscValidLogicalCollectiveInt(da,x,2); 222 PetscValidLogicalCollectiveInt(da,y,3); 223 PetscValidLogicalCollectiveInt(da,z,4); 224 dd->xol = x; 225 dd->yol = y; 226 dd->zol = z; 227 PetscFunctionReturn(0); 228 } 229 230 231 /*@ 232 DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition. 233 234 Not collective 235 236 Input Parameters: 237 . da - The DMDA 238 239 Output Parameters: 240 + Nsub - Number of local subdomains created upon decomposition 241 242 Level: intermediate 243 244 .keywords: distributed array, domain decomposition 245 .seealso: DMDACreateDomainDecomposition(), DMDASetNumLocalSubDomains(), DMDA 246 @*/ 247 PetscErrorCode DMDAGetNumLocalSubDomains(DM da,PetscInt *Nsub) 248 { 249 DM_DA *dd = (DM_DA*)da->data; 250 251 PetscFunctionBegin; 252 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 253 if (Nsub) *Nsub = dd->Nsub; 254 PetscFunctionReturn(0); 255 } 256 257 /*@ 258 DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition. 259 260 Not collective 261 262 Input Parameters: 263 + da - The DMDA 264 - Nsub - The number of local subdomains requested 265 266 Level: intermediate 267 268 .keywords: distributed array, domain decomposition 269 .seealso: DMDACreateDomainDecomposition(), DMDAGetNumLocalSubDomains(), DMDA 270 @*/ 271 PetscErrorCode DMDASetNumLocalSubDomains(DM da,PetscInt Nsub) 272 { 273 DM_DA *dd = (DM_DA*)da->data; 274 275 PetscFunctionBegin; 276 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 277 PetscValidLogicalCollectiveInt(da,Nsub,2); 278 dd->Nsub = Nsub; 279 PetscFunctionReturn(0); 280 } 281 282 /*@ 283 DMDASetOffset - Sets the index offset of the DA. 284 285 Collective on DA 286 287 Input Parameter: 288 + da - The DMDA 289 . xo - The offset in the x direction 290 . yo - The offset in the y direction 291 - zo - The offset in the z direction 292 293 Level: intermediate 294 295 Notes: 296 This is used primarily to overlap a computation on a local DA with that on a global DA without 297 changing boundary conditions or subdomain features that depend upon the global offsets. 298 299 .keywords: distributed array, degrees of freedom 300 .seealso: DMDAGetOffset(), DMDAVecGetArray() 301 @*/ 302 PetscErrorCode DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po) 303 { 304 PetscErrorCode ierr; 305 DM_DA *dd = (DM_DA*)da->data; 306 307 PetscFunctionBegin; 308 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 309 PetscValidLogicalCollectiveInt(da,xo,2); 310 PetscValidLogicalCollectiveInt(da,yo,3); 311 PetscValidLogicalCollectiveInt(da,zo,4); 312 PetscValidLogicalCollectiveInt(da,Mo,5); 313 PetscValidLogicalCollectiveInt(da,No,6); 314 PetscValidLogicalCollectiveInt(da,Po,7); 315 dd->xo = xo; 316 dd->yo = yo; 317 dd->zo = zo; 318 dd->Mo = Mo; 319 dd->No = No; 320 dd->Po = Po; 321 322 if (da->coordinateDM) { 323 ierr = DMDASetOffset(da->coordinateDM,xo,yo,zo,Mo,No,Po);CHKERRQ(ierr); 324 } 325 PetscFunctionReturn(0); 326 } 327 328 /*@ 329 DMDAGetOffset - Gets the index offset of the DA. 330 331 Not collective 332 333 Input Parameter: 334 . da - The DMDA 335 336 Output Parameters: 337 + xo - The offset in the x direction 338 . yo - The offset in the y direction 339 . zo - The offset in the z direction 340 . Mo - The global size in the x direction 341 . No - The global size in the y direction 342 - Po - The global size in the z direction 343 344 Level: intermediate 345 346 .keywords: distributed array, degrees of freedom 347 .seealso: DMDASetOffset(), DMDAVecGetArray() 348 @*/ 349 PetscErrorCode DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po) 350 { 351 DM_DA *dd = (DM_DA*)da->data; 352 353 PetscFunctionBegin; 354 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 355 if (xo) *xo = dd->xo; 356 if (yo) *yo = dd->yo; 357 if (zo) *zo = dd->zo; 358 if (Mo) *Mo = dd->Mo; 359 if (No) *No = dd->No; 360 if (Po) *Po = dd->Po; 361 PetscFunctionReturn(0); 362 } 363 364 /*@ 365 DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain DM. 366 367 Not collective 368 369 Input Parameter: 370 . da - The DMDA 371 372 Output Parameters: 373 + xs - The start of the region in x 374 . ys - The start of the region in y 375 . zs - The start of the region in z 376 . xs - The size of the region in x 377 . ys - The size of the region in y 378 . zs - The size of the region in z 379 380 Level: intermediate 381 382 .keywords: distributed array, degrees of freedom 383 .seealso: DMDAGetOffset(), DMDAVecGetArray() 384 @*/ 385 PetscErrorCode DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm) 386 { 387 DM_DA *dd = (DM_DA*)da->data; 388 389 PetscFunctionBegin; 390 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 391 if (xs) *xs = dd->nonxs; 392 if (ys) *ys = dd->nonys; 393 if (zs) *zs = dd->nonzs; 394 if (xm) *xm = dd->nonxm; 395 if (ym) *ym = dd->nonym; 396 if (zm) *zm = dd->nonzm; 397 PetscFunctionReturn(0); 398 } 399 400 401 /*@ 402 DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain DM. 403 404 Collective on DA 405 406 Input Parameter: 407 + da - The DMDA 408 . xs - The start of the region in x 409 . ys - The start of the region in y 410 . zs - The start of the region in z 411 . xs - The size of the region in x 412 . ys - The size of the region in y 413 . zs - The size of the region in z 414 415 Level: intermediate 416 417 .keywords: distributed array, degrees of freedom 418 .seealso: DMDAGetOffset(), DMDAVecGetArray() 419 @*/ 420 PetscErrorCode DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm) 421 { 422 DM_DA *dd = (DM_DA*)da->data; 423 424 PetscFunctionBegin; 425 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 426 PetscValidLogicalCollectiveInt(da,xs,2); 427 PetscValidLogicalCollectiveInt(da,ys,3); 428 PetscValidLogicalCollectiveInt(da,zs,4); 429 PetscValidLogicalCollectiveInt(da,xm,5); 430 PetscValidLogicalCollectiveInt(da,ym,6); 431 PetscValidLogicalCollectiveInt(da,zm,7); 432 dd->nonxs = xs; 433 dd->nonys = ys; 434 dd->nonzs = zs; 435 dd->nonxm = xm; 436 dd->nonym = ym; 437 dd->nonzm = zm; 438 439 PetscFunctionReturn(0); 440 } 441 442 /*@ 443 DMDASetStencilType - Sets the type of the communication stencil 444 445 Logically Collective on DMDA 446 447 Input Parameter: 448 + da - The DMDA 449 - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 450 451 Level: intermediate 452 453 .keywords: distributed array, stencil 454 .seealso: DMDACreate(), DMDestroy(), DMDA 455 @*/ 456 PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) 457 { 458 DM_DA *dd = (DM_DA*)da->data; 459 460 PetscFunctionBegin; 461 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 462 PetscValidLogicalCollectiveEnum(da,stype,2); 463 if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 464 dd->stencil_type = stype; 465 PetscFunctionReturn(0); 466 } 467 468 /*@ 469 DMDAGetStencilType - Gets the type of the communication stencil 470 471 Not collective 472 473 Input Parameter: 474 . da - The DMDA 475 476 Output Parameter: 477 . stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 478 479 Level: intermediate 480 481 .keywords: distributed array, stencil 482 .seealso: DMDACreate(), DMDestroy(), DMDA 483 @*/ 484 PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype) 485 { 486 DM_DA *dd = (DM_DA*)da->data; 487 488 PetscFunctionBegin; 489 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 490 PetscValidPointer(stype,2); 491 *stype = dd->stencil_type; 492 PetscFunctionReturn(0); 493 } 494 495 /*@ 496 DMDASetStencilWidth - Sets the width of the communication stencil 497 498 Logically Collective on DMDA 499 500 Input Parameter: 501 + da - The DMDA 502 - width - The stencil width 503 504 Level: intermediate 505 506 .keywords: distributed array, stencil 507 .seealso: DMDACreate(), DMDestroy(), DMDA 508 @*/ 509 PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) 510 { 511 DM_DA *dd = (DM_DA*)da->data; 512 513 PetscFunctionBegin; 514 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 515 PetscValidLogicalCollectiveInt(da,width,2); 516 if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 517 dd->s = width; 518 PetscFunctionReturn(0); 519 } 520 521 /*@ 522 DMDAGetStencilWidth - Gets the width of the communication stencil 523 524 Not collective 525 526 Input Parameter: 527 . da - The DMDA 528 529 Output Parameter: 530 . width - The stencil width 531 532 Level: intermediate 533 534 .keywords: distributed array, stencil 535 .seealso: DMDACreate(), DMDestroy(), DMDA 536 @*/ 537 PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width) 538 { 539 DM_DA *dd = (DM_DA *) da->data; 540 541 PetscFunctionBegin; 542 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 543 PetscValidPointer(width,2); 544 *width = dd->s; 545 PetscFunctionReturn(0); 546 } 547 548 static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[]) 549 { 550 PetscInt i,sum; 551 552 PetscFunctionBegin; 553 if (M < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set"); 554 for (i=sum=0; i<m; i++) sum += lx[i]; 555 if (sum != M) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M); 556 PetscFunctionReturn(0); 557 } 558 559 /*@ 560 DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 561 562 Logically Collective on DMDA 563 564 Input Parameter: 565 + da - The DMDA 566 . lx - array containing number of nodes in the X direction on each process, or NULL. If non-null, must be of length da->m 567 . ly - array containing number of nodes in the Y direction on each process, or NULL. If non-null, must be of length da->n 568 - lz - array containing number of nodes in the Z direction on each process, or NULL. If non-null, must be of length da->p. 569 570 Level: intermediate 571 572 Note: these numbers are NOT multiplied by the number of dof per node. 573 574 .keywords: distributed array 575 .seealso: DMDACreate(), DMDestroy(), DMDA 576 @*/ 577 PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 578 { 579 PetscErrorCode ierr; 580 DM_DA *dd = (DM_DA*)da->data; 581 582 PetscFunctionBegin; 583 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 584 if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()"); 585 if (lx) { 586 if (dd->m < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 587 ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr); 588 if (!dd->lx) { 589 ierr = PetscMalloc1(dd->m, &dd->lx);CHKERRQ(ierr); 590 } 591 ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr); 592 } 593 if (ly) { 594 if (dd->n < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 595 ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr); 596 if (!dd->ly) { 597 ierr = PetscMalloc1(dd->n, &dd->ly);CHKERRQ(ierr); 598 } 599 ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr); 600 } 601 if (lz) { 602 if (dd->p < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 603 ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr); 604 if (!dd->lz) { 605 ierr = PetscMalloc1(dd->p, &dd->lz);CHKERRQ(ierr); 606 } 607 ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr); 608 } 609 PetscFunctionReturn(0); 610 } 611 612 /*@ 613 DMDASetInterpolationType - Sets the type of interpolation that will be 614 returned by DMCreateInterpolation() 615 616 Logically Collective on DMDA 617 618 Input Parameter: 619 + da - initial distributed array 620 . ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms 621 622 Level: intermediate 623 624 Notes: 625 you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation() 626 627 .keywords: distributed array, interpolation 628 629 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType 630 @*/ 631 PetscErrorCode DMDASetInterpolationType(DM da,DMDAInterpolationType ctype) 632 { 633 DM_DA *dd = (DM_DA*)da->data; 634 635 PetscFunctionBegin; 636 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 637 PetscValidLogicalCollectiveEnum(da,ctype,2); 638 dd->interptype = ctype; 639 PetscFunctionReturn(0); 640 } 641 642 /*@ 643 DMDAGetInterpolationType - Gets the type of interpolation that will be 644 used by DMCreateInterpolation() 645 646 Not Collective 647 648 Input Parameter: 649 . da - distributed array 650 651 Output Parameter: 652 . ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms) 653 654 Level: intermediate 655 656 .keywords: distributed array, interpolation 657 658 .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation() 659 @*/ 660 PetscErrorCode DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype) 661 { 662 DM_DA *dd = (DM_DA*)da->data; 663 664 PetscFunctionBegin; 665 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 666 PetscValidPointer(ctype,2); 667 *ctype = dd->interptype; 668 PetscFunctionReturn(0); 669 } 670 671 /*@C 672 DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 673 processes neighbors. 674 675 Not Collective 676 677 Input Parameter: 678 . da - the DMDA object 679 680 Output Parameters: 681 . ranks - the neighbors ranks, stored with the x index increasing most rapidly. 682 this process itself is in the list 683 684 Notes: 685 In 2d the array is of length 9, in 3d of length 27 686 Not supported in 1d 687 Do not free the array, it is freed when the DMDA is destroyed. 688 689 Fortran Notes: 690 In fortran you must pass in an array of the appropriate length. 691 692 Level: intermediate 693 694 @*/ 695 PetscErrorCode DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[]) 696 { 697 DM_DA *dd = (DM_DA*)da->data; 698 699 PetscFunctionBegin; 700 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 701 *ranks = dd->neighbors; 702 PetscFunctionReturn(0); 703 } 704 705 /*@C 706 DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 707 708 Not Collective 709 710 Input Parameter: 711 . da - the DMDA object 712 713 Output Parameter: 714 + lx - ownership along x direction (optional) 715 . ly - ownership along y direction (optional) 716 - lz - ownership along z direction (optional) 717 718 Level: intermediate 719 720 Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d() 721 722 In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and 723 eighth arguments from DMDAGetInfo() 724 725 In C you should not free these arrays, nor change the values in them. They will only have valid values while the 726 DMDA they came from still exists (has not been destroyed). 727 728 These numbers are NOT multiplied by the number of dof per node. 729 730 Not available from Fortran 731 732 .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges() 733 @*/ 734 PetscErrorCode DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[]) 735 { 736 DM_DA *dd = (DM_DA*)da->data; 737 738 PetscFunctionBegin; 739 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 740 if (lx) *lx = dd->lx; 741 if (ly) *ly = dd->ly; 742 if (lz) *lz = dd->lz; 743 PetscFunctionReturn(0); 744 } 745 746 /*@ 747 DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined 748 749 Logically Collective on DMDA 750 751 Input Parameters: 752 + da - the DMDA object 753 . refine_x - ratio of fine grid to coarse in x direction (2 by default) 754 . refine_y - ratio of fine grid to coarse in y direction (2 by default) 755 - refine_z - ratio of fine grid to coarse in z direction (2 by default) 756 757 Options Database: 758 + -da_refine_x - refinement ratio in x direction 759 . -da_refine_y - refinement ratio in y direction 760 - -da_refine_z - refinement ratio in z direction 761 762 Level: intermediate 763 764 Notes: 765 Pass PETSC_IGNORE to leave a value unchanged 766 767 .seealso: DMRefine(), DMDAGetRefinementFactor() 768 @*/ 769 PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z) 770 { 771 DM_DA *dd = (DM_DA*)da->data; 772 773 PetscFunctionBegin; 774 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 775 PetscValidLogicalCollectiveInt(da,refine_x,2); 776 PetscValidLogicalCollectiveInt(da,refine_y,3); 777 PetscValidLogicalCollectiveInt(da,refine_z,4); 778 779 if (refine_x > 0) dd->refine_x = refine_x; 780 if (refine_y > 0) dd->refine_y = refine_y; 781 if (refine_z > 0) dd->refine_z = refine_z; 782 PetscFunctionReturn(0); 783 } 784 785 /*@C 786 DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined 787 788 Not Collective 789 790 Input Parameter: 791 . da - the DMDA object 792 793 Output Parameters: 794 + refine_x - ratio of fine grid to coarse in x direction (2 by default) 795 . refine_y - ratio of fine grid to coarse in y direction (2 by default) 796 - refine_z - ratio of fine grid to coarse in z direction (2 by default) 797 798 Level: intermediate 799 800 Notes: 801 Pass NULL for values you do not need 802 803 .seealso: DMRefine(), DMDASetRefinementFactor() 804 @*/ 805 PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z) 806 { 807 DM_DA *dd = (DM_DA*)da->data; 808 809 PetscFunctionBegin; 810 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 811 if (refine_x) *refine_x = dd->refine_x; 812 if (refine_y) *refine_y = dd->refine_y; 813 if (refine_z) *refine_z = dd->refine_z; 814 PetscFunctionReturn(0); 815 } 816 817 /*@C 818 DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix. 819 820 Logically Collective on DMDA 821 822 Input Parameters: 823 + da - the DMDA object 824 - f - the function that allocates the matrix for that specific DMDA 825 826 Level: developer 827 828 Notes: 829 See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for 830 the diagonal and off-diagonal blocks of the matrix 831 832 Not supported from Fortran 833 834 .seealso: DMCreateMatrix(), DMDASetBlockFills() 835 @*/ 836 PetscErrorCode DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, Mat*)) 837 { 838 PetscFunctionBegin; 839 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 840 da->ops->creatematrix = f; 841 PetscFunctionReturn(0); 842 } 843 844 /* 845 Creates "balanced" ownership ranges after refinement, constrained by the need for the 846 fine grid boundaries to fall within one stencil width of the coarse partition. 847 848 Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 849 */ 850 static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[]) 851 { 852 PetscInt i,totalc = 0,remaining,startc = 0,startf = 0; 853 PetscErrorCode ierr; 854 855 PetscFunctionBegin; 856 if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 857 if (ratio == 1) { 858 ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr); 859 PetscFunctionReturn(0); 860 } 861 for (i=0; i<m; i++) totalc += lc[i]; 862 remaining = (!periodic) + ratio * (totalc - (!periodic)); 863 for (i=0; i<m; i++) { 864 PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 865 if (i == m-1) lf[i] = want; 866 else { 867 const PetscInt nextc = startc + lc[i]; 868 /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one 869 * coarse stencil width of the first coarse node in the next subdomain. */ 870 while ((startf+want)/ratio < nextc - stencil_width) want++; 871 /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one 872 * coarse stencil width of the last coarse node in the current subdomain. */ 873 while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--; 874 /* Make sure all constraints are satisfied */ 875 if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width) 876 || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range"); 877 } 878 lf[i] = want; 879 startc += lc[i]; 880 startf += lf[i]; 881 remaining -= lf[i]; 882 } 883 PetscFunctionReturn(0); 884 } 885 886 /* 887 Creates "balanced" ownership ranges after coarsening, constrained by the need for the 888 fine grid boundaries to fall within one stencil width of the coarse partition. 889 890 Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 891 */ 892 static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[]) 893 { 894 PetscInt i,totalf,remaining,startc,startf; 895 PetscErrorCode ierr; 896 897 PetscFunctionBegin; 898 if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 899 if (ratio == 1) { 900 ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr); 901 PetscFunctionReturn(0); 902 } 903 for (i=0,totalf=0; i<m; i++) totalf += lf[i]; 904 remaining = (!periodic) + (totalf - (!periodic)) / ratio; 905 for (i=0,startc=0,startf=0; i<m; i++) { 906 PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 907 if (i == m-1) lc[i] = want; 908 else { 909 const PetscInt nextf = startf+lf[i]; 910 /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 911 * node is within one stencil width. */ 912 while (nextf/ratio < startc+want-stencil_width) want--; 913 /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 914 * fine node is within one stencil width. */ 915 while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++; 916 if (want < 0 || want > remaining 917 || (nextf/ratio < startc+want-stencil_width) || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range"); 918 } 919 lc[i] = want; 920 startc += lc[i]; 921 startf += lf[i]; 922 remaining -= lc[i]; 923 } 924 PetscFunctionReturn(0); 925 } 926 927 PetscErrorCode DMRefine_DA(DM da,MPI_Comm comm,DM *daref) 928 { 929 PetscErrorCode ierr; 930 PetscInt M,N,P,i,dim; 931 DM da2; 932 DM_DA *dd = (DM_DA*)da->data,*dd2; 933 934 PetscFunctionBegin; 935 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 936 PetscValidPointer(daref,3); 937 938 ierr = DMGetDimension(da, &dim);CHKERRQ(ierr); 939 if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 940 M = dd->refine_x*dd->M; 941 } else { 942 M = 1 + dd->refine_x*(dd->M - 1); 943 } 944 if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 945 if (dim > 1) { 946 N = dd->refine_y*dd->N; 947 } else { 948 N = 1; 949 } 950 } else { 951 N = 1 + dd->refine_y*(dd->N - 1); 952 } 953 if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 954 if (dim > 2) { 955 P = dd->refine_z*dd->P; 956 } else { 957 P = 1; 958 } 959 } else { 960 P = 1 + dd->refine_z*(dd->P - 1); 961 } 962 ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr); 963 ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr); 964 ierr = DMSetDimension(da2,dim);CHKERRQ(ierr); 965 ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr); 966 ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr); 967 ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); 968 ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr); 969 ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr); 970 ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr); 971 if (dim == 3) { 972 PetscInt *lx,*ly,*lz; 973 ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr); 974 ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 975 ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 976 ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 977 ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr); 978 ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 979 } else if (dim == 2) { 980 PetscInt *lx,*ly; 981 ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr); 982 ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 983 ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 984 ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr); 985 ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 986 } else if (dim == 1) { 987 PetscInt *lx; 988 ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr); 989 ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 990 ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr); 991 ierr = PetscFree(lx);CHKERRQ(ierr); 992 } 993 dd2 = (DM_DA*)da2->data; 994 995 /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 996 da2->ops->creatematrix = da->ops->creatematrix; 997 /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */ 998 da2->ops->getcoloring = da->ops->getcoloring; 999 dd2->interptype = dd->interptype; 1000 1001 /* copy fill information if given */ 1002 if (dd->dfill) { 1003 ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr); 1004 ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 1005 } 1006 if (dd->ofill) { 1007 ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr); 1008 ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 1009 } 1010 /* copy the refine information */ 1011 dd2->coarsen_x = dd2->refine_x = dd->refine_x; 1012 dd2->coarsen_y = dd2->refine_y = dd->refine_y; 1013 dd2->coarsen_z = dd2->refine_z = dd->refine_z; 1014 1015 if (dd->refine_z_hier) { 1016 if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_z_hier_n) { 1017 dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown + 1]; 1018 } 1019 if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_z_hier_n) { 1020 dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown]; 1021 } 1022 dd2->refine_z_hier_n = dd->refine_z_hier_n; 1023 ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr); 1024 ierr = PetscMemcpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1025 } 1026 if (dd->refine_y_hier) { 1027 if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_y_hier_n) { 1028 dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown + 1]; 1029 } 1030 if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_y_hier_n) { 1031 dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown]; 1032 } 1033 dd2->refine_y_hier_n = dd->refine_y_hier_n; 1034 ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr); 1035 ierr = PetscMemcpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1036 } 1037 if (dd->refine_x_hier) { 1038 if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_x_hier_n) { 1039 dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown + 1]; 1040 } 1041 if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_x_hier_n) { 1042 dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown]; 1043 } 1044 dd2->refine_x_hier_n = dd->refine_x_hier_n; 1045 ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr); 1046 ierr = PetscMemcpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1047 } 1048 1049 1050 /* copy vector type information */ 1051 ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr); 1052 1053 dd2->lf = dd->lf; 1054 dd2->lj = dd->lj; 1055 1056 da2->leveldown = da->leveldown; 1057 da2->levelup = da->levelup + 1; 1058 1059 ierr = DMSetUp(da2);CHKERRQ(ierr); 1060 1061 /* interpolate coordinates if they are set on the coarse grid */ 1062 if (da->coordinates) { 1063 DM cdaf,cdac; 1064 Vec coordsc,coordsf; 1065 Mat II; 1066 1067 ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr); 1068 ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr); 1069 ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr); 1070 /* force creation of the coordinate vector */ 1071 ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 1072 ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr); 1073 ierr = DMCreateInterpolation(cdac,cdaf,&II,NULL);CHKERRQ(ierr); 1074 ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr); 1075 ierr = MatDestroy(&II);CHKERRQ(ierr); 1076 } 1077 1078 for (i=0; i<da->bs; i++) { 1079 const char *fieldname; 1080 ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 1081 ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 1082 } 1083 1084 *daref = da2; 1085 PetscFunctionReturn(0); 1086 } 1087 1088 1089 PetscErrorCode DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref) 1090 { 1091 PetscErrorCode ierr; 1092 PetscInt M,N,P,i,dim; 1093 DM da2; 1094 DM_DA *dd = (DM_DA*)da->data,*dd2; 1095 1096 PetscFunctionBegin; 1097 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 1098 PetscValidPointer(daref,3); 1099 1100 ierr = DMGetDimension(da, &dim);CHKERRQ(ierr); 1101 if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1102 M = dd->M / dd->coarsen_x; 1103 } else { 1104 M = 1 + (dd->M - 1) / dd->coarsen_x; 1105 } 1106 if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1107 if (dim > 1) { 1108 N = dd->N / dd->coarsen_y; 1109 } else { 1110 N = 1; 1111 } 1112 } else { 1113 N = 1 + (dd->N - 1) / dd->coarsen_y; 1114 } 1115 if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 1116 if (dim > 2) { 1117 P = dd->P / dd->coarsen_z; 1118 } else { 1119 P = 1; 1120 } 1121 } else { 1122 P = 1 + (dd->P - 1) / dd->coarsen_z; 1123 } 1124 ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr); 1125 ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr); 1126 ierr = DMSetDimension(da2,dim);CHKERRQ(ierr); 1127 ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr); 1128 ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr); 1129 ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); 1130 ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr); 1131 ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr); 1132 ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr); 1133 if (dim == 3) { 1134 PetscInt *lx,*ly,*lz; 1135 ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr); 1136 ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 1137 ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 1138 ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 1139 ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr); 1140 ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 1141 } else if (dim == 2) { 1142 PetscInt *lx,*ly; 1143 ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr); 1144 ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 1145 ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 1146 ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr); 1147 ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 1148 } else if (dim == 1) { 1149 PetscInt *lx; 1150 ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr); 1151 ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 1152 ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr); 1153 ierr = PetscFree(lx);CHKERRQ(ierr); 1154 } 1155 dd2 = (DM_DA*)da2->data; 1156 1157 /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */ 1158 /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */ 1159 da2->ops->creatematrix = da->ops->creatematrix; 1160 da2->ops->getcoloring = da->ops->getcoloring; 1161 dd2->interptype = dd->interptype; 1162 1163 /* copy fill information if given */ 1164 if (dd->dfill) { 1165 ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr); 1166 ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 1167 } 1168 if (dd->ofill) { 1169 ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr); 1170 ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 1171 } 1172 /* copy the refine information */ 1173 dd2->coarsen_x = dd2->refine_x = dd->coarsen_x; 1174 dd2->coarsen_y = dd2->refine_y = dd->coarsen_y; 1175 dd2->coarsen_z = dd2->refine_z = dd->coarsen_z; 1176 1177 if (dd->refine_z_hier) { 1178 if (da->levelup - da->leveldown -1 > -1 && da->levelup - da->leveldown - 1< dd->refine_z_hier_n) { 1179 dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown - 1]; 1180 } 1181 if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_z_hier_n) { 1182 dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown - 2]; 1183 } 1184 dd2->refine_z_hier_n = dd->refine_z_hier_n; 1185 ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr); 1186 ierr = PetscMemcpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1187 } 1188 if (dd->refine_y_hier) { 1189 if (da->levelup - da->leveldown - 1 > -1 && da->levelup - da->leveldown - 1< dd->refine_y_hier_n) { 1190 dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown - 1]; 1191 } 1192 if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_y_hier_n) { 1193 dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown - 2]; 1194 } 1195 dd2->refine_y_hier_n = dd->refine_y_hier_n; 1196 ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr); 1197 ierr = PetscMemcpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1198 } 1199 if (dd->refine_x_hier) { 1200 if (da->levelup - da->leveldown - 1 > -1 && da->levelup - da->leveldown - 1 < dd->refine_x_hier_n) { 1201 dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown - 1]; 1202 } 1203 if (da->levelup - da->leveldown - 2 > -1 && da->levelup - da->leveldown - 2 < dd->refine_x_hier_n) { 1204 dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown - 2]; 1205 } 1206 dd2->refine_x_hier_n = dd->refine_x_hier_n; 1207 ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr); 1208 ierr = PetscMemcpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n*sizeof(PetscInt));CHKERRQ(ierr); 1209 } 1210 1211 /* copy vector type information */ 1212 ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr); 1213 1214 dd2->lf = dd->lf; 1215 dd2->lj = dd->lj; 1216 1217 da2->leveldown = da->leveldown + 1; 1218 da2->levelup = da->levelup; 1219 1220 ierr = DMSetUp(da2);CHKERRQ(ierr); 1221 1222 /* inject coordinates if they are set on the fine grid */ 1223 if (da->coordinates) { 1224 DM cdaf,cdac; 1225 Vec coordsc,coordsf; 1226 Mat inject; 1227 VecScatter vscat; 1228 1229 ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr); 1230 ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr); 1231 ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr); 1232 /* force creation of the coordinate vector */ 1233 ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 1234 ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr); 1235 1236 ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr); 1237 ierr = MatScatterGetVecScatter(inject,&vscat);CHKERRQ(ierr); 1238 ierr = VecScatterBegin(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1239 ierr = VecScatterEnd(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1240 ierr = MatDestroy(&inject);CHKERRQ(ierr); 1241 } 1242 1243 for (i=0; i<da->bs; i++) { 1244 const char *fieldname; 1245 ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 1246 ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 1247 } 1248 1249 *daref = da2; 1250 PetscFunctionReturn(0); 1251 } 1252 1253 PetscErrorCode DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[]) 1254 { 1255 PetscErrorCode ierr; 1256 PetscInt i,n,*refx,*refy,*refz; 1257 1258 PetscFunctionBegin; 1259 PetscValidHeaderSpecific(da,DM_CLASSID,1); 1260 if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 1261 if (nlevels == 0) PetscFunctionReturn(0); 1262 PetscValidPointer(daf,3); 1263 1264 /* Get refinement factors, defaults taken from the coarse DMDA */ 1265 ierr = PetscMalloc3(nlevels,&refx,nlevels,&refy,nlevels,&refz);CHKERRQ(ierr); 1266 for (i=0; i<nlevels; i++) { 1267 ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr); 1268 } 1269 n = nlevels; 1270 ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr); 1271 n = nlevels; 1272 ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr); 1273 n = nlevels; 1274 ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr); 1275 1276 ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr); 1277 ierr = DMRefine(da,PetscObjectComm((PetscObject)da),&daf[0]);CHKERRQ(ierr); 1278 for (i=1; i<nlevels; i++) { 1279 ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr); 1280 ierr = DMRefine(daf[i-1],PetscObjectComm((PetscObject)da),&daf[i]);CHKERRQ(ierr); 1281 } 1282 ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr); 1283 PetscFunctionReturn(0); 1284 } 1285 1286 PetscErrorCode DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[]) 1287 { 1288 PetscErrorCode ierr; 1289 PetscInt i; 1290 1291 PetscFunctionBegin; 1292 PetscValidHeaderSpecific(da,DM_CLASSID,1); 1293 if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 1294 if (nlevels == 0) PetscFunctionReturn(0); 1295 PetscValidPointer(dac,3); 1296 ierr = DMCoarsen(da,PetscObjectComm((PetscObject)da),&dac[0]);CHKERRQ(ierr); 1297 for (i=1; i<nlevels; i++) { 1298 ierr = DMCoarsen(dac[i-1],PetscObjectComm((PetscObject)da),&dac[i]);CHKERRQ(ierr); 1299 } 1300 PetscFunctionReturn(0); 1301 } 1302 1303 #include <petscgll.h> 1304 1305 PetscErrorCode DMDASetGLLCoordinates_1d(DM dm,PetscGLL *gll) 1306 { 1307 PetscErrorCode ierr; 1308 PetscInt i,j,n = gll->n,xs,xn,q; 1309 PetscScalar *xx; 1310 PetscReal h; 1311 Vec x; 1312 DM_DA *da = (DM_DA*)dm->data; 1313 1314 PetscFunctionBegin; 1315 if (da->bx != DM_BOUNDARY_PERIODIC) { 1316 ierr = DMDAGetInfo(dm,NULL,&q,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 1317 q = (q-1)/(n-1); /* number of spectral elements */ 1318 h = 2.0/q; 1319 ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr); 1320 xs = xs/(n-1); 1321 xn = xn/(n-1); 1322 ierr = DMDASetUniformCoordinates(dm,-1.,1.,0.,0.,0.,0.);CHKERRQ(ierr); 1323 ierr = DMGetCoordinates(dm,&x);CHKERRQ(ierr); 1324 ierr = DMDAVecGetArray(dm,x,&xx);CHKERRQ(ierr); 1325 1326 /* loop over local spectral elements */ 1327 for (j=xs; j<xs+xn; j++) { 1328 /* 1329 Except for the first process, each process starts on the second GLL point of the first element on that process 1330 */ 1331 for (i= (j == xs && xs > 0)? 1 : 0; i<n; i++) { 1332 xx[j*(n-1) + i] = -1.0 + h*j + h*(gll->nodes[i]+1.0)/2.; 1333 } 1334 } 1335 ierr = DMDAVecRestoreArray(dm,x,&xx);CHKERRQ(ierr); 1336 } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not yet implemented for periodic"); 1337 PetscFunctionReturn(0); 1338 } 1339 1340 /*@ 1341 1342 DMDASetGLLCoordinates - Sets the global coordinates from -1 to 1 to the GLL points of as many GLL elements that fit the number of grid points 1343 1344 Collective on DM 1345 1346 Input Parameters: 1347 + da - the DMDA object 1348 - gll - the GLL object 1349 1350 Notes: 1351 the parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points 1352 on each process much be divisible by the number of GLL elements needed per process. This depends on whether the DM is 1353 periodic or not. 1354 1355 Level: advanced 1356 1357 .seealso: DMDACreate(), PetscGLLCreate(), DMGetCoordinates() 1358 @*/ 1359 PetscErrorCode DMDASetGLLCoordinates(DM da,PetscGLL *gll) 1360 { 1361 PetscErrorCode ierr; 1362 1363 PetscFunctionBegin; 1364 if (da->dim == 1) { 1365 ierr = DMDASetGLLCoordinates_1d(da,gll);CHKERRQ(ierr); 1366 } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not yet implemented for 2 or 3d"); 1367 PetscFunctionReturn(0); 1368 } 1369 1370 PETSC_INTERN PetscErrorCode DMGetCompatibility_DA(DM da1,DM dm2,PetscBool *compatible,PetscBool *set) 1371 { 1372 PetscErrorCode ierr; 1373 DM_DA *dd1 = (DM_DA*)da1->data,*dd2; 1374 DM da2; 1375 DMType dmtype2; 1376 PetscBool isda,compatibleLocal; 1377 PetscInt i; 1378 1379 PetscFunctionBegin; 1380 if (!da1->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da1),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() must be called on first DM before DMGetCompatibility()"); 1381 ierr = DMGetType(dm2,&dmtype2);CHKERRQ(ierr); 1382 ierr = PetscStrcmp(dmtype2,DMDA,&isda);CHKERRQ(ierr); 1383 if (isda) { 1384 da2 = dm2; 1385 dd2 = (DM_DA*)da2->data; 1386 if (!da2->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da2),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() must be called on second DM before DMGetCompatibility()"); 1387 compatibleLocal = (PetscBool)(da1->dim == da2->dim); 1388 if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->s == dd2->s)); /* Stencil width */ 1389 /* Global size ranks Boundary type */ 1390 if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->M == dd2->M) && (dd1->m == dd2->m) && (dd1->bx == dd2->bx)); 1391 if (compatibleLocal && da1->dim > 1) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->N == dd2->N) && (dd1->n == dd2->n) && (dd1->by == dd2->by)); 1392 if (compatibleLocal && da1->dim > 2) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->P == dd2->P) && (dd1->p == dd2->p) && (dd1->bz == dd2->bz)); 1393 if (compatibleLocal) { 1394 for (i=0; i<dd1->m; ++i) { 1395 compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lx[i] == dd2->lx[i])); /* Local size */ 1396 } 1397 } 1398 if (compatibleLocal && da1->dim > 1) { 1399 for (i=0; i<dd1->n; ++i) { 1400 compatibleLocal = (PetscBool)(compatibleLocal && (dd1->ly[i] == dd2->ly[i])); 1401 } 1402 } 1403 if (compatibleLocal && da1->dim > 2) { 1404 for (i=0; i<dd1->p; ++i) { 1405 compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lz[i] == dd2->lz[i])); 1406 } 1407 } 1408 *compatible = compatibleLocal; 1409 *set = PETSC_TRUE; 1410 } else { 1411 /* Decline to determine compatibility with other DM types */ 1412 *set = PETSC_FALSE; 1413 } 1414 PetscFunctionReturn(0); 1415 } 1416