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