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