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