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