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