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