1 2 /* 3 The PC (preconditioner) interface routines, callable by users. 4 */ 5 #include <petsc-private/pcimpl.h> /*I "petscksp.h" I*/ 6 #include <petscdm.h> 7 8 /* Logging support */ 9 PetscClassId PC_CLASSID; 10 PetscLogEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft; 11 PetscLogEvent PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks, PC_ApplyOnMproc; 12 13 #undef __FUNCT__ 14 #define __FUNCT__ "PCGetDefaultType_Private" 15 PetscErrorCode PCGetDefaultType_Private(PC pc,const char *type[]) 16 { 17 PetscErrorCode ierr; 18 PetscMPIInt size; 19 PetscBool flg1,flg2,set,flg3; 20 21 PetscFunctionBegin; 22 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 23 if (pc->pmat) { 24 PetscErrorCode (*f)(Mat,MatReuse,Mat*); 25 ierr = PetscObjectQueryFunction((PetscObject)pc->pmat,"MatGetDiagonalBlock_C",&f);CHKERRQ(ierr); 26 if (size == 1) { 27 ierr = MatGetFactorAvailable(pc->pmat,"petsc",MAT_FACTOR_ICC,&flg1);CHKERRQ(ierr); 28 ierr = MatGetFactorAvailable(pc->pmat,"petsc",MAT_FACTOR_ILU,&flg2);CHKERRQ(ierr); 29 ierr = MatIsSymmetricKnown(pc->pmat,&set,&flg3);CHKERRQ(ierr); 30 if (flg1 && (!flg2 || (set && flg3))) { 31 *type = PCICC; 32 } else if (flg2) { 33 *type = PCILU; 34 } else if (f) { /* likely is a parallel matrix run on one processor */ 35 *type = PCBJACOBI; 36 } else { 37 *type = PCNONE; 38 } 39 } else { 40 if (f) { 41 *type = PCBJACOBI; 42 } else { 43 *type = PCNONE; 44 } 45 } 46 } else { 47 if (size == 1) { 48 *type = PCILU; 49 } else { 50 *type = PCBJACOBI; 51 } 52 } 53 PetscFunctionReturn(0); 54 } 55 56 #undef __FUNCT__ 57 #define __FUNCT__ "PCReset" 58 /*@ 59 PCReset - Resets a PC context to the pcsetupcalled = 0 state and removes any allocated Vecs and Mats 60 61 Collective on PC 62 63 Input Parameter: 64 . pc - the preconditioner context 65 66 Level: developer 67 68 Notes: This allows a PC to be reused for a different sized linear system but using the same options that have been previously set in the PC 69 70 .keywords: PC, destroy 71 72 .seealso: PCCreate(), PCSetUp() 73 @*/ 74 PetscErrorCode PCReset(PC pc) 75 { 76 PetscErrorCode ierr; 77 78 PetscFunctionBegin; 79 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 80 if (pc->ops->reset) { 81 ierr = (*pc->ops->reset)(pc);CHKERRQ(ierr); 82 } 83 ierr = VecDestroy(&pc->diagonalscaleright);CHKERRQ(ierr); 84 ierr = VecDestroy(&pc->diagonalscaleleft);CHKERRQ(ierr); 85 ierr = MatDestroy(&pc->pmat);CHKERRQ(ierr); 86 ierr = MatDestroy(&pc->mat);CHKERRQ(ierr); 87 88 pc->setupcalled = 0; 89 PetscFunctionReturn(0); 90 } 91 92 #undef __FUNCT__ 93 #define __FUNCT__ "PCDestroy" 94 /*@ 95 PCDestroy - Destroys PC context that was created with PCCreate(). 96 97 Collective on PC 98 99 Input Parameter: 100 . pc - the preconditioner context 101 102 Level: developer 103 104 .keywords: PC, destroy 105 106 .seealso: PCCreate(), PCSetUp() 107 @*/ 108 PetscErrorCode PCDestroy(PC *pc) 109 { 110 PetscErrorCode ierr; 111 112 PetscFunctionBegin; 113 if (!*pc) PetscFunctionReturn(0); 114 PetscValidHeaderSpecific((*pc),PC_CLASSID,1); 115 if (--((PetscObject)(*pc))->refct > 0) {*pc = 0; PetscFunctionReturn(0);} 116 117 ierr = PCReset(*pc);CHKERRQ(ierr); 118 119 /* if memory was published with SAWs then destroy it */ 120 ierr = PetscObjectSAWsViewOff((PetscObject)*pc);CHKERRQ(ierr); 121 if ((*pc)->ops->destroy) {ierr = (*(*pc)->ops->destroy)((*pc));CHKERRQ(ierr);} 122 ierr = DMDestroy(&(*pc)->dm);CHKERRQ(ierr); 123 ierr = PetscHeaderDestroy(pc);CHKERRQ(ierr); 124 PetscFunctionReturn(0); 125 } 126 127 #undef __FUNCT__ 128 #define __FUNCT__ "PCGetDiagonalScale" 129 /*@C 130 PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right 131 scaling as needed by certain time-stepping codes. 132 133 Logically Collective on PC 134 135 Input Parameter: 136 . pc - the preconditioner context 137 138 Output Parameter: 139 . flag - PETSC_TRUE if it applies the scaling 140 141 Level: developer 142 143 Notes: If this returns PETSC_TRUE then the system solved via the Krylov method is 144 $ D M A D^{-1} y = D M b for left preconditioning or 145 $ D A M D^{-1} z = D b for right preconditioning 146 147 .keywords: PC 148 149 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCSetDiagonalScale() 150 @*/ 151 PetscErrorCode PCGetDiagonalScale(PC pc,PetscBool *flag) 152 { 153 PetscFunctionBegin; 154 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 155 PetscValidPointer(flag,2); 156 *flag = pc->diagonalscale; 157 PetscFunctionReturn(0); 158 } 159 160 #undef __FUNCT__ 161 #define __FUNCT__ "PCSetDiagonalScale" 162 /*@ 163 PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right 164 scaling as needed by certain time-stepping codes. 165 166 Logically Collective on PC 167 168 Input Parameters: 169 + pc - the preconditioner context 170 - s - scaling vector 171 172 Level: intermediate 173 174 Notes: The system solved via the Krylov method is 175 $ D M A D^{-1} y = D M b for left preconditioning or 176 $ D A M D^{-1} z = D b for right preconditioning 177 178 PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}. 179 180 .keywords: PC 181 182 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCGetDiagonalScale() 183 @*/ 184 PetscErrorCode PCSetDiagonalScale(PC pc,Vec s) 185 { 186 PetscErrorCode ierr; 187 188 PetscFunctionBegin; 189 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 190 PetscValidHeaderSpecific(s,VEC_CLASSID,2); 191 pc->diagonalscale = PETSC_TRUE; 192 193 ierr = PetscObjectReference((PetscObject)s);CHKERRQ(ierr); 194 ierr = VecDestroy(&pc->diagonalscaleleft);CHKERRQ(ierr); 195 196 pc->diagonalscaleleft = s; 197 198 ierr = VecDuplicate(s,&pc->diagonalscaleright);CHKERRQ(ierr); 199 ierr = VecCopy(s,pc->diagonalscaleright);CHKERRQ(ierr); 200 ierr = VecReciprocal(pc->diagonalscaleright);CHKERRQ(ierr); 201 PetscFunctionReturn(0); 202 } 203 204 #undef __FUNCT__ 205 #define __FUNCT__ "PCDiagonalScaleLeft" 206 /*@ 207 PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes. 208 209 Logically Collective on PC 210 211 Input Parameters: 212 + pc - the preconditioner context 213 . in - input vector 214 + out - scaled vector (maybe the same as in) 215 216 Level: intermediate 217 218 Notes: The system solved via the Krylov method is 219 $ D M A D^{-1} y = D M b for left preconditioning or 220 $ D A M D^{-1} z = D b for right preconditioning 221 222 PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}. 223 224 If diagonal scaling is turned off and in is not out then in is copied to out 225 226 .keywords: PC 227 228 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale() 229 @*/ 230 PetscErrorCode PCDiagonalScaleLeft(PC pc,Vec in,Vec out) 231 { 232 PetscErrorCode ierr; 233 234 PetscFunctionBegin; 235 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 236 PetscValidHeaderSpecific(in,VEC_CLASSID,2); 237 PetscValidHeaderSpecific(out,VEC_CLASSID,3); 238 if (pc->diagonalscale) { 239 ierr = VecPointwiseMult(out,pc->diagonalscaleleft,in);CHKERRQ(ierr); 240 } else if (in != out) { 241 ierr = VecCopy(in,out);CHKERRQ(ierr); 242 } 243 PetscFunctionReturn(0); 244 } 245 246 #undef __FUNCT__ 247 #define __FUNCT__ "PCDiagonalScaleRight" 248 /*@ 249 PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes. 250 251 Logically Collective on PC 252 253 Input Parameters: 254 + pc - the preconditioner context 255 . in - input vector 256 + out - scaled vector (maybe the same as in) 257 258 Level: intermediate 259 260 Notes: The system solved via the Krylov method is 261 $ D M A D^{-1} y = D M b for left preconditioning or 262 $ D A M D^{-1} z = D b for right preconditioning 263 264 PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}. 265 266 If diagonal scaling is turned off and in is not out then in is copied to out 267 268 .keywords: PC 269 270 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale() 271 @*/ 272 PetscErrorCode PCDiagonalScaleRight(PC pc,Vec in,Vec out) 273 { 274 PetscErrorCode ierr; 275 276 PetscFunctionBegin; 277 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 278 PetscValidHeaderSpecific(in,VEC_CLASSID,2); 279 PetscValidHeaderSpecific(out,VEC_CLASSID,3); 280 if (pc->diagonalscale) { 281 ierr = VecPointwiseMult(out,pc->diagonalscaleright,in);CHKERRQ(ierr); 282 } else if (in != out) { 283 ierr = VecCopy(in,out);CHKERRQ(ierr); 284 } 285 PetscFunctionReturn(0); 286 } 287 288 #undef __FUNCT__ 289 #define __FUNCT__ "PCSetUseAmat" 290 /*@ 291 PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the 292 operator during the preconditioning process it applies the Amat provided to TSSetRHSJacobian(), 293 TSSetIJacobian(), SNESSetJacobian(), KSPSetOperator() or PCSetOperator() not the Pmat. 294 295 Logically Collective on PC 296 297 Input Parameters: 298 + pc - the preconditioner context 299 - flg - PETSC_TRUE to use the Amat, PETSC_FALSE to use the Pmat (default is false) 300 301 Options Database Key: 302 . -pc_use_amat <true,false> 303 304 Notes: 305 For the common case in which the linear system matrix and the matrix used to construct the 306 preconditioner are identical, this routine is does nothing. 307 308 Level: intermediate 309 310 .seealso: PCGetUseAmat(), PCBJACOBI, PGMG, PCFIELDSPLIT, PCCOMPOSITE 311 @*/ 312 PetscErrorCode PCSetUseAmat(PC pc,PetscBool flg) 313 { 314 PetscFunctionBegin; 315 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 316 pc->useAmat = flg; 317 PetscFunctionReturn(0); 318 } 319 320 #undef __FUNCT__ 321 #define __FUNCT__ "PCGetUseAmat" 322 /*@ 323 PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the 324 operator during the preconditioning process it applies the Amat provided to TSSetRHSJacobian(), 325 TSSetIJacobian(), SNESSetJacobian(), KSPSetOperator() or PCSetOperator() not the Pmat. 326 327 Logically Collective on PC 328 329 Input Parameter: 330 . pc - the preconditioner context 331 332 Output Parameter: 333 . flg - PETSC_TRUE to use the Amat, PETSC_FALSE to use the Pmat (default is false) 334 335 Notes: 336 For the common case in which the linear system matrix and the matrix used to construct the 337 preconditioner are identical, this routine is does nothing. 338 339 Level: intermediate 340 341 .seealso: PCSetUseAmat(), PCBJACOBI, PGMG, PCFIELDSPLIT, PCCOMPOSITE 342 @*/ 343 PetscErrorCode PCGetUseAmat(PC pc,PetscBool *flg) 344 { 345 PetscFunctionBegin; 346 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 347 *flg = pc->useAmat; 348 PetscFunctionReturn(0); 349 } 350 351 #undef __FUNCT__ 352 #define __FUNCT__ "PCCreate" 353 /*@ 354 PCCreate - Creates a preconditioner context. 355 356 Collective on MPI_Comm 357 358 Input Parameter: 359 . comm - MPI communicator 360 361 Output Parameter: 362 . pc - location to put the preconditioner context 363 364 Notes: 365 The default preconditioner for sparse matrices is PCILU or PCICC with 0 fill on one process and block Jacobi with PCILU or ICC 366 in parallel. For dense matrices it is always PCNONE. 367 368 Level: developer 369 370 .keywords: PC, create, context 371 372 .seealso: PCSetUp(), PCApply(), PCDestroy() 373 @*/ 374 PetscErrorCode PCCreate(MPI_Comm comm,PC *newpc) 375 { 376 PC pc; 377 PetscErrorCode ierr; 378 379 PetscFunctionBegin; 380 PetscValidPointer(newpc,1); 381 *newpc = 0; 382 ierr = PCInitializePackage();CHKERRQ(ierr); 383 384 ierr = PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_CLASSID,"PC","Preconditioner","PC",comm,PCDestroy,PCView);CHKERRQ(ierr); 385 386 pc->mat = 0; 387 pc->pmat = 0; 388 pc->setupcalled = 0; 389 pc->setfromoptionscalled = 0; 390 pc->data = 0; 391 pc->diagonalscale = PETSC_FALSE; 392 pc->diagonalscaleleft = 0; 393 pc->diagonalscaleright = 0; 394 395 pc->modifysubmatrices = 0; 396 pc->modifysubmatricesP = 0; 397 398 *newpc = pc; 399 PetscFunctionReturn(0); 400 401 } 402 403 /* -------------------------------------------------------------------------------*/ 404 405 #undef __FUNCT__ 406 #define __FUNCT__ "PCApply" 407 /*@ 408 PCApply - Applies the preconditioner to a vector. 409 410 Collective on PC and Vec 411 412 Input Parameters: 413 + pc - the preconditioner context 414 - x - input vector 415 416 Output Parameter: 417 . y - output vector 418 419 Level: developer 420 421 .keywords: PC, apply 422 423 .seealso: PCApplyTranspose(), PCApplyBAorAB() 424 @*/ 425 PetscErrorCode PCApply(PC pc,Vec x,Vec y) 426 { 427 PetscErrorCode ierr; 428 429 PetscFunctionBegin; 430 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 431 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 432 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 433 if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors"); 434 ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr); 435 if (pc->setupcalled < 2) { 436 ierr = PCSetUp(pc);CHKERRQ(ierr); 437 } 438 if (!pc->ops->apply) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply"); 439 ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr); 440 ierr = (*pc->ops->apply)(pc,x,y);CHKERRQ(ierr); 441 ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr); 442 ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr); 443 PetscFunctionReturn(0); 444 } 445 446 #undef __FUNCT__ 447 #define __FUNCT__ "PCApplySymmetricLeft" 448 /*@ 449 PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector. 450 451 Collective on PC and Vec 452 453 Input Parameters: 454 + pc - the preconditioner context 455 - x - input vector 456 457 Output Parameter: 458 . y - output vector 459 460 Notes: 461 Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners. 462 463 Level: developer 464 465 .keywords: PC, apply, symmetric, left 466 467 .seealso: PCApply(), PCApplySymmetricRight() 468 @*/ 469 PetscErrorCode PCApplySymmetricLeft(PC pc,Vec x,Vec y) 470 { 471 PetscErrorCode ierr; 472 473 PetscFunctionBegin; 474 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 475 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 476 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 477 if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors"); 478 ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr); 479 if (pc->setupcalled < 2) { 480 ierr = PCSetUp(pc);CHKERRQ(ierr); 481 } 482 if (!pc->ops->applysymmetricleft) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have left symmetric apply"); 483 ierr = PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr); 484 ierr = (*pc->ops->applysymmetricleft)(pc,x,y);CHKERRQ(ierr); 485 ierr = PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr); 486 ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr); 487 PetscFunctionReturn(0); 488 } 489 490 #undef __FUNCT__ 491 #define __FUNCT__ "PCApplySymmetricRight" 492 /*@ 493 PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector. 494 495 Collective on PC and Vec 496 497 Input Parameters: 498 + pc - the preconditioner context 499 - x - input vector 500 501 Output Parameter: 502 . y - output vector 503 504 Level: developer 505 506 Notes: 507 Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners. 508 509 .keywords: PC, apply, symmetric, right 510 511 .seealso: PCApply(), PCApplySymmetricLeft() 512 @*/ 513 PetscErrorCode PCApplySymmetricRight(PC pc,Vec x,Vec y) 514 { 515 PetscErrorCode ierr; 516 517 PetscFunctionBegin; 518 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 519 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 520 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 521 if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors"); 522 ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr); 523 if (pc->setupcalled < 2) { 524 ierr = PCSetUp(pc);CHKERRQ(ierr); 525 } 526 if (!pc->ops->applysymmetricright) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have left symmetric apply"); 527 ierr = PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr); 528 ierr = (*pc->ops->applysymmetricright)(pc,x,y);CHKERRQ(ierr); 529 ierr = PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr); 530 ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr); 531 PetscFunctionReturn(0); 532 } 533 534 #undef __FUNCT__ 535 #define __FUNCT__ "PCApplyTranspose" 536 /*@ 537 PCApplyTranspose - Applies the transpose of preconditioner to a vector. 538 539 Collective on PC and Vec 540 541 Input Parameters: 542 + pc - the preconditioner context 543 - x - input vector 544 545 Output Parameter: 546 . y - output vector 547 548 Notes: For complex numbers this applies the non-Hermitian transpose. 549 550 Developer Notes: We need to implement a PCApplyHermitianTranspose() 551 552 Level: developer 553 554 .keywords: PC, apply, transpose 555 556 .seealso: PCApply(), PCApplyBAorAB(), PCApplyBAorABTranspose(), PCApplyTransposeExists() 557 @*/ 558 PetscErrorCode PCApplyTranspose(PC pc,Vec x,Vec y) 559 { 560 PetscErrorCode ierr; 561 562 PetscFunctionBegin; 563 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 564 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 565 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 566 if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors"); 567 ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr); 568 if (pc->setupcalled < 2) { 569 ierr = PCSetUp(pc);CHKERRQ(ierr); 570 } 571 if (!pc->ops->applytranspose) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply transpose"); 572 ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr); 573 ierr = (*pc->ops->applytranspose)(pc,x,y);CHKERRQ(ierr); 574 ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr); 575 ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr); 576 PetscFunctionReturn(0); 577 } 578 579 #undef __FUNCT__ 580 #define __FUNCT__ "PCApplyTransposeExists" 581 /*@ 582 PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation 583 584 Collective on PC and Vec 585 586 Input Parameters: 587 . pc - the preconditioner context 588 589 Output Parameter: 590 . flg - PETSC_TRUE if a transpose operation is defined 591 592 Level: developer 593 594 .keywords: PC, apply, transpose 595 596 .seealso: PCApplyTranspose() 597 @*/ 598 PetscErrorCode PCApplyTransposeExists(PC pc,PetscBool *flg) 599 { 600 PetscFunctionBegin; 601 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 602 PetscValidPointer(flg,2); 603 if (pc->ops->applytranspose) *flg = PETSC_TRUE; 604 else *flg = PETSC_FALSE; 605 PetscFunctionReturn(0); 606 } 607 608 #undef __FUNCT__ 609 #define __FUNCT__ "PCApplyBAorAB" 610 /*@ 611 PCApplyBAorAB - Applies the preconditioner and operator to a vector. y = B*A*x or y = A*B*x. 612 613 Collective on PC and Vec 614 615 Input Parameters: 616 + pc - the preconditioner context 617 . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC 618 . x - input vector 619 - work - work vector 620 621 Output Parameter: 622 . y - output vector 623 624 Level: developer 625 626 Notes: If the PC has had PCSetDiagonalScale() set then D M A D^{-1} for left preconditioning or D A M D^{-1} is actually applied. Note that the 627 specific KSPSolve() method must also be written to handle the post-solve "correction" for the diagonal scaling. 628 629 .keywords: PC, apply, operator 630 631 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose() 632 @*/ 633 PetscErrorCode PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work) 634 { 635 PetscErrorCode ierr; 636 637 PetscFunctionBegin; 638 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 639 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 640 PetscValidHeaderSpecific(y,VEC_CLASSID,4); 641 PetscValidHeaderSpecific(work,VEC_CLASSID,5); 642 if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors"); 643 ierr = VecValidValues(x,3,PETSC_TRUE);CHKERRQ(ierr); 644 if (side != PC_LEFT && side != PC_SYMMETRIC && side != PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Side must be right, left, or symmetric"); 645 if (pc->diagonalscale && side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot include diagonal scaling with symmetric preconditioner application"); 646 647 if (pc->setupcalled < 2) { 648 ierr = PCSetUp(pc);CHKERRQ(ierr); 649 } 650 651 if (pc->diagonalscale) { 652 if (pc->ops->applyBA) { 653 Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */ 654 ierr = VecDuplicate(x,&work2);CHKERRQ(ierr); 655 ierr = PCDiagonalScaleRight(pc,x,work2);CHKERRQ(ierr); 656 ierr = (*pc->ops->applyBA)(pc,side,work2,y,work);CHKERRQ(ierr); 657 ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr); 658 ierr = VecDestroy(&work2);CHKERRQ(ierr); 659 } else if (side == PC_RIGHT) { 660 ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr); 661 ierr = PCApply(pc,y,work);CHKERRQ(ierr); 662 ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr); 663 ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr); 664 } else if (side == PC_LEFT) { 665 ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr); 666 ierr = MatMult(pc->mat,y,work);CHKERRQ(ierr); 667 ierr = PCApply(pc,work,y);CHKERRQ(ierr); 668 ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr); 669 } else if (side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot provide diagonal scaling with symmetric application of preconditioner"); 670 } else { 671 if (pc->ops->applyBA) { 672 ierr = (*pc->ops->applyBA)(pc,side,x,y,work);CHKERRQ(ierr); 673 } else if (side == PC_RIGHT) { 674 ierr = PCApply(pc,x,work);CHKERRQ(ierr); 675 ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr); 676 } else if (side == PC_LEFT) { 677 ierr = MatMult(pc->mat,x,work);CHKERRQ(ierr); 678 ierr = PCApply(pc,work,y);CHKERRQ(ierr); 679 } else if (side == PC_SYMMETRIC) { 680 /* There's an extra copy here; maybe should provide 2 work vectors instead? */ 681 ierr = PCApplySymmetricRight(pc,x,work);CHKERRQ(ierr); 682 ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr); 683 ierr = VecCopy(y,work);CHKERRQ(ierr); 684 ierr = PCApplySymmetricLeft(pc,work,y);CHKERRQ(ierr); 685 } 686 } 687 ierr = VecValidValues(y,4,PETSC_FALSE);CHKERRQ(ierr); 688 PetscFunctionReturn(0); 689 } 690 691 #undef __FUNCT__ 692 #define __FUNCT__ "PCApplyBAorABTranspose" 693 /*@ 694 PCApplyBAorABTranspose - Applies the transpose of the preconditioner 695 and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning, 696 NOT tr(B*A) = tr(A)*tr(B). 697 698 Collective on PC and Vec 699 700 Input Parameters: 701 + pc - the preconditioner context 702 . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC 703 . x - input vector 704 - work - work vector 705 706 Output Parameter: 707 . y - output vector 708 709 710 Notes: this routine is used internally so that the same Krylov code can be used to solve A x = b and A' x = b, with a preconditioner 711 defined by B'. This is why this has the funny form that it computes tr(B) * tr(A) 712 713 Level: developer 714 715 .keywords: PC, apply, operator, transpose 716 717 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB() 718 @*/ 719 PetscErrorCode PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work) 720 { 721 PetscErrorCode ierr; 722 723 PetscFunctionBegin; 724 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 725 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 726 PetscValidHeaderSpecific(y,VEC_CLASSID,4); 727 PetscValidHeaderSpecific(work,VEC_CLASSID,5); 728 if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors"); 729 ierr = VecValidValues(x,3,PETSC_TRUE);CHKERRQ(ierr); 730 if (pc->ops->applyBAtranspose) { 731 ierr = (*pc->ops->applyBAtranspose)(pc,side,x,y,work);CHKERRQ(ierr); 732 ierr = VecValidValues(y,4,PETSC_FALSE);CHKERRQ(ierr); 733 PetscFunctionReturn(0); 734 } 735 if (side != PC_LEFT && side != PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left"); 736 737 if (pc->setupcalled < 2) { 738 ierr = PCSetUp(pc);CHKERRQ(ierr); 739 } 740 741 if (side == PC_RIGHT) { 742 ierr = PCApplyTranspose(pc,x,work);CHKERRQ(ierr); 743 ierr = MatMultTranspose(pc->mat,work,y);CHKERRQ(ierr); 744 } else if (side == PC_LEFT) { 745 ierr = MatMultTranspose(pc->mat,x,work);CHKERRQ(ierr); 746 ierr = PCApplyTranspose(pc,work,y);CHKERRQ(ierr); 747 } 748 /* add support for PC_SYMMETRIC */ 749 ierr = VecValidValues(y,4,PETSC_FALSE);CHKERRQ(ierr); 750 PetscFunctionReturn(0); 751 } 752 753 /* -------------------------------------------------------------------------------*/ 754 755 #undef __FUNCT__ 756 #define __FUNCT__ "PCApplyRichardsonExists" 757 /*@ 758 PCApplyRichardsonExists - Determines whether a particular preconditioner has a 759 built-in fast application of Richardson's method. 760 761 Not Collective 762 763 Input Parameter: 764 . pc - the preconditioner 765 766 Output Parameter: 767 . exists - PETSC_TRUE or PETSC_FALSE 768 769 Level: developer 770 771 .keywords: PC, apply, Richardson, exists 772 773 .seealso: PCApplyRichardson() 774 @*/ 775 PetscErrorCode PCApplyRichardsonExists(PC pc,PetscBool *exists) 776 { 777 PetscFunctionBegin; 778 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 779 PetscValidIntPointer(exists,2); 780 if (pc->ops->applyrichardson) *exists = PETSC_TRUE; 781 else *exists = PETSC_FALSE; 782 PetscFunctionReturn(0); 783 } 784 785 #undef __FUNCT__ 786 #define __FUNCT__ "PCApplyRichardson" 787 /*@ 788 PCApplyRichardson - Applies several steps of Richardson iteration with 789 the particular preconditioner. This routine is usually used by the 790 Krylov solvers and not the application code directly. 791 792 Collective on PC 793 794 Input Parameters: 795 + pc - the preconditioner context 796 . b - the right hand side 797 . w - one work vector 798 . rtol - relative decrease in residual norm convergence criteria 799 . abstol - absolute residual norm convergence criteria 800 . dtol - divergence residual norm increase criteria 801 . its - the number of iterations to apply. 802 - guesszero - if the input x contains nonzero initial guess 803 804 Output Parameter: 805 + outits - number of iterations actually used (for SOR this always equals its) 806 . reason - the reason the apply terminated 807 - y - the solution (also contains initial guess if guesszero is PETSC_FALSE 808 809 Notes: 810 Most preconditioners do not support this function. Use the command 811 PCApplyRichardsonExists() to determine if one does. 812 813 Except for the multigrid PC this routine ignores the convergence tolerances 814 and always runs for the number of iterations 815 816 Level: developer 817 818 .keywords: PC, apply, Richardson 819 820 .seealso: PCApplyRichardsonExists() 821 @*/ 822 PetscErrorCode PCApplyRichardson(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason) 823 { 824 PetscErrorCode ierr; 825 826 PetscFunctionBegin; 827 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 828 PetscValidHeaderSpecific(b,VEC_CLASSID,2); 829 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 830 PetscValidHeaderSpecific(w,VEC_CLASSID,4); 831 if (b == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"b and y must be different vectors"); 832 if (pc->setupcalled < 2) { 833 ierr = PCSetUp(pc);CHKERRQ(ierr); 834 } 835 if (!pc->ops->applyrichardson) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply richardson"); 836 ierr = (*pc->ops->applyrichardson)(pc,b,y,w,rtol,abstol,dtol,its,guesszero,outits,reason);CHKERRQ(ierr); 837 PetscFunctionReturn(0); 838 } 839 840 /* 841 a setupcall of 0 indicates never setup, 842 1 needs to be resetup, 843 2 does not need any changes. 844 */ 845 #undef __FUNCT__ 846 #define __FUNCT__ "PCSetUp" 847 /*@ 848 PCSetUp - Prepares for the use of a preconditioner. 849 850 Collective on PC 851 852 Input Parameter: 853 . pc - the preconditioner context 854 855 Level: developer 856 857 .keywords: PC, setup 858 859 .seealso: PCCreate(), PCApply(), PCDestroy() 860 @*/ 861 PetscErrorCode PCSetUp(PC pc) 862 { 863 PetscErrorCode ierr; 864 const char *def; 865 866 PetscFunctionBegin; 867 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 868 if (!pc->mat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first"); 869 870 if (pc->setupcalled > 1) { 871 ierr = PetscInfo(pc,"Setting PC with identical preconditioner\n");CHKERRQ(ierr); 872 PetscFunctionReturn(0); 873 } else if (!pc->setupcalled) { 874 ierr = PetscInfo(pc,"Setting up new PC\n");CHKERRQ(ierr); 875 } else if (pc->flag == SAME_NONZERO_PATTERN) { 876 ierr = PetscInfo(pc,"Setting up PC with same nonzero pattern\n");CHKERRQ(ierr); 877 } else { 878 ierr = PetscInfo(pc,"Setting up PC with different nonzero pattern\n");CHKERRQ(ierr); 879 } 880 881 if (!((PetscObject)pc)->type_name) { 882 ierr = PCGetDefaultType_Private(pc,&def);CHKERRQ(ierr); 883 ierr = PCSetType(pc,def);CHKERRQ(ierr); 884 } 885 886 ierr = PetscLogEventBegin(PC_SetUp,pc,0,0,0);CHKERRQ(ierr); 887 if (pc->ops->setup) { 888 ierr = (*pc->ops->setup)(pc);CHKERRQ(ierr); 889 } 890 pc->setupcalled = 2; 891 892 ierr = PetscLogEventEnd(PC_SetUp,pc,0,0,0);CHKERRQ(ierr); 893 PetscFunctionReturn(0); 894 } 895 896 #undef __FUNCT__ 897 #define __FUNCT__ "PCSetUpOnBlocks" 898 /*@ 899 PCSetUpOnBlocks - Sets up the preconditioner for each block in 900 the block Jacobi, block Gauss-Seidel, and overlapping Schwarz 901 methods. 902 903 Collective on PC 904 905 Input Parameters: 906 . pc - the preconditioner context 907 908 Level: developer 909 910 .keywords: PC, setup, blocks 911 912 .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp() 913 @*/ 914 PetscErrorCode PCSetUpOnBlocks(PC pc) 915 { 916 PetscErrorCode ierr; 917 918 PetscFunctionBegin; 919 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 920 if (!pc->ops->setuponblocks) PetscFunctionReturn(0); 921 ierr = PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr); 922 ierr = (*pc->ops->setuponblocks)(pc);CHKERRQ(ierr); 923 ierr = PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr); 924 PetscFunctionReturn(0); 925 } 926 927 #undef __FUNCT__ 928 #define __FUNCT__ "PCSetModifySubMatrices" 929 /*@C 930 PCSetModifySubMatrices - Sets a user-defined routine for modifying the 931 submatrices that arise within certain subdomain-based preconditioners. 932 The basic submatrices are extracted from the preconditioner matrix as 933 usual; the user can then alter these (for example, to set different boundary 934 conditions for each submatrix) before they are used for the local solves. 935 936 Logically Collective on PC 937 938 Input Parameters: 939 + pc - the preconditioner context 940 . func - routine for modifying the submatrices 941 - ctx - optional user-defined context (may be null) 942 943 Calling sequence of func: 944 $ func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx); 945 946 . row - an array of index sets that contain the global row numbers 947 that comprise each local submatrix 948 . col - an array of index sets that contain the global column numbers 949 that comprise each local submatrix 950 . submat - array of local submatrices 951 - ctx - optional user-defined context for private data for the 952 user-defined func routine (may be null) 953 954 Notes: 955 PCSetModifySubMatrices() MUST be called before KSPSetUp() and 956 KSPSolve(). 957 958 A routine set by PCSetModifySubMatrices() is currently called within 959 the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM) 960 preconditioners. All other preconditioners ignore this routine. 961 962 Level: advanced 963 964 .keywords: PC, set, modify, submatrices 965 966 .seealso: PCModifySubMatrices(), PCASMGetSubMatrices() 967 @*/ 968 PetscErrorCode PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx) 969 { 970 PetscFunctionBegin; 971 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 972 pc->modifysubmatrices = func; 973 pc->modifysubmatricesP = ctx; 974 PetscFunctionReturn(0); 975 } 976 977 #undef __FUNCT__ 978 #define __FUNCT__ "PCModifySubMatrices" 979 /*@C 980 PCModifySubMatrices - Calls an optional user-defined routine within 981 certain preconditioners if one has been set with PCSetModifySubMarices(). 982 983 Collective on PC 984 985 Input Parameters: 986 + pc - the preconditioner context 987 . nsub - the number of local submatrices 988 . row - an array of index sets that contain the global row numbers 989 that comprise each local submatrix 990 . col - an array of index sets that contain the global column numbers 991 that comprise each local submatrix 992 . submat - array of local submatrices 993 - ctx - optional user-defined context for private data for the 994 user-defined routine (may be null) 995 996 Output Parameter: 997 . submat - array of local submatrices (the entries of which may 998 have been modified) 999 1000 Notes: 1001 The user should NOT generally call this routine, as it will 1002 automatically be called within certain preconditioners (currently 1003 block Jacobi, additive Schwarz) if set. 1004 1005 The basic submatrices are extracted from the preconditioner matrix 1006 as usual; the user can then alter these (for example, to set different 1007 boundary conditions for each submatrix) before they are used for the 1008 local solves. 1009 1010 Level: developer 1011 1012 .keywords: PC, modify, submatrices 1013 1014 .seealso: PCSetModifySubMatrices() 1015 @*/ 1016 PetscErrorCode PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx) 1017 { 1018 PetscErrorCode ierr; 1019 1020 PetscFunctionBegin; 1021 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1022 if (!pc->modifysubmatrices) PetscFunctionReturn(0); 1023 ierr = PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr); 1024 ierr = (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);CHKERRQ(ierr); 1025 ierr = PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr); 1026 PetscFunctionReturn(0); 1027 } 1028 1029 #undef __FUNCT__ 1030 #define __FUNCT__ "PCSetOperators" 1031 /*@ 1032 PCSetOperators - Sets the matrix associated with the linear system and 1033 a (possibly) different one associated with the preconditioner. 1034 1035 Logically Collective on PC and Mat 1036 1037 Input Parameters: 1038 + pc - the preconditioner context 1039 . Amat - the matrix that defines the linear system 1040 . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat. 1041 - flag - flag indicating information about the preconditioner matrix structure 1042 during successive linear solves. This flag is ignored the first time a 1043 linear system is solved, and thus is irrelevant when solving just one linear 1044 system. 1045 1046 Notes: 1047 The flag can be used to eliminate unnecessary work in the preconditioner 1048 during the repeated solution of linear systems of the same size. The 1049 available options are 1050 + SAME_PRECONDITIONER - 1051 Pmat is identical during successive linear solves. 1052 This option is intended for folks who are using 1053 different Amat and Pmat matrices and wish to reuse the 1054 same preconditioner matrix. For example, this option 1055 saves work by not recomputing incomplete factorization 1056 for ILU/ICC preconditioners. 1057 . SAME_NONZERO_PATTERN - 1058 Pmat has the same nonzero structure during 1059 successive linear solves. 1060 - DIFFERENT_NONZERO_PATTERN - 1061 Pmat does not have the same nonzero structure. 1062 1063 Passing a NULL for Amat or Pmat removes the matrix that is currently used. 1064 1065 If you wish to replace either Amat or Pmat but leave the other one untouched then 1066 first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference() 1067 on it and then pass it back in in your call to KSPSetOperators(). 1068 1069 Caution: 1070 If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion 1071 and does not check the structure of the matrix. If you erroneously 1072 claim that the structure is the same when it actually is not, the new 1073 preconditioner will not function correctly. Thus, use this optimization 1074 feature carefully! 1075 1076 If in doubt about whether your preconditioner matrix has changed 1077 structure or not, use the flag DIFFERENT_NONZERO_PATTERN. 1078 1079 More Notes about Repeated Solution of Linear Systems: 1080 PETSc does NOT reset the matrix entries of either Amat or Pmat 1081 to zero after a linear solve; the user is completely responsible for 1082 matrix assembly. See the routine MatZeroEntries() if desiring to 1083 zero all elements of a matrix. 1084 1085 Level: intermediate 1086 1087 .keywords: PC, set, operators, matrix, linear system 1088 1089 .seealso: PCGetOperators(), MatZeroEntries() 1090 @*/ 1091 PetscErrorCode PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag) 1092 { 1093 PetscErrorCode ierr; 1094 PetscInt m1,n1,m2,n2; 1095 1096 PetscFunctionBegin; 1097 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1098 if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 1099 if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_CLASSID,3); 1100 if (Amat) PetscCheckSameComm(pc,1,Amat,2); 1101 if (Pmat) PetscCheckSameComm(pc,1,Pmat,3); 1102 if (pc->setupcalled && Amat && Pmat) { 1103 ierr = MatGetLocalSize(Amat,&m1,&n1);CHKERRQ(ierr); 1104 ierr = MatGetLocalSize(pc->mat,&m2,&n2);CHKERRQ(ierr); 1105 if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot change local size of Amat after use old sizes %D %D new sizes %D %D",m2,n2,m1,n1); 1106 ierr = MatGetLocalSize(Pmat,&m1,&n1);CHKERRQ(ierr); 1107 ierr = MatGetLocalSize(pc->pmat,&m2,&n2);CHKERRQ(ierr); 1108 if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot change local size of Pmat after use old sizes %D %D new sizes %D %D",m2,n2,m1,n1); 1109 } 1110 1111 /* reference first in case the matrices are the same */ 1112 if (Amat) {ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);} 1113 ierr = MatDestroy(&pc->mat);CHKERRQ(ierr); 1114 if (Pmat) {ierr = PetscObjectReference((PetscObject)Pmat);CHKERRQ(ierr);} 1115 ierr = MatDestroy(&pc->pmat);CHKERRQ(ierr); 1116 pc->mat = Amat; 1117 pc->pmat = Pmat; 1118 1119 if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) { 1120 pc->setupcalled = 1; 1121 } 1122 pc->flag = flag; 1123 PetscFunctionReturn(0); 1124 } 1125 1126 #undef __FUNCT__ 1127 #define __FUNCT__ "PCGetOperators" 1128 /*@C 1129 PCGetOperators - Gets the matrix associated with the linear system and 1130 possibly a different one associated with the preconditioner. 1131 1132 Not collective, though parallel Mats are returned if the PC is parallel 1133 1134 Input Parameter: 1135 . pc - the preconditioner context 1136 1137 Output Parameters: 1138 + Amat - the matrix defining the linear system 1139 . Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat. 1140 - flag - flag indicating information about the preconditioner 1141 matrix structure. See PCSetOperators() for details. 1142 1143 Level: intermediate 1144 1145 Notes: Does not increase the reference count of the matrices, so you should not destroy them 1146 1147 Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators 1148 are created in PC and returned to the user. In this case, if both operators 1149 mat and pmat are requested, two DIFFERENT operators will be returned. If 1150 only one is requested both operators in the PC will be the same (i.e. as 1151 if one had called KSP/PCSetOperators() with the same argument for both Mats). 1152 The user must set the sizes of the returned matrices and their type etc just 1153 as if the user created them with MatCreate(). For example, 1154 1155 $ KSP/PCGetOperators(ksp/pc,&Amat,NULL,NULL); is equivalent to 1156 $ set size, type, etc of Amat 1157 1158 $ MatCreate(comm,&mat); 1159 $ KSP/PCSetOperators(ksp/pc,Amat,Amat,SAME_NONZERO_PATTERN); 1160 $ PetscObjectDereference((PetscObject)mat); 1161 $ set size, type, etc of Amat 1162 1163 and 1164 1165 $ KSP/PCGetOperators(ksp/pc,&Amat,&Pmat,NULL); is equivalent to 1166 $ set size, type, etc of Amat and Pmat 1167 1168 $ MatCreate(comm,&Amat); 1169 $ MatCreate(comm,&Pmat); 1170 $ KSP/PCSetOperators(ksp/pc,Amat,Pmat,SAME_NONZERO_PATTERN); 1171 $ PetscObjectDereference((PetscObject)Amat); 1172 $ PetscObjectDereference((PetscObject)Pmat); 1173 $ set size, type, etc of Amat and Pmat 1174 1175 The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy 1176 of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely 1177 managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look 1178 at this is when you create a SNES you do not NEED to create a KSP and attach it to 1179 the SNES object (the SNES object manages it for you). Similarly when you create a KSP 1180 you do not need to attach a PC to it (the KSP object manages the PC object for you). 1181 Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when 1182 it can be created for you? 1183 1184 1185 .keywords: PC, get, operators, matrix, linear system 1186 1187 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet() 1188 @*/ 1189 PetscErrorCode PCGetOperators(PC pc,Mat *Amat,Mat *Pmat,MatStructure *flag) 1190 { 1191 PetscErrorCode ierr; 1192 1193 PetscFunctionBegin; 1194 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1195 if (Amat) { 1196 if (!pc->mat) { 1197 if (pc->pmat && !Pmat) { /* Apmat has been set, but user did not request it, so use for Amat */ 1198 pc->mat = pc->pmat; 1199 ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr); 1200 } else { /* both Amat and Pmat are empty */ 1201 ierr = MatCreate(PetscObjectComm((PetscObject)pc),&pc->mat);CHKERRQ(ierr); 1202 if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */ 1203 pc->pmat = pc->mat; 1204 ierr = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr); 1205 } 1206 } 1207 } 1208 *Amat = pc->mat; 1209 } 1210 if (Pmat) { 1211 if (!pc->pmat) { 1212 if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */ 1213 pc->pmat = pc->mat; 1214 ierr = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr); 1215 } else { 1216 ierr = MatCreate(PetscObjectComm((PetscObject)pc),&pc->pmat);CHKERRQ(ierr); 1217 if (!Amat) { /* user did NOT request Amat, so make same as Pmat */ 1218 pc->mat = pc->pmat; 1219 ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr); 1220 } 1221 } 1222 } 1223 *Pmat = pc->pmat; 1224 } 1225 if (flag) *flag = pc->flag; 1226 PetscFunctionReturn(0); 1227 } 1228 1229 #undef __FUNCT__ 1230 #define __FUNCT__ "PCGetOperatorsSet" 1231 /*@C 1232 PCGetOperatorsSet - Determines if the matrix associated with the linear system and 1233 possibly a different one associated with the preconditioner have been set in the PC. 1234 1235 Not collective, though the results on all processes should be the same 1236 1237 Input Parameter: 1238 . pc - the preconditioner context 1239 1240 Output Parameters: 1241 + mat - the matrix associated with the linear system was set 1242 - pmat - matrix associated with the preconditioner was set, usually the same 1243 1244 Level: intermediate 1245 1246 .keywords: PC, get, operators, matrix, linear system 1247 1248 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators() 1249 @*/ 1250 PetscErrorCode PCGetOperatorsSet(PC pc,PetscBool *mat,PetscBool *pmat) 1251 { 1252 PetscFunctionBegin; 1253 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1254 if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE; 1255 if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE; 1256 PetscFunctionReturn(0); 1257 } 1258 1259 #undef __FUNCT__ 1260 #define __FUNCT__ "PCFactorGetMatrix" 1261 /*@ 1262 PCFactorGetMatrix - Gets the factored matrix from the 1263 preconditioner context. This routine is valid only for the LU, 1264 incomplete LU, Cholesky, and incomplete Cholesky methods. 1265 1266 Not Collective on PC though Mat is parallel if PC is parallel 1267 1268 Input Parameters: 1269 . pc - the preconditioner context 1270 1271 Output parameters: 1272 . mat - the factored matrix 1273 1274 Level: advanced 1275 1276 Notes: Does not increase the reference count for the matrix so DO NOT destroy it 1277 1278 .keywords: PC, get, factored, matrix 1279 @*/ 1280 PetscErrorCode PCFactorGetMatrix(PC pc,Mat *mat) 1281 { 1282 PetscErrorCode ierr; 1283 1284 PetscFunctionBegin; 1285 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1286 PetscValidPointer(mat,2); 1287 if (pc->ops->getfactoredmatrix) { 1288 ierr = (*pc->ops->getfactoredmatrix)(pc,mat);CHKERRQ(ierr); 1289 } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC type does not support getting factor matrix"); 1290 PetscFunctionReturn(0); 1291 } 1292 1293 #undef __FUNCT__ 1294 #define __FUNCT__ "PCSetOptionsPrefix" 1295 /*@C 1296 PCSetOptionsPrefix - Sets the prefix used for searching for all 1297 PC options in the database. 1298 1299 Logically Collective on PC 1300 1301 Input Parameters: 1302 + pc - the preconditioner context 1303 - prefix - the prefix string to prepend to all PC option requests 1304 1305 Notes: 1306 A hyphen (-) must NOT be given at the beginning of the prefix name. 1307 The first character of all runtime options is AUTOMATICALLY the 1308 hyphen. 1309 1310 Level: advanced 1311 1312 .keywords: PC, set, options, prefix, database 1313 1314 .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix() 1315 @*/ 1316 PetscErrorCode PCSetOptionsPrefix(PC pc,const char prefix[]) 1317 { 1318 PetscErrorCode ierr; 1319 1320 PetscFunctionBegin; 1321 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1322 ierr = PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr); 1323 PetscFunctionReturn(0); 1324 } 1325 1326 #undef __FUNCT__ 1327 #define __FUNCT__ "PCAppendOptionsPrefix" 1328 /*@C 1329 PCAppendOptionsPrefix - Appends to the prefix used for searching for all 1330 PC options in the database. 1331 1332 Logically Collective on PC 1333 1334 Input Parameters: 1335 + pc - the preconditioner context 1336 - prefix - the prefix string to prepend to all PC option requests 1337 1338 Notes: 1339 A hyphen (-) must NOT be given at the beginning of the prefix name. 1340 The first character of all runtime options is AUTOMATICALLY the 1341 hyphen. 1342 1343 Level: advanced 1344 1345 .keywords: PC, append, options, prefix, database 1346 1347 .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix() 1348 @*/ 1349 PetscErrorCode PCAppendOptionsPrefix(PC pc,const char prefix[]) 1350 { 1351 PetscErrorCode ierr; 1352 1353 PetscFunctionBegin; 1354 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1355 ierr = PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr); 1356 PetscFunctionReturn(0); 1357 } 1358 1359 #undef __FUNCT__ 1360 #define __FUNCT__ "PCGetOptionsPrefix" 1361 /*@C 1362 PCGetOptionsPrefix - Gets the prefix used for searching for all 1363 PC options in the database. 1364 1365 Not Collective 1366 1367 Input Parameters: 1368 . pc - the preconditioner context 1369 1370 Output Parameters: 1371 . prefix - pointer to the prefix string used, is returned 1372 1373 Notes: On the fortran side, the user should pass in a string 'prifix' of 1374 sufficient length to hold the prefix. 1375 1376 Level: advanced 1377 1378 .keywords: PC, get, options, prefix, database 1379 1380 .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix() 1381 @*/ 1382 PetscErrorCode PCGetOptionsPrefix(PC pc,const char *prefix[]) 1383 { 1384 PetscErrorCode ierr; 1385 1386 PetscFunctionBegin; 1387 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1388 PetscValidPointer(prefix,2); 1389 ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr); 1390 PetscFunctionReturn(0); 1391 } 1392 1393 #undef __FUNCT__ 1394 #define __FUNCT__ "PCPreSolve" 1395 /*@ 1396 PCPreSolve - Optional pre-solve phase, intended for any 1397 preconditioner-specific actions that must be performed before 1398 the iterative solve itself. 1399 1400 Collective on PC 1401 1402 Input Parameters: 1403 + pc - the preconditioner context 1404 - ksp - the Krylov subspace context 1405 1406 Level: developer 1407 1408 Sample of Usage: 1409 .vb 1410 PCPreSolve(pc,ksp); 1411 KSPSolve(ksp,b,x); 1412 PCPostSolve(pc,ksp); 1413 .ve 1414 1415 Notes: 1416 The pre-solve phase is distinct from the PCSetUp() phase. 1417 1418 KSPSolve() calls this directly, so is rarely called by the user. 1419 1420 .keywords: PC, pre-solve 1421 1422 .seealso: PCPostSolve() 1423 @*/ 1424 PetscErrorCode PCPreSolve(PC pc,KSP ksp) 1425 { 1426 PetscErrorCode ierr; 1427 Vec x,rhs; 1428 1429 PetscFunctionBegin; 1430 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1431 PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 1432 pc->presolvedone++; 1433 if (pc->presolvedone > 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot embed PCPreSolve() more than twice"); 1434 ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr); 1435 ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr); 1436 1437 if (pc->ops->presolve) { 1438 ierr = (*pc->ops->presolve)(pc,ksp,rhs,x);CHKERRQ(ierr); 1439 } 1440 PetscFunctionReturn(0); 1441 } 1442 1443 #undef __FUNCT__ 1444 #define __FUNCT__ "PCPostSolve" 1445 /*@ 1446 PCPostSolve - Optional post-solve phase, intended for any 1447 preconditioner-specific actions that must be performed after 1448 the iterative solve itself. 1449 1450 Collective on PC 1451 1452 Input Parameters: 1453 + pc - the preconditioner context 1454 - ksp - the Krylov subspace context 1455 1456 Sample of Usage: 1457 .vb 1458 PCPreSolve(pc,ksp); 1459 KSPSolve(ksp,b,x); 1460 PCPostSolve(pc,ksp); 1461 .ve 1462 1463 Note: 1464 KSPSolve() calls this routine directly, so it is rarely called by the user. 1465 1466 Level: developer 1467 1468 .keywords: PC, post-solve 1469 1470 .seealso: PCPreSolve(), KSPSolve() 1471 @*/ 1472 PetscErrorCode PCPostSolve(PC pc,KSP ksp) 1473 { 1474 PetscErrorCode ierr; 1475 Vec x,rhs; 1476 1477 PetscFunctionBegin; 1478 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1479 PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 1480 pc->presolvedone--; 1481 ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr); 1482 ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr); 1483 if (pc->ops->postsolve) { 1484 ierr = (*pc->ops->postsolve)(pc,ksp,rhs,x);CHKERRQ(ierr); 1485 } 1486 PetscFunctionReturn(0); 1487 } 1488 1489 #undef __FUNCT__ 1490 #define __FUNCT__ "PCLoad" 1491 /*@C 1492 PCLoad - Loads a PC that has been stored in binary with PCView(). 1493 1494 Collective on PetscViewer 1495 1496 Input Parameters: 1497 + newdm - the newly loaded PC, this needs to have been created with PCCreate() or 1498 some related function before a call to PCLoad(). 1499 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() 1500 1501 Level: intermediate 1502 1503 Notes: 1504 The type is determined by the data in the file, any type set into the PC before this call is ignored. 1505 1506 Notes for advanced users: 1507 Most users should not need to know the details of the binary storage 1508 format, since PCLoad() and PCView() completely hide these details. 1509 But for anyone who's interested, the standard binary matrix storage 1510 format is 1511 .vb 1512 has not yet been determined 1513 .ve 1514 1515 .seealso: PetscViewerBinaryOpen(), PCView(), MatLoad(), VecLoad() 1516 @*/ 1517 PetscErrorCode PCLoad(PC newdm, PetscViewer viewer) 1518 { 1519 PetscErrorCode ierr; 1520 PetscBool isbinary; 1521 PetscInt classid; 1522 char type[256]; 1523 1524 PetscFunctionBegin; 1525 PetscValidHeaderSpecific(newdm,PC_CLASSID,1); 1526 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 1527 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 1528 if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 1529 1530 ierr = PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);CHKERRQ(ierr); 1531 if (classid != PC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not PC next in file"); 1532 ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr); 1533 ierr = PCSetType(newdm, type);CHKERRQ(ierr); 1534 if (newdm->ops->load) { 1535 ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr); 1536 } 1537 PetscFunctionReturn(0); 1538 } 1539 1540 #include <petscdraw.h> 1541 #if defined(PETSC_HAVE_SAWS) 1542 #include <petscviewersaws.h> 1543 #endif 1544 #undef __FUNCT__ 1545 #define __FUNCT__ "PCView" 1546 /*@C 1547 PCView - Prints the PC data structure. 1548 1549 Collective on PC 1550 1551 Input Parameters: 1552 + PC - the PC context 1553 - viewer - optional visualization context 1554 1555 Note: 1556 The available visualization contexts include 1557 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 1558 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 1559 output where only the first processor opens 1560 the file. All other processors send their 1561 data to the first processor to print. 1562 1563 The user can open an alternative visualization contexts with 1564 PetscViewerASCIIOpen() (output to a specified file). 1565 1566 Level: developer 1567 1568 .keywords: PC, view 1569 1570 .seealso: KSPView(), PetscViewerASCIIOpen() 1571 @*/ 1572 PetscErrorCode PCView(PC pc,PetscViewer viewer) 1573 { 1574 PCType cstr; 1575 PetscErrorCode ierr; 1576 PetscBool iascii,isstring,isbinary,isdraw; 1577 PetscViewerFormat format; 1578 #if defined(PETSC_HAVE_SAWS) 1579 PetscBool isams; 1580 #endif 1581 1582 PetscFunctionBegin; 1583 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1584 if (!viewer) { 1585 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr); 1586 } 1587 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 1588 PetscCheckSameComm(pc,1,viewer,2); 1589 1590 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1591 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 1592 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 1593 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 1594 #if defined(PETSC_HAVE_SAWS) 1595 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&isams);CHKERRQ(ierr); 1596 #endif 1597 1598 if (iascii) { 1599 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1600 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)pc,viewer);CHKERRQ(ierr); 1601 if (!pc->setupcalled) { 1602 ierr = PetscViewerASCIIPrintf(viewer," PC has not been set up so information may be incomplete\n");CHKERRQ(ierr); 1603 } 1604 if (pc->ops->view) { 1605 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1606 ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr); 1607 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1608 } 1609 if (pc->mat) { 1610 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 1611 if (pc->pmat == pc->mat) { 1612 ierr = PetscViewerASCIIPrintf(viewer," linear system matrix = precond matrix:\n");CHKERRQ(ierr); 1613 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1614 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 1615 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1616 } else { 1617 if (pc->pmat) { 1618 ierr = PetscViewerASCIIPrintf(viewer," linear system matrix followed by preconditioner matrix:\n");CHKERRQ(ierr); 1619 } else { 1620 ierr = PetscViewerASCIIPrintf(viewer," linear system matrix:\n");CHKERRQ(ierr); 1621 } 1622 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1623 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 1624 if (pc->pmat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);} 1625 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1626 } 1627 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1628 } 1629 } else if (isstring) { 1630 ierr = PCGetType(pc,&cstr);CHKERRQ(ierr); 1631 ierr = PetscViewerStringSPrintf(viewer," %-7.7s",cstr);CHKERRQ(ierr); 1632 if (pc->ops->view) {ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);} 1633 } else if (isbinary) { 1634 PetscInt classid = PC_FILE_CLASSID; 1635 MPI_Comm comm; 1636 PetscMPIInt rank; 1637 char type[256]; 1638 1639 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1640 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1641 if (!rank) { 1642 ierr = PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 1643 ierr = PetscStrncpy(type,((PetscObject)pc)->type_name,256);CHKERRQ(ierr); 1644 ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 1645 } 1646 if (pc->ops->view) { 1647 ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr); 1648 } 1649 } else if (isdraw) { 1650 PetscDraw draw; 1651 char str[25]; 1652 PetscReal x,y,bottom,h; 1653 PetscInt n; 1654 1655 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1656 ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 1657 if (pc->mat) { 1658 ierr = MatGetSize(pc->mat,&n,NULL);CHKERRQ(ierr); 1659 ierr = PetscSNPrintf(str,25,"PC: %s (%D)",((PetscObject)pc)->type_name,n);CHKERRQ(ierr); 1660 } else { 1661 ierr = PetscSNPrintf(str,25,"PC: %s",((PetscObject)pc)->type_name);CHKERRQ(ierr); 1662 } 1663 ierr = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr); 1664 bottom = y - h; 1665 ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 1666 if (pc->ops->view) { 1667 ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr); 1668 } 1669 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 1670 #if defined(PETSC_HAVE_SAWS) 1671 } else if (isams) { 1672 PetscMPIInt rank; 1673 1674 ierr = PetscObjectName((PetscObject)pc);CHKERRQ(ierr); 1675 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 1676 if (!((PetscObject)pc)->amsmem && !rank) { 1677 ierr = PetscObjectViewSAWs((PetscObject)pc,viewer);CHKERRQ(ierr); 1678 } 1679 if (pc->mat) {ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);} 1680 if (pc->pmat && pc->pmat != pc->mat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);} 1681 #endif 1682 } 1683 PetscFunctionReturn(0); 1684 } 1685 1686 1687 #undef __FUNCT__ 1688 #define __FUNCT__ "PCSetInitialGuessNonzero" 1689 /*@ 1690 PCSetInitialGuessNonzero - Tells the iterative solver that the 1691 initial guess is nonzero; otherwise PC assumes the initial guess 1692 is to be zero (and thus zeros it out before solving). 1693 1694 Logically Collective on PC 1695 1696 Input Parameters: 1697 + pc - iterative context obtained from PCCreate() 1698 - flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero 1699 1700 Level: Developer 1701 1702 Notes: 1703 This is a weird function. Since PC's are linear operators on the right hand side they 1704 CANNOT use an initial guess. This function is for the "pass-through" preconditioners 1705 PCKSP and PCREDUNDANT and causes the inner KSP object to use the nonzero 1706 initial guess. Not currently working for PCREDUNDANT, that has to be rewritten to use KSP. 1707 1708 1709 .keywords: PC, set, initial guess, nonzero 1710 1711 .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll() 1712 @*/ 1713 PetscErrorCode PCSetInitialGuessNonzero(PC pc,PetscBool flg) 1714 { 1715 PetscFunctionBegin; 1716 PetscValidLogicalCollectiveBool(pc,flg,2); 1717 pc->nonzero_guess = flg; 1718 PetscFunctionReturn(0); 1719 } 1720 1721 #undef __FUNCT__ 1722 #define __FUNCT__ "PCRegister" 1723 /*@C 1724 PCRegister - Adds a method to the preconditioner package. 1725 1726 Not collective 1727 1728 Input Parameters: 1729 + name_solver - name of a new user-defined solver 1730 - routine_create - routine to create method context 1731 1732 Notes: 1733 PCRegister() may be called multiple times to add several user-defined preconditioners. 1734 1735 Sample usage: 1736 .vb 1737 PCRegister("my_solver", MySolverCreate); 1738 .ve 1739 1740 Then, your solver can be chosen with the procedural interface via 1741 $ PCSetType(pc,"my_solver") 1742 or at runtime via the option 1743 $ -pc_type my_solver 1744 1745 Level: advanced 1746 1747 .keywords: PC, register 1748 1749 .seealso: PCRegisterAll(), PCRegisterDestroy() 1750 @*/ 1751 PetscErrorCode PCRegister(const char sname[],PetscErrorCode (*function)(PC)) 1752 { 1753 PetscErrorCode ierr; 1754 1755 PetscFunctionBegin; 1756 ierr = PetscFunctionListAdd(&PCList,sname,function);CHKERRQ(ierr); 1757 PetscFunctionReturn(0); 1758 } 1759 1760 #undef __FUNCT__ 1761 #define __FUNCT__ "PCComputeExplicitOperator" 1762 /*@ 1763 PCComputeExplicitOperator - Computes the explicit preconditioned operator. 1764 1765 Collective on PC 1766 1767 Input Parameter: 1768 . pc - the preconditioner object 1769 1770 Output Parameter: 1771 . mat - the explict preconditioned operator 1772 1773 Notes: 1774 This computation is done by applying the operators to columns of the 1775 identity matrix. 1776 1777 Currently, this routine uses a dense matrix format when 1 processor 1778 is used and a sparse format otherwise. This routine is costly in general, 1779 and is recommended for use only with relatively small systems. 1780 1781 Level: advanced 1782 1783 .keywords: PC, compute, explicit, operator 1784 1785 .seealso: KSPComputeExplicitOperator() 1786 1787 @*/ 1788 PetscErrorCode PCComputeExplicitOperator(PC pc,Mat *mat) 1789 { 1790 Vec in,out; 1791 PetscErrorCode ierr; 1792 PetscInt i,M,m,*rows,start,end; 1793 PetscMPIInt size; 1794 MPI_Comm comm; 1795 PetscScalar *array,one = 1.0; 1796 1797 PetscFunctionBegin; 1798 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1799 PetscValidPointer(mat,2); 1800 1801 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1802 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1803 1804 if (!pc->pmat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call"); 1805 ierr = MatGetVecs(pc->pmat,&in,0);CHKERRQ(ierr); 1806 ierr = VecDuplicate(in,&out);CHKERRQ(ierr); 1807 ierr = VecGetOwnershipRange(in,&start,&end);CHKERRQ(ierr); 1808 ierr = VecGetSize(in,&M);CHKERRQ(ierr); 1809 ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr); 1810 ierr = PetscMalloc1((m+1),&rows);CHKERRQ(ierr); 1811 for (i=0; i<m; i++) rows[i] = start + i; 1812 1813 ierr = MatCreate(comm,mat);CHKERRQ(ierr); 1814 ierr = MatSetSizes(*mat,m,m,M,M);CHKERRQ(ierr); 1815 if (size == 1) { 1816 ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 1817 ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 1818 } else { 1819 ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 1820 ierr = MatMPIAIJSetPreallocation(*mat,0,NULL,0,NULL);CHKERRQ(ierr); 1821 } 1822 1823 for (i=0; i<M; i++) { 1824 1825 ierr = VecSet(in,0.0);CHKERRQ(ierr); 1826 ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 1827 ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 1828 ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 1829 1830 /* should fix, allowing user to choose side */ 1831 ierr = PCApply(pc,in,out);CHKERRQ(ierr); 1832 1833 ierr = VecGetArray(out,&array);CHKERRQ(ierr); 1834 ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 1835 ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 1836 1837 } 1838 ierr = PetscFree(rows);CHKERRQ(ierr); 1839 ierr = VecDestroy(&out);CHKERRQ(ierr); 1840 ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1841 ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1842 PetscFunctionReturn(0); 1843 } 1844 1845 #undef __FUNCT__ 1846 #define __FUNCT__ "PCSetCoordinates" 1847 /*@ 1848 PCSetCoordinates - sets the coordinates of all the nodes on the local process 1849 1850 Collective on PC 1851 1852 Input Parameters: 1853 + pc - the solver context 1854 . dim - the dimension of the coordinates 1, 2, or 3 1855 - coords - the coordinates 1856 1857 Level: intermediate 1858 1859 Notes: coords is an array of the 3D coordinates for the nodes on 1860 the local processor. So if there are 108 equation on a processor 1861 for a displacement finite element discretization of elasticity (so 1862 that there are 36 = 108/3 nodes) then the array must have 108 1863 double precision values (ie, 3 * 36). These x y z coordinates 1864 should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x, 1865 ... , N-1.z ]. 1866 1867 .seealso: MatSetNearNullSpace 1868 @*/ 1869 PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords) 1870 { 1871 PetscErrorCode ierr; 1872 1873 PetscFunctionBegin; 1874 ierr = PetscTryMethod(pc,"PCSetCoordinates_C",(PC,PetscInt,PetscInt,PetscReal*),(pc,dim,nloc,coords));CHKERRQ(ierr); 1875 PetscFunctionReturn(0); 1876 } 1877