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