xref: /petsc/src/ksp/pc/interface/precon.c (revision 7530f0189ecda25d023fdcb28178244332f3efb5)
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;
9bded88eaSPierre Jolivet PetscLogEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_MatApply, PC_ApplyCoarse, PC_ApplySymmetricLeft;
10011ea8aeSBarry Smith PetscLogEvent PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks;
11ab83eea4SMatthew G. Knepley PetscInt      PetscMGLevelId;
120316ec64SBarry Smith PetscLogStage PCMPIStage;
134b9ad928SBarry Smith 
PCGetDefaultType_Private(PC pc,const char * type[])14b4f8a55aSStefano Zampini PETSC_INTERN PetscErrorCode PCGetDefaultType_Private(PC pc, const char *type[])
15d71ae5a4SJacob Faibussowitsch {
162e70eadcSBarry Smith   PetscMPIInt size;
17b4f8a55aSStefano Zampini   PetscBool   hasopblock, hasopsolve, flg1, flg2, set, flg3, isnormal;
18566e8bf2SBarry Smith 
19566e8bf2SBarry Smith   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
21566e8bf2SBarry Smith   if (pc->pmat) {
22b4f8a55aSStefano Zampini     PetscCall(MatHasOperation(pc->pmat, MATOP_GET_DIAGONAL_BLOCK, &hasopblock));
23b4f8a55aSStefano Zampini     PetscCall(MatHasOperation(pc->pmat, MATOP_SOLVE, &hasopsolve));
24566e8bf2SBarry Smith     if (size == 1) {
259566063dSJacob Faibussowitsch       PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ICC, &flg1));
269566063dSJacob Faibussowitsch       PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ILU, &flg2));
279566063dSJacob Faibussowitsch       PetscCall(MatIsSymmetricKnown(pc->pmat, &set, &flg3));
28da74c943SJose E. Roman       PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isnormal, MATNORMAL, MATNORMALHERMITIAN, NULL));
29ba8a7149SBarry Smith       if (flg1 && (!flg2 || (set && flg3))) {
30566e8bf2SBarry Smith         *type = PCICC;
31566e8bf2SBarry Smith       } else if (flg2) {
32566e8bf2SBarry Smith         *type = PCILU;
33da74c943SJose E. Roman       } else if (isnormal) {
34da74c943SJose E. Roman         *type = PCNONE;
35b4f8a55aSStefano Zampini       } else if (hasopblock) { /* likely is a parallel matrix run on one processor */
36566e8bf2SBarry Smith         *type = PCBJACOBI;
37b4f8a55aSStefano Zampini       } else if (hasopsolve) {
38b4f8a55aSStefano Zampini         *type = PCMAT;
39566e8bf2SBarry Smith       } else {
40566e8bf2SBarry Smith         *type = PCNONE;
41566e8bf2SBarry Smith       }
42566e8bf2SBarry Smith     } else {
43b4f8a55aSStefano Zampini       if (hasopblock) {
44566e8bf2SBarry Smith         *type = PCBJACOBI;
45b4f8a55aSStefano Zampini       } else if (hasopsolve) {
46b4f8a55aSStefano Zampini         *type = PCMAT;
47566e8bf2SBarry Smith       } else {
48566e8bf2SBarry Smith         *type = PCNONE;
49566e8bf2SBarry Smith       }
50566e8bf2SBarry Smith     }
51b4f8a55aSStefano Zampini   } else *type = NULL;
523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53566e8bf2SBarry Smith }
544b9ad928SBarry Smith 
55fe57bb1aSStefano Zampini /* do not log solves, setup and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */
56fe57bb1aSStefano Zampini PETSC_EXTERN PetscLogEvent KSP_Solve, KSP_SetUp;
57fe57bb1aSStefano Zampini 
PCLogEventsDeactivatePush(void)58fe57bb1aSStefano Zampini static PetscErrorCode PCLogEventsDeactivatePush(void)
59fe57bb1aSStefano Zampini {
60fe57bb1aSStefano Zampini   PetscFunctionBegin;
61fe57bb1aSStefano Zampini   PetscCall(KSPInitializePackage());
62fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePush(KSP_Solve));
63fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePush(KSP_SetUp));
64fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePush(PC_Apply));
65fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePush(PC_SetUp));
66fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePush(PC_SetUpOnBlocks));
67fe57bb1aSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
68fe57bb1aSStefano Zampini }
69fe57bb1aSStefano Zampini 
PCLogEventsDeactivatePop(void)70fe57bb1aSStefano Zampini static PetscErrorCode PCLogEventsDeactivatePop(void)
71fe57bb1aSStefano Zampini {
72fe57bb1aSStefano Zampini   PetscFunctionBegin;
73fe57bb1aSStefano Zampini   PetscCall(KSPInitializePackage());
74fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePop(KSP_Solve));
75fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePop(KSP_SetUp));
76fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePop(PC_Apply));
77fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePop(PC_SetUp));
78fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePop(PC_SetUpOnBlocks));
79fe57bb1aSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
80fe57bb1aSStefano Zampini }
81fe57bb1aSStefano Zampini 
8288584be7SBarry Smith /*@
830b4b7b1cSBarry Smith   PCReset - Resets a `PC` context to the state it was in before `PCSetUp()` was called, and removes any allocated `Vec` and `Mat` from its data structure
8488584be7SBarry Smith 
85c3339decSBarry Smith   Collective
8688584be7SBarry Smith 
8788584be7SBarry Smith   Input Parameter:
880b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
8988584be7SBarry Smith 
9088584be7SBarry Smith   Level: developer
9188584be7SBarry Smith 
920b4b7b1cSBarry Smith   Notes:
930b4b7b1cSBarry Smith   Any options set, including those set with `KSPSetFromOptions()` remain.
940b4b7b1cSBarry Smith 
95562efe2eSBarry 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`
9688584be7SBarry Smith 
97562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`
9888584be7SBarry Smith @*/
PCReset(PC pc)99d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset(PC pc)
100d71ae5a4SJacob Faibussowitsch {
10188584be7SBarry Smith   PetscFunctionBegin;
10288584be7SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
103dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, reset);
1049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleright));
1059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleleft));
1069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->pmat));
1079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->mat));
1082fa5cd67SKarl Rupp 
109371d2eb7SMartin Diehl   pc->setupcalled = PETSC_FALSE;
1103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11188584be7SBarry Smith }
11288584be7SBarry Smith 
1130764c050SBarry Smith /*@
114f1580f4eSBarry Smith   PCDestroy - Destroys `PC` context that was created with `PCCreate()`.
1154b9ad928SBarry Smith 
116c3339decSBarry Smith   Collective
1174b9ad928SBarry Smith 
1184b9ad928SBarry Smith   Input Parameter:
1190b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
1204b9ad928SBarry Smith 
1214b9ad928SBarry Smith   Level: developer
1224b9ad928SBarry Smith 
123562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`
1244b9ad928SBarry Smith @*/
PCDestroy(PC * pc)125d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy(PC *pc)
126d71ae5a4SJacob Faibussowitsch {
1274b9ad928SBarry Smith   PetscFunctionBegin;
1283ba16761SJacob Faibussowitsch   if (!*pc) PetscFunctionReturn(PETSC_SUCCESS);
129f4f49eeaSPierre Jolivet   PetscValidHeaderSpecific(*pc, PC_CLASSID, 1);
130f4f49eeaSPierre Jolivet   if (--((PetscObject)*pc)->refct > 0) {
1319371c9d4SSatish Balay     *pc = NULL;
1323ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1339371c9d4SSatish Balay   }
1344b9ad928SBarry Smith 
1359566063dSJacob Faibussowitsch   PetscCall(PCReset(*pc));
136241cb3abSLisandro Dalcin 
137e04113cfSBarry Smith   /* if memory was published with SAWs then destroy it */
1389566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*pc));
139f4f49eeaSPierre Jolivet   PetscTryTypeMethod(*pc, destroy);
1409566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&(*pc)->dm));
1419566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(pc));
1423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1434b9ad928SBarry Smith }
1444b9ad928SBarry Smith 
145cc4c1da9SBarry Smith /*@
146c5eb9154SBarry Smith   PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right
1474b9ad928SBarry Smith   scaling as needed by certain time-stepping codes.
1484b9ad928SBarry Smith 
149c3339decSBarry Smith   Logically Collective
1504b9ad928SBarry Smith 
1514b9ad928SBarry Smith   Input Parameter:
1520b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
1534b9ad928SBarry Smith 
1544b9ad928SBarry Smith   Output Parameter:
155f1580f4eSBarry Smith . flag - `PETSC_TRUE` if it applies the scaling
1564b9ad928SBarry Smith 
1574b9ad928SBarry Smith   Level: developer
1584b9ad928SBarry Smith 
159f1580f4eSBarry Smith   Note:
160562efe2eSBarry Smith   If this returns `PETSC_TRUE` then the system solved via the Krylov method is, for left and right preconditioning,
1614b9ad928SBarry Smith 
162562efe2eSBarry Smith   $$
163562efe2eSBarry Smith   \begin{align*}
164562efe2eSBarry Smith   D M A D^{-1} y = D M b  \\
165562efe2eSBarry Smith   D A M D^{-1} z = D b.
166562efe2eSBarry Smith   \end{align*}
167562efe2eSBarry Smith   $$
168562efe2eSBarry Smith 
169562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCSetDiagonalScale()`
1704b9ad928SBarry Smith @*/
PCGetDiagonalScale(PC pc,PetscBool * flag)171d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetDiagonalScale(PC pc, PetscBool *flag)
172d71ae5a4SJacob Faibussowitsch {
1734b9ad928SBarry Smith   PetscFunctionBegin;
1740700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1754f572ea9SToby Isaac   PetscAssertPointer(flag, 2);
1764b9ad928SBarry Smith   *flag = pc->diagonalscale;
1773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1784b9ad928SBarry Smith }
1794b9ad928SBarry Smith 
1804b9ad928SBarry Smith /*@
181c5eb9154SBarry Smith   PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right
1824b9ad928SBarry Smith   scaling as needed by certain time-stepping codes.
1834b9ad928SBarry Smith 
184c3339decSBarry Smith   Logically Collective
1854b9ad928SBarry Smith 
1864b9ad928SBarry Smith   Input Parameters:
1870b4b7b1cSBarry Smith + pc - the `PC` preconditioner context
1884b9ad928SBarry Smith - s  - scaling vector
1894b9ad928SBarry Smith 
1904b9ad928SBarry Smith   Level: intermediate
1914b9ad928SBarry Smith 
19295452b02SPatrick Sanan   Notes:
193562efe2eSBarry Smith   The system solved via the Krylov method is, for left and right preconditioning,
194562efe2eSBarry Smith   $$
195562efe2eSBarry Smith   \begin{align*}
196562efe2eSBarry Smith   D M A D^{-1} y = D M b \\
197562efe2eSBarry Smith   D A M D^{-1} z = D b.
198562efe2eSBarry Smith   \end{align*}
199562efe2eSBarry Smith   $$
2004b9ad928SBarry Smith 
201562efe2eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.
2024b9ad928SBarry Smith 
203562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCGetDiagonalScale()`
2044b9ad928SBarry Smith @*/
PCSetDiagonalScale(PC pc,Vec s)205d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetDiagonalScale(PC pc, Vec s)
206d71ae5a4SJacob Faibussowitsch {
2074b9ad928SBarry Smith   PetscFunctionBegin;
2080700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2090700a824SBarry Smith   PetscValidHeaderSpecific(s, VEC_CLASSID, 2);
2104b9ad928SBarry Smith   pc->diagonalscale = PETSC_TRUE;
2112fa5cd67SKarl Rupp 
2129566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)s));
2139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleleft));
2142fa5cd67SKarl Rupp 
2154b9ad928SBarry Smith   pc->diagonalscaleleft = s;
2162fa5cd67SKarl Rupp 
2179566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(s, &pc->diagonalscaleright));
2189566063dSJacob Faibussowitsch   PetscCall(VecCopy(s, pc->diagonalscaleright));
2199566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(pc->diagonalscaleright));
2203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2214b9ad928SBarry Smith }
2224b9ad928SBarry Smith 
223bc08b0f1SBarry Smith /*@
224c5eb9154SBarry Smith   PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes.
2254b9ad928SBarry Smith 
226c3339decSBarry Smith   Logically Collective
2274b9ad928SBarry Smith 
2284b9ad928SBarry Smith   Input Parameters:
2290b4b7b1cSBarry Smith + pc  - the `PC` preconditioner context
2304b9ad928SBarry Smith . in  - input vector
231a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in)
2324b9ad928SBarry Smith 
2334b9ad928SBarry Smith   Level: intermediate
2344b9ad928SBarry Smith 
23595452b02SPatrick Sanan   Notes:
236562efe2eSBarry Smith   The system solved via the Krylov method is, for left and right preconditioning,
2374b9ad928SBarry Smith 
238562efe2eSBarry Smith   $$
239562efe2eSBarry Smith   \begin{align*}
240562efe2eSBarry Smith   D M A D^{-1} y = D M b  \\
241562efe2eSBarry Smith   D A M D^{-1} z = D b.
242562efe2eSBarry Smith   \end{align*}
243562efe2eSBarry Smith   $$
2444b9ad928SBarry Smith 
245562efe2eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.
2464b9ad928SBarry Smith 
247562efe2eSBarry Smith   If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out`
248562efe2eSBarry Smith 
2490241b274SPierre Jolivet .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCSetDiagonalScale()`, `PCDiagonalScaleRight()`, `MatDiagonalScale()`
2504b9ad928SBarry Smith @*/
PCDiagonalScaleLeft(PC pc,Vec in,Vec out)251d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleLeft(PC pc, Vec in, Vec out)
252d71ae5a4SJacob Faibussowitsch {
2534b9ad928SBarry Smith   PetscFunctionBegin;
2540700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2550700a824SBarry Smith   PetscValidHeaderSpecific(in, VEC_CLASSID, 2);
2560700a824SBarry Smith   PetscValidHeaderSpecific(out, VEC_CLASSID, 3);
2574b9ad928SBarry Smith   if (pc->diagonalscale) {
2589566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(out, pc->diagonalscaleleft, in));
2594b9ad928SBarry Smith   } else if (in != out) {
2609566063dSJacob Faibussowitsch     PetscCall(VecCopy(in, out));
2614b9ad928SBarry Smith   }
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2634b9ad928SBarry Smith }
2644b9ad928SBarry Smith 
265bc08b0f1SBarry Smith /*@
2664b9ad928SBarry Smith   PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
2674b9ad928SBarry Smith 
268c3339decSBarry Smith   Logically Collective
2694b9ad928SBarry Smith 
2704b9ad928SBarry Smith   Input Parameters:
2710b4b7b1cSBarry Smith + pc  - the `PC` preconditioner context
2724b9ad928SBarry Smith . in  - input vector
273a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in)
2744b9ad928SBarry Smith 
2754b9ad928SBarry Smith   Level: intermediate
2764b9ad928SBarry Smith 
27795452b02SPatrick Sanan   Notes:
278562efe2eSBarry Smith   The system solved via the Krylov method is, for left and right preconditioning,
2794b9ad928SBarry Smith 
280562efe2eSBarry Smith   $$
281562efe2eSBarry Smith   \begin{align*}
282562efe2eSBarry Smith   D M A D^{-1} y = D M b  \\
283562efe2eSBarry Smith   D A M D^{-1} z = D b.
284a3524536SPierre Jolivet   \end{align*}
285562efe2eSBarry Smith   $$
2864b9ad928SBarry Smith 
287562efe2eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.
2884b9ad928SBarry Smith 
289562efe2eSBarry Smith   If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out`
290562efe2eSBarry Smith 
2910241b274SPierre Jolivet .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCSetDiagonalScale()`, `MatDiagonalScale()`
2924b9ad928SBarry Smith @*/
PCDiagonalScaleRight(PC pc,Vec in,Vec out)293d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleRight(PC pc, Vec in, Vec out)
294d71ae5a4SJacob Faibussowitsch {
2954b9ad928SBarry Smith   PetscFunctionBegin;
2960700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2970700a824SBarry Smith   PetscValidHeaderSpecific(in, VEC_CLASSID, 2);
2980700a824SBarry Smith   PetscValidHeaderSpecific(out, VEC_CLASSID, 3);
2994b9ad928SBarry Smith   if (pc->diagonalscale) {
3009566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(out, pc->diagonalscaleright, in));
3014b9ad928SBarry Smith   } else if (in != out) {
3029566063dSJacob Faibussowitsch     PetscCall(VecCopy(in, out));
3034b9ad928SBarry Smith   }
3043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3054b9ad928SBarry Smith }
3064b9ad928SBarry Smith 
30749517cdeSBarry Smith /*@
30849517cdeSBarry Smith   PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the
309f1580f4eSBarry Smith   operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
310f1580f4eSBarry Smith   `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.
31149517cdeSBarry Smith 
312c3339decSBarry Smith   Logically Collective
31349517cdeSBarry Smith 
31449517cdeSBarry Smith   Input Parameters:
3150b4b7b1cSBarry Smith + pc  - the `PC` preconditioner context
316f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)
31749517cdeSBarry Smith 
31849517cdeSBarry Smith   Options Database Key:
319562efe2eSBarry Smith . -pc_use_amat <true,false> - use the amat argument to `KSPSetOperators()` or `PCSetOperators()` to apply the operator
32049517cdeSBarry Smith 
32120f4b53cSBarry Smith   Level: intermediate
32220f4b53cSBarry Smith 
323f1580f4eSBarry Smith   Note:
32449517cdeSBarry Smith   For the common case in which the linear system matrix and the matrix used to construct the
325562efe2eSBarry Smith   preconditioner are identical, this routine has no affect.
32649517cdeSBarry Smith 
3270241b274SPierre Jolivet .seealso: [](ch_ksp), `PC`, `PCGetUseAmat()`, `PCBJACOBI`, `PCMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`,
328562efe2eSBarry Smith           `KSPSetOperators()`, `PCSetOperators()`
32949517cdeSBarry Smith @*/
PCSetUseAmat(PC pc,PetscBool flg)330d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUseAmat(PC pc, PetscBool flg)
331d71ae5a4SJacob Faibussowitsch {
33249517cdeSBarry Smith   PetscFunctionBegin;
33349517cdeSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
33449517cdeSBarry Smith   pc->useAmat = flg;
3353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33649517cdeSBarry Smith }
33749517cdeSBarry Smith 
338422a814eSBarry Smith /*@
339562efe2eSBarry Smith   PCSetErrorIfFailure - Causes `PC` to generate an error if a floating point exception, for example a zero pivot, is detected.
340422a814eSBarry Smith 
341c3339decSBarry Smith   Logically Collective
342422a814eSBarry Smith 
343422a814eSBarry Smith   Input Parameters:
344562efe2eSBarry Smith + pc  - iterative context obtained from `PCCreate()`
345f1580f4eSBarry Smith - flg - `PETSC_TRUE` indicates you want the error generated
346422a814eSBarry Smith 
347422a814eSBarry Smith   Level: advanced
348422a814eSBarry Smith 
349422a814eSBarry Smith   Notes:
350f1580f4eSBarry Smith   Normally PETSc continues if a linear solver fails due to a failed setup of a preconditioner, you can call `KSPGetConvergedReason()` after a `KSPSolve()`
351422a814eSBarry 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
352422a814eSBarry Smith   detected.
353422a814eSBarry Smith 
354562efe2eSBarry Smith   This is propagated into `KSP`s used by this `PC`, which then propagate it into `PC`s used by those `KSP`s
355422a814eSBarry Smith 
356562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPSetErrorIfNotConverged()`, `PCGetInitialGuessNonzero()`, `PCSetInitialGuessKnoll()`, `PCGetInitialGuessKnoll()`
357422a814eSBarry Smith @*/
PCSetErrorIfFailure(PC pc,PetscBool flg)358d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetErrorIfFailure(PC pc, PetscBool flg)
359d71ae5a4SJacob Faibussowitsch {
360422a814eSBarry Smith   PetscFunctionBegin;
361422a814eSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
362422a814eSBarry Smith   PetscValidLogicalCollectiveBool(pc, flg, 2);
363422a814eSBarry Smith   pc->erroriffailure = flg;
3643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
365422a814eSBarry Smith }
366422a814eSBarry Smith 
36749517cdeSBarry Smith /*@
36849517cdeSBarry Smith   PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the
369f1580f4eSBarry Smith   operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
370f1580f4eSBarry Smith   `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.
37149517cdeSBarry Smith 
372c3339decSBarry Smith   Logically Collective
37349517cdeSBarry Smith 
37449517cdeSBarry Smith   Input Parameter:
3750b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
37649517cdeSBarry Smith 
37749517cdeSBarry Smith   Output Parameter:
378f1580f4eSBarry Smith . flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)
37949517cdeSBarry Smith 
38020f4b53cSBarry Smith   Level: intermediate
38120f4b53cSBarry Smith 
382f1580f4eSBarry Smith   Note:
38349517cdeSBarry Smith   For the common case in which the linear system matrix and the matrix used to construct the
38449517cdeSBarry Smith   preconditioner are identical, this routine is does nothing.
38549517cdeSBarry Smith 
3860241b274SPierre Jolivet .seealso: [](ch_ksp), `PC`, `PCSetUseAmat()`, `PCBJACOBI`, `PCMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`
38749517cdeSBarry Smith @*/
PCGetUseAmat(PC pc,PetscBool * flg)388d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetUseAmat(PC pc, PetscBool *flg)
389d71ae5a4SJacob Faibussowitsch {
39049517cdeSBarry Smith   PetscFunctionBegin;
39149517cdeSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
39249517cdeSBarry Smith   *flg = pc->useAmat;
3933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39449517cdeSBarry Smith }
39549517cdeSBarry Smith 
396f39d8e23SSatish Balay /*@
3973821be0aSBarry Smith   PCSetKSPNestLevel - sets the amount of nesting the `KSP` that contains this `PC` has
3983821be0aSBarry Smith 
3993821be0aSBarry Smith   Collective
4003821be0aSBarry Smith 
4013821be0aSBarry Smith   Input Parameters:
4023821be0aSBarry Smith + pc    - the `PC`
4033821be0aSBarry Smith - level - the nest level
4043821be0aSBarry Smith 
4053821be0aSBarry Smith   Level: developer
4063821be0aSBarry Smith 
4073821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCGetKSPNestLevel()`, `KSPSetNestLevel()`
4083821be0aSBarry Smith @*/
PCSetKSPNestLevel(PC pc,PetscInt level)4093821be0aSBarry Smith PetscErrorCode PCSetKSPNestLevel(PC pc, PetscInt level)
4103821be0aSBarry Smith {
4113821be0aSBarry Smith   PetscFunctionBegin;
4127a99bfcaSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4137a99bfcaSBarry Smith   PetscValidLogicalCollectiveInt(pc, level, 2);
4143821be0aSBarry Smith   pc->kspnestlevel = level;
4153821be0aSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
4163821be0aSBarry Smith }
4173821be0aSBarry Smith 
4183821be0aSBarry Smith /*@
4193821be0aSBarry Smith   PCGetKSPNestLevel - gets the amount of nesting the `KSP` that contains this `PC` has
4203821be0aSBarry Smith 
4217a99bfcaSBarry Smith   Not Collective
4223821be0aSBarry Smith 
4233821be0aSBarry Smith   Input Parameter:
4243821be0aSBarry Smith . pc - the `PC`
4253821be0aSBarry Smith 
4263821be0aSBarry Smith   Output Parameter:
4273821be0aSBarry Smith . level - the nest level
4283821be0aSBarry Smith 
4293821be0aSBarry Smith   Level: developer
4303821be0aSBarry Smith 
4313821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `KSPGetNestLevel()`
4323821be0aSBarry Smith @*/
PCGetKSPNestLevel(PC pc,PetscInt * level)4333821be0aSBarry Smith PetscErrorCode PCGetKSPNestLevel(PC pc, PetscInt *level)
4343821be0aSBarry Smith {
4353821be0aSBarry Smith   PetscFunctionBegin;
4367a99bfcaSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4377a99bfcaSBarry Smith   PetscAssertPointer(level, 2);
4383821be0aSBarry Smith   *level = pc->kspnestlevel;
4393821be0aSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
4403821be0aSBarry Smith }
4413821be0aSBarry Smith 
4423821be0aSBarry Smith /*@
443f1580f4eSBarry Smith   PCCreate - Creates a preconditioner context, `PC`
4444b9ad928SBarry Smith 
445d083f849SBarry Smith   Collective
4464b9ad928SBarry Smith 
4474b9ad928SBarry Smith   Input Parameter:
4484b9ad928SBarry Smith . comm - MPI communicator
4494b9ad928SBarry Smith 
4504b9ad928SBarry Smith   Output Parameter:
4510b4b7b1cSBarry Smith . newpc - location to put the `PC` preconditioner context
4524b9ad928SBarry Smith 
45320f4b53cSBarry Smith   Level: developer
45420f4b53cSBarry Smith 
4550b4b7b1cSBarry Smith   Notes:
4560b4b7b1cSBarry Smith   This is rarely called directly by users since `KSP` manages the `PC` objects it uses. Use `KSPGetPC()` to access the `PC` used by a `KSP`.
4570b4b7b1cSBarry Smith 
4580b4b7b1cSBarry Smith   Use `PCSetType()` or `PCSetFromOptions()` with the option `-pc_type pctype` to set the `PCType` for this `PC`
4590b4b7b1cSBarry Smith 
4600b4b7b1cSBarry Smith   The default preconditioner type `PCType` for sparse matrices is `PCILU` or `PCICC` with 0 fill on one process and block Jacobi (`PCBJACOBI`) with `PCILU` or `PCICC`
461f1580f4eSBarry Smith   in parallel. For dense matrices it is always `PCNONE`.
4624b9ad928SBarry Smith 
4630b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCType`, `PCSetType`, `PCSetUp()`, `PCApply()`, `PCDestroy()`, `KSP`, `KSPGetPC()`
4644b9ad928SBarry Smith @*/
PCCreate(MPI_Comm comm,PC * newpc)465d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreate(MPI_Comm comm, PC *newpc)
466d71ae5a4SJacob Faibussowitsch {
4674b9ad928SBarry Smith   PC pc;
4684b9ad928SBarry Smith 
4694b9ad928SBarry Smith   PetscFunctionBegin;
4704f572ea9SToby Isaac   PetscAssertPointer(newpc, 2);
4719566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
4724b9ad928SBarry Smith 
4739566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(pc, PC_CLASSID, "PC", "Preconditioner", "PC", comm, PCDestroy, PCView));
4740a545947SLisandro Dalcin   pc->mat                  = NULL;
4750a545947SLisandro Dalcin   pc->pmat                 = NULL;
476371d2eb7SMartin Diehl   pc->setupcalled          = PETSC_FALSE;
4770c24e6a1SHong Zhang   pc->setfromoptionscalled = 0;
4780a545947SLisandro Dalcin   pc->data                 = NULL;
4794b9ad928SBarry Smith   pc->diagonalscale        = PETSC_FALSE;
4800a545947SLisandro Dalcin   pc->diagonalscaleleft    = NULL;
4810a545947SLisandro Dalcin   pc->diagonalscaleright   = NULL;
4824b9ad928SBarry Smith 
4830a545947SLisandro Dalcin   pc->modifysubmatrices  = NULL;
4840a545947SLisandro Dalcin   pc->modifysubmatricesP = NULL;
4852fa5cd67SKarl Rupp 
486a7d21a52SLisandro Dalcin   *newpc = pc;
4873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4884b9ad928SBarry Smith }
4894b9ad928SBarry Smith 
4904b9ad928SBarry Smith /*@
4914b9ad928SBarry Smith   PCApply - Applies the preconditioner to a vector.
4924b9ad928SBarry Smith 
493c3339decSBarry Smith   Collective
4944b9ad928SBarry Smith 
4954b9ad928SBarry Smith   Input Parameters:
4960b4b7b1cSBarry Smith + pc - the `PC` preconditioner context
4974b9ad928SBarry Smith - x  - input vector
4984b9ad928SBarry Smith 
4994b9ad928SBarry Smith   Output Parameter:
5004b9ad928SBarry Smith . y - output vector
5014b9ad928SBarry Smith 
5024b9ad928SBarry Smith   Level: developer
5034b9ad928SBarry Smith 
504562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`, `PCApplyBAorAB()`
5054b9ad928SBarry Smith @*/
PCApply(PC pc,Vec x,Vec y)506d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApply(PC pc, Vec x, Vec y)
507d71ae5a4SJacob Faibussowitsch {
508106e20bfSBarry Smith   PetscInt m, n, mv, nv;
5094b9ad928SBarry Smith 
5104b9ad928SBarry Smith   PetscFunctionBegin;
5110700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
5120700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
5130700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
5147827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
515e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
516540bdfbdSVaclav Hapla   /* use pmat to check vector sizes since for KSPLSQR the pmat may be of a different size than mat */
5179566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &m, &n));
5189566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x, &mv));
5199566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(y, &nv));
520540bdfbdSVaclav Hapla   /* check pmat * y = x is feasible */
52163a3b9bcSJacob 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);
52263a3b9bcSJacob 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);
5239566063dSJacob Faibussowitsch   PetscCall(VecSetErrorIfLocked(y, 3));
524106e20bfSBarry Smith 
5259566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
5269566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
5279566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
528dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, apply, x, y);
5299566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
530e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
5319566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
5323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5334b9ad928SBarry Smith }
5344b9ad928SBarry Smith 
PCMatApplyTranspose_Private(PC pc,Mat X,Mat Y,PetscBool transpose)5354dbf25a8SPierre Jolivet static PetscErrorCode PCMatApplyTranspose_Private(PC pc, Mat X, Mat Y, PetscBool transpose)
536d71ae5a4SJacob Faibussowitsch {
537c677e75fSPierre Jolivet   Mat       A;
538c677e75fSPierre Jolivet   Vec       cy, cx;
539bd82155bSPierre Jolivet   PetscInt  m1, M1, m2, M2, n1, N1, n2, N2, m3, M3, n3, N3;
540c677e75fSPierre Jolivet   PetscBool match;
541c677e75fSPierre Jolivet 
542c677e75fSPierre Jolivet   PetscFunctionBegin;
543c677e75fSPierre Jolivet   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
544e2b200f6SPierre Jolivet   PetscValidHeaderSpecific(X, MAT_CLASSID, 2);
545e2b200f6SPierre Jolivet   PetscValidHeaderSpecific(Y, MAT_CLASSID, 3);
546e2b200f6SPierre Jolivet   PetscCheckSameComm(pc, 1, X, 2);
547e2b200f6SPierre Jolivet   PetscCheckSameComm(pc, 1, Y, 3);
5487827d75bSBarry Smith   PetscCheck(Y != X, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "Y and X must be different matrices");
5499566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &A));
5509566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m3, &n3));
5519566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m2, &n2));
5529566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &m1, &n1));
5539566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M3, &N3));
5549566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &M2, &N2));
5559566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Y, &M1, &N1));
55663a3b9bcSJacob 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);
55763a3b9bcSJacob 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);
55863a3b9bcSJacob 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);
5599566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &match, MATSEQDENSE, MATMPIDENSE, ""));
56028b400f6SJacob Faibussowitsch   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of output vectors not stored in a dense Mat");
5619566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
56228b400f6SJacob Faibussowitsch   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of input vectors not stored in a dense Mat");
5639566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
5644dbf25a8SPierre Jolivet   if (!transpose && pc->ops->matapply) {
5659566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0));
566dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, matapply, X, Y);
5679566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0));
5684dbf25a8SPierre Jolivet   } else if (transpose && pc->ops->matapplytranspose) {
5694dbf25a8SPierre Jolivet     PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0));
5704dbf25a8SPierre Jolivet     PetscUseTypeMethod(pc, matapplytranspose, X, Y);
5714dbf25a8SPierre Jolivet     PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0));
572c677e75fSPierre Jolivet   } else {
5739566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "PC type %s applying column by column\n", ((PetscObject)pc)->type_name));
574bd82155bSPierre Jolivet     for (n1 = 0; n1 < N1; ++n1) {
5759566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecRead(X, n1, &cx));
5769566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecWrite(Y, n1, &cy));
5774dbf25a8SPierre Jolivet       if (!transpose) PetscCall(PCApply(pc, cx, cy));
5784dbf25a8SPierre Jolivet       else PetscCall(PCApplyTranspose(pc, cx, cy));
5799566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecWrite(Y, n1, &cy));
5809566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecRead(X, n1, &cx));
581c677e75fSPierre Jolivet     }
582c677e75fSPierre Jolivet   }
5833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
584c677e75fSPierre Jolivet }
585c677e75fSPierre Jolivet 
586c677e75fSPierre Jolivet /*@
5874dbf25a8SPierre Jolivet   PCMatApply - Applies the preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApply()`, `Y` and `X` must be different matrices.
5884dbf25a8SPierre Jolivet 
5894dbf25a8SPierre Jolivet   Collective
5904dbf25a8SPierre Jolivet 
5914dbf25a8SPierre Jolivet   Input Parameters:
5924dbf25a8SPierre Jolivet + pc - the `PC` preconditioner context
5934dbf25a8SPierre Jolivet - X  - block of input vectors
5944dbf25a8SPierre Jolivet 
5954dbf25a8SPierre Jolivet   Output Parameter:
5964dbf25a8SPierre Jolivet . Y - block of output vectors
5974dbf25a8SPierre Jolivet 
5984dbf25a8SPierre Jolivet   Level: developer
5994dbf25a8SPierre Jolivet 
6004dbf25a8SPierre Jolivet .seealso: [](ch_ksp), `PC`, `PCApply()`, `KSPMatSolve()`
6014dbf25a8SPierre Jolivet @*/
PCMatApply(PC pc,Mat X,Mat Y)6024dbf25a8SPierre Jolivet PetscErrorCode PCMatApply(PC pc, Mat X, Mat Y)
6034dbf25a8SPierre Jolivet {
6044dbf25a8SPierre Jolivet   PetscFunctionBegin;
6054dbf25a8SPierre Jolivet   PetscCall(PCMatApplyTranspose_Private(pc, X, Y, PETSC_FALSE));
6064dbf25a8SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
6074dbf25a8SPierre Jolivet }
6084dbf25a8SPierre Jolivet 
6094dbf25a8SPierre Jolivet /*@
6104dbf25a8SPierre Jolivet   PCMatApplyTranspose - Applies the transpose of preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApplyTranspose()`, `Y` and `X` must be different matrices.
6114dbf25a8SPierre Jolivet 
6124dbf25a8SPierre Jolivet   Collective
6134dbf25a8SPierre Jolivet 
6144dbf25a8SPierre Jolivet   Input Parameters:
6154dbf25a8SPierre Jolivet + pc - the `PC` preconditioner context
6164dbf25a8SPierre Jolivet - X  - block of input vectors
6174dbf25a8SPierre Jolivet 
6184dbf25a8SPierre Jolivet   Output Parameter:
6194dbf25a8SPierre Jolivet . Y - block of output vectors
6204dbf25a8SPierre Jolivet 
6214dbf25a8SPierre Jolivet   Level: developer
6224dbf25a8SPierre Jolivet 
6234dbf25a8SPierre Jolivet .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`, `KSPMatSolveTranspose()`
6244dbf25a8SPierre Jolivet @*/
PCMatApplyTranspose(PC pc,Mat X,Mat Y)6254dbf25a8SPierre Jolivet PetscErrorCode PCMatApplyTranspose(PC pc, Mat X, Mat Y)
6264dbf25a8SPierre Jolivet {
6274dbf25a8SPierre Jolivet   PetscFunctionBegin;
6284dbf25a8SPierre Jolivet   PetscCall(PCMatApplyTranspose_Private(pc, X, Y, PETSC_TRUE));
6294dbf25a8SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
6304dbf25a8SPierre Jolivet }
6314dbf25a8SPierre Jolivet 
6324dbf25a8SPierre Jolivet /*@
6334b9ad928SBarry Smith   PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
6344b9ad928SBarry Smith 
635c3339decSBarry Smith   Collective
6364b9ad928SBarry Smith 
6374b9ad928SBarry Smith   Input Parameters:
6380b4b7b1cSBarry Smith + pc - the `PC` preconditioner context
6394b9ad928SBarry Smith - x  - input vector
6404b9ad928SBarry Smith 
6414b9ad928SBarry Smith   Output Parameter:
6424b9ad928SBarry Smith . y - output vector
6434b9ad928SBarry Smith 
64420f4b53cSBarry Smith   Level: developer
64520f4b53cSBarry Smith 
646f1580f4eSBarry Smith   Note:
647f1580f4eSBarry Smith   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
6484b9ad928SBarry Smith 
649562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricRight()`
6504b9ad928SBarry Smith @*/
PCApplySymmetricLeft(PC pc,Vec x,Vec y)651d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricLeft(PC pc, Vec x, Vec y)
652d71ae5a4SJacob Faibussowitsch {
6534b9ad928SBarry Smith   PetscFunctionBegin;
6540700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6550700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6560700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6577827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
658e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
6599566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
6609566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
6619566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplySymmetricLeft, pc, x, y, 0));
662dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applysymmetricleft, x, y);
6639566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplySymmetricLeft, pc, x, y, 0));
6649566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
665e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
6663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6674b9ad928SBarry Smith }
6684b9ad928SBarry Smith 
6694b9ad928SBarry Smith /*@
6704b9ad928SBarry Smith   PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
6714b9ad928SBarry Smith 
672c3339decSBarry Smith   Collective
6734b9ad928SBarry Smith 
6744b9ad928SBarry Smith   Input Parameters:
6750b4b7b1cSBarry Smith + pc - the `PC` preconditioner context
6764b9ad928SBarry Smith - x  - input vector
6774b9ad928SBarry Smith 
6784b9ad928SBarry Smith   Output Parameter:
6794b9ad928SBarry Smith . y - output vector
6804b9ad928SBarry Smith 
6814b9ad928SBarry Smith   Level: developer
6824b9ad928SBarry Smith 
683f1580f4eSBarry Smith   Note:
684f1580f4eSBarry Smith   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
6854b9ad928SBarry Smith 
686562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricLeft()`
6874b9ad928SBarry Smith @*/
PCApplySymmetricRight(PC pc,Vec x,Vec y)688d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricRight(PC pc, Vec x, Vec y)
689d71ae5a4SJacob Faibussowitsch {
6904b9ad928SBarry Smith   PetscFunctionBegin;
6910700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6920700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6930700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6947827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
695e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
6969566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
6979566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
6989566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplySymmetricRight, pc, x, y, 0));
699dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applysymmetricright, x, y);
7009566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplySymmetricRight, pc, x, y, 0));
7019566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
702e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
7033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7044b9ad928SBarry Smith }
7054b9ad928SBarry Smith 
7064b9ad928SBarry Smith /*@
7074b9ad928SBarry Smith   PCApplyTranspose - Applies the transpose of preconditioner to a vector.
7084b9ad928SBarry Smith 
709c3339decSBarry Smith   Collective
7104b9ad928SBarry Smith 
7114b9ad928SBarry Smith   Input Parameters:
7120b4b7b1cSBarry Smith + pc - the `PC` preconditioner context
7134b9ad928SBarry Smith - x  - input vector
7144b9ad928SBarry Smith 
7154b9ad928SBarry Smith   Output Parameter:
7164b9ad928SBarry Smith . y - output vector
7174b9ad928SBarry Smith 
71820f4b53cSBarry Smith   Level: developer
71920f4b53cSBarry Smith 
720f1580f4eSBarry Smith   Note:
72195452b02SPatrick Sanan   For complex numbers this applies the non-Hermitian transpose.
7224c97465dSBarry Smith 
723562efe2eSBarry Smith   Developer Note:
724f1580f4eSBarry Smith   We need to implement a `PCApplyHermitianTranspose()`
7254c97465dSBarry Smith 
726562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyBAorAB()`, `PCApplyBAorABTranspose()`, `PCApplyTransposeExists()`
7274b9ad928SBarry Smith @*/
PCApplyTranspose(PC pc,Vec x,Vec y)728d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTranspose(PC pc, Vec x, Vec y)
729d71ae5a4SJacob Faibussowitsch {
7304b9ad928SBarry Smith   PetscFunctionBegin;
7310700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
7320700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
7330700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
7347827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
735e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
7369566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
7379566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
7389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
739dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applytranspose, x, y);
7409566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
7419566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
742e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
7433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7444b9ad928SBarry Smith }
7454b9ad928SBarry Smith 
7464b9ad928SBarry Smith /*@
7479cc050e5SLisandro Dalcin   PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
7484b9ad928SBarry Smith 
749c3339decSBarry Smith   Collective
7504b9ad928SBarry Smith 
751f1580f4eSBarry Smith   Input Parameter:
7520b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
7534b9ad928SBarry Smith 
7544b9ad928SBarry Smith   Output Parameter:
755f1580f4eSBarry Smith . flg - `PETSC_TRUE` if a transpose operation is defined
7564b9ad928SBarry Smith 
7574b9ad928SBarry Smith   Level: developer
7584b9ad928SBarry Smith 
759562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`
7604b9ad928SBarry Smith @*/
PCApplyTransposeExists(PC pc,PetscBool * flg)761d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTransposeExists(PC pc, PetscBool *flg)
762d71ae5a4SJacob Faibussowitsch {
7634b9ad928SBarry Smith   PetscFunctionBegin;
7640700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
7654f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
7666ce5e81cSLisandro Dalcin   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
7676ce5e81cSLisandro Dalcin   else *flg = PETSC_FALSE;
7683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7694b9ad928SBarry Smith }
7704b9ad928SBarry Smith 
7714b9ad928SBarry Smith /*@
772562efe2eSBarry Smith   PCApplyBAorAB - Applies the preconditioner and operator to a vector. $y = B*A*x $ or $ y = A*B*x$.
7734b9ad928SBarry Smith 
774c3339decSBarry Smith   Collective
7754b9ad928SBarry Smith 
7764b9ad928SBarry Smith   Input Parameters:
7770b4b7b1cSBarry Smith + pc   - the `PC` preconditioner context
778f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
7794b9ad928SBarry Smith . x    - input vector
7804b9ad928SBarry Smith - work - work vector
7814b9ad928SBarry Smith 
7824b9ad928SBarry Smith   Output Parameter:
7834b9ad928SBarry Smith . y - output vector
7844b9ad928SBarry Smith 
7854b9ad928SBarry Smith   Level: developer
7864b9ad928SBarry Smith 
787f1580f4eSBarry Smith   Note:
788562efe2eSBarry 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.
789562efe2eSBarry Smith   The specific `KSPSolve()` method must also be written to handle the post-solve "correction" for the diagonal scaling.
79013e0d083SBarry Smith 
791562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorABTranspose()`
7924b9ad928SBarry Smith @*/
PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)793d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorAB(PC pc, PCSide side, Vec x, Vec y, Vec work)
794d71ae5a4SJacob Faibussowitsch {
7954b9ad928SBarry Smith   PetscFunctionBegin;
7960700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
797a69c7061SStefano Zampini   PetscValidLogicalCollectiveEnum(pc, side, 2);
7980700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
7990700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
8000700a824SBarry Smith   PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
801a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, x, 3);
802a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, y, 4);
803a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, work, 5);
8047827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
8057827d75bSBarry 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");
8067827d75bSBarry Smith   PetscCheck(!pc->diagonalscale || side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot include diagonal scaling with symmetric preconditioner application");
807e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
8084b9ad928SBarry Smith 
8099566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
8104b9ad928SBarry Smith   if (pc->diagonalscale) {
8114b9ad928SBarry Smith     if (pc->ops->applyBA) {
8124b9ad928SBarry Smith       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
8139566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &work2));
8149566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, work2));
815dbbe0bcdSBarry Smith       PetscUseTypeMethod(pc, applyBA, side, work2, y, work);
8169566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
8179566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&work2));
8184b9ad928SBarry Smith     } else if (side == PC_RIGHT) {
8199566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, y));
8209566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, y, work));
8219566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
8229566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
8234b9ad928SBarry Smith     } else if (side == PC_LEFT) {
8249566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, y));
8259566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, y, work));
8269566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, work, y));
8279566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
8287827d75bSBarry Smith     } else PetscCheck(side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot provide diagonal scaling with symmetric application of preconditioner");
8294b9ad928SBarry Smith   } else {
8304b9ad928SBarry Smith     if (pc->ops->applyBA) {
831dbbe0bcdSBarry Smith       PetscUseTypeMethod(pc, applyBA, side, x, y, work);
8324b9ad928SBarry Smith     } else if (side == PC_RIGHT) {
8339566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, x, work));
8349566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
8354b9ad928SBarry Smith     } else if (side == PC_LEFT) {
8369566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, x, work));
8379566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, work, y));
8384b9ad928SBarry Smith     } else if (side == PC_SYMMETRIC) {
8394b9ad928SBarry Smith       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
8409566063dSJacob Faibussowitsch       PetscCall(PCApplySymmetricRight(pc, x, work));
8419566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
8429566063dSJacob Faibussowitsch       PetscCall(VecCopy(y, work));
8439566063dSJacob Faibussowitsch       PetscCall(PCApplySymmetricLeft(pc, work, y));
8444b9ad928SBarry Smith     }
8454b9ad928SBarry Smith   }
846e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8484b9ad928SBarry Smith }
8494b9ad928SBarry Smith 
8504b9ad928SBarry Smith /*@
8514b9ad928SBarry Smith   PCApplyBAorABTranspose - Applies the transpose of the preconditioner
852562efe2eSBarry Smith   and operator to a vector. That is, applies $B^T * A^T$ with left preconditioning,
853562efe2eSBarry Smith   NOT $(B*A)^T = A^T*B^T$.
8544b9ad928SBarry Smith 
855c3339decSBarry Smith   Collective
8564b9ad928SBarry Smith 
8574b9ad928SBarry Smith   Input Parameters:
8580b4b7b1cSBarry Smith + pc   - the `PC` preconditioner context
859f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
8604b9ad928SBarry Smith . x    - input vector
8614b9ad928SBarry Smith - work - work vector
8624b9ad928SBarry Smith 
8634b9ad928SBarry Smith   Output Parameter:
8644b9ad928SBarry Smith . y - output vector
8654b9ad928SBarry Smith 
86620f4b53cSBarry Smith   Level: developer
86720f4b53cSBarry Smith 
868f1580f4eSBarry Smith   Note:
869562efe2eSBarry 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
870562efe2eSBarry Smith   defined by $B^T$. This is why this has the funny form that it computes $B^T * A^T$
8719b3038f0SBarry Smith 
872562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorAB()`
8734b9ad928SBarry Smith @*/
PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)874d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorABTranspose(PC pc, PCSide side, Vec x, Vec y, Vec work)
875d71ae5a4SJacob Faibussowitsch {
8764b9ad928SBarry Smith   PetscFunctionBegin;
8770700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8780700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
8790700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
8800700a824SBarry Smith   PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
8817827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
882e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
8834b9ad928SBarry Smith   if (pc->ops->applyBAtranspose) {
884dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, applyBAtranspose, side, x, y, work);
885e0f629ddSJacob Faibussowitsch     if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8863ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
8874b9ad928SBarry Smith   }
8887827d75bSBarry Smith   PetscCheck(side == PC_LEFT || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right or left");
8894b9ad928SBarry Smith 
8909566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
8914b9ad928SBarry Smith   if (side == PC_RIGHT) {
8929566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(pc, x, work));
8939566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->mat, work, y));
8944b9ad928SBarry Smith   } else if (side == PC_LEFT) {
8959566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->mat, x, work));
8969566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(pc, work, y));
8974b9ad928SBarry Smith   }
8984b9ad928SBarry Smith   /* add support for PC_SYMMETRIC */
899e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
9003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9014b9ad928SBarry Smith }
9024b9ad928SBarry Smith 
9034b9ad928SBarry Smith /*@
9044b9ad928SBarry Smith   PCApplyRichardsonExists - Determines whether a particular preconditioner has a
9054b9ad928SBarry Smith   built-in fast application of Richardson's method.
9064b9ad928SBarry Smith 
9074b9ad928SBarry Smith   Not Collective
9084b9ad928SBarry Smith 
9094b9ad928SBarry Smith   Input Parameter:
9104b9ad928SBarry Smith . pc - the preconditioner
9114b9ad928SBarry Smith 
9124b9ad928SBarry Smith   Output Parameter:
913f1580f4eSBarry Smith . exists - `PETSC_TRUE` or `PETSC_FALSE`
9144b9ad928SBarry Smith 
9154b9ad928SBarry Smith   Level: developer
9164b9ad928SBarry Smith 
91739b1ba1bSPierre Jolivet .seealso: [](ch_ksp), `PC`, `KSPRICHARDSON`, `PCApplyRichardson()`
9184b9ad928SBarry Smith @*/
PCApplyRichardsonExists(PC pc,PetscBool * exists)919d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyRichardsonExists(PC pc, PetscBool *exists)
920d71ae5a4SJacob Faibussowitsch {
9214b9ad928SBarry Smith   PetscFunctionBegin;
9220700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9234f572ea9SToby Isaac   PetscAssertPointer(exists, 2);
9244b9ad928SBarry Smith   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
9254b9ad928SBarry Smith   else *exists = PETSC_FALSE;
9263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9274b9ad928SBarry Smith }
9284b9ad928SBarry Smith 
9294b9ad928SBarry Smith /*@
9304b9ad928SBarry Smith   PCApplyRichardson - Applies several steps of Richardson iteration with
9314b9ad928SBarry Smith   the particular preconditioner. This routine is usually used by the
9324b9ad928SBarry Smith   Krylov solvers and not the application code directly.
9334b9ad928SBarry Smith 
934c3339decSBarry Smith   Collective
9354b9ad928SBarry Smith 
9364b9ad928SBarry Smith   Input Parameters:
9370b4b7b1cSBarry Smith + pc        - the `PC` preconditioner context
938dd8e379bSPierre Jolivet . b         - the right-hand side
9394b9ad928SBarry Smith . w         - one work vector
9404b9ad928SBarry Smith . rtol      - relative decrease in residual norm convergence criteria
94170441072SBarry Smith . abstol    - absolute residual norm convergence criteria
9424b9ad928SBarry Smith . dtol      - divergence residual norm increase criteria
9437319c654SBarry Smith . its       - the number of iterations to apply.
9447319c654SBarry Smith - guesszero - if the input x contains nonzero initial guess
9454b9ad928SBarry Smith 
946d8d19677SJose E. Roman   Output Parameters:
9474d0a8057SBarry Smith + outits - number of iterations actually used (for SOR this always equals its)
9484d0a8057SBarry Smith . reason - the reason the apply terminated
949f1580f4eSBarry Smith - y      - the solution (also contains initial guess if guesszero is `PETSC_FALSE`
9504b9ad928SBarry Smith 
95120f4b53cSBarry Smith   Level: developer
95220f4b53cSBarry Smith 
9534b9ad928SBarry Smith   Notes:
9544b9ad928SBarry Smith   Most preconditioners do not support this function. Use the command
955f1580f4eSBarry Smith   `PCApplyRichardsonExists()` to determine if one does.
9564b9ad928SBarry Smith 
957f1580f4eSBarry Smith   Except for the `PCMG` this routine ignores the convergence tolerances
9584b9ad928SBarry Smith   and always runs for the number of iterations
9594b9ad928SBarry Smith 
960562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyRichardsonExists()`
9614b9ad928SBarry Smith @*/
PCApplyRichardson(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt * outits,PCRichardsonConvergedReason * reason)962d71ae5a4SJacob 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)
963d71ae5a4SJacob Faibussowitsch {
9644b9ad928SBarry Smith   PetscFunctionBegin;
9650700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9660700a824SBarry Smith   PetscValidHeaderSpecific(b, VEC_CLASSID, 2);
9670700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
9680700a824SBarry Smith   PetscValidHeaderSpecific(w, VEC_CLASSID, 4);
9697827d75bSBarry Smith   PetscCheck(b != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "b and y must be different vectors");
9709566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
971dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applyrichardson, b, y, w, rtol, abstol, dtol, its, guesszero, outits, reason);
9723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9734b9ad928SBarry Smith }
9744b9ad928SBarry Smith 
975422a814eSBarry Smith /*@
976f1580f4eSBarry Smith   PCSetFailedReason - Sets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
9771b2b9847SBarry Smith 
978c3339decSBarry Smith   Logically Collective
9791b2b9847SBarry Smith 
980d8d19677SJose E. Roman   Input Parameters:
9810b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
9820b4b7b1cSBarry Smith - reason - the reason it failed
9831b2b9847SBarry Smith 
9841b2b9847SBarry Smith   Level: advanced
9851b2b9847SBarry Smith 
986562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCFailedReason`
9871b2b9847SBarry Smith @*/
PCSetFailedReason(PC pc,PCFailedReason reason)988d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFailedReason(PC pc, PCFailedReason reason)
989d71ae5a4SJacob Faibussowitsch {
9901b2b9847SBarry Smith   PetscFunctionBegin;
9916479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9921b2b9847SBarry Smith   pc->failedreason = reason;
9933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9941b2b9847SBarry Smith }
9951b2b9847SBarry Smith 
9961b2b9847SBarry Smith /*@
997f1580f4eSBarry Smith   PCGetFailedReason - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
998422a814eSBarry Smith 
999f1580f4eSBarry Smith   Not Collective
10001b2b9847SBarry Smith 
10011b2b9847SBarry Smith   Input Parameter:
10020b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
10031b2b9847SBarry Smith 
10041b2b9847SBarry Smith   Output Parameter:
10051b2b9847SBarry Smith . reason - the reason it failed
10061b2b9847SBarry Smith 
100720f4b53cSBarry Smith   Level: advanced
100820f4b53cSBarry Smith 
1009f1580f4eSBarry Smith   Note:
10100b4b7b1cSBarry Smith   After a call to `KSPCheckDot()` or  `KSPCheckNorm()` inside a `KSPSolve()` or a call to `PCReduceFailedReason()`
10110b4b7b1cSBarry Smith   this is the maximum reason over all MPI processes in the `PC` communicator and hence logically collective.
10129a6c2652SBarry Smith   Otherwise it returns the local value.
10131b2b9847SBarry Smith 
10140b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCSetFailedReason()`, `PCFailedReason`
10151b2b9847SBarry Smith @*/
PCGetFailedReason(PC pc,PCFailedReason * reason)10169a6c2652SBarry Smith PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason)
1017d71ae5a4SJacob Faibussowitsch {
10181b2b9847SBarry Smith   PetscFunctionBegin;
10196479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1020371d2eb7SMartin Diehl   *reason = pc->failedreason;
10213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10221b2b9847SBarry Smith }
1023422a814eSBarry Smith 
10246479e835SStefano Zampini /*@
10256479e835SStefano Zampini   PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC`
10266479e835SStefano Zampini 
10276479e835SStefano Zampini   Collective
10286479e835SStefano Zampini 
10296479e835SStefano Zampini   Input Parameter:
10300b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
10316479e835SStefano Zampini 
10326479e835SStefano Zampini   Level: advanced
10336479e835SStefano Zampini 
10346479e835SStefano Zampini   Note:
1035562efe2eSBarry Smith   Different MPI processes may have different reasons or no reason, see `PCGetFailedReason()`. This routine
10366479e835SStefano Zampini   makes them have a common value (failure if any MPI process had a failure).
10376479e835SStefano Zampini 
10380b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`, `PCFailedReason`
10396479e835SStefano Zampini @*/
PCReduceFailedReason(PC pc)10406479e835SStefano Zampini PetscErrorCode PCReduceFailedReason(PC pc)
10416479e835SStefano Zampini {
10426479e835SStefano Zampini   PetscInt buf;
10436479e835SStefano Zampini 
10446479e835SStefano Zampini   PetscFunctionBegin;
10456479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
10466479e835SStefano Zampini   buf = (PetscInt)pc->failedreason;
1047462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
10486479e835SStefano Zampini   pc->failedreason = (PCFailedReason)buf;
10496479e835SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
10506479e835SStefano Zampini }
10516479e835SStefano Zampini 
10524b9ad928SBarry Smith /*
10534b9ad928SBarry Smith       a setupcall of 0 indicates never setup,
105423ee1639SBarry Smith                      1 indicates has been previously setup
1055422a814eSBarry Smith                     -1 indicates a PCSetUp() was attempted and failed
10564b9ad928SBarry Smith */
10574b9ad928SBarry Smith /*@
10580b4b7b1cSBarry Smith   PCSetUp - Prepares for the use of a preconditioner. Performs all the one-time operations needed before the preconditioner
10590b4b7b1cSBarry Smith   can be used with `PCApply()`
10604b9ad928SBarry Smith 
1061c3339decSBarry Smith   Collective
10624b9ad928SBarry Smith 
10634b9ad928SBarry Smith   Input Parameter:
10640b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
10654b9ad928SBarry Smith 
10664b9ad928SBarry Smith   Level: developer
10674b9ad928SBarry Smith 
10680b4b7b1cSBarry Smith   Notes:
10690b4b7b1cSBarry Smith   For example, for `PCLU` this will compute the factorization.
10700b4b7b1cSBarry Smith 
10710b4b7b1cSBarry Smith   This is called automatically by `KSPSetUp()` or `PCApply()` so rarely needs to be called directly.
10720b4b7b1cSBarry Smith 
10730b4b7b1cSBarry Smith   For nested preconditioners, such as `PCFIELDSPLIT` or `PCBJACOBI` this may not finish the construction of the preconditioner
10740b4b7b1cSBarry Smith   on the inner levels, the routine `PCSetUpOnBlocks()` may compute more of the preconditioner in those situations.
10750b4b7b1cSBarry Smith 
10760b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `KSPSetUp()`, `PCSetUpOnBlocks()`
10774b9ad928SBarry Smith @*/
PCSetUp(PC pc)1078d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp(PC pc)
1079d71ae5a4SJacob Faibussowitsch {
1080566e8bf2SBarry Smith   const char      *def;
1081fc9b51d3SBarry Smith   PetscObjectState matstate, matnonzerostate;
10824b9ad928SBarry Smith 
10834b9ad928SBarry Smith   PetscFunctionBegin;
10840700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
108528b400f6SJacob Faibussowitsch   PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first");
10866ce5e81cSLisandro Dalcin 
108723ee1639SBarry Smith   if (pc->setupcalled && pc->reusepreconditioner) {
10889566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n"));
10893ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10904b9ad928SBarry Smith   }
10914b9ad928SBarry Smith 
10929566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate));
10939566063dSJacob Faibussowitsch   PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate));
109423ee1639SBarry Smith   if (!pc->setupcalled) {
10955e62d202SMark Adams     //PetscCall(PetscInfo(pc, "Setting up PC for first time\n"));
109623ee1639SBarry Smith     pc->flag = DIFFERENT_NONZERO_PATTERN;
10979df67409SStefano Zampini   } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS);
10989df67409SStefano Zampini   else {
10999df67409SStefano Zampini     if (matnonzerostate != pc->matnonzerostate) {
11009566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n"));
110123ee1639SBarry Smith       pc->flag = DIFFERENT_NONZERO_PATTERN;
110223ee1639SBarry Smith     } else {
11035e62d202SMark Adams       //PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n"));
110423ee1639SBarry Smith       pc->flag = SAME_NONZERO_PATTERN;
110523ee1639SBarry Smith     }
110623ee1639SBarry Smith   }
110723ee1639SBarry Smith   pc->matstate        = matstate;
1108fc9b51d3SBarry Smith   pc->matnonzerostate = matnonzerostate;
110923ee1639SBarry Smith 
11107adad957SLisandro Dalcin   if (!((PetscObject)pc)->type_name) {
11119566063dSJacob Faibussowitsch     PetscCall(PCGetDefaultType_Private(pc, &def));
11129566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, def));
1113566e8bf2SBarry Smith   }
11144b9ad928SBarry Smith 
11159566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
11169566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure));
11179566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0));
11184b9ad928SBarry Smith   if (pc->ops->setup) {
1119fe57bb1aSStefano Zampini     PetscCall(PCLogEventsDeactivatePush());
1120dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, setup);
1121fe57bb1aSStefano Zampini     PetscCall(PCLogEventsDeactivatePop());
11224b9ad928SBarry Smith   }
11239566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0));
1124ab76df0cSBarry Smith   if (pc->postsetup) PetscCall((*pc->postsetup)(pc));
1125371d2eb7SMartin Diehl   if (!pc->setupcalled) pc->setupcalled = PETSC_TRUE;
11263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11274b9ad928SBarry Smith }
11284b9ad928SBarry Smith 
11294b9ad928SBarry Smith /*@
11304b9ad928SBarry Smith   PCSetUpOnBlocks - Sets up the preconditioner for each block in
113173716367SStefano Zampini   the block Jacobi, overlapping Schwarz, and fieldsplit methods.
11324b9ad928SBarry Smith 
1133c3339decSBarry Smith   Collective
11344b9ad928SBarry Smith 
1135f1580f4eSBarry Smith   Input Parameter:
11360b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
11374b9ad928SBarry Smith 
11384b9ad928SBarry Smith   Level: developer
11394b9ad928SBarry Smith 
114073716367SStefano Zampini   Notes:
114173716367SStefano Zampini   For nested preconditioners such as `PCBJACOBI`, `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is
1142f1580f4eSBarry Smith   called on the outer `PC`, this routine ensures it is called.
1143f1580f4eSBarry Smith 
114473716367SStefano Zampini   It calls `PCSetUp()` if not yet called.
114573716367SStefano Zampini 
1146562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()`
11474b9ad928SBarry Smith @*/
PCSetUpOnBlocks(PC pc)1148d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUpOnBlocks(PC pc)
1149d71ae5a4SJacob Faibussowitsch {
11504b9ad928SBarry Smith   PetscFunctionBegin;
11510700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
115273716367SStefano Zampini   if (!pc->setupcalled) PetscCall(PCSetUp(pc)); /* "if" to prevent -info extra prints */
11533ba16761SJacob Faibussowitsch   if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS);
115421604f62SStefano Zampini   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
11559566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0));
1156fe57bb1aSStefano Zampini   PetscCall(PCLogEventsDeactivatePush());
1157dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, setuponblocks);
1158fe57bb1aSStefano Zampini   PetscCall(PCLogEventsDeactivatePop());
11599566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0));
11603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11614b9ad928SBarry Smith }
11624b9ad928SBarry Smith 
11634b9ad928SBarry Smith /*@C
11644b9ad928SBarry Smith   PCSetModifySubMatrices - Sets a user-defined routine for modifying the
116504c3f3b8SBarry Smith   submatrices that arise within certain subdomain-based preconditioners such as `PCASM`
11664b9ad928SBarry Smith 
1167c3339decSBarry Smith   Logically Collective
11684b9ad928SBarry Smith 
11694b9ad928SBarry Smith   Input Parameters:
11700b4b7b1cSBarry Smith + pc   - the `PC` preconditioner context
11714d4d2bdcSBarry Smith . func - routine for modifying the submatrices, see `PCModifySubMatricesFn`
117204c3f3b8SBarry Smith - ctx  - optional user-defined context (may be `NULL`)
11734b9ad928SBarry Smith 
117420f4b53cSBarry Smith   Level: advanced
117520f4b53cSBarry Smith 
11764b9ad928SBarry Smith   Notes:
11777addb90fSBarry Smith   The basic submatrices are extracted from the matrix used to construct the preconditioner as
117804c3f3b8SBarry Smith   usual; the user can then alter these (for example, to set different boundary
117904c3f3b8SBarry Smith   conditions for each submatrix) before they are used for the local solves.
118004c3f3b8SBarry Smith 
1181f1580f4eSBarry Smith   `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and
1182f1580f4eSBarry Smith   `KSPSolve()`.
11834b9ad928SBarry Smith 
1184f1580f4eSBarry Smith   A routine set by `PCSetModifySubMatrices()` is currently called within
1185*80e22091SPierre Jolivet   `PCBJACOBI`, `PCASM`, `PCGASM`, and `PCHPDDM`.
1186*80e22091SPierre Jolivet   All other preconditioners ignore this routine.
11874b9ad928SBarry Smith 
11884d4d2bdcSBarry Smith .seealso: [](ch_ksp), `PC`, `PCModifySubMatricesFn`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()`
11894b9ad928SBarry Smith @*/
PCSetModifySubMatrices(PC pc,PCModifySubMatricesFn * func,PetscCtx ctx)11902a8381b2SBarry Smith PetscErrorCode PCSetModifySubMatrices(PC pc, PCModifySubMatricesFn *func, PetscCtx ctx)
1191d71ae5a4SJacob Faibussowitsch {
11924b9ad928SBarry Smith   PetscFunctionBegin;
11930700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11944b9ad928SBarry Smith   pc->modifysubmatrices  = func;
11954b9ad928SBarry Smith   pc->modifysubmatricesP = ctx;
11963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11974b9ad928SBarry Smith }
11984b9ad928SBarry Smith 
11994b9ad928SBarry Smith /*@C
12004b9ad928SBarry Smith   PCModifySubMatrices - Calls an optional user-defined routine within
1201f1580f4eSBarry Smith   certain preconditioners if one has been set with `PCSetModifySubMatrices()`.
12024b9ad928SBarry Smith 
1203c3339decSBarry Smith   Collective
12044b9ad928SBarry Smith 
12054b9ad928SBarry Smith   Input Parameters:
12060b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
12074b9ad928SBarry Smith . nsub   - the number of local submatrices
12084b9ad928SBarry Smith . row    - an array of index sets that contain the global row numbers
12094b9ad928SBarry Smith          that comprise each local submatrix
12104b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
12114b9ad928SBarry Smith          that comprise each local submatrix
12124b9ad928SBarry Smith . submat - array of local submatrices
12134b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
1214562efe2eSBarry Smith          user-defined routine (may be `NULL`)
12154b9ad928SBarry Smith 
12164b9ad928SBarry Smith   Output Parameter:
12174b9ad928SBarry Smith . submat - array of local submatrices (the entries of which may
12184b9ad928SBarry Smith             have been modified)
12194b9ad928SBarry Smith 
122020f4b53cSBarry Smith   Level: developer
122120f4b53cSBarry Smith 
122204c3f3b8SBarry Smith   Note:
12234b9ad928SBarry Smith   The user should NOT generally call this routine, as it will
122404c3f3b8SBarry Smith   automatically be called within certain preconditioners.
12254b9ad928SBarry Smith 
12264d4d2bdcSBarry Smith .seealso: [](ch_ksp), `PC`, `PCModifySubMatricesFn`, `PCSetModifySubMatrices()`
12274b9ad928SBarry Smith @*/
PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],PetscCtx ctx)12282a8381b2SBarry Smith PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], PetscCtx ctx)
1229d71ae5a4SJacob Faibussowitsch {
12304b9ad928SBarry Smith   PetscFunctionBegin;
12310700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12323ba16761SJacob Faibussowitsch   if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS);
12339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0));
12349566063dSJacob Faibussowitsch   PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx));
12359566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0));
12363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12374b9ad928SBarry Smith }
12384b9ad928SBarry Smith 
12394b9ad928SBarry Smith /*@
12404b9ad928SBarry Smith   PCSetOperators - Sets the matrix associated with the linear system and
12411c02f50aSBarry Smith   a (possibly) different one from which the preconditioner will be constructed.
12424b9ad928SBarry Smith 
1243c3339decSBarry Smith   Logically Collective
12444b9ad928SBarry Smith 
12454b9ad928SBarry Smith   Input Parameters:
12460b4b7b1cSBarry Smith + pc   - the `PC` preconditioner context
1247e5d3d808SBarry Smith . Amat - the matrix that defines the linear system
1248d1e9a80fSBarry Smith - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
12494b9ad928SBarry Smith 
12501c02f50aSBarry Smith   Level: advanced
1251189c0b8aSBarry Smith 
125220f4b53cSBarry Smith   Notes:
12531c02f50aSBarry Smith   Using this routine directly is rarely needed, the preferred, and equivalent, usage is `KSPSetOperators()`.
12541c02f50aSBarry Smith 
125520f4b53cSBarry Smith   Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used.
125620f4b53cSBarry Smith 
125720f4b53cSBarry Smith   If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then
1258f1580f4eSBarry Smith   first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()`
1259f1580f4eSBarry Smith   on it and then pass it back in in your call to `KSPSetOperators()`.
1260189c0b8aSBarry Smith 
12614b9ad928SBarry Smith   More Notes about Repeated Solution of Linear Systems:
126220f4b53cSBarry Smith   PETSc does NOT reset the matrix entries of either `Amat` or `Pmat`
12634b9ad928SBarry Smith   to zero after a linear solve; the user is completely responsible for
1264f1580f4eSBarry Smith   matrix assembly.  See the routine `MatZeroEntries()` if desiring to
12654b9ad928SBarry Smith   zero all elements of a matrix.
12664b9ad928SBarry Smith 
1267562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`
12684b9ad928SBarry Smith  @*/
PCSetOperators(PC pc,Mat Amat,Mat Pmat)1269d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat)
1270d71ae5a4SJacob Faibussowitsch {
12713b3f6333SBarry Smith   PetscInt m1, n1, m2, n2;
12724b9ad928SBarry Smith 
12734b9ad928SBarry Smith   PetscFunctionBegin;
12740700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12750700a824SBarry Smith   if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
12760700a824SBarry Smith   if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3);
12773fc8bf9cSBarry Smith   if (Amat) PetscCheckSameComm(pc, 1, Amat, 2);
12783fc8bf9cSBarry Smith   if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3);
127931641f1aSBarry Smith   if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) {
12809566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Amat, &m1, &n1));
12819566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->mat, &m2, &n2));
128263a3b9bcSJacob 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);
12839566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Pmat, &m1, &n1));
12849566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2));
128563a3b9bcSJacob 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);
12863b3f6333SBarry Smith   }
12874b9ad928SBarry Smith 
128823ee1639SBarry Smith   if (Pmat != pc->pmat) {
128923ee1639SBarry Smith     /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
129023ee1639SBarry Smith     pc->matnonzerostate = -1;
129123ee1639SBarry Smith     pc->matstate        = -1;
129223ee1639SBarry Smith   }
129323ee1639SBarry Smith 
1294906ed7ccSBarry Smith   /* reference first in case the matrices are the same */
12959566063dSJacob Faibussowitsch   if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat));
12969566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->mat));
12979566063dSJacob Faibussowitsch   if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat));
12989566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->pmat));
12994b9ad928SBarry Smith   pc->mat  = Amat;
13004b9ad928SBarry Smith   pc->pmat = Pmat;
13013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1302e56f5c9eSBarry Smith }
1303e56f5c9eSBarry Smith 
130423ee1639SBarry Smith /*@
13050b4b7b1cSBarry Smith   PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner `PC` has changed.
130623ee1639SBarry Smith 
1307c3339decSBarry Smith   Logically Collective
130823ee1639SBarry Smith 
130923ee1639SBarry Smith   Input Parameters:
13100b4b7b1cSBarry Smith + pc   - the `PC` preconditioner context
1311f1580f4eSBarry Smith - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
131223ee1639SBarry Smith 
13132b26979fSBarry Smith   Level: intermediate
13142b26979fSBarry Smith 
1315f1580f4eSBarry Smith   Note:
1316f1580f4eSBarry Smith   Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option
1317f1580f4eSBarry Smith   prevents this.
1318f1580f4eSBarry Smith 
1319562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()`
132023ee1639SBarry Smith  @*/
PCSetReusePreconditioner(PC pc,PetscBool flag)1321d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag)
1322d71ae5a4SJacob Faibussowitsch {
132323ee1639SBarry Smith   PetscFunctionBegin;
132423ee1639SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1325f9177522SStefano Zampini   PetscValidLogicalCollectiveBool(pc, flag, 2);
132623ee1639SBarry Smith   pc->reusepreconditioner = flag;
1327bbbcb081SMark Adams   PetscTryMethod(pc, "PCSetReusePreconditioner_C", (PC, PetscBool), (pc, flag));
13283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13294b9ad928SBarry Smith }
13304b9ad928SBarry Smith 
1331c60c7ad4SBarry Smith /*@
1332f1580f4eSBarry Smith   PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed.
1333c60c7ad4SBarry Smith 
1334bba28a21SBarry Smith   Not Collective
1335c60c7ad4SBarry Smith 
1336c60c7ad4SBarry Smith   Input Parameter:
13370b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
1338c60c7ad4SBarry Smith 
1339c60c7ad4SBarry Smith   Output Parameter:
1340f1580f4eSBarry Smith . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
1341c60c7ad4SBarry Smith 
1342d0418729SSatish Balay   Level: intermediate
1343d0418729SSatish Balay 
1344562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()`
1345c60c7ad4SBarry Smith  @*/
PCGetReusePreconditioner(PC pc,PetscBool * flag)1346d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag)
1347d71ae5a4SJacob Faibussowitsch {
1348c60c7ad4SBarry Smith   PetscFunctionBegin;
1349c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
13504f572ea9SToby Isaac   PetscAssertPointer(flag, 2);
1351c60c7ad4SBarry Smith   *flag = pc->reusepreconditioner;
13523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1353c60c7ad4SBarry Smith }
1354c60c7ad4SBarry Smith 
1355487a658cSBarry Smith /*@
13564b9ad928SBarry Smith   PCGetOperators - Gets the matrix associated with the linear system and
13570b4b7b1cSBarry Smith   possibly a different one which is used to construct the preconditioner.
13584b9ad928SBarry Smith 
1359562efe2eSBarry Smith   Not Collective, though parallel `Mat`s are returned if `pc` is parallel
13604b9ad928SBarry Smith 
13614b9ad928SBarry Smith   Input Parameter:
13620b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
13634b9ad928SBarry Smith 
13644b9ad928SBarry Smith   Output Parameters:
1365e5d3d808SBarry Smith + Amat - the matrix defining the linear system
136623ee1639SBarry Smith - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
13674b9ad928SBarry Smith 
13684b9ad928SBarry Smith   Level: intermediate
13694b9ad928SBarry Smith 
1370f1580f4eSBarry Smith   Note:
137195452b02SPatrick Sanan   Does not increase the reference count of the matrices, so you should not destroy them
1372298cc208SBarry Smith 
1373f1580f4eSBarry Smith   Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators
1374f1580f4eSBarry Smith   are created in `PC` and returned to the user. In this case, if both operators
137573950996SBarry Smith   mat and pmat are requested, two DIFFERENT operators will be returned. If
137673950996SBarry Smith   only one is requested both operators in the PC will be the same (i.e. as
1377f1580f4eSBarry Smith   if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats).
137873950996SBarry Smith   The user must set the sizes of the returned matrices and their type etc just
1379f1580f4eSBarry Smith   as if the user created them with `MatCreate()`. For example,
138073950996SBarry Smith 
1381f1580f4eSBarry Smith .vb
1382f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1383f1580f4eSBarry Smith            set size, type, etc of Amat
138473950996SBarry Smith 
1385f1580f4eSBarry Smith          MatCreate(comm,&mat);
1386f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Amat);
1387f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)mat);
1388f1580f4eSBarry Smith            set size, type, etc of Amat
1389f1580f4eSBarry Smith .ve
139073950996SBarry Smith 
139173950996SBarry Smith   and
139273950996SBarry Smith 
1393f1580f4eSBarry Smith .vb
1394f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1395f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
139673950996SBarry Smith 
1397f1580f4eSBarry Smith          MatCreate(comm,&Amat);
1398f1580f4eSBarry Smith          MatCreate(comm,&Pmat);
1399f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1400f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Amat);
1401f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Pmat);
1402f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
1403f1580f4eSBarry Smith .ve
140473950996SBarry Smith 
1405f1580f4eSBarry Smith   The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy
1406b8d709abSRichard Tran Mills   of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely
1407f1580f4eSBarry Smith   managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look
1408f1580f4eSBarry Smith   at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to
1409f1580f4eSBarry Smith   the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP
1410f1580f4eSBarry Smith   you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you).
1411f1580f4eSBarry Smith   Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when
141273950996SBarry Smith   it can be created for you?
141373950996SBarry Smith 
1414562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()`
14154b9ad928SBarry Smith @*/
PCGetOperators(PC pc,Mat * Amat,Mat * Pmat)1416d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat)
1417d71ae5a4SJacob Faibussowitsch {
14184b9ad928SBarry Smith   PetscFunctionBegin;
14190700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1420e5d3d808SBarry Smith   if (Amat) {
142173950996SBarry Smith     if (!pc->mat) {
14229fca8c71SStefano Zampini       if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */
14239a4708feSJed Brown         pc->mat = pc->pmat;
14249566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->mat));
1425e5d3d808SBarry Smith       } else { /* both Amat and Pmat are empty */
14269566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat));
1427e5d3d808SBarry Smith         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
142873950996SBarry Smith           pc->pmat = pc->mat;
14299566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->pmat));
143073950996SBarry Smith         }
143173950996SBarry Smith       }
14329a4708feSJed Brown     }
1433e5d3d808SBarry Smith     *Amat = pc->mat;
143473950996SBarry Smith   }
1435e5d3d808SBarry Smith   if (Pmat) {
143673950996SBarry Smith     if (!pc->pmat) {
1437e5d3d808SBarry Smith       if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */
14389a4708feSJed Brown         pc->pmat = pc->mat;
14399566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->pmat));
14409a4708feSJed Brown       } else {
14419566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat));
1442e5d3d808SBarry Smith         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
144373950996SBarry Smith           pc->mat = pc->pmat;
14449566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->mat));
144573950996SBarry Smith         }
144673950996SBarry Smith       }
14479a4708feSJed Brown     }
1448e5d3d808SBarry Smith     *Pmat = pc->pmat;
144973950996SBarry Smith   }
14503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14514b9ad928SBarry Smith }
14524b9ad928SBarry Smith 
14535d83a8b1SBarry Smith /*@
1454906ed7ccSBarry Smith   PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1455f1580f4eSBarry Smith   possibly a different one associated with the preconditioner have been set in the `PC`.
1456906ed7ccSBarry Smith 
145720f4b53cSBarry Smith   Not Collective, though the results on all processes should be the same
1458906ed7ccSBarry Smith 
1459906ed7ccSBarry Smith   Input Parameter:
14600b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
1461906ed7ccSBarry Smith 
1462906ed7ccSBarry Smith   Output Parameters:
1463906ed7ccSBarry Smith + mat  - the matrix associated with the linear system was set
1464906ed7ccSBarry Smith - pmat - matrix associated with the preconditioner was set, usually the same
1465906ed7ccSBarry Smith 
1466906ed7ccSBarry Smith   Level: intermediate
1467906ed7ccSBarry Smith 
1468562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`
1469906ed7ccSBarry Smith @*/
PCGetOperatorsSet(PC pc,PetscBool * mat,PetscBool * pmat)1470d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat)
1471d71ae5a4SJacob Faibussowitsch {
1472906ed7ccSBarry Smith   PetscFunctionBegin;
14730700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1474906ed7ccSBarry Smith   if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1475906ed7ccSBarry Smith   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
14763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1477906ed7ccSBarry Smith }
1478906ed7ccSBarry Smith 
1479f39d8e23SSatish Balay /*@
1480a4fd02acSBarry Smith   PCFactorGetMatrix - Gets the factored matrix from the
1481f1580f4eSBarry Smith   preconditioner context.  This routine is valid only for the `PCLU`,
1482f1580f4eSBarry Smith   `PCILU`, `PCCHOLESKY`, and `PCICC` methods.
14834b9ad928SBarry Smith 
1484562efe2eSBarry Smith   Not Collective though `mat` is parallel if `pc` is parallel
14854b9ad928SBarry Smith 
1486f1580f4eSBarry Smith   Input Parameter:
14870b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
14884b9ad928SBarry Smith 
1489feefa0e1SJacob Faibussowitsch   Output Parameters:
14904b9ad928SBarry Smith . mat - the factored matrix
14914b9ad928SBarry Smith 
14924b9ad928SBarry Smith   Level: advanced
14934b9ad928SBarry Smith 
1494f1580f4eSBarry Smith   Note:
1495562efe2eSBarry Smith   Does not increase the reference count for `mat` so DO NOT destroy it
14969405f653SBarry Smith 
1497562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC`
14984b9ad928SBarry Smith @*/
PCFactorGetMatrix(PC pc,Mat * mat)1499d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat)
1500d71ae5a4SJacob Faibussowitsch {
15014b9ad928SBarry Smith   PetscFunctionBegin;
15020700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15034f572ea9SToby Isaac   PetscAssertPointer(mat, 2);
15047def7855SStefano Zampini   PetscCall(PCFactorSetUpMatSolverType(pc));
1505dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, getfactoredmatrix, mat);
15063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15074b9ad928SBarry Smith }
15084b9ad928SBarry Smith 
1509cc4c1da9SBarry Smith /*@
15104b9ad928SBarry Smith   PCSetOptionsPrefix - Sets the prefix used for searching for all
1511f1580f4eSBarry Smith   `PC` options in the database.
15124b9ad928SBarry Smith 
1513c3339decSBarry Smith   Logically Collective
15144b9ad928SBarry Smith 
15154b9ad928SBarry Smith   Input Parameters:
15160b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
1517f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
15184b9ad928SBarry Smith 
1519f1580f4eSBarry Smith   Note:
15204b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
15214b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
15224b9ad928SBarry Smith   hyphen.
15234b9ad928SBarry Smith 
15244b9ad928SBarry Smith   Level: advanced
15254b9ad928SBarry Smith 
1526562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()`
15274b9ad928SBarry Smith @*/
PCSetOptionsPrefix(PC pc,const char prefix[])1528d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[])
1529d71ae5a4SJacob Faibussowitsch {
15304b9ad928SBarry Smith   PetscFunctionBegin;
15310700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15329566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix));
15333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15344b9ad928SBarry Smith }
15354b9ad928SBarry Smith 
1536cc4c1da9SBarry Smith /*@
15374b9ad928SBarry Smith   PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1538f1580f4eSBarry Smith   `PC` options in the database.
15394b9ad928SBarry Smith 
1540c3339decSBarry Smith   Logically Collective
15414b9ad928SBarry Smith 
15424b9ad928SBarry Smith   Input Parameters:
15430b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
1544f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
15454b9ad928SBarry Smith 
1546f1580f4eSBarry Smith   Note:
15474b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
15484b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
15494b9ad928SBarry Smith   hyphen.
15504b9ad928SBarry Smith 
15514b9ad928SBarry Smith   Level: advanced
15524b9ad928SBarry Smith 
1553562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()`
15544b9ad928SBarry Smith @*/
PCAppendOptionsPrefix(PC pc,const char prefix[])1555d71ae5a4SJacob Faibussowitsch PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[])
1556d71ae5a4SJacob Faibussowitsch {
15574b9ad928SBarry Smith   PetscFunctionBegin;
15580700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15599566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix));
15603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15614b9ad928SBarry Smith }
15624b9ad928SBarry Smith 
1563cc4c1da9SBarry Smith /*@
15644b9ad928SBarry Smith   PCGetOptionsPrefix - Gets the prefix used for searching for all
15654b9ad928SBarry Smith   PC options in the database.
15664b9ad928SBarry Smith 
15674b9ad928SBarry Smith   Not Collective
15684b9ad928SBarry Smith 
1569f1580f4eSBarry Smith   Input Parameter:
15700b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
15714b9ad928SBarry Smith 
1572f1580f4eSBarry Smith   Output Parameter:
15734b9ad928SBarry Smith . prefix - pointer to the prefix string used, is returned
15744b9ad928SBarry Smith 
15754b9ad928SBarry Smith   Level: advanced
15764b9ad928SBarry Smith 
1577562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()`
15784b9ad928SBarry Smith @*/
PCGetOptionsPrefix(PC pc,const char * prefix[])1579d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[])
1580d71ae5a4SJacob Faibussowitsch {
15814b9ad928SBarry Smith   PetscFunctionBegin;
15820700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15834f572ea9SToby Isaac   PetscAssertPointer(prefix, 2);
15849566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix));
15853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15864b9ad928SBarry Smith }
15874b9ad928SBarry Smith 
15888066bbecSBarry Smith /*
1589dd8e379bSPierre Jolivet    Indicates the right-hand side will be changed by KSPSolve(), this occurs for a few
15908066bbecSBarry Smith   preconditioners including BDDC and Eisentat that transform the equations before applying
15918066bbecSBarry Smith   the Krylov methods
15928066bbecSBarry Smith */
PCPreSolveChangeRHS(PC pc,PetscBool * change)1593d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change)
1594d71ae5a4SJacob Faibussowitsch {
1595a06fd7f2SStefano Zampini   PetscFunctionBegin;
1596a06fd7f2SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15974f572ea9SToby Isaac   PetscAssertPointer(change, 2);
1598a06fd7f2SStefano Zampini   *change = PETSC_FALSE;
1599cac4c232SBarry Smith   PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change));
16003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1601a06fd7f2SStefano Zampini }
1602a06fd7f2SStefano Zampini 
16034b9ad928SBarry Smith /*@
16040b4b7b1cSBarry Smith   PCPreSolve - Optional pre-solve phase, intended for any preconditioner-specific actions that must be performed before
16050b4b7b1cSBarry Smith   the iterative solve itself. Used in conjunction with `PCPostSolve()`
16064b9ad928SBarry Smith 
1607c3339decSBarry Smith   Collective
16084b9ad928SBarry Smith 
16094b9ad928SBarry Smith   Input Parameters:
16100b4b7b1cSBarry Smith + pc  - the `PC` preconditioner context
16114b9ad928SBarry Smith - ksp - the Krylov subspace context
16124b9ad928SBarry Smith 
16134b9ad928SBarry Smith   Level: developer
16144b9ad928SBarry Smith 
16154b9ad928SBarry Smith   Notes:
1616f1580f4eSBarry Smith   `KSPSolve()` calls this directly, so is rarely called by the user.
16174b9ad928SBarry Smith 
16180b4b7b1cSBarry Smith   Certain preconditioners, such as the `PCType` of `PCEISENSTAT`, change the formulation of the linear system to be solved iteratively.
16190b4b7b1cSBarry Smith   This function performs that transformation. `PCPostSolve()` then transforms the system back to its original form after the solve.
16200b4b7b1cSBarry Smith   `PCPostSolve()` also transforms the resulting solution of the transformed system to the solution of the original problem.
16210b4b7b1cSBarry Smith 
1622ab76df0cSBarry Smith   `KSPSetPostSolve()` provides an alternative way to provide such transformations.
16230b4b7b1cSBarry Smith 
1624ab76df0cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCPostSolve()`, `KSP`, `PCSetPostSetUp()`, `KSPSetPreSolve()`, `KSPSetPostSolve()`
16254b9ad928SBarry Smith @*/
PCPreSolve(PC pc,KSP ksp)1626d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPreSolve(PC pc, KSP ksp)
1627d71ae5a4SJacob Faibussowitsch {
16284b9ad928SBarry Smith   Vec x, rhs;
16294b9ad928SBarry Smith 
16304b9ad928SBarry Smith   PetscFunctionBegin;
16310700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16320700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1633d9a03883SBarry Smith   pc->presolvedone++;
16347827d75bSBarry Smith   PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice");
16359566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
16369566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
1637ab76df0cSBarry Smith   PetscTryTypeMethod(pc, presolve, ksp, rhs, x);
16383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16394b9ad928SBarry Smith }
16404b9ad928SBarry Smith 
1641f560b561SHong Zhang /*@C
1642ab76df0cSBarry Smith   PCSetPostSetUp - Sets function called at the end of `PCSetUp()` to adjust the computed preconditioner
1643f560b561SHong Zhang 
1644c3339decSBarry Smith   Logically Collective
1645f560b561SHong Zhang 
1646f560b561SHong Zhang   Input Parameters:
1647f560b561SHong Zhang + pc        - the preconditioner object
1648ab76df0cSBarry Smith - postsetup - the function to call after `PCSetUp()`
1649f560b561SHong Zhang 
1650ab76df0cSBarry Smith   Calling sequence of `postsetup`:
1651ab76df0cSBarry Smith . pc - the `PC` context
1652f560b561SHong Zhang 
1653f560b561SHong Zhang   Level: developer
1654f560b561SHong Zhang 
1655ab76df0cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`
1656f560b561SHong Zhang @*/
PCSetPostSetUp(PC pc,PetscErrorCode (* postsetup)(PC pc))1657ab76df0cSBarry Smith PetscErrorCode PCSetPostSetUp(PC pc, PetscErrorCode (*postsetup)(PC pc))
1658d71ae5a4SJacob Faibussowitsch {
1659f560b561SHong Zhang   PetscFunctionBegin;
1660f560b561SHong Zhang   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1661ab76df0cSBarry Smith   pc->postsetup = postsetup;
16623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1663f560b561SHong Zhang }
1664f560b561SHong Zhang 
16654b9ad928SBarry Smith /*@
16664b9ad928SBarry Smith   PCPostSolve - Optional post-solve phase, intended for any
16674b9ad928SBarry Smith   preconditioner-specific actions that must be performed after
16684b9ad928SBarry Smith   the iterative solve itself.
16694b9ad928SBarry Smith 
1670c3339decSBarry Smith   Collective
16714b9ad928SBarry Smith 
16724b9ad928SBarry Smith   Input Parameters:
16730b4b7b1cSBarry Smith + pc  - the `PC` preconditioner context
16740b4b7b1cSBarry Smith - ksp - the `KSP` Krylov subspace context
16754b9ad928SBarry Smith 
1676feefa0e1SJacob Faibussowitsch   Example Usage:
16774b9ad928SBarry Smith .vb
16784b9ad928SBarry Smith     PCPreSolve(pc,ksp);
167923ce1328SBarry Smith     KSPSolve(ksp,b,x);
16804b9ad928SBarry Smith     PCPostSolve(pc,ksp);
16814b9ad928SBarry Smith .ve
16824b9ad928SBarry Smith 
1683562efe2eSBarry Smith   Level: developer
1684562efe2eSBarry Smith 
16854b9ad928SBarry Smith   Note:
1686f1580f4eSBarry Smith   `KSPSolve()` calls this routine directly, so it is rarely called by the user.
16874b9ad928SBarry Smith 
16881a02b73cSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPSetPostSolve()`, `KSPSetPreSolve()`, `PCPreSolve()`, `KSPSolve()`
16894b9ad928SBarry Smith @*/
PCPostSolve(PC pc,KSP ksp)1690d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPostSolve(PC pc, KSP ksp)
1691d71ae5a4SJacob Faibussowitsch {
16924b9ad928SBarry Smith   Vec x, rhs;
16934b9ad928SBarry Smith 
16944b9ad928SBarry Smith   PetscFunctionBegin;
16950700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16960700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1697d9a03883SBarry Smith   pc->presolvedone--;
16989566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
16999566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
1700dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, postsolve, ksp, rhs, x);
17013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17024b9ad928SBarry Smith }
17034b9ad928SBarry Smith 
1704ffeef943SBarry Smith /*@
1705f1580f4eSBarry Smith   PCLoad - Loads a `PC` that has been stored in binary  with `PCView()`.
170655849f57SBarry Smith 
1707c3339decSBarry Smith   Collective
170855849f57SBarry Smith 
170955849f57SBarry Smith   Input Parameters:
1710f1580f4eSBarry Smith + newdm  - the newly loaded `PC`, this needs to have been created with `PCCreate()` or
1711f1580f4eSBarry Smith            some related function before a call to `PCLoad()`.
17120b4b7b1cSBarry Smith - viewer - binary file viewer `PETSCVIEWERBINARY`, obtained from `PetscViewerBinaryOpen()`
171355849f57SBarry Smith 
171455849f57SBarry Smith   Level: intermediate
171555849f57SBarry Smith 
1716f1580f4eSBarry Smith   Note:
1717562efe2eSBarry Smith   The type is determined by the data in the file, any `PCType` set into the `PC` before this call is ignored.
171855849f57SBarry Smith 
17190b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()`, `PETSCVIEWERBINARY`
172055849f57SBarry Smith @*/
PCLoad(PC newdm,PetscViewer viewer)1721d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLoad(PC newdm, PetscViewer viewer)
1722d71ae5a4SJacob Faibussowitsch {
172355849f57SBarry Smith   PetscBool isbinary;
1724060da220SMatthew G. Knepley   PetscInt  classid;
172555849f57SBarry Smith   char      type[256];
172655849f57SBarry Smith 
172755849f57SBarry Smith   PetscFunctionBegin;
172855849f57SBarry Smith   PetscValidHeaderSpecific(newdm, PC_CLASSID, 1);
172955849f57SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
17309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
173128b400f6SJacob Faibussowitsch   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
173255849f57SBarry Smith 
17339566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
17347827d75bSBarry Smith   PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file");
17359566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
17369566063dSJacob Faibussowitsch   PetscCall(PCSetType(newdm, type));
1737dbbe0bcdSBarry Smith   PetscTryTypeMethod(newdm, load, viewer);
17383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
173955849f57SBarry Smith }
174055849f57SBarry Smith 
17419804daf3SBarry Smith #include <petscdraw.h>
1742e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1743e04113cfSBarry Smith   #include <petscviewersaws.h>
17440acecf5bSBarry Smith #endif
1745fe2efc57SMark 
1746ffeef943SBarry Smith /*@
17470b4b7b1cSBarry Smith   PCViewFromOptions - View (print or provide information about) the `PC`, based on options in the options database
1748fe2efc57SMark 
1749c3339decSBarry Smith   Collective
1750fe2efc57SMark 
1751fe2efc57SMark   Input Parameters:
175220f4b53cSBarry Smith + A    - the `PC` context
1753f1580f4eSBarry Smith . obj  - Optional object that provides the options prefix
17540b4b7b1cSBarry Smith - name - command line option name
1755fe2efc57SMark 
17560b4b7b1cSBarry Smith   Level: developer
1757f1580f4eSBarry Smith 
1758562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()`
1759fe2efc57SMark @*/
PCViewFromOptions(PC A,PetscObject obj,const char name[])1760d71ae5a4SJacob Faibussowitsch PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[])
1761d71ae5a4SJacob Faibussowitsch {
1762fe2efc57SMark   PetscFunctionBegin;
1763fe2efc57SMark   PetscValidHeaderSpecific(A, PC_CLASSID, 1);
17649566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
17653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1766fe2efc57SMark }
1767fe2efc57SMark 
1768ffeef943SBarry Smith /*@
1769f1580f4eSBarry Smith   PCView - Prints information about the `PC`
17704b9ad928SBarry Smith 
1771c3339decSBarry Smith   Collective
17724b9ad928SBarry Smith 
17734b9ad928SBarry Smith   Input Parameters:
17740b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
17750b4b7b1cSBarry Smith - viewer - optional `PetscViewer` visualization context
17764b9ad928SBarry Smith 
17770b4b7b1cSBarry Smith   Level: intermediate
177820f4b53cSBarry Smith 
1779f1580f4eSBarry Smith   Notes:
17804b9ad928SBarry Smith   The available visualization contexts include
1781f1580f4eSBarry Smith +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
1782f1580f4eSBarry Smith -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
17834b9ad928SBarry Smith   output where only the first processor opens
17844b9ad928SBarry Smith   the file. All other processors send their
17854b9ad928SBarry Smith   data to the first processor to print.
17864b9ad928SBarry Smith 
17874b9ad928SBarry Smith   The user can open an alternative visualization contexts with
1788f1580f4eSBarry Smith   `PetscViewerASCIIOpen()` (output to a specified file).
17894b9ad928SBarry Smith 
17900b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewer`, `PetscViewerType`, `KSPView()`, `PetscViewerASCIIOpen()`
17914b9ad928SBarry Smith @*/
PCView(PC pc,PetscViewer viewer)1792d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView(PC pc, PetscViewer viewer)
1793d71ae5a4SJacob Faibussowitsch {
179419fd82e9SBarry Smith   PCType            cstr;
17956cd81132SPierre Jolivet   PetscViewerFormat format;
17969f196a02SMartin Diehl   PetscBool         isascii, isstring, isbinary, isdraw, pop = PETSC_FALSE;
1797e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1798536b137fSBarry Smith   PetscBool issaws;
17990acecf5bSBarry Smith #endif
18004b9ad928SBarry Smith 
18014b9ad928SBarry Smith   PetscFunctionBegin;
18020700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
180348a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
18040700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1805c9780b6fSBarry Smith   PetscCheckSameComm(pc, 1, viewer, 2);
18064b9ad928SBarry Smith 
18079f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
18089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
18099566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
18109566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1811e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
18129566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
18130acecf5bSBarry Smith #endif
1814219b1045SBarry Smith 
18159f196a02SMartin Diehl   if (isascii) {
18169566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer));
181748a46eb9SPierre Jolivet     if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, "  PC has not been set up so information may be incomplete\n"));
18189566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1819dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18209566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
1821834dbeb0SBarry Smith     if (pc->mat) {
18226cd81132SPierre Jolivet       PetscCall(PetscViewerGetFormat(viewer, &format));
18236cd81132SPierre Jolivet       if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
18249566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
18256cd81132SPierre Jolivet         pop = PETSC_TRUE;
18266cd81132SPierre Jolivet       }
18274b9ad928SBarry Smith       if (pc->pmat == pc->mat) {
1828ecf3d421SBarry Smith         PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix, which is also used to construct the preconditioner:\n"));
18299566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
18309566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
18319566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
18324b9ad928SBarry Smith       } else {
1833834dbeb0SBarry Smith         if (pc->pmat) {
1834ecf3d421SBarry Smith           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix, followed by the matrix used to construct the preconditioner:\n"));
18354b9ad928SBarry Smith         } else {
18369566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix:\n"));
18374b9ad928SBarry Smith         }
18389566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
18399566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
18409566063dSJacob Faibussowitsch         if (pc->pmat) PetscCall(MatView(pc->pmat, viewer));
18419566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
18424b9ad928SBarry Smith       }
18436cd81132SPierre Jolivet       if (pop) PetscCall(PetscViewerPopFormat(viewer));
18444b9ad928SBarry Smith     }
18454b9ad928SBarry Smith   } else if (isstring) {
18469566063dSJacob Faibussowitsch     PetscCall(PCGetType(pc, &cstr));
18479566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr));
1848dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18499566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18509566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18515b0b0462SBarry Smith   } else if (isbinary) {
185255849f57SBarry Smith     PetscInt    classid = PC_FILE_CLASSID;
185355849f57SBarry Smith     MPI_Comm    comm;
185455849f57SBarry Smith     PetscMPIInt rank;
185555849f57SBarry Smith     char        type[256];
185655849f57SBarry Smith 
18579566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
18589566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1859dd400576SPatrick Sanan     if (rank == 0) {
18609566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
18619566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256));
18629566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
186355849f57SBarry Smith     }
1864dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
1865219b1045SBarry Smith   } else if (isdraw) {
1866219b1045SBarry Smith     PetscDraw draw;
1867d9884438SBarry Smith     char      str[25];
186889fd9fafSBarry Smith     PetscReal x, y, bottom, h;
1869d9884438SBarry Smith     PetscInt  n;
1870219b1045SBarry Smith 
18719566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
18729566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
18731d840656SPeter Brune     if (pc->mat) {
18749566063dSJacob Faibussowitsch       PetscCall(MatGetSize(pc->mat, &n, NULL));
187563a3b9bcSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n));
18761d840656SPeter Brune     } else {
18779566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name));
18781d840656SPeter Brune     }
18799566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
188089fd9fafSBarry Smith     bottom = y - h;
18819566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
1882dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18839566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
1884e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1885536b137fSBarry Smith   } else if (issaws) {
1886d45a07a7SBarry Smith     PetscMPIInt rank;
1887d45a07a7SBarry Smith 
18889566063dSJacob Faibussowitsch     PetscCall(PetscObjectName((PetscObject)pc));
18899566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
189048a46eb9SPierre Jolivet     if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer));
18919566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18929566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18930acecf5bSBarry Smith #endif
18944b9ad928SBarry Smith   }
18953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18964b9ad928SBarry Smith }
18974b9ad928SBarry Smith 
18984b9ad928SBarry Smith /*@C
18990b4b7b1cSBarry Smith   PCRegister -  Adds a method (`PCType`) to the PETSc preconditioner package.
19001c84c290SBarry Smith 
1901cc4c1da9SBarry Smith   Not collective. No Fortran Support
19021c84c290SBarry Smith 
19031c84c290SBarry Smith   Input Parameters:
190420f4b53cSBarry Smith + sname    - name of a new user-defined solver
19050b4b7b1cSBarry Smith - function - routine to create the method context which will be stored in a `PC` when `PCSetType()` is called
19061c84c290SBarry Smith 
1907feefa0e1SJacob Faibussowitsch   Example Usage:
19081c84c290SBarry Smith .vb
1909bdf89e91SBarry Smith    PCRegister("my_solver", MySolverCreate);
19101c84c290SBarry Smith .ve
19111c84c290SBarry Smith 
19121c84c290SBarry Smith   Then, your solver can be chosen with the procedural interface via
1913b44f4de4SBarry Smith .vb
1914b44f4de4SBarry Smith   PCSetType(pc, "my_solver")
1915b44f4de4SBarry Smith .ve
19161c84c290SBarry Smith   or at runtime via the option
1917b44f4de4SBarry Smith .vb
1918b44f4de4SBarry Smith   -pc_type my_solver
1919b44f4de4SBarry Smith .ve
19204b9ad928SBarry Smith 
19214b9ad928SBarry Smith   Level: advanced
19221c84c290SBarry Smith 
192320f4b53cSBarry Smith   Note:
19240b4b7b1cSBarry Smith   A simpler alternative to using `PCRegister()` for an application specific preconditioner is to use a `PC` of `PCType` `PCSHELL` and
19250b4b7b1cSBarry Smith   provide your customizations with `PCShellSetContext()` and `PCShellSetApply()`
19260b4b7b1cSBarry Smith 
192720f4b53cSBarry Smith   `PCRegister()` may be called multiple times to add several user-defined preconditioners.
192820f4b53cSBarry Smith 
19290b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCType`, `PCRegisterAll()`, `PCSetType()`, `PCShellSetContext()`, `PCShellSetApply()`, `PCSHELL`
19304b9ad928SBarry Smith @*/
PCRegister(const char sname[],PetscErrorCode (* function)(PC))1931d71ae5a4SJacob Faibussowitsch PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC))
1932d71ae5a4SJacob Faibussowitsch {
19334b9ad928SBarry Smith   PetscFunctionBegin;
19349566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
19359566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCList, sname, function));
19363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19374b9ad928SBarry Smith }
19384b9ad928SBarry Smith 
MatMult_PC(Mat A,Vec X,Vec Y)1939d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y)
1940d71ae5a4SJacob Faibussowitsch {
1941186a3e20SStefano Zampini   PC pc;
1942186a3e20SStefano Zampini 
1943186a3e20SStefano Zampini   PetscFunctionBegin;
19449566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &pc));
19459566063dSJacob Faibussowitsch   PetscCall(PCApply(pc, X, Y));
19463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1947186a3e20SStefano Zampini }
1948186a3e20SStefano Zampini 
19495d83a8b1SBarry Smith /*@
19500b4b7b1cSBarry Smith   PCComputeOperator - Computes the explicit preconditioned operator as a matrix `Mat`.
19514b9ad928SBarry Smith 
1952c3339decSBarry Smith   Collective
19534b9ad928SBarry Smith 
1954d8d19677SJose E. Roman   Input Parameters:
19550b4b7b1cSBarry Smith + pc      - the `PC` preconditioner object
1956562efe2eSBarry Smith - mattype - the `MatType` to be used for the operator
19574b9ad928SBarry Smith 
19584b9ad928SBarry Smith   Output Parameter:
1959a5b23f4aSJose E. Roman . mat - the explicit preconditioned operator
19604b9ad928SBarry Smith 
196120f4b53cSBarry Smith   Level: advanced
196220f4b53cSBarry Smith 
1963f1580f4eSBarry Smith   Note:
1964186a3e20SStefano Zampini   This computation is done by applying the operators to columns of the identity matrix.
1965186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
1966562efe2eSBarry Smith   Currently, this routine uses a dense matrix format when `mattype` == `NULL`
19674b9ad928SBarry Smith 
19680b4b7b1cSBarry Smith   Developer Note:
19690b4b7b1cSBarry Smith   This should be called `PCCreateExplicitOperator()`
19700b4b7b1cSBarry Smith 
1971562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPComputeOperator()`, `MatType`
19724b9ad928SBarry Smith @*/
PCComputeOperator(PC pc,MatType mattype,Mat * mat)1973d71ae5a4SJacob Faibussowitsch PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat)
1974d71ae5a4SJacob Faibussowitsch {
1975186a3e20SStefano Zampini   PetscInt N, M, m, n;
1976186a3e20SStefano Zampini   Mat      A, Apc;
19774b9ad928SBarry Smith 
19784b9ad928SBarry Smith   PetscFunctionBegin;
19790700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
19804f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
19819566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, &A, NULL));
19829566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
19839566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
19849566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc));
198557d50842SBarry Smith   PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (PetscErrorCodeFn *)MatMult_PC));
19869566063dSJacob Faibussowitsch   PetscCall(MatComputeOperator(Apc, mattype, mat));
19879566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Apc));
19883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19894b9ad928SBarry Smith }
19904b9ad928SBarry Smith 
19916c237d78SBarry Smith /*@
19920b4b7b1cSBarry Smith   PCSetCoordinates - sets the coordinates of all the nodes (degrees of freedom in the vector) on the local process
19936c237d78SBarry Smith 
1994c3339decSBarry Smith   Collective
19956c237d78SBarry Smith 
19966c237d78SBarry Smith   Input Parameters:
19970b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
19986c237d78SBarry Smith . dim    - the dimension of the coordinates 1, 2, or 3
199914893cbeSStefano Zampini . nloc   - the blocked size of the coordinates array
200014893cbeSStefano Zampini - coords - the coordinates array
20016c237d78SBarry Smith 
20026c237d78SBarry Smith   Level: intermediate
20036c237d78SBarry Smith 
20040b4b7b1cSBarry Smith   Notes:
200520f4b53cSBarry Smith   `coords` is an array of the dim coordinates for the nodes on
200620f4b53cSBarry Smith   the local processor, of size `dim`*`nloc`.
20076e415bd2SNuno Nobre   If there are 108 equations (dofs) on a processor
20086e415bd2SNuno Nobre   for a 3d displacement finite element discretization of elasticity (so
200914893cbeSStefano Zampini   that there are nloc = 36 = 108/3 nodes) then the array must have 108
20106c237d78SBarry Smith   double precision values (ie, 3 * 36).  These x y z coordinates
20116c237d78SBarry Smith   should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
20126c237d78SBarry Smith   ... , N-1.z ].
20136c237d78SBarry Smith 
20140b4b7b1cSBarry Smith   The information provided here can be used by some preconditioners, such as `PCGAMG`, to produce a better preconditioner.
20150b4b7b1cSBarry Smith   See also  `MatSetNearNullSpace()`.
20160b4b7b1cSBarry Smith 
2017562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `MatSetNearNullSpace()`
20186c237d78SBarry Smith @*/
PCSetCoordinates(PC pc,PetscInt dim,PetscInt nloc,PetscReal coords[])2019d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[])
2020d71ae5a4SJacob Faibussowitsch {
20216c237d78SBarry Smith   PetscFunctionBegin;
202232954da3SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
202332954da3SStefano Zampini   PetscValidLogicalCollectiveInt(pc, dim, 2);
202422794d57SStefano Zampini   PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords));
20253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20266c237d78SBarry Smith }
2027fd2dd295SFande Kong 
2028fd2dd295SFande Kong /*@
2029fd2dd295SFande Kong   PCGetInterpolations - Gets interpolation matrices for all levels (except level 0)
2030fd2dd295SFande Kong 
2031c3339decSBarry Smith   Logically Collective
2032fd2dd295SFande Kong 
2033d8d19677SJose E. Roman   Input Parameter:
2034d8d19677SJose E. Roman . pc - the precondition context
2035fd2dd295SFande Kong 
2036d8d19677SJose E. Roman   Output Parameters:
2037d8d19677SJose E. Roman + num_levels     - the number of levels
2038562efe2eSBarry Smith - interpolations - the interpolation matrices (size of `num_levels`-1)
2039fd2dd295SFande Kong 
2040fd2dd295SFande Kong   Level: advanced
2041fd2dd295SFande Kong 
2042562efe2eSBarry Smith   Developer Note:
2043f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
2044fd2dd295SFande Kong 
2045562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()`
2046fd2dd295SFande Kong @*/
PCGetInterpolations(PC pc,PetscInt * num_levels,Mat * interpolations[])2047d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[])
2048d71ae5a4SJacob Faibussowitsch {
2049fd2dd295SFande Kong   PetscFunctionBegin;
2050fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20514f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20524f572ea9SToby Isaac   PetscAssertPointer(interpolations, 3);
2053cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations));
20543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2055fd2dd295SFande Kong }
2056fd2dd295SFande Kong 
2057fd2dd295SFande Kong /*@
2058fd2dd295SFande Kong   PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level)
2059fd2dd295SFande Kong 
2060c3339decSBarry Smith   Logically Collective
2061fd2dd295SFande Kong 
2062d8d19677SJose E. Roman   Input Parameter:
2063d8d19677SJose E. Roman . pc - the precondition context
2064fd2dd295SFande Kong 
2065d8d19677SJose E. Roman   Output Parameters:
2066d8d19677SJose E. Roman + num_levels      - the number of levels
2067562efe2eSBarry Smith - coarseOperators - the coarse operator matrices (size of `num_levels`-1)
2068fd2dd295SFande Kong 
2069fd2dd295SFande Kong   Level: advanced
2070fd2dd295SFande Kong 
2071562efe2eSBarry Smith   Developer Note:
2072f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
2073fd2dd295SFande Kong 
2074562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()`
2075fd2dd295SFande Kong @*/
PCGetCoarseOperators(PC pc,PetscInt * num_levels,Mat * coarseOperators[])2076d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
2077d71ae5a4SJacob Faibussowitsch {
2078fd2dd295SFande Kong   PetscFunctionBegin;
2079fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20804f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20814f572ea9SToby Isaac   PetscAssertPointer(coarseOperators, 3);
2082cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators));
20833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2084fd2dd295SFande Kong }
2085