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