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