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