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