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