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