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