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