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 ierr = PCSetUp(pc);CHKERRQ(ierr); 479 if (!pc->ops->apply) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply"); 480 ierr = VecLockPush(x);CHKERRQ(ierr); 481 ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr); 482 ierr = (*pc->ops->apply)(pc,x,y);CHKERRQ(ierr); 483 ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr); 484 if (pc->erroriffailure) {ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);} 485 ierr = VecLockPop(x);CHKERRQ(ierr); 486 PetscFunctionReturn(0); 487 } 488 489 #undef __FUNCT__ 490 #define __FUNCT__ "PCApplySymmetricLeft" 491 /*@ 492 PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector. 493 494 Collective on PC and Vec 495 496 Input Parameters: 497 + pc - the preconditioner context 498 - x - input vector 499 500 Output Parameter: 501 . y - output vector 502 503 Notes: 504 Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners. 505 506 Level: developer 507 508 .keywords: PC, apply, symmetric, left 509 510 .seealso: PCApply(), PCApplySymmetricRight() 511 @*/ 512 PetscErrorCode PCApplySymmetricLeft(PC pc,Vec x,Vec y) 513 { 514 PetscErrorCode ierr; 515 516 PetscFunctionBegin; 517 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 518 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 519 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 520 if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors"); 521 if (pc->erroriffailure) {ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);} 522 ierr = PCSetUp(pc);CHKERRQ(ierr); 523 if (!pc->ops->applysymmetricleft) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have left symmetric apply"); 524 ierr = VecLockPush(x);CHKERRQ(ierr); 525 ierr = PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr); 526 ierr = (*pc->ops->applysymmetricleft)(pc,x,y);CHKERRQ(ierr); 527 ierr = PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr); 528 ierr = VecLockPop(x);CHKERRQ(ierr); 529 if (pc->erroriffailure) {ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);} 530 PetscFunctionReturn(0); 531 } 532 533 #undef __FUNCT__ 534 #define __FUNCT__ "PCApplySymmetricRight" 535 /*@ 536 PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector. 537 538 Collective on PC and Vec 539 540 Input Parameters: 541 + pc - the preconditioner context 542 - x - input vector 543 544 Output Parameter: 545 . y - output vector 546 547 Level: developer 548 549 Notes: 550 Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners. 551 552 .keywords: PC, apply, symmetric, right 553 554 .seealso: PCApply(), PCApplySymmetricLeft() 555 @*/ 556 PetscErrorCode PCApplySymmetricRight(PC pc,Vec x,Vec y) 557 { 558 PetscErrorCode ierr; 559 560 PetscFunctionBegin; 561 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 562 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 563 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 564 if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors"); 565 if (pc->erroriffailure) {ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);} 566 ierr = PCSetUp(pc);CHKERRQ(ierr); 567 if (!pc->ops->applysymmetricright) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have left symmetric apply"); 568 ierr = VecLockPush(x);CHKERRQ(ierr); 569 ierr = PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr); 570 ierr = (*pc->ops->applysymmetricright)(pc,x,y);CHKERRQ(ierr); 571 ierr = PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr); 572 ierr = VecLockPop(x);CHKERRQ(ierr); 573 if (pc->erroriffailure) {ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);} 574 PetscFunctionReturn(0); 575 } 576 577 #undef __FUNCT__ 578 #define __FUNCT__ "PCApplyTranspose" 579 /*@ 580 PCApplyTranspose - Applies the transpose of preconditioner to a vector. 581 582 Collective on PC and Vec 583 584 Input Parameters: 585 + pc - the preconditioner context 586 - x - input vector 587 588 Output Parameter: 589 . y - output vector 590 591 Notes: For complex numbers this applies the non-Hermitian transpose. 592 593 Developer Notes: We need to implement a PCApplyHermitianTranspose() 594 595 Level: developer 596 597 .keywords: PC, apply, transpose 598 599 .seealso: PCApply(), PCApplyBAorAB(), PCApplyBAorABTranspose(), PCApplyTransposeExists() 600 @*/ 601 PetscErrorCode PCApplyTranspose(PC pc,Vec x,Vec y) 602 { 603 PetscErrorCode ierr; 604 605 PetscFunctionBegin; 606 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 607 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 608 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 609 if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors"); 610 if (pc->erroriffailure) {ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr);} 611 ierr = PCSetUp(pc);CHKERRQ(ierr); 612 if (!pc->ops->applytranspose) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply transpose"); 613 ierr = VecLockPush(x);CHKERRQ(ierr); 614 ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr); 615 ierr = (*pc->ops->applytranspose)(pc,x,y);CHKERRQ(ierr); 616 ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr); 617 ierr = VecLockPop(x);CHKERRQ(ierr); 618 if (pc->erroriffailure) {ierr = VecValidValues(y,3,PETSC_FALSE);CHKERRQ(ierr);} 619 PetscFunctionReturn(0); 620 } 621 622 #undef __FUNCT__ 623 #define __FUNCT__ "PCApplyTransposeExists" 624 /*@ 625 PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation 626 627 Collective on PC and Vec 628 629 Input Parameters: 630 . pc - the preconditioner context 631 632 Output Parameter: 633 . flg - PETSC_TRUE if a transpose operation is defined 634 635 Level: developer 636 637 .keywords: PC, apply, transpose 638 639 .seealso: PCApplyTranspose() 640 @*/ 641 PetscErrorCode PCApplyTransposeExists(PC pc,PetscBool *flg) 642 { 643 PetscFunctionBegin; 644 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 645 PetscValidPointer(flg,2); 646 if (pc->ops->applytranspose) *flg = PETSC_TRUE; 647 else *flg = PETSC_FALSE; 648 PetscFunctionReturn(0); 649 } 650 651 #undef __FUNCT__ 652 #define __FUNCT__ "PCApplyBAorAB" 653 /*@ 654 PCApplyBAorAB - Applies the preconditioner and operator to a vector. y = B*A*x or y = A*B*x. 655 656 Collective on PC and Vec 657 658 Input Parameters: 659 + pc - the preconditioner context 660 . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC 661 . x - input vector 662 - work - work vector 663 664 Output Parameter: 665 . y - output vector 666 667 Level: developer 668 669 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 670 specific KSPSolve() method must also be written to handle the post-solve "correction" for the diagonal scaling. 671 672 .keywords: PC, apply, operator 673 674 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose() 675 @*/ 676 PetscErrorCode PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work) 677 { 678 PetscErrorCode ierr; 679 680 PetscFunctionBegin; 681 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 682 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 683 PetscValidHeaderSpecific(y,VEC_CLASSID,4); 684 PetscValidHeaderSpecific(work,VEC_CLASSID,5); 685 if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors"); 686 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"); 687 if (pc->diagonalscale && side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot include diagonal scaling with symmetric preconditioner application"); 688 if (pc->erroriffailure) {ierr = VecValidValues(x,3,PETSC_TRUE);CHKERRQ(ierr);} 689 690 ierr = PCSetUp(pc);CHKERRQ(ierr); 691 if (pc->diagonalscale) { 692 if (pc->ops->applyBA) { 693 Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */ 694 ierr = VecDuplicate(x,&work2);CHKERRQ(ierr); 695 ierr = PCDiagonalScaleRight(pc,x,work2);CHKERRQ(ierr); 696 ierr = (*pc->ops->applyBA)(pc,side,work2,y,work);CHKERRQ(ierr); 697 ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr); 698 ierr = VecDestroy(&work2);CHKERRQ(ierr); 699 } else if (side == PC_RIGHT) { 700 ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr); 701 ierr = PCApply(pc,y,work);CHKERRQ(ierr); 702 ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr); 703 ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr); 704 } else if (side == PC_LEFT) { 705 ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr); 706 ierr = MatMult(pc->mat,y,work);CHKERRQ(ierr); 707 ierr = PCApply(pc,work,y);CHKERRQ(ierr); 708 ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr); 709 } else if (side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot provide diagonal scaling with symmetric application of preconditioner"); 710 } else { 711 if (pc->ops->applyBA) { 712 ierr = (*pc->ops->applyBA)(pc,side,x,y,work);CHKERRQ(ierr); 713 } else if (side == PC_RIGHT) { 714 ierr = PCApply(pc,x,work);CHKERRQ(ierr); 715 ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr); 716 } else if (side == PC_LEFT) { 717 ierr = MatMult(pc->mat,x,work);CHKERRQ(ierr); 718 ierr = PCApply(pc,work,y);CHKERRQ(ierr); 719 } else if (side == PC_SYMMETRIC) { 720 /* There's an extra copy here; maybe should provide 2 work vectors instead? */ 721 ierr = PCApplySymmetricRight(pc,x,work);CHKERRQ(ierr); 722 ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr); 723 ierr = VecCopy(y,work);CHKERRQ(ierr); 724 ierr = PCApplySymmetricLeft(pc,work,y);CHKERRQ(ierr); 725 } 726 } 727 if (pc->erroriffailure) {ierr = VecValidValues(y,4,PETSC_FALSE);CHKERRQ(ierr);} 728 PetscFunctionReturn(0); 729 } 730 731 #undef __FUNCT__ 732 #define __FUNCT__ "PCApplyBAorABTranspose" 733 /*@ 734 PCApplyBAorABTranspose - Applies the transpose of the preconditioner 735 and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning, 736 NOT tr(B*A) = tr(A)*tr(B). 737 738 Collective on PC and Vec 739 740 Input Parameters: 741 + pc - the preconditioner context 742 . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC 743 . x - input vector 744 - work - work vector 745 746 Output Parameter: 747 . y - output vector 748 749 750 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 751 defined by B'. This is why this has the funny form that it computes tr(B) * tr(A) 752 753 Level: developer 754 755 .keywords: PC, apply, operator, transpose 756 757 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB() 758 @*/ 759 PetscErrorCode PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work) 760 { 761 PetscErrorCode ierr; 762 763 PetscFunctionBegin; 764 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 765 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 766 PetscValidHeaderSpecific(y,VEC_CLASSID,4); 767 PetscValidHeaderSpecific(work,VEC_CLASSID,5); 768 if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors"); 769 if (pc->erroriffailure) {ierr = VecValidValues(x,3,PETSC_TRUE);CHKERRQ(ierr);} 770 if (pc->ops->applyBAtranspose) { 771 ierr = (*pc->ops->applyBAtranspose)(pc,side,x,y,work);CHKERRQ(ierr); 772 if (pc->erroriffailure) {ierr = VecValidValues(y,4,PETSC_FALSE);CHKERRQ(ierr);} 773 PetscFunctionReturn(0); 774 } 775 if (side != PC_LEFT && side != PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left"); 776 777 ierr = PCSetUp(pc);CHKERRQ(ierr); 778 if (side == PC_RIGHT) { 779 ierr = PCApplyTranspose(pc,x,work);CHKERRQ(ierr); 780 ierr = MatMultTranspose(pc->mat,work,y);CHKERRQ(ierr); 781 } else if (side == PC_LEFT) { 782 ierr = MatMultTranspose(pc->mat,x,work);CHKERRQ(ierr); 783 ierr = PCApplyTranspose(pc,work,y);CHKERRQ(ierr); 784 } 785 /* add support for PC_SYMMETRIC */ 786 if (pc->erroriffailure) {ierr = VecValidValues(y,4,PETSC_FALSE);CHKERRQ(ierr);} 787 PetscFunctionReturn(0); 788 } 789 790 /* -------------------------------------------------------------------------------*/ 791 792 #undef __FUNCT__ 793 #define __FUNCT__ "PCApplyRichardsonExists" 794 /*@ 795 PCApplyRichardsonExists - Determines whether a particular preconditioner has a 796 built-in fast application of Richardson's method. 797 798 Not Collective 799 800 Input Parameter: 801 . pc - the preconditioner 802 803 Output Parameter: 804 . exists - PETSC_TRUE or PETSC_FALSE 805 806 Level: developer 807 808 .keywords: PC, apply, Richardson, exists 809 810 .seealso: PCApplyRichardson() 811 @*/ 812 PetscErrorCode PCApplyRichardsonExists(PC pc,PetscBool *exists) 813 { 814 PetscFunctionBegin; 815 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 816 PetscValidIntPointer(exists,2); 817 if (pc->ops->applyrichardson) *exists = PETSC_TRUE; 818 else *exists = PETSC_FALSE; 819 PetscFunctionReturn(0); 820 } 821 822 #undef __FUNCT__ 823 #define __FUNCT__ "PCApplyRichardson" 824 /*@ 825 PCApplyRichardson - Applies several steps of Richardson iteration with 826 the particular preconditioner. This routine is usually used by the 827 Krylov solvers and not the application code directly. 828 829 Collective on PC 830 831 Input Parameters: 832 + pc - the preconditioner context 833 . b - the right hand side 834 . w - one work vector 835 . rtol - relative decrease in residual norm convergence criteria 836 . abstol - absolute residual norm convergence criteria 837 . dtol - divergence residual norm increase criteria 838 . its - the number of iterations to apply. 839 - guesszero - if the input x contains nonzero initial guess 840 841 Output Parameter: 842 + outits - number of iterations actually used (for SOR this always equals its) 843 . reason - the reason the apply terminated 844 - y - the solution (also contains initial guess if guesszero is PETSC_FALSE 845 846 Notes: 847 Most preconditioners do not support this function. Use the command 848 PCApplyRichardsonExists() to determine if one does. 849 850 Except for the multigrid PC this routine ignores the convergence tolerances 851 and always runs for the number of iterations 852 853 Level: developer 854 855 .keywords: PC, apply, Richardson 856 857 .seealso: PCApplyRichardsonExists() 858 @*/ 859 PetscErrorCode PCApplyRichardson(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason) 860 { 861 PetscErrorCode ierr; 862 863 PetscFunctionBegin; 864 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 865 PetscValidHeaderSpecific(b,VEC_CLASSID,2); 866 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 867 PetscValidHeaderSpecific(w,VEC_CLASSID,4); 868 if (b == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"b and y must be different vectors"); 869 ierr = PCSetUp(pc);CHKERRQ(ierr); 870 if (!pc->ops->applyrichardson) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply richardson"); 871 ierr = (*pc->ops->applyrichardson)(pc,b,y,w,rtol,abstol,dtol,its,guesszero,outits,reason);CHKERRQ(ierr); 872 PetscFunctionReturn(0); 873 } 874 875 #undef __FUNCT__ 876 #define __FUNCT__ "PCGetSetUpFailedReason" 877 /*@ 878 PCGetSetUpFailedReason - Gets the reason a PCSetUp() failed or 0 if it did not fail 879 880 Logically Collective on PC 881 882 Input Parameter: 883 . pc - the preconditioner context 884 885 Output Parameter: 886 . reason - the reason it failed, currently only -1 887 888 Level: advanced 889 890 .keywords: PC, setup 891 892 .seealso: PCCreate(), PCApply(), PCDestroy() 893 @*/ 894 PetscErrorCode PCGetSetUpFailedReason(PC pc,PCFailedReason *reason) 895 { 896 PetscFunctionBegin; 897 if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled; 898 else *reason = pc->failedreason; 899 PetscFunctionReturn(0); 900 } 901 902 903 /* 904 a setupcall of 0 indicates never setup, 905 1 indicates has been previously setup 906 -1 indicates a PCSetUp() was attempted and failed 907 */ 908 #undef __FUNCT__ 909 #define __FUNCT__ "PCSetUp" 910 /*@ 911 PCSetUp - Prepares for the use of a preconditioner. 912 913 Collective on PC 914 915 Input Parameter: 916 . pc - the preconditioner context 917 918 Level: developer 919 920 .keywords: PC, setup 921 922 .seealso: PCCreate(), PCApply(), PCDestroy() 923 @*/ 924 PetscErrorCode PCSetUp(PC pc) 925 { 926 PetscErrorCode ierr; 927 const char *def; 928 PetscObjectState matstate, matnonzerostate; 929 930 PetscFunctionBegin; 931 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 932 if (!pc->mat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first"); 933 934 if (pc->setupcalled && pc->reusepreconditioner) { 935 ierr = PetscInfo(pc,"Leaving PC with identical preconditioner since reuse preconditioner is set\n");CHKERRQ(ierr); 936 PetscFunctionReturn(0); 937 } 938 939 ierr = PetscObjectStateGet((PetscObject)pc->pmat,&matstate);CHKERRQ(ierr); 940 ierr = MatGetNonzeroState(pc->pmat,&matnonzerostate);CHKERRQ(ierr); 941 if (!pc->setupcalled) { 942 ierr = PetscInfo(pc,"Setting up PC for first time\n");CHKERRQ(ierr); 943 pc->flag = DIFFERENT_NONZERO_PATTERN; 944 } else if (matstate == pc->matstate) { 945 ierr = PetscInfo(pc,"Leaving PC with identical preconditioner since operator is unchanged\n");CHKERRQ(ierr); 946 PetscFunctionReturn(0); 947 } else { 948 if (matnonzerostate > pc->matnonzerostate) { 949 ierr = PetscInfo(pc,"Setting up PC with different nonzero pattern\n");CHKERRQ(ierr); 950 pc->flag = DIFFERENT_NONZERO_PATTERN; 951 } else { 952 ierr = PetscInfo(pc,"Setting up PC with same nonzero pattern\n");CHKERRQ(ierr); 953 pc->flag = SAME_NONZERO_PATTERN; 954 } 955 } 956 pc->matstate = matstate; 957 pc->matnonzerostate = matnonzerostate; 958 959 if (!((PetscObject)pc)->type_name) { 960 ierr = PCGetDefaultType_Private(pc,&def);CHKERRQ(ierr); 961 ierr = PCSetType(pc,def);CHKERRQ(ierr); 962 } 963 964 ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr); 965 ierr = MatSetErrorIfFailure(pc->mat,pc->erroriffailure);CHKERRQ(ierr); 966 ierr = PetscLogEventBegin(PC_SetUp,pc,0,0,0);CHKERRQ(ierr); 967 if (pc->ops->setup) { 968 ierr = (*pc->ops->setup)(pc);CHKERRQ(ierr); 969 } 970 ierr = PetscLogEventEnd(PC_SetUp,pc,0,0,0);CHKERRQ(ierr); 971 if (!pc->setupcalled) pc->setupcalled = 1; 972 PetscFunctionReturn(0); 973 } 974 975 #undef __FUNCT__ 976 #define __FUNCT__ "PCSetUpOnBlocks" 977 /*@ 978 PCSetUpOnBlocks - Sets up the preconditioner for each block in 979 the block Jacobi, block Gauss-Seidel, and overlapping Schwarz 980 methods. 981 982 Collective on PC 983 984 Input Parameters: 985 . pc - the preconditioner context 986 987 Level: developer 988 989 .keywords: PC, setup, blocks 990 991 .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp() 992 @*/ 993 PetscErrorCode PCSetUpOnBlocks(PC pc) 994 { 995 PetscErrorCode ierr; 996 997 PetscFunctionBegin; 998 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 999 if (!pc->ops->setuponblocks) PetscFunctionReturn(0); 1000 ierr = PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr); 1001 ierr = (*pc->ops->setuponblocks)(pc);CHKERRQ(ierr); 1002 ierr = PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr); 1003 PetscFunctionReturn(0); 1004 } 1005 1006 #undef __FUNCT__ 1007 #define __FUNCT__ "PCSetModifySubMatrices" 1008 /*@C 1009 PCSetModifySubMatrices - Sets a user-defined routine for modifying the 1010 submatrices that arise within certain subdomain-based preconditioners. 1011 The basic submatrices are extracted from the preconditioner matrix as 1012 usual; the user can then alter these (for example, to set different boundary 1013 conditions for each submatrix) before they are used for the local solves. 1014 1015 Logically Collective on PC 1016 1017 Input Parameters: 1018 + pc - the preconditioner context 1019 . func - routine for modifying the submatrices 1020 - ctx - optional user-defined context (may be null) 1021 1022 Calling sequence of func: 1023 $ func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx); 1024 1025 . row - an array of index sets that contain the global row numbers 1026 that comprise each local submatrix 1027 . col - an array of index sets that contain the global column numbers 1028 that comprise each local submatrix 1029 . submat - array of local submatrices 1030 - ctx - optional user-defined context for private data for the 1031 user-defined func routine (may be null) 1032 1033 Notes: 1034 PCSetModifySubMatrices() MUST be called before KSPSetUp() and 1035 KSPSolve(). 1036 1037 A routine set by PCSetModifySubMatrices() is currently called within 1038 the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM) 1039 preconditioners. All other preconditioners ignore this routine. 1040 1041 Level: advanced 1042 1043 .keywords: PC, set, modify, submatrices 1044 1045 .seealso: PCModifySubMatrices(), PCASMGetSubMatrices() 1046 @*/ 1047 PetscErrorCode PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx) 1048 { 1049 PetscFunctionBegin; 1050 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1051 pc->modifysubmatrices = func; 1052 pc->modifysubmatricesP = ctx; 1053 PetscFunctionReturn(0); 1054 } 1055 1056 #undef __FUNCT__ 1057 #define __FUNCT__ "PCModifySubMatrices" 1058 /*@C 1059 PCModifySubMatrices - Calls an optional user-defined routine within 1060 certain preconditioners if one has been set with PCSetModifySubMarices(). 1061 1062 Collective on PC 1063 1064 Input Parameters: 1065 + pc - the preconditioner context 1066 . nsub - the number of local submatrices 1067 . row - an array of index sets that contain the global row numbers 1068 that comprise each local submatrix 1069 . col - an array of index sets that contain the global column numbers 1070 that comprise each local submatrix 1071 . submat - array of local submatrices 1072 - ctx - optional user-defined context for private data for the 1073 user-defined routine (may be null) 1074 1075 Output Parameter: 1076 . submat - array of local submatrices (the entries of which may 1077 have been modified) 1078 1079 Notes: 1080 The user should NOT generally call this routine, as it will 1081 automatically be called within certain preconditioners (currently 1082 block Jacobi, additive Schwarz) if set. 1083 1084 The basic submatrices are extracted from the preconditioner matrix 1085 as usual; the user can then alter these (for example, to set different 1086 boundary conditions for each submatrix) before they are used for the 1087 local solves. 1088 1089 Level: developer 1090 1091 .keywords: PC, modify, submatrices 1092 1093 .seealso: PCSetModifySubMatrices() 1094 @*/ 1095 PetscErrorCode PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx) 1096 { 1097 PetscErrorCode ierr; 1098 1099 PetscFunctionBegin; 1100 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1101 if (!pc->modifysubmatrices) PetscFunctionReturn(0); 1102 ierr = PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr); 1103 ierr = (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);CHKERRQ(ierr); 1104 ierr = PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr); 1105 PetscFunctionReturn(0); 1106 } 1107 1108 #undef __FUNCT__ 1109 #define __FUNCT__ "PCSetOperators" 1110 /*@ 1111 PCSetOperators - Sets the matrix associated with the linear system and 1112 a (possibly) different one associated with the preconditioner. 1113 1114 Logically Collective on PC and Mat 1115 1116 Input Parameters: 1117 + pc - the preconditioner context 1118 . Amat - the matrix that defines the linear system 1119 - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat. 1120 1121 Notes: 1122 Passing a NULL for Amat or Pmat removes the matrix that is currently used. 1123 1124 If you wish to replace either Amat or Pmat but leave the other one untouched then 1125 first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference() 1126 on it and then pass it back in in your call to KSPSetOperators(). 1127 1128 More Notes about Repeated Solution of Linear Systems: 1129 PETSc does NOT reset the matrix entries of either Amat or Pmat 1130 to zero after a linear solve; the user is completely responsible for 1131 matrix assembly. See the routine MatZeroEntries() if desiring to 1132 zero all elements of a matrix. 1133 1134 Level: intermediate 1135 1136 .keywords: PC, set, operators, matrix, linear system 1137 1138 .seealso: PCGetOperators(), MatZeroEntries() 1139 @*/ 1140 PetscErrorCode PCSetOperators(PC pc,Mat Amat,Mat Pmat) 1141 { 1142 PetscErrorCode ierr; 1143 PetscInt m1,n1,m2,n2; 1144 1145 PetscFunctionBegin; 1146 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1147 if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 1148 if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_CLASSID,3); 1149 if (Amat) PetscCheckSameComm(pc,1,Amat,2); 1150 if (Pmat) PetscCheckSameComm(pc,1,Pmat,3); 1151 if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) { 1152 ierr = MatGetLocalSize(Amat,&m1,&n1);CHKERRQ(ierr); 1153 ierr = MatGetLocalSize(pc->mat,&m2,&n2);CHKERRQ(ierr); 1154 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); 1155 ierr = MatGetLocalSize(Pmat,&m1,&n1);CHKERRQ(ierr); 1156 ierr = MatGetLocalSize(pc->pmat,&m2,&n2);CHKERRQ(ierr); 1157 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); 1158 } 1159 1160 if (Pmat != pc->pmat) { 1161 /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */ 1162 pc->matnonzerostate = -1; 1163 pc->matstate = -1; 1164 } 1165 1166 /* reference first in case the matrices are the same */ 1167 if (Amat) {ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);} 1168 ierr = MatDestroy(&pc->mat);CHKERRQ(ierr); 1169 if (Pmat) {ierr = PetscObjectReference((PetscObject)Pmat);CHKERRQ(ierr);} 1170 ierr = MatDestroy(&pc->pmat);CHKERRQ(ierr); 1171 pc->mat = Amat; 1172 pc->pmat = Pmat; 1173 PetscFunctionReturn(0); 1174 } 1175 1176 #undef __FUNCT__ 1177 #define __FUNCT__ "PCSetReusePreconditioner" 1178 /*@ 1179 PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed. 1180 1181 Logically Collective on PC 1182 1183 Input Parameters: 1184 + pc - the preconditioner context 1185 - flag - PETSC_TRUE do not compute a new preconditioner, PETSC_FALSE do compute a new preconditioner 1186 1187 Level: intermediate 1188 1189 .seealso: PCGetOperators(), MatZeroEntries(), PCGetReusePreconditioner(), KSPSetReusePreconditioner() 1190 @*/ 1191 PetscErrorCode PCSetReusePreconditioner(PC pc,PetscBool flag) 1192 { 1193 PetscFunctionBegin; 1194 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1195 pc->reusepreconditioner = flag; 1196 PetscFunctionReturn(0); 1197 } 1198 1199 #undef __FUNCT__ 1200 #define __FUNCT__ "PCGetReusePreconditioner" 1201 /*@ 1202 PCGetReusePreconditioner - Determines if the PC reuses the current preconditioner even if the operator in the preconditioner has changed. 1203 1204 Not Collective 1205 1206 Input Parameter: 1207 . pc - the preconditioner context 1208 1209 Output Parameter: 1210 . flag - PETSC_TRUE do not compute a new preconditioner, PETSC_FALSE do compute a new preconditioner 1211 1212 Level: intermediate 1213 1214 .seealso: PCGetOperators(), MatZeroEntries(), PCSetReusePreconditioner() 1215 @*/ 1216 PetscErrorCode PCGetReusePreconditioner(PC pc,PetscBool *flag) 1217 { 1218 PetscFunctionBegin; 1219 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1220 *flag = pc->reusepreconditioner; 1221 PetscFunctionReturn(0); 1222 } 1223 1224 #undef __FUNCT__ 1225 #define __FUNCT__ "PCGetOperators" 1226 /*@C 1227 PCGetOperators - Gets the matrix associated with the linear system and 1228 possibly a different one associated with the preconditioner. 1229 1230 Not collective, though parallel Mats are returned if the PC is parallel 1231 1232 Input Parameter: 1233 . pc - the preconditioner context 1234 1235 Output Parameters: 1236 + Amat - the matrix defining the linear system 1237 - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat. 1238 1239 Level: intermediate 1240 1241 Notes: Does not increase the reference count of the matrices, so you should not destroy them 1242 1243 Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators 1244 are created in PC and returned to the user. In this case, if both operators 1245 mat and pmat are requested, two DIFFERENT operators will be returned. If 1246 only one is requested both operators in the PC will be the same (i.e. as 1247 if one had called KSP/PCSetOperators() with the same argument for both Mats). 1248 The user must set the sizes of the returned matrices and their type etc just 1249 as if the user created them with MatCreate(). For example, 1250 1251 $ KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to 1252 $ set size, type, etc of Amat 1253 1254 $ MatCreate(comm,&mat); 1255 $ KSP/PCSetOperators(ksp/pc,Amat,Amat); 1256 $ PetscObjectDereference((PetscObject)mat); 1257 $ set size, type, etc of Amat 1258 1259 and 1260 1261 $ KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to 1262 $ set size, type, etc of Amat and Pmat 1263 1264 $ MatCreate(comm,&Amat); 1265 $ MatCreate(comm,&Pmat); 1266 $ KSP/PCSetOperators(ksp/pc,Amat,Pmat); 1267 $ PetscObjectDereference((PetscObject)Amat); 1268 $ PetscObjectDereference((PetscObject)Pmat); 1269 $ set size, type, etc of Amat and Pmat 1270 1271 The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy 1272 of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely 1273 managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look 1274 at this is when you create a SNES you do not NEED to create a KSP and attach it to 1275 the SNES object (the SNES object manages it for you). Similarly when you create a KSP 1276 you do not need to attach a PC to it (the KSP object manages the PC object for you). 1277 Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when 1278 it can be created for you? 1279 1280 1281 .keywords: PC, get, operators, matrix, linear system 1282 1283 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet() 1284 @*/ 1285 PetscErrorCode PCGetOperators(PC pc,Mat *Amat,Mat *Pmat) 1286 { 1287 PetscErrorCode ierr; 1288 1289 PetscFunctionBegin; 1290 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1291 if (Amat) { 1292 if (!pc->mat) { 1293 if (pc->pmat && !Pmat) { /* Apmat has been set, but user did not request it, so use for Amat */ 1294 pc->mat = pc->pmat; 1295 ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr); 1296 } else { /* both Amat and Pmat are empty */ 1297 ierr = MatCreate(PetscObjectComm((PetscObject)pc),&pc->mat);CHKERRQ(ierr); 1298 if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */ 1299 pc->pmat = pc->mat; 1300 ierr = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr); 1301 } 1302 } 1303 } 1304 *Amat = pc->mat; 1305 } 1306 if (Pmat) { 1307 if (!pc->pmat) { 1308 if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */ 1309 pc->pmat = pc->mat; 1310 ierr = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr); 1311 } else { 1312 ierr = MatCreate(PetscObjectComm((PetscObject)pc),&pc->pmat);CHKERRQ(ierr); 1313 if (!Amat) { /* user did NOT request Amat, so make same as Pmat */ 1314 pc->mat = pc->pmat; 1315 ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr); 1316 } 1317 } 1318 } 1319 *Pmat = pc->pmat; 1320 } 1321 PetscFunctionReturn(0); 1322 } 1323 1324 #undef __FUNCT__ 1325 #define __FUNCT__ "PCGetOperatorsSet" 1326 /*@C 1327 PCGetOperatorsSet - Determines if the matrix associated with the linear system and 1328 possibly a different one associated with the preconditioner have been set in the PC. 1329 1330 Not collective, though the results on all processes should be the same 1331 1332 Input Parameter: 1333 . pc - the preconditioner context 1334 1335 Output Parameters: 1336 + mat - the matrix associated with the linear system was set 1337 - pmat - matrix associated with the preconditioner was set, usually the same 1338 1339 Level: intermediate 1340 1341 .keywords: PC, get, operators, matrix, linear system 1342 1343 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators() 1344 @*/ 1345 PetscErrorCode PCGetOperatorsSet(PC pc,PetscBool *mat,PetscBool *pmat) 1346 { 1347 PetscFunctionBegin; 1348 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1349 if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE; 1350 if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE; 1351 PetscFunctionReturn(0); 1352 } 1353 1354 #undef __FUNCT__ 1355 #define __FUNCT__ "PCFactorGetMatrix" 1356 /*@ 1357 PCFactorGetMatrix - Gets the factored matrix from the 1358 preconditioner context. This routine is valid only for the LU, 1359 incomplete LU, Cholesky, and incomplete Cholesky methods. 1360 1361 Not Collective on PC though Mat is parallel if PC is parallel 1362 1363 Input Parameters: 1364 . pc - the preconditioner context 1365 1366 Output parameters: 1367 . mat - the factored matrix 1368 1369 Level: advanced 1370 1371 Notes: Does not increase the reference count for the matrix so DO NOT destroy it 1372 1373 .keywords: PC, get, factored, matrix 1374 @*/ 1375 PetscErrorCode PCFactorGetMatrix(PC pc,Mat *mat) 1376 { 1377 PetscErrorCode ierr; 1378 1379 PetscFunctionBegin; 1380 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1381 PetscValidPointer(mat,2); 1382 if (pc->ops->getfactoredmatrix) { 1383 ierr = (*pc->ops->getfactoredmatrix)(pc,mat);CHKERRQ(ierr); 1384 } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC type does not support getting factor matrix"); 1385 PetscFunctionReturn(0); 1386 } 1387 1388 #undef __FUNCT__ 1389 #define __FUNCT__ "PCSetOptionsPrefix" 1390 /*@C 1391 PCSetOptionsPrefix - Sets the prefix used for searching for all 1392 PC options in the database. 1393 1394 Logically Collective on PC 1395 1396 Input Parameters: 1397 + pc - the preconditioner context 1398 - prefix - the prefix string to prepend to all PC option requests 1399 1400 Notes: 1401 A hyphen (-) must NOT be given at the beginning of the prefix name. 1402 The first character of all runtime options is AUTOMATICALLY the 1403 hyphen. 1404 1405 Level: advanced 1406 1407 .keywords: PC, set, options, prefix, database 1408 1409 .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix() 1410 @*/ 1411 PetscErrorCode PCSetOptionsPrefix(PC pc,const char prefix[]) 1412 { 1413 PetscErrorCode ierr; 1414 1415 PetscFunctionBegin; 1416 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1417 ierr = PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr); 1418 PetscFunctionReturn(0); 1419 } 1420 1421 #undef __FUNCT__ 1422 #define __FUNCT__ "PCAppendOptionsPrefix" 1423 /*@C 1424 PCAppendOptionsPrefix - Appends to the prefix used for searching for all 1425 PC options in the database. 1426 1427 Logically Collective on PC 1428 1429 Input Parameters: 1430 + pc - the preconditioner context 1431 - prefix - the prefix string to prepend to all PC option requests 1432 1433 Notes: 1434 A hyphen (-) must NOT be given at the beginning of the prefix name. 1435 The first character of all runtime options is AUTOMATICALLY the 1436 hyphen. 1437 1438 Level: advanced 1439 1440 .keywords: PC, append, options, prefix, database 1441 1442 .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix() 1443 @*/ 1444 PetscErrorCode PCAppendOptionsPrefix(PC pc,const char prefix[]) 1445 { 1446 PetscErrorCode ierr; 1447 1448 PetscFunctionBegin; 1449 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1450 ierr = PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr); 1451 PetscFunctionReturn(0); 1452 } 1453 1454 #undef __FUNCT__ 1455 #define __FUNCT__ "PCGetOptionsPrefix" 1456 /*@C 1457 PCGetOptionsPrefix - Gets the prefix used for searching for all 1458 PC options in the database. 1459 1460 Not Collective 1461 1462 Input Parameters: 1463 . pc - the preconditioner context 1464 1465 Output Parameters: 1466 . prefix - pointer to the prefix string used, is returned 1467 1468 Notes: On the fortran side, the user should pass in a string 'prifix' of 1469 sufficient length to hold the prefix. 1470 1471 Level: advanced 1472 1473 .keywords: PC, get, options, prefix, database 1474 1475 .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix() 1476 @*/ 1477 PetscErrorCode PCGetOptionsPrefix(PC pc,const char *prefix[]) 1478 { 1479 PetscErrorCode ierr; 1480 1481 PetscFunctionBegin; 1482 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1483 PetscValidPointer(prefix,2); 1484 ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr); 1485 PetscFunctionReturn(0); 1486 } 1487 1488 #undef __FUNCT__ 1489 #define __FUNCT__ "PCPreSolve" 1490 /*@ 1491 PCPreSolve - Optional pre-solve phase, intended for any 1492 preconditioner-specific actions that must be performed before 1493 the iterative solve itself. 1494 1495 Collective on PC 1496 1497 Input Parameters: 1498 + pc - the preconditioner context 1499 - ksp - the Krylov subspace context 1500 1501 Level: developer 1502 1503 Sample of Usage: 1504 .vb 1505 PCPreSolve(pc,ksp); 1506 KSPSolve(ksp,b,x); 1507 PCPostSolve(pc,ksp); 1508 .ve 1509 1510 Notes: 1511 The pre-solve phase is distinct from the PCSetUp() phase. 1512 1513 KSPSolve() calls this directly, so is rarely called by the user. 1514 1515 .keywords: PC, pre-solve 1516 1517 .seealso: PCPostSolve() 1518 @*/ 1519 PetscErrorCode PCPreSolve(PC pc,KSP ksp) 1520 { 1521 PetscErrorCode ierr; 1522 Vec x,rhs; 1523 1524 PetscFunctionBegin; 1525 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1526 PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 1527 pc->presolvedone++; 1528 if (pc->presolvedone > 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot embed PCPreSolve() more than twice"); 1529 ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr); 1530 ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr); 1531 1532 if (pc->ops->presolve) { 1533 ierr = (*pc->ops->presolve)(pc,ksp,rhs,x);CHKERRQ(ierr); 1534 } 1535 PetscFunctionReturn(0); 1536 } 1537 1538 #undef __FUNCT__ 1539 #define __FUNCT__ "PCPostSolve" 1540 /*@ 1541 PCPostSolve - Optional post-solve phase, intended for any 1542 preconditioner-specific actions that must be performed after 1543 the iterative solve itself. 1544 1545 Collective on PC 1546 1547 Input Parameters: 1548 + pc - the preconditioner context 1549 - ksp - the Krylov subspace context 1550 1551 Sample of Usage: 1552 .vb 1553 PCPreSolve(pc,ksp); 1554 KSPSolve(ksp,b,x); 1555 PCPostSolve(pc,ksp); 1556 .ve 1557 1558 Note: 1559 KSPSolve() calls this routine directly, so it is rarely called by the user. 1560 1561 Level: developer 1562 1563 .keywords: PC, post-solve 1564 1565 .seealso: PCPreSolve(), KSPSolve() 1566 @*/ 1567 PetscErrorCode PCPostSolve(PC pc,KSP ksp) 1568 { 1569 PetscErrorCode ierr; 1570 Vec x,rhs; 1571 1572 PetscFunctionBegin; 1573 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1574 PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 1575 pc->presolvedone--; 1576 ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr); 1577 ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr); 1578 if (pc->ops->postsolve) { 1579 ierr = (*pc->ops->postsolve)(pc,ksp,rhs,x);CHKERRQ(ierr); 1580 } 1581 PetscFunctionReturn(0); 1582 } 1583 1584 #undef __FUNCT__ 1585 #define __FUNCT__ "PCLoad" 1586 /*@C 1587 PCLoad - Loads a PC that has been stored in binary with PCView(). 1588 1589 Collective on PetscViewer 1590 1591 Input Parameters: 1592 + newdm - the newly loaded PC, this needs to have been created with PCCreate() or 1593 some related function before a call to PCLoad(). 1594 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() 1595 1596 Level: intermediate 1597 1598 Notes: 1599 The type is determined by the data in the file, any type set into the PC before this call is ignored. 1600 1601 Notes for advanced users: 1602 Most users should not need to know the details of the binary storage 1603 format, since PCLoad() and PCView() completely hide these details. 1604 But for anyone who's interested, the standard binary matrix storage 1605 format is 1606 .vb 1607 has not yet been determined 1608 .ve 1609 1610 .seealso: PetscViewerBinaryOpen(), PCView(), MatLoad(), VecLoad() 1611 @*/ 1612 PetscErrorCode PCLoad(PC newdm, PetscViewer viewer) 1613 { 1614 PetscErrorCode ierr; 1615 PetscBool isbinary; 1616 PetscInt classid; 1617 char type[256]; 1618 1619 PetscFunctionBegin; 1620 PetscValidHeaderSpecific(newdm,PC_CLASSID,1); 1621 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 1622 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 1623 if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 1624 1625 ierr = PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);CHKERRQ(ierr); 1626 if (classid != PC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not PC next in file"); 1627 ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr); 1628 ierr = PCSetType(newdm, type);CHKERRQ(ierr); 1629 if (newdm->ops->load) { 1630 ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr); 1631 } 1632 PetscFunctionReturn(0); 1633 } 1634 1635 #include <petscdraw.h> 1636 #if defined(PETSC_HAVE_SAWS) 1637 #include <petscviewersaws.h> 1638 #endif 1639 #undef __FUNCT__ 1640 #define __FUNCT__ "PCView" 1641 /*@C 1642 PCView - Prints the PC data structure. 1643 1644 Collective on PC 1645 1646 Input Parameters: 1647 + PC - the PC context 1648 - viewer - optional visualization context 1649 1650 Note: 1651 The available visualization contexts include 1652 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 1653 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 1654 output where only the first processor opens 1655 the file. All other processors send their 1656 data to the first processor to print. 1657 1658 The user can open an alternative visualization contexts with 1659 PetscViewerASCIIOpen() (output to a specified file). 1660 1661 Level: developer 1662 1663 .keywords: PC, view 1664 1665 .seealso: KSPView(), PetscViewerASCIIOpen() 1666 @*/ 1667 PetscErrorCode PCView(PC pc,PetscViewer viewer) 1668 { 1669 PCType cstr; 1670 PetscErrorCode ierr; 1671 PetscBool iascii,isstring,isbinary,isdraw; 1672 PetscViewerFormat format; 1673 #if defined(PETSC_HAVE_SAWS) 1674 PetscBool issaws; 1675 #endif 1676 1677 PetscFunctionBegin; 1678 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1679 if (!viewer) { 1680 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);CHKERRQ(ierr); 1681 } 1682 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 1683 PetscCheckSameComm(pc,1,viewer,2); 1684 1685 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1686 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 1687 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 1688 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 1689 #if defined(PETSC_HAVE_SAWS) 1690 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);CHKERRQ(ierr); 1691 #endif 1692 1693 if (iascii) { 1694 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1695 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)pc,viewer);CHKERRQ(ierr); 1696 if (!pc->setupcalled) { 1697 ierr = PetscViewerASCIIPrintf(viewer," PC has not been set up so information may be incomplete\n");CHKERRQ(ierr); 1698 } 1699 if (pc->ops->view) { 1700 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1701 ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr); 1702 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1703 } 1704 if (pc->mat) { 1705 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 1706 if (pc->pmat == pc->mat) { 1707 ierr = PetscViewerASCIIPrintf(viewer," linear system matrix = precond matrix:\n");CHKERRQ(ierr); 1708 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1709 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 1710 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1711 } else { 1712 if (pc->pmat) { 1713 ierr = PetscViewerASCIIPrintf(viewer," linear system matrix followed by preconditioner matrix:\n");CHKERRQ(ierr); 1714 } else { 1715 ierr = PetscViewerASCIIPrintf(viewer," linear system matrix:\n");CHKERRQ(ierr); 1716 } 1717 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1718 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 1719 if (pc->pmat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);} 1720 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1721 } 1722 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1723 } 1724 } else if (isstring) { 1725 ierr = PCGetType(pc,&cstr);CHKERRQ(ierr); 1726 ierr = PetscViewerStringSPrintf(viewer," %-7.7s",cstr);CHKERRQ(ierr); 1727 if (pc->ops->view) {ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);} 1728 } else if (isbinary) { 1729 PetscInt classid = PC_FILE_CLASSID; 1730 MPI_Comm comm; 1731 PetscMPIInt rank; 1732 char type[256]; 1733 1734 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1735 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1736 if (!rank) { 1737 ierr = PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 1738 ierr = PetscStrncpy(type,((PetscObject)pc)->type_name,256);CHKERRQ(ierr); 1739 ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 1740 } 1741 if (pc->ops->view) { 1742 ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr); 1743 } 1744 } else if (isdraw) { 1745 PetscDraw draw; 1746 char str[25]; 1747 PetscReal x,y,bottom,h; 1748 PetscInt n; 1749 1750 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1751 ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 1752 if (pc->mat) { 1753 ierr = MatGetSize(pc->mat,&n,NULL);CHKERRQ(ierr); 1754 ierr = PetscSNPrintf(str,25,"PC: %s (%D)",((PetscObject)pc)->type_name,n);CHKERRQ(ierr); 1755 } else { 1756 ierr = PetscSNPrintf(str,25,"PC: %s",((PetscObject)pc)->type_name);CHKERRQ(ierr); 1757 } 1758 ierr = PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr); 1759 bottom = y - h; 1760 ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 1761 if (pc->ops->view) { 1762 ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr); 1763 } 1764 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 1765 #if defined(PETSC_HAVE_SAWS) 1766 } else if (issaws) { 1767 PetscMPIInt rank; 1768 1769 ierr = PetscObjectName((PetscObject)pc);CHKERRQ(ierr); 1770 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 1771 if (!((PetscObject)pc)->amsmem && !rank) { 1772 ierr = PetscObjectViewSAWs((PetscObject)pc,viewer);CHKERRQ(ierr); 1773 } 1774 if (pc->mat) {ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);} 1775 if (pc->pmat && pc->pmat != pc->mat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);} 1776 #endif 1777 } 1778 PetscFunctionReturn(0); 1779 } 1780 1781 #undef __FUNCT__ 1782 #define __FUNCT__ "PCRegister" 1783 /*@C 1784 PCRegister - Adds a method to the preconditioner package. 1785 1786 Not collective 1787 1788 Input Parameters: 1789 + name_solver - name of a new user-defined solver 1790 - routine_create - routine to create method context 1791 1792 Notes: 1793 PCRegister() may be called multiple times to add several user-defined preconditioners. 1794 1795 Sample usage: 1796 .vb 1797 PCRegister("my_solver", MySolverCreate); 1798 .ve 1799 1800 Then, your solver can be chosen with the procedural interface via 1801 $ PCSetType(pc,"my_solver") 1802 or at runtime via the option 1803 $ -pc_type my_solver 1804 1805 Level: advanced 1806 1807 .keywords: PC, register 1808 1809 .seealso: PCRegisterAll(), PCRegisterDestroy() 1810 @*/ 1811 PetscErrorCode PCRegister(const char sname[],PetscErrorCode (*function)(PC)) 1812 { 1813 PetscErrorCode ierr; 1814 1815 PetscFunctionBegin; 1816 ierr = PetscFunctionListAdd(&PCList,sname,function);CHKERRQ(ierr); 1817 PetscFunctionReturn(0); 1818 } 1819 1820 #undef __FUNCT__ 1821 #define __FUNCT__ "PCComputeExplicitOperator" 1822 /*@ 1823 PCComputeExplicitOperator - Computes the explicit preconditioned operator. 1824 1825 Collective on PC 1826 1827 Input Parameter: 1828 . pc - the preconditioner object 1829 1830 Output Parameter: 1831 . mat - the explict preconditioned operator 1832 1833 Notes: 1834 This computation is done by applying the operators to columns of the 1835 identity matrix. 1836 1837 Currently, this routine uses a dense matrix format when 1 processor 1838 is used and a sparse format otherwise. This routine is costly in general, 1839 and is recommended for use only with relatively small systems. 1840 1841 Level: advanced 1842 1843 .keywords: PC, compute, explicit, operator 1844 1845 .seealso: KSPComputeExplicitOperator() 1846 1847 @*/ 1848 PetscErrorCode PCComputeExplicitOperator(PC pc,Mat *mat) 1849 { 1850 Vec in,out; 1851 PetscErrorCode ierr; 1852 PetscInt i,M,m,*rows,start,end; 1853 PetscMPIInt size; 1854 MPI_Comm comm; 1855 PetscScalar *array,one = 1.0; 1856 1857 PetscFunctionBegin; 1858 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1859 PetscValidPointer(mat,2); 1860 1861 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1862 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1863 1864 if (!pc->pmat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call"); 1865 ierr = MatCreateVecs(pc->pmat,&in,0);CHKERRQ(ierr); 1866 ierr = VecDuplicate(in,&out);CHKERRQ(ierr); 1867 ierr = VecGetOwnershipRange(in,&start,&end);CHKERRQ(ierr); 1868 ierr = VecGetSize(in,&M);CHKERRQ(ierr); 1869 ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr); 1870 ierr = PetscMalloc1(m+1,&rows);CHKERRQ(ierr); 1871 for (i=0; i<m; i++) rows[i] = start + i; 1872 1873 ierr = MatCreate(comm,mat);CHKERRQ(ierr); 1874 ierr = MatSetSizes(*mat,m,m,M,M);CHKERRQ(ierr); 1875 if (size == 1) { 1876 ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 1877 ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 1878 } else { 1879 ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 1880 ierr = MatMPIAIJSetPreallocation(*mat,0,NULL,0,NULL);CHKERRQ(ierr); 1881 } 1882 ierr = MatSetOption(*mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 1883 1884 for (i=0; i<M; i++) { 1885 1886 ierr = VecSet(in,0.0);CHKERRQ(ierr); 1887 ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 1888 ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 1889 ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 1890 1891 /* should fix, allowing user to choose side */ 1892 ierr = PCApply(pc,in,out);CHKERRQ(ierr); 1893 1894 ierr = VecGetArray(out,&array);CHKERRQ(ierr); 1895 ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 1896 ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 1897 1898 } 1899 ierr = PetscFree(rows);CHKERRQ(ierr); 1900 ierr = VecDestroy(&out);CHKERRQ(ierr); 1901 ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1902 ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1903 PetscFunctionReturn(0); 1904 } 1905 1906 #undef __FUNCT__ 1907 #define __FUNCT__ "PCSetCoordinates" 1908 /*@ 1909 PCSetCoordinates - sets the coordinates of all the nodes on the local process 1910 1911 Collective on PC 1912 1913 Input Parameters: 1914 + pc - the solver context 1915 . dim - the dimension of the coordinates 1, 2, or 3 1916 - coords - the coordinates 1917 1918 Level: intermediate 1919 1920 Notes: coords is an array of the 3D coordinates for the nodes on 1921 the local processor. So if there are 108 equation on a processor 1922 for a displacement finite element discretization of elasticity (so 1923 that there are 36 = 108/3 nodes) then the array must have 108 1924 double precision values (ie, 3 * 36). These x y z coordinates 1925 should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x, 1926 ... , N-1.z ]. 1927 1928 .seealso: MatSetNearNullSpace 1929 @*/ 1930 PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords) 1931 { 1932 PetscErrorCode ierr; 1933 1934 PetscFunctionBegin; 1935 ierr = PetscTryMethod(pc,"PCSetCoordinates_C",(PC,PetscInt,PetscInt,PetscReal*),(pc,dim,nloc,coords));CHKERRQ(ierr); 1936 PetscFunctionReturn(0); 1937 } 1938