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