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