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