xref: /petsc/src/ksp/pc/interface/precon.c (revision 4dbf25a8fa98e38799e7b47dcb2d8a9309975f41)
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 
14b4f8a55aSStefano Zampini PETSC_INTERN PetscErrorCode PCGetDefaultType_Private(PC pc, const char *type[])
15d71ae5a4SJacob Faibussowitsch {
162e70eadcSBarry Smith   PetscMPIInt size;
17b4f8a55aSStefano Zampini   PetscBool   hasopblock, hasopsolve, flg1, flg2, set, flg3, isnormal;
18566e8bf2SBarry Smith 
19566e8bf2SBarry Smith   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
21566e8bf2SBarry Smith   if (pc->pmat) {
22b4f8a55aSStefano Zampini     PetscCall(MatHasOperation(pc->pmat, MATOP_GET_DIAGONAL_BLOCK, &hasopblock));
23b4f8a55aSStefano Zampini     PetscCall(MatHasOperation(pc->pmat, MATOP_SOLVE, &hasopsolve));
24566e8bf2SBarry Smith     if (size == 1) {
259566063dSJacob Faibussowitsch       PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ICC, &flg1));
269566063dSJacob Faibussowitsch       PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ILU, &flg2));
279566063dSJacob Faibussowitsch       PetscCall(MatIsSymmetricKnown(pc->pmat, &set, &flg3));
28da74c943SJose E. Roman       PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isnormal, MATNORMAL, MATNORMALHERMITIAN, NULL));
29ba8a7149SBarry Smith       if (flg1 && (!flg2 || (set && flg3))) {
30566e8bf2SBarry Smith         *type = PCICC;
31566e8bf2SBarry Smith       } else if (flg2) {
32566e8bf2SBarry Smith         *type = PCILU;
33da74c943SJose E. Roman       } else if (isnormal) {
34da74c943SJose E. Roman         *type = PCNONE;
35b4f8a55aSStefano Zampini       } else if (hasopblock) { /* likely is a parallel matrix run on one processor */
36566e8bf2SBarry Smith         *type = PCBJACOBI;
37b4f8a55aSStefano Zampini       } else if (hasopsolve) {
38b4f8a55aSStefano Zampini         *type = PCMAT;
39566e8bf2SBarry Smith       } else {
40566e8bf2SBarry Smith         *type = PCNONE;
41566e8bf2SBarry Smith       }
42566e8bf2SBarry Smith     } else {
43b4f8a55aSStefano Zampini       if (hasopblock) {
44566e8bf2SBarry Smith         *type = PCBJACOBI;
45b4f8a55aSStefano Zampini       } else if (hasopsolve) {
46b4f8a55aSStefano Zampini         *type = PCMAT;
47566e8bf2SBarry Smith       } else {
48566e8bf2SBarry Smith         *type = PCNONE;
49566e8bf2SBarry Smith       }
50566e8bf2SBarry Smith     }
51b4f8a55aSStefano Zampini   } else *type = NULL;
523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53566e8bf2SBarry Smith }
544b9ad928SBarry Smith 
55fe57bb1aSStefano Zampini /* do not log solves, setup and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */
56fe57bb1aSStefano Zampini PETSC_EXTERN PetscLogEvent KSP_Solve, KSP_SetUp;
57fe57bb1aSStefano Zampini 
58fe57bb1aSStefano Zampini static PetscErrorCode PCLogEventsDeactivatePush(void)
59fe57bb1aSStefano Zampini {
60fe57bb1aSStefano Zampini   PetscFunctionBegin;
61fe57bb1aSStefano Zampini   PetscCall(KSPInitializePackage());
62fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePush(KSP_Solve));
63fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePush(KSP_SetUp));
64fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePush(PC_Apply));
65fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePush(PC_SetUp));
66fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePush(PC_SetUpOnBlocks));
67fe57bb1aSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
68fe57bb1aSStefano Zampini }
69fe57bb1aSStefano Zampini 
70fe57bb1aSStefano Zampini static PetscErrorCode PCLogEventsDeactivatePop(void)
71fe57bb1aSStefano Zampini {
72fe57bb1aSStefano Zampini   PetscFunctionBegin;
73fe57bb1aSStefano Zampini   PetscCall(KSPInitializePackage());
74fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePop(KSP_Solve));
75fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePop(KSP_SetUp));
76fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePop(PC_Apply));
77fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePop(PC_SetUp));
78fe57bb1aSStefano Zampini   PetscCall(PetscLogEventDeactivatePop(PC_SetUpOnBlocks));
79fe57bb1aSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
80fe57bb1aSStefano Zampini }
81fe57bb1aSStefano Zampini 
8288584be7SBarry Smith /*@
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 @*/
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 
1090ce6c5a2SBarry Smith   pc->setupcalled = 0;
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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
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;
4764b9ad928SBarry Smith   pc->setupcalled          = 0;
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 @*/
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 
535*4dbf25a8SPierre 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));
564*4dbf25a8SPierre 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));
568*4dbf25a8SPierre Jolivet   } else if (transpose && pc->ops->matapplytranspose) {
569*4dbf25a8SPierre Jolivet     PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0));
570*4dbf25a8SPierre Jolivet     PetscUseTypeMethod(pc, matapplytranspose, X, Y);
571*4dbf25a8SPierre 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));
577*4dbf25a8SPierre Jolivet       if (!transpose) PetscCall(PCApply(pc, cx, cy));
578*4dbf25a8SPierre 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 /*@
587*4dbf25a8SPierre Jolivet   PCMatApply - Applies the preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApply()`, `Y` and `X` must be different matrices.
588*4dbf25a8SPierre Jolivet 
589*4dbf25a8SPierre Jolivet   Collective
590*4dbf25a8SPierre Jolivet 
591*4dbf25a8SPierre Jolivet   Input Parameters:
592*4dbf25a8SPierre Jolivet + pc - the `PC` preconditioner context
593*4dbf25a8SPierre Jolivet - X  - block of input vectors
594*4dbf25a8SPierre Jolivet 
595*4dbf25a8SPierre Jolivet   Output Parameter:
596*4dbf25a8SPierre Jolivet . Y - block of output vectors
597*4dbf25a8SPierre Jolivet 
598*4dbf25a8SPierre Jolivet   Level: developer
599*4dbf25a8SPierre Jolivet 
600*4dbf25a8SPierre Jolivet .seealso: [](ch_ksp), `PC`, `PCApply()`, `KSPMatSolve()`
601*4dbf25a8SPierre Jolivet @*/
602*4dbf25a8SPierre Jolivet PetscErrorCode PCMatApply(PC pc, Mat X, Mat Y)
603*4dbf25a8SPierre Jolivet {
604*4dbf25a8SPierre Jolivet   PetscFunctionBegin;
605*4dbf25a8SPierre Jolivet   PetscCall(PCMatApplyTranspose_Private(pc, X, Y, PETSC_FALSE));
606*4dbf25a8SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
607*4dbf25a8SPierre Jolivet }
608*4dbf25a8SPierre Jolivet 
609*4dbf25a8SPierre Jolivet /*@
610*4dbf25a8SPierre Jolivet   PCMatApplyTranspose - Applies the transpose of preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApplyTranspose()`, `Y` and `X` must be different matrices.
611*4dbf25a8SPierre Jolivet 
612*4dbf25a8SPierre Jolivet   Collective
613*4dbf25a8SPierre Jolivet 
614*4dbf25a8SPierre Jolivet   Input Parameters:
615*4dbf25a8SPierre Jolivet + pc - the `PC` preconditioner context
616*4dbf25a8SPierre Jolivet - X  - block of input vectors
617*4dbf25a8SPierre Jolivet 
618*4dbf25a8SPierre Jolivet   Output Parameter:
619*4dbf25a8SPierre Jolivet . Y - block of output vectors
620*4dbf25a8SPierre Jolivet 
621*4dbf25a8SPierre Jolivet   Level: developer
622*4dbf25a8SPierre Jolivet 
623*4dbf25a8SPierre Jolivet .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`, `KSPMatSolveTranspose()`
624*4dbf25a8SPierre Jolivet @*/
625*4dbf25a8SPierre Jolivet PetscErrorCode PCMatApplyTranspose(PC pc, Mat X, Mat Y)
626*4dbf25a8SPierre Jolivet {
627*4dbf25a8SPierre Jolivet   PetscFunctionBegin;
628*4dbf25a8SPierre Jolivet   PetscCall(PCMatApplyTranspose_Private(pc, X, Y, PETSC_TRUE));
629*4dbf25a8SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
630*4dbf25a8SPierre Jolivet }
631*4dbf25a8SPierre Jolivet 
632*4dbf25a8SPierre 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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
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 @*/
10169a6c2652SBarry Smith PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason)
1017d71ae5a4SJacob Faibussowitsch {
10181b2b9847SBarry Smith   PetscFunctionBegin;
10196479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
10201b2b9847SBarry Smith   if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled;
10211b2b9847SBarry Smith   else *reason = pc->failedreason;
10223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10231b2b9847SBarry Smith }
1024422a814eSBarry Smith 
10256479e835SStefano Zampini /*@
10266479e835SStefano Zampini   PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC`
10276479e835SStefano Zampini 
10286479e835SStefano Zampini   Collective
10296479e835SStefano Zampini 
10306479e835SStefano Zampini   Input Parameter:
10310b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
10326479e835SStefano Zampini 
10336479e835SStefano Zampini   Level: advanced
10346479e835SStefano Zampini 
10356479e835SStefano Zampini   Note:
1036562efe2eSBarry Smith   Different MPI processes may have different reasons or no reason, see `PCGetFailedReason()`. This routine
10376479e835SStefano Zampini   makes them have a common value (failure if any MPI process had a failure).
10386479e835SStefano Zampini 
10390b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`, `PCFailedReason`
10406479e835SStefano Zampini @*/
10416479e835SStefano Zampini PetscErrorCode PCReduceFailedReason(PC pc)
10426479e835SStefano Zampini {
10436479e835SStefano Zampini   PetscInt buf;
10446479e835SStefano Zampini 
10456479e835SStefano Zampini   PetscFunctionBegin;
10466479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
10476479e835SStefano Zampini   buf = (PetscInt)pc->failedreason;
1048462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
10496479e835SStefano Zampini   pc->failedreason = (PCFailedReason)buf;
10506479e835SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
10516479e835SStefano Zampini }
10526479e835SStefano Zampini 
10534b9ad928SBarry Smith /*
10544b9ad928SBarry Smith       a setupcall of 0 indicates never setup,
105523ee1639SBarry Smith                      1 indicates has been previously setup
1056422a814eSBarry Smith                     -1 indicates a PCSetUp() was attempted and failed
10574b9ad928SBarry Smith */
10584b9ad928SBarry Smith /*@
10590b4b7b1cSBarry Smith   PCSetUp - Prepares for the use of a preconditioner. Performs all the one-time operations needed before the preconditioner
10600b4b7b1cSBarry Smith   can be used with `PCApply()`
10614b9ad928SBarry Smith 
1062c3339decSBarry Smith   Collective
10634b9ad928SBarry Smith 
10644b9ad928SBarry Smith   Input Parameter:
10650b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
10664b9ad928SBarry Smith 
10674b9ad928SBarry Smith   Level: developer
10684b9ad928SBarry Smith 
10690b4b7b1cSBarry Smith   Notes:
10700b4b7b1cSBarry Smith   For example, for `PCLU` this will compute the factorization.
10710b4b7b1cSBarry Smith 
10720b4b7b1cSBarry Smith   This is called automatically by `KSPSetUp()` or `PCApply()` so rarely needs to be called directly.
10730b4b7b1cSBarry Smith 
10740b4b7b1cSBarry Smith   For nested preconditioners, such as `PCFIELDSPLIT` or `PCBJACOBI` this may not finish the construction of the preconditioner
10750b4b7b1cSBarry Smith   on the inner levels, the routine `PCSetUpOnBlocks()` may compute more of the preconditioner in those situations.
10760b4b7b1cSBarry Smith 
10770b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `KSPSetUp()`, `PCSetUpOnBlocks()`
10784b9ad928SBarry Smith @*/
1079d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp(PC pc)
1080d71ae5a4SJacob Faibussowitsch {
1081566e8bf2SBarry Smith   const char      *def;
1082fc9b51d3SBarry Smith   PetscObjectState matstate, matnonzerostate;
10834b9ad928SBarry Smith 
10844b9ad928SBarry Smith   PetscFunctionBegin;
10850700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
108628b400f6SJacob Faibussowitsch   PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first");
10876ce5e81cSLisandro Dalcin 
108823ee1639SBarry Smith   if (pc->setupcalled && pc->reusepreconditioner) {
10899566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n"));
10903ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10914b9ad928SBarry Smith   }
10924b9ad928SBarry Smith 
10939566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate));
10949566063dSJacob Faibussowitsch   PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate));
109523ee1639SBarry Smith   if (!pc->setupcalled) {
10965e62d202SMark Adams     //PetscCall(PetscInfo(pc, "Setting up PC for first time\n"));
109723ee1639SBarry Smith     pc->flag = DIFFERENT_NONZERO_PATTERN;
10989df67409SStefano Zampini   } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS);
10999df67409SStefano Zampini   else {
11009df67409SStefano Zampini     if (matnonzerostate != pc->matnonzerostate) {
11019566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n"));
110223ee1639SBarry Smith       pc->flag = DIFFERENT_NONZERO_PATTERN;
110323ee1639SBarry Smith     } else {
11045e62d202SMark Adams       //PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n"));
110523ee1639SBarry Smith       pc->flag = SAME_NONZERO_PATTERN;
110623ee1639SBarry Smith     }
110723ee1639SBarry Smith   }
110823ee1639SBarry Smith   pc->matstate        = matstate;
1109fc9b51d3SBarry Smith   pc->matnonzerostate = matnonzerostate;
111023ee1639SBarry Smith 
11117adad957SLisandro Dalcin   if (!((PetscObject)pc)->type_name) {
11129566063dSJacob Faibussowitsch     PetscCall(PCGetDefaultType_Private(pc, &def));
11139566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, def));
1114566e8bf2SBarry Smith   }
11154b9ad928SBarry Smith 
11169566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
11179566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure));
11189566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0));
11194b9ad928SBarry Smith   if (pc->ops->setup) {
1120fe57bb1aSStefano Zampini     PetscCall(PCLogEventsDeactivatePush());
1121dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, setup);
1122fe57bb1aSStefano Zampini     PetscCall(PCLogEventsDeactivatePop());
11234b9ad928SBarry Smith   }
11249566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0));
1125422a814eSBarry Smith   if (!pc->setupcalled) pc->setupcalled = 1;
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 @*/
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
11714b9ad928SBarry Smith . func - routine for modifying the submatrices
117204c3f3b8SBarry Smith - ctx  - optional user-defined context (may be `NULL`)
11734b9ad928SBarry Smith 
117420f4b53cSBarry Smith   Calling sequence of `func`:
11750b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
117604c3f3b8SBarry Smith . nsub   - number of index sets
117720f4b53cSBarry Smith . row    - an array of index sets that contain the global row numbers
11784b9ad928SBarry Smith          that comprise each local submatrix
11794b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
11804b9ad928SBarry Smith          that comprise each local submatrix
11814b9ad928SBarry Smith . submat - array of local submatrices
11824b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
118304c3f3b8SBarry Smith          user-defined func routine (may be `NULL`)
11844b9ad928SBarry Smith 
118520f4b53cSBarry Smith   Level: advanced
118620f4b53cSBarry Smith 
11874b9ad928SBarry Smith   Notes:
118804c3f3b8SBarry Smith   The basic submatrices are extracted from the preconditioner matrix as
118904c3f3b8SBarry Smith   usual; the user can then alter these (for example, to set different boundary
119004c3f3b8SBarry Smith   conditions for each submatrix) before they are used for the local solves.
119104c3f3b8SBarry Smith 
1192f1580f4eSBarry Smith   `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and
1193f1580f4eSBarry Smith   `KSPSolve()`.
11944b9ad928SBarry Smith 
1195f1580f4eSBarry Smith   A routine set by `PCSetModifySubMatrices()` is currently called within
1196f1580f4eSBarry Smith   the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`)
11974b9ad928SBarry Smith   preconditioners.  All other preconditioners ignore this routine.
11984b9ad928SBarry Smith 
1199562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()`
12004b9ad928SBarry Smith @*/
120104c3f3b8SBarry Smith PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx), void *ctx)
1202d71ae5a4SJacob Faibussowitsch {
12034b9ad928SBarry Smith   PetscFunctionBegin;
12040700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12054b9ad928SBarry Smith   pc->modifysubmatrices  = func;
12064b9ad928SBarry Smith   pc->modifysubmatricesP = ctx;
12073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12084b9ad928SBarry Smith }
12094b9ad928SBarry Smith 
12104b9ad928SBarry Smith /*@C
12114b9ad928SBarry Smith   PCModifySubMatrices - Calls an optional user-defined routine within
1212f1580f4eSBarry Smith   certain preconditioners if one has been set with `PCSetModifySubMatrices()`.
12134b9ad928SBarry Smith 
1214c3339decSBarry Smith   Collective
12154b9ad928SBarry Smith 
12164b9ad928SBarry Smith   Input Parameters:
12170b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
12184b9ad928SBarry Smith . nsub   - the number of local submatrices
12194b9ad928SBarry Smith . row    - an array of index sets that contain the global row numbers
12204b9ad928SBarry Smith          that comprise each local submatrix
12214b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
12224b9ad928SBarry Smith          that comprise each local submatrix
12234b9ad928SBarry Smith . submat - array of local submatrices
12244b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
1225562efe2eSBarry Smith          user-defined routine (may be `NULL`)
12264b9ad928SBarry Smith 
12274b9ad928SBarry Smith   Output Parameter:
12284b9ad928SBarry Smith . submat - array of local submatrices (the entries of which may
12294b9ad928SBarry Smith             have been modified)
12304b9ad928SBarry Smith 
123120f4b53cSBarry Smith   Level: developer
123220f4b53cSBarry Smith 
123304c3f3b8SBarry Smith   Note:
12344b9ad928SBarry Smith   The user should NOT generally call this routine, as it will
123504c3f3b8SBarry Smith   automatically be called within certain preconditioners.
12364b9ad928SBarry Smith 
1237562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetModifySubMatrices()`
12384b9ad928SBarry Smith @*/
1239d71ae5a4SJacob Faibussowitsch PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx)
1240d71ae5a4SJacob Faibussowitsch {
12414b9ad928SBarry Smith   PetscFunctionBegin;
12420700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12433ba16761SJacob Faibussowitsch   if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS);
12449566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0));
12459566063dSJacob Faibussowitsch   PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx));
12469566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0));
12473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12484b9ad928SBarry Smith }
12494b9ad928SBarry Smith 
12504b9ad928SBarry Smith /*@
12514b9ad928SBarry Smith   PCSetOperators - Sets the matrix associated with the linear system and
12524b9ad928SBarry Smith   a (possibly) different one associated with the preconditioner.
12534b9ad928SBarry Smith 
1254c3339decSBarry Smith   Logically Collective
12554b9ad928SBarry Smith 
12564b9ad928SBarry Smith   Input Parameters:
12570b4b7b1cSBarry Smith + pc   - the `PC` preconditioner context
1258e5d3d808SBarry Smith . Amat - the matrix that defines the linear system
1259d1e9a80fSBarry Smith - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
12604b9ad928SBarry Smith 
126120f4b53cSBarry Smith   Level: intermediate
1262189c0b8aSBarry Smith 
126320f4b53cSBarry Smith   Notes:
126420f4b53cSBarry Smith   Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used.
126520f4b53cSBarry Smith 
126620f4b53cSBarry Smith   If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then
1267f1580f4eSBarry Smith   first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()`
1268f1580f4eSBarry Smith   on it and then pass it back in in your call to `KSPSetOperators()`.
1269189c0b8aSBarry Smith 
12704b9ad928SBarry Smith   More Notes about Repeated Solution of Linear Systems:
127120f4b53cSBarry Smith   PETSc does NOT reset the matrix entries of either `Amat` or `Pmat`
12724b9ad928SBarry Smith   to zero after a linear solve; the user is completely responsible for
1273f1580f4eSBarry Smith   matrix assembly.  See the routine `MatZeroEntries()` if desiring to
12744b9ad928SBarry Smith   zero all elements of a matrix.
12754b9ad928SBarry Smith 
1276562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`
12774b9ad928SBarry Smith  @*/
1278d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat)
1279d71ae5a4SJacob Faibussowitsch {
12803b3f6333SBarry Smith   PetscInt m1, n1, m2, n2;
12814b9ad928SBarry Smith 
12824b9ad928SBarry Smith   PetscFunctionBegin;
12830700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12840700a824SBarry Smith   if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
12850700a824SBarry Smith   if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3);
12863fc8bf9cSBarry Smith   if (Amat) PetscCheckSameComm(pc, 1, Amat, 2);
12873fc8bf9cSBarry Smith   if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3);
128831641f1aSBarry Smith   if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) {
12899566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Amat, &m1, &n1));
12909566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->mat, &m2, &n2));
129163a3b9bcSJacob 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);
12929566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Pmat, &m1, &n1));
12939566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2));
129463a3b9bcSJacob 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);
12953b3f6333SBarry Smith   }
12964b9ad928SBarry Smith 
129723ee1639SBarry Smith   if (Pmat != pc->pmat) {
129823ee1639SBarry Smith     /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
129923ee1639SBarry Smith     pc->matnonzerostate = -1;
130023ee1639SBarry Smith     pc->matstate        = -1;
130123ee1639SBarry Smith   }
130223ee1639SBarry Smith 
1303906ed7ccSBarry Smith   /* reference first in case the matrices are the same */
13049566063dSJacob Faibussowitsch   if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat));
13059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->mat));
13069566063dSJacob Faibussowitsch   if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat));
13079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->pmat));
13084b9ad928SBarry Smith   pc->mat  = Amat;
13094b9ad928SBarry Smith   pc->pmat = Pmat;
13103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1311e56f5c9eSBarry Smith }
1312e56f5c9eSBarry Smith 
131323ee1639SBarry Smith /*@
13140b4b7b1cSBarry Smith   PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner `PC` has changed.
131523ee1639SBarry Smith 
1316c3339decSBarry Smith   Logically Collective
131723ee1639SBarry Smith 
131823ee1639SBarry Smith   Input Parameters:
13190b4b7b1cSBarry Smith + pc   - the `PC` preconditioner context
1320f1580f4eSBarry Smith - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
132123ee1639SBarry Smith 
13222b26979fSBarry Smith   Level: intermediate
13232b26979fSBarry Smith 
1324f1580f4eSBarry Smith   Note:
1325f1580f4eSBarry Smith   Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option
1326f1580f4eSBarry Smith   prevents this.
1327f1580f4eSBarry Smith 
1328562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()`
132923ee1639SBarry Smith  @*/
1330d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag)
1331d71ae5a4SJacob Faibussowitsch {
133223ee1639SBarry Smith   PetscFunctionBegin;
133323ee1639SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1334f9177522SStefano Zampini   PetscValidLogicalCollectiveBool(pc, flag, 2);
133523ee1639SBarry Smith   pc->reusepreconditioner = flag;
13363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13374b9ad928SBarry Smith }
13384b9ad928SBarry Smith 
1339c60c7ad4SBarry Smith /*@
1340f1580f4eSBarry Smith   PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed.
1341c60c7ad4SBarry Smith 
1342bba28a21SBarry Smith   Not Collective
1343c60c7ad4SBarry Smith 
1344c60c7ad4SBarry Smith   Input Parameter:
13450b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
1346c60c7ad4SBarry Smith 
1347c60c7ad4SBarry Smith   Output Parameter:
1348f1580f4eSBarry Smith . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
1349c60c7ad4SBarry Smith 
1350d0418729SSatish Balay   Level: intermediate
1351d0418729SSatish Balay 
1352562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()`
1353c60c7ad4SBarry Smith  @*/
1354d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag)
1355d71ae5a4SJacob Faibussowitsch {
1356c60c7ad4SBarry Smith   PetscFunctionBegin;
1357c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
13584f572ea9SToby Isaac   PetscAssertPointer(flag, 2);
1359c60c7ad4SBarry Smith   *flag = pc->reusepreconditioner;
13603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1361c60c7ad4SBarry Smith }
1362c60c7ad4SBarry Smith 
1363487a658cSBarry Smith /*@
13644b9ad928SBarry Smith   PCGetOperators - Gets the matrix associated with the linear system and
13650b4b7b1cSBarry Smith   possibly a different one which is used to construct the preconditioner.
13664b9ad928SBarry Smith 
1367562efe2eSBarry Smith   Not Collective, though parallel `Mat`s are returned if `pc` is parallel
13684b9ad928SBarry Smith 
13694b9ad928SBarry Smith   Input Parameter:
13700b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
13714b9ad928SBarry Smith 
13724b9ad928SBarry Smith   Output Parameters:
1373e5d3d808SBarry Smith + Amat - the matrix defining the linear system
137423ee1639SBarry Smith - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
13754b9ad928SBarry Smith 
13764b9ad928SBarry Smith   Level: intermediate
13774b9ad928SBarry Smith 
1378f1580f4eSBarry Smith   Note:
137995452b02SPatrick Sanan   Does not increase the reference count of the matrices, so you should not destroy them
1380298cc208SBarry Smith 
1381f1580f4eSBarry Smith   Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators
1382f1580f4eSBarry Smith   are created in `PC` and returned to the user. In this case, if both operators
138373950996SBarry Smith   mat and pmat are requested, two DIFFERENT operators will be returned. If
138473950996SBarry Smith   only one is requested both operators in the PC will be the same (i.e. as
1385f1580f4eSBarry Smith   if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats).
138673950996SBarry Smith   The user must set the sizes of the returned matrices and their type etc just
1387f1580f4eSBarry Smith   as if the user created them with `MatCreate()`. For example,
138873950996SBarry Smith 
1389f1580f4eSBarry Smith .vb
1390f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1391f1580f4eSBarry Smith            set size, type, etc of Amat
139273950996SBarry Smith 
1393f1580f4eSBarry Smith          MatCreate(comm,&mat);
1394f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Amat);
1395f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)mat);
1396f1580f4eSBarry Smith            set size, type, etc of Amat
1397f1580f4eSBarry Smith .ve
139873950996SBarry Smith 
139973950996SBarry Smith   and
140073950996SBarry Smith 
1401f1580f4eSBarry Smith .vb
1402f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1403f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
140473950996SBarry Smith 
1405f1580f4eSBarry Smith          MatCreate(comm,&Amat);
1406f1580f4eSBarry Smith          MatCreate(comm,&Pmat);
1407f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1408f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Amat);
1409f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Pmat);
1410f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
1411f1580f4eSBarry Smith .ve
141273950996SBarry Smith 
1413f1580f4eSBarry Smith   The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy
1414b8d709abSRichard Tran Mills   of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely
1415f1580f4eSBarry Smith   managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look
1416f1580f4eSBarry Smith   at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to
1417f1580f4eSBarry Smith   the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP
1418f1580f4eSBarry Smith   you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you).
1419f1580f4eSBarry Smith   Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when
142073950996SBarry Smith   it can be created for you?
142173950996SBarry Smith 
1422562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()`
14234b9ad928SBarry Smith @*/
1424d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat)
1425d71ae5a4SJacob Faibussowitsch {
14264b9ad928SBarry Smith   PetscFunctionBegin;
14270700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1428e5d3d808SBarry Smith   if (Amat) {
142973950996SBarry Smith     if (!pc->mat) {
14309fca8c71SStefano Zampini       if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */
14319a4708feSJed Brown         pc->mat = pc->pmat;
14329566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->mat));
1433e5d3d808SBarry Smith       } else { /* both Amat and Pmat are empty */
14349566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat));
1435e5d3d808SBarry Smith         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
143673950996SBarry Smith           pc->pmat = pc->mat;
14379566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->pmat));
143873950996SBarry Smith         }
143973950996SBarry Smith       }
14409a4708feSJed Brown     }
1441e5d3d808SBarry Smith     *Amat = pc->mat;
144273950996SBarry Smith   }
1443e5d3d808SBarry Smith   if (Pmat) {
144473950996SBarry Smith     if (!pc->pmat) {
1445e5d3d808SBarry Smith       if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */
14469a4708feSJed Brown         pc->pmat = pc->mat;
14479566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->pmat));
14489a4708feSJed Brown       } else {
14499566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat));
1450e5d3d808SBarry Smith         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
145173950996SBarry Smith           pc->mat = pc->pmat;
14529566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->mat));
145373950996SBarry Smith         }
145473950996SBarry Smith       }
14559a4708feSJed Brown     }
1456e5d3d808SBarry Smith     *Pmat = pc->pmat;
145773950996SBarry Smith   }
14583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14594b9ad928SBarry Smith }
14604b9ad928SBarry Smith 
14615d83a8b1SBarry Smith /*@
1462906ed7ccSBarry Smith   PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1463f1580f4eSBarry Smith   possibly a different one associated with the preconditioner have been set in the `PC`.
1464906ed7ccSBarry Smith 
146520f4b53cSBarry Smith   Not Collective, though the results on all processes should be the same
1466906ed7ccSBarry Smith 
1467906ed7ccSBarry Smith   Input Parameter:
14680b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
1469906ed7ccSBarry Smith 
1470906ed7ccSBarry Smith   Output Parameters:
1471906ed7ccSBarry Smith + mat  - the matrix associated with the linear system was set
1472906ed7ccSBarry Smith - pmat - matrix associated with the preconditioner was set, usually the same
1473906ed7ccSBarry Smith 
1474906ed7ccSBarry Smith   Level: intermediate
1475906ed7ccSBarry Smith 
1476562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`
1477906ed7ccSBarry Smith @*/
1478d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat)
1479d71ae5a4SJacob Faibussowitsch {
1480906ed7ccSBarry Smith   PetscFunctionBegin;
14810700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1482906ed7ccSBarry Smith   if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1483906ed7ccSBarry Smith   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
14843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1485906ed7ccSBarry Smith }
1486906ed7ccSBarry Smith 
1487f39d8e23SSatish Balay /*@
1488a4fd02acSBarry Smith   PCFactorGetMatrix - Gets the factored matrix from the
1489f1580f4eSBarry Smith   preconditioner context.  This routine is valid only for the `PCLU`,
1490f1580f4eSBarry Smith   `PCILU`, `PCCHOLESKY`, and `PCICC` methods.
14914b9ad928SBarry Smith 
1492562efe2eSBarry Smith   Not Collective though `mat` is parallel if `pc` is parallel
14934b9ad928SBarry Smith 
1494f1580f4eSBarry Smith   Input Parameter:
14950b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
14964b9ad928SBarry Smith 
1497feefa0e1SJacob Faibussowitsch   Output Parameters:
14984b9ad928SBarry Smith . mat - the factored matrix
14994b9ad928SBarry Smith 
15004b9ad928SBarry Smith   Level: advanced
15014b9ad928SBarry Smith 
1502f1580f4eSBarry Smith   Note:
1503562efe2eSBarry Smith   Does not increase the reference count for `mat` so DO NOT destroy it
15049405f653SBarry Smith 
1505562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC`
15064b9ad928SBarry Smith @*/
1507d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat)
1508d71ae5a4SJacob Faibussowitsch {
15094b9ad928SBarry Smith   PetscFunctionBegin;
15100700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15114f572ea9SToby Isaac   PetscAssertPointer(mat, 2);
15127def7855SStefano Zampini   PetscCall(PCFactorSetUpMatSolverType(pc));
1513dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, getfactoredmatrix, mat);
15143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15154b9ad928SBarry Smith }
15164b9ad928SBarry Smith 
1517cc4c1da9SBarry Smith /*@
15184b9ad928SBarry Smith   PCSetOptionsPrefix - Sets the prefix used for searching for all
1519f1580f4eSBarry Smith   `PC` options in the database.
15204b9ad928SBarry Smith 
1521c3339decSBarry Smith   Logically Collective
15224b9ad928SBarry Smith 
15234b9ad928SBarry Smith   Input Parameters:
15240b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
1525f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
15264b9ad928SBarry Smith 
1527f1580f4eSBarry Smith   Note:
15284b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
15294b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
15304b9ad928SBarry Smith   hyphen.
15314b9ad928SBarry Smith 
15324b9ad928SBarry Smith   Level: advanced
15334b9ad928SBarry Smith 
1534562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()`
15354b9ad928SBarry Smith @*/
1536d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[])
1537d71ae5a4SJacob Faibussowitsch {
15384b9ad928SBarry Smith   PetscFunctionBegin;
15390700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15409566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix));
15413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15424b9ad928SBarry Smith }
15434b9ad928SBarry Smith 
1544cc4c1da9SBarry Smith /*@
15454b9ad928SBarry Smith   PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1546f1580f4eSBarry Smith   `PC` options in the database.
15474b9ad928SBarry Smith 
1548c3339decSBarry Smith   Logically Collective
15494b9ad928SBarry Smith 
15504b9ad928SBarry Smith   Input Parameters:
15510b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
1552f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
15534b9ad928SBarry Smith 
1554f1580f4eSBarry Smith   Note:
15554b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
15564b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
15574b9ad928SBarry Smith   hyphen.
15584b9ad928SBarry Smith 
15594b9ad928SBarry Smith   Level: advanced
15604b9ad928SBarry Smith 
1561562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()`
15624b9ad928SBarry Smith @*/
1563d71ae5a4SJacob Faibussowitsch PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[])
1564d71ae5a4SJacob Faibussowitsch {
15654b9ad928SBarry Smith   PetscFunctionBegin;
15660700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15679566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix));
15683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15694b9ad928SBarry Smith }
15704b9ad928SBarry Smith 
1571cc4c1da9SBarry Smith /*@
15724b9ad928SBarry Smith   PCGetOptionsPrefix - Gets the prefix used for searching for all
15734b9ad928SBarry Smith   PC options in the database.
15744b9ad928SBarry Smith 
15754b9ad928SBarry Smith   Not Collective
15764b9ad928SBarry Smith 
1577f1580f4eSBarry Smith   Input Parameter:
15780b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
15794b9ad928SBarry Smith 
1580f1580f4eSBarry Smith   Output Parameter:
15814b9ad928SBarry Smith . prefix - pointer to the prefix string used, is returned
15824b9ad928SBarry Smith 
15834b9ad928SBarry Smith   Level: advanced
15844b9ad928SBarry Smith 
1585562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()`
15864b9ad928SBarry Smith @*/
1587d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[])
1588d71ae5a4SJacob Faibussowitsch {
15894b9ad928SBarry Smith   PetscFunctionBegin;
15900700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15914f572ea9SToby Isaac   PetscAssertPointer(prefix, 2);
15929566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix));
15933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15944b9ad928SBarry Smith }
15954b9ad928SBarry Smith 
15968066bbecSBarry Smith /*
1597dd8e379bSPierre Jolivet    Indicates the right-hand side will be changed by KSPSolve(), this occurs for a few
15988066bbecSBarry Smith   preconditioners including BDDC and Eisentat that transform the equations before applying
15998066bbecSBarry Smith   the Krylov methods
16008066bbecSBarry Smith */
1601d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change)
1602d71ae5a4SJacob Faibussowitsch {
1603a06fd7f2SStefano Zampini   PetscFunctionBegin;
1604a06fd7f2SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16054f572ea9SToby Isaac   PetscAssertPointer(change, 2);
1606a06fd7f2SStefano Zampini   *change = PETSC_FALSE;
1607cac4c232SBarry Smith   PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change));
16083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1609a06fd7f2SStefano Zampini }
1610a06fd7f2SStefano Zampini 
16114b9ad928SBarry Smith /*@
16120b4b7b1cSBarry Smith   PCPreSolve - Optional pre-solve phase, intended for any preconditioner-specific actions that must be performed before
16130b4b7b1cSBarry Smith   the iterative solve itself. Used in conjunction with `PCPostSolve()`
16144b9ad928SBarry Smith 
1615c3339decSBarry Smith   Collective
16164b9ad928SBarry Smith 
16174b9ad928SBarry Smith   Input Parameters:
16180b4b7b1cSBarry Smith + pc  - the `PC` preconditioner context
16194b9ad928SBarry Smith - ksp - the Krylov subspace context
16204b9ad928SBarry Smith 
16214b9ad928SBarry Smith   Level: developer
16224b9ad928SBarry Smith 
16234b9ad928SBarry Smith   Notes:
1624f1580f4eSBarry Smith   `KSPSolve()` calls this directly, so is rarely called by the user.
16254b9ad928SBarry Smith 
16260b4b7b1cSBarry Smith   Certain preconditioners, such as the `PCType` of `PCEISENSTAT`, change the formulation of the linear system to be solved iteratively.
16270b4b7b1cSBarry Smith   This function performs that transformation. `PCPostSolve()` then transforms the system back to its original form after the solve.
16280b4b7b1cSBarry Smith   `PCPostSolve()` also transforms the resulting solution of the transformed system to the solution of the original problem.
16290b4b7b1cSBarry Smith 
163054c05997SPierre Jolivet   `KSPSetPreSolve()` and `KSPSetPostSolve()` provide an alternative way to provide such transformations.
16310b4b7b1cSBarry Smith 
16320b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCPostSolve()`, `KSP`, `PCSetPreSolve()`, `KSPSetPreSolve()`, `KSPSetPostSolve()`
16334b9ad928SBarry Smith @*/
1634d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPreSolve(PC pc, KSP ksp)
1635d71ae5a4SJacob Faibussowitsch {
16364b9ad928SBarry Smith   Vec x, rhs;
16374b9ad928SBarry Smith 
16384b9ad928SBarry Smith   PetscFunctionBegin;
16390700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16400700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1641d9a03883SBarry Smith   pc->presolvedone++;
16427827d75bSBarry Smith   PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice");
16439566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
16449566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
16454b9ad928SBarry Smith 
1646dbbe0bcdSBarry Smith   if (pc->ops->presolve) PetscUseTypeMethod(pc, presolve, ksp, rhs, x);
1647f4f49eeaSPierre Jolivet   else if (pc->presolve) PetscCall(pc->presolve(pc, ksp));
16483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16494b9ad928SBarry Smith }
16504b9ad928SBarry Smith 
1651f560b561SHong Zhang /*@C
1652f1580f4eSBarry Smith   PCSetPreSolve - Sets function used by `PCPreSolve()` which is intended for any
1653f560b561SHong Zhang   preconditioner-specific actions that must be performed before
1654f560b561SHong Zhang   the iterative solve itself.
1655f560b561SHong Zhang 
1656c3339decSBarry Smith   Logically Collective
1657f560b561SHong Zhang 
1658f560b561SHong Zhang   Input Parameters:
1659f560b561SHong Zhang + pc       - the preconditioner object
1660f560b561SHong Zhang - presolve - the function to call before the solve
1661f560b561SHong Zhang 
166220f4b53cSBarry Smith   Calling sequence of `presolve`:
166320f4b53cSBarry Smith + pc  - the `PC` context
166420f4b53cSBarry Smith - ksp - the `KSP` context
1665f560b561SHong Zhang 
1666f560b561SHong Zhang   Level: developer
1667f560b561SHong Zhang 
1668562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCPreSolve()`
1669f560b561SHong Zhang @*/
167004c3f3b8SBarry Smith PetscErrorCode PCSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC pc, KSP ksp))
1671d71ae5a4SJacob Faibussowitsch {
1672f560b561SHong Zhang   PetscFunctionBegin;
1673f560b561SHong Zhang   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1674f560b561SHong Zhang   pc->presolve = presolve;
16753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1676f560b561SHong Zhang }
1677f560b561SHong Zhang 
16784b9ad928SBarry Smith /*@
16794b9ad928SBarry Smith   PCPostSolve - Optional post-solve phase, intended for any
16804b9ad928SBarry Smith   preconditioner-specific actions that must be performed after
16814b9ad928SBarry Smith   the iterative solve itself.
16824b9ad928SBarry Smith 
1683c3339decSBarry Smith   Collective
16844b9ad928SBarry Smith 
16854b9ad928SBarry Smith   Input Parameters:
16860b4b7b1cSBarry Smith + pc  - the `PC` preconditioner context
16870b4b7b1cSBarry Smith - ksp - the `KSP` Krylov subspace context
16884b9ad928SBarry Smith 
1689feefa0e1SJacob Faibussowitsch   Example Usage:
16904b9ad928SBarry Smith .vb
16914b9ad928SBarry Smith     PCPreSolve(pc,ksp);
169223ce1328SBarry Smith     KSPSolve(ksp,b,x);
16934b9ad928SBarry Smith     PCPostSolve(pc,ksp);
16944b9ad928SBarry Smith .ve
16954b9ad928SBarry Smith 
1696562efe2eSBarry Smith   Level: developer
1697562efe2eSBarry Smith 
16984b9ad928SBarry Smith   Note:
1699f1580f4eSBarry Smith   `KSPSolve()` calls this routine directly, so it is rarely called by the user.
17004b9ad928SBarry Smith 
170154c05997SPierre Jolivet .seealso: [](ch_ksp), `PC`, `PCSetPreSolve()`, `KSPSetPostSolve()`, `KSPSetPreSolve()`, `PCPreSolve()`, `KSPSolve()`
17024b9ad928SBarry Smith @*/
1703d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPostSolve(PC pc, KSP ksp)
1704d71ae5a4SJacob Faibussowitsch {
17054b9ad928SBarry Smith   Vec x, rhs;
17064b9ad928SBarry Smith 
17074b9ad928SBarry Smith   PetscFunctionBegin;
17080700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
17090700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1710d9a03883SBarry Smith   pc->presolvedone--;
17119566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
17129566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
1713dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, postsolve, ksp, rhs, x);
17143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17154b9ad928SBarry Smith }
17164b9ad928SBarry Smith 
1717ffeef943SBarry Smith /*@
1718f1580f4eSBarry Smith   PCLoad - Loads a `PC` that has been stored in binary  with `PCView()`.
171955849f57SBarry Smith 
1720c3339decSBarry Smith   Collective
172155849f57SBarry Smith 
172255849f57SBarry Smith   Input Parameters:
1723f1580f4eSBarry Smith + newdm  - the newly loaded `PC`, this needs to have been created with `PCCreate()` or
1724f1580f4eSBarry Smith            some related function before a call to `PCLoad()`.
17250b4b7b1cSBarry Smith - viewer - binary file viewer `PETSCVIEWERBINARY`, obtained from `PetscViewerBinaryOpen()`
172655849f57SBarry Smith 
172755849f57SBarry Smith   Level: intermediate
172855849f57SBarry Smith 
1729f1580f4eSBarry Smith   Note:
1730562efe2eSBarry Smith   The type is determined by the data in the file, any `PCType` set into the `PC` before this call is ignored.
173155849f57SBarry Smith 
17320b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()`, `PETSCVIEWERBINARY`
173355849f57SBarry Smith @*/
1734d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLoad(PC newdm, PetscViewer viewer)
1735d71ae5a4SJacob Faibussowitsch {
173655849f57SBarry Smith   PetscBool isbinary;
1737060da220SMatthew G. Knepley   PetscInt  classid;
173855849f57SBarry Smith   char      type[256];
173955849f57SBarry Smith 
174055849f57SBarry Smith   PetscFunctionBegin;
174155849f57SBarry Smith   PetscValidHeaderSpecific(newdm, PC_CLASSID, 1);
174255849f57SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
17439566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
174428b400f6SJacob Faibussowitsch   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
174555849f57SBarry Smith 
17469566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
17477827d75bSBarry Smith   PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file");
17489566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
17499566063dSJacob Faibussowitsch   PetscCall(PCSetType(newdm, type));
1750dbbe0bcdSBarry Smith   PetscTryTypeMethod(newdm, load, viewer);
17513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
175255849f57SBarry Smith }
175355849f57SBarry Smith 
17549804daf3SBarry Smith #include <petscdraw.h>
1755e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1756e04113cfSBarry Smith   #include <petscviewersaws.h>
17570acecf5bSBarry Smith #endif
1758fe2efc57SMark 
1759ffeef943SBarry Smith /*@
17600b4b7b1cSBarry Smith   PCViewFromOptions - View (print or provide information about) the `PC`, based on options in the options database
1761fe2efc57SMark 
1762c3339decSBarry Smith   Collective
1763fe2efc57SMark 
1764fe2efc57SMark   Input Parameters:
176520f4b53cSBarry Smith + A    - the `PC` context
1766f1580f4eSBarry Smith . obj  - Optional object that provides the options prefix
17670b4b7b1cSBarry Smith - name - command line option name
1768fe2efc57SMark 
17690b4b7b1cSBarry Smith   Level: developer
1770f1580f4eSBarry Smith 
1771562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()`
1772fe2efc57SMark @*/
1773d71ae5a4SJacob Faibussowitsch PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[])
1774d71ae5a4SJacob Faibussowitsch {
1775fe2efc57SMark   PetscFunctionBegin;
1776fe2efc57SMark   PetscValidHeaderSpecific(A, PC_CLASSID, 1);
17779566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
17783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1779fe2efc57SMark }
1780fe2efc57SMark 
1781ffeef943SBarry Smith /*@
1782f1580f4eSBarry Smith   PCView - Prints information about the `PC`
17834b9ad928SBarry Smith 
1784c3339decSBarry Smith   Collective
17854b9ad928SBarry Smith 
17864b9ad928SBarry Smith   Input Parameters:
17870b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
17880b4b7b1cSBarry Smith - viewer - optional `PetscViewer` visualization context
17894b9ad928SBarry Smith 
17900b4b7b1cSBarry Smith   Level: intermediate
179120f4b53cSBarry Smith 
1792f1580f4eSBarry Smith   Notes:
17934b9ad928SBarry Smith   The available visualization contexts include
1794f1580f4eSBarry Smith +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
1795f1580f4eSBarry Smith -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
17964b9ad928SBarry Smith   output where only the first processor opens
17974b9ad928SBarry Smith   the file. All other processors send their
17984b9ad928SBarry Smith   data to the first processor to print.
17994b9ad928SBarry Smith 
18004b9ad928SBarry Smith   The user can open an alternative visualization contexts with
1801f1580f4eSBarry Smith   `PetscViewerASCIIOpen()` (output to a specified file).
18024b9ad928SBarry Smith 
18030b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewer`, `PetscViewerType`, `KSPView()`, `PetscViewerASCIIOpen()`
18044b9ad928SBarry Smith @*/
1805d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView(PC pc, PetscViewer viewer)
1806d71ae5a4SJacob Faibussowitsch {
180719fd82e9SBarry Smith   PCType            cstr;
18086cd81132SPierre Jolivet   PetscViewerFormat format;
18096cd81132SPierre Jolivet   PetscBool         iascii, isstring, isbinary, isdraw, pop = PETSC_FALSE;
1810e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1811536b137fSBarry Smith   PetscBool issaws;
18120acecf5bSBarry Smith #endif
18134b9ad928SBarry Smith 
18144b9ad928SBarry Smith   PetscFunctionBegin;
18150700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
181648a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
18170700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1818c9780b6fSBarry Smith   PetscCheckSameComm(pc, 1, viewer, 2);
18194b9ad928SBarry Smith 
18209566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
18219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
18229566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
18239566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1824e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
18259566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
18260acecf5bSBarry Smith #endif
1827219b1045SBarry Smith 
182832077d6dSBarry Smith   if (iascii) {
18299566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer));
183048a46eb9SPierre Jolivet     if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, "  PC has not been set up so information may be incomplete\n"));
18319566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1832dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
1834834dbeb0SBarry Smith     if (pc->mat) {
18356cd81132SPierre Jolivet       PetscCall(PetscViewerGetFormat(viewer, &format));
18366cd81132SPierre Jolivet       if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
18379566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
18386cd81132SPierre Jolivet         pop = PETSC_TRUE;
18396cd81132SPierre Jolivet       }
18404b9ad928SBarry Smith       if (pc->pmat == pc->mat) {
18419566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix = precond matrix:\n"));
18429566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
18439566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
18449566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
18454b9ad928SBarry Smith       } else {
1846834dbeb0SBarry Smith         if (pc->pmat) {
18479566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix followed by preconditioner matrix:\n"));
18484b9ad928SBarry Smith         } else {
18499566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix:\n"));
18504b9ad928SBarry Smith         }
18519566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
18529566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
18539566063dSJacob Faibussowitsch         if (pc->pmat) PetscCall(MatView(pc->pmat, viewer));
18549566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
18554b9ad928SBarry Smith       }
18566cd81132SPierre Jolivet       if (pop) PetscCall(PetscViewerPopFormat(viewer));
18574b9ad928SBarry Smith     }
18584b9ad928SBarry Smith   } else if (isstring) {
18599566063dSJacob Faibussowitsch     PetscCall(PCGetType(pc, &cstr));
18609566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr));
1861dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18629566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18639566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18645b0b0462SBarry Smith   } else if (isbinary) {
186555849f57SBarry Smith     PetscInt    classid = PC_FILE_CLASSID;
186655849f57SBarry Smith     MPI_Comm    comm;
186755849f57SBarry Smith     PetscMPIInt rank;
186855849f57SBarry Smith     char        type[256];
186955849f57SBarry Smith 
18709566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
18719566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1872dd400576SPatrick Sanan     if (rank == 0) {
18739566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
18749566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256));
18759566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
187655849f57SBarry Smith     }
1877dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
1878219b1045SBarry Smith   } else if (isdraw) {
1879219b1045SBarry Smith     PetscDraw draw;
1880d9884438SBarry Smith     char      str[25];
188189fd9fafSBarry Smith     PetscReal x, y, bottom, h;
1882d9884438SBarry Smith     PetscInt  n;
1883219b1045SBarry Smith 
18849566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
18859566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
18861d840656SPeter Brune     if (pc->mat) {
18879566063dSJacob Faibussowitsch       PetscCall(MatGetSize(pc->mat, &n, NULL));
188863a3b9bcSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n));
18891d840656SPeter Brune     } else {
18909566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name));
18911d840656SPeter Brune     }
18929566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
189389fd9fafSBarry Smith     bottom = y - h;
18949566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
1895dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18969566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
1897e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1898536b137fSBarry Smith   } else if (issaws) {
1899d45a07a7SBarry Smith     PetscMPIInt rank;
1900d45a07a7SBarry Smith 
19019566063dSJacob Faibussowitsch     PetscCall(PetscObjectName((PetscObject)pc));
19029566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
190348a46eb9SPierre Jolivet     if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer));
19049566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
19059566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
19060acecf5bSBarry Smith #endif
19074b9ad928SBarry Smith   }
19083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19094b9ad928SBarry Smith }
19104b9ad928SBarry Smith 
19114b9ad928SBarry Smith /*@C
19120b4b7b1cSBarry Smith   PCRegister -  Adds a method (`PCType`) to the PETSc preconditioner package.
19131c84c290SBarry Smith 
1914cc4c1da9SBarry Smith   Not collective. No Fortran Support
19151c84c290SBarry Smith 
19161c84c290SBarry Smith   Input Parameters:
191720f4b53cSBarry Smith + sname    - name of a new user-defined solver
19180b4b7b1cSBarry Smith - function - routine to create the method context which will be stored in a `PC` when `PCSetType()` is called
19191c84c290SBarry Smith 
1920feefa0e1SJacob Faibussowitsch   Example Usage:
19211c84c290SBarry Smith .vb
1922bdf89e91SBarry Smith    PCRegister("my_solver", MySolverCreate);
19231c84c290SBarry Smith .ve
19241c84c290SBarry Smith 
19251c84c290SBarry Smith   Then, your solver can be chosen with the procedural interface via
19261c84c290SBarry Smith $     PCSetType(pc, "my_solver")
19271c84c290SBarry Smith   or at runtime via the option
19281c84c290SBarry Smith $     -pc_type my_solver
19294b9ad928SBarry Smith 
19304b9ad928SBarry Smith   Level: advanced
19311c84c290SBarry Smith 
193220f4b53cSBarry Smith   Note:
19330b4b7b1cSBarry Smith   A simpler alternative to using `PCRegister()` for an application specific preconditioner is to use a `PC` of `PCType` `PCSHELL` and
19340b4b7b1cSBarry Smith   provide your customizations with `PCShellSetContext()` and `PCShellSetApply()`
19350b4b7b1cSBarry Smith 
193620f4b53cSBarry Smith   `PCRegister()` may be called multiple times to add several user-defined preconditioners.
193720f4b53cSBarry Smith 
19380b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCType`, `PCRegisterAll()`, `PCSetType()`, `PCShellSetContext()`, `PCShellSetApply()`, `PCSHELL`
19394b9ad928SBarry Smith @*/
1940d71ae5a4SJacob Faibussowitsch PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC))
1941d71ae5a4SJacob Faibussowitsch {
19424b9ad928SBarry Smith   PetscFunctionBegin;
19439566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
19449566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCList, sname, function));
19453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19464b9ad928SBarry Smith }
19474b9ad928SBarry Smith 
1948d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y)
1949d71ae5a4SJacob Faibussowitsch {
1950186a3e20SStefano Zampini   PC pc;
1951186a3e20SStefano Zampini 
1952186a3e20SStefano Zampini   PetscFunctionBegin;
19539566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &pc));
19549566063dSJacob Faibussowitsch   PetscCall(PCApply(pc, X, Y));
19553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1956186a3e20SStefano Zampini }
1957186a3e20SStefano Zampini 
19585d83a8b1SBarry Smith /*@
19590b4b7b1cSBarry Smith   PCComputeOperator - Computes the explicit preconditioned operator as a matrix `Mat`.
19604b9ad928SBarry Smith 
1961c3339decSBarry Smith   Collective
19624b9ad928SBarry Smith 
1963d8d19677SJose E. Roman   Input Parameters:
19640b4b7b1cSBarry Smith + pc      - the `PC` preconditioner object
1965562efe2eSBarry Smith - mattype - the `MatType` to be used for the operator
19664b9ad928SBarry Smith 
19674b9ad928SBarry Smith   Output Parameter:
1968a5b23f4aSJose E. Roman . mat - the explicit preconditioned operator
19694b9ad928SBarry Smith 
197020f4b53cSBarry Smith   Level: advanced
197120f4b53cSBarry Smith 
1972f1580f4eSBarry Smith   Note:
1973186a3e20SStefano Zampini   This computation is done by applying the operators to columns of the identity matrix.
1974186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
1975562efe2eSBarry Smith   Currently, this routine uses a dense matrix format when `mattype` == `NULL`
19764b9ad928SBarry Smith 
19770b4b7b1cSBarry Smith   Developer Note:
19780b4b7b1cSBarry Smith   This should be called `PCCreateExplicitOperator()`
19790b4b7b1cSBarry Smith 
1980562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPComputeOperator()`, `MatType`
19814b9ad928SBarry Smith @*/
1982d71ae5a4SJacob Faibussowitsch PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat)
1983d71ae5a4SJacob Faibussowitsch {
1984186a3e20SStefano Zampini   PetscInt N, M, m, n;
1985186a3e20SStefano Zampini   Mat      A, Apc;
19864b9ad928SBarry Smith 
19874b9ad928SBarry Smith   PetscFunctionBegin;
19880700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
19894f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
19909566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, &A, NULL));
19919566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
19929566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
19939566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc));
19949566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (void (*)(void))MatMult_PC));
19959566063dSJacob Faibussowitsch   PetscCall(MatComputeOperator(Apc, mattype, mat));
19969566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Apc));
19973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19984b9ad928SBarry Smith }
19994b9ad928SBarry Smith 
20006c237d78SBarry Smith /*@
20010b4b7b1cSBarry Smith   PCSetCoordinates - sets the coordinates of all the nodes (degrees of freedom in the vector) on the local process
20026c237d78SBarry Smith 
2003c3339decSBarry Smith   Collective
20046c237d78SBarry Smith 
20056c237d78SBarry Smith   Input Parameters:
20060b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
20076c237d78SBarry Smith . dim    - the dimension of the coordinates 1, 2, or 3
200814893cbeSStefano Zampini . nloc   - the blocked size of the coordinates array
200914893cbeSStefano Zampini - coords - the coordinates array
20106c237d78SBarry Smith 
20116c237d78SBarry Smith   Level: intermediate
20126c237d78SBarry Smith 
20130b4b7b1cSBarry Smith   Notes:
201420f4b53cSBarry Smith   `coords` is an array of the dim coordinates for the nodes on
201520f4b53cSBarry Smith   the local processor, of size `dim`*`nloc`.
20166e415bd2SNuno Nobre   If there are 108 equations (dofs) on a processor
20176e415bd2SNuno Nobre   for a 3d displacement finite element discretization of elasticity (so
201814893cbeSStefano Zampini   that there are nloc = 36 = 108/3 nodes) then the array must have 108
20196c237d78SBarry Smith   double precision values (ie, 3 * 36).  These x y z coordinates
20206c237d78SBarry Smith   should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
20216c237d78SBarry Smith   ... , N-1.z ].
20226c237d78SBarry Smith 
20230b4b7b1cSBarry Smith   The information provided here can be used by some preconditioners, such as `PCGAMG`, to produce a better preconditioner.
20240b4b7b1cSBarry Smith   See also  `MatSetNearNullSpace()`.
20250b4b7b1cSBarry Smith 
2026562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `MatSetNearNullSpace()`
20276c237d78SBarry Smith @*/
2028d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[])
2029d71ae5a4SJacob Faibussowitsch {
20306c237d78SBarry Smith   PetscFunctionBegin;
203132954da3SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
203232954da3SStefano Zampini   PetscValidLogicalCollectiveInt(pc, dim, 2);
203322794d57SStefano Zampini   PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords));
20343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20356c237d78SBarry Smith }
2036fd2dd295SFande Kong 
2037fd2dd295SFande Kong /*@
2038fd2dd295SFande Kong   PCGetInterpolations - Gets interpolation matrices for all levels (except level 0)
2039fd2dd295SFande Kong 
2040c3339decSBarry Smith   Logically Collective
2041fd2dd295SFande Kong 
2042d8d19677SJose E. Roman   Input Parameter:
2043d8d19677SJose E. Roman . pc - the precondition context
2044fd2dd295SFande Kong 
2045d8d19677SJose E. Roman   Output Parameters:
2046d8d19677SJose E. Roman + num_levels     - the number of levels
2047562efe2eSBarry Smith - interpolations - the interpolation matrices (size of `num_levels`-1)
2048fd2dd295SFande Kong 
2049fd2dd295SFande Kong   Level: advanced
2050fd2dd295SFande Kong 
2051562efe2eSBarry Smith   Developer Note:
2052f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
2053fd2dd295SFande Kong 
2054562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()`
2055fd2dd295SFande Kong @*/
2056d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[])
2057d71ae5a4SJacob Faibussowitsch {
2058fd2dd295SFande Kong   PetscFunctionBegin;
2059fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20604f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20614f572ea9SToby Isaac   PetscAssertPointer(interpolations, 3);
2062cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations));
20633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2064fd2dd295SFande Kong }
2065fd2dd295SFande Kong 
2066fd2dd295SFande Kong /*@
2067fd2dd295SFande Kong   PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level)
2068fd2dd295SFande Kong 
2069c3339decSBarry Smith   Logically Collective
2070fd2dd295SFande Kong 
2071d8d19677SJose E. Roman   Input Parameter:
2072d8d19677SJose E. Roman . pc - the precondition context
2073fd2dd295SFande Kong 
2074d8d19677SJose E. Roman   Output Parameters:
2075d8d19677SJose E. Roman + num_levels      - the number of levels
2076562efe2eSBarry Smith - coarseOperators - the coarse operator matrices (size of `num_levels`-1)
2077fd2dd295SFande Kong 
2078fd2dd295SFande Kong   Level: advanced
2079fd2dd295SFande Kong 
2080562efe2eSBarry Smith   Developer Note:
2081f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
2082fd2dd295SFande Kong 
2083562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()`
2084fd2dd295SFande Kong @*/
2085d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
2086d71ae5a4SJacob Faibussowitsch {
2087fd2dd295SFande Kong   PetscFunctionBegin;
2088fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20894f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20904f572ea9SToby Isaac   PetscAssertPointer(coarseOperators, 3);
2091cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators));
20923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2093fd2dd295SFande Kong }
2094