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 55fe57bb1aSStefano Zampini /* do not log solves, setup and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */ 56fe57bb1aSStefano Zampini PETSC_EXTERN PetscLogEvent KSP_Solve, KSP_SetUp; 57fe57bb1aSStefano Zampini 58fe57bb1aSStefano Zampini static PetscErrorCode PCLogEventsDeactivatePush(void) 59fe57bb1aSStefano Zampini { 60fe57bb1aSStefano Zampini PetscFunctionBegin; 61fe57bb1aSStefano Zampini PetscCall(KSPInitializePackage()); 62fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePush(KSP_Solve)); 63fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePush(KSP_SetUp)); 64fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePush(PC_Apply)); 65fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePush(PC_SetUp)); 66fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePush(PC_SetUpOnBlocks)); 67fe57bb1aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 68fe57bb1aSStefano Zampini } 69fe57bb1aSStefano Zampini 70fe57bb1aSStefano Zampini static PetscErrorCode PCLogEventsDeactivatePop(void) 71fe57bb1aSStefano Zampini { 72fe57bb1aSStefano Zampini PetscFunctionBegin; 73fe57bb1aSStefano Zampini PetscCall(KSPInitializePackage()); 74fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePop(KSP_Solve)); 75fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePop(KSP_SetUp)); 76fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePop(PC_Apply)); 77fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePop(PC_SetUp)); 78fe57bb1aSStefano Zampini PetscCall(PetscLogEventDeactivatePop(PC_SetUpOnBlocks)); 79fe57bb1aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 80fe57bb1aSStefano Zampini } 81fe57bb1aSStefano 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 958f1580f4eSBarry Smith Not Collective 9591b2b9847SBarry Smith 9601b2b9847SBarry Smith Input Parameter: 9611b2b9847SBarry Smith . pc - the preconditioner context 9621b2b9847SBarry Smith 9631b2b9847SBarry Smith Output Parameter: 9641b2b9847SBarry Smith . reason - the reason it failed 9651b2b9847SBarry Smith 96620f4b53cSBarry Smith Level: advanced 96720f4b53cSBarry Smith 968f1580f4eSBarry Smith Note: 9699a6c2652SBarry Smith After call `KSPCheckDot()` or `KSPCheckNorm()` inside a `KSPSolve()` or a call to `PCReduceFailedReason()` 9709a6c2652SBarry Smith this is the maximum over reason over all ranks in the `PC` communicator and hence logically collective. 9719a6c2652SBarry Smith Otherwise it returns the local value. 9721b2b9847SBarry Smith 9739a6c2652SBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCSetFailedReason()` 9741b2b9847SBarry Smith @*/ 9759a6c2652SBarry Smith PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason) 976d71ae5a4SJacob Faibussowitsch { 9771b2b9847SBarry Smith PetscFunctionBegin; 9786479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9791b2b9847SBarry Smith if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled; 9801b2b9847SBarry Smith else *reason = pc->failedreason; 9813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9821b2b9847SBarry Smith } 983422a814eSBarry Smith 9846479e835SStefano Zampini /*@ 9856479e835SStefano Zampini PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC` 9866479e835SStefano Zampini 9876479e835SStefano Zampini Collective 9886479e835SStefano Zampini 9896479e835SStefano Zampini Input Parameter: 9906479e835SStefano Zampini . pc - the preconditioner context 9916479e835SStefano Zampini 9926479e835SStefano Zampini Level: advanced 9936479e835SStefano Zampini 9946479e835SStefano Zampini Note: 995562efe2eSBarry Smith Different MPI processes may have different reasons or no reason, see `PCGetFailedReason()`. This routine 9966479e835SStefano Zampini makes them have a common value (failure if any MPI process had a failure). 9976479e835SStefano Zampini 998562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()` 9996479e835SStefano Zampini @*/ 10006479e835SStefano Zampini PetscErrorCode PCReduceFailedReason(PC pc) 10016479e835SStefano Zampini { 10026479e835SStefano Zampini PetscInt buf; 10036479e835SStefano Zampini 10046479e835SStefano Zampini PetscFunctionBegin; 10056479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 10066479e835SStefano Zampini buf = (PetscInt)pc->failedreason; 1007462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc))); 10086479e835SStefano Zampini pc->failedreason = (PCFailedReason)buf; 10096479e835SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 10106479e835SStefano Zampini } 10116479e835SStefano Zampini 10124b9ad928SBarry Smith /* 10134b9ad928SBarry Smith a setupcall of 0 indicates never setup, 101423ee1639SBarry Smith 1 indicates has been previously setup 1015422a814eSBarry Smith -1 indicates a PCSetUp() was attempted and failed 10164b9ad928SBarry Smith */ 10174b9ad928SBarry Smith /*@ 10184b9ad928SBarry Smith PCSetUp - Prepares for the use of a preconditioner. 10194b9ad928SBarry Smith 1020c3339decSBarry Smith Collective 10214b9ad928SBarry Smith 10224b9ad928SBarry Smith Input Parameter: 10234b9ad928SBarry Smith . pc - the preconditioner context 10244b9ad928SBarry Smith 10254b9ad928SBarry Smith Level: developer 10264b9ad928SBarry Smith 1027562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()` 10284b9ad928SBarry Smith @*/ 1029d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp(PC pc) 1030d71ae5a4SJacob Faibussowitsch { 1031566e8bf2SBarry Smith const char *def; 1032fc9b51d3SBarry Smith PetscObjectState matstate, matnonzerostate; 10334b9ad928SBarry Smith 10344b9ad928SBarry Smith PetscFunctionBegin; 10350700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 103628b400f6SJacob Faibussowitsch PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first"); 10376ce5e81cSLisandro Dalcin 103823ee1639SBarry Smith if (pc->setupcalled && pc->reusepreconditioner) { 10399566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n")); 10403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10414b9ad928SBarry Smith } 10424b9ad928SBarry Smith 10439566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate)); 10449566063dSJacob Faibussowitsch PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate)); 104523ee1639SBarry Smith if (!pc->setupcalled) { 10465e62d202SMark Adams //PetscCall(PetscInfo(pc, "Setting up PC for first time\n")); 104723ee1639SBarry Smith pc->flag = DIFFERENT_NONZERO_PATTERN; 10489df67409SStefano Zampini } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS); 10499df67409SStefano Zampini else { 10509df67409SStefano Zampini if (matnonzerostate != pc->matnonzerostate) { 10519566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n")); 105223ee1639SBarry Smith pc->flag = DIFFERENT_NONZERO_PATTERN; 105323ee1639SBarry Smith } else { 10545e62d202SMark Adams //PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n")); 105523ee1639SBarry Smith pc->flag = SAME_NONZERO_PATTERN; 105623ee1639SBarry Smith } 105723ee1639SBarry Smith } 105823ee1639SBarry Smith pc->matstate = matstate; 1059fc9b51d3SBarry Smith pc->matnonzerostate = matnonzerostate; 106023ee1639SBarry Smith 10617adad957SLisandro Dalcin if (!((PetscObject)pc)->type_name) { 10629566063dSJacob Faibussowitsch PetscCall(PCGetDefaultType_Private(pc, &def)); 10639566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, def)); 1064566e8bf2SBarry Smith } 10654b9ad928SBarry Smith 10669566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure)); 10679566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure)); 10689566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0)); 10694b9ad928SBarry Smith if (pc->ops->setup) { 1070fe57bb1aSStefano Zampini PetscCall(PCLogEventsDeactivatePush()); 1071dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, setup); 1072fe57bb1aSStefano Zampini PetscCall(PCLogEventsDeactivatePop()); 10734b9ad928SBarry Smith } 10749566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0)); 1075422a814eSBarry Smith if (!pc->setupcalled) pc->setupcalled = 1; 10763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10774b9ad928SBarry Smith } 10784b9ad928SBarry Smith 10794b9ad928SBarry Smith /*@ 10804b9ad928SBarry Smith PCSetUpOnBlocks - Sets up the preconditioner for each block in 108173716367SStefano Zampini the block Jacobi, overlapping Schwarz, and fieldsplit methods. 10824b9ad928SBarry Smith 1083c3339decSBarry Smith Collective 10844b9ad928SBarry Smith 1085f1580f4eSBarry Smith Input Parameter: 10864b9ad928SBarry Smith . pc - the preconditioner context 10874b9ad928SBarry Smith 10884b9ad928SBarry Smith Level: developer 10894b9ad928SBarry Smith 109073716367SStefano Zampini Notes: 109173716367SStefano Zampini For nested preconditioners such as `PCBJACOBI`, `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is 1092f1580f4eSBarry Smith called on the outer `PC`, this routine ensures it is called. 1093f1580f4eSBarry Smith 109473716367SStefano Zampini It calls `PCSetUp()` if not yet called. 109573716367SStefano Zampini 1096562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()` 10974b9ad928SBarry Smith @*/ 1098d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUpOnBlocks(PC pc) 1099d71ae5a4SJacob Faibussowitsch { 11004b9ad928SBarry Smith PetscFunctionBegin; 11010700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 110273716367SStefano Zampini if (!pc->setupcalled) PetscCall(PCSetUp(pc)); /* "if" to prevent -info extra prints */ 11033ba16761SJacob Faibussowitsch if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS); 11049566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0)); 1105fe57bb1aSStefano Zampini PetscCall(PCLogEventsDeactivatePush()); 1106dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, setuponblocks); 1107fe57bb1aSStefano Zampini PetscCall(PCLogEventsDeactivatePop()); 11089566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0)); 11093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11104b9ad928SBarry Smith } 11114b9ad928SBarry Smith 11124b9ad928SBarry Smith /*@C 11134b9ad928SBarry Smith PCSetModifySubMatrices - Sets a user-defined routine for modifying the 111404c3f3b8SBarry Smith submatrices that arise within certain subdomain-based preconditioners such as `PCASM` 11154b9ad928SBarry Smith 1116c3339decSBarry Smith Logically Collective 11174b9ad928SBarry Smith 11184b9ad928SBarry Smith Input Parameters: 11194b9ad928SBarry Smith + pc - the preconditioner context 11204b9ad928SBarry Smith . func - routine for modifying the submatrices 112104c3f3b8SBarry Smith - ctx - optional user-defined context (may be `NULL`) 11224b9ad928SBarry Smith 112320f4b53cSBarry Smith Calling sequence of `func`: 112420f4b53cSBarry Smith + pc - the preconditioner context 112504c3f3b8SBarry Smith . nsub - number of index sets 112620f4b53cSBarry Smith . row - an array of index sets that contain the global row numbers 11274b9ad928SBarry Smith that comprise each local submatrix 11284b9ad928SBarry Smith . col - an array of index sets that contain the global column numbers 11294b9ad928SBarry Smith that comprise each local submatrix 11304b9ad928SBarry Smith . submat - array of local submatrices 11314b9ad928SBarry Smith - ctx - optional user-defined context for private data for the 113204c3f3b8SBarry Smith user-defined func routine (may be `NULL`) 11334b9ad928SBarry Smith 113420f4b53cSBarry Smith Level: advanced 113520f4b53cSBarry Smith 11364b9ad928SBarry Smith Notes: 113704c3f3b8SBarry Smith The basic submatrices are extracted from the preconditioner matrix as 113804c3f3b8SBarry Smith usual; the user can then alter these (for example, to set different boundary 113904c3f3b8SBarry Smith conditions for each submatrix) before they are used for the local solves. 114004c3f3b8SBarry Smith 1141f1580f4eSBarry Smith `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and 1142f1580f4eSBarry Smith `KSPSolve()`. 11434b9ad928SBarry Smith 1144f1580f4eSBarry Smith A routine set by `PCSetModifySubMatrices()` is currently called within 1145f1580f4eSBarry Smith the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`) 11464b9ad928SBarry Smith preconditioners. All other preconditioners ignore this routine. 11474b9ad928SBarry Smith 1148562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()` 11494b9ad928SBarry Smith @*/ 115004c3f3b8SBarry Smith PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx), void *ctx) 1151d71ae5a4SJacob Faibussowitsch { 11524b9ad928SBarry Smith PetscFunctionBegin; 11530700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 11544b9ad928SBarry Smith pc->modifysubmatrices = func; 11554b9ad928SBarry Smith pc->modifysubmatricesP = ctx; 11563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11574b9ad928SBarry Smith } 11584b9ad928SBarry Smith 11594b9ad928SBarry Smith /*@C 11604b9ad928SBarry Smith PCModifySubMatrices - Calls an optional user-defined routine within 1161f1580f4eSBarry Smith certain preconditioners if one has been set with `PCSetModifySubMatrices()`. 11624b9ad928SBarry Smith 1163c3339decSBarry Smith Collective 11644b9ad928SBarry Smith 11654b9ad928SBarry Smith Input Parameters: 11664b9ad928SBarry Smith + pc - the preconditioner context 11674b9ad928SBarry Smith . nsub - the number of local submatrices 11684b9ad928SBarry Smith . row - an array of index sets that contain the global row numbers 11694b9ad928SBarry Smith that comprise each local submatrix 11704b9ad928SBarry Smith . col - an array of index sets that contain the global column numbers 11714b9ad928SBarry Smith that comprise each local submatrix 11724b9ad928SBarry Smith . submat - array of local submatrices 11734b9ad928SBarry Smith - ctx - optional user-defined context for private data for the 1174562efe2eSBarry Smith user-defined routine (may be `NULL`) 11754b9ad928SBarry Smith 11764b9ad928SBarry Smith Output Parameter: 11774b9ad928SBarry Smith . submat - array of local submatrices (the entries of which may 11784b9ad928SBarry Smith have been modified) 11794b9ad928SBarry Smith 118020f4b53cSBarry Smith Level: developer 118120f4b53cSBarry Smith 118204c3f3b8SBarry Smith Note: 11834b9ad928SBarry Smith The user should NOT generally call this routine, as it will 118404c3f3b8SBarry Smith automatically be called within certain preconditioners. 11854b9ad928SBarry Smith 1186562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetModifySubMatrices()` 11874b9ad928SBarry Smith @*/ 1188d71ae5a4SJacob Faibussowitsch PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx) 1189d71ae5a4SJacob Faibussowitsch { 11904b9ad928SBarry Smith PetscFunctionBegin; 11910700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 11923ba16761SJacob Faibussowitsch if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS); 11939566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0)); 11949566063dSJacob Faibussowitsch PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx)); 11959566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0)); 11963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11974b9ad928SBarry Smith } 11984b9ad928SBarry Smith 11994b9ad928SBarry Smith /*@ 12004b9ad928SBarry Smith PCSetOperators - Sets the matrix associated with the linear system and 12014b9ad928SBarry Smith a (possibly) different one associated with the preconditioner. 12024b9ad928SBarry Smith 1203c3339decSBarry Smith Logically Collective 12044b9ad928SBarry Smith 12054b9ad928SBarry Smith Input Parameters: 12064b9ad928SBarry Smith + pc - the preconditioner context 1207e5d3d808SBarry Smith . Amat - the matrix that defines the linear system 1208d1e9a80fSBarry Smith - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat. 12094b9ad928SBarry Smith 121020f4b53cSBarry Smith Level: intermediate 1211189c0b8aSBarry Smith 121220f4b53cSBarry Smith Notes: 121320f4b53cSBarry Smith Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used. 121420f4b53cSBarry Smith 121520f4b53cSBarry Smith If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then 1216f1580f4eSBarry Smith first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()` 1217f1580f4eSBarry Smith on it and then pass it back in in your call to `KSPSetOperators()`. 1218189c0b8aSBarry Smith 12194b9ad928SBarry Smith More Notes about Repeated Solution of Linear Systems: 122020f4b53cSBarry Smith PETSc does NOT reset the matrix entries of either `Amat` or `Pmat` 12214b9ad928SBarry Smith to zero after a linear solve; the user is completely responsible for 1222f1580f4eSBarry Smith matrix assembly. See the routine `MatZeroEntries()` if desiring to 12234b9ad928SBarry Smith zero all elements of a matrix. 12244b9ad928SBarry Smith 1225562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()` 12264b9ad928SBarry Smith @*/ 1227d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat) 1228d71ae5a4SJacob Faibussowitsch { 12293b3f6333SBarry Smith PetscInt m1, n1, m2, n2; 12304b9ad928SBarry Smith 12314b9ad928SBarry Smith PetscFunctionBegin; 12320700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 12330700a824SBarry Smith if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 12340700a824SBarry Smith if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3); 12353fc8bf9cSBarry Smith if (Amat) PetscCheckSameComm(pc, 1, Amat, 2); 12363fc8bf9cSBarry Smith if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3); 123731641f1aSBarry Smith if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) { 12389566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Amat, &m1, &n1)); 12399566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->mat, &m2, &n2)); 124063a3b9bcSJacob 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); 12419566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Pmat, &m1, &n1)); 12429566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2)); 124363a3b9bcSJacob 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); 12443b3f6333SBarry Smith } 12454b9ad928SBarry Smith 124623ee1639SBarry Smith if (Pmat != pc->pmat) { 124723ee1639SBarry Smith /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */ 124823ee1639SBarry Smith pc->matnonzerostate = -1; 124923ee1639SBarry Smith pc->matstate = -1; 125023ee1639SBarry Smith } 125123ee1639SBarry Smith 1252906ed7ccSBarry Smith /* reference first in case the matrices are the same */ 12539566063dSJacob Faibussowitsch if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat)); 12549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->mat)); 12559566063dSJacob Faibussowitsch if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat)); 12569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->pmat)); 12574b9ad928SBarry Smith pc->mat = Amat; 12584b9ad928SBarry Smith pc->pmat = Pmat; 12593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1260e56f5c9eSBarry Smith } 1261e56f5c9eSBarry Smith 126223ee1639SBarry Smith /*@ 126323ee1639SBarry Smith PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed. 126423ee1639SBarry Smith 1265c3339decSBarry Smith Logically Collective 126623ee1639SBarry Smith 126723ee1639SBarry Smith Input Parameters: 126823ee1639SBarry Smith + pc - the preconditioner context 1269f1580f4eSBarry Smith - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner 127023ee1639SBarry Smith 12712b26979fSBarry Smith Level: intermediate 12722b26979fSBarry Smith 1273f1580f4eSBarry Smith Note: 1274f1580f4eSBarry Smith Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option 1275f1580f4eSBarry Smith prevents this. 1276f1580f4eSBarry Smith 1277562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()` 127823ee1639SBarry Smith @*/ 1279d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag) 1280d71ae5a4SJacob Faibussowitsch { 128123ee1639SBarry Smith PetscFunctionBegin; 128223ee1639SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1283f9177522SStefano Zampini PetscValidLogicalCollectiveBool(pc, flag, 2); 128423ee1639SBarry Smith pc->reusepreconditioner = flag; 12853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12864b9ad928SBarry Smith } 12874b9ad928SBarry Smith 1288c60c7ad4SBarry Smith /*@ 1289f1580f4eSBarry Smith PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed. 1290c60c7ad4SBarry Smith 1291bba28a21SBarry Smith Not Collective 1292c60c7ad4SBarry Smith 1293c60c7ad4SBarry Smith Input Parameter: 1294c60c7ad4SBarry Smith . pc - the preconditioner context 1295c60c7ad4SBarry Smith 1296c60c7ad4SBarry Smith Output Parameter: 1297f1580f4eSBarry Smith . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner 1298c60c7ad4SBarry Smith 1299d0418729SSatish Balay Level: intermediate 1300d0418729SSatish Balay 1301562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()` 1302c60c7ad4SBarry Smith @*/ 1303d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag) 1304d71ae5a4SJacob Faibussowitsch { 1305c60c7ad4SBarry Smith PetscFunctionBegin; 1306c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 13074f572ea9SToby Isaac PetscAssertPointer(flag, 2); 1308c60c7ad4SBarry Smith *flag = pc->reusepreconditioner; 13093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1310c60c7ad4SBarry Smith } 1311c60c7ad4SBarry Smith 1312487a658cSBarry Smith /*@ 13134b9ad928SBarry Smith PCGetOperators - Gets the matrix associated with the linear system and 13144b9ad928SBarry Smith possibly a different one associated with the preconditioner. 13154b9ad928SBarry Smith 1316562efe2eSBarry Smith Not Collective, though parallel `Mat`s are returned if `pc` is parallel 13174b9ad928SBarry Smith 13184b9ad928SBarry Smith Input Parameter: 13194b9ad928SBarry Smith . pc - the preconditioner context 13204b9ad928SBarry Smith 13214b9ad928SBarry Smith Output Parameters: 1322e5d3d808SBarry Smith + Amat - the matrix defining the linear system 132323ee1639SBarry Smith - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat. 13244b9ad928SBarry Smith 13254b9ad928SBarry Smith Level: intermediate 13264b9ad928SBarry Smith 1327f1580f4eSBarry Smith Note: 132895452b02SPatrick Sanan Does not increase the reference count of the matrices, so you should not destroy them 1329298cc208SBarry Smith 1330f1580f4eSBarry Smith Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators 1331f1580f4eSBarry Smith are created in `PC` and returned to the user. In this case, if both operators 133273950996SBarry Smith mat and pmat are requested, two DIFFERENT operators will be returned. If 133373950996SBarry Smith only one is requested both operators in the PC will be the same (i.e. as 1334f1580f4eSBarry Smith if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats). 133573950996SBarry Smith The user must set the sizes of the returned matrices and their type etc just 1336f1580f4eSBarry Smith as if the user created them with `MatCreate()`. For example, 133773950996SBarry Smith 1338f1580f4eSBarry Smith .vb 1339f1580f4eSBarry Smith KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to 1340f1580f4eSBarry Smith set size, type, etc of Amat 134173950996SBarry Smith 1342f1580f4eSBarry Smith MatCreate(comm,&mat); 1343f1580f4eSBarry Smith KSP/PCSetOperators(ksp/pc,Amat,Amat); 1344f1580f4eSBarry Smith PetscObjectDereference((PetscObject)mat); 1345f1580f4eSBarry Smith set size, type, etc of Amat 1346f1580f4eSBarry Smith .ve 134773950996SBarry Smith 134873950996SBarry Smith and 134973950996SBarry Smith 1350f1580f4eSBarry Smith .vb 1351f1580f4eSBarry Smith KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to 1352f1580f4eSBarry Smith set size, type, etc of Amat and Pmat 135373950996SBarry Smith 1354f1580f4eSBarry Smith MatCreate(comm,&Amat); 1355f1580f4eSBarry Smith MatCreate(comm,&Pmat); 1356f1580f4eSBarry Smith KSP/PCSetOperators(ksp/pc,Amat,Pmat); 1357f1580f4eSBarry Smith PetscObjectDereference((PetscObject)Amat); 1358f1580f4eSBarry Smith PetscObjectDereference((PetscObject)Pmat); 1359f1580f4eSBarry Smith set size, type, etc of Amat and Pmat 1360f1580f4eSBarry Smith .ve 136173950996SBarry Smith 1362f1580f4eSBarry Smith The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy 1363b8d709abSRichard Tran Mills of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely 1364f1580f4eSBarry Smith managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look 1365f1580f4eSBarry Smith at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to 1366f1580f4eSBarry Smith the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP 1367f1580f4eSBarry Smith you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you). 1368f1580f4eSBarry Smith Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when 136973950996SBarry Smith it can be created for you? 137073950996SBarry Smith 1371562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()` 13724b9ad928SBarry Smith @*/ 1373d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat) 1374d71ae5a4SJacob Faibussowitsch { 13754b9ad928SBarry Smith PetscFunctionBegin; 13760700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1377e5d3d808SBarry Smith if (Amat) { 137873950996SBarry Smith if (!pc->mat) { 13799fca8c71SStefano Zampini if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */ 13809a4708feSJed Brown pc->mat = pc->pmat; 13819566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->mat)); 1382e5d3d808SBarry Smith } else { /* both Amat and Pmat are empty */ 13839566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat)); 1384e5d3d808SBarry Smith if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */ 138573950996SBarry Smith pc->pmat = pc->mat; 13869566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->pmat)); 138773950996SBarry Smith } 138873950996SBarry Smith } 13899a4708feSJed Brown } 1390e5d3d808SBarry Smith *Amat = pc->mat; 139173950996SBarry Smith } 1392e5d3d808SBarry Smith if (Pmat) { 139373950996SBarry Smith if (!pc->pmat) { 1394e5d3d808SBarry Smith if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */ 13959a4708feSJed Brown pc->pmat = pc->mat; 13969566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->pmat)); 13979a4708feSJed Brown } else { 13989566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat)); 1399e5d3d808SBarry Smith if (!Amat) { /* user did NOT request Amat, so make same as Pmat */ 140073950996SBarry Smith pc->mat = pc->pmat; 14019566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->mat)); 140273950996SBarry Smith } 140373950996SBarry Smith } 14049a4708feSJed Brown } 1405e5d3d808SBarry Smith *Pmat = pc->pmat; 140673950996SBarry Smith } 14073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14084b9ad928SBarry Smith } 14094b9ad928SBarry Smith 14105d83a8b1SBarry Smith /*@ 1411906ed7ccSBarry Smith PCGetOperatorsSet - Determines if the matrix associated with the linear system and 1412f1580f4eSBarry Smith possibly a different one associated with the preconditioner have been set in the `PC`. 1413906ed7ccSBarry Smith 141420f4b53cSBarry Smith Not Collective, though the results on all processes should be the same 1415906ed7ccSBarry Smith 1416906ed7ccSBarry Smith Input Parameter: 1417906ed7ccSBarry Smith . pc - the preconditioner context 1418906ed7ccSBarry Smith 1419906ed7ccSBarry Smith Output Parameters: 1420906ed7ccSBarry Smith + mat - the matrix associated with the linear system was set 1421906ed7ccSBarry Smith - pmat - matrix associated with the preconditioner was set, usually the same 1422906ed7ccSBarry Smith 1423906ed7ccSBarry Smith Level: intermediate 1424906ed7ccSBarry Smith 1425562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()` 1426906ed7ccSBarry Smith @*/ 1427d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat) 1428d71ae5a4SJacob Faibussowitsch { 1429906ed7ccSBarry Smith PetscFunctionBegin; 14300700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1431906ed7ccSBarry Smith if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE; 1432906ed7ccSBarry Smith if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE; 14333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1434906ed7ccSBarry Smith } 1435906ed7ccSBarry Smith 1436f39d8e23SSatish Balay /*@ 1437a4fd02acSBarry Smith PCFactorGetMatrix - Gets the factored matrix from the 1438f1580f4eSBarry Smith preconditioner context. This routine is valid only for the `PCLU`, 1439f1580f4eSBarry Smith `PCILU`, `PCCHOLESKY`, and `PCICC` methods. 14404b9ad928SBarry Smith 1441562efe2eSBarry Smith Not Collective though `mat` is parallel if `pc` is parallel 14424b9ad928SBarry Smith 1443f1580f4eSBarry Smith Input Parameter: 14444b9ad928SBarry Smith . pc - the preconditioner context 14454b9ad928SBarry Smith 1446feefa0e1SJacob Faibussowitsch Output Parameters: 14474b9ad928SBarry Smith . mat - the factored matrix 14484b9ad928SBarry Smith 14494b9ad928SBarry Smith Level: advanced 14504b9ad928SBarry Smith 1451f1580f4eSBarry Smith Note: 1452562efe2eSBarry Smith Does not increase the reference count for `mat` so DO NOT destroy it 14539405f653SBarry Smith 1454562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC` 14554b9ad928SBarry Smith @*/ 1456d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat) 1457d71ae5a4SJacob Faibussowitsch { 14584b9ad928SBarry Smith PetscFunctionBegin; 14590700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 14604f572ea9SToby Isaac PetscAssertPointer(mat, 2); 14617def7855SStefano Zampini PetscCall(PCFactorSetUpMatSolverType(pc)); 1462dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, getfactoredmatrix, mat); 14633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14644b9ad928SBarry Smith } 14654b9ad928SBarry Smith 1466cc4c1da9SBarry Smith /*@ 14674b9ad928SBarry Smith PCSetOptionsPrefix - Sets the prefix used for searching for all 1468f1580f4eSBarry Smith `PC` options in the database. 14694b9ad928SBarry Smith 1470c3339decSBarry Smith Logically Collective 14714b9ad928SBarry Smith 14724b9ad928SBarry Smith Input Parameters: 14734b9ad928SBarry Smith + pc - the preconditioner context 1474f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests 14754b9ad928SBarry Smith 1476f1580f4eSBarry Smith Note: 14774b9ad928SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 14784b9ad928SBarry Smith The first character of all runtime options is AUTOMATICALLY the 14794b9ad928SBarry Smith hyphen. 14804b9ad928SBarry Smith 14814b9ad928SBarry Smith Level: advanced 14824b9ad928SBarry Smith 1483562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()` 14844b9ad928SBarry Smith @*/ 1485d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[]) 1486d71ae5a4SJacob Faibussowitsch { 14874b9ad928SBarry Smith PetscFunctionBegin; 14880700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 14899566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix)); 14903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14914b9ad928SBarry Smith } 14924b9ad928SBarry Smith 1493cc4c1da9SBarry Smith /*@ 14944b9ad928SBarry Smith PCAppendOptionsPrefix - Appends to the prefix used for searching for all 1495f1580f4eSBarry Smith `PC` options in the database. 14964b9ad928SBarry Smith 1497c3339decSBarry Smith Logically Collective 14984b9ad928SBarry Smith 14994b9ad928SBarry Smith Input Parameters: 15004b9ad928SBarry Smith + pc - the preconditioner context 1501f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests 15024b9ad928SBarry Smith 1503f1580f4eSBarry Smith Note: 15044b9ad928SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 15054b9ad928SBarry Smith The first character of all runtime options is AUTOMATICALLY the 15064b9ad928SBarry Smith hyphen. 15074b9ad928SBarry Smith 15084b9ad928SBarry Smith Level: advanced 15094b9ad928SBarry Smith 1510562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()` 15114b9ad928SBarry Smith @*/ 1512d71ae5a4SJacob Faibussowitsch PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[]) 1513d71ae5a4SJacob Faibussowitsch { 15144b9ad928SBarry Smith PetscFunctionBegin; 15150700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15169566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix)); 15173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15184b9ad928SBarry Smith } 15194b9ad928SBarry Smith 1520cc4c1da9SBarry Smith /*@ 15214b9ad928SBarry Smith PCGetOptionsPrefix - Gets the prefix used for searching for all 15224b9ad928SBarry Smith PC options in the database. 15234b9ad928SBarry Smith 15244b9ad928SBarry Smith Not Collective 15254b9ad928SBarry Smith 1526f1580f4eSBarry Smith Input Parameter: 15274b9ad928SBarry Smith . pc - the preconditioner context 15284b9ad928SBarry Smith 1529f1580f4eSBarry Smith Output Parameter: 15304b9ad928SBarry Smith . prefix - pointer to the prefix string used, is returned 15314b9ad928SBarry Smith 15324b9ad928SBarry Smith Level: advanced 15334b9ad928SBarry Smith 1534562efe2eSBarry Smith Fortran Note: 1535562efe2eSBarry Smith The user should pass in a string `prefix` of 1536562efe2eSBarry Smith sufficient length to hold the prefix. 1537562efe2eSBarry Smith 1538562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()` 15394b9ad928SBarry Smith @*/ 1540d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[]) 1541d71ae5a4SJacob Faibussowitsch { 15424b9ad928SBarry Smith PetscFunctionBegin; 15430700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15444f572ea9SToby Isaac PetscAssertPointer(prefix, 2); 15459566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix)); 15463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15474b9ad928SBarry Smith } 15484b9ad928SBarry Smith 15498066bbecSBarry Smith /* 1550dd8e379bSPierre Jolivet Indicates the right-hand side will be changed by KSPSolve(), this occurs for a few 15518066bbecSBarry Smith preconditioners including BDDC and Eisentat that transform the equations before applying 15528066bbecSBarry Smith the Krylov methods 15538066bbecSBarry Smith */ 1554d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change) 1555d71ae5a4SJacob Faibussowitsch { 1556a06fd7f2SStefano Zampini PetscFunctionBegin; 1557a06fd7f2SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15584f572ea9SToby Isaac PetscAssertPointer(change, 2); 1559a06fd7f2SStefano Zampini *change = PETSC_FALSE; 1560cac4c232SBarry Smith PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change)); 15613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1562a06fd7f2SStefano Zampini } 1563a06fd7f2SStefano Zampini 15644b9ad928SBarry Smith /*@ 15654b9ad928SBarry Smith PCPreSolve - Optional pre-solve phase, intended for any 15664b9ad928SBarry Smith preconditioner-specific actions that must be performed before 15674b9ad928SBarry Smith the iterative solve itself. 15684b9ad928SBarry Smith 1569c3339decSBarry Smith Collective 15704b9ad928SBarry Smith 15714b9ad928SBarry Smith Input Parameters: 15724b9ad928SBarry Smith + pc - the preconditioner context 15734b9ad928SBarry Smith - ksp - the Krylov subspace context 15744b9ad928SBarry Smith 15754b9ad928SBarry Smith Level: developer 15764b9ad928SBarry Smith 1577feefa0e1SJacob Faibussowitsch Example Usage: 15784b9ad928SBarry Smith .vb 15794b9ad928SBarry Smith PCPreSolve(pc,ksp); 158023ce1328SBarry Smith KSPSolve(ksp,b,x); 15814b9ad928SBarry Smith PCPostSolve(pc,ksp); 15824b9ad928SBarry Smith .ve 15834b9ad928SBarry Smith 15844b9ad928SBarry Smith Notes: 1585f1580f4eSBarry Smith The pre-solve phase is distinct from the `PCSetUp()` phase. 15864b9ad928SBarry Smith 1587f1580f4eSBarry Smith `KSPSolve()` calls this directly, so is rarely called by the user. 15884b9ad928SBarry Smith 1589562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCPostSolve()` 15904b9ad928SBarry Smith @*/ 1591d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPreSolve(PC pc, KSP ksp) 1592d71ae5a4SJacob Faibussowitsch { 15934b9ad928SBarry Smith Vec x, rhs; 15944b9ad928SBarry Smith 15954b9ad928SBarry Smith PetscFunctionBegin; 15960700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15970700a824SBarry Smith PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2); 1598d9a03883SBarry Smith pc->presolvedone++; 15997827d75bSBarry Smith PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice"); 16009566063dSJacob Faibussowitsch PetscCall(KSPGetSolution(ksp, &x)); 16019566063dSJacob Faibussowitsch PetscCall(KSPGetRhs(ksp, &rhs)); 16024b9ad928SBarry Smith 1603dbbe0bcdSBarry Smith if (pc->ops->presolve) PetscUseTypeMethod(pc, presolve, ksp, rhs, x); 1604f4f49eeaSPierre Jolivet else if (pc->presolve) PetscCall(pc->presolve(pc, ksp)); 16053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16064b9ad928SBarry Smith } 16074b9ad928SBarry Smith 1608f560b561SHong Zhang /*@C 1609f1580f4eSBarry Smith PCSetPreSolve - Sets function used by `PCPreSolve()` which is intended for any 1610f560b561SHong Zhang preconditioner-specific actions that must be performed before 1611f560b561SHong Zhang the iterative solve itself. 1612f560b561SHong Zhang 1613c3339decSBarry Smith Logically Collective 1614f560b561SHong Zhang 1615f560b561SHong Zhang Input Parameters: 1616f560b561SHong Zhang + pc - the preconditioner object 1617f560b561SHong Zhang - presolve - the function to call before the solve 1618f560b561SHong Zhang 161920f4b53cSBarry Smith Calling sequence of `presolve`: 162020f4b53cSBarry Smith + pc - the `PC` context 162120f4b53cSBarry Smith - ksp - the `KSP` context 1622f560b561SHong Zhang 1623f560b561SHong Zhang Level: developer 1624f560b561SHong Zhang 1625562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCPreSolve()` 1626f560b561SHong Zhang @*/ 162704c3f3b8SBarry Smith PetscErrorCode PCSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC pc, KSP ksp)) 1628d71ae5a4SJacob Faibussowitsch { 1629f560b561SHong Zhang PetscFunctionBegin; 1630f560b561SHong Zhang PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1631f560b561SHong Zhang pc->presolve = presolve; 16323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1633f560b561SHong Zhang } 1634f560b561SHong Zhang 16354b9ad928SBarry Smith /*@ 16364b9ad928SBarry Smith PCPostSolve - Optional post-solve phase, intended for any 16374b9ad928SBarry Smith preconditioner-specific actions that must be performed after 16384b9ad928SBarry Smith the iterative solve itself. 16394b9ad928SBarry Smith 1640c3339decSBarry Smith Collective 16414b9ad928SBarry Smith 16424b9ad928SBarry Smith Input Parameters: 16434b9ad928SBarry Smith + pc - the preconditioner context 16444b9ad928SBarry Smith - ksp - the Krylov subspace context 16454b9ad928SBarry Smith 1646feefa0e1SJacob Faibussowitsch Example Usage: 16474b9ad928SBarry Smith .vb 16484b9ad928SBarry Smith PCPreSolve(pc,ksp); 164923ce1328SBarry Smith KSPSolve(ksp,b,x); 16504b9ad928SBarry Smith PCPostSolve(pc,ksp); 16514b9ad928SBarry Smith .ve 16524b9ad928SBarry Smith 1653562efe2eSBarry Smith Level: developer 1654562efe2eSBarry Smith 16554b9ad928SBarry Smith Note: 1656f1580f4eSBarry Smith `KSPSolve()` calls this routine directly, so it is rarely called by the user. 16574b9ad928SBarry Smith 1658562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetPostSolve()`, `PCSetPresolve()`, `PCPreSolve()`, `KSPSolve()` 16594b9ad928SBarry Smith @*/ 1660d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPostSolve(PC pc, KSP ksp) 1661d71ae5a4SJacob Faibussowitsch { 16624b9ad928SBarry Smith Vec x, rhs; 16634b9ad928SBarry Smith 16644b9ad928SBarry Smith PetscFunctionBegin; 16650700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16660700a824SBarry Smith PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2); 1667d9a03883SBarry Smith pc->presolvedone--; 16689566063dSJacob Faibussowitsch PetscCall(KSPGetSolution(ksp, &x)); 16699566063dSJacob Faibussowitsch PetscCall(KSPGetRhs(ksp, &rhs)); 1670dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, postsolve, ksp, rhs, x); 16713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16724b9ad928SBarry Smith } 16734b9ad928SBarry Smith 1674ffeef943SBarry Smith /*@ 1675f1580f4eSBarry Smith PCLoad - Loads a `PC` that has been stored in binary with `PCView()`. 167655849f57SBarry Smith 1677c3339decSBarry Smith Collective 167855849f57SBarry Smith 167955849f57SBarry Smith Input Parameters: 1680f1580f4eSBarry Smith + newdm - the newly loaded `PC`, this needs to have been created with `PCCreate()` or 1681f1580f4eSBarry Smith some related function before a call to `PCLoad()`. 1682f1580f4eSBarry Smith - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()` 168355849f57SBarry Smith 168455849f57SBarry Smith Level: intermediate 168555849f57SBarry Smith 1686f1580f4eSBarry Smith Note: 1687562efe2eSBarry Smith The type is determined by the data in the file, any `PCType` set into the `PC` before this call is ignored. 168855849f57SBarry Smith 1689562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()` 169055849f57SBarry Smith @*/ 1691d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLoad(PC newdm, PetscViewer viewer) 1692d71ae5a4SJacob Faibussowitsch { 169355849f57SBarry Smith PetscBool isbinary; 1694060da220SMatthew G. Knepley PetscInt classid; 169555849f57SBarry Smith char type[256]; 169655849f57SBarry Smith 169755849f57SBarry Smith PetscFunctionBegin; 169855849f57SBarry Smith PetscValidHeaderSpecific(newdm, PC_CLASSID, 1); 169955849f57SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 17009566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 170128b400f6SJacob Faibussowitsch PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 170255849f57SBarry Smith 17039566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT)); 17047827d75bSBarry Smith PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file"); 17059566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR)); 17069566063dSJacob Faibussowitsch PetscCall(PCSetType(newdm, type)); 1707dbbe0bcdSBarry Smith PetscTryTypeMethod(newdm, load, viewer); 17083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 170955849f57SBarry Smith } 171055849f57SBarry Smith 17119804daf3SBarry Smith #include <petscdraw.h> 1712e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1713e04113cfSBarry Smith #include <petscviewersaws.h> 17140acecf5bSBarry Smith #endif 1715fe2efc57SMark 1716ffeef943SBarry Smith /*@ 1717562efe2eSBarry Smith PCViewFromOptions - View from the `PC` based on options in the options database 1718fe2efc57SMark 1719c3339decSBarry Smith Collective 1720fe2efc57SMark 1721fe2efc57SMark Input Parameters: 172220f4b53cSBarry Smith + A - the `PC` context 1723f1580f4eSBarry Smith . obj - Optional object that provides the options prefix 1724736c3998SJose E. Roman - name - command line option 1725fe2efc57SMark 1726fe2efc57SMark Level: intermediate 1727f1580f4eSBarry Smith 1728562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()` 1729fe2efc57SMark @*/ 1730d71ae5a4SJacob Faibussowitsch PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[]) 1731d71ae5a4SJacob Faibussowitsch { 1732fe2efc57SMark PetscFunctionBegin; 1733fe2efc57SMark PetscValidHeaderSpecific(A, PC_CLASSID, 1); 17349566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 17353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1736fe2efc57SMark } 1737fe2efc57SMark 1738ffeef943SBarry Smith /*@ 1739f1580f4eSBarry Smith PCView - Prints information about the `PC` 17404b9ad928SBarry Smith 1741c3339decSBarry Smith Collective 17424b9ad928SBarry Smith 17434b9ad928SBarry Smith Input Parameters: 1744feefa0e1SJacob Faibussowitsch + pc - the `PC` context 17454b9ad928SBarry Smith - viewer - optional visualization context 17464b9ad928SBarry Smith 174720f4b53cSBarry Smith Level: developer 174820f4b53cSBarry Smith 1749f1580f4eSBarry Smith Notes: 17504b9ad928SBarry Smith The available visualization contexts include 1751f1580f4eSBarry Smith + `PETSC_VIEWER_STDOUT_SELF` - standard output (default) 1752f1580f4eSBarry Smith - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard 17534b9ad928SBarry Smith output where only the first processor opens 17544b9ad928SBarry Smith the file. All other processors send their 17554b9ad928SBarry Smith data to the first processor to print. 17564b9ad928SBarry Smith 17574b9ad928SBarry Smith The user can open an alternative visualization contexts with 1758f1580f4eSBarry Smith `PetscViewerASCIIOpen()` (output to a specified file). 17594b9ad928SBarry Smith 1760562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewer`, `KSPView()`, `PetscViewerASCIIOpen()` 17614b9ad928SBarry Smith @*/ 1762d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView(PC pc, PetscViewer viewer) 1763d71ae5a4SJacob Faibussowitsch { 176419fd82e9SBarry Smith PCType cstr; 17656cd81132SPierre Jolivet PetscViewerFormat format; 17666cd81132SPierre Jolivet PetscBool iascii, isstring, isbinary, isdraw, pop = PETSC_FALSE; 1767e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1768536b137fSBarry Smith PetscBool issaws; 17690acecf5bSBarry Smith #endif 17704b9ad928SBarry Smith 17714b9ad928SBarry Smith PetscFunctionBegin; 17720700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 177348a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer)); 17740700a824SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1775c9780b6fSBarry Smith PetscCheckSameComm(pc, 1, viewer, 2); 17764b9ad928SBarry Smith 17779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 17789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 17799566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 17809566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 1781e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 17829566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws)); 17830acecf5bSBarry Smith #endif 1784219b1045SBarry Smith 178532077d6dSBarry Smith if (iascii) { 17869566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer)); 178748a46eb9SPierre Jolivet if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, " PC has not been set up so information may be incomplete\n")); 17889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1789dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 17909566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1791834dbeb0SBarry Smith if (pc->mat) { 17926cd81132SPierre Jolivet PetscCall(PetscViewerGetFormat(viewer, &format)); 17936cd81132SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 17949566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 17956cd81132SPierre Jolivet pop = PETSC_TRUE; 17966cd81132SPierre Jolivet } 17974b9ad928SBarry Smith if (pc->pmat == pc->mat) { 17989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix = precond matrix:\n")); 17999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 18009566063dSJacob Faibussowitsch PetscCall(MatView(pc->mat, viewer)); 18019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 18024b9ad928SBarry Smith } else { 1803834dbeb0SBarry Smith if (pc->pmat) { 18049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix followed by preconditioner matrix:\n")); 18054b9ad928SBarry Smith } else { 18069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix:\n")); 18074b9ad928SBarry Smith } 18089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 18099566063dSJacob Faibussowitsch PetscCall(MatView(pc->mat, viewer)); 18109566063dSJacob Faibussowitsch if (pc->pmat) PetscCall(MatView(pc->pmat, viewer)); 18119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 18124b9ad928SBarry Smith } 18136cd81132SPierre Jolivet if (pop) PetscCall(PetscViewerPopFormat(viewer)); 18144b9ad928SBarry Smith } 18154b9ad928SBarry Smith } else if (isstring) { 18169566063dSJacob Faibussowitsch PetscCall(PCGetType(pc, &cstr)); 18179566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr)); 1818dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 18199566063dSJacob Faibussowitsch if (pc->mat) PetscCall(MatView(pc->mat, viewer)); 18209566063dSJacob Faibussowitsch if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer)); 18215b0b0462SBarry Smith } else if (isbinary) { 182255849f57SBarry Smith PetscInt classid = PC_FILE_CLASSID; 182355849f57SBarry Smith MPI_Comm comm; 182455849f57SBarry Smith PetscMPIInt rank; 182555849f57SBarry Smith char type[256]; 182655849f57SBarry Smith 18279566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 18289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1829dd400576SPatrick Sanan if (rank == 0) { 18309566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT)); 18319566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256)); 18329566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR)); 183355849f57SBarry Smith } 1834dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 1835219b1045SBarry Smith } else if (isdraw) { 1836219b1045SBarry Smith PetscDraw draw; 1837d9884438SBarry Smith char str[25]; 183889fd9fafSBarry Smith PetscReal x, y, bottom, h; 1839d9884438SBarry Smith PetscInt n; 1840219b1045SBarry Smith 18419566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 18429566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 18431d840656SPeter Brune if (pc->mat) { 18449566063dSJacob Faibussowitsch PetscCall(MatGetSize(pc->mat, &n, NULL)); 184563a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n)); 18461d840656SPeter Brune } else { 18479566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name)); 18481d840656SPeter Brune } 18499566063dSJacob Faibussowitsch PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h)); 185089fd9fafSBarry Smith bottom = y - h; 18519566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom)); 1852dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 18539566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 1854e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1855536b137fSBarry Smith } else if (issaws) { 1856d45a07a7SBarry Smith PetscMPIInt rank; 1857d45a07a7SBarry Smith 18589566063dSJacob Faibussowitsch PetscCall(PetscObjectName((PetscObject)pc)); 18599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 186048a46eb9SPierre Jolivet if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer)); 18619566063dSJacob Faibussowitsch if (pc->mat) PetscCall(MatView(pc->mat, viewer)); 18629566063dSJacob Faibussowitsch if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer)); 18630acecf5bSBarry Smith #endif 18644b9ad928SBarry Smith } 18653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18664b9ad928SBarry Smith } 18674b9ad928SBarry Smith 18684b9ad928SBarry Smith /*@C 1869562efe2eSBarry Smith PCRegister - Adds a method (`PCType`) to the preconditioner package. 18701c84c290SBarry Smith 1871cc4c1da9SBarry Smith Not collective. No Fortran Support 18721c84c290SBarry Smith 18731c84c290SBarry Smith Input Parameters: 187420f4b53cSBarry Smith + sname - name of a new user-defined solver 187520f4b53cSBarry Smith - function - routine to create method context 18761c84c290SBarry Smith 1877feefa0e1SJacob Faibussowitsch Example Usage: 18781c84c290SBarry Smith .vb 1879bdf89e91SBarry Smith PCRegister("my_solver", MySolverCreate); 18801c84c290SBarry Smith .ve 18811c84c290SBarry Smith 18821c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 18831c84c290SBarry Smith $ PCSetType(pc, "my_solver") 18841c84c290SBarry Smith or at runtime via the option 18851c84c290SBarry Smith $ -pc_type my_solver 18864b9ad928SBarry Smith 18874b9ad928SBarry Smith Level: advanced 18881c84c290SBarry Smith 188920f4b53cSBarry Smith Note: 189020f4b53cSBarry Smith `PCRegister()` may be called multiple times to add several user-defined preconditioners. 189120f4b53cSBarry Smith 1892562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCType`, `PCRegisterAll()` 18934b9ad928SBarry Smith @*/ 1894d71ae5a4SJacob Faibussowitsch PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC)) 1895d71ae5a4SJacob Faibussowitsch { 18964b9ad928SBarry Smith PetscFunctionBegin; 18979566063dSJacob Faibussowitsch PetscCall(PCInitializePackage()); 18989566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCList, sname, function)); 18993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19004b9ad928SBarry Smith } 19014b9ad928SBarry Smith 1902d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y) 1903d71ae5a4SJacob Faibussowitsch { 1904186a3e20SStefano Zampini PC pc; 1905186a3e20SStefano Zampini 1906186a3e20SStefano Zampini PetscFunctionBegin; 19079566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &pc)); 19089566063dSJacob Faibussowitsch PetscCall(PCApply(pc, X, Y)); 19093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1910186a3e20SStefano Zampini } 1911186a3e20SStefano Zampini 19125d83a8b1SBarry Smith /*@ 19130bacdadaSStefano Zampini PCComputeOperator - Computes the explicit preconditioned operator. 19144b9ad928SBarry Smith 1915c3339decSBarry Smith Collective 19164b9ad928SBarry Smith 1917d8d19677SJose E. Roman Input Parameters: 1918186a3e20SStefano Zampini + pc - the preconditioner object 1919562efe2eSBarry Smith - mattype - the `MatType` to be used for the operator 19204b9ad928SBarry Smith 19214b9ad928SBarry Smith Output Parameter: 1922a5b23f4aSJose E. Roman . mat - the explicit preconditioned operator 19234b9ad928SBarry Smith 192420f4b53cSBarry Smith Level: advanced 192520f4b53cSBarry Smith 1926f1580f4eSBarry Smith Note: 1927186a3e20SStefano Zampini This computation is done by applying the operators to columns of the identity matrix. 1928186a3e20SStefano Zampini This routine is costly in general, and is recommended for use only with relatively small systems. 1929562efe2eSBarry Smith Currently, this routine uses a dense matrix format when `mattype` == `NULL` 19304b9ad928SBarry Smith 1931562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPComputeOperator()`, `MatType` 19324b9ad928SBarry Smith @*/ 1933d71ae5a4SJacob Faibussowitsch PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat) 1934d71ae5a4SJacob Faibussowitsch { 1935186a3e20SStefano Zampini PetscInt N, M, m, n; 1936186a3e20SStefano Zampini Mat A, Apc; 19374b9ad928SBarry Smith 19384b9ad928SBarry Smith PetscFunctionBegin; 19390700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 19404f572ea9SToby Isaac PetscAssertPointer(mat, 3); 19419566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, &A, NULL)); 19429566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 19439566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 19449566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc)); 19459566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (void (*)(void))MatMult_PC)); 19469566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(Apc, mattype, mat)); 19479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Apc)); 19483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19494b9ad928SBarry Smith } 19504b9ad928SBarry Smith 19516c237d78SBarry Smith /*@ 19526c237d78SBarry Smith PCSetCoordinates - sets the coordinates of all the nodes on the local process 19536c237d78SBarry Smith 1954c3339decSBarry Smith Collective 19556c237d78SBarry Smith 19566c237d78SBarry Smith Input Parameters: 19576c237d78SBarry Smith + pc - the solver context 19586c237d78SBarry Smith . dim - the dimension of the coordinates 1, 2, or 3 195914893cbeSStefano Zampini . nloc - the blocked size of the coordinates array 196014893cbeSStefano Zampini - coords - the coordinates array 19616c237d78SBarry Smith 19626c237d78SBarry Smith Level: intermediate 19636c237d78SBarry Smith 1964f1580f4eSBarry Smith Note: 196520f4b53cSBarry Smith `coords` is an array of the dim coordinates for the nodes on 196620f4b53cSBarry Smith the local processor, of size `dim`*`nloc`. 1967*6e415bd2SNuno Nobre If there are 108 equations (dofs) on a processor 1968*6e415bd2SNuno Nobre for a 3d displacement finite element discretization of elasticity (so 196914893cbeSStefano Zampini that there are nloc = 36 = 108/3 nodes) then the array must have 108 19706c237d78SBarry Smith double precision values (ie, 3 * 36). These x y z coordinates 19716c237d78SBarry Smith should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x, 19726c237d78SBarry Smith ... , N-1.z ]. 19736c237d78SBarry Smith 1974562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `MatSetNearNullSpace()` 19756c237d78SBarry Smith @*/ 1976d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[]) 1977d71ae5a4SJacob Faibussowitsch { 19786c237d78SBarry Smith PetscFunctionBegin; 197932954da3SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 198032954da3SStefano Zampini PetscValidLogicalCollectiveInt(pc, dim, 2); 198122794d57SStefano Zampini PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords)); 19823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19836c237d78SBarry Smith } 1984fd2dd295SFande Kong 1985fd2dd295SFande Kong /*@ 1986fd2dd295SFande Kong PCGetInterpolations - Gets interpolation matrices for all levels (except level 0) 1987fd2dd295SFande Kong 1988c3339decSBarry Smith Logically Collective 1989fd2dd295SFande Kong 1990d8d19677SJose E. Roman Input Parameter: 1991d8d19677SJose E. Roman . pc - the precondition context 1992fd2dd295SFande Kong 1993d8d19677SJose E. Roman Output Parameters: 1994d8d19677SJose E. Roman + num_levels - the number of levels 1995562efe2eSBarry Smith - interpolations - the interpolation matrices (size of `num_levels`-1) 1996fd2dd295SFande Kong 1997fd2dd295SFande Kong Level: advanced 1998fd2dd295SFande Kong 1999562efe2eSBarry Smith Developer Note: 2000f1580f4eSBarry Smith Why is this here instead of in `PCMG` etc? 2001fd2dd295SFande Kong 2002562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()` 2003fd2dd295SFande Kong @*/ 2004d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[]) 2005d71ae5a4SJacob Faibussowitsch { 2006fd2dd295SFande Kong PetscFunctionBegin; 2007fd2dd295SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 20084f572ea9SToby Isaac PetscAssertPointer(num_levels, 2); 20094f572ea9SToby Isaac PetscAssertPointer(interpolations, 3); 2010cac4c232SBarry Smith PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations)); 20113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2012fd2dd295SFande Kong } 2013fd2dd295SFande Kong 2014fd2dd295SFande Kong /*@ 2015fd2dd295SFande Kong PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level) 2016fd2dd295SFande Kong 2017c3339decSBarry Smith Logically Collective 2018fd2dd295SFande Kong 2019d8d19677SJose E. Roman Input Parameter: 2020d8d19677SJose E. Roman . pc - the precondition context 2021fd2dd295SFande Kong 2022d8d19677SJose E. Roman Output Parameters: 2023d8d19677SJose E. Roman + num_levels - the number of levels 2024562efe2eSBarry Smith - coarseOperators - the coarse operator matrices (size of `num_levels`-1) 2025fd2dd295SFande Kong 2026fd2dd295SFande Kong Level: advanced 2027fd2dd295SFande Kong 2028562efe2eSBarry Smith Developer Note: 2029f1580f4eSBarry Smith Why is this here instead of in `PCMG` etc? 2030fd2dd295SFande Kong 2031562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()` 2032fd2dd295SFande Kong @*/ 2033d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[]) 2034d71ae5a4SJacob Faibussowitsch { 2035fd2dd295SFande Kong PetscFunctionBegin; 2036fd2dd295SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 20374f572ea9SToby Isaac PetscAssertPointer(num_levels, 2); 20384f572ea9SToby Isaac PetscAssertPointer(coarseOperators, 3); 2039cac4c232SBarry Smith PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators)); 20403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2041fd2dd295SFande Kong } 2042