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