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,1); 85 PetscValidHeaderSpecific(nullsp,MAT_NULLSPACE_COOKIE,2); 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,1); 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,1); 160 PetscValidPointer(flag,2); 161 *flag = pc->diagonalscale; 162 PetscFunctionReturn(0); 163 } 164 165 #undef __FUNCT__ 166 #define __FUNCT__ "PCDiagonalScaleSet" 167 /*@ 168 PCDiagonalScaleSet - Indicates the left scaling to use to apply an additional left and right 169 scaling as needed by certain time-stepping codes. 170 171 Collective on PC 172 173 Input Parameters: 174 + pc - the preconditioner context 175 - s - scaling vector 176 177 Level: intermediate 178 179 Notes: The system solved via the Krylov method is 180 $ D M A D^{-1} y = D M b for left preconditioning or 181 $ D A M D^{-1} z = D b for right preconditioning 182 183 PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}. 184 185 .keywords: PC 186 187 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScale() 188 @*/ 189 int PCDiagonalScaleSet(PC pc,Vec s) 190 { 191 int ierr; 192 193 PetscFunctionBegin; 194 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 195 PetscValidHeaderSpecific(s,VEC_COOKIE,2); 196 pc->diagonalscale = PETSC_TRUE; 197 if (pc->diagonalscaleleft) { 198 ierr = VecDestroy(pc->diagonalscaleleft);CHKERRQ(ierr); 199 } 200 pc->diagonalscaleleft = s; 201 ierr = PetscObjectReference((PetscObject)s);CHKERRQ(ierr); 202 if (!pc->diagonalscaleright) { 203 ierr = VecDuplicate(s,&pc->diagonalscaleright);CHKERRQ(ierr); 204 } 205 ierr = VecCopy(s,pc->diagonalscaleright);CHKERRQ(ierr); 206 ierr = VecReciprocal(pc->diagonalscaleright);CHKERRQ(ierr); 207 PetscFunctionReturn(0); 208 } 209 210 #undef __FUNCT__ 211 #define __FUNCT__ "PCDiagonalScaleLeft" 212 /*@C 213 PCDiagonalScaleLeft - Indicates the left scaling to use to apply an additional left and right 214 scaling as needed by certain time-stepping codes. 215 216 Collective on PC 217 218 Input Parameters: 219 + pc - the preconditioner context 220 . in - input vector 221 + out - scaled vector (maybe the same as in) 222 223 Level: intermediate 224 225 Notes: The system solved via the Krylov method is 226 $ D M A D^{-1} y = D M b for left preconditioning or 227 $ D A M D^{-1} z = D b for right preconditioning 228 229 PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}. 230 231 If diagonal scaling is turned off and in is not out then in is copied to out 232 233 .keywords: PC 234 235 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale() 236 @*/ 237 int PCDiagonalScaleLeft(PC pc,Vec in,Vec out) 238 { 239 int ierr; 240 241 PetscFunctionBegin; 242 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 243 PetscValidHeaderSpecific(in,VEC_COOKIE,2); 244 PetscValidHeaderSpecific(out,VEC_COOKIE,3); 245 if (pc->diagonalscale) { 246 ierr = VecPointwiseMult(pc->diagonalscaleleft,in,out);CHKERRQ(ierr); 247 } else if (in != out) { 248 ierr = VecCopy(in,out);CHKERRQ(ierr); 249 } 250 PetscFunctionReturn(0); 251 } 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "PCDiagonalScaleRight" 255 /*@C 256 PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes. 257 258 Collective on PC 259 260 Input Parameters: 261 + pc - the preconditioner context 262 . in - input vector 263 + out - scaled vector (maybe the same as in) 264 265 Level: intermediate 266 267 Notes: The system solved via the Krylov method is 268 $ D M A D^{-1} y = D M b for left preconditioning or 269 $ D A M D^{-1} z = D b for right preconditioning 270 271 PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}. 272 273 If diagonal scaling is turned off and in is not out then in is copied to out 274 275 .keywords: PC 276 277 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale() 278 @*/ 279 int PCDiagonalScaleRight(PC pc,Vec in,Vec out) 280 { 281 int ierr; 282 283 PetscFunctionBegin; 284 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 285 PetscValidHeaderSpecific(in,VEC_COOKIE,2); 286 PetscValidHeaderSpecific(out,VEC_COOKIE,3); 287 if (pc->diagonalscale) { 288 ierr = VecPointwiseMult(pc->diagonalscaleright,in,out);CHKERRQ(ierr); 289 } else if (in != out) { 290 ierr = VecCopy(in,out);CHKERRQ(ierr); 291 } 292 PetscFunctionReturn(0); 293 } 294 295 #undef __FUNCT__ 296 #define __FUNCT__ "PCPublish_Petsc" 297 static int PCPublish_Petsc(PetscObject obj) 298 { 299 #if defined(PETSC_HAVE_AMS) 300 PC v = (PC) obj; 301 int ierr; 302 #endif 303 304 PetscFunctionBegin; 305 306 #if defined(PETSC_HAVE_AMS) 307 /* if it is already published then return */ 308 if (v->amem >=0) PetscFunctionReturn(0); 309 310 ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr); 311 ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr); 312 #endif 313 314 PetscFunctionReturn(0); 315 } 316 317 #undef __FUNCT__ 318 #define __FUNCT__ "PCCreate" 319 /*@C 320 PCCreate - Creates a preconditioner context. 321 322 Collective on MPI_Comm 323 324 Input Parameter: 325 . comm - MPI communicator 326 327 Output Parameter: 328 . pc - location to put the preconditioner context 329 330 Notes: 331 The default preconditioner on one processor is PCILU with 0 fill on more 332 then one it is PCBJACOBI with ILU() on each processor. 333 334 Level: developer 335 336 .keywords: PC, create, context 337 338 .seealso: PCSetUp(), PCApply(), PCDestroy() 339 @*/ 340 int PCCreate(MPI_Comm comm,PC *newpc) 341 { 342 PC pc; 343 int ierr; 344 345 PetscFunctionBegin; 346 PetscValidPointer(newpc,1) 347 *newpc = 0; 348 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 349 ierr = PCInitializePackage(PETSC_NULL); CHKERRQ(ierr); 350 #endif 351 352 PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_COOKIE,-1,"PC",comm,PCDestroy,PCView); 353 PetscLogObjectCreate(pc); 354 pc->bops->publish = PCPublish_Petsc; 355 pc->vec = 0; 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->vec) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Vector must be set first");} 870 if (!pc->mat) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");} 871 872 if (!pc->type_name) { 873 ierr = PCGetDefaultType_Private(pc,&def);CHKERRQ(ierr); 874 ierr = PCSetType(pc,def);CHKERRQ(ierr); 875 } 876 877 if (pc->ops->setup) { 878 ierr = (*pc->ops->setup)(pc);CHKERRQ(ierr); 879 } 880 pc->setupcalled = 2; 881 if (pc->nullsp) { 882 PetscTruth test; 883 ierr = PetscOptionsHasName(pc->prefix,"-pc_test_null_space",&test);CHKERRQ(ierr); 884 if (test) { 885 ierr = MatNullSpaceTest(pc->nullsp,pc->mat);CHKERRQ(ierr); 886 } 887 } 888 ierr = PetscLogEventEnd(PC_SetUp,pc,0,0,0);CHKERRQ(ierr); 889 PetscFunctionReturn(0); 890 } 891 892 #undef __FUNCT__ 893 #define __FUNCT__ "PCSetUpOnBlocks" 894 /*@ 895 PCSetUpOnBlocks - Sets up the preconditioner for each block in 896 the block Jacobi, block Gauss-Seidel, and overlapping Schwarz 897 methods. 898 899 Collective on PC 900 901 Input Parameters: 902 . pc - the preconditioner context 903 904 Level: developer 905 906 .keywords: PC, setup, blocks 907 908 .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp() 909 @*/ 910 int PCSetUpOnBlocks(PC pc) 911 { 912 int ierr; 913 914 PetscFunctionBegin; 915 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 916 if (!pc->ops->setuponblocks) PetscFunctionReturn(0); 917 ierr = PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr); 918 ierr = (*pc->ops->setuponblocks)(pc);CHKERRQ(ierr); 919 ierr = PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr); 920 PetscFunctionReturn(0); 921 } 922 923 #undef __FUNCT__ 924 #define __FUNCT__ "PCSetModifySubMatrices" 925 /*@C 926 PCSetModifySubMatrices - Sets a user-defined routine for modifying the 927 submatrices that arise within certain subdomain-based preconditioners. 928 The basic submatrices are extracted from the preconditioner matrix as 929 usual; the user can then alter these (for example, to set different boundary 930 conditions for each submatrix) before they are used for the local solves. 931 932 Collective on PC 933 934 Input Parameters: 935 + pc - the preconditioner context 936 . func - routine for modifying the submatrices 937 - ctx - optional user-defined context (may be null) 938 939 Calling sequence of func: 940 $ func (PC pc,int nsub,IS *row,IS *col,Mat *submat,void *ctx); 941 942 . row - an array of index sets that contain the global row numbers 943 that comprise each local submatrix 944 . col - an array of index sets that contain the global column numbers 945 that comprise each local submatrix 946 . submat - array of local submatrices 947 - ctx - optional user-defined context for private data for the 948 user-defined func routine (may be null) 949 950 Notes: 951 PCSetModifySubMatrices() MUST be called before KSPSetUp() and 952 KSPSolve(). 953 954 A routine set by PCSetModifySubMatrices() is currently called within 955 the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM) 956 preconditioners. All other preconditioners ignore this routine. 957 958 Level: advanced 959 960 .keywords: PC, set, modify, submatrices 961 962 .seealso: PCModifySubMatrices() 963 @*/ 964 int PCSetModifySubMatrices(PC pc,int(*func)(PC,int,const IS[],const IS[],Mat[],void*),void *ctx) 965 { 966 PetscFunctionBegin; 967 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 968 pc->modifysubmatrices = func; 969 pc->modifysubmatricesP = ctx; 970 PetscFunctionReturn(0); 971 } 972 973 #undef __FUNCT__ 974 #define __FUNCT__ "PCModifySubMatrices" 975 /*@C 976 PCModifySubMatrices - Calls an optional user-defined routine within 977 certain preconditioners if one has been set with PCSetModifySubMarices(). 978 979 Collective on PC 980 981 Input Parameters: 982 + pc - the preconditioner context 983 . nsub - the number of local submatrices 984 . row - an array of index sets that contain the global row numbers 985 that comprise each local submatrix 986 . col - an array of index sets that contain the global column numbers 987 that comprise each local submatrix 988 . submat - array of local submatrices 989 - ctx - optional user-defined context for private data for the 990 user-defined routine (may be null) 991 992 Output Parameter: 993 . submat - array of local submatrices (the entries of which may 994 have been modified) 995 996 Notes: 997 The user should NOT generally call this routine, as it will 998 automatically be called within certain preconditioners (currently 999 block Jacobi, additive Schwarz) if set. 1000 1001 The basic submatrices are extracted from the preconditioner matrix 1002 as usual; the user can then alter these (for example, to set different 1003 boundary conditions for each submatrix) before they are used for the 1004 local solves. 1005 1006 Level: developer 1007 1008 .keywords: PC, modify, submatrices 1009 1010 .seealso: PCSetModifySubMatrices() 1011 @*/ 1012 int PCModifySubMatrices(PC pc,int nsub,const IS row[],const IS col[],Mat submat[],void *ctx) 1013 { 1014 int ierr; 1015 1016 PetscFunctionBegin; 1017 if (!pc->modifysubmatrices) PetscFunctionReturn(0); 1018 ierr = PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr); 1019 ierr = (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);CHKERRQ(ierr); 1020 ierr = PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr); 1021 PetscFunctionReturn(0); 1022 } 1023 1024 #undef __FUNCT__ 1025 #define __FUNCT__ "PCSetOperators" 1026 /*@ 1027 PCSetOperators - Sets the matrix associated with the linear system and 1028 a (possibly) different one associated with the preconditioner. 1029 1030 Collective on PC and Mat 1031 1032 Input Parameters: 1033 + pc - the preconditioner context 1034 . Amat - the matrix associated with the linear system 1035 . Pmat - the matrix to be used in constructing the preconditioner, usually 1036 the same as Amat. 1037 - flag - flag indicating information about the preconditioner matrix structure 1038 during successive linear solves. This flag is ignored the first time a 1039 linear system is solved, and thus is irrelevant when solving just one linear 1040 system. 1041 1042 Notes: 1043 The flag can be used to eliminate unnecessary work in the preconditioner 1044 during the repeated solution of linear systems of the same size. The 1045 available options are 1046 + SAME_PRECONDITIONER - 1047 Pmat is identical during successive linear solves. 1048 This option is intended for folks who are using 1049 different Amat and Pmat matrices and wish to reuse the 1050 same preconditioner matrix. For example, this option 1051 saves work by not recomputing incomplete factorization 1052 for ILU/ICC preconditioners. 1053 . SAME_NONZERO_PATTERN - 1054 Pmat has the same nonzero structure during 1055 successive linear solves. 1056 - DIFFERENT_NONZERO_PATTERN - 1057 Pmat does not have the same nonzero structure. 1058 1059 Caution: 1060 If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion 1061 and does not check the structure of the matrix. If you erroneously 1062 claim that the structure is the same when it actually is not, the new 1063 preconditioner will not function correctly. Thus, use this optimization 1064 feature carefully! 1065 1066 If in doubt about whether your preconditioner matrix has changed 1067 structure or not, use the flag DIFFERENT_NONZERO_PATTERN. 1068 1069 More Notes about Repeated Solution of Linear Systems: 1070 PETSc does NOT reset the matrix entries of either Amat or Pmat 1071 to zero after a linear solve; the user is completely responsible for 1072 matrix assembly. See the routine MatZeroEntries() if desiring to 1073 zero all elements of a matrix. 1074 1075 Level: intermediate 1076 1077 .keywords: PC, set, operators, matrix, linear system 1078 1079 .seealso: PCGetOperators(), MatZeroEntries() 1080 @*/ 1081 int PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag) 1082 { 1083 int ierr; 1084 PetscTruth isbjacobi,isrowbs; 1085 1086 PetscFunctionBegin; 1087 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1088 PetscValidHeaderSpecific(Amat,MAT_COOKIE,2); 1089 PetscValidHeaderSpecific(Pmat,MAT_COOKIE,3); 1090 1091 /* 1092 BlockSolve95 cannot use default BJacobi preconditioning 1093 */ 1094 ierr = PetscTypeCompare((PetscObject)Amat,MATMPIROWBS,&isrowbs);CHKERRQ(ierr); 1095 if (isrowbs) { 1096 ierr = PetscTypeCompare((PetscObject)pc,PCBJACOBI,&isbjacobi);CHKERRQ(ierr); 1097 if (isbjacobi) { 1098 ierr = PCSetType(pc,PCILU);CHKERRQ(ierr); 1099 PetscLogInfo(pc,"PCSetOperators:Switching default PC to PCILU since BS95 doesn't support PCBJACOBI\n"); 1100 } 1101 } 1102 1103 pc->mat = Amat; 1104 pc->pmat = Pmat; 1105 if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) { 1106 pc->setupcalled = 1; 1107 } 1108 pc->flag = flag; 1109 PetscFunctionReturn(0); 1110 } 1111 1112 #undef __FUNCT__ 1113 #define __FUNCT__ "PCGetOperators" 1114 /*@C 1115 PCGetOperators - Gets the matrix associated with the linear system and 1116 possibly a different one associated with the preconditioner. 1117 1118 Not collective, though parallel Mats are returned if the PC is parallel 1119 1120 Input Parameter: 1121 . pc - the preconditioner context 1122 1123 Output Parameters: 1124 + mat - the matrix associated with the linear system 1125 . pmat - matrix associated with the preconditioner, usually the same 1126 as mat. 1127 - flag - flag indicating information about the preconditioner 1128 matrix structure. See PCSetOperators() for details. 1129 1130 Level: intermediate 1131 1132 .keywords: PC, get, operators, matrix, linear system 1133 1134 .seealso: PCSetOperators() 1135 @*/ 1136 int PCGetOperators(PC pc,Mat *mat,Mat *pmat,MatStructure *flag) 1137 { 1138 PetscFunctionBegin; 1139 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1140 if (mat) *mat = pc->mat; 1141 if (pmat) *pmat = pc->pmat; 1142 if (flag) *flag = pc->flag; 1143 PetscFunctionReturn(0); 1144 } 1145 1146 #undef __FUNCT__ 1147 #define __FUNCT__ "PCSetVector" 1148 /*@ 1149 PCSetVector - Sets a vector associated with the preconditioner. 1150 1151 Collective on PC and Vec 1152 1153 Input Parameters: 1154 + pc - the preconditioner context 1155 - vec - the vector 1156 1157 Notes: 1158 The vector must be set so that the preconditioner knows what type 1159 of vector to allocate if necessary. 1160 1161 Level: intermediate 1162 1163 .keywords: PC, set, vector 1164 1165 .seealso: PCGetVector() 1166 1167 @*/ 1168 int PCSetVector(PC pc,Vec vec) 1169 { 1170 PetscFunctionBegin; 1171 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1172 PetscValidHeaderSpecific(vec,VEC_COOKIE,2); 1173 PetscCheckSameComm(pc,vec); 1174 pc->vec = vec; 1175 PetscFunctionReturn(0); 1176 } 1177 1178 #undef __FUNCT__ 1179 #define __FUNCT__ "PCGetVector" 1180 /*@ 1181 PCGetVector - Gets a vector associated with the preconditioner; if the 1182 vector was not get set it will return a 0 pointer. 1183 1184 Not collective, but vector is shared by all processors that share the PC 1185 1186 Input Parameter: 1187 . pc - the preconditioner context 1188 1189 Output Parameter: 1190 . vec - the vector 1191 1192 Level: intermediate 1193 1194 .keywords: PC, get, vector 1195 1196 .seealso: PCSetVector() 1197 1198 @*/ 1199 int PCGetVector(PC pc,Vec *vec) 1200 { 1201 PetscFunctionBegin; 1202 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1203 *vec = pc->vec; 1204 PetscFunctionReturn(0); 1205 } 1206 1207 #undef __FUNCT__ 1208 #define __FUNCT__ "PCGetFactoredMatrix" 1209 /*@C 1210 PCGetFactoredMatrix - Gets the factored matrix from the 1211 preconditioner context. This routine is valid only for the LU, 1212 incomplete LU, Cholesky, and incomplete Cholesky methods. 1213 1214 Not Collective on PC though Mat is parallel if PC is parallel 1215 1216 Input Parameters: 1217 . pc - the preconditioner context 1218 1219 Output parameters: 1220 . mat - the factored matrix 1221 1222 Level: advanced 1223 1224 .keywords: PC, get, factored, matrix 1225 @*/ 1226 int PCGetFactoredMatrix(PC pc,Mat *mat) 1227 { 1228 int ierr; 1229 1230 PetscFunctionBegin; 1231 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1232 PetscValidPointer(mat,2); 1233 if (pc->ops->getfactoredmatrix) { 1234 ierr = (*pc->ops->getfactoredmatrix)(pc,mat);CHKERRQ(ierr); 1235 } 1236 PetscFunctionReturn(0); 1237 } 1238 1239 #undef __FUNCT__ 1240 #define __FUNCT__ "PCSetOptionsPrefix" 1241 /*@C 1242 PCSetOptionsPrefix - Sets the prefix used for searching for all 1243 PC options in the database. 1244 1245 Collective on PC 1246 1247 Input Parameters: 1248 + pc - the preconditioner context 1249 - prefix - the prefix string to prepend to all PC option requests 1250 1251 Notes: 1252 A hyphen (-) must NOT be given at the beginning of the prefix name. 1253 The first character of all runtime options is AUTOMATICALLY the 1254 hyphen. 1255 1256 Level: advanced 1257 1258 .keywords: PC, set, options, prefix, database 1259 1260 .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix() 1261 @*/ 1262 int PCSetOptionsPrefix(PC pc,const char prefix[]) 1263 { 1264 int ierr; 1265 1266 PetscFunctionBegin; 1267 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1268 ierr = PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr); 1269 PetscFunctionReturn(0); 1270 } 1271 1272 #undef __FUNCT__ 1273 #define __FUNCT__ "PCAppendOptionsPrefix" 1274 /*@C 1275 PCAppendOptionsPrefix - Appends to the prefix used for searching for all 1276 PC options in the database. 1277 1278 Collective on PC 1279 1280 Input Parameters: 1281 + pc - the preconditioner context 1282 - prefix - the prefix string to prepend to all PC option requests 1283 1284 Notes: 1285 A hyphen (-) must NOT be given at the beginning of the prefix name. 1286 The first character of all runtime options is AUTOMATICALLY the 1287 hyphen. 1288 1289 Level: advanced 1290 1291 .keywords: PC, append, options, prefix, database 1292 1293 .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix() 1294 @*/ 1295 int PCAppendOptionsPrefix(PC pc,const char prefix[]) 1296 { 1297 int ierr; 1298 1299 PetscFunctionBegin; 1300 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1301 ierr = PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr); 1302 PetscFunctionReturn(0); 1303 } 1304 1305 #undef __FUNCT__ 1306 #define __FUNCT__ "PCGetOptionsPrefix" 1307 /*@C 1308 PCGetOptionsPrefix - Gets the prefix used for searching for all 1309 PC options in the database. 1310 1311 Not Collective 1312 1313 Input Parameters: 1314 . pc - the preconditioner context 1315 1316 Output Parameters: 1317 . prefix - pointer to the prefix string used, is returned 1318 1319 Notes: On the fortran side, the user should pass in a string 'prifix' of 1320 sufficient length to hold the prefix. 1321 1322 Level: advanced 1323 1324 .keywords: PC, get, options, prefix, database 1325 1326 .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix() 1327 @*/ 1328 int PCGetOptionsPrefix(PC pc,char *prefix[]) 1329 { 1330 int ierr; 1331 1332 PetscFunctionBegin; 1333 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1334 PetscValidPointer(prefix,2); 1335 ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr); 1336 PetscFunctionReturn(0); 1337 } 1338 1339 #undef __FUNCT__ 1340 #define __FUNCT__ "PCPreSolve" 1341 /*@ 1342 PCPreSolve - Optional pre-solve phase, intended for any 1343 preconditioner-specific actions that must be performed before 1344 the iterative solve itself. 1345 1346 Collective on PC 1347 1348 Input Parameters: 1349 + pc - the preconditioner context 1350 - ksp - the Krylov subspace context 1351 1352 Level: developer 1353 1354 Sample of Usage: 1355 .vb 1356 PCPreSolve(pc,ksp); 1357 KSPSolve(ksp,its); 1358 PCPostSolve(pc,ksp); 1359 .ve 1360 1361 Notes: 1362 The pre-solve phase is distinct from the PCSetUp() phase. 1363 1364 KSPSolve() calls this directly, so is rarely called by the user. 1365 1366 .keywords: PC, pre-solve 1367 1368 .seealso: PCPostSolve() 1369 @*/ 1370 int PCPreSolve(PC pc,KSP ksp) 1371 { 1372 int ierr; 1373 Vec x,rhs; 1374 Mat A,B; 1375 1376 PetscFunctionBegin; 1377 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1378 PetscValidHeaderSpecific(ksp,KSP_COOKIE,2); 1379 ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr); 1380 ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr); 1381 /* 1382 Scale the system and have the matrices use the scaled form 1383 only if the two matrices are actually the same (and hence 1384 have the same scaling 1385 */ 1386 ierr = PCGetOperators(pc,&A,&B,PETSC_NULL);CHKERRQ(ierr); 1387 if (A == B) { 1388 ierr = MatScaleSystem(pc->mat,x,rhs);CHKERRQ(ierr); 1389 ierr = MatUseScaledForm(pc->mat,PETSC_TRUE);CHKERRQ(ierr); 1390 } 1391 1392 if (pc->ops->presolve) { 1393 ierr = (*pc->ops->presolve)(pc,ksp,x,rhs);CHKERRQ(ierr); 1394 } 1395 PetscFunctionReturn(0); 1396 } 1397 1398 #undef __FUNCT__ 1399 #define __FUNCT__ "PCPostSolve" 1400 /*@ 1401 PCPostSolve - Optional post-solve phase, intended for any 1402 preconditioner-specific actions that must be performed after 1403 the iterative solve itself. 1404 1405 Collective on PC 1406 1407 Input Parameters: 1408 + pc - the preconditioner context 1409 - ksp - the Krylov subspace context 1410 1411 Sample of Usage: 1412 .vb 1413 PCPreSolve(pc,ksp); 1414 KSPSolve(ksp,its); 1415 PCPostSolve(pc,ksp); 1416 .ve 1417 1418 Note: 1419 KSPSolve() calls this routine directly, so it is rarely called by the user. 1420 1421 Level: developer 1422 1423 .keywords: PC, post-solve 1424 1425 .seealso: PCPreSolve(), KSPSolve() 1426 @*/ 1427 int PCPostSolve(PC pc,KSP ksp) 1428 { 1429 int ierr; 1430 Vec x,rhs; 1431 Mat A,B; 1432 1433 PetscFunctionBegin; 1434 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1435 PetscValidHeaderSpecific(ksp,KSP_COOKIE,2); 1436 ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr); 1437 ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr); 1438 if (pc->ops->postsolve) { 1439 ierr = (*pc->ops->postsolve)(pc,ksp,x,rhs);CHKERRQ(ierr); 1440 } 1441 1442 /* 1443 Scale the system and have the matrices use the scaled form 1444 only if the two matrices are actually the same (and hence 1445 have the same scaling 1446 */ 1447 ierr = PCGetOperators(pc,&A,&B,PETSC_NULL);CHKERRQ(ierr); 1448 if (A == B) { 1449 ierr = MatUnScaleSystem(pc->mat,x,rhs);CHKERRQ(ierr); 1450 ierr = MatUseScaledForm(pc->mat,PETSC_FALSE);CHKERRQ(ierr); 1451 } 1452 PetscFunctionReturn(0); 1453 } 1454 1455 #undef __FUNCT__ 1456 #define __FUNCT__ "PCView" 1457 /*@C 1458 PCView - Prints the PC data structure. 1459 1460 Collective on PC 1461 1462 Input Parameters: 1463 + PC - the PC context 1464 - viewer - optional visualization context 1465 1466 Note: 1467 The available visualization contexts include 1468 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 1469 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 1470 output where only the first processor opens 1471 the file. All other processors send their 1472 data to the first processor to print. 1473 1474 The user can open an alternative visualization contexts with 1475 PetscViewerASCIIOpen() (output to a specified file). 1476 1477 Level: developer 1478 1479 .keywords: PC, view 1480 1481 .seealso: KSPView(), PetscViewerASCIIOpen() 1482 @*/ 1483 int PCView(PC pc,PetscViewer viewer) 1484 { 1485 PCType cstr; 1486 int ierr; 1487 PetscTruth mat_exists,isascii,isstring; 1488 PetscViewerFormat format; 1489 1490 PetscFunctionBegin; 1491 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1492 if (!viewer) viewer = PETSC_VIEWER_STDOUT_(pc->comm); 1493 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 1494 PetscCheckSameComm(pc,viewer); 1495 1496 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 1497 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 1498 if (isascii) { 1499 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1500 if (pc->prefix) { 1501 ierr = PetscViewerASCIIPrintf(viewer,"PC Object:(%s)\n",pc->prefix);CHKERRQ(ierr); 1502 } else { 1503 ierr = PetscViewerASCIIPrintf(viewer,"PC Object:\n");CHKERRQ(ierr); 1504 } 1505 ierr = PCGetType(pc,&cstr);CHKERRQ(ierr); 1506 if (cstr) { 1507 ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",cstr);CHKERRQ(ierr); 1508 } else { 1509 ierr = PetscViewerASCIIPrintf(viewer," type: not yet set\n");CHKERRQ(ierr); 1510 } 1511 if (pc->ops->view) { 1512 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1513 ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr); 1514 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1515 } 1516 ierr = PetscObjectExists((PetscObject)pc->mat,&mat_exists);CHKERRQ(ierr); 1517 if (mat_exists) { 1518 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 1519 if (pc->pmat == pc->mat) { 1520 ierr = PetscViewerASCIIPrintf(viewer," linear system matrix = precond matrix:\n");CHKERRQ(ierr); 1521 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1522 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 1523 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1524 } else { 1525 ierr = PetscObjectExists((PetscObject)pc->pmat,&mat_exists);CHKERRQ(ierr); 1526 if (mat_exists) { 1527 ierr = PetscViewerASCIIPrintf(viewer," linear system matrix followed by preconditioner matrix:\n");CHKERRQ(ierr); 1528 } else { 1529 ierr = PetscViewerASCIIPrintf(viewer," linear system matrix:\n");CHKERRQ(ierr); 1530 } 1531 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1532 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 1533 if (mat_exists) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);} 1534 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1535 } 1536 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1537 } 1538 } else if (isstring) { 1539 ierr = PCGetType(pc,&cstr);CHKERRQ(ierr); 1540 ierr = PetscViewerStringSPrintf(viewer," %-7.7s",cstr);CHKERRQ(ierr); 1541 if (pc->ops->view) {ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);} 1542 } else { 1543 SETERRQ1(1,"Viewer type %s not supported by PC",((PetscObject)viewer)->type_name); 1544 } 1545 PetscFunctionReturn(0); 1546 } 1547 1548 #undef __FUNCT__ 1549 #define __FUNCT__ "PCRegister" 1550 /*@C 1551 PCRegister - See PCRegisterDynamic() 1552 1553 Level: advanced 1554 @*/ 1555 int PCRegister(const char sname[],const char path[],const char name[],int (*function)(PC)) 1556 { 1557 int ierr; 1558 char fullname[256]; 1559 1560 PetscFunctionBegin; 1561 1562 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 1563 ierr = PetscFListAdd(&PCList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 1564 PetscFunctionReturn(0); 1565 } 1566 1567 #undef __FUNCT__ 1568 #define __FUNCT__ "PCComputeExplicitOperator" 1569 /*@ 1570 PCComputeExplicitOperator - Computes the explicit preconditioned operator. 1571 1572 Collective on PC 1573 1574 Input Parameter: 1575 . pc - the preconditioner object 1576 1577 Output Parameter: 1578 . mat - the explict preconditioned operator 1579 1580 Notes: 1581 This computation is done by applying the operators to columns of the 1582 identity matrix. 1583 1584 Currently, this routine uses a dense matrix format when 1 processor 1585 is used and a sparse format otherwise. This routine is costly in general, 1586 and is recommended for use only with relatively small systems. 1587 1588 Level: advanced 1589 1590 .keywords: PC, compute, explicit, operator 1591 1592 @*/ 1593 int PCComputeExplicitOperator(PC pc,Mat *mat) 1594 { 1595 Vec in,out; 1596 int ierr,i,M,m,size,*rows,start,end; 1597 MPI_Comm comm; 1598 PetscScalar *array,zero = 0.0,one = 1.0; 1599 1600 PetscFunctionBegin; 1601 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1602 PetscValidPointer(mat,2); 1603 1604 comm = pc->comm; 1605 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1606 1607 ierr = PCGetVector(pc,&in);CHKERRQ(ierr); 1608 ierr = VecDuplicate(in,&out);CHKERRQ(ierr); 1609 ierr = VecGetOwnershipRange(in,&start,&end);CHKERRQ(ierr); 1610 ierr = VecGetSize(in,&M);CHKERRQ(ierr); 1611 ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr); 1612 ierr = PetscMalloc((m+1)*sizeof(int),&rows);CHKERRQ(ierr); 1613 for (i=0; i<m; i++) {rows[i] = start + i;} 1614 1615 if (size == 1) { 1616 ierr = MatCreateSeqDense(comm,M,M,PETSC_NULL,mat);CHKERRQ(ierr); 1617 } else { 1618 ierr = MatCreateMPIAIJ(comm,m,m,M,M,0,0,0,0,mat);CHKERRQ(ierr); 1619 } 1620 1621 for (i=0; i<M; i++) { 1622 1623 ierr = VecSet(&zero,in);CHKERRQ(ierr); 1624 ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 1625 ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 1626 ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 1627 1628 /* should fix, allowing user to choose side */ 1629 ierr = PCApply(pc,in,out,PC_LEFT);CHKERRQ(ierr); 1630 1631 ierr = VecGetArray(out,&array);CHKERRQ(ierr); 1632 ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 1633 ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 1634 1635 } 1636 ierr = PetscFree(rows);CHKERRQ(ierr); 1637 ierr = VecDestroy(out);CHKERRQ(ierr); 1638 ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1639 ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1640 PetscFunctionReturn(0); 1641 } 1642 1643