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 5588584be7SBarry Smith /*@ 56562efe2eSBarry Smith PCReset - Resets a `PC` context to the pcsetupcalled = 0 state and removes any allocated `Vec`s and `Mat`s 5788584be7SBarry Smith 58c3339decSBarry Smith Collective 5988584be7SBarry Smith 6088584be7SBarry Smith Input Parameter: 6188584be7SBarry Smith . pc - the preconditioner context 6288584be7SBarry Smith 6388584be7SBarry Smith Level: developer 6488584be7SBarry Smith 65f1580f4eSBarry Smith Note: 66562efe2eSBarry 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` 6788584be7SBarry Smith 68562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()` 6988584be7SBarry Smith @*/ 70d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset(PC pc) 71d71ae5a4SJacob Faibussowitsch { 7288584be7SBarry Smith PetscFunctionBegin; 7388584be7SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 74dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, reset); 759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc->diagonalscaleright)); 769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc->diagonalscaleleft)); 779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->pmat)); 789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->mat)); 792fa5cd67SKarl Rupp 800ce6c5a2SBarry Smith pc->setupcalled = 0; 813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8288584be7SBarry Smith } 8388584be7SBarry Smith 841fb7b255SJunchao Zhang /*@C 85f1580f4eSBarry Smith PCDestroy - Destroys `PC` context that was created with `PCCreate()`. 864b9ad928SBarry Smith 87c3339decSBarry Smith Collective 884b9ad928SBarry Smith 894b9ad928SBarry Smith Input Parameter: 904b9ad928SBarry Smith . pc - the preconditioner context 914b9ad928SBarry Smith 924b9ad928SBarry Smith Level: developer 934b9ad928SBarry Smith 94562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()` 954b9ad928SBarry Smith @*/ 96d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy(PC *pc) 97d71ae5a4SJacob Faibussowitsch { 984b9ad928SBarry Smith PetscFunctionBegin; 993ba16761SJacob Faibussowitsch if (!*pc) PetscFunctionReturn(PETSC_SUCCESS); 100f4f49eeaSPierre Jolivet PetscValidHeaderSpecific(*pc, PC_CLASSID, 1); 101f4f49eeaSPierre Jolivet if (--((PetscObject)*pc)->refct > 0) { 1029371c9d4SSatish Balay *pc = NULL; 1033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1049371c9d4SSatish Balay } 1054b9ad928SBarry Smith 1069566063dSJacob Faibussowitsch PetscCall(PCReset(*pc)); 107241cb3abSLisandro Dalcin 108e04113cfSBarry Smith /* if memory was published with SAWs then destroy it */ 1099566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsViewOff((PetscObject)*pc)); 110f4f49eeaSPierre Jolivet PetscTryTypeMethod(*pc, destroy); 1119566063dSJacob Faibussowitsch PetscCall(DMDestroy(&(*pc)->dm)); 1129566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(pc)); 1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1144b9ad928SBarry Smith } 1154b9ad928SBarry Smith 116cc4c1da9SBarry Smith /*@ 117c5eb9154SBarry Smith PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right 1184b9ad928SBarry Smith scaling as needed by certain time-stepping codes. 1194b9ad928SBarry Smith 120c3339decSBarry Smith Logically Collective 1214b9ad928SBarry Smith 1224b9ad928SBarry Smith Input Parameter: 1234b9ad928SBarry Smith . pc - the preconditioner context 1244b9ad928SBarry Smith 1254b9ad928SBarry Smith Output Parameter: 126f1580f4eSBarry Smith . flag - `PETSC_TRUE` if it applies the scaling 1274b9ad928SBarry Smith 1284b9ad928SBarry Smith Level: developer 1294b9ad928SBarry Smith 130f1580f4eSBarry Smith Note: 131562efe2eSBarry Smith If this returns `PETSC_TRUE` then the system solved via the Krylov method is, for left and right preconditioning, 1324b9ad928SBarry Smith 133562efe2eSBarry Smith $$ 134562efe2eSBarry Smith \begin{align*} 135562efe2eSBarry Smith D M A D^{-1} y = D M b \\ 136562efe2eSBarry Smith D A M D^{-1} z = D b. 137562efe2eSBarry Smith \end{align*} 138562efe2eSBarry Smith $$ 139562efe2eSBarry Smith 140562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCSetDiagonalScale()` 1414b9ad928SBarry Smith @*/ 142d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetDiagonalScale(PC pc, PetscBool *flag) 143d71ae5a4SJacob Faibussowitsch { 1444b9ad928SBarry Smith PetscFunctionBegin; 1450700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1464f572ea9SToby Isaac PetscAssertPointer(flag, 2); 1474b9ad928SBarry Smith *flag = pc->diagonalscale; 1483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1494b9ad928SBarry Smith } 1504b9ad928SBarry Smith 1514b9ad928SBarry Smith /*@ 152c5eb9154SBarry Smith PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right 1534b9ad928SBarry Smith scaling as needed by certain time-stepping codes. 1544b9ad928SBarry Smith 155c3339decSBarry Smith Logically Collective 1564b9ad928SBarry Smith 1574b9ad928SBarry Smith Input Parameters: 1584b9ad928SBarry Smith + pc - the preconditioner context 1594b9ad928SBarry Smith - s - scaling vector 1604b9ad928SBarry Smith 1614b9ad928SBarry Smith Level: intermediate 1624b9ad928SBarry Smith 16395452b02SPatrick Sanan Notes: 164562efe2eSBarry Smith The system solved via the Krylov method is, for left and right preconditioning, 165562efe2eSBarry Smith $$ 166562efe2eSBarry Smith \begin{align*} 167562efe2eSBarry Smith D M A D^{-1} y = D M b \\ 168562efe2eSBarry Smith D A M D^{-1} z = D b. 169562efe2eSBarry Smith \end{align*} 170562efe2eSBarry Smith $$ 1714b9ad928SBarry Smith 172562efe2eSBarry Smith `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$. 1734b9ad928SBarry Smith 174562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCGetDiagonalScale()` 1754b9ad928SBarry Smith @*/ 176d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetDiagonalScale(PC pc, Vec s) 177d71ae5a4SJacob Faibussowitsch { 1784b9ad928SBarry Smith PetscFunctionBegin; 1790700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1800700a824SBarry Smith PetscValidHeaderSpecific(s, VEC_CLASSID, 2); 1814b9ad928SBarry Smith pc->diagonalscale = PETSC_TRUE; 1822fa5cd67SKarl Rupp 1839566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)s)); 1849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc->diagonalscaleleft)); 1852fa5cd67SKarl Rupp 1864b9ad928SBarry Smith pc->diagonalscaleleft = s; 1872fa5cd67SKarl Rupp 1889566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s, &pc->diagonalscaleright)); 1899566063dSJacob Faibussowitsch PetscCall(VecCopy(s, pc->diagonalscaleright)); 1909566063dSJacob Faibussowitsch PetscCall(VecReciprocal(pc->diagonalscaleright)); 1913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1924b9ad928SBarry Smith } 1934b9ad928SBarry Smith 194bc08b0f1SBarry Smith /*@ 195c5eb9154SBarry Smith PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes. 1964b9ad928SBarry Smith 197c3339decSBarry Smith Logically Collective 1984b9ad928SBarry Smith 1994b9ad928SBarry Smith Input Parameters: 2004b9ad928SBarry Smith + pc - the preconditioner context 2014b9ad928SBarry Smith . in - input vector 202a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in) 2034b9ad928SBarry Smith 2044b9ad928SBarry Smith Level: intermediate 2054b9ad928SBarry Smith 20695452b02SPatrick Sanan Notes: 207562efe2eSBarry Smith The system solved via the Krylov method is, for left and right preconditioning, 2084b9ad928SBarry Smith 209562efe2eSBarry Smith $$ 210562efe2eSBarry Smith \begin{align*} 211562efe2eSBarry Smith D M A D^{-1} y = D M b \\ 212562efe2eSBarry Smith D A M D^{-1} z = D b. 213562efe2eSBarry Smith \end{align*} 214562efe2eSBarry Smith $$ 2154b9ad928SBarry Smith 216562efe2eSBarry Smith `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$. 2174b9ad928SBarry Smith 218562efe2eSBarry Smith If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out` 219562efe2eSBarry Smith 220562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCSetDiagonalScale()`, `PCDiagonalScaleRight()`, `PCDiagonalScale()` 2214b9ad928SBarry Smith @*/ 222d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleLeft(PC pc, Vec in, Vec out) 223d71ae5a4SJacob Faibussowitsch { 2244b9ad928SBarry Smith PetscFunctionBegin; 2250700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2260700a824SBarry Smith PetscValidHeaderSpecific(in, VEC_CLASSID, 2); 2270700a824SBarry Smith PetscValidHeaderSpecific(out, VEC_CLASSID, 3); 2284b9ad928SBarry Smith if (pc->diagonalscale) { 2299566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(out, pc->diagonalscaleleft, in)); 2304b9ad928SBarry Smith } else if (in != out) { 2319566063dSJacob Faibussowitsch PetscCall(VecCopy(in, out)); 2324b9ad928SBarry Smith } 2333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2344b9ad928SBarry Smith } 2354b9ad928SBarry Smith 236bc08b0f1SBarry Smith /*@ 2374b9ad928SBarry Smith PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes. 2384b9ad928SBarry Smith 239c3339decSBarry Smith Logically Collective 2404b9ad928SBarry Smith 2414b9ad928SBarry Smith Input Parameters: 2424b9ad928SBarry Smith + pc - the preconditioner context 2434b9ad928SBarry Smith . in - input vector 244a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in) 2454b9ad928SBarry Smith 2464b9ad928SBarry Smith Level: intermediate 2474b9ad928SBarry Smith 24895452b02SPatrick Sanan Notes: 249562efe2eSBarry Smith The system solved via the Krylov method is, for left and right preconditioning, 2504b9ad928SBarry Smith 251562efe2eSBarry Smith $$ 252562efe2eSBarry Smith \begin{align*} 253562efe2eSBarry Smith D M A D^{-1} y = D M b \\ 254562efe2eSBarry Smith D A M D^{-1} z = D b. 255a3524536SPierre Jolivet \end{align*} 256562efe2eSBarry Smith $$ 2574b9ad928SBarry Smith 258562efe2eSBarry Smith `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$. 2594b9ad928SBarry Smith 260562efe2eSBarry Smith If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out` 261562efe2eSBarry Smith 262562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCSetDiagonalScale()`, `PCDiagonalScale()` 2634b9ad928SBarry Smith @*/ 264d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleRight(PC pc, Vec in, Vec out) 265d71ae5a4SJacob Faibussowitsch { 2664b9ad928SBarry Smith PetscFunctionBegin; 2670700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2680700a824SBarry Smith PetscValidHeaderSpecific(in, VEC_CLASSID, 2); 2690700a824SBarry Smith PetscValidHeaderSpecific(out, VEC_CLASSID, 3); 2704b9ad928SBarry Smith if (pc->diagonalscale) { 2719566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(out, pc->diagonalscaleright, in)); 2724b9ad928SBarry Smith } else if (in != out) { 2739566063dSJacob Faibussowitsch PetscCall(VecCopy(in, out)); 2744b9ad928SBarry Smith } 2753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2764b9ad928SBarry Smith } 2774b9ad928SBarry Smith 27849517cdeSBarry Smith /*@ 27949517cdeSBarry Smith PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the 280f1580f4eSBarry Smith operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`, 281f1580f4eSBarry Smith `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat. 28249517cdeSBarry Smith 283c3339decSBarry Smith Logically Collective 28449517cdeSBarry Smith 28549517cdeSBarry Smith Input Parameters: 28649517cdeSBarry Smith + pc - the preconditioner context 287f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false) 28849517cdeSBarry Smith 28949517cdeSBarry Smith Options Database Key: 290562efe2eSBarry Smith . -pc_use_amat <true,false> - use the amat argument to `KSPSetOperators()` or `PCSetOperators()` to apply the operator 29149517cdeSBarry Smith 29220f4b53cSBarry Smith Level: intermediate 29320f4b53cSBarry Smith 294f1580f4eSBarry Smith Note: 29549517cdeSBarry Smith For the common case in which the linear system matrix and the matrix used to construct the 296562efe2eSBarry Smith preconditioner are identical, this routine has no affect. 29749517cdeSBarry Smith 298562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetUseAmat()`, `PCBJACOBI`, `PGMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`, 299562efe2eSBarry Smith `KSPSetOperators()`, `PCSetOperators()` 30049517cdeSBarry Smith @*/ 301d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUseAmat(PC pc, PetscBool flg) 302d71ae5a4SJacob Faibussowitsch { 30349517cdeSBarry Smith PetscFunctionBegin; 30449517cdeSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 30549517cdeSBarry Smith pc->useAmat = flg; 3063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 30749517cdeSBarry Smith } 30849517cdeSBarry Smith 309422a814eSBarry Smith /*@ 310562efe2eSBarry Smith PCSetErrorIfFailure - Causes `PC` to generate an error if a floating point exception, for example a zero pivot, is detected. 311422a814eSBarry Smith 312c3339decSBarry Smith Logically Collective 313422a814eSBarry Smith 314422a814eSBarry Smith Input Parameters: 315562efe2eSBarry Smith + pc - iterative context obtained from `PCCreate()` 316f1580f4eSBarry Smith - flg - `PETSC_TRUE` indicates you want the error generated 317422a814eSBarry Smith 318422a814eSBarry Smith Level: advanced 319422a814eSBarry Smith 320422a814eSBarry Smith Notes: 321f1580f4eSBarry Smith Normally PETSc continues if a linear solver fails due to a failed setup of a preconditioner, you can call `KSPGetConvergedReason()` after a `KSPSolve()` 322422a814eSBarry 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 323422a814eSBarry Smith detected. 324422a814eSBarry Smith 325562efe2eSBarry Smith This is propagated into `KSP`s used by this `PC`, which then propagate it into `PC`s used by those `KSP`s 326422a814eSBarry Smith 327562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPSetErrorIfNotConverged()`, `PCGetInitialGuessNonzero()`, `PCSetInitialGuessKnoll()`, `PCGetInitialGuessKnoll()` 328422a814eSBarry Smith @*/ 329d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetErrorIfFailure(PC pc, PetscBool flg) 330d71ae5a4SJacob Faibussowitsch { 331422a814eSBarry Smith PetscFunctionBegin; 332422a814eSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 333422a814eSBarry Smith PetscValidLogicalCollectiveBool(pc, flg, 2); 334422a814eSBarry Smith pc->erroriffailure = flg; 3353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 336422a814eSBarry Smith } 337422a814eSBarry Smith 33849517cdeSBarry Smith /*@ 33949517cdeSBarry Smith PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the 340f1580f4eSBarry Smith operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`, 341f1580f4eSBarry Smith `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat. 34249517cdeSBarry Smith 343c3339decSBarry Smith Logically Collective 34449517cdeSBarry Smith 34549517cdeSBarry Smith Input Parameter: 34649517cdeSBarry Smith . pc - the preconditioner context 34749517cdeSBarry Smith 34849517cdeSBarry Smith Output Parameter: 349f1580f4eSBarry Smith . flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false) 35049517cdeSBarry Smith 35120f4b53cSBarry Smith Level: intermediate 35220f4b53cSBarry Smith 353f1580f4eSBarry Smith Note: 35449517cdeSBarry Smith For the common case in which the linear system matrix and the matrix used to construct the 35549517cdeSBarry Smith preconditioner are identical, this routine is does nothing. 35649517cdeSBarry Smith 357562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUseAmat()`, `PCBJACOBI`, `PGMG`, `PCFIELDSPLIT`, `PCCOMPOSITE` 35849517cdeSBarry Smith @*/ 359d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetUseAmat(PC pc, PetscBool *flg) 360d71ae5a4SJacob Faibussowitsch { 36149517cdeSBarry Smith PetscFunctionBegin; 36249517cdeSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 36349517cdeSBarry Smith *flg = pc->useAmat; 3643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 36549517cdeSBarry Smith } 36649517cdeSBarry Smith 367f39d8e23SSatish Balay /*@ 3683821be0aSBarry Smith PCSetKSPNestLevel - sets the amount of nesting the `KSP` that contains this `PC` has 3693821be0aSBarry Smith 3703821be0aSBarry Smith Collective 3713821be0aSBarry Smith 3723821be0aSBarry Smith Input Parameters: 3733821be0aSBarry Smith + pc - the `PC` 3743821be0aSBarry Smith - level - the nest level 3753821be0aSBarry Smith 3763821be0aSBarry Smith Level: developer 3773821be0aSBarry Smith 3783821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCGetKSPNestLevel()`, `KSPSetNestLevel()` 3793821be0aSBarry Smith @*/ 3803821be0aSBarry Smith PetscErrorCode PCSetKSPNestLevel(PC pc, PetscInt level) 3813821be0aSBarry Smith { 3823821be0aSBarry Smith PetscFunctionBegin; 3837a99bfcaSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 3847a99bfcaSBarry Smith PetscValidLogicalCollectiveInt(pc, level, 2); 3853821be0aSBarry Smith pc->kspnestlevel = level; 3863821be0aSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 3873821be0aSBarry Smith } 3883821be0aSBarry Smith 3893821be0aSBarry Smith /*@ 3903821be0aSBarry Smith PCGetKSPNestLevel - gets the amount of nesting the `KSP` that contains this `PC` has 3913821be0aSBarry Smith 3927a99bfcaSBarry Smith Not Collective 3933821be0aSBarry Smith 3943821be0aSBarry Smith Input Parameter: 3953821be0aSBarry Smith . pc - the `PC` 3963821be0aSBarry Smith 3973821be0aSBarry Smith Output Parameter: 3983821be0aSBarry Smith . level - the nest level 3993821be0aSBarry Smith 4003821be0aSBarry Smith Level: developer 4013821be0aSBarry Smith 4023821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `KSPGetNestLevel()` 4033821be0aSBarry Smith @*/ 4043821be0aSBarry Smith PetscErrorCode PCGetKSPNestLevel(PC pc, PetscInt *level) 4053821be0aSBarry Smith { 4063821be0aSBarry Smith PetscFunctionBegin; 4077a99bfcaSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 4087a99bfcaSBarry Smith PetscAssertPointer(level, 2); 4093821be0aSBarry Smith *level = pc->kspnestlevel; 4103821be0aSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 4113821be0aSBarry Smith } 4123821be0aSBarry Smith 4133821be0aSBarry Smith /*@ 414f1580f4eSBarry Smith PCCreate - Creates a preconditioner context, `PC` 4154b9ad928SBarry Smith 416d083f849SBarry Smith Collective 4174b9ad928SBarry Smith 4184b9ad928SBarry Smith Input Parameter: 4194b9ad928SBarry Smith . comm - MPI communicator 4204b9ad928SBarry Smith 4214b9ad928SBarry Smith Output Parameter: 4227a99bfcaSBarry Smith . newpc - location to put the preconditioner context 4234b9ad928SBarry Smith 42420f4b53cSBarry Smith Level: developer 42520f4b53cSBarry Smith 426f1580f4eSBarry Smith Note: 427f1580f4eSBarry 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` 428f1580f4eSBarry Smith in parallel. For dense matrices it is always `PCNONE`. 4294b9ad928SBarry Smith 430562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCApply()`, `PCDestroy()` 4314b9ad928SBarry Smith @*/ 432d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreate(MPI_Comm comm, PC *newpc) 433d71ae5a4SJacob Faibussowitsch { 4344b9ad928SBarry Smith PC pc; 4354b9ad928SBarry Smith 4364b9ad928SBarry Smith PetscFunctionBegin; 4374f572ea9SToby Isaac PetscAssertPointer(newpc, 2); 4380a545947SLisandro Dalcin *newpc = NULL; 4399566063dSJacob Faibussowitsch PetscCall(PCInitializePackage()); 4404b9ad928SBarry Smith 4419566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(pc, PC_CLASSID, "PC", "Preconditioner", "PC", comm, PCDestroy, PCView)); 442a7d21a52SLisandro Dalcin 4430a545947SLisandro Dalcin pc->mat = NULL; 4440a545947SLisandro Dalcin pc->pmat = NULL; 4454b9ad928SBarry Smith pc->setupcalled = 0; 4460c24e6a1SHong Zhang pc->setfromoptionscalled = 0; 4470a545947SLisandro Dalcin pc->data = NULL; 4484b9ad928SBarry Smith pc->diagonalscale = PETSC_FALSE; 4490a545947SLisandro Dalcin pc->diagonalscaleleft = NULL; 4500a545947SLisandro Dalcin pc->diagonalscaleright = NULL; 4514b9ad928SBarry Smith 4520a545947SLisandro Dalcin pc->modifysubmatrices = NULL; 4530a545947SLisandro Dalcin pc->modifysubmatricesP = NULL; 4542fa5cd67SKarl Rupp 455a7d21a52SLisandro Dalcin *newpc = pc; 4563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4574b9ad928SBarry Smith } 4584b9ad928SBarry Smith 4594b9ad928SBarry Smith /*@ 4604b9ad928SBarry Smith PCApply - Applies the preconditioner to a vector. 4614b9ad928SBarry Smith 462c3339decSBarry Smith Collective 4634b9ad928SBarry Smith 4644b9ad928SBarry Smith Input Parameters: 4654b9ad928SBarry Smith + pc - the preconditioner context 4664b9ad928SBarry Smith - x - input vector 4674b9ad928SBarry Smith 4684b9ad928SBarry Smith Output Parameter: 4694b9ad928SBarry Smith . y - output vector 4704b9ad928SBarry Smith 4714b9ad928SBarry Smith Level: developer 4724b9ad928SBarry Smith 473562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`, `PCApplyBAorAB()` 4744b9ad928SBarry Smith @*/ 475d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApply(PC pc, Vec x, Vec y) 476d71ae5a4SJacob Faibussowitsch { 477106e20bfSBarry Smith PetscInt m, n, mv, nv; 4784b9ad928SBarry Smith 4794b9ad928SBarry Smith PetscFunctionBegin; 4800700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 4810700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 4820700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 4837827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 484e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 485540bdfbdSVaclav Hapla /* use pmat to check vector sizes since for KSPLSQR the pmat may be of a different size than mat */ 4869566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &m, &n)); 4879566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &mv)); 4889566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &nv)); 489540bdfbdSVaclav Hapla /* check pmat * y = x is feasible */ 49063a3b9bcSJacob 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); 49163a3b9bcSJacob 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); 4929566063dSJacob Faibussowitsch PetscCall(VecSetErrorIfLocked(y, 3)); 493106e20bfSBarry Smith 4949566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 4959566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(x)); 4969566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0)); 497dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, apply, x, y); 4989566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0)); 499e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 5009566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(x)); 5013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5024b9ad928SBarry Smith } 5034b9ad928SBarry Smith 5044b9ad928SBarry Smith /*@ 505562efe2eSBarry Smith PCMatApply - Applies the preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApply()`, `Y` and `X` must be different matrices. 506c41ea70eSPierre Jolivet 507c3339decSBarry Smith Collective 508c677e75fSPierre Jolivet 509c677e75fSPierre Jolivet Input Parameters: 510c677e75fSPierre Jolivet + pc - the preconditioner context 511c677e75fSPierre Jolivet - X - block of input vectors 512c677e75fSPierre Jolivet 513c677e75fSPierre Jolivet Output Parameter: 514c677e75fSPierre Jolivet . Y - block of output vectors 515c677e75fSPierre Jolivet 516c677e75fSPierre Jolivet Level: developer 517c677e75fSPierre Jolivet 518562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `KSPMatSolve()` 519c677e75fSPierre Jolivet @*/ 520d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMatApply(PC pc, Mat X, Mat Y) 521d71ae5a4SJacob Faibussowitsch { 522c677e75fSPierre Jolivet Mat A; 523c677e75fSPierre Jolivet Vec cy, cx; 524bd82155bSPierre Jolivet PetscInt m1, M1, m2, M2, n1, N1, n2, N2, m3, M3, n3, N3; 525c677e75fSPierre Jolivet PetscBool match; 526c677e75fSPierre Jolivet 527c677e75fSPierre Jolivet PetscFunctionBegin; 528c677e75fSPierre Jolivet PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 529e2b200f6SPierre Jolivet PetscValidHeaderSpecific(X, MAT_CLASSID, 2); 530e2b200f6SPierre Jolivet PetscValidHeaderSpecific(Y, MAT_CLASSID, 3); 531e2b200f6SPierre Jolivet PetscCheckSameComm(pc, 1, X, 2); 532e2b200f6SPierre Jolivet PetscCheckSameComm(pc, 1, Y, 3); 5337827d75bSBarry Smith PetscCheck(Y != X, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "Y and X must be different matrices"); 5349566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, NULL, &A)); 5359566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m3, &n3)); 5369566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(X, &m2, &n2)); 5379566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Y, &m1, &n1)); 5389566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M3, &N3)); 5399566063dSJacob Faibussowitsch PetscCall(MatGetSize(X, &M2, &N2)); 5409566063dSJacob Faibussowitsch PetscCall(MatGetSize(Y, &M1, &N1)); 54163a3b9bcSJacob 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); 54263a3b9bcSJacob 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); 54363a3b9bcSJacob 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); 5449566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &match, MATSEQDENSE, MATMPIDENSE, "")); 54528b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of output vectors not stored in a dense Mat"); 5469566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, "")); 54728b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of input vectors not stored in a dense Mat"); 5489566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 549c677e75fSPierre Jolivet if (pc->ops->matapply) { 5509566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0)); 551dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, matapply, X, Y); 5529566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0)); 553c677e75fSPierre Jolivet } else { 5549566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "PC type %s applying column by column\n", ((PetscObject)pc)->type_name)); 555bd82155bSPierre Jolivet for (n1 = 0; n1 < N1; ++n1) { 5569566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(X, n1, &cx)); 5579566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Y, n1, &cy)); 5589566063dSJacob Faibussowitsch PetscCall(PCApply(pc, cx, cy)); 5599566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Y, n1, &cy)); 5609566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(X, n1, &cx)); 561c677e75fSPierre Jolivet } 562c677e75fSPierre Jolivet } 5633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 564c677e75fSPierre Jolivet } 565c677e75fSPierre Jolivet 566c677e75fSPierre Jolivet /*@ 5674b9ad928SBarry Smith PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector. 5684b9ad928SBarry Smith 569c3339decSBarry Smith Collective 5704b9ad928SBarry Smith 5714b9ad928SBarry Smith Input Parameters: 5724b9ad928SBarry Smith + pc - the preconditioner context 5734b9ad928SBarry Smith - x - input vector 5744b9ad928SBarry Smith 5754b9ad928SBarry Smith Output Parameter: 5764b9ad928SBarry Smith . y - output vector 5774b9ad928SBarry Smith 57820f4b53cSBarry Smith Level: developer 57920f4b53cSBarry Smith 580f1580f4eSBarry Smith Note: 581f1580f4eSBarry Smith Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners. 5824b9ad928SBarry Smith 583562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricRight()` 5844b9ad928SBarry Smith @*/ 585d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricLeft(PC pc, Vec x, Vec y) 586d71ae5a4SJacob Faibussowitsch { 5874b9ad928SBarry Smith PetscFunctionBegin; 5880700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 5890700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 5900700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 5917827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 592e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 5939566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 5949566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(x)); 5959566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplySymmetricLeft, pc, x, y, 0)); 596dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applysymmetricleft, x, y); 5979566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplySymmetricLeft, pc, x, y, 0)); 5989566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(x)); 599e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 6003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6014b9ad928SBarry Smith } 6024b9ad928SBarry Smith 6034b9ad928SBarry Smith /*@ 6044b9ad928SBarry Smith PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector. 6054b9ad928SBarry Smith 606c3339decSBarry Smith Collective 6074b9ad928SBarry Smith 6084b9ad928SBarry Smith Input Parameters: 6094b9ad928SBarry Smith + pc - the preconditioner context 6104b9ad928SBarry Smith - x - input vector 6114b9ad928SBarry Smith 6124b9ad928SBarry Smith Output Parameter: 6134b9ad928SBarry Smith . y - output vector 6144b9ad928SBarry Smith 6154b9ad928SBarry Smith Level: developer 6164b9ad928SBarry Smith 617f1580f4eSBarry Smith Note: 618f1580f4eSBarry Smith Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners. 6194b9ad928SBarry Smith 620562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricLeft()` 6214b9ad928SBarry Smith @*/ 622d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricRight(PC pc, Vec x, Vec y) 623d71ae5a4SJacob Faibussowitsch { 6244b9ad928SBarry Smith PetscFunctionBegin; 6250700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 6260700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 6270700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 6287827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 629e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 6309566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6319566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(x)); 6329566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplySymmetricRight, pc, x, y, 0)); 633dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applysymmetricright, x, y); 6349566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplySymmetricRight, pc, x, y, 0)); 6359566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(x)); 636e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 6373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6384b9ad928SBarry Smith } 6394b9ad928SBarry Smith 6404b9ad928SBarry Smith /*@ 6414b9ad928SBarry Smith PCApplyTranspose - Applies the transpose of preconditioner to a vector. 6424b9ad928SBarry Smith 643c3339decSBarry Smith Collective 6444b9ad928SBarry Smith 6454b9ad928SBarry Smith Input Parameters: 6464b9ad928SBarry Smith + pc - the preconditioner context 6474b9ad928SBarry Smith - x - input vector 6484b9ad928SBarry Smith 6494b9ad928SBarry Smith Output Parameter: 6504b9ad928SBarry Smith . y - output vector 6514b9ad928SBarry Smith 65220f4b53cSBarry Smith Level: developer 65320f4b53cSBarry Smith 654f1580f4eSBarry Smith Note: 65595452b02SPatrick Sanan For complex numbers this applies the non-Hermitian transpose. 6564c97465dSBarry Smith 657562efe2eSBarry Smith Developer Note: 658f1580f4eSBarry Smith We need to implement a `PCApplyHermitianTranspose()` 6594c97465dSBarry Smith 660562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyBAorAB()`, `PCApplyBAorABTranspose()`, `PCApplyTransposeExists()` 6614b9ad928SBarry Smith @*/ 662d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTranspose(PC pc, Vec x, Vec y) 663d71ae5a4SJacob Faibussowitsch { 6644b9ad928SBarry Smith PetscFunctionBegin; 6650700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 6660700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 6670700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 6687827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 669e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 6709566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6719566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(x)); 6729566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0)); 673dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applytranspose, x, y); 6749566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0)); 6759566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(x)); 676e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 6773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6784b9ad928SBarry Smith } 6794b9ad928SBarry Smith 6804b9ad928SBarry Smith /*@ 6819cc050e5SLisandro Dalcin PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation 6824b9ad928SBarry Smith 683c3339decSBarry Smith Collective 6844b9ad928SBarry Smith 685f1580f4eSBarry Smith Input Parameter: 6864b9ad928SBarry Smith . pc - the preconditioner context 6874b9ad928SBarry Smith 6884b9ad928SBarry Smith Output Parameter: 689f1580f4eSBarry Smith . flg - `PETSC_TRUE` if a transpose operation is defined 6904b9ad928SBarry Smith 6914b9ad928SBarry Smith Level: developer 6924b9ad928SBarry Smith 693562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()` 6944b9ad928SBarry Smith @*/ 695d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTransposeExists(PC pc, PetscBool *flg) 696d71ae5a4SJacob Faibussowitsch { 6974b9ad928SBarry Smith PetscFunctionBegin; 6980700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 6994f572ea9SToby Isaac PetscAssertPointer(flg, 2); 7006ce5e81cSLisandro Dalcin if (pc->ops->applytranspose) *flg = PETSC_TRUE; 7016ce5e81cSLisandro Dalcin else *flg = PETSC_FALSE; 7023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7034b9ad928SBarry Smith } 7044b9ad928SBarry Smith 7054b9ad928SBarry Smith /*@ 706562efe2eSBarry Smith PCApplyBAorAB - Applies the preconditioner and operator to a vector. $y = B*A*x $ or $ y = A*B*x$. 7074b9ad928SBarry Smith 708c3339decSBarry Smith Collective 7094b9ad928SBarry Smith 7104b9ad928SBarry Smith Input Parameters: 7114b9ad928SBarry Smith + pc - the preconditioner context 712f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC` 7134b9ad928SBarry Smith . x - input vector 7144b9ad928SBarry Smith - work - work vector 7154b9ad928SBarry Smith 7164b9ad928SBarry Smith Output Parameter: 7174b9ad928SBarry Smith . y - output vector 7184b9ad928SBarry Smith 7194b9ad928SBarry Smith Level: developer 7204b9ad928SBarry Smith 721f1580f4eSBarry Smith Note: 722562efe2eSBarry 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. 723562efe2eSBarry Smith The specific `KSPSolve()` method must also be written to handle the post-solve "correction" for the diagonal scaling. 72413e0d083SBarry Smith 725562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorABTranspose()` 7264b9ad928SBarry Smith @*/ 727d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorAB(PC pc, PCSide side, Vec x, Vec y, Vec work) 728d71ae5a4SJacob Faibussowitsch { 7294b9ad928SBarry Smith PetscFunctionBegin; 7300700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 731a69c7061SStefano Zampini PetscValidLogicalCollectiveEnum(pc, side, 2); 7320700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 7330700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 4); 7340700a824SBarry Smith PetscValidHeaderSpecific(work, VEC_CLASSID, 5); 735a69c7061SStefano Zampini PetscCheckSameComm(pc, 1, x, 3); 736a69c7061SStefano Zampini PetscCheckSameComm(pc, 1, y, 4); 737a69c7061SStefano Zampini PetscCheckSameComm(pc, 1, work, 5); 7387827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 7397827d75bSBarry 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"); 7407827d75bSBarry Smith PetscCheck(!pc->diagonalscale || side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot include diagonal scaling with symmetric preconditioner application"); 741e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE)); 7424b9ad928SBarry Smith 7439566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 7444b9ad928SBarry Smith if (pc->diagonalscale) { 7454b9ad928SBarry Smith if (pc->ops->applyBA) { 7464b9ad928SBarry Smith Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */ 7479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &work2)); 7489566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleRight(pc, x, work2)); 749dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applyBA, side, work2, y, work); 7509566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleLeft(pc, y, y)); 7519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work2)); 7524b9ad928SBarry Smith } else if (side == PC_RIGHT) { 7539566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleRight(pc, x, y)); 7549566063dSJacob Faibussowitsch PetscCall(PCApply(pc, y, work)); 7559566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, work, y)); 7569566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleLeft(pc, y, y)); 7574b9ad928SBarry Smith } else if (side == PC_LEFT) { 7589566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleRight(pc, x, y)); 7599566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, y, work)); 7609566063dSJacob Faibussowitsch PetscCall(PCApply(pc, work, y)); 7619566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleLeft(pc, y, y)); 7627827d75bSBarry Smith } else PetscCheck(side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot provide diagonal scaling with symmetric application of preconditioner"); 7634b9ad928SBarry Smith } else { 7644b9ad928SBarry Smith if (pc->ops->applyBA) { 765dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applyBA, side, x, y, work); 7664b9ad928SBarry Smith } else if (side == PC_RIGHT) { 7679566063dSJacob Faibussowitsch PetscCall(PCApply(pc, x, work)); 7689566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, work, y)); 7694b9ad928SBarry Smith } else if (side == PC_LEFT) { 7709566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, x, work)); 7719566063dSJacob Faibussowitsch PetscCall(PCApply(pc, work, y)); 7724b9ad928SBarry Smith } else if (side == PC_SYMMETRIC) { 7734b9ad928SBarry Smith /* There's an extra copy here; maybe should provide 2 work vectors instead? */ 7749566063dSJacob Faibussowitsch PetscCall(PCApplySymmetricRight(pc, x, work)); 7759566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, work, y)); 7769566063dSJacob Faibussowitsch PetscCall(VecCopy(y, work)); 7779566063dSJacob Faibussowitsch PetscCall(PCApplySymmetricLeft(pc, work, y)); 7784b9ad928SBarry Smith } 7794b9ad928SBarry Smith } 780e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE)); 7813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7824b9ad928SBarry Smith } 7834b9ad928SBarry Smith 7844b9ad928SBarry Smith /*@ 7854b9ad928SBarry Smith PCApplyBAorABTranspose - Applies the transpose of the preconditioner 786562efe2eSBarry Smith and operator to a vector. That is, applies $B^T * A^T$ with left preconditioning, 787562efe2eSBarry Smith NOT $(B*A)^T = A^T*B^T$. 7884b9ad928SBarry Smith 789c3339decSBarry Smith Collective 7904b9ad928SBarry Smith 7914b9ad928SBarry Smith Input Parameters: 7924b9ad928SBarry Smith + pc - the preconditioner context 793f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC` 7944b9ad928SBarry Smith . x - input vector 7954b9ad928SBarry Smith - work - work vector 7964b9ad928SBarry Smith 7974b9ad928SBarry Smith Output Parameter: 7984b9ad928SBarry Smith . y - output vector 7994b9ad928SBarry Smith 80020f4b53cSBarry Smith Level: developer 80120f4b53cSBarry Smith 802f1580f4eSBarry Smith Note: 803562efe2eSBarry 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 804562efe2eSBarry Smith defined by $B^T$. This is why this has the funny form that it computes $B^T * A^T$ 8059b3038f0SBarry Smith 806562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorAB()` 8074b9ad928SBarry Smith @*/ 808d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorABTranspose(PC pc, PCSide side, Vec x, Vec y, Vec work) 809d71ae5a4SJacob Faibussowitsch { 8104b9ad928SBarry Smith PetscFunctionBegin; 8110700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8120700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 8130700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 4); 8140700a824SBarry Smith PetscValidHeaderSpecific(work, VEC_CLASSID, 5); 8157827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 816e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE)); 8174b9ad928SBarry Smith if (pc->ops->applyBAtranspose) { 818dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applyBAtranspose, side, x, y, work); 819e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE)); 8203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8214b9ad928SBarry Smith } 8227827d75bSBarry Smith PetscCheck(side == PC_LEFT || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right or left"); 8234b9ad928SBarry Smith 8249566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 8254b9ad928SBarry Smith if (side == PC_RIGHT) { 8269566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose(pc, x, work)); 8279566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc->mat, work, y)); 8284b9ad928SBarry Smith } else if (side == PC_LEFT) { 8299566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc->mat, x, work)); 8309566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose(pc, work, y)); 8314b9ad928SBarry Smith } 8324b9ad928SBarry Smith /* add support for PC_SYMMETRIC */ 833e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE)); 8343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8354b9ad928SBarry Smith } 8364b9ad928SBarry Smith 8374b9ad928SBarry Smith /*@ 8384b9ad928SBarry Smith PCApplyRichardsonExists - Determines whether a particular preconditioner has a 8394b9ad928SBarry Smith built-in fast application of Richardson's method. 8404b9ad928SBarry Smith 8414b9ad928SBarry Smith Not Collective 8424b9ad928SBarry Smith 8434b9ad928SBarry Smith Input Parameter: 8444b9ad928SBarry Smith . pc - the preconditioner 8454b9ad928SBarry Smith 8464b9ad928SBarry Smith Output Parameter: 847f1580f4eSBarry Smith . exists - `PETSC_TRUE` or `PETSC_FALSE` 8484b9ad928SBarry Smith 8494b9ad928SBarry Smith Level: developer 8504b9ad928SBarry Smith 85139b1ba1bSPierre Jolivet .seealso: [](ch_ksp), `PC`, `KSPRICHARDSON`, `PCApplyRichardson()` 8524b9ad928SBarry Smith @*/ 853d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyRichardsonExists(PC pc, PetscBool *exists) 854d71ae5a4SJacob Faibussowitsch { 8554b9ad928SBarry Smith PetscFunctionBegin; 8560700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8574f572ea9SToby Isaac PetscAssertPointer(exists, 2); 8584b9ad928SBarry Smith if (pc->ops->applyrichardson) *exists = PETSC_TRUE; 8594b9ad928SBarry Smith else *exists = PETSC_FALSE; 8603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8614b9ad928SBarry Smith } 8624b9ad928SBarry Smith 8634b9ad928SBarry Smith /*@ 8644b9ad928SBarry Smith PCApplyRichardson - Applies several steps of Richardson iteration with 8654b9ad928SBarry Smith the particular preconditioner. This routine is usually used by the 8664b9ad928SBarry Smith Krylov solvers and not the application code directly. 8674b9ad928SBarry Smith 868c3339decSBarry Smith Collective 8694b9ad928SBarry Smith 8704b9ad928SBarry Smith Input Parameters: 8714b9ad928SBarry Smith + pc - the preconditioner context 872dd8e379bSPierre Jolivet . b - the right-hand side 8734b9ad928SBarry Smith . w - one work vector 8744b9ad928SBarry Smith . rtol - relative decrease in residual norm convergence criteria 87570441072SBarry Smith . abstol - absolute residual norm convergence criteria 8764b9ad928SBarry Smith . dtol - divergence residual norm increase criteria 8777319c654SBarry Smith . its - the number of iterations to apply. 8787319c654SBarry Smith - guesszero - if the input x contains nonzero initial guess 8794b9ad928SBarry Smith 880d8d19677SJose E. Roman Output Parameters: 8814d0a8057SBarry Smith + outits - number of iterations actually used (for SOR this always equals its) 8824d0a8057SBarry Smith . reason - the reason the apply terminated 883f1580f4eSBarry Smith - y - the solution (also contains initial guess if guesszero is `PETSC_FALSE` 8844b9ad928SBarry Smith 88520f4b53cSBarry Smith Level: developer 88620f4b53cSBarry Smith 8874b9ad928SBarry Smith Notes: 8884b9ad928SBarry Smith Most preconditioners do not support this function. Use the command 889f1580f4eSBarry Smith `PCApplyRichardsonExists()` to determine if one does. 8904b9ad928SBarry Smith 891f1580f4eSBarry Smith Except for the `PCMG` this routine ignores the convergence tolerances 8924b9ad928SBarry Smith and always runs for the number of iterations 8934b9ad928SBarry Smith 894562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyRichardsonExists()` 8954b9ad928SBarry Smith @*/ 896d71ae5a4SJacob 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) 897d71ae5a4SJacob Faibussowitsch { 8984b9ad928SBarry Smith PetscFunctionBegin; 8990700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9000700a824SBarry Smith PetscValidHeaderSpecific(b, VEC_CLASSID, 2); 9010700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 9020700a824SBarry Smith PetscValidHeaderSpecific(w, VEC_CLASSID, 4); 9037827d75bSBarry Smith PetscCheck(b != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "b and y must be different vectors"); 9049566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 905dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applyrichardson, b, y, w, rtol, abstol, dtol, its, guesszero, outits, reason); 9063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9074b9ad928SBarry Smith } 9084b9ad928SBarry Smith 909422a814eSBarry Smith /*@ 910f1580f4eSBarry Smith PCSetFailedReason - Sets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail 9111b2b9847SBarry Smith 912c3339decSBarry Smith Logically Collective 9131b2b9847SBarry Smith 914d8d19677SJose E. Roman Input Parameters: 9151b2b9847SBarry Smith + pc - the preconditioner context 9161b2b9847SBarry Smith - reason - the reason it failedx 9171b2b9847SBarry Smith 9181b2b9847SBarry Smith Level: advanced 9191b2b9847SBarry Smith 920562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCFailedReason` 9211b2b9847SBarry Smith @*/ 922d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFailedReason(PC pc, PCFailedReason reason) 923d71ae5a4SJacob Faibussowitsch { 9241b2b9847SBarry Smith PetscFunctionBegin; 9256479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9261b2b9847SBarry Smith pc->failedreason = reason; 9273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9281b2b9847SBarry Smith } 9291b2b9847SBarry Smith 9301b2b9847SBarry Smith /*@ 931f1580f4eSBarry Smith PCGetFailedReason - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail 932422a814eSBarry Smith 933c3339decSBarry Smith Logically Collective 934422a814eSBarry Smith 935422a814eSBarry Smith Input Parameter: 936422a814eSBarry Smith . pc - the preconditioner context 937422a814eSBarry Smith 938422a814eSBarry Smith Output Parameter: 9391b2b9847SBarry Smith . reason - the reason it failed 940422a814eSBarry Smith 941422a814eSBarry Smith Level: advanced 942422a814eSBarry Smith 943f1580f4eSBarry Smith Note: 944f1580f4eSBarry Smith This is the maximum over reason over all ranks in the PC communicator. It is only valid after 9456479e835SStefano Zampini a call `KSPCheckDot()` or `KSPCheckNorm()` inside a `KSPSolve()` or `PCReduceFailedReason()`. 9466479e835SStefano Zampini It is not valid immediately after a `PCSetUp()` or `PCApply()`, then use `PCGetFailedReasonRank()` 9471b2b9847SBarry Smith 948562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, ``PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReasonRank()`, `PCSetFailedReason()` 949422a814eSBarry Smith @*/ 950d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason) 951d71ae5a4SJacob Faibussowitsch { 952422a814eSBarry Smith PetscFunctionBegin; 9536479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9549754764cSHong Zhang if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled; 9550ae0b32dSHong Zhang else *reason = pc->failedreason; 9563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 957422a814eSBarry Smith } 958422a814eSBarry Smith 9591b2b9847SBarry Smith /*@ 960f1580f4eSBarry Smith PCGetFailedReasonRank - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail on this MPI rank 9611b2b9847SBarry Smith 962f1580f4eSBarry Smith Not Collective 9631b2b9847SBarry Smith 9641b2b9847SBarry Smith Input Parameter: 9651b2b9847SBarry Smith . pc - the preconditioner context 9661b2b9847SBarry Smith 9671b2b9847SBarry Smith Output Parameter: 9681b2b9847SBarry Smith . reason - the reason it failed 9691b2b9847SBarry Smith 97020f4b53cSBarry Smith Level: advanced 97120f4b53cSBarry Smith 972f1580f4eSBarry Smith Note: 973562efe2eSBarry Smith Different processes may have different reasons or no reason, see `PCGetFailedReason()` 9741b2b9847SBarry Smith 975562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`, `PCReduceFailedReason()` 9761b2b9847SBarry Smith @*/ 977d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReasonRank(PC pc, PCFailedReason *reason) 978d71ae5a4SJacob Faibussowitsch { 9791b2b9847SBarry Smith PetscFunctionBegin; 9806479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9811b2b9847SBarry Smith if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled; 9821b2b9847SBarry Smith else *reason = pc->failedreason; 9833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9841b2b9847SBarry Smith } 985422a814eSBarry Smith 9866479e835SStefano Zampini /*@ 9876479e835SStefano Zampini PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC` 9886479e835SStefano Zampini 9896479e835SStefano Zampini Collective 9906479e835SStefano Zampini 9916479e835SStefano Zampini Input Parameter: 9926479e835SStefano Zampini . pc - the preconditioner context 9936479e835SStefano Zampini 9946479e835SStefano Zampini Level: advanced 9956479e835SStefano Zampini 9966479e835SStefano Zampini Note: 997562efe2eSBarry Smith Different MPI processes may have different reasons or no reason, see `PCGetFailedReason()`. This routine 9986479e835SStefano Zampini makes them have a common value (failure if any MPI process had a failure). 9996479e835SStefano Zampini 1000562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()` 10016479e835SStefano Zampini @*/ 10026479e835SStefano Zampini PetscErrorCode PCReduceFailedReason(PC pc) 10036479e835SStefano Zampini { 10046479e835SStefano Zampini PetscInt buf; 10056479e835SStefano Zampini 10066479e835SStefano Zampini PetscFunctionBegin; 10076479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 10086479e835SStefano Zampini buf = (PetscInt)pc->failedreason; 10096479e835SStefano Zampini PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc))); 10106479e835SStefano Zampini pc->failedreason = (PCFailedReason)buf; 10116479e835SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 10126479e835SStefano Zampini } 10136479e835SStefano Zampini 1014c00cb57fSBarry Smith /* Next line needed to deactivate KSP_Solve logging */ 1015c00cb57fSBarry Smith #include <petsc/private/kspimpl.h> 1016c00cb57fSBarry Smith 10174b9ad928SBarry Smith /* 10184b9ad928SBarry Smith a setupcall of 0 indicates never setup, 101923ee1639SBarry Smith 1 indicates has been previously setup 1020422a814eSBarry Smith -1 indicates a PCSetUp() was attempted and failed 10214b9ad928SBarry Smith */ 10224b9ad928SBarry Smith /*@ 10234b9ad928SBarry Smith PCSetUp - Prepares for the use of a preconditioner. 10244b9ad928SBarry Smith 1025c3339decSBarry Smith Collective 10264b9ad928SBarry Smith 10274b9ad928SBarry Smith Input Parameter: 10284b9ad928SBarry Smith . pc - the preconditioner context 10294b9ad928SBarry Smith 10304b9ad928SBarry Smith Level: developer 10314b9ad928SBarry Smith 1032562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()` 10334b9ad928SBarry Smith @*/ 1034d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp(PC pc) 1035d71ae5a4SJacob Faibussowitsch { 1036566e8bf2SBarry Smith const char *def; 1037fc9b51d3SBarry Smith PetscObjectState matstate, matnonzerostate; 10384b9ad928SBarry Smith 10394b9ad928SBarry Smith PetscFunctionBegin; 10400700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 104128b400f6SJacob Faibussowitsch PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first"); 10426ce5e81cSLisandro Dalcin 104323ee1639SBarry Smith if (pc->setupcalled && pc->reusepreconditioner) { 10449566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n")); 10453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10464b9ad928SBarry Smith } 10474b9ad928SBarry Smith 10489566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate)); 10499566063dSJacob Faibussowitsch PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate)); 105023ee1639SBarry Smith if (!pc->setupcalled) { 10515e62d202SMark Adams //PetscCall(PetscInfo(pc, "Setting up PC for first time\n")); 105223ee1639SBarry Smith pc->flag = DIFFERENT_NONZERO_PATTERN; 10539df67409SStefano Zampini } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS); 10549df67409SStefano Zampini else { 10559df67409SStefano Zampini if (matnonzerostate != pc->matnonzerostate) { 10569566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n")); 105723ee1639SBarry Smith pc->flag = DIFFERENT_NONZERO_PATTERN; 105823ee1639SBarry Smith } else { 10595e62d202SMark Adams //PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n")); 106023ee1639SBarry Smith pc->flag = SAME_NONZERO_PATTERN; 106123ee1639SBarry Smith } 106223ee1639SBarry Smith } 106323ee1639SBarry Smith pc->matstate = matstate; 1064fc9b51d3SBarry Smith pc->matnonzerostate = matnonzerostate; 106523ee1639SBarry Smith 10667adad957SLisandro Dalcin if (!((PetscObject)pc)->type_name) { 10679566063dSJacob Faibussowitsch PetscCall(PCGetDefaultType_Private(pc, &def)); 10689566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, def)); 1069566e8bf2SBarry Smith } 10704b9ad928SBarry Smith 10719566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure)); 10729566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure)); 10739566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0)); 10744b9ad928SBarry Smith if (pc->ops->setup) { 1075c00cb57fSBarry Smith /* do not log solves and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */ 10769566063dSJacob Faibussowitsch PetscCall(KSPInitializePackage()); 10779566063dSJacob Faibussowitsch PetscCall(PetscLogEventDeactivatePush(KSP_Solve)); 10789566063dSJacob Faibussowitsch PetscCall(PetscLogEventDeactivatePush(PC_Apply)); 1079dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, setup); 10809566063dSJacob Faibussowitsch PetscCall(PetscLogEventDeactivatePop(KSP_Solve)); 10819566063dSJacob Faibussowitsch PetscCall(PetscLogEventDeactivatePop(PC_Apply)); 10824b9ad928SBarry Smith } 10839566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0)); 1084422a814eSBarry Smith if (!pc->setupcalled) pc->setupcalled = 1; 10853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10864b9ad928SBarry Smith } 10874b9ad928SBarry Smith 10884b9ad928SBarry Smith /*@ 10894b9ad928SBarry Smith PCSetUpOnBlocks - Sets up the preconditioner for each block in 10904b9ad928SBarry Smith the block Jacobi, block Gauss-Seidel, and overlapping Schwarz 10914b9ad928SBarry Smith methods. 10924b9ad928SBarry Smith 1093c3339decSBarry Smith Collective 10944b9ad928SBarry Smith 1095f1580f4eSBarry Smith Input Parameter: 10964b9ad928SBarry Smith . pc - the preconditioner context 10974b9ad928SBarry Smith 10984b9ad928SBarry Smith Level: developer 10994b9ad928SBarry Smith 1100f1580f4eSBarry Smith Note: 1101f1580f4eSBarry Smith For nested preconditioners such as `PCBJACOBI` `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is 1102f1580f4eSBarry Smith called on the outer `PC`, this routine ensures it is called. 1103f1580f4eSBarry Smith 1104562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()` 11054b9ad928SBarry Smith @*/ 1106d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUpOnBlocks(PC pc) 1107d71ae5a4SJacob Faibussowitsch { 11084b9ad928SBarry Smith PetscFunctionBegin; 11090700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 11103ba16761SJacob Faibussowitsch if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS); 11119566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0)); 1112dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, setuponblocks); 11139566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0)); 11143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11154b9ad928SBarry Smith } 11164b9ad928SBarry Smith 11174b9ad928SBarry Smith /*@C 11184b9ad928SBarry Smith PCSetModifySubMatrices - Sets a user-defined routine for modifying the 111904c3f3b8SBarry Smith submatrices that arise within certain subdomain-based preconditioners such as `PCASM` 11204b9ad928SBarry Smith 1121c3339decSBarry Smith Logically Collective 11224b9ad928SBarry Smith 11234b9ad928SBarry Smith Input Parameters: 11244b9ad928SBarry Smith + pc - the preconditioner context 11254b9ad928SBarry Smith . func - routine for modifying the submatrices 112604c3f3b8SBarry Smith - ctx - optional user-defined context (may be `NULL`) 11274b9ad928SBarry Smith 112820f4b53cSBarry Smith Calling sequence of `func`: 112920f4b53cSBarry Smith + pc - the preconditioner context 113004c3f3b8SBarry Smith . nsub - number of index sets 113120f4b53cSBarry Smith . row - an array of index sets that contain the global row numbers 11324b9ad928SBarry Smith that comprise each local submatrix 11334b9ad928SBarry Smith . col - an array of index sets that contain the global column numbers 11344b9ad928SBarry Smith that comprise each local submatrix 11354b9ad928SBarry Smith . submat - array of local submatrices 11364b9ad928SBarry Smith - ctx - optional user-defined context for private data for the 113704c3f3b8SBarry Smith user-defined func routine (may be `NULL`) 11384b9ad928SBarry Smith 113920f4b53cSBarry Smith Level: advanced 114020f4b53cSBarry Smith 11414b9ad928SBarry Smith Notes: 114204c3f3b8SBarry Smith The basic submatrices are extracted from the preconditioner matrix as 114304c3f3b8SBarry Smith usual; the user can then alter these (for example, to set different boundary 114404c3f3b8SBarry Smith conditions for each submatrix) before they are used for the local solves. 114504c3f3b8SBarry Smith 1146f1580f4eSBarry Smith `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and 1147f1580f4eSBarry Smith `KSPSolve()`. 11484b9ad928SBarry Smith 1149f1580f4eSBarry Smith A routine set by `PCSetModifySubMatrices()` is currently called within 1150f1580f4eSBarry Smith the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`) 11514b9ad928SBarry Smith preconditioners. All other preconditioners ignore this routine. 11524b9ad928SBarry Smith 1153562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()` 11544b9ad928SBarry Smith @*/ 115504c3f3b8SBarry Smith PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx), void *ctx) 1156d71ae5a4SJacob Faibussowitsch { 11574b9ad928SBarry Smith PetscFunctionBegin; 11580700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 11594b9ad928SBarry Smith pc->modifysubmatrices = func; 11604b9ad928SBarry Smith pc->modifysubmatricesP = ctx; 11613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11624b9ad928SBarry Smith } 11634b9ad928SBarry Smith 11644b9ad928SBarry Smith /*@C 11654b9ad928SBarry Smith PCModifySubMatrices - Calls an optional user-defined routine within 1166f1580f4eSBarry Smith certain preconditioners if one has been set with `PCSetModifySubMatrices()`. 11674b9ad928SBarry Smith 1168c3339decSBarry Smith Collective 11694b9ad928SBarry Smith 11704b9ad928SBarry Smith Input Parameters: 11714b9ad928SBarry Smith + pc - the preconditioner context 11724b9ad928SBarry Smith . nsub - the number of local submatrices 11734b9ad928SBarry Smith . row - an array of index sets that contain the global row numbers 11744b9ad928SBarry Smith that comprise each local submatrix 11754b9ad928SBarry Smith . col - an array of index sets that contain the global column numbers 11764b9ad928SBarry Smith that comprise each local submatrix 11774b9ad928SBarry Smith . submat - array of local submatrices 11784b9ad928SBarry Smith - ctx - optional user-defined context for private data for the 1179562efe2eSBarry Smith user-defined routine (may be `NULL`) 11804b9ad928SBarry Smith 11814b9ad928SBarry Smith Output Parameter: 11824b9ad928SBarry Smith . submat - array of local submatrices (the entries of which may 11834b9ad928SBarry Smith have been modified) 11844b9ad928SBarry Smith 118520f4b53cSBarry Smith Level: developer 118620f4b53cSBarry Smith 118704c3f3b8SBarry Smith Note: 11884b9ad928SBarry Smith The user should NOT generally call this routine, as it will 118904c3f3b8SBarry Smith automatically be called within certain preconditioners. 11904b9ad928SBarry Smith 1191562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetModifySubMatrices()` 11924b9ad928SBarry Smith @*/ 1193d71ae5a4SJacob Faibussowitsch PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx) 1194d71ae5a4SJacob Faibussowitsch { 11954b9ad928SBarry Smith PetscFunctionBegin; 11960700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 11973ba16761SJacob Faibussowitsch if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS); 11989566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0)); 11999566063dSJacob Faibussowitsch PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx)); 12009566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0)); 12013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12024b9ad928SBarry Smith } 12034b9ad928SBarry Smith 12044b9ad928SBarry Smith /*@ 12054b9ad928SBarry Smith PCSetOperators - Sets the matrix associated with the linear system and 12064b9ad928SBarry Smith a (possibly) different one associated with the preconditioner. 12074b9ad928SBarry Smith 1208c3339decSBarry Smith Logically Collective 12094b9ad928SBarry Smith 12104b9ad928SBarry Smith Input Parameters: 12114b9ad928SBarry Smith + pc - the preconditioner context 1212e5d3d808SBarry Smith . Amat - the matrix that defines the linear system 1213d1e9a80fSBarry Smith - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat. 12144b9ad928SBarry Smith 121520f4b53cSBarry Smith Level: intermediate 1216189c0b8aSBarry Smith 121720f4b53cSBarry Smith Notes: 121820f4b53cSBarry Smith Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used. 121920f4b53cSBarry Smith 122020f4b53cSBarry Smith If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then 1221f1580f4eSBarry Smith first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()` 1222f1580f4eSBarry Smith on it and then pass it back in in your call to `KSPSetOperators()`. 1223189c0b8aSBarry Smith 12244b9ad928SBarry Smith More Notes about Repeated Solution of Linear Systems: 122520f4b53cSBarry Smith PETSc does NOT reset the matrix entries of either `Amat` or `Pmat` 12264b9ad928SBarry Smith to zero after a linear solve; the user is completely responsible for 1227f1580f4eSBarry Smith matrix assembly. See the routine `MatZeroEntries()` if desiring to 12284b9ad928SBarry Smith zero all elements of a matrix. 12294b9ad928SBarry Smith 1230562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()` 12314b9ad928SBarry Smith @*/ 1232d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat) 1233d71ae5a4SJacob Faibussowitsch { 12343b3f6333SBarry Smith PetscInt m1, n1, m2, n2; 12354b9ad928SBarry Smith 12364b9ad928SBarry Smith PetscFunctionBegin; 12370700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 12380700a824SBarry Smith if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 12390700a824SBarry Smith if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3); 12403fc8bf9cSBarry Smith if (Amat) PetscCheckSameComm(pc, 1, Amat, 2); 12413fc8bf9cSBarry Smith if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3); 124231641f1aSBarry Smith if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) { 12439566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Amat, &m1, &n1)); 12449566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->mat, &m2, &n2)); 124563a3b9bcSJacob 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); 12469566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Pmat, &m1, &n1)); 12479566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2)); 124863a3b9bcSJacob 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); 12493b3f6333SBarry Smith } 12504b9ad928SBarry Smith 125123ee1639SBarry Smith if (Pmat != pc->pmat) { 125223ee1639SBarry Smith /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */ 125323ee1639SBarry Smith pc->matnonzerostate = -1; 125423ee1639SBarry Smith pc->matstate = -1; 125523ee1639SBarry Smith } 125623ee1639SBarry Smith 1257906ed7ccSBarry Smith /* reference first in case the matrices are the same */ 12589566063dSJacob Faibussowitsch if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat)); 12599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->mat)); 12609566063dSJacob Faibussowitsch if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat)); 12619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->pmat)); 12624b9ad928SBarry Smith pc->mat = Amat; 12634b9ad928SBarry Smith pc->pmat = Pmat; 12643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1265e56f5c9eSBarry Smith } 1266e56f5c9eSBarry Smith 126723ee1639SBarry Smith /*@ 126823ee1639SBarry Smith PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed. 126923ee1639SBarry Smith 1270c3339decSBarry Smith Logically Collective 127123ee1639SBarry Smith 127223ee1639SBarry Smith Input Parameters: 127323ee1639SBarry Smith + pc - the preconditioner context 1274f1580f4eSBarry Smith - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner 127523ee1639SBarry Smith 12762b26979fSBarry Smith Level: intermediate 12772b26979fSBarry Smith 1278f1580f4eSBarry Smith Note: 1279f1580f4eSBarry Smith Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option 1280f1580f4eSBarry Smith prevents this. 1281f1580f4eSBarry Smith 1282562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()` 128323ee1639SBarry Smith @*/ 1284d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag) 1285d71ae5a4SJacob Faibussowitsch { 128623ee1639SBarry Smith PetscFunctionBegin; 128723ee1639SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1288f9177522SStefano Zampini PetscValidLogicalCollectiveBool(pc, flag, 2); 128923ee1639SBarry Smith pc->reusepreconditioner = flag; 12903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12914b9ad928SBarry Smith } 12924b9ad928SBarry Smith 1293c60c7ad4SBarry Smith /*@ 1294f1580f4eSBarry Smith PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed. 1295c60c7ad4SBarry Smith 1296bba28a21SBarry Smith Not Collective 1297c60c7ad4SBarry Smith 1298c60c7ad4SBarry Smith Input Parameter: 1299c60c7ad4SBarry Smith . pc - the preconditioner context 1300c60c7ad4SBarry Smith 1301c60c7ad4SBarry Smith Output Parameter: 1302f1580f4eSBarry Smith . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner 1303c60c7ad4SBarry Smith 1304d0418729SSatish Balay Level: intermediate 1305d0418729SSatish Balay 1306562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()` 1307c60c7ad4SBarry Smith @*/ 1308d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag) 1309d71ae5a4SJacob Faibussowitsch { 1310c60c7ad4SBarry Smith PetscFunctionBegin; 1311c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 13124f572ea9SToby Isaac PetscAssertPointer(flag, 2); 1313c60c7ad4SBarry Smith *flag = pc->reusepreconditioner; 13143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1315c60c7ad4SBarry Smith } 1316c60c7ad4SBarry Smith 1317487a658cSBarry Smith /*@ 13184b9ad928SBarry Smith PCGetOperators - Gets the matrix associated with the linear system and 13194b9ad928SBarry Smith possibly a different one associated with the preconditioner. 13204b9ad928SBarry Smith 1321562efe2eSBarry Smith Not Collective, though parallel `Mat`s are returned if `pc` is parallel 13224b9ad928SBarry Smith 13234b9ad928SBarry Smith Input Parameter: 13244b9ad928SBarry Smith . pc - the preconditioner context 13254b9ad928SBarry Smith 13264b9ad928SBarry Smith Output Parameters: 1327e5d3d808SBarry Smith + Amat - the matrix defining the linear system 132823ee1639SBarry Smith - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat. 13294b9ad928SBarry Smith 13304b9ad928SBarry Smith Level: intermediate 13314b9ad928SBarry Smith 1332f1580f4eSBarry Smith Note: 133395452b02SPatrick Sanan Does not increase the reference count of the matrices, so you should not destroy them 1334298cc208SBarry Smith 1335f1580f4eSBarry Smith Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators 1336f1580f4eSBarry Smith are created in `PC` and returned to the user. In this case, if both operators 133773950996SBarry Smith mat and pmat are requested, two DIFFERENT operators will be returned. If 133873950996SBarry Smith only one is requested both operators in the PC will be the same (i.e. as 1339f1580f4eSBarry Smith if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats). 134073950996SBarry Smith The user must set the sizes of the returned matrices and their type etc just 1341f1580f4eSBarry Smith as if the user created them with `MatCreate()`. For example, 134273950996SBarry Smith 1343f1580f4eSBarry Smith .vb 1344f1580f4eSBarry Smith KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to 1345f1580f4eSBarry Smith set size, type, etc of Amat 134673950996SBarry Smith 1347f1580f4eSBarry Smith MatCreate(comm,&mat); 1348f1580f4eSBarry Smith KSP/PCSetOperators(ksp/pc,Amat,Amat); 1349f1580f4eSBarry Smith PetscObjectDereference((PetscObject)mat); 1350f1580f4eSBarry Smith set size, type, etc of Amat 1351f1580f4eSBarry Smith .ve 135273950996SBarry Smith 135373950996SBarry Smith and 135473950996SBarry Smith 1355f1580f4eSBarry Smith .vb 1356f1580f4eSBarry Smith KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to 1357f1580f4eSBarry Smith set size, type, etc of Amat and Pmat 135873950996SBarry Smith 1359f1580f4eSBarry Smith MatCreate(comm,&Amat); 1360f1580f4eSBarry Smith MatCreate(comm,&Pmat); 1361f1580f4eSBarry Smith KSP/PCSetOperators(ksp/pc,Amat,Pmat); 1362f1580f4eSBarry Smith PetscObjectDereference((PetscObject)Amat); 1363f1580f4eSBarry Smith PetscObjectDereference((PetscObject)Pmat); 1364f1580f4eSBarry Smith set size, type, etc of Amat and Pmat 1365f1580f4eSBarry Smith .ve 136673950996SBarry Smith 1367f1580f4eSBarry Smith The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy 1368b8d709abSRichard Tran Mills of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely 1369f1580f4eSBarry Smith managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look 1370f1580f4eSBarry Smith at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to 1371f1580f4eSBarry Smith the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP 1372f1580f4eSBarry Smith you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you). 1373f1580f4eSBarry Smith Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when 137473950996SBarry Smith it can be created for you? 137573950996SBarry Smith 1376562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()` 13774b9ad928SBarry Smith @*/ 1378d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat) 1379d71ae5a4SJacob Faibussowitsch { 13804b9ad928SBarry Smith PetscFunctionBegin; 13810700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1382e5d3d808SBarry Smith if (Amat) { 138373950996SBarry Smith if (!pc->mat) { 13849fca8c71SStefano Zampini if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */ 13859a4708feSJed Brown pc->mat = pc->pmat; 13869566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->mat)); 1387e5d3d808SBarry Smith } else { /* both Amat and Pmat are empty */ 13889566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat)); 1389e5d3d808SBarry Smith if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */ 139073950996SBarry Smith pc->pmat = pc->mat; 13919566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->pmat)); 139273950996SBarry Smith } 139373950996SBarry Smith } 13949a4708feSJed Brown } 1395e5d3d808SBarry Smith *Amat = pc->mat; 139673950996SBarry Smith } 1397e5d3d808SBarry Smith if (Pmat) { 139873950996SBarry Smith if (!pc->pmat) { 1399e5d3d808SBarry Smith if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */ 14009a4708feSJed Brown pc->pmat = pc->mat; 14019566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->pmat)); 14029a4708feSJed Brown } else { 14039566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat)); 1404e5d3d808SBarry Smith if (!Amat) { /* user did NOT request Amat, so make same as Pmat */ 140573950996SBarry Smith pc->mat = pc->pmat; 14069566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->mat)); 140773950996SBarry Smith } 140873950996SBarry Smith } 14099a4708feSJed Brown } 1410e5d3d808SBarry Smith *Pmat = pc->pmat; 141173950996SBarry Smith } 14123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14134b9ad928SBarry Smith } 14144b9ad928SBarry Smith 1415906ed7ccSBarry Smith /*@C 1416906ed7ccSBarry Smith PCGetOperatorsSet - Determines if the matrix associated with the linear system and 1417f1580f4eSBarry Smith possibly a different one associated with the preconditioner have been set in the `PC`. 1418906ed7ccSBarry Smith 141920f4b53cSBarry Smith Not Collective, though the results on all processes should be the same 1420906ed7ccSBarry Smith 1421906ed7ccSBarry Smith Input Parameter: 1422906ed7ccSBarry Smith . pc - the preconditioner context 1423906ed7ccSBarry Smith 1424906ed7ccSBarry Smith Output Parameters: 1425906ed7ccSBarry Smith + mat - the matrix associated with the linear system was set 1426906ed7ccSBarry Smith - pmat - matrix associated with the preconditioner was set, usually the same 1427906ed7ccSBarry Smith 1428906ed7ccSBarry Smith Level: intermediate 1429906ed7ccSBarry Smith 1430562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()` 1431906ed7ccSBarry Smith @*/ 1432d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat) 1433d71ae5a4SJacob Faibussowitsch { 1434906ed7ccSBarry Smith PetscFunctionBegin; 14350700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1436906ed7ccSBarry Smith if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE; 1437906ed7ccSBarry Smith if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE; 14383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1439906ed7ccSBarry Smith } 1440906ed7ccSBarry Smith 1441f39d8e23SSatish Balay /*@ 1442a4fd02acSBarry Smith PCFactorGetMatrix - Gets the factored matrix from the 1443f1580f4eSBarry Smith preconditioner context. This routine is valid only for the `PCLU`, 1444f1580f4eSBarry Smith `PCILU`, `PCCHOLESKY`, and `PCICC` methods. 14454b9ad928SBarry Smith 1446562efe2eSBarry Smith Not Collective though `mat` is parallel if `pc` is parallel 14474b9ad928SBarry Smith 1448f1580f4eSBarry Smith Input Parameter: 14494b9ad928SBarry Smith . pc - the preconditioner context 14504b9ad928SBarry Smith 1451feefa0e1SJacob Faibussowitsch Output Parameters: 14524b9ad928SBarry Smith . mat - the factored matrix 14534b9ad928SBarry Smith 14544b9ad928SBarry Smith Level: advanced 14554b9ad928SBarry Smith 1456f1580f4eSBarry Smith Note: 1457562efe2eSBarry Smith Does not increase the reference count for `mat` so DO NOT destroy it 14589405f653SBarry Smith 1459562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC` 14604b9ad928SBarry Smith @*/ 1461d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat) 1462d71ae5a4SJacob Faibussowitsch { 14634b9ad928SBarry Smith PetscFunctionBegin; 14640700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 14654f572ea9SToby Isaac PetscAssertPointer(mat, 2); 14667def7855SStefano Zampini PetscCall(PCFactorSetUpMatSolverType(pc)); 1467dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, getfactoredmatrix, mat); 14683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14694b9ad928SBarry Smith } 14704b9ad928SBarry Smith 1471cc4c1da9SBarry Smith /*@ 14724b9ad928SBarry Smith PCSetOptionsPrefix - Sets the prefix used for searching for all 1473f1580f4eSBarry Smith `PC` options in the database. 14744b9ad928SBarry Smith 1475c3339decSBarry Smith Logically Collective 14764b9ad928SBarry Smith 14774b9ad928SBarry Smith Input Parameters: 14784b9ad928SBarry Smith + pc - the preconditioner context 1479f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests 14804b9ad928SBarry Smith 1481f1580f4eSBarry Smith Note: 14824b9ad928SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 14834b9ad928SBarry Smith The first character of all runtime options is AUTOMATICALLY the 14844b9ad928SBarry Smith hyphen. 14854b9ad928SBarry Smith 14864b9ad928SBarry Smith Level: advanced 14874b9ad928SBarry Smith 1488562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()` 14894b9ad928SBarry Smith @*/ 1490d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[]) 1491d71ae5a4SJacob Faibussowitsch { 14924b9ad928SBarry Smith PetscFunctionBegin; 14930700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 14949566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix)); 14953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14964b9ad928SBarry Smith } 14974b9ad928SBarry Smith 1498cc4c1da9SBarry Smith /*@ 14994b9ad928SBarry Smith PCAppendOptionsPrefix - Appends to the prefix used for searching for all 1500f1580f4eSBarry Smith `PC` options in the database. 15014b9ad928SBarry Smith 1502c3339decSBarry Smith Logically Collective 15034b9ad928SBarry Smith 15044b9ad928SBarry Smith Input Parameters: 15054b9ad928SBarry Smith + pc - the preconditioner context 1506f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests 15074b9ad928SBarry Smith 1508f1580f4eSBarry Smith Note: 15094b9ad928SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 15104b9ad928SBarry Smith The first character of all runtime options is AUTOMATICALLY the 15114b9ad928SBarry Smith hyphen. 15124b9ad928SBarry Smith 15134b9ad928SBarry Smith Level: advanced 15144b9ad928SBarry Smith 1515562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()` 15164b9ad928SBarry Smith @*/ 1517d71ae5a4SJacob Faibussowitsch PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[]) 1518d71ae5a4SJacob Faibussowitsch { 15194b9ad928SBarry Smith PetscFunctionBegin; 15200700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15219566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix)); 15223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15234b9ad928SBarry Smith } 15244b9ad928SBarry Smith 1525cc4c1da9SBarry Smith /*@ 15264b9ad928SBarry Smith PCGetOptionsPrefix - Gets the prefix used for searching for all 15274b9ad928SBarry Smith PC options in the database. 15284b9ad928SBarry Smith 15294b9ad928SBarry Smith Not Collective 15304b9ad928SBarry Smith 1531f1580f4eSBarry Smith Input Parameter: 15324b9ad928SBarry Smith . pc - the preconditioner context 15334b9ad928SBarry Smith 1534f1580f4eSBarry Smith Output Parameter: 15354b9ad928SBarry Smith . prefix - pointer to the prefix string used, is returned 15364b9ad928SBarry Smith 15374b9ad928SBarry Smith Level: advanced 15384b9ad928SBarry Smith 1539562efe2eSBarry Smith Fortran Note: 1540562efe2eSBarry Smith The user should pass in a string `prefix` of 1541562efe2eSBarry Smith sufficient length to hold the prefix. 1542562efe2eSBarry Smith 1543562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()` 15444b9ad928SBarry Smith @*/ 1545d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[]) 1546d71ae5a4SJacob Faibussowitsch { 15474b9ad928SBarry Smith PetscFunctionBegin; 15480700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15494f572ea9SToby Isaac PetscAssertPointer(prefix, 2); 15509566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix)); 15513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15524b9ad928SBarry Smith } 15534b9ad928SBarry Smith 15548066bbecSBarry Smith /* 1555dd8e379bSPierre Jolivet Indicates the right-hand side will be changed by KSPSolve(), this occurs for a few 15568066bbecSBarry Smith preconditioners including BDDC and Eisentat that transform the equations before applying 15578066bbecSBarry Smith the Krylov methods 15588066bbecSBarry Smith */ 1559d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change) 1560d71ae5a4SJacob Faibussowitsch { 1561a06fd7f2SStefano Zampini PetscFunctionBegin; 1562a06fd7f2SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15634f572ea9SToby Isaac PetscAssertPointer(change, 2); 1564a06fd7f2SStefano Zampini *change = PETSC_FALSE; 1565cac4c232SBarry Smith PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change)); 15663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1567a06fd7f2SStefano Zampini } 1568a06fd7f2SStefano Zampini 15694b9ad928SBarry Smith /*@ 15704b9ad928SBarry Smith PCPreSolve - Optional pre-solve phase, intended for any 15714b9ad928SBarry Smith preconditioner-specific actions that must be performed before 15724b9ad928SBarry Smith the iterative solve itself. 15734b9ad928SBarry Smith 1574c3339decSBarry Smith Collective 15754b9ad928SBarry Smith 15764b9ad928SBarry Smith Input Parameters: 15774b9ad928SBarry Smith + pc - the preconditioner context 15784b9ad928SBarry Smith - ksp - the Krylov subspace context 15794b9ad928SBarry Smith 15804b9ad928SBarry Smith Level: developer 15814b9ad928SBarry Smith 1582feefa0e1SJacob Faibussowitsch Example Usage: 15834b9ad928SBarry Smith .vb 15844b9ad928SBarry Smith PCPreSolve(pc,ksp); 158523ce1328SBarry Smith KSPSolve(ksp,b,x); 15864b9ad928SBarry Smith PCPostSolve(pc,ksp); 15874b9ad928SBarry Smith .ve 15884b9ad928SBarry Smith 15894b9ad928SBarry Smith Notes: 1590f1580f4eSBarry Smith The pre-solve phase is distinct from the `PCSetUp()` phase. 15914b9ad928SBarry Smith 1592f1580f4eSBarry Smith `KSPSolve()` calls this directly, so is rarely called by the user. 15934b9ad928SBarry Smith 1594562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCPostSolve()` 15954b9ad928SBarry Smith @*/ 1596d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPreSolve(PC pc, KSP ksp) 1597d71ae5a4SJacob Faibussowitsch { 15984b9ad928SBarry Smith Vec x, rhs; 15994b9ad928SBarry Smith 16004b9ad928SBarry Smith PetscFunctionBegin; 16010700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16020700a824SBarry Smith PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2); 1603d9a03883SBarry Smith pc->presolvedone++; 16047827d75bSBarry Smith PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice"); 16059566063dSJacob Faibussowitsch PetscCall(KSPGetSolution(ksp, &x)); 16069566063dSJacob Faibussowitsch PetscCall(KSPGetRhs(ksp, &rhs)); 16074b9ad928SBarry Smith 1608dbbe0bcdSBarry Smith if (pc->ops->presolve) PetscUseTypeMethod(pc, presolve, ksp, rhs, x); 1609f4f49eeaSPierre Jolivet else if (pc->presolve) PetscCall(pc->presolve(pc, ksp)); 16103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16114b9ad928SBarry Smith } 16124b9ad928SBarry Smith 1613f560b561SHong Zhang /*@C 1614f1580f4eSBarry Smith PCSetPreSolve - Sets function used by `PCPreSolve()` which is intended for any 1615f560b561SHong Zhang preconditioner-specific actions that must be performed before 1616f560b561SHong Zhang the iterative solve itself. 1617f560b561SHong Zhang 1618c3339decSBarry Smith Logically Collective 1619f560b561SHong Zhang 1620f560b561SHong Zhang Input Parameters: 1621f560b561SHong Zhang + pc - the preconditioner object 1622f560b561SHong Zhang - presolve - the function to call before the solve 1623f560b561SHong Zhang 162420f4b53cSBarry Smith Calling sequence of `presolve`: 162520f4b53cSBarry Smith + pc - the `PC` context 162620f4b53cSBarry Smith - ksp - the `KSP` context 1627f560b561SHong Zhang 1628f560b561SHong Zhang Level: developer 1629f560b561SHong Zhang 1630562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCPreSolve()` 1631f560b561SHong Zhang @*/ 163204c3f3b8SBarry Smith PetscErrorCode PCSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC pc, KSP ksp)) 1633d71ae5a4SJacob Faibussowitsch { 1634f560b561SHong Zhang PetscFunctionBegin; 1635f560b561SHong Zhang PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1636f560b561SHong Zhang pc->presolve = presolve; 16373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1638f560b561SHong Zhang } 1639f560b561SHong Zhang 16404b9ad928SBarry Smith /*@ 16414b9ad928SBarry Smith PCPostSolve - Optional post-solve phase, intended for any 16424b9ad928SBarry Smith preconditioner-specific actions that must be performed after 16434b9ad928SBarry Smith the iterative solve itself. 16444b9ad928SBarry Smith 1645c3339decSBarry Smith Collective 16464b9ad928SBarry Smith 16474b9ad928SBarry Smith Input Parameters: 16484b9ad928SBarry Smith + pc - the preconditioner context 16494b9ad928SBarry Smith - ksp - the Krylov subspace context 16504b9ad928SBarry Smith 1651feefa0e1SJacob Faibussowitsch Example Usage: 16524b9ad928SBarry Smith .vb 16534b9ad928SBarry Smith PCPreSolve(pc,ksp); 165423ce1328SBarry Smith KSPSolve(ksp,b,x); 16554b9ad928SBarry Smith PCPostSolve(pc,ksp); 16564b9ad928SBarry Smith .ve 16574b9ad928SBarry Smith 1658562efe2eSBarry Smith Level: developer 1659562efe2eSBarry Smith 16604b9ad928SBarry Smith Note: 1661f1580f4eSBarry Smith `KSPSolve()` calls this routine directly, so it is rarely called by the user. 16624b9ad928SBarry Smith 1663562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetPostSolve()`, `PCSetPresolve()`, `PCPreSolve()`, `KSPSolve()` 16644b9ad928SBarry Smith @*/ 1665d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPostSolve(PC pc, KSP ksp) 1666d71ae5a4SJacob Faibussowitsch { 16674b9ad928SBarry Smith Vec x, rhs; 16684b9ad928SBarry Smith 16694b9ad928SBarry Smith PetscFunctionBegin; 16700700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16710700a824SBarry Smith PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2); 1672d9a03883SBarry Smith pc->presolvedone--; 16739566063dSJacob Faibussowitsch PetscCall(KSPGetSolution(ksp, &x)); 16749566063dSJacob Faibussowitsch PetscCall(KSPGetRhs(ksp, &rhs)); 1675dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, postsolve, ksp, rhs, x); 16763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16774b9ad928SBarry Smith } 16784b9ad928SBarry Smith 1679*ffeef943SBarry Smith /*@ 1680f1580f4eSBarry Smith PCLoad - Loads a `PC` that has been stored in binary with `PCView()`. 168155849f57SBarry Smith 1682c3339decSBarry Smith Collective 168355849f57SBarry Smith 168455849f57SBarry Smith Input Parameters: 1685f1580f4eSBarry Smith + newdm - the newly loaded `PC`, this needs to have been created with `PCCreate()` or 1686f1580f4eSBarry Smith some related function before a call to `PCLoad()`. 1687f1580f4eSBarry Smith - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()` 168855849f57SBarry Smith 168955849f57SBarry Smith Level: intermediate 169055849f57SBarry Smith 1691f1580f4eSBarry Smith Note: 1692562efe2eSBarry Smith The type is determined by the data in the file, any `PCType` set into the `PC` before this call is ignored. 169355849f57SBarry Smith 1694562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()` 169555849f57SBarry Smith @*/ 1696d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLoad(PC newdm, PetscViewer viewer) 1697d71ae5a4SJacob Faibussowitsch { 169855849f57SBarry Smith PetscBool isbinary; 1699060da220SMatthew G. Knepley PetscInt classid; 170055849f57SBarry Smith char type[256]; 170155849f57SBarry Smith 170255849f57SBarry Smith PetscFunctionBegin; 170355849f57SBarry Smith PetscValidHeaderSpecific(newdm, PC_CLASSID, 1); 170455849f57SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 17059566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 170628b400f6SJacob Faibussowitsch PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 170755849f57SBarry Smith 17089566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT)); 17097827d75bSBarry Smith PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file"); 17109566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR)); 17119566063dSJacob Faibussowitsch PetscCall(PCSetType(newdm, type)); 1712dbbe0bcdSBarry Smith PetscTryTypeMethod(newdm, load, viewer); 17133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 171455849f57SBarry Smith } 171555849f57SBarry Smith 17169804daf3SBarry Smith #include <petscdraw.h> 1717e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1718e04113cfSBarry Smith #include <petscviewersaws.h> 17190acecf5bSBarry Smith #endif 1720fe2efc57SMark 1721*ffeef943SBarry Smith /*@ 1722562efe2eSBarry Smith PCViewFromOptions - View from the `PC` based on options in the options database 1723fe2efc57SMark 1724c3339decSBarry Smith Collective 1725fe2efc57SMark 1726fe2efc57SMark Input Parameters: 172720f4b53cSBarry Smith + A - the `PC` context 1728f1580f4eSBarry Smith . obj - Optional object that provides the options prefix 1729736c3998SJose E. Roman - name - command line option 1730fe2efc57SMark 1731fe2efc57SMark Level: intermediate 1732f1580f4eSBarry Smith 1733562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()` 1734fe2efc57SMark @*/ 1735d71ae5a4SJacob Faibussowitsch PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[]) 1736d71ae5a4SJacob Faibussowitsch { 1737fe2efc57SMark PetscFunctionBegin; 1738fe2efc57SMark PetscValidHeaderSpecific(A, PC_CLASSID, 1); 17399566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 17403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1741fe2efc57SMark } 1742fe2efc57SMark 1743*ffeef943SBarry Smith /*@ 1744f1580f4eSBarry Smith PCView - Prints information about the `PC` 17454b9ad928SBarry Smith 1746c3339decSBarry Smith Collective 17474b9ad928SBarry Smith 17484b9ad928SBarry Smith Input Parameters: 1749feefa0e1SJacob Faibussowitsch + pc - the `PC` context 17504b9ad928SBarry Smith - viewer - optional visualization context 17514b9ad928SBarry Smith 175220f4b53cSBarry Smith Level: developer 175320f4b53cSBarry Smith 1754f1580f4eSBarry Smith Notes: 17554b9ad928SBarry Smith The available visualization contexts include 1756f1580f4eSBarry Smith + `PETSC_VIEWER_STDOUT_SELF` - standard output (default) 1757f1580f4eSBarry Smith - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard 17584b9ad928SBarry Smith output where only the first processor opens 17594b9ad928SBarry Smith the file. All other processors send their 17604b9ad928SBarry Smith data to the first processor to print. 17614b9ad928SBarry Smith 17624b9ad928SBarry Smith The user can open an alternative visualization contexts with 1763f1580f4eSBarry Smith `PetscViewerASCIIOpen()` (output to a specified file). 17644b9ad928SBarry Smith 1765562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewer`, `KSPView()`, `PetscViewerASCIIOpen()` 17664b9ad928SBarry Smith @*/ 1767d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView(PC pc, PetscViewer viewer) 1768d71ae5a4SJacob Faibussowitsch { 176919fd82e9SBarry Smith PCType cstr; 17706cd81132SPierre Jolivet PetscViewerFormat format; 17716cd81132SPierre Jolivet PetscBool iascii, isstring, isbinary, isdraw, pop = PETSC_FALSE; 1772e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1773536b137fSBarry Smith PetscBool issaws; 17740acecf5bSBarry Smith #endif 17754b9ad928SBarry Smith 17764b9ad928SBarry Smith PetscFunctionBegin; 17770700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 177848a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer)); 17790700a824SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1780c9780b6fSBarry Smith PetscCheckSameComm(pc, 1, viewer, 2); 17814b9ad928SBarry Smith 17829566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 17839566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 17849566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 17859566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 1786e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 17879566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws)); 17880acecf5bSBarry Smith #endif 1789219b1045SBarry Smith 179032077d6dSBarry Smith if (iascii) { 17919566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer)); 179248a46eb9SPierre Jolivet if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, " PC has not been set up so information may be incomplete\n")); 17939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1794dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 17959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1796834dbeb0SBarry Smith if (pc->mat) { 17976cd81132SPierre Jolivet PetscCall(PetscViewerGetFormat(viewer, &format)); 17986cd81132SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 17999566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 18006cd81132SPierre Jolivet pop = PETSC_TRUE; 18016cd81132SPierre Jolivet } 18024b9ad928SBarry Smith if (pc->pmat == pc->mat) { 18039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix = precond matrix:\n")); 18049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 18059566063dSJacob Faibussowitsch PetscCall(MatView(pc->mat, viewer)); 18069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 18074b9ad928SBarry Smith } else { 1808834dbeb0SBarry Smith if (pc->pmat) { 18099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix followed by preconditioner matrix:\n")); 18104b9ad928SBarry Smith } else { 18119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix:\n")); 18124b9ad928SBarry Smith } 18139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 18149566063dSJacob Faibussowitsch PetscCall(MatView(pc->mat, viewer)); 18159566063dSJacob Faibussowitsch if (pc->pmat) PetscCall(MatView(pc->pmat, viewer)); 18169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 18174b9ad928SBarry Smith } 18186cd81132SPierre Jolivet if (pop) PetscCall(PetscViewerPopFormat(viewer)); 18194b9ad928SBarry Smith } 18204b9ad928SBarry Smith } else if (isstring) { 18219566063dSJacob Faibussowitsch PetscCall(PCGetType(pc, &cstr)); 18229566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr)); 1823dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 18249566063dSJacob Faibussowitsch if (pc->mat) PetscCall(MatView(pc->mat, viewer)); 18259566063dSJacob Faibussowitsch if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer)); 18265b0b0462SBarry Smith } else if (isbinary) { 182755849f57SBarry Smith PetscInt classid = PC_FILE_CLASSID; 182855849f57SBarry Smith MPI_Comm comm; 182955849f57SBarry Smith PetscMPIInt rank; 183055849f57SBarry Smith char type[256]; 183155849f57SBarry Smith 18329566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 18339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1834dd400576SPatrick Sanan if (rank == 0) { 18359566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT)); 18369566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256)); 18379566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR)); 183855849f57SBarry Smith } 1839dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 1840219b1045SBarry Smith } else if (isdraw) { 1841219b1045SBarry Smith PetscDraw draw; 1842d9884438SBarry Smith char str[25]; 184389fd9fafSBarry Smith PetscReal x, y, bottom, h; 1844d9884438SBarry Smith PetscInt n; 1845219b1045SBarry Smith 18469566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 18479566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 18481d840656SPeter Brune if (pc->mat) { 18499566063dSJacob Faibussowitsch PetscCall(MatGetSize(pc->mat, &n, NULL)); 185063a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n)); 18511d840656SPeter Brune } else { 18529566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name)); 18531d840656SPeter Brune } 18549566063dSJacob Faibussowitsch PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h)); 185589fd9fafSBarry Smith bottom = y - h; 18569566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom)); 1857dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 18589566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 1859e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1860536b137fSBarry Smith } else if (issaws) { 1861d45a07a7SBarry Smith PetscMPIInt rank; 1862d45a07a7SBarry Smith 18639566063dSJacob Faibussowitsch PetscCall(PetscObjectName((PetscObject)pc)); 18649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 186548a46eb9SPierre Jolivet if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer)); 18669566063dSJacob Faibussowitsch if (pc->mat) PetscCall(MatView(pc->mat, viewer)); 18679566063dSJacob Faibussowitsch if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer)); 18680acecf5bSBarry Smith #endif 18694b9ad928SBarry Smith } 18703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18714b9ad928SBarry Smith } 18724b9ad928SBarry Smith 18734b9ad928SBarry Smith /*@C 1874562efe2eSBarry Smith PCRegister - Adds a method (`PCType`) to the preconditioner package. 18751c84c290SBarry Smith 1876cc4c1da9SBarry Smith Not collective. No Fortran Support 18771c84c290SBarry Smith 18781c84c290SBarry Smith Input Parameters: 187920f4b53cSBarry Smith + sname - name of a new user-defined solver 188020f4b53cSBarry Smith - function - routine to create method context 18811c84c290SBarry Smith 1882feefa0e1SJacob Faibussowitsch Example Usage: 18831c84c290SBarry Smith .vb 1884bdf89e91SBarry Smith PCRegister("my_solver", MySolverCreate); 18851c84c290SBarry Smith .ve 18861c84c290SBarry Smith 18871c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 18881c84c290SBarry Smith $ PCSetType(pc, "my_solver") 18891c84c290SBarry Smith or at runtime via the option 18901c84c290SBarry Smith $ -pc_type my_solver 18914b9ad928SBarry Smith 18924b9ad928SBarry Smith Level: advanced 18931c84c290SBarry Smith 189420f4b53cSBarry Smith Note: 189520f4b53cSBarry Smith `PCRegister()` may be called multiple times to add several user-defined preconditioners. 189620f4b53cSBarry Smith 1897562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCType`, `PCRegisterAll()` 18984b9ad928SBarry Smith @*/ 1899d71ae5a4SJacob Faibussowitsch PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC)) 1900d71ae5a4SJacob Faibussowitsch { 19014b9ad928SBarry Smith PetscFunctionBegin; 19029566063dSJacob Faibussowitsch PetscCall(PCInitializePackage()); 19039566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCList, sname, function)); 19043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19054b9ad928SBarry Smith } 19064b9ad928SBarry Smith 1907d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y) 1908d71ae5a4SJacob Faibussowitsch { 1909186a3e20SStefano Zampini PC pc; 1910186a3e20SStefano Zampini 1911186a3e20SStefano Zampini PetscFunctionBegin; 19129566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &pc)); 19139566063dSJacob Faibussowitsch PetscCall(PCApply(pc, X, Y)); 19143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1915186a3e20SStefano Zampini } 1916186a3e20SStefano Zampini 1917feefa0e1SJacob Faibussowitsch /*@C 19180bacdadaSStefano Zampini PCComputeOperator - Computes the explicit preconditioned operator. 19194b9ad928SBarry Smith 1920c3339decSBarry Smith Collective 19214b9ad928SBarry Smith 1922d8d19677SJose E. Roman Input Parameters: 1923186a3e20SStefano Zampini + pc - the preconditioner object 1924562efe2eSBarry Smith - mattype - the `MatType` to be used for the operator 19254b9ad928SBarry Smith 19264b9ad928SBarry Smith Output Parameter: 1927a5b23f4aSJose E. Roman . mat - the explicit preconditioned operator 19284b9ad928SBarry Smith 192920f4b53cSBarry Smith Level: advanced 193020f4b53cSBarry Smith 1931f1580f4eSBarry Smith Note: 1932186a3e20SStefano Zampini This computation is done by applying the operators to columns of the identity matrix. 1933186a3e20SStefano Zampini This routine is costly in general, and is recommended for use only with relatively small systems. 1934562efe2eSBarry Smith Currently, this routine uses a dense matrix format when `mattype` == `NULL` 19354b9ad928SBarry Smith 1936562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPComputeOperator()`, `MatType` 19374b9ad928SBarry Smith @*/ 1938d71ae5a4SJacob Faibussowitsch PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat) 1939d71ae5a4SJacob Faibussowitsch { 1940186a3e20SStefano Zampini PetscInt N, M, m, n; 1941186a3e20SStefano Zampini Mat A, Apc; 19424b9ad928SBarry Smith 19434b9ad928SBarry Smith PetscFunctionBegin; 19440700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 19454f572ea9SToby Isaac PetscAssertPointer(mat, 3); 19469566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, &A, NULL)); 19479566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 19489566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 19499566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc)); 19509566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (void (*)(void))MatMult_PC)); 19519566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(Apc, mattype, mat)); 19529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Apc)); 19533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19544b9ad928SBarry Smith } 19554b9ad928SBarry Smith 19566c237d78SBarry Smith /*@ 19576c237d78SBarry Smith PCSetCoordinates - sets the coordinates of all the nodes on the local process 19586c237d78SBarry Smith 1959c3339decSBarry Smith Collective 19606c237d78SBarry Smith 19616c237d78SBarry Smith Input Parameters: 19626c237d78SBarry Smith + pc - the solver context 19636c237d78SBarry Smith . dim - the dimension of the coordinates 1, 2, or 3 196414893cbeSStefano Zampini . nloc - the blocked size of the coordinates array 196514893cbeSStefano Zampini - coords - the coordinates array 19666c237d78SBarry Smith 19676c237d78SBarry Smith Level: intermediate 19686c237d78SBarry Smith 1969f1580f4eSBarry Smith Note: 197020f4b53cSBarry Smith `coords` is an array of the dim coordinates for the nodes on 197120f4b53cSBarry Smith the local processor, of size `dim`*`nloc`. 197214893cbeSStefano Zampini If there are 108 equation on a processor 19736c237d78SBarry Smith for a displacement finite element discretization of elasticity (so 197414893cbeSStefano Zampini that there are nloc = 36 = 108/3 nodes) then the array must have 108 19756c237d78SBarry Smith double precision values (ie, 3 * 36). These x y z coordinates 19766c237d78SBarry Smith should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x, 19776c237d78SBarry Smith ... , N-1.z ]. 19786c237d78SBarry Smith 1979562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `MatSetNearNullSpace()` 19806c237d78SBarry Smith @*/ 1981d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[]) 1982d71ae5a4SJacob Faibussowitsch { 19836c237d78SBarry Smith PetscFunctionBegin; 198432954da3SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 198532954da3SStefano Zampini PetscValidLogicalCollectiveInt(pc, dim, 2); 198622794d57SStefano Zampini PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords)); 19873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19886c237d78SBarry Smith } 1989fd2dd295SFande Kong 1990fd2dd295SFande Kong /*@ 1991fd2dd295SFande Kong PCGetInterpolations - Gets interpolation matrices for all levels (except level 0) 1992fd2dd295SFande Kong 1993c3339decSBarry Smith Logically Collective 1994fd2dd295SFande Kong 1995d8d19677SJose E. Roman Input Parameter: 1996d8d19677SJose E. Roman . pc - the precondition context 1997fd2dd295SFande Kong 1998d8d19677SJose E. Roman Output Parameters: 1999d8d19677SJose E. Roman + num_levels - the number of levels 2000562efe2eSBarry Smith - interpolations - the interpolation matrices (size of `num_levels`-1) 2001fd2dd295SFande Kong 2002fd2dd295SFande Kong Level: advanced 2003fd2dd295SFande Kong 2004562efe2eSBarry Smith Developer Note: 2005f1580f4eSBarry Smith Why is this here instead of in `PCMG` etc? 2006fd2dd295SFande Kong 2007562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()` 2008fd2dd295SFande Kong @*/ 2009d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[]) 2010d71ae5a4SJacob Faibussowitsch { 2011fd2dd295SFande Kong PetscFunctionBegin; 2012fd2dd295SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 20134f572ea9SToby Isaac PetscAssertPointer(num_levels, 2); 20144f572ea9SToby Isaac PetscAssertPointer(interpolations, 3); 2015cac4c232SBarry Smith PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations)); 20163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2017fd2dd295SFande Kong } 2018fd2dd295SFande Kong 2019fd2dd295SFande Kong /*@ 2020fd2dd295SFande Kong PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level) 2021fd2dd295SFande Kong 2022c3339decSBarry Smith Logically Collective 2023fd2dd295SFande Kong 2024d8d19677SJose E. Roman Input Parameter: 2025d8d19677SJose E. Roman . pc - the precondition context 2026fd2dd295SFande Kong 2027d8d19677SJose E. Roman Output Parameters: 2028d8d19677SJose E. Roman + num_levels - the number of levels 2029562efe2eSBarry Smith - coarseOperators - the coarse operator matrices (size of `num_levels`-1) 2030fd2dd295SFande Kong 2031fd2dd295SFande Kong Level: advanced 2032fd2dd295SFande Kong 2033562efe2eSBarry Smith Developer Note: 2034f1580f4eSBarry Smith Why is this here instead of in `PCMG` etc? 2035fd2dd295SFande Kong 2036562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()` 2037fd2dd295SFande Kong @*/ 2038d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[]) 2039d71ae5a4SJacob Faibussowitsch { 2040fd2dd295SFande Kong PetscFunctionBegin; 2041fd2dd295SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 20424f572ea9SToby Isaac PetscAssertPointer(num_levels, 2); 20434f572ea9SToby Isaac PetscAssertPointer(coarseOperators, 3); 2044cac4c232SBarry Smith PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators)); 20453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2046fd2dd295SFande Kong } 2047