1 #include <private/daimpl.h> /*I "petscdmda.h" I*/ 2 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "DMDASetDim" 6 /*@ 7 DMDASetDim - Sets the dimension 8 9 Collective on DMDA 10 11 Input Parameters: 12 + da - the DMDA 13 - dim - the dimension (or PETSC_DECIDE) 14 15 Level: intermediate 16 17 .seealso: DaGetDim(), DMDASetSizes() 18 @*/ 19 PetscErrorCode DMDASetDim(DM da, PetscInt dim) 20 { 21 DM_DA *dd = (DM_DA*)da->data; 22 23 PetscFunctionBegin; 24 PetscValidHeaderSpecific(da, DM_CLASSID, 1); 25 if (dd->dim > 0 && dim != dd->dim) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change DMDA dim from %D after it was set to %D",dd->dim,dim); 26 dd->dim = dim; 27 PetscFunctionReturn(0); 28 } 29 30 #undef __FUNCT__ 31 #define __FUNCT__ "DMDASetSizes" 32 /*@ 33 DMDASetSizes - Sets the global sizes 34 35 Logically Collective on DMDA 36 37 Input Parameters: 38 + da - the DMDA 39 . M - the global X size (or PETSC_DECIDE) 40 . N - the global Y size (or PETSC_DECIDE) 41 - P - the global Z size (or PETSC_DECIDE) 42 43 Level: intermediate 44 45 .seealso: DMDAGetSize(), PetscSplitOwnership() 46 @*/ 47 PetscErrorCode DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P) 48 { 49 DM_DA *dd = (DM_DA*)da->data; 50 51 PetscFunctionBegin; 52 PetscValidHeaderSpecific(da, DM_CLASSID, 1); 53 PetscValidLogicalCollectiveInt(da,M,2); 54 PetscValidLogicalCollectiveInt(da,N,3); 55 PetscValidLogicalCollectiveInt(da,P,4); 56 57 dd->M = M; 58 dd->N = N; 59 dd->P = P; 60 PetscFunctionReturn(0); 61 } 62 63 #undef __FUNCT__ 64 #define __FUNCT__ "DMDASetNumProcs" 65 /*@ 66 DMDASetNumProcs - Sets the number of processes in each dimension 67 68 Logically Collective on DMDA 69 70 Input Parameters: 71 + da - the DMDA 72 . m - the number of X procs (or PETSC_DECIDE) 73 . n - the number of Y procs (or PETSC_DECIDE) 74 - p - the number of Z procs (or PETSC_DECIDE) 75 76 Level: intermediate 77 78 .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership() 79 @*/ 80 PetscErrorCode DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p) 81 { 82 DM_DA *dd = (DM_DA*)da->data; 83 84 PetscFunctionBegin; 85 PetscValidHeaderSpecific(da, DM_CLASSID, 1); 86 PetscValidLogicalCollectiveInt(da,m,2); 87 PetscValidLogicalCollectiveInt(da,n,3); 88 PetscValidLogicalCollectiveInt(da,p,4); 89 dd->m = m; 90 dd->n = n; 91 dd->p = p; 92 PetscFunctionReturn(0); 93 } 94 95 #undef __FUNCT__ 96 #define __FUNCT__ "DMDASetBoundaryType" 97 /*@ 98 DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries. 99 100 Not collective 101 102 Input Parameter: 103 + da - The DMDA 104 - bx,by,bz - One of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC 105 106 Level: intermediate 107 108 .keywords: distributed array, periodicity 109 .seealso: DMDACreate(), DMDestroy(), DMDA, DMDABoundaryType 110 @*/ 111 PetscErrorCode DMDASetBoundaryType(DM da,DMDABoundaryType bx,DMDABoundaryType by,DMDABoundaryType bz) 112 { 113 DM_DA *dd = (DM_DA*)da->data; 114 115 PetscFunctionBegin; 116 PetscValidHeaderSpecific(da,DM_CLASSID,1); 117 PetscValidLogicalCollectiveEnum(da,bx,2); 118 PetscValidLogicalCollectiveEnum(da,by,3); 119 PetscValidLogicalCollectiveEnum(da,bz,4); 120 dd->bx = bx; 121 dd->by = by; 122 dd->bz = bz; 123 PetscFunctionReturn(0); 124 } 125 126 #undef __FUNCT__ 127 #define __FUNCT__ "DMDASetDof" 128 /*@ 129 DMDASetDof - Sets the number of degrees of freedom per vertex 130 131 Not collective 132 133 Input Parameter: 134 + da - The DMDA 135 - dof - Number of degrees of freedom 136 137 Level: intermediate 138 139 .keywords: distributed array, degrees of freedom 140 .seealso: DMDACreate(), DMDestroy(), DMDA 141 @*/ 142 PetscErrorCode DMDASetDof(DM da, PetscInt dof) 143 { 144 DM_DA *dd = (DM_DA*)da->data; 145 146 PetscFunctionBegin; 147 PetscValidHeaderSpecific(da,DM_CLASSID,1); 148 PetscValidLogicalCollectiveInt(da,dof,2); 149 dd->w = dof; 150 da->bs = dof; 151 PetscFunctionReturn(0); 152 } 153 154 #undef __FUNCT__ 155 #define __FUNCT__ "DMDASetStencilType" 156 /*@ 157 DMDASetStencilType - Sets the type of the communication stencil 158 159 Logically Collective on DMDA 160 161 Input Parameter: 162 + da - The DMDA 163 - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR. 164 165 Level: intermediate 166 167 .keywords: distributed array, stencil 168 .seealso: DMDACreate(), DMDestroy(), DMDA 169 @*/ 170 PetscErrorCode DMDASetStencilType(DM da, DMDAStencilType stype) 171 { 172 DM_DA *dd = (DM_DA*)da->data; 173 174 PetscFunctionBegin; 175 PetscValidHeaderSpecific(da,DM_CLASSID,1); 176 PetscValidLogicalCollectiveEnum(da,stype,2); 177 dd->stencil_type = stype; 178 PetscFunctionReturn(0); 179 } 180 181 #undef __FUNCT__ 182 #define __FUNCT__ "DMDASetStencilWidth" 183 /*@ 184 DMDASetStencilWidth - Sets the width of the communication stencil 185 186 Logically Collective on DMDA 187 188 Input Parameter: 189 + da - The DMDA 190 - width - The stencil width 191 192 Level: intermediate 193 194 .keywords: distributed array, stencil 195 .seealso: DMDACreate(), DMDestroy(), DMDA 196 @*/ 197 PetscErrorCode DMDASetStencilWidth(DM da, PetscInt width) 198 { 199 DM_DA *dd = (DM_DA*)da->data; 200 201 PetscFunctionBegin; 202 PetscValidHeaderSpecific(da,DM_CLASSID,1); 203 PetscValidLogicalCollectiveInt(da,width,2); 204 dd->s = width; 205 PetscFunctionReturn(0); 206 } 207 208 #undef __FUNCT__ 209 #define __FUNCT__ "DMDACheckOwnershipRanges_Private" 210 static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[]) 211 { 212 PetscInt i,sum; 213 214 PetscFunctionBegin; 215 if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set"); 216 for (i=sum=0; i<m; i++) sum += lx[i]; 217 if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M); 218 PetscFunctionReturn(0); 219 } 220 221 #undef __FUNCT__ 222 #define __FUNCT__ "DMDASetOwnershipRanges" 223 /*@ 224 DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process 225 226 Logically Collective on DMDA 227 228 Input Parameter: 229 + da - The DMDA 230 . lx - array containing number of nodes in the X direction on each process, or PETSC_NULL. If non-null, must be of length da->m 231 . ly - array containing number of nodes in the Y direction on each process, or PETSC_NULL. If non-null, must be of length da->n 232 - lz - array containing number of nodes in the Z direction on each process, or PETSC_NULL. If non-null, must be of length da->p. 233 234 Level: intermediate 235 236 .keywords: distributed array 237 .seealso: DMDACreate(), DMDestroy(), DMDA 238 @*/ 239 PetscErrorCode DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[]) 240 { 241 PetscErrorCode ierr; 242 DM_DA *dd = (DM_DA*)da->data; 243 244 PetscFunctionBegin; 245 PetscValidHeaderSpecific(da,DM_CLASSID,1); 246 if (lx) { 247 if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 248 ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr); 249 if (!dd->lx) { 250 ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr); 251 } 252 ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr); 253 } 254 if (ly) { 255 if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 256 ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr); 257 if (!dd->ly) { 258 ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr); 259 } 260 ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr); 261 } 262 if (lz) { 263 if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs"); 264 ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr); 265 if (!dd->lz) { 266 ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr); 267 } 268 ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr); 269 } 270 PetscFunctionReturn(0); 271 } 272 273 #undef __FUNCT__ 274 #define __FUNCT__ "DMDASetInterpolationType" 275 /*@ 276 DMDASetInterpolationType - Sets the type of interpolation that will be 277 returned by DMCreateInterpolation() 278 279 Logically Collective on DMDA 280 281 Input Parameter: 282 + da - initial distributed array 283 . ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms 284 285 Level: intermediate 286 287 Notes: you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation() 288 289 .keywords: distributed array, interpolation 290 291 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType 292 @*/ 293 PetscErrorCode DMDASetInterpolationType(DM da,DMDAInterpolationType ctype) 294 { 295 DM_DA *dd = (DM_DA*)da->data; 296 297 PetscFunctionBegin; 298 PetscValidHeaderSpecific(da,DM_CLASSID,1); 299 PetscValidLogicalCollectiveEnum(da,ctype,2); 300 dd->interptype = ctype; 301 PetscFunctionReturn(0); 302 } 303 304 #undef __FUNCT__ 305 #define __FUNCT__ "DMDAGetInterpolationType" 306 /*@ 307 DMDAGetInterpolationType - Gets the type of interpolation that will be 308 used by DMCreateInterpolation() 309 310 Not Collective 311 312 Input Parameter: 313 . da - distributed array 314 315 Output Parameter: 316 . ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms) 317 318 Level: intermediate 319 320 .keywords: distributed array, interpolation 321 322 .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation() 323 @*/ 324 PetscErrorCode DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype) 325 { 326 DM_DA *dd = (DM_DA*)da->data; 327 328 PetscFunctionBegin; 329 PetscValidHeaderSpecific(da,DM_CLASSID,1); 330 PetscValidPointer(ctype,2); 331 *ctype = dd->interptype; 332 PetscFunctionReturn(0); 333 } 334 335 #undef __FUNCT__ 336 #define __FUNCT__ "DMDAGetNeighbors" 337 /*@C 338 DMDAGetNeighbors - Gets an array containing the MPI rank of all the current 339 processes neighbors. 340 341 Not Collective 342 343 Input Parameter: 344 . da - the DMDA object 345 346 Output Parameters: 347 . ranks - the neighbors ranks, stored with the x index increasing most rapidly. 348 this process itself is in the list 349 350 Notes: In 2d the array is of length 9, in 3d of length 27 351 Not supported in 1d 352 Do not free the array, it is freed when the DMDA is destroyed. 353 354 Fortran Notes: In fortran you must pass in an array of the appropriate length. 355 356 Level: intermediate 357 358 @*/ 359 PetscErrorCode DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[]) 360 { 361 DM_DA *dd = (DM_DA*)da->data; 362 PetscFunctionBegin; 363 PetscValidHeaderSpecific(da,DM_CLASSID,1); 364 *ranks = dd->neighbors; 365 PetscFunctionReturn(0); 366 } 367 368 #undef __FUNCT__ 369 #define __FUNCT__ "DMDAGetOwnershipRanges" 370 /*@C 371 DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process 372 373 Not Collective 374 375 Input Parameter: 376 . da - the DMDA object 377 378 Output Parameter: 379 + lx - ownership along x direction (optional) 380 . ly - ownership along y direction (optional) 381 - lz - ownership along z direction (optional) 382 383 Level: intermediate 384 385 Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d() 386 387 In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and 388 eighth arguments from DMDAGetInfo() 389 390 In C you should not free these arrays, nor change the values in them. They will only have valid values while the 391 DMDA they came from still exists (has not been destroyed). 392 393 .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges() 394 @*/ 395 PetscErrorCode DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[]) 396 { 397 DM_DA *dd = (DM_DA*)da->data; 398 399 PetscFunctionBegin; 400 PetscValidHeaderSpecific(da,DM_CLASSID,1); 401 if (lx) *lx = dd->lx; 402 if (ly) *ly = dd->ly; 403 if (lz) *lz = dd->lz; 404 PetscFunctionReturn(0); 405 } 406 407 #undef __FUNCT__ 408 #define __FUNCT__ "DMDASetRefinementFactor" 409 /*@ 410 DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined 411 412 Logically Collective on DMDA 413 414 Input Parameters: 415 + da - the DMDA object 416 . refine_x - ratio of fine grid to coarse in x direction (2 by default) 417 . refine_y - ratio of fine grid to coarse in y direction (2 by default) 418 - refine_z - ratio of fine grid to coarse in z direction (2 by default) 419 420 Options Database: 421 + -da_refine_x - refinement ratio in x direction 422 . -da_refine_y - refinement ratio in y direction 423 - -da_refine_z - refinement ratio in z direction 424 425 Level: intermediate 426 427 Notes: Pass PETSC_IGNORE to leave a value unchanged 428 429 .seealso: DMRefine(), DMDAGetRefinementFactor() 430 @*/ 431 PetscErrorCode DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z) 432 { 433 DM_DA *dd = (DM_DA*)da->data; 434 435 PetscFunctionBegin; 436 PetscValidHeaderSpecific(da,DM_CLASSID,1); 437 PetscValidLogicalCollectiveInt(da,refine_x,2); 438 PetscValidLogicalCollectiveInt(da,refine_y,3); 439 PetscValidLogicalCollectiveInt(da,refine_z,4); 440 441 if (refine_x > 0) dd->refine_x = refine_x; 442 if (refine_y > 0) dd->refine_y = refine_y; 443 if (refine_z > 0) dd->refine_z = refine_z; 444 PetscFunctionReturn(0); 445 } 446 447 #undef __FUNCT__ 448 #define __FUNCT__ "DMDAGetRefinementFactor" 449 /*@C 450 DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined 451 452 Not Collective 453 454 Input Parameter: 455 . da - the DMDA object 456 457 Output Parameters: 458 + refine_x - ratio of fine grid to coarse in x direction (2 by default) 459 . refine_y - ratio of fine grid to coarse in y direction (2 by default) 460 - refine_z - ratio of fine grid to coarse in z direction (2 by default) 461 462 Level: intermediate 463 464 Notes: Pass PETSC_NULL for values you do not need 465 466 .seealso: DMRefine(), DMDASetRefinementFactor() 467 @*/ 468 PetscErrorCode DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z) 469 { 470 DM_DA *dd = (DM_DA*)da->data; 471 472 PetscFunctionBegin; 473 PetscValidHeaderSpecific(da,DM_CLASSID,1); 474 if (refine_x) *refine_x = dd->refine_x; 475 if (refine_y) *refine_y = dd->refine_y; 476 if (refine_z) *refine_z = dd->refine_z; 477 PetscFunctionReturn(0); 478 } 479 480 #undef __FUNCT__ 481 #define __FUNCT__ "DMDASetGetMatrix" 482 /*@C 483 DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix. 484 485 Logically Collective on DMDA 486 487 Input Parameters: 488 + da - the DMDA object 489 - f - the function that allocates the matrix for that specific DMDA 490 491 Level: developer 492 493 Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for 494 the diagonal and off-diagonal blocks of the matrix 495 496 .seealso: DMCreateMatrix(), DMDASetBlockFills() 497 @*/ 498 PetscErrorCode DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, const MatType,Mat*)) 499 { 500 PetscFunctionBegin; 501 PetscValidHeaderSpecific(da,DM_CLASSID,1); 502 da->ops->creatematrix = f; 503 PetscFunctionReturn(0); 504 } 505 506 #undef __FUNCT__ 507 #define __FUNCT__ "DMDARefineOwnershipRanges" 508 /* 509 Creates "balanced" ownership ranges after refinement, constrained by the need for the 510 fine grid boundaries to fall within one stencil width of the coarse partition. 511 512 Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 513 */ 514 static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[]) 515 { 516 PetscInt i,totalc = 0,remaining,startc = 0,startf = 0; 517 PetscErrorCode ierr; 518 519 PetscFunctionBegin; 520 if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 521 if (ratio == 1) { 522 ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr); 523 PetscFunctionReturn(0); 524 } 525 for (i=0; i<m; i++) totalc += lc[i]; 526 remaining = (!periodic) + ratio * (totalc - (!periodic)); 527 for (i=0; i<m; i++) { 528 PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 529 if (i == m-1) lf[i] = want; 530 else { 531 const PetscInt nextc = startc + lc[i]; 532 /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one 533 * coarse stencil width of the first coarse node in the next subdomain. */ 534 while ((startf+want)/ratio < nextc - stencil_width) want++; 535 /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one 536 * coarse stencil width of the last coarse node in the current subdomain. */ 537 while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--; 538 /* Make sure all constraints are satisfied */ 539 if (want < 0 || want > remaining 540 || ((startf+want)/ratio < nextc - stencil_width) 541 || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) 542 SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range"); 543 } 544 lf[i] = want; 545 startc += lc[i]; 546 startf += lf[i]; 547 remaining -= lf[i]; 548 } 549 PetscFunctionReturn(0); 550 } 551 552 #undef __FUNCT__ 553 #define __FUNCT__ "DMDACoarsenOwnershipRanges" 554 /* 555 Creates "balanced" ownership ranges after coarsening, constrained by the need for the 556 fine grid boundaries to fall within one stencil width of the coarse partition. 557 558 Uses a greedy algorithm to handle non-ideal layouts, could probably do something better. 559 */ 560 static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[]) 561 { 562 PetscInt i,totalf,remaining,startc,startf; 563 PetscErrorCode ierr; 564 565 PetscFunctionBegin; 566 if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio); 567 if (ratio == 1) { 568 ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr); 569 PetscFunctionReturn(0); 570 } 571 for (i=0,totalf=0; i<m; i++) totalf += lf[i]; 572 remaining = (!periodic) + (totalf - (!periodic)) / ratio; 573 for (i=0,startc=0,startf=0; i<m; i++) { 574 PetscInt want = remaining/(m-i) + !!(remaining%(m-i)); 575 if (i == m-1) lc[i] = want; 576 else { 577 const PetscInt nextf = startf+lf[i]; 578 /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine 579 * node is within one stencil width. */ 580 while (nextf/ratio < startc+want-stencil_width) want--; 581 /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last 582 * fine node is within one stencil width. */ 583 while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++; 584 if (want < 0 || want > remaining 585 || (nextf/ratio < startc+want-stencil_width) 586 || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) 587 SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range"); 588 } 589 lc[i] = want; 590 startc += lc[i]; 591 startf += lf[i]; 592 remaining -= lc[i]; 593 } 594 PetscFunctionReturn(0); 595 } 596 597 #undef __FUNCT__ 598 #define __FUNCT__ "DMRefine_DA" 599 PetscErrorCode DMRefine_DA(DM da,MPI_Comm comm,DM *daref) 600 { 601 PetscErrorCode ierr; 602 PetscInt M,N,P,i; 603 DM da2; 604 DM_DA *dd = (DM_DA*)da->data,*dd2; 605 606 PetscFunctionBegin; 607 PetscValidHeaderSpecific(da,DM_CLASSID,1); 608 PetscValidPointer(daref,3); 609 610 if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 611 M = dd->refine_x*dd->M; 612 } else { 613 M = 1 + dd->refine_x*(dd->M - 1); 614 } 615 if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 616 N = dd->refine_y*dd->N; 617 } else { 618 N = 1 + dd->refine_y*(dd->N - 1); 619 } 620 if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 621 P = dd->refine_z*dd->P; 622 } else { 623 P = 1 + dd->refine_z*(dd->P - 1); 624 } 625 if (dd->dim == 3) { 626 PetscInt *lx,*ly,*lz; 627 ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr); 628 ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 629 ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 630 ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 631 ierr = DMDACreate3d(((PetscObject)da)->comm,dd->bx,dd->by,dd->bz,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr); 632 ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 633 } else if (dd->dim == 2) { 634 PetscInt *lx,*ly; 635 ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr); 636 ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 637 ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 638 ierr = DMDACreate2d(((PetscObject)da)->comm,dd->bx,dd->by,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr); 639 ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 640 } else if (dd->dim == 1) { 641 PetscInt *lx; 642 ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); 643 ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 644 ierr = DMDACreate1d(((PetscObject)da)->comm,dd->bx,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr); 645 ierr = PetscFree(lx);CHKERRQ(ierr); 646 } 647 dd2 = (DM_DA*)da2->data; 648 649 /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 650 da2->ops->creatematrix = da->ops->creatematrix; 651 /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */ 652 da2->ops->getcoloring = da->ops->getcoloring; 653 dd2->interptype = dd->interptype; 654 655 /* copy fill information if given */ 656 if (dd->dfill) { 657 ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); 658 ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 659 } 660 if (dd->ofill) { 661 ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); 662 ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 663 } 664 /* copy the refine information */ 665 dd2->refine_x = dd->refine_x; 666 dd2->refine_y = dd->refine_y; 667 dd2->refine_z = dd->refine_z; 668 669 /* copy vector type information */ 670 ierr = PetscFree(da2->vectype);CHKERRQ(ierr); 671 ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr); 672 673 dd2->lf = dd->lf; 674 dd2->lj = dd->lj; 675 676 /* interpolate coordinates if they are set on the coarse grid */ 677 if (dd->coordinates) { 678 DM cdaf,cdac; 679 Vec coordsc,coordsf; 680 Mat II; 681 682 ierr = DMDAGetCoordinateDA(da,&cdac);CHKERRQ(ierr); 683 ierr = DMDAGetCoordinates(da,&coordsc);CHKERRQ(ierr); 684 ierr = DMDAGetCoordinateDA(da2,&cdaf);CHKERRQ(ierr); 685 /* force creation of the coordinate vector */ 686 ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 687 ierr = DMDAGetCoordinates(da2,&coordsf);CHKERRQ(ierr); 688 ierr = DMCreateInterpolation(cdac,cdaf,&II,PETSC_NULL);CHKERRQ(ierr); 689 ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr); 690 ierr = MatDestroy(&II);CHKERRQ(ierr); 691 } 692 for (i=0; i<da->bs; i++) { 693 const char *fieldname; 694 ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); 695 ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); 696 } 697 *daref = da2; 698 PetscFunctionReturn(0); 699 } 700 701 #undef __FUNCT__ 702 #define __FUNCT__ "DMCoarsen_DA" 703 PetscErrorCode DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref) 704 { 705 PetscErrorCode ierr; 706 PetscInt M,N,P; 707 DM da2; 708 DM_DA *dd = (DM_DA*)da->data,*dd2; 709 710 PetscFunctionBegin; 711 PetscValidHeaderSpecific(da,DM_CLASSID,1); 712 PetscValidPointer(daref,3); 713 714 if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 715 M = dd->M / dd->refine_x; 716 } else { 717 M = 1 + (dd->M - 1) / dd->refine_x; 718 } 719 if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 720 N = dd->N / dd->refine_y; 721 } else { 722 N = 1 + (dd->N - 1) / dd->refine_y; 723 } 724 if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ 725 P = dd->P / dd->refine_z; 726 } else { 727 P = 1 + (dd->P - 1) / dd->refine_z; 728 } 729 if (dd->dim == 3) { 730 PetscInt *lx,*ly,*lz; 731 ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr); 732 ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 733 ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 734 ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr); 735 ierr = DMDACreate3d(((PetscObject)da)->comm,dd->bx,dd->by,dd->bz,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr); 736 ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); 737 } else if (dd->dim == 2) { 738 PetscInt *lx,*ly; 739 ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr); 740 ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 741 ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr); 742 ierr = DMDACreate2d(((PetscObject)da)->comm,dd->bx,dd->by,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr); 743 ierr = PetscFree2(lx,ly);CHKERRQ(ierr); 744 } else if (dd->dim == 1) { 745 PetscInt *lx; 746 ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); 747 ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr); 748 ierr = DMDACreate1d(((PetscObject)da)->comm,dd->bx,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr); 749 ierr = PetscFree(lx);CHKERRQ(ierr); 750 } 751 dd2 = (DM_DA*)da2->data; 752 753 /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */ 754 /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */ 755 da2->ops->creatematrix = da->ops->creatematrix; 756 da2->ops->getcoloring = da->ops->getcoloring; 757 dd2->interptype = dd->interptype; 758 759 /* copy fill information if given */ 760 if (dd->dfill) { 761 ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); 762 ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 763 } 764 if (dd->ofill) { 765 ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); 766 ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); 767 } 768 /* copy the refine information */ 769 dd2->refine_x = dd->refine_x; 770 dd2->refine_y = dd->refine_y; 771 dd2->refine_z = dd->refine_z; 772 773 /* copy vector type information */ 774 ierr = PetscFree(da2->vectype);CHKERRQ(ierr); 775 ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr); 776 777 dd2->lf = dd->lf; 778 dd2->lj = dd->lj; 779 780 /* inject coordinates if they are set on the fine grid */ 781 if (dd->coordinates) { 782 DM cdaf,cdac; 783 Vec coordsc,coordsf; 784 VecScatter inject; 785 786 ierr = DMDAGetCoordinateDA(da,&cdaf);CHKERRQ(ierr); 787 ierr = DMDAGetCoordinates(da,&coordsf);CHKERRQ(ierr); 788 ierr = DMDAGetCoordinateDA(da2,&cdac);CHKERRQ(ierr); 789 /* force creation of the coordinate vector */ 790 ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 791 ierr = DMDAGetCoordinates(da2,&coordsc);CHKERRQ(ierr); 792 793 ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr); 794 ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 795 ierr = VecScatterEnd(inject ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 796 ierr = VecScatterDestroy(&inject);CHKERRQ(ierr); 797 } 798 *daref = da2; 799 PetscFunctionReturn(0); 800 } 801 802 #undef __FUNCT__ 803 #define __FUNCT__ "DMRefineHierarchy_DA" 804 PetscErrorCode DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[]) 805 { 806 PetscErrorCode ierr; 807 PetscInt i,n,*refx,*refy,*refz; 808 809 PetscFunctionBegin; 810 PetscValidHeaderSpecific(da,DM_CLASSID,1); 811 if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 812 if (nlevels == 0) PetscFunctionReturn(0); 813 PetscValidPointer(daf,3); 814 815 /* Get refinement factors, defaults taken from the coarse DMDA */ 816 ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr); 817 for (i=0; i<nlevels; i++) { 818 ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr); 819 } 820 n = nlevels; 821 ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr); 822 n = nlevels; 823 ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr); 824 n = nlevels; 825 ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr); 826 827 ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr); 828 ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr); 829 for (i=1; i<nlevels; i++) { 830 ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr); 831 ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr); 832 } 833 ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr); 834 PetscFunctionReturn(0); 835 } 836 837 #undef __FUNCT__ 838 #define __FUNCT__ "DMCoarsenHierarchy_DA" 839 PetscErrorCode DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[]) 840 { 841 PetscErrorCode ierr; 842 PetscInt i; 843 844 PetscFunctionBegin; 845 PetscValidHeaderSpecific(da,DM_CLASSID,1); 846 if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 847 if (nlevels == 0) PetscFunctionReturn(0); 848 PetscValidPointer(dac,3); 849 ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr); 850 for (i=1; i<nlevels; i++) { 851 ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr); 852 } 853 PetscFunctionReturn(0); 854 } 855