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