14b9ad928SBarry Smith /* 24b9ad928SBarry Smith The PC (preconditioner) interface routines, callable by users. 34b9ad928SBarry Smith */ 4af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscksp.h" I*/ 51e25c274SJed Brown #include <petscdm.h> 64b9ad928SBarry Smith 74b9ad928SBarry Smith /* Logging support */ 80700a824SBarry Smith PetscClassId PC_CLASSID; 9c677e75fSPierre Jolivet PetscLogEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_MatApply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft; 10011ea8aeSBarry Smith PetscLogEvent PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks; 11ab83eea4SMatthew G. Knepley PetscInt PetscMGLevelId; 120316ec64SBarry Smith PetscLogStage PCMPIStage; 134b9ad928SBarry Smith 14b4f8a55aSStefano Zampini PETSC_INTERN PetscErrorCode PCGetDefaultType_Private(PC pc, const char *type[]) 15d71ae5a4SJacob Faibussowitsch { 162e70eadcSBarry Smith PetscMPIInt size; 17b4f8a55aSStefano Zampini PetscBool hasopblock, hasopsolve, flg1, flg2, set, flg3, isnormal; 18566e8bf2SBarry Smith 19566e8bf2SBarry Smith PetscFunctionBegin; 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 21566e8bf2SBarry Smith if (pc->pmat) { 22b4f8a55aSStefano Zampini PetscCall(MatHasOperation(pc->pmat, MATOP_GET_DIAGONAL_BLOCK, &hasopblock)); 23b4f8a55aSStefano Zampini PetscCall(MatHasOperation(pc->pmat, MATOP_SOLVE, &hasopsolve)); 24566e8bf2SBarry Smith if (size == 1) { 259566063dSJacob Faibussowitsch PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ICC, &flg1)); 269566063dSJacob Faibussowitsch PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ILU, &flg2)); 279566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(pc->pmat, &set, &flg3)); 28da74c943SJose E. Roman PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isnormal, MATNORMAL, MATNORMALHERMITIAN, NULL)); 29ba8a7149SBarry Smith if (flg1 && (!flg2 || (set && flg3))) { 30566e8bf2SBarry Smith *type = PCICC; 31566e8bf2SBarry Smith } else if (flg2) { 32566e8bf2SBarry Smith *type = PCILU; 33da74c943SJose E. Roman } else if (isnormal) { 34da74c943SJose E. Roman *type = PCNONE; 35b4f8a55aSStefano Zampini } else if (hasopblock) { /* likely is a parallel matrix run on one processor */ 36566e8bf2SBarry Smith *type = PCBJACOBI; 37b4f8a55aSStefano Zampini } else if (hasopsolve) { 38b4f8a55aSStefano Zampini *type = PCMAT; 39566e8bf2SBarry Smith } else { 40566e8bf2SBarry Smith *type = PCNONE; 41566e8bf2SBarry Smith } 42566e8bf2SBarry Smith } else { 43b4f8a55aSStefano Zampini if (hasopblock) { 44566e8bf2SBarry Smith *type = PCBJACOBI; 45b4f8a55aSStefano Zampini } else if (hasopsolve) { 46b4f8a55aSStefano Zampini *type = PCMAT; 47566e8bf2SBarry Smith } else { 48566e8bf2SBarry Smith *type = PCNONE; 49566e8bf2SBarry Smith } 50566e8bf2SBarry Smith } 51b4f8a55aSStefano Zampini } else *type = NULL; 523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53566e8bf2SBarry Smith } 544b9ad928SBarry Smith 55*fe57bb1aSStefano Zampini /* do not log solves, setup and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */ 56*fe57bb1aSStefano Zampini PETSC_EXTERN PetscLogEvent KSP_Solve, KSP_SetUp; 57*fe57bb1aSStefano Zampini 58*fe57bb1aSStefano Zampini static PetscErrorCode PCLogEventsDeactivatePush(void) 59*fe57bb1aSStefano Zampini { 60*fe57bb1aSStefano Zampini PetscFunctionBegin; 61*fe57bb1aSStefano Zampini PetscCall(KSPInitializePackage()); 62*fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePush(KSP_Solve)); 63*fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePush(KSP_SetUp)); 64*fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePush(PC_Apply)); 65*fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePush(PC_SetUp)); 66*fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePush(PC_SetUpOnBlocks)); 67*fe57bb1aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 68*fe57bb1aSStefano Zampini } 69*fe57bb1aSStefano Zampini 70*fe57bb1aSStefano Zampini static PetscErrorCode PCLogEventsDeactivatePop(void) 71*fe57bb1aSStefano Zampini { 72*fe57bb1aSStefano Zampini PetscFunctionBegin; 73*fe57bb1aSStefano Zampini PetscCall(KSPInitializePackage()); 74*fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePop(KSP_Solve)); 75*fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePop(KSP_SetUp)); 76*fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePop(PC_Apply)); 77*fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePop(PC_SetUp)); 78*fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePop(PC_SetUpOnBlocks)); 79*fe57bb1aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 80*fe57bb1aSStefano Zampini } 81*fe57bb1aSStefano Zampini 8288584be7SBarry Smith /*@ 83562efe2eSBarry Smith PCReset - Resets a `PC` context to the pcsetupcalled = 0 state and removes any allocated `Vec`s and `Mat`s 8488584be7SBarry Smith 85c3339decSBarry Smith Collective 8688584be7SBarry Smith 8788584be7SBarry Smith Input Parameter: 8888584be7SBarry Smith . pc - the preconditioner context 8988584be7SBarry Smith 9088584be7SBarry Smith Level: developer 9188584be7SBarry Smith 92f1580f4eSBarry Smith Note: 93562efe2eSBarry Smith This allows a `PC` to be reused for a different sized linear system but using the same options that have been previously set in `pc` 9488584be7SBarry Smith 95562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()` 9688584be7SBarry Smith @*/ 97d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset(PC pc) 98d71ae5a4SJacob Faibussowitsch { 9988584be7SBarry Smith PetscFunctionBegin; 10088584be7SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 101dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, reset); 1029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc->diagonalscaleright)); 1039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc->diagonalscaleleft)); 1049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->pmat)); 1059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->mat)); 1062fa5cd67SKarl Rupp 1070ce6c5a2SBarry Smith pc->setupcalled = 0; 1083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10988584be7SBarry Smith } 11088584be7SBarry Smith 1110764c050SBarry Smith /*@ 112f1580f4eSBarry Smith PCDestroy - Destroys `PC` context that was created with `PCCreate()`. 1134b9ad928SBarry Smith 114c3339decSBarry Smith Collective 1154b9ad928SBarry Smith 1164b9ad928SBarry Smith Input Parameter: 1174b9ad928SBarry Smith . pc - the preconditioner context 1184b9ad928SBarry Smith 1194b9ad928SBarry Smith Level: developer 1204b9ad928SBarry Smith 121562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()` 1224b9ad928SBarry Smith @*/ 123d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy(PC *pc) 124d71ae5a4SJacob Faibussowitsch { 1254b9ad928SBarry Smith PetscFunctionBegin; 1263ba16761SJacob Faibussowitsch if (!*pc) PetscFunctionReturn(PETSC_SUCCESS); 127f4f49eeaSPierre Jolivet PetscValidHeaderSpecific(*pc, PC_CLASSID, 1); 128f4f49eeaSPierre Jolivet if (--((PetscObject)*pc)->refct > 0) { 1299371c9d4SSatish Balay *pc = NULL; 1303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1319371c9d4SSatish Balay } 1324b9ad928SBarry Smith 1339566063dSJacob Faibussowitsch PetscCall(PCReset(*pc)); 134241cb3abSLisandro Dalcin 135e04113cfSBarry Smith /* if memory was published with SAWs then destroy it */ 1369566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsViewOff((PetscObject)*pc)); 137f4f49eeaSPierre Jolivet PetscTryTypeMethod(*pc, destroy); 1389566063dSJacob Faibussowitsch PetscCall(DMDestroy(&(*pc)->dm)); 1399566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(pc)); 1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1414b9ad928SBarry Smith } 1424b9ad928SBarry Smith 143cc4c1da9SBarry Smith /*@ 144c5eb9154SBarry Smith PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right 1454b9ad928SBarry Smith scaling as needed by certain time-stepping codes. 1464b9ad928SBarry Smith 147c3339decSBarry Smith Logically Collective 1484b9ad928SBarry Smith 1494b9ad928SBarry Smith Input Parameter: 1504b9ad928SBarry Smith . pc - the preconditioner context 1514b9ad928SBarry Smith 1524b9ad928SBarry Smith Output Parameter: 153f1580f4eSBarry Smith . flag - `PETSC_TRUE` if it applies the scaling 1544b9ad928SBarry Smith 1554b9ad928SBarry Smith Level: developer 1564b9ad928SBarry Smith 157f1580f4eSBarry Smith Note: 158562efe2eSBarry Smith If this returns `PETSC_TRUE` then the system solved via the Krylov method is, for left and right preconditioning, 1594b9ad928SBarry Smith 160562efe2eSBarry Smith $$ 161562efe2eSBarry Smith \begin{align*} 162562efe2eSBarry Smith D M A D^{-1} y = D M b \\ 163562efe2eSBarry Smith D A M D^{-1} z = D b. 164562efe2eSBarry Smith \end{align*} 165562efe2eSBarry Smith $$ 166562efe2eSBarry Smith 167562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCSetDiagonalScale()` 1684b9ad928SBarry Smith @*/ 169d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetDiagonalScale(PC pc, PetscBool *flag) 170d71ae5a4SJacob Faibussowitsch { 1714b9ad928SBarry Smith PetscFunctionBegin; 1720700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1734f572ea9SToby Isaac PetscAssertPointer(flag, 2); 1744b9ad928SBarry Smith *flag = pc->diagonalscale; 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1764b9ad928SBarry Smith } 1774b9ad928SBarry Smith 1784b9ad928SBarry Smith /*@ 179c5eb9154SBarry Smith PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right 1804b9ad928SBarry Smith scaling as needed by certain time-stepping codes. 1814b9ad928SBarry Smith 182c3339decSBarry Smith Logically Collective 1834b9ad928SBarry Smith 1844b9ad928SBarry Smith Input Parameters: 1854b9ad928SBarry Smith + pc - the preconditioner context 1864b9ad928SBarry Smith - s - scaling vector 1874b9ad928SBarry Smith 1884b9ad928SBarry Smith Level: intermediate 1894b9ad928SBarry Smith 19095452b02SPatrick Sanan Notes: 191562efe2eSBarry Smith The system solved via the Krylov method is, for left and right preconditioning, 192562efe2eSBarry Smith $$ 193562efe2eSBarry Smith \begin{align*} 194562efe2eSBarry Smith D M A D^{-1} y = D M b \\ 195562efe2eSBarry Smith D A M D^{-1} z = D b. 196562efe2eSBarry Smith \end{align*} 197562efe2eSBarry Smith $$ 1984b9ad928SBarry Smith 199562efe2eSBarry Smith `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$. 2004b9ad928SBarry Smith 201562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCGetDiagonalScale()` 2024b9ad928SBarry Smith @*/ 203d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetDiagonalScale(PC pc, Vec s) 204d71ae5a4SJacob Faibussowitsch { 2054b9ad928SBarry Smith PetscFunctionBegin; 2060700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2070700a824SBarry Smith PetscValidHeaderSpecific(s, VEC_CLASSID, 2); 2084b9ad928SBarry Smith pc->diagonalscale = PETSC_TRUE; 2092fa5cd67SKarl Rupp 2109566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)s)); 2119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc->diagonalscaleleft)); 2122fa5cd67SKarl Rupp 2134b9ad928SBarry Smith pc->diagonalscaleleft = s; 2142fa5cd67SKarl Rupp 2159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s, &pc->diagonalscaleright)); 2169566063dSJacob Faibussowitsch PetscCall(VecCopy(s, pc->diagonalscaleright)); 2179566063dSJacob Faibussowitsch PetscCall(VecReciprocal(pc->diagonalscaleright)); 2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2194b9ad928SBarry Smith } 2204b9ad928SBarry Smith 221bc08b0f1SBarry Smith /*@ 222c5eb9154SBarry Smith PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes. 2234b9ad928SBarry Smith 224c3339decSBarry Smith Logically Collective 2254b9ad928SBarry Smith 2264b9ad928SBarry Smith Input Parameters: 2274b9ad928SBarry Smith + pc - the preconditioner context 2284b9ad928SBarry Smith . in - input vector 229a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in) 2304b9ad928SBarry Smith 2314b9ad928SBarry Smith Level: intermediate 2324b9ad928SBarry Smith 23395452b02SPatrick Sanan Notes: 234562efe2eSBarry Smith The system solved via the Krylov method is, for left and right preconditioning, 2354b9ad928SBarry Smith 236562efe2eSBarry Smith $$ 237562efe2eSBarry Smith \begin{align*} 238562efe2eSBarry Smith D M A D^{-1} y = D M b \\ 239562efe2eSBarry Smith D A M D^{-1} z = D b. 240562efe2eSBarry Smith \end{align*} 241562efe2eSBarry Smith $$ 2424b9ad928SBarry Smith 243562efe2eSBarry Smith `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$. 2444b9ad928SBarry Smith 245562efe2eSBarry Smith If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out` 246562efe2eSBarry Smith 2470241b274SPierre Jolivet .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCSetDiagonalScale()`, `PCDiagonalScaleRight()`, `MatDiagonalScale()` 2484b9ad928SBarry Smith @*/ 249d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleLeft(PC pc, Vec in, Vec out) 250d71ae5a4SJacob Faibussowitsch { 2514b9ad928SBarry Smith PetscFunctionBegin; 2520700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2530700a824SBarry Smith PetscValidHeaderSpecific(in, VEC_CLASSID, 2); 2540700a824SBarry Smith PetscValidHeaderSpecific(out, VEC_CLASSID, 3); 2554b9ad928SBarry Smith if (pc->diagonalscale) { 2569566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(out, pc->diagonalscaleleft, in)); 2574b9ad928SBarry Smith } else if (in != out) { 2589566063dSJacob Faibussowitsch PetscCall(VecCopy(in, out)); 2594b9ad928SBarry Smith } 2603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2614b9ad928SBarry Smith } 2624b9ad928SBarry Smith 263bc08b0f1SBarry Smith /*@ 2644b9ad928SBarry Smith PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes. 2654b9ad928SBarry Smith 266c3339decSBarry Smith Logically Collective 2674b9ad928SBarry Smith 2684b9ad928SBarry Smith Input Parameters: 2694b9ad928SBarry Smith + pc - the preconditioner context 2704b9ad928SBarry Smith . in - input vector 271a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in) 2724b9ad928SBarry Smith 2734b9ad928SBarry Smith Level: intermediate 2744b9ad928SBarry Smith 27595452b02SPatrick Sanan Notes: 276562efe2eSBarry Smith The system solved via the Krylov method is, for left and right preconditioning, 2774b9ad928SBarry Smith 278562efe2eSBarry Smith $$ 279562efe2eSBarry Smith \begin{align*} 280562efe2eSBarry Smith D M A D^{-1} y = D M b \\ 281562efe2eSBarry Smith D A M D^{-1} z = D b. 282a3524536SPierre Jolivet \end{align*} 283562efe2eSBarry Smith $$ 2844b9ad928SBarry Smith 285562efe2eSBarry Smith `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$. 2864b9ad928SBarry Smith 287562efe2eSBarry Smith If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out` 288562efe2eSBarry Smith 2890241b274SPierre Jolivet .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCSetDiagonalScale()`, `MatDiagonalScale()` 2904b9ad928SBarry Smith @*/ 291d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleRight(PC pc, Vec in, Vec out) 292d71ae5a4SJacob Faibussowitsch { 2934b9ad928SBarry Smith PetscFunctionBegin; 2940700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2950700a824SBarry Smith PetscValidHeaderSpecific(in, VEC_CLASSID, 2); 2960700a824SBarry Smith PetscValidHeaderSpecific(out, VEC_CLASSID, 3); 2974b9ad928SBarry Smith if (pc->diagonalscale) { 2989566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(out, pc->diagonalscaleright, in)); 2994b9ad928SBarry Smith } else if (in != out) { 3009566063dSJacob Faibussowitsch PetscCall(VecCopy(in, out)); 3014b9ad928SBarry Smith } 3023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3034b9ad928SBarry Smith } 3044b9ad928SBarry Smith 30549517cdeSBarry Smith /*@ 30649517cdeSBarry Smith PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the 307f1580f4eSBarry Smith operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`, 308f1580f4eSBarry Smith `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat. 30949517cdeSBarry Smith 310c3339decSBarry Smith Logically Collective 31149517cdeSBarry Smith 31249517cdeSBarry Smith Input Parameters: 31349517cdeSBarry Smith + pc - the preconditioner context 314f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false) 31549517cdeSBarry Smith 31649517cdeSBarry Smith Options Database Key: 317562efe2eSBarry Smith . -pc_use_amat <true,false> - use the amat argument to `KSPSetOperators()` or `PCSetOperators()` to apply the operator 31849517cdeSBarry Smith 31920f4b53cSBarry Smith Level: intermediate 32020f4b53cSBarry Smith 321f1580f4eSBarry Smith Note: 32249517cdeSBarry Smith For the common case in which the linear system matrix and the matrix used to construct the 323562efe2eSBarry Smith preconditioner are identical, this routine has no affect. 32449517cdeSBarry Smith 3250241b274SPierre Jolivet .seealso: [](ch_ksp), `PC`, `PCGetUseAmat()`, `PCBJACOBI`, `PCMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`, 326562efe2eSBarry Smith `KSPSetOperators()`, `PCSetOperators()` 32749517cdeSBarry Smith @*/ 328d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUseAmat(PC pc, PetscBool flg) 329d71ae5a4SJacob Faibussowitsch { 33049517cdeSBarry Smith PetscFunctionBegin; 33149517cdeSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 33249517cdeSBarry Smith pc->useAmat = flg; 3333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33449517cdeSBarry Smith } 33549517cdeSBarry Smith 336422a814eSBarry Smith /*@ 337562efe2eSBarry Smith PCSetErrorIfFailure - Causes `PC` to generate an error if a floating point exception, for example a zero pivot, is detected. 338422a814eSBarry Smith 339c3339decSBarry Smith Logically Collective 340422a814eSBarry Smith 341422a814eSBarry Smith Input Parameters: 342562efe2eSBarry Smith + pc - iterative context obtained from `PCCreate()` 343f1580f4eSBarry Smith - flg - `PETSC_TRUE` indicates you want the error generated 344422a814eSBarry Smith 345422a814eSBarry Smith Level: advanced 346422a814eSBarry Smith 347422a814eSBarry Smith Notes: 348f1580f4eSBarry Smith Normally PETSc continues if a linear solver fails due to a failed setup of a preconditioner, you can call `KSPGetConvergedReason()` after a `KSPSolve()` 349422a814eSBarry Smith 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 350422a814eSBarry Smith detected. 351422a814eSBarry Smith 352562efe2eSBarry Smith This is propagated into `KSP`s used by this `PC`, which then propagate it into `PC`s used by those `KSP`s 353422a814eSBarry Smith 354562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPSetErrorIfNotConverged()`, `PCGetInitialGuessNonzero()`, `PCSetInitialGuessKnoll()`, `PCGetInitialGuessKnoll()` 355422a814eSBarry Smith @*/ 356d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetErrorIfFailure(PC pc, PetscBool flg) 357d71ae5a4SJacob Faibussowitsch { 358422a814eSBarry Smith PetscFunctionBegin; 359422a814eSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 360422a814eSBarry Smith PetscValidLogicalCollectiveBool(pc, flg, 2); 361422a814eSBarry Smith pc->erroriffailure = flg; 3623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 363422a814eSBarry Smith } 364422a814eSBarry Smith 36549517cdeSBarry Smith /*@ 36649517cdeSBarry Smith PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the 367f1580f4eSBarry Smith operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`, 368f1580f4eSBarry Smith `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat. 36949517cdeSBarry Smith 370c3339decSBarry Smith Logically Collective 37149517cdeSBarry Smith 37249517cdeSBarry Smith Input Parameter: 37349517cdeSBarry Smith . pc - the preconditioner context 37449517cdeSBarry Smith 37549517cdeSBarry Smith Output Parameter: 376f1580f4eSBarry Smith . flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false) 37749517cdeSBarry Smith 37820f4b53cSBarry Smith Level: intermediate 37920f4b53cSBarry Smith 380f1580f4eSBarry Smith Note: 38149517cdeSBarry Smith For the common case in which the linear system matrix and the matrix used to construct the 38249517cdeSBarry Smith preconditioner are identical, this routine is does nothing. 38349517cdeSBarry Smith 3840241b274SPierre Jolivet .seealso: [](ch_ksp), `PC`, `PCSetUseAmat()`, `PCBJACOBI`, `PCMG`, `PCFIELDSPLIT`, `PCCOMPOSITE` 38549517cdeSBarry Smith @*/ 386d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetUseAmat(PC pc, PetscBool *flg) 387d71ae5a4SJacob Faibussowitsch { 38849517cdeSBarry Smith PetscFunctionBegin; 38949517cdeSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 39049517cdeSBarry Smith *flg = pc->useAmat; 3913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39249517cdeSBarry Smith } 39349517cdeSBarry Smith 394f39d8e23SSatish Balay /*@ 3953821be0aSBarry Smith PCSetKSPNestLevel - sets the amount of nesting the `KSP` that contains this `PC` has 3963821be0aSBarry Smith 3973821be0aSBarry Smith Collective 3983821be0aSBarry Smith 3993821be0aSBarry Smith Input Parameters: 4003821be0aSBarry Smith + pc - the `PC` 4013821be0aSBarry Smith - level - the nest level 4023821be0aSBarry Smith 4033821be0aSBarry Smith Level: developer 4043821be0aSBarry Smith 4053821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCGetKSPNestLevel()`, `KSPSetNestLevel()` 4063821be0aSBarry Smith @*/ 4073821be0aSBarry Smith PetscErrorCode PCSetKSPNestLevel(PC pc, PetscInt level) 4083821be0aSBarry Smith { 4093821be0aSBarry Smith PetscFunctionBegin; 4107a99bfcaSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 4117a99bfcaSBarry Smith PetscValidLogicalCollectiveInt(pc, level, 2); 4123821be0aSBarry Smith pc->kspnestlevel = level; 4133821be0aSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 4143821be0aSBarry Smith } 4153821be0aSBarry Smith 4163821be0aSBarry Smith /*@ 4173821be0aSBarry Smith PCGetKSPNestLevel - gets the amount of nesting the `KSP` that contains this `PC` has 4183821be0aSBarry Smith 4197a99bfcaSBarry Smith Not Collective 4203821be0aSBarry Smith 4213821be0aSBarry Smith Input Parameter: 4223821be0aSBarry Smith . pc - the `PC` 4233821be0aSBarry Smith 4243821be0aSBarry Smith Output Parameter: 4253821be0aSBarry Smith . level - the nest level 4263821be0aSBarry Smith 4273821be0aSBarry Smith Level: developer 4283821be0aSBarry Smith 4293821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `KSPGetNestLevel()` 4303821be0aSBarry Smith @*/ 4313821be0aSBarry Smith PetscErrorCode PCGetKSPNestLevel(PC pc, PetscInt *level) 4323821be0aSBarry Smith { 4333821be0aSBarry Smith PetscFunctionBegin; 4347a99bfcaSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 4357a99bfcaSBarry Smith PetscAssertPointer(level, 2); 4363821be0aSBarry Smith *level = pc->kspnestlevel; 4373821be0aSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 4383821be0aSBarry Smith } 4393821be0aSBarry Smith 4403821be0aSBarry Smith /*@ 441f1580f4eSBarry Smith PCCreate - Creates a preconditioner context, `PC` 4424b9ad928SBarry Smith 443d083f849SBarry Smith Collective 4444b9ad928SBarry Smith 4454b9ad928SBarry Smith Input Parameter: 4464b9ad928SBarry Smith . comm - MPI communicator 4474b9ad928SBarry Smith 4484b9ad928SBarry Smith Output Parameter: 4497a99bfcaSBarry Smith . newpc - location to put the preconditioner context 4504b9ad928SBarry Smith 45120f4b53cSBarry Smith Level: developer 45220f4b53cSBarry Smith 453f1580f4eSBarry Smith Note: 454f1580f4eSBarry Smith The default preconditioner for sparse matrices is `PCILU` or `PCICC` with 0 fill on one process and block Jacobi (`PCBJACOBI`) with `PCILU` or `PCICC` 455f1580f4eSBarry Smith in parallel. For dense matrices it is always `PCNONE`. 4564b9ad928SBarry Smith 457562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCApply()`, `PCDestroy()` 4584b9ad928SBarry Smith @*/ 459d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreate(MPI_Comm comm, PC *newpc) 460d71ae5a4SJacob Faibussowitsch { 4614b9ad928SBarry Smith PC pc; 4624b9ad928SBarry Smith 4634b9ad928SBarry Smith PetscFunctionBegin; 4644f572ea9SToby Isaac PetscAssertPointer(newpc, 2); 4659566063dSJacob Faibussowitsch PetscCall(PCInitializePackage()); 4664b9ad928SBarry Smith 4679566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(pc, PC_CLASSID, "PC", "Preconditioner", "PC", comm, PCDestroy, PCView)); 4680a545947SLisandro Dalcin pc->mat = NULL; 4690a545947SLisandro Dalcin pc->pmat = NULL; 4704b9ad928SBarry Smith pc->setupcalled = 0; 4710c24e6a1SHong Zhang pc->setfromoptionscalled = 0; 4720a545947SLisandro Dalcin pc->data = NULL; 4734b9ad928SBarry Smith pc->diagonalscale = PETSC_FALSE; 4740a545947SLisandro Dalcin pc->diagonalscaleleft = NULL; 4750a545947SLisandro Dalcin pc->diagonalscaleright = NULL; 4764b9ad928SBarry Smith 4770a545947SLisandro Dalcin pc->modifysubmatrices = NULL; 4780a545947SLisandro Dalcin pc->modifysubmatricesP = NULL; 4792fa5cd67SKarl Rupp 480a7d21a52SLisandro Dalcin *newpc = pc; 4813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4824b9ad928SBarry Smith } 4834b9ad928SBarry Smith 4844b9ad928SBarry Smith /*@ 4854b9ad928SBarry Smith PCApply - Applies the preconditioner to a vector. 4864b9ad928SBarry Smith 487c3339decSBarry Smith Collective 4884b9ad928SBarry Smith 4894b9ad928SBarry Smith Input Parameters: 4904b9ad928SBarry Smith + pc - the preconditioner context 4914b9ad928SBarry Smith - x - input vector 4924b9ad928SBarry Smith 4934b9ad928SBarry Smith Output Parameter: 4944b9ad928SBarry Smith . y - output vector 4954b9ad928SBarry Smith 4964b9ad928SBarry Smith Level: developer 4974b9ad928SBarry Smith 498562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`, `PCApplyBAorAB()` 4994b9ad928SBarry Smith @*/ 500d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApply(PC pc, Vec x, Vec y) 501d71ae5a4SJacob Faibussowitsch { 502106e20bfSBarry Smith PetscInt m, n, mv, nv; 5034b9ad928SBarry Smith 5044b9ad928SBarry Smith PetscFunctionBegin; 5050700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 5060700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 5070700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 5087827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 509e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 510540bdfbdSVaclav Hapla /* use pmat to check vector sizes since for KSPLSQR the pmat may be of a different size than mat */ 5119566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &m, &n)); 5129566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &mv)); 5139566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &nv)); 514540bdfbdSVaclav Hapla /* check pmat * y = x is feasible */ 51563a3b9bcSJacob Faibussowitsch PetscCheck(mv == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Preconditioner number of local rows %" PetscInt_FMT " does not equal input vector size %" PetscInt_FMT, m, mv); 51663a3b9bcSJacob Faibussowitsch PetscCheck(nv == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Preconditioner number of local columns %" PetscInt_FMT " does not equal output vector size %" PetscInt_FMT, n, nv); 5179566063dSJacob Faibussowitsch PetscCall(VecSetErrorIfLocked(y, 3)); 518106e20bfSBarry Smith 5199566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 5209566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(x)); 5219566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0)); 522dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, apply, x, y); 5239566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0)); 524e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 5259566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(x)); 5263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5274b9ad928SBarry Smith } 5284b9ad928SBarry Smith 5294b9ad928SBarry Smith /*@ 530562efe2eSBarry Smith PCMatApply - Applies the preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApply()`, `Y` and `X` must be different matrices. 531c41ea70eSPierre Jolivet 532c3339decSBarry Smith Collective 533c677e75fSPierre Jolivet 534c677e75fSPierre Jolivet Input Parameters: 535c677e75fSPierre Jolivet + pc - the preconditioner context 536c677e75fSPierre Jolivet - X - block of input vectors 537c677e75fSPierre Jolivet 538c677e75fSPierre Jolivet Output Parameter: 539c677e75fSPierre Jolivet . Y - block of output vectors 540c677e75fSPierre Jolivet 541c677e75fSPierre Jolivet Level: developer 542c677e75fSPierre Jolivet 543562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `KSPMatSolve()` 544c677e75fSPierre Jolivet @*/ 545d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMatApply(PC pc, Mat X, Mat Y) 546d71ae5a4SJacob Faibussowitsch { 547c677e75fSPierre Jolivet Mat A; 548c677e75fSPierre Jolivet Vec cy, cx; 549bd82155bSPierre Jolivet PetscInt m1, M1, m2, M2, n1, N1, n2, N2, m3, M3, n3, N3; 550c677e75fSPierre Jolivet PetscBool match; 551c677e75fSPierre Jolivet 552c677e75fSPierre Jolivet PetscFunctionBegin; 553c677e75fSPierre Jolivet PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 554e2b200f6SPierre Jolivet PetscValidHeaderSpecific(X, MAT_CLASSID, 2); 555e2b200f6SPierre Jolivet PetscValidHeaderSpecific(Y, MAT_CLASSID, 3); 556e2b200f6SPierre Jolivet PetscCheckSameComm(pc, 1, X, 2); 557e2b200f6SPierre Jolivet PetscCheckSameComm(pc, 1, Y, 3); 5587827d75bSBarry Smith PetscCheck(Y != X, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "Y and X must be different matrices"); 5599566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, NULL, &A)); 5609566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m3, &n3)); 5619566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(X, &m2, &n2)); 5629566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Y, &m1, &n1)); 5639566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M3, &N3)); 5649566063dSJacob Faibussowitsch PetscCall(MatGetSize(X, &M2, &N2)); 5659566063dSJacob Faibussowitsch PetscCall(MatGetSize(Y, &M1, &N1)); 56663a3b9bcSJacob Faibussowitsch PetscCheck(n1 == n2 && N1 == N2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible number of columns between block of input vectors (n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and block of output vectors (n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")", n2, N2, n1, N1); 56763a3b9bcSJacob Faibussowitsch PetscCheck(m2 == m3 && M2 == M3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout between block of input vectors (m,M) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and Pmat (m,M)x(n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")x(%" PetscInt_FMT ",%" PetscInt_FMT ")", m2, M2, m3, M3, n3, N3); 56863a3b9bcSJacob Faibussowitsch PetscCheck(m1 == n3 && M1 == N3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout between block of output vectors (m,M) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and Pmat (m,M)x(n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")x(%" PetscInt_FMT ",%" PetscInt_FMT ")", m1, M1, m3, M3, n3, N3); 5699566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &match, MATSEQDENSE, MATMPIDENSE, "")); 57028b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of output vectors not stored in a dense Mat"); 5719566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, "")); 57228b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of input vectors not stored in a dense Mat"); 5739566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 574c677e75fSPierre Jolivet if (pc->ops->matapply) { 5759566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0)); 576dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, matapply, X, Y); 5779566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0)); 578c677e75fSPierre Jolivet } else { 5799566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "PC type %s applying column by column\n", ((PetscObject)pc)->type_name)); 580bd82155bSPierre Jolivet for (n1 = 0; n1 < N1; ++n1) { 5819566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(X, n1, &cx)); 5829566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Y, n1, &cy)); 5839566063dSJacob Faibussowitsch PetscCall(PCApply(pc, cx, cy)); 5849566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Y, n1, &cy)); 5859566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(X, n1, &cx)); 586c677e75fSPierre Jolivet } 587c677e75fSPierre Jolivet } 5883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 589c677e75fSPierre Jolivet } 590c677e75fSPierre Jolivet 591c677e75fSPierre Jolivet /*@ 5924b9ad928SBarry Smith PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector. 5934b9ad928SBarry Smith 594c3339decSBarry Smith Collective 5954b9ad928SBarry Smith 5964b9ad928SBarry Smith Input Parameters: 5974b9ad928SBarry Smith + pc - the preconditioner context 5984b9ad928SBarry Smith - x - input vector 5994b9ad928SBarry Smith 6004b9ad928SBarry Smith Output Parameter: 6014b9ad928SBarry Smith . y - output vector 6024b9ad928SBarry Smith 60320f4b53cSBarry Smith Level: developer 60420f4b53cSBarry Smith 605f1580f4eSBarry Smith Note: 606f1580f4eSBarry Smith Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners. 6074b9ad928SBarry Smith 608562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricRight()` 6094b9ad928SBarry Smith @*/ 610d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricLeft(PC pc, Vec x, Vec y) 611d71ae5a4SJacob Faibussowitsch { 6124b9ad928SBarry Smith PetscFunctionBegin; 6130700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 6140700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 6150700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 6167827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 617e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 6189566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6199566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(x)); 6209566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplySymmetricLeft, pc, x, y, 0)); 621dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applysymmetricleft, x, y); 6229566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplySymmetricLeft, pc, x, y, 0)); 6239566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(x)); 624e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 6253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6264b9ad928SBarry Smith } 6274b9ad928SBarry Smith 6284b9ad928SBarry Smith /*@ 6294b9ad928SBarry Smith PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector. 6304b9ad928SBarry Smith 631c3339decSBarry Smith Collective 6324b9ad928SBarry Smith 6334b9ad928SBarry Smith Input Parameters: 6344b9ad928SBarry Smith + pc - the preconditioner context 6354b9ad928SBarry Smith - x - input vector 6364b9ad928SBarry Smith 6374b9ad928SBarry Smith Output Parameter: 6384b9ad928SBarry Smith . y - output vector 6394b9ad928SBarry Smith 6404b9ad928SBarry Smith Level: developer 6414b9ad928SBarry Smith 642f1580f4eSBarry Smith Note: 643f1580f4eSBarry Smith Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners. 6444b9ad928SBarry Smith 645562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricLeft()` 6464b9ad928SBarry Smith @*/ 647d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricRight(PC pc, Vec x, Vec y) 648d71ae5a4SJacob Faibussowitsch { 6494b9ad928SBarry Smith PetscFunctionBegin; 6500700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 6510700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 6520700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 6537827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 654e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 6559566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6569566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(x)); 6579566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplySymmetricRight, pc, x, y, 0)); 658dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applysymmetricright, x, y); 6599566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplySymmetricRight, pc, x, y, 0)); 6609566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(x)); 661e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 6623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6634b9ad928SBarry Smith } 6644b9ad928SBarry Smith 6654b9ad928SBarry Smith /*@ 6664b9ad928SBarry Smith PCApplyTranspose - Applies the transpose of preconditioner to a vector. 6674b9ad928SBarry Smith 668c3339decSBarry Smith Collective 6694b9ad928SBarry Smith 6704b9ad928SBarry Smith Input Parameters: 6714b9ad928SBarry Smith + pc - the preconditioner context 6724b9ad928SBarry Smith - x - input vector 6734b9ad928SBarry Smith 6744b9ad928SBarry Smith Output Parameter: 6754b9ad928SBarry Smith . y - output vector 6764b9ad928SBarry Smith 67720f4b53cSBarry Smith Level: developer 67820f4b53cSBarry Smith 679f1580f4eSBarry Smith Note: 68095452b02SPatrick Sanan For complex numbers this applies the non-Hermitian transpose. 6814c97465dSBarry Smith 682562efe2eSBarry Smith Developer Note: 683f1580f4eSBarry Smith We need to implement a `PCApplyHermitianTranspose()` 6844c97465dSBarry Smith 685562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyBAorAB()`, `PCApplyBAorABTranspose()`, `PCApplyTransposeExists()` 6864b9ad928SBarry Smith @*/ 687d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTranspose(PC pc, Vec x, Vec y) 688d71ae5a4SJacob Faibussowitsch { 6894b9ad928SBarry Smith PetscFunctionBegin; 6900700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 6910700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 6920700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 6937827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 694e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 6959566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6969566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(x)); 6979566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0)); 698dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applytranspose, x, y); 6999566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0)); 7009566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(x)); 701e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 7023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7034b9ad928SBarry Smith } 7044b9ad928SBarry Smith 7054b9ad928SBarry Smith /*@ 7069cc050e5SLisandro Dalcin PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation 7074b9ad928SBarry Smith 708c3339decSBarry Smith Collective 7094b9ad928SBarry Smith 710f1580f4eSBarry Smith Input Parameter: 7114b9ad928SBarry Smith . pc - the preconditioner context 7124b9ad928SBarry Smith 7134b9ad928SBarry Smith Output Parameter: 714f1580f4eSBarry Smith . flg - `PETSC_TRUE` if a transpose operation is defined 7154b9ad928SBarry Smith 7164b9ad928SBarry Smith Level: developer 7174b9ad928SBarry Smith 718562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()` 7194b9ad928SBarry Smith @*/ 720d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTransposeExists(PC pc, PetscBool *flg) 721d71ae5a4SJacob Faibussowitsch { 7224b9ad928SBarry Smith PetscFunctionBegin; 7230700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 7244f572ea9SToby Isaac PetscAssertPointer(flg, 2); 7256ce5e81cSLisandro Dalcin if (pc->ops->applytranspose) *flg = PETSC_TRUE; 7266ce5e81cSLisandro Dalcin else *flg = PETSC_FALSE; 7273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7284b9ad928SBarry Smith } 7294b9ad928SBarry Smith 7304b9ad928SBarry Smith /*@ 731562efe2eSBarry Smith PCApplyBAorAB - Applies the preconditioner and operator to a vector. $y = B*A*x $ or $ y = A*B*x$. 7324b9ad928SBarry Smith 733c3339decSBarry Smith Collective 7344b9ad928SBarry Smith 7354b9ad928SBarry Smith Input Parameters: 7364b9ad928SBarry Smith + pc - the preconditioner context 737f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC` 7384b9ad928SBarry Smith . x - input vector 7394b9ad928SBarry Smith - work - work vector 7404b9ad928SBarry Smith 7414b9ad928SBarry Smith Output Parameter: 7424b9ad928SBarry Smith . y - output vector 7434b9ad928SBarry Smith 7444b9ad928SBarry Smith Level: developer 7454b9ad928SBarry Smith 746f1580f4eSBarry Smith Note: 747562efe2eSBarry Smith 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. 748562efe2eSBarry Smith The specific `KSPSolve()` method must also be written to handle the post-solve "correction" for the diagonal scaling. 74913e0d083SBarry Smith 750562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorABTranspose()` 7514b9ad928SBarry Smith @*/ 752d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorAB(PC pc, PCSide side, Vec x, Vec y, Vec work) 753d71ae5a4SJacob Faibussowitsch { 7544b9ad928SBarry Smith PetscFunctionBegin; 7550700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 756a69c7061SStefano Zampini PetscValidLogicalCollectiveEnum(pc, side, 2); 7570700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 7580700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 4); 7590700a824SBarry Smith PetscValidHeaderSpecific(work, VEC_CLASSID, 5); 760a69c7061SStefano Zampini PetscCheckSameComm(pc, 1, x, 3); 761a69c7061SStefano Zampini PetscCheckSameComm(pc, 1, y, 4); 762a69c7061SStefano Zampini PetscCheckSameComm(pc, 1, work, 5); 7637827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 7647827d75bSBarry Smith PetscCheck(side == PC_LEFT || side == PC_SYMMETRIC || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right, left, or symmetric"); 7657827d75bSBarry Smith PetscCheck(!pc->diagonalscale || side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot include diagonal scaling with symmetric preconditioner application"); 766e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE)); 7674b9ad928SBarry Smith 7689566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 7694b9ad928SBarry Smith if (pc->diagonalscale) { 7704b9ad928SBarry Smith if (pc->ops->applyBA) { 7714b9ad928SBarry Smith Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */ 7729566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &work2)); 7739566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleRight(pc, x, work2)); 774dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applyBA, side, work2, y, work); 7759566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleLeft(pc, y, y)); 7769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work2)); 7774b9ad928SBarry Smith } else if (side == PC_RIGHT) { 7789566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleRight(pc, x, y)); 7799566063dSJacob Faibussowitsch PetscCall(PCApply(pc, y, work)); 7809566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, work, y)); 7819566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleLeft(pc, y, y)); 7824b9ad928SBarry Smith } else if (side == PC_LEFT) { 7839566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleRight(pc, x, y)); 7849566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, y, work)); 7859566063dSJacob Faibussowitsch PetscCall(PCApply(pc, work, y)); 7869566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleLeft(pc, y, y)); 7877827d75bSBarry Smith } else PetscCheck(side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot provide diagonal scaling with symmetric application of preconditioner"); 7884b9ad928SBarry Smith } else { 7894b9ad928SBarry Smith if (pc->ops->applyBA) { 790dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applyBA, side, x, y, work); 7914b9ad928SBarry Smith } else if (side == PC_RIGHT) { 7929566063dSJacob Faibussowitsch PetscCall(PCApply(pc, x, work)); 7939566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, work, y)); 7944b9ad928SBarry Smith } else if (side == PC_LEFT) { 7959566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, x, work)); 7969566063dSJacob Faibussowitsch PetscCall(PCApply(pc, work, y)); 7974b9ad928SBarry Smith } else if (side == PC_SYMMETRIC) { 7984b9ad928SBarry Smith /* There's an extra copy here; maybe should provide 2 work vectors instead? */ 7999566063dSJacob Faibussowitsch PetscCall(PCApplySymmetricRight(pc, x, work)); 8009566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, work, y)); 8019566063dSJacob Faibussowitsch PetscCall(VecCopy(y, work)); 8029566063dSJacob Faibussowitsch PetscCall(PCApplySymmetricLeft(pc, work, y)); 8034b9ad928SBarry Smith } 8044b9ad928SBarry Smith } 805e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE)); 8063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8074b9ad928SBarry Smith } 8084b9ad928SBarry Smith 8094b9ad928SBarry Smith /*@ 8104b9ad928SBarry Smith PCApplyBAorABTranspose - Applies the transpose of the preconditioner 811562efe2eSBarry Smith and operator to a vector. That is, applies $B^T * A^T$ with left preconditioning, 812562efe2eSBarry Smith NOT $(B*A)^T = A^T*B^T$. 8134b9ad928SBarry Smith 814c3339decSBarry Smith Collective 8154b9ad928SBarry Smith 8164b9ad928SBarry Smith Input Parameters: 8174b9ad928SBarry Smith + pc - the preconditioner context 818f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC` 8194b9ad928SBarry Smith . x - input vector 8204b9ad928SBarry Smith - work - work vector 8214b9ad928SBarry Smith 8224b9ad928SBarry Smith Output Parameter: 8234b9ad928SBarry Smith . y - output vector 8244b9ad928SBarry Smith 82520f4b53cSBarry Smith Level: developer 82620f4b53cSBarry Smith 827f1580f4eSBarry Smith Note: 828562efe2eSBarry Smith This routine is used internally so that the same Krylov code can be used to solve $A x = b$ and $A^T x = b$, with a preconditioner 829562efe2eSBarry Smith defined by $B^T$. This is why this has the funny form that it computes $B^T * A^T$ 8309b3038f0SBarry Smith 831562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorAB()` 8324b9ad928SBarry Smith @*/ 833d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorABTranspose(PC pc, PCSide side, Vec x, Vec y, Vec work) 834d71ae5a4SJacob Faibussowitsch { 8354b9ad928SBarry Smith PetscFunctionBegin; 8360700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8370700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 8380700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 4); 8390700a824SBarry Smith PetscValidHeaderSpecific(work, VEC_CLASSID, 5); 8407827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 841e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE)); 8424b9ad928SBarry Smith if (pc->ops->applyBAtranspose) { 843dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applyBAtranspose, side, x, y, work); 844e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE)); 8453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8464b9ad928SBarry Smith } 8477827d75bSBarry Smith PetscCheck(side == PC_LEFT || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right or left"); 8484b9ad928SBarry Smith 8499566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 8504b9ad928SBarry Smith if (side == PC_RIGHT) { 8519566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose(pc, x, work)); 8529566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc->mat, work, y)); 8534b9ad928SBarry Smith } else if (side == PC_LEFT) { 8549566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc->mat, x, work)); 8559566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose(pc, work, y)); 8564b9ad928SBarry Smith } 8574b9ad928SBarry Smith /* add support for PC_SYMMETRIC */ 858e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE)); 8593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8604b9ad928SBarry Smith } 8614b9ad928SBarry Smith 8624b9ad928SBarry Smith /*@ 8634b9ad928SBarry Smith PCApplyRichardsonExists - Determines whether a particular preconditioner has a 8644b9ad928SBarry Smith built-in fast application of Richardson's method. 8654b9ad928SBarry Smith 8664b9ad928SBarry Smith Not Collective 8674b9ad928SBarry Smith 8684b9ad928SBarry Smith Input Parameter: 8694b9ad928SBarry Smith . pc - the preconditioner 8704b9ad928SBarry Smith 8714b9ad928SBarry Smith Output Parameter: 872f1580f4eSBarry Smith . exists - `PETSC_TRUE` or `PETSC_FALSE` 8734b9ad928SBarry Smith 8744b9ad928SBarry Smith Level: developer 8754b9ad928SBarry Smith 87639b1ba1bSPierre Jolivet .seealso: [](ch_ksp), `PC`, `KSPRICHARDSON`, `PCApplyRichardson()` 8774b9ad928SBarry Smith @*/ 878d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyRichardsonExists(PC pc, PetscBool *exists) 879d71ae5a4SJacob Faibussowitsch { 8804b9ad928SBarry Smith PetscFunctionBegin; 8810700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8824f572ea9SToby Isaac PetscAssertPointer(exists, 2); 8834b9ad928SBarry Smith if (pc->ops->applyrichardson) *exists = PETSC_TRUE; 8844b9ad928SBarry Smith else *exists = PETSC_FALSE; 8853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8864b9ad928SBarry Smith } 8874b9ad928SBarry Smith 8884b9ad928SBarry Smith /*@ 8894b9ad928SBarry Smith PCApplyRichardson - Applies several steps of Richardson iteration with 8904b9ad928SBarry Smith the particular preconditioner. This routine is usually used by the 8914b9ad928SBarry Smith Krylov solvers and not the application code directly. 8924b9ad928SBarry Smith 893c3339decSBarry Smith Collective 8944b9ad928SBarry Smith 8954b9ad928SBarry Smith Input Parameters: 8964b9ad928SBarry Smith + pc - the preconditioner context 897dd8e379bSPierre Jolivet . b - the right-hand side 8984b9ad928SBarry Smith . w - one work vector 8994b9ad928SBarry Smith . rtol - relative decrease in residual norm convergence criteria 90070441072SBarry Smith . abstol - absolute residual norm convergence criteria 9014b9ad928SBarry Smith . dtol - divergence residual norm increase criteria 9027319c654SBarry Smith . its - the number of iterations to apply. 9037319c654SBarry Smith - guesszero - if the input x contains nonzero initial guess 9044b9ad928SBarry Smith 905d8d19677SJose E. Roman Output Parameters: 9064d0a8057SBarry Smith + outits - number of iterations actually used (for SOR this always equals its) 9074d0a8057SBarry Smith . reason - the reason the apply terminated 908f1580f4eSBarry Smith - y - the solution (also contains initial guess if guesszero is `PETSC_FALSE` 9094b9ad928SBarry Smith 91020f4b53cSBarry Smith Level: developer 91120f4b53cSBarry Smith 9124b9ad928SBarry Smith Notes: 9134b9ad928SBarry Smith Most preconditioners do not support this function. Use the command 914f1580f4eSBarry Smith `PCApplyRichardsonExists()` to determine if one does. 9154b9ad928SBarry Smith 916f1580f4eSBarry Smith Except for the `PCMG` this routine ignores the convergence tolerances 9174b9ad928SBarry Smith and always runs for the number of iterations 9184b9ad928SBarry Smith 919562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyRichardsonExists()` 9204b9ad928SBarry Smith @*/ 921d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyRichardson(PC pc, Vec b, Vec y, Vec w, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt its, PetscBool guesszero, PetscInt *outits, PCRichardsonConvergedReason *reason) 922d71ae5a4SJacob Faibussowitsch { 9234b9ad928SBarry Smith PetscFunctionBegin; 9240700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9250700a824SBarry Smith PetscValidHeaderSpecific(b, VEC_CLASSID, 2); 9260700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 9270700a824SBarry Smith PetscValidHeaderSpecific(w, VEC_CLASSID, 4); 9287827d75bSBarry Smith PetscCheck(b != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "b and y must be different vectors"); 9299566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 930dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applyrichardson, b, y, w, rtol, abstol, dtol, its, guesszero, outits, reason); 9313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9324b9ad928SBarry Smith } 9334b9ad928SBarry Smith 934422a814eSBarry Smith /*@ 935f1580f4eSBarry Smith PCSetFailedReason - Sets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail 9361b2b9847SBarry Smith 937c3339decSBarry Smith Logically Collective 9381b2b9847SBarry Smith 939d8d19677SJose E. Roman Input Parameters: 9401b2b9847SBarry Smith + pc - the preconditioner context 9411b2b9847SBarry Smith - reason - the reason it failedx 9421b2b9847SBarry Smith 9431b2b9847SBarry Smith Level: advanced 9441b2b9847SBarry Smith 945562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCFailedReason` 9461b2b9847SBarry Smith @*/ 947d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFailedReason(PC pc, PCFailedReason reason) 948d71ae5a4SJacob Faibussowitsch { 9491b2b9847SBarry Smith PetscFunctionBegin; 9506479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9511b2b9847SBarry Smith pc->failedreason = reason; 9523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9531b2b9847SBarry Smith } 9541b2b9847SBarry Smith 9551b2b9847SBarry Smith /*@ 956f1580f4eSBarry Smith PCGetFailedReason - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail 957422a814eSBarry Smith 958c3339decSBarry Smith Logically Collective 959422a814eSBarry Smith 960422a814eSBarry Smith Input Parameter: 961422a814eSBarry Smith . pc - the preconditioner context 962422a814eSBarry Smith 963422a814eSBarry Smith Output Parameter: 9641b2b9847SBarry Smith . reason - the reason it failed 965422a814eSBarry Smith 966422a814eSBarry Smith Level: advanced 967422a814eSBarry Smith 968f1580f4eSBarry Smith Note: 969f1580f4eSBarry Smith This is the maximum over reason over all ranks in the PC communicator. It is only valid after 9706479e835SStefano Zampini a call `KSPCheckDot()` or `KSPCheckNorm()` inside a `KSPSolve()` or `PCReduceFailedReason()`. 9716479e835SStefano Zampini It is not valid immediately after a `PCSetUp()` or `PCApply()`, then use `PCGetFailedReasonRank()` 9721b2b9847SBarry Smith 973a94f484eSPierre Jolivet .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReasonRank()`, `PCSetFailedReason()` 974422a814eSBarry Smith @*/ 975d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason) 976d71ae5a4SJacob Faibussowitsch { 977422a814eSBarry Smith PetscFunctionBegin; 9786479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9799754764cSHong Zhang if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled; 9800ae0b32dSHong Zhang else *reason = pc->failedreason; 9813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 982422a814eSBarry Smith } 983422a814eSBarry Smith 9841b2b9847SBarry Smith /*@ 985f1580f4eSBarry Smith PCGetFailedReasonRank - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail on this MPI rank 9861b2b9847SBarry Smith 987f1580f4eSBarry Smith Not Collective 9881b2b9847SBarry Smith 9891b2b9847SBarry Smith Input Parameter: 9901b2b9847SBarry Smith . pc - the preconditioner context 9911b2b9847SBarry Smith 9921b2b9847SBarry Smith Output Parameter: 9931b2b9847SBarry Smith . reason - the reason it failed 9941b2b9847SBarry Smith 99520f4b53cSBarry Smith Level: advanced 99620f4b53cSBarry Smith 997f1580f4eSBarry Smith Note: 998562efe2eSBarry Smith Different processes may have different reasons or no reason, see `PCGetFailedReason()` 9991b2b9847SBarry Smith 1000562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`, `PCReduceFailedReason()` 10011b2b9847SBarry Smith @*/ 1002d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReasonRank(PC pc, PCFailedReason *reason) 1003d71ae5a4SJacob Faibussowitsch { 10041b2b9847SBarry Smith PetscFunctionBegin; 10056479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 10061b2b9847SBarry Smith if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled; 10071b2b9847SBarry Smith else *reason = pc->failedreason; 10083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10091b2b9847SBarry Smith } 1010422a814eSBarry Smith 10116479e835SStefano Zampini /*@ 10126479e835SStefano Zampini PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC` 10136479e835SStefano Zampini 10146479e835SStefano Zampini Collective 10156479e835SStefano Zampini 10166479e835SStefano Zampini Input Parameter: 10176479e835SStefano Zampini . pc - the preconditioner context 10186479e835SStefano Zampini 10196479e835SStefano Zampini Level: advanced 10206479e835SStefano Zampini 10216479e835SStefano Zampini Note: 1022562efe2eSBarry Smith Different MPI processes may have different reasons or no reason, see `PCGetFailedReason()`. This routine 10236479e835SStefano Zampini makes them have a common value (failure if any MPI process had a failure). 10246479e835SStefano Zampini 1025562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()` 10266479e835SStefano Zampini @*/ 10276479e835SStefano Zampini PetscErrorCode PCReduceFailedReason(PC pc) 10286479e835SStefano Zampini { 10296479e835SStefano Zampini PetscInt buf; 10306479e835SStefano Zampini 10316479e835SStefano Zampini PetscFunctionBegin; 10326479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 10336479e835SStefano Zampini buf = (PetscInt)pc->failedreason; 10346479e835SStefano Zampini PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc))); 10356479e835SStefano Zampini pc->failedreason = (PCFailedReason)buf; 10366479e835SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 10376479e835SStefano Zampini } 10386479e835SStefano Zampini 10394b9ad928SBarry Smith /* 10404b9ad928SBarry Smith a setupcall of 0 indicates never setup, 104123ee1639SBarry Smith 1 indicates has been previously setup 1042422a814eSBarry Smith -1 indicates a PCSetUp() was attempted and failed 10434b9ad928SBarry Smith */ 10444b9ad928SBarry Smith /*@ 10454b9ad928SBarry Smith PCSetUp - Prepares for the use of a preconditioner. 10464b9ad928SBarry Smith 1047c3339decSBarry Smith Collective 10484b9ad928SBarry Smith 10494b9ad928SBarry Smith Input Parameter: 10504b9ad928SBarry Smith . pc - the preconditioner context 10514b9ad928SBarry Smith 10524b9ad928SBarry Smith Level: developer 10534b9ad928SBarry Smith 1054562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()` 10554b9ad928SBarry Smith @*/ 1056d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp(PC pc) 1057d71ae5a4SJacob Faibussowitsch { 1058566e8bf2SBarry Smith const char *def; 1059fc9b51d3SBarry Smith PetscObjectState matstate, matnonzerostate; 10604b9ad928SBarry Smith 10614b9ad928SBarry Smith PetscFunctionBegin; 10620700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 106328b400f6SJacob Faibussowitsch PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first"); 10646ce5e81cSLisandro Dalcin 106523ee1639SBarry Smith if (pc->setupcalled && pc->reusepreconditioner) { 10669566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n")); 10673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10684b9ad928SBarry Smith } 10694b9ad928SBarry Smith 10709566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate)); 10719566063dSJacob Faibussowitsch PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate)); 107223ee1639SBarry Smith if (!pc->setupcalled) { 10735e62d202SMark Adams //PetscCall(PetscInfo(pc, "Setting up PC for first time\n")); 107423ee1639SBarry Smith pc->flag = DIFFERENT_NONZERO_PATTERN; 10759df67409SStefano Zampini } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS); 10769df67409SStefano Zampini else { 10779df67409SStefano Zampini if (matnonzerostate != pc->matnonzerostate) { 10789566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n")); 107923ee1639SBarry Smith pc->flag = DIFFERENT_NONZERO_PATTERN; 108023ee1639SBarry Smith } else { 10815e62d202SMark Adams //PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n")); 108223ee1639SBarry Smith pc->flag = SAME_NONZERO_PATTERN; 108323ee1639SBarry Smith } 108423ee1639SBarry Smith } 108523ee1639SBarry Smith pc->matstate = matstate; 1086fc9b51d3SBarry Smith pc->matnonzerostate = matnonzerostate; 108723ee1639SBarry Smith 10887adad957SLisandro Dalcin if (!((PetscObject)pc)->type_name) { 10899566063dSJacob Faibussowitsch PetscCall(PCGetDefaultType_Private(pc, &def)); 10909566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, def)); 1091566e8bf2SBarry Smith } 10924b9ad928SBarry Smith 10939566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure)); 10949566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure)); 10959566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0)); 10964b9ad928SBarry Smith if (pc->ops->setup) { 1097*fe57bb1aSStefano Zampini PetscCall(PCLogEventsDeactivatePush()); 1098dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, setup); 1099*fe57bb1aSStefano Zampini PetscCall(PCLogEventsDeactivatePop()); 11004b9ad928SBarry Smith } 11019566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0)); 1102422a814eSBarry Smith if (!pc->setupcalled) pc->setupcalled = 1; 11033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11044b9ad928SBarry Smith } 11054b9ad928SBarry Smith 11064b9ad928SBarry Smith /*@ 11074b9ad928SBarry Smith PCSetUpOnBlocks - Sets up the preconditioner for each block in 11084b9ad928SBarry Smith the block Jacobi, block Gauss-Seidel, and overlapping Schwarz 11094b9ad928SBarry Smith methods. 11104b9ad928SBarry Smith 1111c3339decSBarry Smith Collective 11124b9ad928SBarry Smith 1113f1580f4eSBarry Smith Input Parameter: 11144b9ad928SBarry Smith . pc - the preconditioner context 11154b9ad928SBarry Smith 11164b9ad928SBarry Smith Level: developer 11174b9ad928SBarry Smith 1118f1580f4eSBarry Smith Note: 1119f1580f4eSBarry Smith For nested preconditioners such as `PCBJACOBI` `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is 1120f1580f4eSBarry Smith called on the outer `PC`, this routine ensures it is called. 1121f1580f4eSBarry Smith 1122562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()` 11234b9ad928SBarry Smith @*/ 1124d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUpOnBlocks(PC pc) 1125d71ae5a4SJacob Faibussowitsch { 11264b9ad928SBarry Smith PetscFunctionBegin; 11270700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 11283ba16761SJacob Faibussowitsch if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS); 11299566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0)); 1130*fe57bb1aSStefano Zampini PetscCall(PCLogEventsDeactivatePush()); 1131dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, setuponblocks); 1132*fe57bb1aSStefano Zampini PetscCall(PCLogEventsDeactivatePop()); 11339566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0)); 11343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11354b9ad928SBarry Smith } 11364b9ad928SBarry Smith 11374b9ad928SBarry Smith /*@C 11384b9ad928SBarry Smith PCSetModifySubMatrices - Sets a user-defined routine for modifying the 113904c3f3b8SBarry Smith submatrices that arise within certain subdomain-based preconditioners such as `PCASM` 11404b9ad928SBarry Smith 1141c3339decSBarry Smith Logically Collective 11424b9ad928SBarry Smith 11434b9ad928SBarry Smith Input Parameters: 11444b9ad928SBarry Smith + pc - the preconditioner context 11454b9ad928SBarry Smith . func - routine for modifying the submatrices 114604c3f3b8SBarry Smith - ctx - optional user-defined context (may be `NULL`) 11474b9ad928SBarry Smith 114820f4b53cSBarry Smith Calling sequence of `func`: 114920f4b53cSBarry Smith + pc - the preconditioner context 115004c3f3b8SBarry Smith . nsub - number of index sets 115120f4b53cSBarry Smith . row - an array of index sets that contain the global row numbers 11524b9ad928SBarry Smith that comprise each local submatrix 11534b9ad928SBarry Smith . col - an array of index sets that contain the global column numbers 11544b9ad928SBarry Smith that comprise each local submatrix 11554b9ad928SBarry Smith . submat - array of local submatrices 11564b9ad928SBarry Smith - ctx - optional user-defined context for private data for the 115704c3f3b8SBarry Smith user-defined func routine (may be `NULL`) 11584b9ad928SBarry Smith 115920f4b53cSBarry Smith Level: advanced 116020f4b53cSBarry Smith 11614b9ad928SBarry Smith Notes: 116204c3f3b8SBarry Smith The basic submatrices are extracted from the preconditioner matrix as 116304c3f3b8SBarry Smith usual; the user can then alter these (for example, to set different boundary 116404c3f3b8SBarry Smith conditions for each submatrix) before they are used for the local solves. 116504c3f3b8SBarry Smith 1166f1580f4eSBarry Smith `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and 1167f1580f4eSBarry Smith `KSPSolve()`. 11684b9ad928SBarry Smith 1169f1580f4eSBarry Smith A routine set by `PCSetModifySubMatrices()` is currently called within 1170f1580f4eSBarry Smith the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`) 11714b9ad928SBarry Smith preconditioners. All other preconditioners ignore this routine. 11724b9ad928SBarry Smith 1173562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()` 11744b9ad928SBarry Smith @*/ 117504c3f3b8SBarry Smith PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx), void *ctx) 1176d71ae5a4SJacob Faibussowitsch { 11774b9ad928SBarry Smith PetscFunctionBegin; 11780700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 11794b9ad928SBarry Smith pc->modifysubmatrices = func; 11804b9ad928SBarry Smith pc->modifysubmatricesP = ctx; 11813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11824b9ad928SBarry Smith } 11834b9ad928SBarry Smith 11844b9ad928SBarry Smith /*@C 11854b9ad928SBarry Smith PCModifySubMatrices - Calls an optional user-defined routine within 1186f1580f4eSBarry Smith certain preconditioners if one has been set with `PCSetModifySubMatrices()`. 11874b9ad928SBarry Smith 1188c3339decSBarry Smith Collective 11894b9ad928SBarry Smith 11904b9ad928SBarry Smith Input Parameters: 11914b9ad928SBarry Smith + pc - the preconditioner context 11924b9ad928SBarry Smith . nsub - the number of local submatrices 11934b9ad928SBarry Smith . row - an array of index sets that contain the global row numbers 11944b9ad928SBarry Smith that comprise each local submatrix 11954b9ad928SBarry Smith . col - an array of index sets that contain the global column numbers 11964b9ad928SBarry Smith that comprise each local submatrix 11974b9ad928SBarry Smith . submat - array of local submatrices 11984b9ad928SBarry Smith - ctx - optional user-defined context for private data for the 1199562efe2eSBarry Smith user-defined routine (may be `NULL`) 12004b9ad928SBarry Smith 12014b9ad928SBarry Smith Output Parameter: 12024b9ad928SBarry Smith . submat - array of local submatrices (the entries of which may 12034b9ad928SBarry Smith have been modified) 12044b9ad928SBarry Smith 120520f4b53cSBarry Smith Level: developer 120620f4b53cSBarry Smith 120704c3f3b8SBarry Smith Note: 12084b9ad928SBarry Smith The user should NOT generally call this routine, as it will 120904c3f3b8SBarry Smith automatically be called within certain preconditioners. 12104b9ad928SBarry Smith 1211562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetModifySubMatrices()` 12124b9ad928SBarry Smith @*/ 1213d71ae5a4SJacob Faibussowitsch PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx) 1214d71ae5a4SJacob Faibussowitsch { 12154b9ad928SBarry Smith PetscFunctionBegin; 12160700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 12173ba16761SJacob Faibussowitsch if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS); 12189566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0)); 12199566063dSJacob Faibussowitsch PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx)); 12209566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0)); 12213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12224b9ad928SBarry Smith } 12234b9ad928SBarry Smith 12244b9ad928SBarry Smith /*@ 12254b9ad928SBarry Smith PCSetOperators - Sets the matrix associated with the linear system and 12264b9ad928SBarry Smith a (possibly) different one associated with the preconditioner. 12274b9ad928SBarry Smith 1228c3339decSBarry Smith Logically Collective 12294b9ad928SBarry Smith 12304b9ad928SBarry Smith Input Parameters: 12314b9ad928SBarry Smith + pc - the preconditioner context 1232e5d3d808SBarry Smith . Amat - the matrix that defines the linear system 1233d1e9a80fSBarry Smith - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat. 12344b9ad928SBarry Smith 123520f4b53cSBarry Smith Level: intermediate 1236189c0b8aSBarry Smith 123720f4b53cSBarry Smith Notes: 123820f4b53cSBarry Smith Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used. 123920f4b53cSBarry Smith 124020f4b53cSBarry Smith If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then 1241f1580f4eSBarry Smith first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()` 1242f1580f4eSBarry Smith on it and then pass it back in in your call to `KSPSetOperators()`. 1243189c0b8aSBarry Smith 12444b9ad928SBarry Smith More Notes about Repeated Solution of Linear Systems: 124520f4b53cSBarry Smith PETSc does NOT reset the matrix entries of either `Amat` or `Pmat` 12464b9ad928SBarry Smith to zero after a linear solve; the user is completely responsible for 1247f1580f4eSBarry Smith matrix assembly. See the routine `MatZeroEntries()` if desiring to 12484b9ad928SBarry Smith zero all elements of a matrix. 12494b9ad928SBarry Smith 1250562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()` 12514b9ad928SBarry Smith @*/ 1252d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat) 1253d71ae5a4SJacob Faibussowitsch { 12543b3f6333SBarry Smith PetscInt m1, n1, m2, n2; 12554b9ad928SBarry Smith 12564b9ad928SBarry Smith PetscFunctionBegin; 12570700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 12580700a824SBarry Smith if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 12590700a824SBarry Smith if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3); 12603fc8bf9cSBarry Smith if (Amat) PetscCheckSameComm(pc, 1, Amat, 2); 12613fc8bf9cSBarry Smith if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3); 126231641f1aSBarry Smith if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) { 12639566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Amat, &m1, &n1)); 12649566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->mat, &m2, &n2)); 126563a3b9bcSJacob Faibussowitsch PetscCheck(m1 == m2 && n1 == n2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot change local size of Amat after use old sizes %" PetscInt_FMT " %" PetscInt_FMT " new sizes %" PetscInt_FMT " %" PetscInt_FMT, m2, n2, m1, n1); 12669566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Pmat, &m1, &n1)); 12679566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2)); 126863a3b9bcSJacob Faibussowitsch PetscCheck(m1 == m2 && n1 == n2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot change local size of Pmat after use old sizes %" PetscInt_FMT " %" PetscInt_FMT " new sizes %" PetscInt_FMT " %" PetscInt_FMT, m2, n2, m1, n1); 12693b3f6333SBarry Smith } 12704b9ad928SBarry Smith 127123ee1639SBarry Smith if (Pmat != pc->pmat) { 127223ee1639SBarry Smith /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */ 127323ee1639SBarry Smith pc->matnonzerostate = -1; 127423ee1639SBarry Smith pc->matstate = -1; 127523ee1639SBarry Smith } 127623ee1639SBarry Smith 1277906ed7ccSBarry Smith /* reference first in case the matrices are the same */ 12789566063dSJacob Faibussowitsch if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat)); 12799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->mat)); 12809566063dSJacob Faibussowitsch if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat)); 12819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->pmat)); 12824b9ad928SBarry Smith pc->mat = Amat; 12834b9ad928SBarry Smith pc->pmat = Pmat; 12843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1285e56f5c9eSBarry Smith } 1286e56f5c9eSBarry Smith 128723ee1639SBarry Smith /*@ 128823ee1639SBarry Smith PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed. 128923ee1639SBarry Smith 1290c3339decSBarry Smith Logically Collective 129123ee1639SBarry Smith 129223ee1639SBarry Smith Input Parameters: 129323ee1639SBarry Smith + pc - the preconditioner context 1294f1580f4eSBarry Smith - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner 129523ee1639SBarry Smith 12962b26979fSBarry Smith Level: intermediate 12972b26979fSBarry Smith 1298f1580f4eSBarry Smith Note: 1299f1580f4eSBarry Smith Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option 1300f1580f4eSBarry Smith prevents this. 1301f1580f4eSBarry Smith 1302562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()` 130323ee1639SBarry Smith @*/ 1304d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag) 1305d71ae5a4SJacob Faibussowitsch { 130623ee1639SBarry Smith PetscFunctionBegin; 130723ee1639SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1308f9177522SStefano Zampini PetscValidLogicalCollectiveBool(pc, flag, 2); 130923ee1639SBarry Smith pc->reusepreconditioner = flag; 13103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13114b9ad928SBarry Smith } 13124b9ad928SBarry Smith 1313c60c7ad4SBarry Smith /*@ 1314f1580f4eSBarry Smith PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed. 1315c60c7ad4SBarry Smith 1316bba28a21SBarry Smith Not Collective 1317c60c7ad4SBarry Smith 1318c60c7ad4SBarry Smith Input Parameter: 1319c60c7ad4SBarry Smith . pc - the preconditioner context 1320c60c7ad4SBarry Smith 1321c60c7ad4SBarry Smith Output Parameter: 1322f1580f4eSBarry Smith . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner 1323c60c7ad4SBarry Smith 1324d0418729SSatish Balay Level: intermediate 1325d0418729SSatish Balay 1326562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()` 1327c60c7ad4SBarry Smith @*/ 1328d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag) 1329d71ae5a4SJacob Faibussowitsch { 1330c60c7ad4SBarry Smith PetscFunctionBegin; 1331c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 13324f572ea9SToby Isaac PetscAssertPointer(flag, 2); 1333c60c7ad4SBarry Smith *flag = pc->reusepreconditioner; 13343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1335c60c7ad4SBarry Smith } 1336c60c7ad4SBarry Smith 1337487a658cSBarry Smith /*@ 13384b9ad928SBarry Smith PCGetOperators - Gets the matrix associated with the linear system and 13394b9ad928SBarry Smith possibly a different one associated with the preconditioner. 13404b9ad928SBarry Smith 1341562efe2eSBarry Smith Not Collective, though parallel `Mat`s are returned if `pc` is parallel 13424b9ad928SBarry Smith 13434b9ad928SBarry Smith Input Parameter: 13444b9ad928SBarry Smith . pc - the preconditioner context 13454b9ad928SBarry Smith 13464b9ad928SBarry Smith Output Parameters: 1347e5d3d808SBarry Smith + Amat - the matrix defining the linear system 134823ee1639SBarry Smith - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat. 13494b9ad928SBarry Smith 13504b9ad928SBarry Smith Level: intermediate 13514b9ad928SBarry Smith 1352f1580f4eSBarry Smith Note: 135395452b02SPatrick Sanan Does not increase the reference count of the matrices, so you should not destroy them 1354298cc208SBarry Smith 1355f1580f4eSBarry Smith Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators 1356f1580f4eSBarry Smith are created in `PC` and returned to the user. In this case, if both operators 135773950996SBarry Smith mat and pmat are requested, two DIFFERENT operators will be returned. If 135873950996SBarry Smith only one is requested both operators in the PC will be the same (i.e. as 1359f1580f4eSBarry Smith if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats). 136073950996SBarry Smith The user must set the sizes of the returned matrices and their type etc just 1361f1580f4eSBarry Smith as if the user created them with `MatCreate()`. For example, 136273950996SBarry Smith 1363f1580f4eSBarry Smith .vb 1364f1580f4eSBarry Smith KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to 1365f1580f4eSBarry Smith set size, type, etc of Amat 136673950996SBarry Smith 1367f1580f4eSBarry Smith MatCreate(comm,&mat); 1368f1580f4eSBarry Smith KSP/PCSetOperators(ksp/pc,Amat,Amat); 1369f1580f4eSBarry Smith PetscObjectDereference((PetscObject)mat); 1370f1580f4eSBarry Smith set size, type, etc of Amat 1371f1580f4eSBarry Smith .ve 137273950996SBarry Smith 137373950996SBarry Smith and 137473950996SBarry Smith 1375f1580f4eSBarry Smith .vb 1376f1580f4eSBarry Smith KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to 1377f1580f4eSBarry Smith set size, type, etc of Amat and Pmat 137873950996SBarry Smith 1379f1580f4eSBarry Smith MatCreate(comm,&Amat); 1380f1580f4eSBarry Smith MatCreate(comm,&Pmat); 1381f1580f4eSBarry Smith KSP/PCSetOperators(ksp/pc,Amat,Pmat); 1382f1580f4eSBarry Smith PetscObjectDereference((PetscObject)Amat); 1383f1580f4eSBarry Smith PetscObjectDereference((PetscObject)Pmat); 1384f1580f4eSBarry Smith set size, type, etc of Amat and Pmat 1385f1580f4eSBarry Smith .ve 138673950996SBarry Smith 1387f1580f4eSBarry Smith The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy 1388b8d709abSRichard Tran Mills of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely 1389f1580f4eSBarry Smith managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look 1390f1580f4eSBarry Smith at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to 1391f1580f4eSBarry Smith the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP 1392f1580f4eSBarry Smith you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you). 1393f1580f4eSBarry Smith Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when 139473950996SBarry Smith it can be created for you? 139573950996SBarry Smith 1396562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()` 13974b9ad928SBarry Smith @*/ 1398d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat) 1399d71ae5a4SJacob Faibussowitsch { 14004b9ad928SBarry Smith PetscFunctionBegin; 14010700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1402e5d3d808SBarry Smith if (Amat) { 140373950996SBarry Smith if (!pc->mat) { 14049fca8c71SStefano Zampini if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */ 14059a4708feSJed Brown pc->mat = pc->pmat; 14069566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->mat)); 1407e5d3d808SBarry Smith } else { /* both Amat and Pmat are empty */ 14089566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat)); 1409e5d3d808SBarry Smith if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */ 141073950996SBarry Smith pc->pmat = pc->mat; 14119566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->pmat)); 141273950996SBarry Smith } 141373950996SBarry Smith } 14149a4708feSJed Brown } 1415e5d3d808SBarry Smith *Amat = pc->mat; 141673950996SBarry Smith } 1417e5d3d808SBarry Smith if (Pmat) { 141873950996SBarry Smith if (!pc->pmat) { 1419e5d3d808SBarry Smith if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */ 14209a4708feSJed Brown pc->pmat = pc->mat; 14219566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->pmat)); 14229a4708feSJed Brown } else { 14239566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat)); 1424e5d3d808SBarry Smith if (!Amat) { /* user did NOT request Amat, so make same as Pmat */ 142573950996SBarry Smith pc->mat = pc->pmat; 14269566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->mat)); 142773950996SBarry Smith } 142873950996SBarry Smith } 14299a4708feSJed Brown } 1430e5d3d808SBarry Smith *Pmat = pc->pmat; 143173950996SBarry Smith } 14323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14334b9ad928SBarry Smith } 14344b9ad928SBarry Smith 14355d83a8b1SBarry Smith /*@ 1436906ed7ccSBarry Smith PCGetOperatorsSet - Determines if the matrix associated with the linear system and 1437f1580f4eSBarry Smith possibly a different one associated with the preconditioner have been set in the `PC`. 1438906ed7ccSBarry Smith 143920f4b53cSBarry Smith Not Collective, though the results on all processes should be the same 1440906ed7ccSBarry Smith 1441906ed7ccSBarry Smith Input Parameter: 1442906ed7ccSBarry Smith . pc - the preconditioner context 1443906ed7ccSBarry Smith 1444906ed7ccSBarry Smith Output Parameters: 1445906ed7ccSBarry Smith + mat - the matrix associated with the linear system was set 1446906ed7ccSBarry Smith - pmat - matrix associated with the preconditioner was set, usually the same 1447906ed7ccSBarry Smith 1448906ed7ccSBarry Smith Level: intermediate 1449906ed7ccSBarry Smith 1450562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()` 1451906ed7ccSBarry Smith @*/ 1452d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat) 1453d71ae5a4SJacob Faibussowitsch { 1454906ed7ccSBarry Smith PetscFunctionBegin; 14550700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1456906ed7ccSBarry Smith if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE; 1457906ed7ccSBarry Smith if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE; 14583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1459906ed7ccSBarry Smith } 1460906ed7ccSBarry Smith 1461f39d8e23SSatish Balay /*@ 1462a4fd02acSBarry Smith PCFactorGetMatrix - Gets the factored matrix from the 1463f1580f4eSBarry Smith preconditioner context. This routine is valid only for the `PCLU`, 1464f1580f4eSBarry Smith `PCILU`, `PCCHOLESKY`, and `PCICC` methods. 14654b9ad928SBarry Smith 1466562efe2eSBarry Smith Not Collective though `mat` is parallel if `pc` is parallel 14674b9ad928SBarry Smith 1468f1580f4eSBarry Smith Input Parameter: 14694b9ad928SBarry Smith . pc - the preconditioner context 14704b9ad928SBarry Smith 1471feefa0e1SJacob Faibussowitsch Output Parameters: 14724b9ad928SBarry Smith . mat - the factored matrix 14734b9ad928SBarry Smith 14744b9ad928SBarry Smith Level: advanced 14754b9ad928SBarry Smith 1476f1580f4eSBarry Smith Note: 1477562efe2eSBarry Smith Does not increase the reference count for `mat` so DO NOT destroy it 14789405f653SBarry Smith 1479562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC` 14804b9ad928SBarry Smith @*/ 1481d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat) 1482d71ae5a4SJacob Faibussowitsch { 14834b9ad928SBarry Smith PetscFunctionBegin; 14840700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 14854f572ea9SToby Isaac PetscAssertPointer(mat, 2); 14867def7855SStefano Zampini PetscCall(PCFactorSetUpMatSolverType(pc)); 1487dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, getfactoredmatrix, mat); 14883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14894b9ad928SBarry Smith } 14904b9ad928SBarry Smith 1491cc4c1da9SBarry Smith /*@ 14924b9ad928SBarry Smith PCSetOptionsPrefix - Sets the prefix used for searching for all 1493f1580f4eSBarry Smith `PC` options in the database. 14944b9ad928SBarry Smith 1495c3339decSBarry Smith Logically Collective 14964b9ad928SBarry Smith 14974b9ad928SBarry Smith Input Parameters: 14984b9ad928SBarry Smith + pc - the preconditioner context 1499f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests 15004b9ad928SBarry Smith 1501f1580f4eSBarry Smith Note: 15024b9ad928SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 15034b9ad928SBarry Smith The first character of all runtime options is AUTOMATICALLY the 15044b9ad928SBarry Smith hyphen. 15054b9ad928SBarry Smith 15064b9ad928SBarry Smith Level: advanced 15074b9ad928SBarry Smith 1508562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()` 15094b9ad928SBarry Smith @*/ 1510d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[]) 1511d71ae5a4SJacob Faibussowitsch { 15124b9ad928SBarry Smith PetscFunctionBegin; 15130700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix)); 15153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15164b9ad928SBarry Smith } 15174b9ad928SBarry Smith 1518cc4c1da9SBarry Smith /*@ 15194b9ad928SBarry Smith PCAppendOptionsPrefix - Appends to the prefix used for searching for all 1520f1580f4eSBarry Smith `PC` options in the database. 15214b9ad928SBarry Smith 1522c3339decSBarry Smith Logically Collective 15234b9ad928SBarry Smith 15244b9ad928SBarry Smith Input Parameters: 15254b9ad928SBarry Smith + pc - the preconditioner context 1526f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests 15274b9ad928SBarry Smith 1528f1580f4eSBarry Smith Note: 15294b9ad928SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 15304b9ad928SBarry Smith The first character of all runtime options is AUTOMATICALLY the 15314b9ad928SBarry Smith hyphen. 15324b9ad928SBarry Smith 15334b9ad928SBarry Smith Level: advanced 15344b9ad928SBarry Smith 1535562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()` 15364b9ad928SBarry Smith @*/ 1537d71ae5a4SJacob Faibussowitsch PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[]) 1538d71ae5a4SJacob Faibussowitsch { 15394b9ad928SBarry Smith PetscFunctionBegin; 15400700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15419566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix)); 15423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15434b9ad928SBarry Smith } 15444b9ad928SBarry Smith 1545cc4c1da9SBarry Smith /*@ 15464b9ad928SBarry Smith PCGetOptionsPrefix - Gets the prefix used for searching for all 15474b9ad928SBarry Smith PC options in the database. 15484b9ad928SBarry Smith 15494b9ad928SBarry Smith Not Collective 15504b9ad928SBarry Smith 1551f1580f4eSBarry Smith Input Parameter: 15524b9ad928SBarry Smith . pc - the preconditioner context 15534b9ad928SBarry Smith 1554f1580f4eSBarry Smith Output Parameter: 15554b9ad928SBarry Smith . prefix - pointer to the prefix string used, is returned 15564b9ad928SBarry Smith 15574b9ad928SBarry Smith Level: advanced 15584b9ad928SBarry Smith 1559562efe2eSBarry Smith Fortran Note: 1560562efe2eSBarry Smith The user should pass in a string `prefix` of 1561562efe2eSBarry Smith sufficient length to hold the prefix. 1562562efe2eSBarry Smith 1563562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()` 15644b9ad928SBarry Smith @*/ 1565d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[]) 1566d71ae5a4SJacob Faibussowitsch { 15674b9ad928SBarry Smith PetscFunctionBegin; 15680700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15694f572ea9SToby Isaac PetscAssertPointer(prefix, 2); 15709566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix)); 15713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15724b9ad928SBarry Smith } 15734b9ad928SBarry Smith 15748066bbecSBarry Smith /* 1575dd8e379bSPierre Jolivet Indicates the right-hand side will be changed by KSPSolve(), this occurs for a few 15768066bbecSBarry Smith preconditioners including BDDC and Eisentat that transform the equations before applying 15778066bbecSBarry Smith the Krylov methods 15788066bbecSBarry Smith */ 1579d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change) 1580d71ae5a4SJacob Faibussowitsch { 1581a06fd7f2SStefano Zampini PetscFunctionBegin; 1582a06fd7f2SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15834f572ea9SToby Isaac PetscAssertPointer(change, 2); 1584a06fd7f2SStefano Zampini *change = PETSC_FALSE; 1585cac4c232SBarry Smith PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change)); 15863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1587a06fd7f2SStefano Zampini } 1588a06fd7f2SStefano Zampini 15894b9ad928SBarry Smith /*@ 15904b9ad928SBarry Smith PCPreSolve - Optional pre-solve phase, intended for any 15914b9ad928SBarry Smith preconditioner-specific actions that must be performed before 15924b9ad928SBarry Smith the iterative solve itself. 15934b9ad928SBarry Smith 1594c3339decSBarry Smith Collective 15954b9ad928SBarry Smith 15964b9ad928SBarry Smith Input Parameters: 15974b9ad928SBarry Smith + pc - the preconditioner context 15984b9ad928SBarry Smith - ksp - the Krylov subspace context 15994b9ad928SBarry Smith 16004b9ad928SBarry Smith Level: developer 16014b9ad928SBarry Smith 1602feefa0e1SJacob Faibussowitsch Example Usage: 16034b9ad928SBarry Smith .vb 16044b9ad928SBarry Smith PCPreSolve(pc,ksp); 160523ce1328SBarry Smith KSPSolve(ksp,b,x); 16064b9ad928SBarry Smith PCPostSolve(pc,ksp); 16074b9ad928SBarry Smith .ve 16084b9ad928SBarry Smith 16094b9ad928SBarry Smith Notes: 1610f1580f4eSBarry Smith The pre-solve phase is distinct from the `PCSetUp()` phase. 16114b9ad928SBarry Smith 1612f1580f4eSBarry Smith `KSPSolve()` calls this directly, so is rarely called by the user. 16134b9ad928SBarry Smith 1614562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCPostSolve()` 16154b9ad928SBarry Smith @*/ 1616d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPreSolve(PC pc, KSP ksp) 1617d71ae5a4SJacob Faibussowitsch { 16184b9ad928SBarry Smith Vec x, rhs; 16194b9ad928SBarry Smith 16204b9ad928SBarry Smith PetscFunctionBegin; 16210700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16220700a824SBarry Smith PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2); 1623d9a03883SBarry Smith pc->presolvedone++; 16247827d75bSBarry Smith PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice"); 16259566063dSJacob Faibussowitsch PetscCall(KSPGetSolution(ksp, &x)); 16269566063dSJacob Faibussowitsch PetscCall(KSPGetRhs(ksp, &rhs)); 16274b9ad928SBarry Smith 1628dbbe0bcdSBarry Smith if (pc->ops->presolve) PetscUseTypeMethod(pc, presolve, ksp, rhs, x); 1629f4f49eeaSPierre Jolivet else if (pc->presolve) PetscCall(pc->presolve(pc, ksp)); 16303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16314b9ad928SBarry Smith } 16324b9ad928SBarry Smith 1633f560b561SHong Zhang /*@C 1634f1580f4eSBarry Smith PCSetPreSolve - Sets function used by `PCPreSolve()` which is intended for any 1635f560b561SHong Zhang preconditioner-specific actions that must be performed before 1636f560b561SHong Zhang the iterative solve itself. 1637f560b561SHong Zhang 1638c3339decSBarry Smith Logically Collective 1639f560b561SHong Zhang 1640f560b561SHong Zhang Input Parameters: 1641f560b561SHong Zhang + pc - the preconditioner object 1642f560b561SHong Zhang - presolve - the function to call before the solve 1643f560b561SHong Zhang 164420f4b53cSBarry Smith Calling sequence of `presolve`: 164520f4b53cSBarry Smith + pc - the `PC` context 164620f4b53cSBarry Smith - ksp - the `KSP` context 1647f560b561SHong Zhang 1648f560b561SHong Zhang Level: developer 1649f560b561SHong Zhang 1650562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCPreSolve()` 1651f560b561SHong Zhang @*/ 165204c3f3b8SBarry Smith PetscErrorCode PCSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC pc, KSP ksp)) 1653d71ae5a4SJacob Faibussowitsch { 1654f560b561SHong Zhang PetscFunctionBegin; 1655f560b561SHong Zhang PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1656f560b561SHong Zhang pc->presolve = presolve; 16573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1658f560b561SHong Zhang } 1659f560b561SHong Zhang 16604b9ad928SBarry Smith /*@ 16614b9ad928SBarry Smith PCPostSolve - Optional post-solve phase, intended for any 16624b9ad928SBarry Smith preconditioner-specific actions that must be performed after 16634b9ad928SBarry Smith the iterative solve itself. 16644b9ad928SBarry Smith 1665c3339decSBarry Smith Collective 16664b9ad928SBarry Smith 16674b9ad928SBarry Smith Input Parameters: 16684b9ad928SBarry Smith + pc - the preconditioner context 16694b9ad928SBarry Smith - ksp - the Krylov subspace context 16704b9ad928SBarry Smith 1671feefa0e1SJacob Faibussowitsch Example Usage: 16724b9ad928SBarry Smith .vb 16734b9ad928SBarry Smith PCPreSolve(pc,ksp); 167423ce1328SBarry Smith KSPSolve(ksp,b,x); 16754b9ad928SBarry Smith PCPostSolve(pc,ksp); 16764b9ad928SBarry Smith .ve 16774b9ad928SBarry Smith 1678562efe2eSBarry Smith Level: developer 1679562efe2eSBarry Smith 16804b9ad928SBarry Smith Note: 1681f1580f4eSBarry Smith `KSPSolve()` calls this routine directly, so it is rarely called by the user. 16824b9ad928SBarry Smith 1683562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetPostSolve()`, `PCSetPresolve()`, `PCPreSolve()`, `KSPSolve()` 16844b9ad928SBarry Smith @*/ 1685d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPostSolve(PC pc, KSP ksp) 1686d71ae5a4SJacob Faibussowitsch { 16874b9ad928SBarry Smith Vec x, rhs; 16884b9ad928SBarry Smith 16894b9ad928SBarry Smith PetscFunctionBegin; 16900700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16910700a824SBarry Smith PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2); 1692d9a03883SBarry Smith pc->presolvedone--; 16939566063dSJacob Faibussowitsch PetscCall(KSPGetSolution(ksp, &x)); 16949566063dSJacob Faibussowitsch PetscCall(KSPGetRhs(ksp, &rhs)); 1695dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, postsolve, ksp, rhs, x); 16963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16974b9ad928SBarry Smith } 16984b9ad928SBarry Smith 1699ffeef943SBarry Smith /*@ 1700f1580f4eSBarry Smith PCLoad - Loads a `PC` that has been stored in binary with `PCView()`. 170155849f57SBarry Smith 1702c3339decSBarry Smith Collective 170355849f57SBarry Smith 170455849f57SBarry Smith Input Parameters: 1705f1580f4eSBarry Smith + newdm - the newly loaded `PC`, this needs to have been created with `PCCreate()` or 1706f1580f4eSBarry Smith some related function before a call to `PCLoad()`. 1707f1580f4eSBarry Smith - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()` 170855849f57SBarry Smith 170955849f57SBarry Smith Level: intermediate 171055849f57SBarry Smith 1711f1580f4eSBarry Smith Note: 1712562efe2eSBarry Smith The type is determined by the data in the file, any `PCType` set into the `PC` before this call is ignored. 171355849f57SBarry Smith 1714562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()` 171555849f57SBarry Smith @*/ 1716d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLoad(PC newdm, PetscViewer viewer) 1717d71ae5a4SJacob Faibussowitsch { 171855849f57SBarry Smith PetscBool isbinary; 1719060da220SMatthew G. Knepley PetscInt classid; 172055849f57SBarry Smith char type[256]; 172155849f57SBarry Smith 172255849f57SBarry Smith PetscFunctionBegin; 172355849f57SBarry Smith PetscValidHeaderSpecific(newdm, PC_CLASSID, 1); 172455849f57SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 17259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 172628b400f6SJacob Faibussowitsch PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 172755849f57SBarry Smith 17289566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT)); 17297827d75bSBarry Smith PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file"); 17309566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR)); 17319566063dSJacob Faibussowitsch PetscCall(PCSetType(newdm, type)); 1732dbbe0bcdSBarry Smith PetscTryTypeMethod(newdm, load, viewer); 17333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173455849f57SBarry Smith } 173555849f57SBarry Smith 17369804daf3SBarry Smith #include <petscdraw.h> 1737e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1738e04113cfSBarry Smith #include <petscviewersaws.h> 17390acecf5bSBarry Smith #endif 1740fe2efc57SMark 1741ffeef943SBarry Smith /*@ 1742562efe2eSBarry Smith PCViewFromOptions - View from the `PC` based on options in the options database 1743fe2efc57SMark 1744c3339decSBarry Smith Collective 1745fe2efc57SMark 1746fe2efc57SMark Input Parameters: 174720f4b53cSBarry Smith + A - the `PC` context 1748f1580f4eSBarry Smith . obj - Optional object that provides the options prefix 1749736c3998SJose E. Roman - name - command line option 1750fe2efc57SMark 1751fe2efc57SMark Level: intermediate 1752f1580f4eSBarry Smith 1753562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()` 1754fe2efc57SMark @*/ 1755d71ae5a4SJacob Faibussowitsch PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[]) 1756d71ae5a4SJacob Faibussowitsch { 1757fe2efc57SMark PetscFunctionBegin; 1758fe2efc57SMark PetscValidHeaderSpecific(A, PC_CLASSID, 1); 17599566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 17603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1761fe2efc57SMark } 1762fe2efc57SMark 1763ffeef943SBarry Smith /*@ 1764f1580f4eSBarry Smith PCView - Prints information about the `PC` 17654b9ad928SBarry Smith 1766c3339decSBarry Smith Collective 17674b9ad928SBarry Smith 17684b9ad928SBarry Smith Input Parameters: 1769feefa0e1SJacob Faibussowitsch + pc - the `PC` context 17704b9ad928SBarry Smith - viewer - optional visualization context 17714b9ad928SBarry Smith 177220f4b53cSBarry Smith Level: developer 177320f4b53cSBarry Smith 1774f1580f4eSBarry Smith Notes: 17754b9ad928SBarry Smith The available visualization contexts include 1776f1580f4eSBarry Smith + `PETSC_VIEWER_STDOUT_SELF` - standard output (default) 1777f1580f4eSBarry Smith - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard 17784b9ad928SBarry Smith output where only the first processor opens 17794b9ad928SBarry Smith the file. All other processors send their 17804b9ad928SBarry Smith data to the first processor to print. 17814b9ad928SBarry Smith 17824b9ad928SBarry Smith The user can open an alternative visualization contexts with 1783f1580f4eSBarry Smith `PetscViewerASCIIOpen()` (output to a specified file). 17844b9ad928SBarry Smith 1785562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewer`, `KSPView()`, `PetscViewerASCIIOpen()` 17864b9ad928SBarry Smith @*/ 1787d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView(PC pc, PetscViewer viewer) 1788d71ae5a4SJacob Faibussowitsch { 178919fd82e9SBarry Smith PCType cstr; 17906cd81132SPierre Jolivet PetscViewerFormat format; 17916cd81132SPierre Jolivet PetscBool iascii, isstring, isbinary, isdraw, pop = PETSC_FALSE; 1792e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1793536b137fSBarry Smith PetscBool issaws; 17940acecf5bSBarry Smith #endif 17954b9ad928SBarry Smith 17964b9ad928SBarry Smith PetscFunctionBegin; 17970700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 179848a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer)); 17990700a824SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1800c9780b6fSBarry Smith PetscCheckSameComm(pc, 1, viewer, 2); 18014b9ad928SBarry Smith 18029566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 18039566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 18049566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 18059566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 1806e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 18079566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws)); 18080acecf5bSBarry Smith #endif 1809219b1045SBarry Smith 181032077d6dSBarry Smith if (iascii) { 18119566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer)); 181248a46eb9SPierre Jolivet if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, " PC has not been set up so information may be incomplete\n")); 18139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1814dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 18159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1816834dbeb0SBarry Smith if (pc->mat) { 18176cd81132SPierre Jolivet PetscCall(PetscViewerGetFormat(viewer, &format)); 18186cd81132SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 18199566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 18206cd81132SPierre Jolivet pop = PETSC_TRUE; 18216cd81132SPierre Jolivet } 18224b9ad928SBarry Smith if (pc->pmat == pc->mat) { 18239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix = precond matrix:\n")); 18249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 18259566063dSJacob Faibussowitsch PetscCall(MatView(pc->mat, viewer)); 18269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 18274b9ad928SBarry Smith } else { 1828834dbeb0SBarry Smith if (pc->pmat) { 18299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix followed by preconditioner matrix:\n")); 18304b9ad928SBarry Smith } else { 18319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix:\n")); 18324b9ad928SBarry Smith } 18339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 18349566063dSJacob Faibussowitsch PetscCall(MatView(pc->mat, viewer)); 18359566063dSJacob Faibussowitsch if (pc->pmat) PetscCall(MatView(pc->pmat, viewer)); 18369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 18374b9ad928SBarry Smith } 18386cd81132SPierre Jolivet if (pop) PetscCall(PetscViewerPopFormat(viewer)); 18394b9ad928SBarry Smith } 18404b9ad928SBarry Smith } else if (isstring) { 18419566063dSJacob Faibussowitsch PetscCall(PCGetType(pc, &cstr)); 18429566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr)); 1843dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 18449566063dSJacob Faibussowitsch if (pc->mat) PetscCall(MatView(pc->mat, viewer)); 18459566063dSJacob Faibussowitsch if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer)); 18465b0b0462SBarry Smith } else if (isbinary) { 184755849f57SBarry Smith PetscInt classid = PC_FILE_CLASSID; 184855849f57SBarry Smith MPI_Comm comm; 184955849f57SBarry Smith PetscMPIInt rank; 185055849f57SBarry Smith char type[256]; 185155849f57SBarry Smith 18529566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 18539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1854dd400576SPatrick Sanan if (rank == 0) { 18559566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT)); 18569566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256)); 18579566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR)); 185855849f57SBarry Smith } 1859dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 1860219b1045SBarry Smith } else if (isdraw) { 1861219b1045SBarry Smith PetscDraw draw; 1862d9884438SBarry Smith char str[25]; 186389fd9fafSBarry Smith PetscReal x, y, bottom, h; 1864d9884438SBarry Smith PetscInt n; 1865219b1045SBarry Smith 18669566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 18679566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 18681d840656SPeter Brune if (pc->mat) { 18699566063dSJacob Faibussowitsch PetscCall(MatGetSize(pc->mat, &n, NULL)); 187063a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n)); 18711d840656SPeter Brune } else { 18729566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name)); 18731d840656SPeter Brune } 18749566063dSJacob Faibussowitsch PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h)); 187589fd9fafSBarry Smith bottom = y - h; 18769566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom)); 1877dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 18789566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 1879e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1880536b137fSBarry Smith } else if (issaws) { 1881d45a07a7SBarry Smith PetscMPIInt rank; 1882d45a07a7SBarry Smith 18839566063dSJacob Faibussowitsch PetscCall(PetscObjectName((PetscObject)pc)); 18849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 188548a46eb9SPierre Jolivet if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer)); 18869566063dSJacob Faibussowitsch if (pc->mat) PetscCall(MatView(pc->mat, viewer)); 18879566063dSJacob Faibussowitsch if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer)); 18880acecf5bSBarry Smith #endif 18894b9ad928SBarry Smith } 18903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18914b9ad928SBarry Smith } 18924b9ad928SBarry Smith 18934b9ad928SBarry Smith /*@C 1894562efe2eSBarry Smith PCRegister - Adds a method (`PCType`) to the preconditioner package. 18951c84c290SBarry Smith 1896cc4c1da9SBarry Smith Not collective. No Fortran Support 18971c84c290SBarry Smith 18981c84c290SBarry Smith Input Parameters: 189920f4b53cSBarry Smith + sname - name of a new user-defined solver 190020f4b53cSBarry Smith - function - routine to create method context 19011c84c290SBarry Smith 1902feefa0e1SJacob Faibussowitsch Example Usage: 19031c84c290SBarry Smith .vb 1904bdf89e91SBarry Smith PCRegister("my_solver", MySolverCreate); 19051c84c290SBarry Smith .ve 19061c84c290SBarry Smith 19071c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 19081c84c290SBarry Smith $ PCSetType(pc, "my_solver") 19091c84c290SBarry Smith or at runtime via the option 19101c84c290SBarry Smith $ -pc_type my_solver 19114b9ad928SBarry Smith 19124b9ad928SBarry Smith Level: advanced 19131c84c290SBarry Smith 191420f4b53cSBarry Smith Note: 191520f4b53cSBarry Smith `PCRegister()` may be called multiple times to add several user-defined preconditioners. 191620f4b53cSBarry Smith 1917562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCType`, `PCRegisterAll()` 19184b9ad928SBarry Smith @*/ 1919d71ae5a4SJacob Faibussowitsch PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC)) 1920d71ae5a4SJacob Faibussowitsch { 19214b9ad928SBarry Smith PetscFunctionBegin; 19229566063dSJacob Faibussowitsch PetscCall(PCInitializePackage()); 19239566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCList, sname, function)); 19243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19254b9ad928SBarry Smith } 19264b9ad928SBarry Smith 1927d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y) 1928d71ae5a4SJacob Faibussowitsch { 1929186a3e20SStefano Zampini PC pc; 1930186a3e20SStefano Zampini 1931186a3e20SStefano Zampini PetscFunctionBegin; 19329566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &pc)); 19339566063dSJacob Faibussowitsch PetscCall(PCApply(pc, X, Y)); 19343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1935186a3e20SStefano Zampini } 1936186a3e20SStefano Zampini 19375d83a8b1SBarry Smith /*@ 19380bacdadaSStefano Zampini PCComputeOperator - Computes the explicit preconditioned operator. 19394b9ad928SBarry Smith 1940c3339decSBarry Smith Collective 19414b9ad928SBarry Smith 1942d8d19677SJose E. Roman Input Parameters: 1943186a3e20SStefano Zampini + pc - the preconditioner object 1944562efe2eSBarry Smith - mattype - the `MatType` to be used for the operator 19454b9ad928SBarry Smith 19464b9ad928SBarry Smith Output Parameter: 1947a5b23f4aSJose E. Roman . mat - the explicit preconditioned operator 19484b9ad928SBarry Smith 194920f4b53cSBarry Smith Level: advanced 195020f4b53cSBarry Smith 1951f1580f4eSBarry Smith Note: 1952186a3e20SStefano Zampini This computation is done by applying the operators to columns of the identity matrix. 1953186a3e20SStefano Zampini This routine is costly in general, and is recommended for use only with relatively small systems. 1954562efe2eSBarry Smith Currently, this routine uses a dense matrix format when `mattype` == `NULL` 19554b9ad928SBarry Smith 1956562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPComputeOperator()`, `MatType` 19574b9ad928SBarry Smith @*/ 1958d71ae5a4SJacob Faibussowitsch PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat) 1959d71ae5a4SJacob Faibussowitsch { 1960186a3e20SStefano Zampini PetscInt N, M, m, n; 1961186a3e20SStefano Zampini Mat A, Apc; 19624b9ad928SBarry Smith 19634b9ad928SBarry Smith PetscFunctionBegin; 19640700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 19654f572ea9SToby Isaac PetscAssertPointer(mat, 3); 19669566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, &A, NULL)); 19679566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 19689566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 19699566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc)); 19709566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (void (*)(void))MatMult_PC)); 19719566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(Apc, mattype, mat)); 19729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Apc)); 19733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19744b9ad928SBarry Smith } 19754b9ad928SBarry Smith 19766c237d78SBarry Smith /*@ 19776c237d78SBarry Smith PCSetCoordinates - sets the coordinates of all the nodes on the local process 19786c237d78SBarry Smith 1979c3339decSBarry Smith Collective 19806c237d78SBarry Smith 19816c237d78SBarry Smith Input Parameters: 19826c237d78SBarry Smith + pc - the solver context 19836c237d78SBarry Smith . dim - the dimension of the coordinates 1, 2, or 3 198414893cbeSStefano Zampini . nloc - the blocked size of the coordinates array 198514893cbeSStefano Zampini - coords - the coordinates array 19866c237d78SBarry Smith 19876c237d78SBarry Smith Level: intermediate 19886c237d78SBarry Smith 1989f1580f4eSBarry Smith Note: 199020f4b53cSBarry Smith `coords` is an array of the dim coordinates for the nodes on 199120f4b53cSBarry Smith the local processor, of size `dim`*`nloc`. 199214893cbeSStefano Zampini If there are 108 equation on a processor 19936c237d78SBarry Smith for a displacement finite element discretization of elasticity (so 199414893cbeSStefano Zampini that there are nloc = 36 = 108/3 nodes) then the array must have 108 19956c237d78SBarry Smith double precision values (ie, 3 * 36). These x y z coordinates 19966c237d78SBarry Smith should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x, 19976c237d78SBarry Smith ... , N-1.z ]. 19986c237d78SBarry Smith 1999562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `MatSetNearNullSpace()` 20006c237d78SBarry Smith @*/ 2001d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[]) 2002d71ae5a4SJacob Faibussowitsch { 20036c237d78SBarry Smith PetscFunctionBegin; 200432954da3SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 200532954da3SStefano Zampini PetscValidLogicalCollectiveInt(pc, dim, 2); 200622794d57SStefano Zampini PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords)); 20073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20086c237d78SBarry Smith } 2009fd2dd295SFande Kong 2010fd2dd295SFande Kong /*@ 2011fd2dd295SFande Kong PCGetInterpolations - Gets interpolation matrices for all levels (except level 0) 2012fd2dd295SFande Kong 2013c3339decSBarry Smith Logically Collective 2014fd2dd295SFande Kong 2015d8d19677SJose E. Roman Input Parameter: 2016d8d19677SJose E. Roman . pc - the precondition context 2017fd2dd295SFande Kong 2018d8d19677SJose E. Roman Output Parameters: 2019d8d19677SJose E. Roman + num_levels - the number of levels 2020562efe2eSBarry Smith - interpolations - the interpolation matrices (size of `num_levels`-1) 2021fd2dd295SFande Kong 2022fd2dd295SFande Kong Level: advanced 2023fd2dd295SFande Kong 2024562efe2eSBarry Smith Developer Note: 2025f1580f4eSBarry Smith Why is this here instead of in `PCMG` etc? 2026fd2dd295SFande Kong 2027562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()` 2028fd2dd295SFande Kong @*/ 2029d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[]) 2030d71ae5a4SJacob Faibussowitsch { 2031fd2dd295SFande Kong PetscFunctionBegin; 2032fd2dd295SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 20334f572ea9SToby Isaac PetscAssertPointer(num_levels, 2); 20344f572ea9SToby Isaac PetscAssertPointer(interpolations, 3); 2035cac4c232SBarry Smith PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations)); 20363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2037fd2dd295SFande Kong } 2038fd2dd295SFande Kong 2039fd2dd295SFande Kong /*@ 2040fd2dd295SFande Kong PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level) 2041fd2dd295SFande Kong 2042c3339decSBarry Smith Logically Collective 2043fd2dd295SFande Kong 2044d8d19677SJose E. Roman Input Parameter: 2045d8d19677SJose E. Roman . pc - the precondition context 2046fd2dd295SFande Kong 2047d8d19677SJose E. Roman Output Parameters: 2048d8d19677SJose E. Roman + num_levels - the number of levels 2049562efe2eSBarry Smith - coarseOperators - the coarse operator matrices (size of `num_levels`-1) 2050fd2dd295SFande Kong 2051fd2dd295SFande Kong Level: advanced 2052fd2dd295SFande Kong 2053562efe2eSBarry Smith Developer Note: 2054f1580f4eSBarry Smith Why is this here instead of in `PCMG` etc? 2055fd2dd295SFande Kong 2056562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()` 2057fd2dd295SFande Kong @*/ 2058d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[]) 2059d71ae5a4SJacob Faibussowitsch { 2060fd2dd295SFande Kong PetscFunctionBegin; 2061fd2dd295SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 20624f572ea9SToby Isaac PetscAssertPointer(num_levels, 2); 20634f572ea9SToby Isaac PetscAssertPointer(coarseOperators, 3); 2064cac4c232SBarry Smith PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators)); 20653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2066fd2dd295SFande Kong } 2067