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