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