xref: /petsc/src/ksp/pc/interface/precon.c (revision bded88ea7b7af78c219a7725f5ce5bb2a37a9d77)
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;
9*bded88eaSPierre 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 
5354b9ad928SBarry Smith /*@
536562efe2eSBarry Smith   PCMatApply - Applies the preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApply()`, `Y` and `X` must be different matrices.
537c41ea70eSPierre Jolivet 
538c3339decSBarry Smith   Collective
539c677e75fSPierre Jolivet 
540c677e75fSPierre Jolivet   Input Parameters:
5410b4b7b1cSBarry Smith + pc - the `PC` preconditioner context
542c677e75fSPierre Jolivet - X  - block of input vectors
543c677e75fSPierre Jolivet 
544c677e75fSPierre Jolivet   Output Parameter:
545c677e75fSPierre Jolivet . Y - block of output vectors
546c677e75fSPierre Jolivet 
547c677e75fSPierre Jolivet   Level: developer
548c677e75fSPierre Jolivet 
549562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `KSPMatSolve()`
550c677e75fSPierre Jolivet @*/
551d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMatApply(PC pc, Mat X, Mat Y)
552d71ae5a4SJacob Faibussowitsch {
553c677e75fSPierre Jolivet   Mat       A;
554c677e75fSPierre Jolivet   Vec       cy, cx;
555bd82155bSPierre Jolivet   PetscInt  m1, M1, m2, M2, n1, N1, n2, N2, m3, M3, n3, N3;
556c677e75fSPierre Jolivet   PetscBool match;
557c677e75fSPierre Jolivet 
558c677e75fSPierre Jolivet   PetscFunctionBegin;
559c677e75fSPierre Jolivet   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
560e2b200f6SPierre Jolivet   PetscValidHeaderSpecific(X, MAT_CLASSID, 2);
561e2b200f6SPierre Jolivet   PetscValidHeaderSpecific(Y, MAT_CLASSID, 3);
562e2b200f6SPierre Jolivet   PetscCheckSameComm(pc, 1, X, 2);
563e2b200f6SPierre Jolivet   PetscCheckSameComm(pc, 1, Y, 3);
5647827d75bSBarry Smith   PetscCheck(Y != X, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "Y and X must be different matrices");
5659566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &A));
5669566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m3, &n3));
5679566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m2, &n2));
5689566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &m1, &n1));
5699566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M3, &N3));
5709566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &M2, &N2));
5719566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Y, &M1, &N1));
57263a3b9bcSJacob 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);
57363a3b9bcSJacob 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);
57463a3b9bcSJacob 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);
5759566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &match, MATSEQDENSE, MATMPIDENSE, ""));
57628b400f6SJacob Faibussowitsch   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of output vectors not stored in a dense Mat");
5779566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
57828b400f6SJacob Faibussowitsch   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of input vectors not stored in a dense Mat");
5799566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
580c677e75fSPierre Jolivet   if (pc->ops->matapply) {
5819566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0));
582dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, matapply, X, Y);
5839566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0));
584c677e75fSPierre Jolivet   } else {
5859566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "PC type %s applying column by column\n", ((PetscObject)pc)->type_name));
586bd82155bSPierre Jolivet     for (n1 = 0; n1 < N1; ++n1) {
5879566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecRead(X, n1, &cx));
5889566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecWrite(Y, n1, &cy));
5899566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, cx, cy));
5909566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecWrite(Y, n1, &cy));
5919566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecRead(X, n1, &cx));
592c677e75fSPierre Jolivet     }
593c677e75fSPierre Jolivet   }
5943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
595c677e75fSPierre Jolivet }
596c677e75fSPierre Jolivet 
597c677e75fSPierre Jolivet /*@
5984b9ad928SBarry Smith   PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
5994b9ad928SBarry Smith 
600c3339decSBarry Smith   Collective
6014b9ad928SBarry Smith 
6024b9ad928SBarry Smith   Input Parameters:
6030b4b7b1cSBarry Smith + pc - the `PC` preconditioner context
6044b9ad928SBarry Smith - x  - input vector
6054b9ad928SBarry Smith 
6064b9ad928SBarry Smith   Output Parameter:
6074b9ad928SBarry Smith . y - output vector
6084b9ad928SBarry Smith 
60920f4b53cSBarry Smith   Level: developer
61020f4b53cSBarry Smith 
611f1580f4eSBarry Smith   Note:
612f1580f4eSBarry Smith   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
6134b9ad928SBarry Smith 
614562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricRight()`
6154b9ad928SBarry Smith @*/
616d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricLeft(PC pc, Vec x, Vec y)
617d71ae5a4SJacob Faibussowitsch {
6184b9ad928SBarry Smith   PetscFunctionBegin;
6190700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6200700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6210700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6227827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
623e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
6249566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
6259566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
6269566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplySymmetricLeft, pc, x, y, 0));
627dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applysymmetricleft, x, y);
6289566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplySymmetricLeft, pc, x, y, 0));
6299566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
630e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
6313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6324b9ad928SBarry Smith }
6334b9ad928SBarry Smith 
6344b9ad928SBarry Smith /*@
6354b9ad928SBarry Smith   PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
6364b9ad928SBarry Smith 
637c3339decSBarry Smith   Collective
6384b9ad928SBarry Smith 
6394b9ad928SBarry Smith   Input Parameters:
6400b4b7b1cSBarry Smith + pc - the `PC` preconditioner context
6414b9ad928SBarry Smith - x  - input vector
6424b9ad928SBarry Smith 
6434b9ad928SBarry Smith   Output Parameter:
6444b9ad928SBarry Smith . y - output vector
6454b9ad928SBarry Smith 
6464b9ad928SBarry Smith   Level: developer
6474b9ad928SBarry Smith 
648f1580f4eSBarry Smith   Note:
649f1580f4eSBarry Smith   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
6504b9ad928SBarry Smith 
651562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricLeft()`
6524b9ad928SBarry Smith @*/
653d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricRight(PC pc, Vec x, Vec y)
654d71ae5a4SJacob Faibussowitsch {
6554b9ad928SBarry Smith   PetscFunctionBegin;
6560700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6570700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6580700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6597827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
660e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
6619566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
6629566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
6639566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplySymmetricRight, pc, x, y, 0));
664dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applysymmetricright, x, y);
6659566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplySymmetricRight, pc, x, y, 0));
6669566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
667e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
6683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6694b9ad928SBarry Smith }
6704b9ad928SBarry Smith 
6714b9ad928SBarry Smith /*@
6724b9ad928SBarry Smith   PCApplyTranspose - Applies the transpose of preconditioner to a vector.
6734b9ad928SBarry Smith 
674c3339decSBarry Smith   Collective
6754b9ad928SBarry Smith 
6764b9ad928SBarry Smith   Input Parameters:
6770b4b7b1cSBarry Smith + pc - the `PC` preconditioner context
6784b9ad928SBarry Smith - x  - input vector
6794b9ad928SBarry Smith 
6804b9ad928SBarry Smith   Output Parameter:
6814b9ad928SBarry Smith . y - output vector
6824b9ad928SBarry Smith 
68320f4b53cSBarry Smith   Level: developer
68420f4b53cSBarry Smith 
685f1580f4eSBarry Smith   Note:
68695452b02SPatrick Sanan   For complex numbers this applies the non-Hermitian transpose.
6874c97465dSBarry Smith 
688562efe2eSBarry Smith   Developer Note:
689f1580f4eSBarry Smith   We need to implement a `PCApplyHermitianTranspose()`
6904c97465dSBarry Smith 
691562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyBAorAB()`, `PCApplyBAorABTranspose()`, `PCApplyTransposeExists()`
6924b9ad928SBarry Smith @*/
693d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTranspose(PC pc, Vec x, Vec y)
694d71ae5a4SJacob Faibussowitsch {
6954b9ad928SBarry Smith   PetscFunctionBegin;
6960700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6970700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6980700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6997827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
700e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
7019566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
7029566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
7039566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
704dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applytranspose, x, y);
7059566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
7069566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
707e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
7083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7094b9ad928SBarry Smith }
7104b9ad928SBarry Smith 
7114b9ad928SBarry Smith /*@
7129cc050e5SLisandro Dalcin   PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
7134b9ad928SBarry Smith 
714c3339decSBarry Smith   Collective
7154b9ad928SBarry Smith 
716f1580f4eSBarry Smith   Input Parameter:
7170b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
7184b9ad928SBarry Smith 
7194b9ad928SBarry Smith   Output Parameter:
720f1580f4eSBarry Smith . flg - `PETSC_TRUE` if a transpose operation is defined
7214b9ad928SBarry Smith 
7224b9ad928SBarry Smith   Level: developer
7234b9ad928SBarry Smith 
724562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`
7254b9ad928SBarry Smith @*/
726d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTransposeExists(PC pc, PetscBool *flg)
727d71ae5a4SJacob Faibussowitsch {
7284b9ad928SBarry Smith   PetscFunctionBegin;
7290700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
7304f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
7316ce5e81cSLisandro Dalcin   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
7326ce5e81cSLisandro Dalcin   else *flg = PETSC_FALSE;
7333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7344b9ad928SBarry Smith }
7354b9ad928SBarry Smith 
7364b9ad928SBarry Smith /*@
737562efe2eSBarry Smith   PCApplyBAorAB - Applies the preconditioner and operator to a vector. $y = B*A*x $ or $ y = A*B*x$.
7384b9ad928SBarry Smith 
739c3339decSBarry Smith   Collective
7404b9ad928SBarry Smith 
7414b9ad928SBarry Smith   Input Parameters:
7420b4b7b1cSBarry Smith + pc   - the `PC` preconditioner context
743f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
7444b9ad928SBarry Smith . x    - input vector
7454b9ad928SBarry Smith - work - work vector
7464b9ad928SBarry Smith 
7474b9ad928SBarry Smith   Output Parameter:
7484b9ad928SBarry Smith . y - output vector
7494b9ad928SBarry Smith 
7504b9ad928SBarry Smith   Level: developer
7514b9ad928SBarry Smith 
752f1580f4eSBarry Smith   Note:
753562efe2eSBarry 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.
754562efe2eSBarry Smith   The specific `KSPSolve()` method must also be written to handle the post-solve "correction" for the diagonal scaling.
75513e0d083SBarry Smith 
756562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorABTranspose()`
7574b9ad928SBarry Smith @*/
758d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorAB(PC pc, PCSide side, Vec x, Vec y, Vec work)
759d71ae5a4SJacob Faibussowitsch {
7604b9ad928SBarry Smith   PetscFunctionBegin;
7610700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
762a69c7061SStefano Zampini   PetscValidLogicalCollectiveEnum(pc, side, 2);
7630700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
7640700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
7650700a824SBarry Smith   PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
766a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, x, 3);
767a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, y, 4);
768a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, work, 5);
7697827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
7707827d75bSBarry 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");
7717827d75bSBarry Smith   PetscCheck(!pc->diagonalscale || side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot include diagonal scaling with symmetric preconditioner application");
772e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
7734b9ad928SBarry Smith 
7749566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
7754b9ad928SBarry Smith   if (pc->diagonalscale) {
7764b9ad928SBarry Smith     if (pc->ops->applyBA) {
7774b9ad928SBarry Smith       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
7789566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &work2));
7799566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, work2));
780dbbe0bcdSBarry Smith       PetscUseTypeMethod(pc, applyBA, side, work2, y, work);
7819566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7829566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&work2));
7834b9ad928SBarry Smith     } else if (side == PC_RIGHT) {
7849566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, y));
7859566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, y, work));
7869566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7879566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7884b9ad928SBarry Smith     } else if (side == PC_LEFT) {
7899566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, y));
7909566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, y, work));
7919566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, work, y));
7929566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7937827d75bSBarry Smith     } else PetscCheck(side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot provide diagonal scaling with symmetric application of preconditioner");
7944b9ad928SBarry Smith   } else {
7954b9ad928SBarry Smith     if (pc->ops->applyBA) {
796dbbe0bcdSBarry Smith       PetscUseTypeMethod(pc, applyBA, side, x, y, work);
7974b9ad928SBarry Smith     } else if (side == PC_RIGHT) {
7989566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, x, work));
7999566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
8004b9ad928SBarry Smith     } else if (side == PC_LEFT) {
8019566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, x, work));
8029566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, work, y));
8034b9ad928SBarry Smith     } else if (side == PC_SYMMETRIC) {
8044b9ad928SBarry Smith       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
8059566063dSJacob Faibussowitsch       PetscCall(PCApplySymmetricRight(pc, x, work));
8069566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
8079566063dSJacob Faibussowitsch       PetscCall(VecCopy(y, work));
8089566063dSJacob Faibussowitsch       PetscCall(PCApplySymmetricLeft(pc, work, y));
8094b9ad928SBarry Smith     }
8104b9ad928SBarry Smith   }
811e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8134b9ad928SBarry Smith }
8144b9ad928SBarry Smith 
8154b9ad928SBarry Smith /*@
8164b9ad928SBarry Smith   PCApplyBAorABTranspose - Applies the transpose of the preconditioner
817562efe2eSBarry Smith   and operator to a vector. That is, applies $B^T * A^T$ with left preconditioning,
818562efe2eSBarry Smith   NOT $(B*A)^T = A^T*B^T$.
8194b9ad928SBarry Smith 
820c3339decSBarry Smith   Collective
8214b9ad928SBarry Smith 
8224b9ad928SBarry Smith   Input Parameters:
8230b4b7b1cSBarry Smith + pc   - the `PC` preconditioner context
824f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
8254b9ad928SBarry Smith . x    - input vector
8264b9ad928SBarry Smith - work - work vector
8274b9ad928SBarry Smith 
8284b9ad928SBarry Smith   Output Parameter:
8294b9ad928SBarry Smith . y - output vector
8304b9ad928SBarry Smith 
83120f4b53cSBarry Smith   Level: developer
83220f4b53cSBarry Smith 
833f1580f4eSBarry Smith   Note:
834562efe2eSBarry 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
835562efe2eSBarry Smith   defined by $B^T$. This is why this has the funny form that it computes $B^T * A^T$
8369b3038f0SBarry Smith 
837562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorAB()`
8384b9ad928SBarry Smith @*/
839d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorABTranspose(PC pc, PCSide side, Vec x, Vec y, Vec work)
840d71ae5a4SJacob Faibussowitsch {
8414b9ad928SBarry Smith   PetscFunctionBegin;
8420700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8430700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
8440700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
8450700a824SBarry Smith   PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
8467827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
847e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
8484b9ad928SBarry Smith   if (pc->ops->applyBAtranspose) {
849dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, applyBAtranspose, side, x, y, work);
850e0f629ddSJacob Faibussowitsch     if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8513ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
8524b9ad928SBarry Smith   }
8537827d75bSBarry Smith   PetscCheck(side == PC_LEFT || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right or left");
8544b9ad928SBarry Smith 
8559566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
8564b9ad928SBarry Smith   if (side == PC_RIGHT) {
8579566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(pc, x, work));
8589566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->mat, work, y));
8594b9ad928SBarry Smith   } else if (side == PC_LEFT) {
8609566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->mat, x, work));
8619566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(pc, work, y));
8624b9ad928SBarry Smith   }
8634b9ad928SBarry Smith   /* add support for PC_SYMMETRIC */
864e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8664b9ad928SBarry Smith }
8674b9ad928SBarry Smith 
8684b9ad928SBarry Smith /*@
8694b9ad928SBarry Smith   PCApplyRichardsonExists - Determines whether a particular preconditioner has a
8704b9ad928SBarry Smith   built-in fast application of Richardson's method.
8714b9ad928SBarry Smith 
8724b9ad928SBarry Smith   Not Collective
8734b9ad928SBarry Smith 
8744b9ad928SBarry Smith   Input Parameter:
8754b9ad928SBarry Smith . pc - the preconditioner
8764b9ad928SBarry Smith 
8774b9ad928SBarry Smith   Output Parameter:
878f1580f4eSBarry Smith . exists - `PETSC_TRUE` or `PETSC_FALSE`
8794b9ad928SBarry Smith 
8804b9ad928SBarry Smith   Level: developer
8814b9ad928SBarry Smith 
88239b1ba1bSPierre Jolivet .seealso: [](ch_ksp), `PC`, `KSPRICHARDSON`, `PCApplyRichardson()`
8834b9ad928SBarry Smith @*/
884d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyRichardsonExists(PC pc, PetscBool *exists)
885d71ae5a4SJacob Faibussowitsch {
8864b9ad928SBarry Smith   PetscFunctionBegin;
8870700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8884f572ea9SToby Isaac   PetscAssertPointer(exists, 2);
8894b9ad928SBarry Smith   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
8904b9ad928SBarry Smith   else *exists = PETSC_FALSE;
8913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8924b9ad928SBarry Smith }
8934b9ad928SBarry Smith 
8944b9ad928SBarry Smith /*@
8954b9ad928SBarry Smith   PCApplyRichardson - Applies several steps of Richardson iteration with
8964b9ad928SBarry Smith   the particular preconditioner. This routine is usually used by the
8974b9ad928SBarry Smith   Krylov solvers and not the application code directly.
8984b9ad928SBarry Smith 
899c3339decSBarry Smith   Collective
9004b9ad928SBarry Smith 
9014b9ad928SBarry Smith   Input Parameters:
9020b4b7b1cSBarry Smith + pc        - the `PC` preconditioner context
903dd8e379bSPierre Jolivet . b         - the right-hand side
9044b9ad928SBarry Smith . w         - one work vector
9054b9ad928SBarry Smith . rtol      - relative decrease in residual norm convergence criteria
90670441072SBarry Smith . abstol    - absolute residual norm convergence criteria
9074b9ad928SBarry Smith . dtol      - divergence residual norm increase criteria
9087319c654SBarry Smith . its       - the number of iterations to apply.
9097319c654SBarry Smith - guesszero - if the input x contains nonzero initial guess
9104b9ad928SBarry Smith 
911d8d19677SJose E. Roman   Output Parameters:
9124d0a8057SBarry Smith + outits - number of iterations actually used (for SOR this always equals its)
9134d0a8057SBarry Smith . reason - the reason the apply terminated
914f1580f4eSBarry Smith - y      - the solution (also contains initial guess if guesszero is `PETSC_FALSE`
9154b9ad928SBarry Smith 
91620f4b53cSBarry Smith   Level: developer
91720f4b53cSBarry Smith 
9184b9ad928SBarry Smith   Notes:
9194b9ad928SBarry Smith   Most preconditioners do not support this function. Use the command
920f1580f4eSBarry Smith   `PCApplyRichardsonExists()` to determine if one does.
9214b9ad928SBarry Smith 
922f1580f4eSBarry Smith   Except for the `PCMG` this routine ignores the convergence tolerances
9234b9ad928SBarry Smith   and always runs for the number of iterations
9244b9ad928SBarry Smith 
925562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyRichardsonExists()`
9264b9ad928SBarry Smith @*/
927d71ae5a4SJacob 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)
928d71ae5a4SJacob Faibussowitsch {
9294b9ad928SBarry Smith   PetscFunctionBegin;
9300700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9310700a824SBarry Smith   PetscValidHeaderSpecific(b, VEC_CLASSID, 2);
9320700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
9330700a824SBarry Smith   PetscValidHeaderSpecific(w, VEC_CLASSID, 4);
9347827d75bSBarry Smith   PetscCheck(b != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "b and y must be different vectors");
9359566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
936dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applyrichardson, b, y, w, rtol, abstol, dtol, its, guesszero, outits, reason);
9373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9384b9ad928SBarry Smith }
9394b9ad928SBarry Smith 
940422a814eSBarry Smith /*@
941f1580f4eSBarry Smith   PCSetFailedReason - Sets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
9421b2b9847SBarry Smith 
943c3339decSBarry Smith   Logically Collective
9441b2b9847SBarry Smith 
945d8d19677SJose E. Roman   Input Parameters:
9460b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
9470b4b7b1cSBarry Smith - reason - the reason it failed
9481b2b9847SBarry Smith 
9491b2b9847SBarry Smith   Level: advanced
9501b2b9847SBarry Smith 
951562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCFailedReason`
9521b2b9847SBarry Smith @*/
953d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFailedReason(PC pc, PCFailedReason reason)
954d71ae5a4SJacob Faibussowitsch {
9551b2b9847SBarry Smith   PetscFunctionBegin;
9566479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9571b2b9847SBarry Smith   pc->failedreason = reason;
9583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9591b2b9847SBarry Smith }
9601b2b9847SBarry Smith 
9611b2b9847SBarry Smith /*@
962f1580f4eSBarry Smith   PCGetFailedReason - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
963422a814eSBarry Smith 
964f1580f4eSBarry Smith   Not Collective
9651b2b9847SBarry Smith 
9661b2b9847SBarry Smith   Input Parameter:
9670b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
9681b2b9847SBarry Smith 
9691b2b9847SBarry Smith   Output Parameter:
9701b2b9847SBarry Smith . reason - the reason it failed
9711b2b9847SBarry Smith 
97220f4b53cSBarry Smith   Level: advanced
97320f4b53cSBarry Smith 
974f1580f4eSBarry Smith   Note:
9750b4b7b1cSBarry Smith   After a call to `KSPCheckDot()` or  `KSPCheckNorm()` inside a `KSPSolve()` or a call to `PCReduceFailedReason()`
9760b4b7b1cSBarry Smith   this is the maximum reason over all MPI processes in the `PC` communicator and hence logically collective.
9779a6c2652SBarry Smith   Otherwise it returns the local value.
9781b2b9847SBarry Smith 
9790b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCSetFailedReason()`, `PCFailedReason`
9801b2b9847SBarry Smith @*/
9819a6c2652SBarry Smith PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason)
982d71ae5a4SJacob Faibussowitsch {
9831b2b9847SBarry Smith   PetscFunctionBegin;
9846479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9851b2b9847SBarry Smith   if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled;
9861b2b9847SBarry Smith   else *reason = pc->failedreason;
9873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9881b2b9847SBarry Smith }
989422a814eSBarry Smith 
9906479e835SStefano Zampini /*@
9916479e835SStefano Zampini   PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC`
9926479e835SStefano Zampini 
9936479e835SStefano Zampini   Collective
9946479e835SStefano Zampini 
9956479e835SStefano Zampini   Input Parameter:
9960b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
9976479e835SStefano Zampini 
9986479e835SStefano Zampini   Level: advanced
9996479e835SStefano Zampini 
10006479e835SStefano Zampini   Note:
1001562efe2eSBarry Smith   Different MPI processes may have different reasons or no reason, see `PCGetFailedReason()`. This routine
10026479e835SStefano Zampini   makes them have a common value (failure if any MPI process had a failure).
10036479e835SStefano Zampini 
10040b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`, `PCFailedReason`
10056479e835SStefano Zampini @*/
10066479e835SStefano Zampini PetscErrorCode PCReduceFailedReason(PC pc)
10076479e835SStefano Zampini {
10086479e835SStefano Zampini   PetscInt buf;
10096479e835SStefano Zampini 
10106479e835SStefano Zampini   PetscFunctionBegin;
10116479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
10126479e835SStefano Zampini   buf = (PetscInt)pc->failedreason;
1013462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
10146479e835SStefano Zampini   pc->failedreason = (PCFailedReason)buf;
10156479e835SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
10166479e835SStefano Zampini }
10176479e835SStefano Zampini 
10184b9ad928SBarry Smith /*
10194b9ad928SBarry Smith       a setupcall of 0 indicates never setup,
102023ee1639SBarry Smith                      1 indicates has been previously setup
1021422a814eSBarry Smith                     -1 indicates a PCSetUp() was attempted and failed
10224b9ad928SBarry Smith */
10234b9ad928SBarry Smith /*@
10240b4b7b1cSBarry Smith   PCSetUp - Prepares for the use of a preconditioner. Performs all the one-time operations needed before the preconditioner
10250b4b7b1cSBarry Smith   can be used with `PCApply()`
10264b9ad928SBarry Smith 
1027c3339decSBarry Smith   Collective
10284b9ad928SBarry Smith 
10294b9ad928SBarry Smith   Input Parameter:
10300b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
10314b9ad928SBarry Smith 
10324b9ad928SBarry Smith   Level: developer
10334b9ad928SBarry Smith 
10340b4b7b1cSBarry Smith   Notes:
10350b4b7b1cSBarry Smith   For example, for `PCLU` this will compute the factorization.
10360b4b7b1cSBarry Smith 
10370b4b7b1cSBarry Smith   This is called automatically by `KSPSetUp()` or `PCApply()` so rarely needs to be called directly.
10380b4b7b1cSBarry Smith 
10390b4b7b1cSBarry Smith   For nested preconditioners, such as `PCFIELDSPLIT` or `PCBJACOBI` this may not finish the construction of the preconditioner
10400b4b7b1cSBarry Smith   on the inner levels, the routine `PCSetUpOnBlocks()` may compute more of the preconditioner in those situations.
10410b4b7b1cSBarry Smith 
10420b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `KSPSetUp()`, `PCSetUpOnBlocks()`
10434b9ad928SBarry Smith @*/
1044d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp(PC pc)
1045d71ae5a4SJacob Faibussowitsch {
1046566e8bf2SBarry Smith   const char      *def;
1047fc9b51d3SBarry Smith   PetscObjectState matstate, matnonzerostate;
10484b9ad928SBarry Smith 
10494b9ad928SBarry Smith   PetscFunctionBegin;
10500700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
105128b400f6SJacob Faibussowitsch   PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first");
10526ce5e81cSLisandro Dalcin 
105323ee1639SBarry Smith   if (pc->setupcalled && pc->reusepreconditioner) {
10549566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n"));
10553ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10564b9ad928SBarry Smith   }
10574b9ad928SBarry Smith 
10589566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate));
10599566063dSJacob Faibussowitsch   PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate));
106023ee1639SBarry Smith   if (!pc->setupcalled) {
10615e62d202SMark Adams     //PetscCall(PetscInfo(pc, "Setting up PC for first time\n"));
106223ee1639SBarry Smith     pc->flag = DIFFERENT_NONZERO_PATTERN;
10639df67409SStefano Zampini   } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS);
10649df67409SStefano Zampini   else {
10659df67409SStefano Zampini     if (matnonzerostate != pc->matnonzerostate) {
10669566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n"));
106723ee1639SBarry Smith       pc->flag = DIFFERENT_NONZERO_PATTERN;
106823ee1639SBarry Smith     } else {
10695e62d202SMark Adams       //PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n"));
107023ee1639SBarry Smith       pc->flag = SAME_NONZERO_PATTERN;
107123ee1639SBarry Smith     }
107223ee1639SBarry Smith   }
107323ee1639SBarry Smith   pc->matstate        = matstate;
1074fc9b51d3SBarry Smith   pc->matnonzerostate = matnonzerostate;
107523ee1639SBarry Smith 
10767adad957SLisandro Dalcin   if (!((PetscObject)pc)->type_name) {
10779566063dSJacob Faibussowitsch     PetscCall(PCGetDefaultType_Private(pc, &def));
10789566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, def));
1079566e8bf2SBarry Smith   }
10804b9ad928SBarry Smith 
10819566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
10829566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure));
10839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0));
10844b9ad928SBarry Smith   if (pc->ops->setup) {
1085fe57bb1aSStefano Zampini     PetscCall(PCLogEventsDeactivatePush());
1086dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, setup);
1087fe57bb1aSStefano Zampini     PetscCall(PCLogEventsDeactivatePop());
10884b9ad928SBarry Smith   }
10899566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0));
1090422a814eSBarry Smith   if (!pc->setupcalled) pc->setupcalled = 1;
10913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10924b9ad928SBarry Smith }
10934b9ad928SBarry Smith 
10944b9ad928SBarry Smith /*@
10954b9ad928SBarry Smith   PCSetUpOnBlocks - Sets up the preconditioner for each block in
109673716367SStefano Zampini   the block Jacobi, overlapping Schwarz, and fieldsplit methods.
10974b9ad928SBarry Smith 
1098c3339decSBarry Smith   Collective
10994b9ad928SBarry Smith 
1100f1580f4eSBarry Smith   Input Parameter:
11010b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
11024b9ad928SBarry Smith 
11034b9ad928SBarry Smith   Level: developer
11044b9ad928SBarry Smith 
110573716367SStefano Zampini   Notes:
110673716367SStefano Zampini   For nested preconditioners such as `PCBJACOBI`, `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is
1107f1580f4eSBarry Smith   called on the outer `PC`, this routine ensures it is called.
1108f1580f4eSBarry Smith 
110973716367SStefano Zampini   It calls `PCSetUp()` if not yet called.
111073716367SStefano Zampini 
1111562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()`
11124b9ad928SBarry Smith @*/
1113d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUpOnBlocks(PC pc)
1114d71ae5a4SJacob Faibussowitsch {
11154b9ad928SBarry Smith   PetscFunctionBegin;
11160700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
111773716367SStefano Zampini   if (!pc->setupcalled) PetscCall(PCSetUp(pc)); /* "if" to prevent -info extra prints */
11183ba16761SJacob Faibussowitsch   if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS);
111921604f62SStefano Zampini   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
11209566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0));
1121fe57bb1aSStefano Zampini   PetscCall(PCLogEventsDeactivatePush());
1122dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, setuponblocks);
1123fe57bb1aSStefano Zampini   PetscCall(PCLogEventsDeactivatePop());
11249566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0));
11253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11264b9ad928SBarry Smith }
11274b9ad928SBarry Smith 
11284b9ad928SBarry Smith /*@C
11294b9ad928SBarry Smith   PCSetModifySubMatrices - Sets a user-defined routine for modifying the
113004c3f3b8SBarry Smith   submatrices that arise within certain subdomain-based preconditioners such as `PCASM`
11314b9ad928SBarry Smith 
1132c3339decSBarry Smith   Logically Collective
11334b9ad928SBarry Smith 
11344b9ad928SBarry Smith   Input Parameters:
11350b4b7b1cSBarry Smith + pc   - the `PC` preconditioner context
11364b9ad928SBarry Smith . func - routine for modifying the submatrices
113704c3f3b8SBarry Smith - ctx  - optional user-defined context (may be `NULL`)
11384b9ad928SBarry Smith 
113920f4b53cSBarry Smith   Calling sequence of `func`:
11400b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
114104c3f3b8SBarry Smith . nsub   - number of index sets
114220f4b53cSBarry Smith . row    - an array of index sets that contain the global row numbers
11434b9ad928SBarry Smith          that comprise each local submatrix
11444b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
11454b9ad928SBarry Smith          that comprise each local submatrix
11464b9ad928SBarry Smith . submat - array of local submatrices
11474b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
114804c3f3b8SBarry Smith          user-defined func routine (may be `NULL`)
11494b9ad928SBarry Smith 
115020f4b53cSBarry Smith   Level: advanced
115120f4b53cSBarry Smith 
11524b9ad928SBarry Smith   Notes:
115304c3f3b8SBarry Smith   The basic submatrices are extracted from the preconditioner matrix as
115404c3f3b8SBarry Smith   usual; the user can then alter these (for example, to set different boundary
115504c3f3b8SBarry Smith   conditions for each submatrix) before they are used for the local solves.
115604c3f3b8SBarry Smith 
1157f1580f4eSBarry Smith   `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and
1158f1580f4eSBarry Smith   `KSPSolve()`.
11594b9ad928SBarry Smith 
1160f1580f4eSBarry Smith   A routine set by `PCSetModifySubMatrices()` is currently called within
1161f1580f4eSBarry Smith   the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`)
11624b9ad928SBarry Smith   preconditioners.  All other preconditioners ignore this routine.
11634b9ad928SBarry Smith 
1164562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()`
11654b9ad928SBarry Smith @*/
116604c3f3b8SBarry Smith PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx), void *ctx)
1167d71ae5a4SJacob Faibussowitsch {
11684b9ad928SBarry Smith   PetscFunctionBegin;
11690700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11704b9ad928SBarry Smith   pc->modifysubmatrices  = func;
11714b9ad928SBarry Smith   pc->modifysubmatricesP = ctx;
11723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11734b9ad928SBarry Smith }
11744b9ad928SBarry Smith 
11754b9ad928SBarry Smith /*@C
11764b9ad928SBarry Smith   PCModifySubMatrices - Calls an optional user-defined routine within
1177f1580f4eSBarry Smith   certain preconditioners if one has been set with `PCSetModifySubMatrices()`.
11784b9ad928SBarry Smith 
1179c3339decSBarry Smith   Collective
11804b9ad928SBarry Smith 
11814b9ad928SBarry Smith   Input Parameters:
11820b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
11834b9ad928SBarry Smith . nsub   - the number of local submatrices
11844b9ad928SBarry Smith . row    - an array of index sets that contain the global row numbers
11854b9ad928SBarry Smith          that comprise each local submatrix
11864b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
11874b9ad928SBarry Smith          that comprise each local submatrix
11884b9ad928SBarry Smith . submat - array of local submatrices
11894b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
1190562efe2eSBarry Smith          user-defined routine (may be `NULL`)
11914b9ad928SBarry Smith 
11924b9ad928SBarry Smith   Output Parameter:
11934b9ad928SBarry Smith . submat - array of local submatrices (the entries of which may
11944b9ad928SBarry Smith             have been modified)
11954b9ad928SBarry Smith 
119620f4b53cSBarry Smith   Level: developer
119720f4b53cSBarry Smith 
119804c3f3b8SBarry Smith   Note:
11994b9ad928SBarry Smith   The user should NOT generally call this routine, as it will
120004c3f3b8SBarry Smith   automatically be called within certain preconditioners.
12014b9ad928SBarry Smith 
1202562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetModifySubMatrices()`
12034b9ad928SBarry Smith @*/
1204d71ae5a4SJacob Faibussowitsch PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx)
1205d71ae5a4SJacob Faibussowitsch {
12064b9ad928SBarry Smith   PetscFunctionBegin;
12070700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12083ba16761SJacob Faibussowitsch   if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS);
12099566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0));
12109566063dSJacob Faibussowitsch   PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx));
12119566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0));
12123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12134b9ad928SBarry Smith }
12144b9ad928SBarry Smith 
12154b9ad928SBarry Smith /*@
12164b9ad928SBarry Smith   PCSetOperators - Sets the matrix associated with the linear system and
12174b9ad928SBarry Smith   a (possibly) different one associated with the preconditioner.
12184b9ad928SBarry Smith 
1219c3339decSBarry Smith   Logically Collective
12204b9ad928SBarry Smith 
12214b9ad928SBarry Smith   Input Parameters:
12220b4b7b1cSBarry Smith + pc   - the `PC` preconditioner context
1223e5d3d808SBarry Smith . Amat - the matrix that defines the linear system
1224d1e9a80fSBarry Smith - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
12254b9ad928SBarry Smith 
122620f4b53cSBarry Smith   Level: intermediate
1227189c0b8aSBarry Smith 
122820f4b53cSBarry Smith   Notes:
122920f4b53cSBarry Smith   Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used.
123020f4b53cSBarry Smith 
123120f4b53cSBarry Smith   If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then
1232f1580f4eSBarry Smith   first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()`
1233f1580f4eSBarry Smith   on it and then pass it back in in your call to `KSPSetOperators()`.
1234189c0b8aSBarry Smith 
12354b9ad928SBarry Smith   More Notes about Repeated Solution of Linear Systems:
123620f4b53cSBarry Smith   PETSc does NOT reset the matrix entries of either `Amat` or `Pmat`
12374b9ad928SBarry Smith   to zero after a linear solve; the user is completely responsible for
1238f1580f4eSBarry Smith   matrix assembly.  See the routine `MatZeroEntries()` if desiring to
12394b9ad928SBarry Smith   zero all elements of a matrix.
12404b9ad928SBarry Smith 
1241562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`
12424b9ad928SBarry Smith  @*/
1243d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat)
1244d71ae5a4SJacob Faibussowitsch {
12453b3f6333SBarry Smith   PetscInt m1, n1, m2, n2;
12464b9ad928SBarry Smith 
12474b9ad928SBarry Smith   PetscFunctionBegin;
12480700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12490700a824SBarry Smith   if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
12500700a824SBarry Smith   if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3);
12513fc8bf9cSBarry Smith   if (Amat) PetscCheckSameComm(pc, 1, Amat, 2);
12523fc8bf9cSBarry Smith   if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3);
125331641f1aSBarry Smith   if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) {
12549566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Amat, &m1, &n1));
12559566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->mat, &m2, &n2));
125663a3b9bcSJacob 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);
12579566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Pmat, &m1, &n1));
12589566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2));
125963a3b9bcSJacob 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);
12603b3f6333SBarry Smith   }
12614b9ad928SBarry Smith 
126223ee1639SBarry Smith   if (Pmat != pc->pmat) {
126323ee1639SBarry Smith     /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
126423ee1639SBarry Smith     pc->matnonzerostate = -1;
126523ee1639SBarry Smith     pc->matstate        = -1;
126623ee1639SBarry Smith   }
126723ee1639SBarry Smith 
1268906ed7ccSBarry Smith   /* reference first in case the matrices are the same */
12699566063dSJacob Faibussowitsch   if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat));
12709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->mat));
12719566063dSJacob Faibussowitsch   if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat));
12729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->pmat));
12734b9ad928SBarry Smith   pc->mat  = Amat;
12744b9ad928SBarry Smith   pc->pmat = Pmat;
12753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1276e56f5c9eSBarry Smith }
1277e56f5c9eSBarry Smith 
127823ee1639SBarry Smith /*@
12790b4b7b1cSBarry Smith   PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner `PC` has changed.
128023ee1639SBarry Smith 
1281c3339decSBarry Smith   Logically Collective
128223ee1639SBarry Smith 
128323ee1639SBarry Smith   Input Parameters:
12840b4b7b1cSBarry Smith + pc   - the `PC` preconditioner context
1285f1580f4eSBarry Smith - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
128623ee1639SBarry Smith 
12872b26979fSBarry Smith   Level: intermediate
12882b26979fSBarry Smith 
1289f1580f4eSBarry Smith   Note:
1290f1580f4eSBarry Smith   Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option
1291f1580f4eSBarry Smith   prevents this.
1292f1580f4eSBarry Smith 
1293562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()`
129423ee1639SBarry Smith  @*/
1295d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag)
1296d71ae5a4SJacob Faibussowitsch {
129723ee1639SBarry Smith   PetscFunctionBegin;
129823ee1639SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1299f9177522SStefano Zampini   PetscValidLogicalCollectiveBool(pc, flag, 2);
130023ee1639SBarry Smith   pc->reusepreconditioner = flag;
13013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13024b9ad928SBarry Smith }
13034b9ad928SBarry Smith 
1304c60c7ad4SBarry Smith /*@
1305f1580f4eSBarry Smith   PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed.
1306c60c7ad4SBarry Smith 
1307bba28a21SBarry Smith   Not Collective
1308c60c7ad4SBarry Smith 
1309c60c7ad4SBarry Smith   Input Parameter:
13100b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
1311c60c7ad4SBarry Smith 
1312c60c7ad4SBarry Smith   Output Parameter:
1313f1580f4eSBarry Smith . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
1314c60c7ad4SBarry Smith 
1315d0418729SSatish Balay   Level: intermediate
1316d0418729SSatish Balay 
1317562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()`
1318c60c7ad4SBarry Smith  @*/
1319d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag)
1320d71ae5a4SJacob Faibussowitsch {
1321c60c7ad4SBarry Smith   PetscFunctionBegin;
1322c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
13234f572ea9SToby Isaac   PetscAssertPointer(flag, 2);
1324c60c7ad4SBarry Smith   *flag = pc->reusepreconditioner;
13253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1326c60c7ad4SBarry Smith }
1327c60c7ad4SBarry Smith 
1328487a658cSBarry Smith /*@
13294b9ad928SBarry Smith   PCGetOperators - Gets the matrix associated with the linear system and
13300b4b7b1cSBarry Smith   possibly a different one which is used to construct the preconditioner.
13314b9ad928SBarry Smith 
1332562efe2eSBarry Smith   Not Collective, though parallel `Mat`s are returned if `pc` is parallel
13334b9ad928SBarry Smith 
13344b9ad928SBarry Smith   Input Parameter:
13350b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
13364b9ad928SBarry Smith 
13374b9ad928SBarry Smith   Output Parameters:
1338e5d3d808SBarry Smith + Amat - the matrix defining the linear system
133923ee1639SBarry Smith - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
13404b9ad928SBarry Smith 
13414b9ad928SBarry Smith   Level: intermediate
13424b9ad928SBarry Smith 
1343f1580f4eSBarry Smith   Note:
134495452b02SPatrick Sanan   Does not increase the reference count of the matrices, so you should not destroy them
1345298cc208SBarry Smith 
1346f1580f4eSBarry Smith   Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators
1347f1580f4eSBarry Smith   are created in `PC` and returned to the user. In this case, if both operators
134873950996SBarry Smith   mat and pmat are requested, two DIFFERENT operators will be returned. If
134973950996SBarry Smith   only one is requested both operators in the PC will be the same (i.e. as
1350f1580f4eSBarry Smith   if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats).
135173950996SBarry Smith   The user must set the sizes of the returned matrices and their type etc just
1352f1580f4eSBarry Smith   as if the user created them with `MatCreate()`. For example,
135373950996SBarry Smith 
1354f1580f4eSBarry Smith .vb
1355f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1356f1580f4eSBarry Smith            set size, type, etc of Amat
135773950996SBarry Smith 
1358f1580f4eSBarry Smith          MatCreate(comm,&mat);
1359f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Amat);
1360f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)mat);
1361f1580f4eSBarry Smith            set size, type, etc of Amat
1362f1580f4eSBarry Smith .ve
136373950996SBarry Smith 
136473950996SBarry Smith   and
136573950996SBarry Smith 
1366f1580f4eSBarry Smith .vb
1367f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1368f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
136973950996SBarry Smith 
1370f1580f4eSBarry Smith          MatCreate(comm,&Amat);
1371f1580f4eSBarry Smith          MatCreate(comm,&Pmat);
1372f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1373f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Amat);
1374f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Pmat);
1375f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
1376f1580f4eSBarry Smith .ve
137773950996SBarry Smith 
1378f1580f4eSBarry Smith   The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy
1379b8d709abSRichard Tran Mills   of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely
1380f1580f4eSBarry Smith   managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look
1381f1580f4eSBarry Smith   at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to
1382f1580f4eSBarry Smith   the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP
1383f1580f4eSBarry Smith   you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you).
1384f1580f4eSBarry Smith   Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when
138573950996SBarry Smith   it can be created for you?
138673950996SBarry Smith 
1387562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()`
13884b9ad928SBarry Smith @*/
1389d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat)
1390d71ae5a4SJacob Faibussowitsch {
13914b9ad928SBarry Smith   PetscFunctionBegin;
13920700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1393e5d3d808SBarry Smith   if (Amat) {
139473950996SBarry Smith     if (!pc->mat) {
13959fca8c71SStefano Zampini       if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */
13969a4708feSJed Brown         pc->mat = pc->pmat;
13979566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->mat));
1398e5d3d808SBarry Smith       } else { /* both Amat and Pmat are empty */
13999566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat));
1400e5d3d808SBarry Smith         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
140173950996SBarry Smith           pc->pmat = pc->mat;
14029566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->pmat));
140373950996SBarry Smith         }
140473950996SBarry Smith       }
14059a4708feSJed Brown     }
1406e5d3d808SBarry Smith     *Amat = pc->mat;
140773950996SBarry Smith   }
1408e5d3d808SBarry Smith   if (Pmat) {
140973950996SBarry Smith     if (!pc->pmat) {
1410e5d3d808SBarry Smith       if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */
14119a4708feSJed Brown         pc->pmat = pc->mat;
14129566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->pmat));
14139a4708feSJed Brown       } else {
14149566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat));
1415e5d3d808SBarry Smith         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
141673950996SBarry Smith           pc->mat = pc->pmat;
14179566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->mat));
141873950996SBarry Smith         }
141973950996SBarry Smith       }
14209a4708feSJed Brown     }
1421e5d3d808SBarry Smith     *Pmat = pc->pmat;
142273950996SBarry Smith   }
14233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14244b9ad928SBarry Smith }
14254b9ad928SBarry Smith 
14265d83a8b1SBarry Smith /*@
1427906ed7ccSBarry Smith   PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1428f1580f4eSBarry Smith   possibly a different one associated with the preconditioner have been set in the `PC`.
1429906ed7ccSBarry Smith 
143020f4b53cSBarry Smith   Not Collective, though the results on all processes should be the same
1431906ed7ccSBarry Smith 
1432906ed7ccSBarry Smith   Input Parameter:
14330b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
1434906ed7ccSBarry Smith 
1435906ed7ccSBarry Smith   Output Parameters:
1436906ed7ccSBarry Smith + mat  - the matrix associated with the linear system was set
1437906ed7ccSBarry Smith - pmat - matrix associated with the preconditioner was set, usually the same
1438906ed7ccSBarry Smith 
1439906ed7ccSBarry Smith   Level: intermediate
1440906ed7ccSBarry Smith 
1441562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`
1442906ed7ccSBarry Smith @*/
1443d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat)
1444d71ae5a4SJacob Faibussowitsch {
1445906ed7ccSBarry Smith   PetscFunctionBegin;
14460700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1447906ed7ccSBarry Smith   if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1448906ed7ccSBarry Smith   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
14493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1450906ed7ccSBarry Smith }
1451906ed7ccSBarry Smith 
1452f39d8e23SSatish Balay /*@
1453a4fd02acSBarry Smith   PCFactorGetMatrix - Gets the factored matrix from the
1454f1580f4eSBarry Smith   preconditioner context.  This routine is valid only for the `PCLU`,
1455f1580f4eSBarry Smith   `PCILU`, `PCCHOLESKY`, and `PCICC` methods.
14564b9ad928SBarry Smith 
1457562efe2eSBarry Smith   Not Collective though `mat` is parallel if `pc` is parallel
14584b9ad928SBarry Smith 
1459f1580f4eSBarry Smith   Input Parameter:
14600b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
14614b9ad928SBarry Smith 
1462feefa0e1SJacob Faibussowitsch   Output Parameters:
14634b9ad928SBarry Smith . mat - the factored matrix
14644b9ad928SBarry Smith 
14654b9ad928SBarry Smith   Level: advanced
14664b9ad928SBarry Smith 
1467f1580f4eSBarry Smith   Note:
1468562efe2eSBarry Smith   Does not increase the reference count for `mat` so DO NOT destroy it
14699405f653SBarry Smith 
1470562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC`
14714b9ad928SBarry Smith @*/
1472d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat)
1473d71ae5a4SJacob Faibussowitsch {
14744b9ad928SBarry Smith   PetscFunctionBegin;
14750700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14764f572ea9SToby Isaac   PetscAssertPointer(mat, 2);
14777def7855SStefano Zampini   PetscCall(PCFactorSetUpMatSolverType(pc));
1478dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, getfactoredmatrix, mat);
14793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14804b9ad928SBarry Smith }
14814b9ad928SBarry Smith 
1482cc4c1da9SBarry Smith /*@
14834b9ad928SBarry Smith   PCSetOptionsPrefix - Sets the prefix used for searching for all
1484f1580f4eSBarry Smith   `PC` options in the database.
14854b9ad928SBarry Smith 
1486c3339decSBarry Smith   Logically Collective
14874b9ad928SBarry Smith 
14884b9ad928SBarry Smith   Input Parameters:
14890b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
1490f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
14914b9ad928SBarry Smith 
1492f1580f4eSBarry Smith   Note:
14934b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
14944b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
14954b9ad928SBarry Smith   hyphen.
14964b9ad928SBarry Smith 
14974b9ad928SBarry Smith   Level: advanced
14984b9ad928SBarry Smith 
1499562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()`
15004b9ad928SBarry Smith @*/
1501d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[])
1502d71ae5a4SJacob Faibussowitsch {
15034b9ad928SBarry Smith   PetscFunctionBegin;
15040700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15059566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix));
15063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15074b9ad928SBarry Smith }
15084b9ad928SBarry Smith 
1509cc4c1da9SBarry Smith /*@
15104b9ad928SBarry Smith   PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1511f1580f4eSBarry Smith   `PC` options in the database.
15124b9ad928SBarry Smith 
1513c3339decSBarry Smith   Logically Collective
15144b9ad928SBarry Smith 
15154b9ad928SBarry Smith   Input Parameters:
15160b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
1517f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
15184b9ad928SBarry Smith 
1519f1580f4eSBarry Smith   Note:
15204b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
15214b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
15224b9ad928SBarry Smith   hyphen.
15234b9ad928SBarry Smith 
15244b9ad928SBarry Smith   Level: advanced
15254b9ad928SBarry Smith 
1526562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()`
15274b9ad928SBarry Smith @*/
1528d71ae5a4SJacob Faibussowitsch PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[])
1529d71ae5a4SJacob Faibussowitsch {
15304b9ad928SBarry Smith   PetscFunctionBegin;
15310700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15329566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix));
15333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15344b9ad928SBarry Smith }
15354b9ad928SBarry Smith 
1536cc4c1da9SBarry Smith /*@
15374b9ad928SBarry Smith   PCGetOptionsPrefix - Gets the prefix used for searching for all
15384b9ad928SBarry Smith   PC options in the database.
15394b9ad928SBarry Smith 
15404b9ad928SBarry Smith   Not Collective
15414b9ad928SBarry Smith 
1542f1580f4eSBarry Smith   Input Parameter:
15430b4b7b1cSBarry Smith . pc - the `PC` preconditioner context
15444b9ad928SBarry Smith 
1545f1580f4eSBarry Smith   Output Parameter:
15464b9ad928SBarry Smith . prefix - pointer to the prefix string used, is returned
15474b9ad928SBarry Smith 
15484b9ad928SBarry Smith   Level: advanced
15494b9ad928SBarry Smith 
1550562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()`
15514b9ad928SBarry Smith @*/
1552d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[])
1553d71ae5a4SJacob Faibussowitsch {
15544b9ad928SBarry Smith   PetscFunctionBegin;
15550700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15564f572ea9SToby Isaac   PetscAssertPointer(prefix, 2);
15579566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix));
15583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15594b9ad928SBarry Smith }
15604b9ad928SBarry Smith 
15618066bbecSBarry Smith /*
1562dd8e379bSPierre Jolivet    Indicates the right-hand side will be changed by KSPSolve(), this occurs for a few
15638066bbecSBarry Smith   preconditioners including BDDC and Eisentat that transform the equations before applying
15648066bbecSBarry Smith   the Krylov methods
15658066bbecSBarry Smith */
1566d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change)
1567d71ae5a4SJacob Faibussowitsch {
1568a06fd7f2SStefano Zampini   PetscFunctionBegin;
1569a06fd7f2SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15704f572ea9SToby Isaac   PetscAssertPointer(change, 2);
1571a06fd7f2SStefano Zampini   *change = PETSC_FALSE;
1572cac4c232SBarry Smith   PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change));
15733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1574a06fd7f2SStefano Zampini }
1575a06fd7f2SStefano Zampini 
15764b9ad928SBarry Smith /*@
15770b4b7b1cSBarry Smith   PCPreSolve - Optional pre-solve phase, intended for any preconditioner-specific actions that must be performed before
15780b4b7b1cSBarry Smith   the iterative solve itself. Used in conjunction with `PCPostSolve()`
15794b9ad928SBarry Smith 
1580c3339decSBarry Smith   Collective
15814b9ad928SBarry Smith 
15824b9ad928SBarry Smith   Input Parameters:
15830b4b7b1cSBarry Smith + pc  - the `PC` preconditioner context
15844b9ad928SBarry Smith - ksp - the Krylov subspace context
15854b9ad928SBarry Smith 
15864b9ad928SBarry Smith   Level: developer
15874b9ad928SBarry Smith 
15884b9ad928SBarry Smith   Notes:
1589f1580f4eSBarry Smith   `KSPSolve()` calls this directly, so is rarely called by the user.
15904b9ad928SBarry Smith 
15910b4b7b1cSBarry Smith   Certain preconditioners, such as the `PCType` of `PCEISENSTAT`, change the formulation of the linear system to be solved iteratively.
15920b4b7b1cSBarry Smith   This function performs that transformation. `PCPostSolve()` then transforms the system back to its original form after the solve.
15930b4b7b1cSBarry Smith   `PCPostSolve()` also transforms the resulting solution of the transformed system to the solution of the original problem.
15940b4b7b1cSBarry Smith 
159554c05997SPierre Jolivet   `KSPSetPreSolve()` and `KSPSetPostSolve()` provide an alternative way to provide such transformations.
15960b4b7b1cSBarry Smith 
15970b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCPostSolve()`, `KSP`, `PCSetPreSolve()`, `KSPSetPreSolve()`, `KSPSetPostSolve()`
15984b9ad928SBarry Smith @*/
1599d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPreSolve(PC pc, KSP ksp)
1600d71ae5a4SJacob Faibussowitsch {
16014b9ad928SBarry Smith   Vec x, rhs;
16024b9ad928SBarry Smith 
16034b9ad928SBarry Smith   PetscFunctionBegin;
16040700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16050700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1606d9a03883SBarry Smith   pc->presolvedone++;
16077827d75bSBarry Smith   PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice");
16089566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
16099566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
16104b9ad928SBarry Smith 
1611dbbe0bcdSBarry Smith   if (pc->ops->presolve) PetscUseTypeMethod(pc, presolve, ksp, rhs, x);
1612f4f49eeaSPierre Jolivet   else if (pc->presolve) PetscCall(pc->presolve(pc, ksp));
16133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16144b9ad928SBarry Smith }
16154b9ad928SBarry Smith 
1616f560b561SHong Zhang /*@C
1617f1580f4eSBarry Smith   PCSetPreSolve - Sets function used by `PCPreSolve()` which is intended for any
1618f560b561SHong Zhang   preconditioner-specific actions that must be performed before
1619f560b561SHong Zhang   the iterative solve itself.
1620f560b561SHong Zhang 
1621c3339decSBarry Smith   Logically Collective
1622f560b561SHong Zhang 
1623f560b561SHong Zhang   Input Parameters:
1624f560b561SHong Zhang + pc       - the preconditioner object
1625f560b561SHong Zhang - presolve - the function to call before the solve
1626f560b561SHong Zhang 
162720f4b53cSBarry Smith   Calling sequence of `presolve`:
162820f4b53cSBarry Smith + pc  - the `PC` context
162920f4b53cSBarry Smith - ksp - the `KSP` context
1630f560b561SHong Zhang 
1631f560b561SHong Zhang   Level: developer
1632f560b561SHong Zhang 
1633562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCPreSolve()`
1634f560b561SHong Zhang @*/
163504c3f3b8SBarry Smith PetscErrorCode PCSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC pc, KSP ksp))
1636d71ae5a4SJacob Faibussowitsch {
1637f560b561SHong Zhang   PetscFunctionBegin;
1638f560b561SHong Zhang   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1639f560b561SHong Zhang   pc->presolve = presolve;
16403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1641f560b561SHong Zhang }
1642f560b561SHong Zhang 
16434b9ad928SBarry Smith /*@
16444b9ad928SBarry Smith   PCPostSolve - Optional post-solve phase, intended for any
16454b9ad928SBarry Smith   preconditioner-specific actions that must be performed after
16464b9ad928SBarry Smith   the iterative solve itself.
16474b9ad928SBarry Smith 
1648c3339decSBarry Smith   Collective
16494b9ad928SBarry Smith 
16504b9ad928SBarry Smith   Input Parameters:
16510b4b7b1cSBarry Smith + pc  - the `PC` preconditioner context
16520b4b7b1cSBarry Smith - ksp - the `KSP` Krylov subspace context
16534b9ad928SBarry Smith 
1654feefa0e1SJacob Faibussowitsch   Example Usage:
16554b9ad928SBarry Smith .vb
16564b9ad928SBarry Smith     PCPreSolve(pc,ksp);
165723ce1328SBarry Smith     KSPSolve(ksp,b,x);
16584b9ad928SBarry Smith     PCPostSolve(pc,ksp);
16594b9ad928SBarry Smith .ve
16604b9ad928SBarry Smith 
1661562efe2eSBarry Smith   Level: developer
1662562efe2eSBarry Smith 
16634b9ad928SBarry Smith   Note:
1664f1580f4eSBarry Smith   `KSPSolve()` calls this routine directly, so it is rarely called by the user.
16654b9ad928SBarry Smith 
166654c05997SPierre Jolivet .seealso: [](ch_ksp), `PC`, `PCSetPreSolve()`, `KSPSetPostSolve()`, `KSPSetPreSolve()`, `PCPreSolve()`, `KSPSolve()`
16674b9ad928SBarry Smith @*/
1668d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPostSolve(PC pc, KSP ksp)
1669d71ae5a4SJacob Faibussowitsch {
16704b9ad928SBarry Smith   Vec x, rhs;
16714b9ad928SBarry Smith 
16724b9ad928SBarry Smith   PetscFunctionBegin;
16730700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16740700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1675d9a03883SBarry Smith   pc->presolvedone--;
16769566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
16779566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
1678dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, postsolve, ksp, rhs, x);
16793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16804b9ad928SBarry Smith }
16814b9ad928SBarry Smith 
1682ffeef943SBarry Smith /*@
1683f1580f4eSBarry Smith   PCLoad - Loads a `PC` that has been stored in binary  with `PCView()`.
168455849f57SBarry Smith 
1685c3339decSBarry Smith   Collective
168655849f57SBarry Smith 
168755849f57SBarry Smith   Input Parameters:
1688f1580f4eSBarry Smith + newdm  - the newly loaded `PC`, this needs to have been created with `PCCreate()` or
1689f1580f4eSBarry Smith            some related function before a call to `PCLoad()`.
16900b4b7b1cSBarry Smith - viewer - binary file viewer `PETSCVIEWERBINARY`, obtained from `PetscViewerBinaryOpen()`
169155849f57SBarry Smith 
169255849f57SBarry Smith   Level: intermediate
169355849f57SBarry Smith 
1694f1580f4eSBarry Smith   Note:
1695562efe2eSBarry Smith   The type is determined by the data in the file, any `PCType` set into the `PC` before this call is ignored.
169655849f57SBarry Smith 
16970b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()`, `PETSCVIEWERBINARY`
169855849f57SBarry Smith @*/
1699d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLoad(PC newdm, PetscViewer viewer)
1700d71ae5a4SJacob Faibussowitsch {
170155849f57SBarry Smith   PetscBool isbinary;
1702060da220SMatthew G. Knepley   PetscInt  classid;
170355849f57SBarry Smith   char      type[256];
170455849f57SBarry Smith 
170555849f57SBarry Smith   PetscFunctionBegin;
170655849f57SBarry Smith   PetscValidHeaderSpecific(newdm, PC_CLASSID, 1);
170755849f57SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
17089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
170928b400f6SJacob Faibussowitsch   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
171055849f57SBarry Smith 
17119566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
17127827d75bSBarry Smith   PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file");
17139566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
17149566063dSJacob Faibussowitsch   PetscCall(PCSetType(newdm, type));
1715dbbe0bcdSBarry Smith   PetscTryTypeMethod(newdm, load, viewer);
17163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
171755849f57SBarry Smith }
171855849f57SBarry Smith 
17199804daf3SBarry Smith #include <petscdraw.h>
1720e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1721e04113cfSBarry Smith   #include <petscviewersaws.h>
17220acecf5bSBarry Smith #endif
1723fe2efc57SMark 
1724ffeef943SBarry Smith /*@
17250b4b7b1cSBarry Smith   PCViewFromOptions - View (print or provide information about) the `PC`, based on options in the options database
1726fe2efc57SMark 
1727c3339decSBarry Smith   Collective
1728fe2efc57SMark 
1729fe2efc57SMark   Input Parameters:
173020f4b53cSBarry Smith + A    - the `PC` context
1731f1580f4eSBarry Smith . obj  - Optional object that provides the options prefix
17320b4b7b1cSBarry Smith - name - command line option name
1733fe2efc57SMark 
17340b4b7b1cSBarry Smith   Level: developer
1735f1580f4eSBarry Smith 
1736562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()`
1737fe2efc57SMark @*/
1738d71ae5a4SJacob Faibussowitsch PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[])
1739d71ae5a4SJacob Faibussowitsch {
1740fe2efc57SMark   PetscFunctionBegin;
1741fe2efc57SMark   PetscValidHeaderSpecific(A, PC_CLASSID, 1);
17429566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
17433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1744fe2efc57SMark }
1745fe2efc57SMark 
1746ffeef943SBarry Smith /*@
1747f1580f4eSBarry Smith   PCView - Prints information about the `PC`
17484b9ad928SBarry Smith 
1749c3339decSBarry Smith   Collective
17504b9ad928SBarry Smith 
17514b9ad928SBarry Smith   Input Parameters:
17520b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
17530b4b7b1cSBarry Smith - viewer - optional `PetscViewer` visualization context
17544b9ad928SBarry Smith 
17550b4b7b1cSBarry Smith   Level: intermediate
175620f4b53cSBarry Smith 
1757f1580f4eSBarry Smith   Notes:
17584b9ad928SBarry Smith   The available visualization contexts include
1759f1580f4eSBarry Smith +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
1760f1580f4eSBarry Smith -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
17614b9ad928SBarry Smith   output where only the first processor opens
17624b9ad928SBarry Smith   the file. All other processors send their
17634b9ad928SBarry Smith   data to the first processor to print.
17644b9ad928SBarry Smith 
17654b9ad928SBarry Smith   The user can open an alternative visualization contexts with
1766f1580f4eSBarry Smith   `PetscViewerASCIIOpen()` (output to a specified file).
17674b9ad928SBarry Smith 
17680b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewer`, `PetscViewerType`, `KSPView()`, `PetscViewerASCIIOpen()`
17694b9ad928SBarry Smith @*/
1770d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView(PC pc, PetscViewer viewer)
1771d71ae5a4SJacob Faibussowitsch {
177219fd82e9SBarry Smith   PCType            cstr;
17736cd81132SPierre Jolivet   PetscViewerFormat format;
17746cd81132SPierre Jolivet   PetscBool         iascii, isstring, isbinary, isdraw, pop = PETSC_FALSE;
1775e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1776536b137fSBarry Smith   PetscBool issaws;
17770acecf5bSBarry Smith #endif
17784b9ad928SBarry Smith 
17794b9ad928SBarry Smith   PetscFunctionBegin;
17800700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
178148a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
17820700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1783c9780b6fSBarry Smith   PetscCheckSameComm(pc, 1, viewer, 2);
17844b9ad928SBarry Smith 
17859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
17869566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
17879566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
17889566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1789e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
17909566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
17910acecf5bSBarry Smith #endif
1792219b1045SBarry Smith 
179332077d6dSBarry Smith   if (iascii) {
17949566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer));
179548a46eb9SPierre Jolivet     if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, "  PC has not been set up so information may be incomplete\n"));
17969566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1797dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
17989566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
1799834dbeb0SBarry Smith     if (pc->mat) {
18006cd81132SPierre Jolivet       PetscCall(PetscViewerGetFormat(viewer, &format));
18016cd81132SPierre Jolivet       if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
18029566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
18036cd81132SPierre Jolivet         pop = PETSC_TRUE;
18046cd81132SPierre Jolivet       }
18054b9ad928SBarry Smith       if (pc->pmat == pc->mat) {
18069566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix = precond matrix:\n"));
18079566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
18089566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
18099566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
18104b9ad928SBarry Smith       } else {
1811834dbeb0SBarry Smith         if (pc->pmat) {
18129566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix followed by preconditioner matrix:\n"));
18134b9ad928SBarry Smith         } else {
18149566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix:\n"));
18154b9ad928SBarry Smith         }
18169566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
18179566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
18189566063dSJacob Faibussowitsch         if (pc->pmat) PetscCall(MatView(pc->pmat, viewer));
18199566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
18204b9ad928SBarry Smith       }
18216cd81132SPierre Jolivet       if (pop) PetscCall(PetscViewerPopFormat(viewer));
18224b9ad928SBarry Smith     }
18234b9ad928SBarry Smith   } else if (isstring) {
18249566063dSJacob Faibussowitsch     PetscCall(PCGetType(pc, &cstr));
18259566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr));
1826dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18279566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18289566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18295b0b0462SBarry Smith   } else if (isbinary) {
183055849f57SBarry Smith     PetscInt    classid = PC_FILE_CLASSID;
183155849f57SBarry Smith     MPI_Comm    comm;
183255849f57SBarry Smith     PetscMPIInt rank;
183355849f57SBarry Smith     char        type[256];
183455849f57SBarry Smith 
18359566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
18369566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1837dd400576SPatrick Sanan     if (rank == 0) {
18389566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
18399566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256));
18409566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
184155849f57SBarry Smith     }
1842dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
1843219b1045SBarry Smith   } else if (isdraw) {
1844219b1045SBarry Smith     PetscDraw draw;
1845d9884438SBarry Smith     char      str[25];
184689fd9fafSBarry Smith     PetscReal x, y, bottom, h;
1847d9884438SBarry Smith     PetscInt  n;
1848219b1045SBarry Smith 
18499566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
18509566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
18511d840656SPeter Brune     if (pc->mat) {
18529566063dSJacob Faibussowitsch       PetscCall(MatGetSize(pc->mat, &n, NULL));
185363a3b9bcSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n));
18541d840656SPeter Brune     } else {
18559566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name));
18561d840656SPeter Brune     }
18579566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
185889fd9fafSBarry Smith     bottom = y - h;
18599566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
1860dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18619566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
1862e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1863536b137fSBarry Smith   } else if (issaws) {
1864d45a07a7SBarry Smith     PetscMPIInt rank;
1865d45a07a7SBarry Smith 
18669566063dSJacob Faibussowitsch     PetscCall(PetscObjectName((PetscObject)pc));
18679566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
186848a46eb9SPierre Jolivet     if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer));
18699566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18709566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18710acecf5bSBarry Smith #endif
18724b9ad928SBarry Smith   }
18733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18744b9ad928SBarry Smith }
18754b9ad928SBarry Smith 
18764b9ad928SBarry Smith /*@C
18770b4b7b1cSBarry Smith   PCRegister -  Adds a method (`PCType`) to the PETSc preconditioner package.
18781c84c290SBarry Smith 
1879cc4c1da9SBarry Smith   Not collective. No Fortran Support
18801c84c290SBarry Smith 
18811c84c290SBarry Smith   Input Parameters:
188220f4b53cSBarry Smith + sname    - name of a new user-defined solver
18830b4b7b1cSBarry Smith - function - routine to create the method context which will be stored in a `PC` when `PCSetType()` is called
18841c84c290SBarry Smith 
1885feefa0e1SJacob Faibussowitsch   Example Usage:
18861c84c290SBarry Smith .vb
1887bdf89e91SBarry Smith    PCRegister("my_solver", MySolverCreate);
18881c84c290SBarry Smith .ve
18891c84c290SBarry Smith 
18901c84c290SBarry Smith   Then, your solver can be chosen with the procedural interface via
18911c84c290SBarry Smith $     PCSetType(pc, "my_solver")
18921c84c290SBarry Smith   or at runtime via the option
18931c84c290SBarry Smith $     -pc_type my_solver
18944b9ad928SBarry Smith 
18954b9ad928SBarry Smith   Level: advanced
18961c84c290SBarry Smith 
189720f4b53cSBarry Smith   Note:
18980b4b7b1cSBarry Smith   A simpler alternative to using `PCRegister()` for an application specific preconditioner is to use a `PC` of `PCType` `PCSHELL` and
18990b4b7b1cSBarry Smith   provide your customizations with `PCShellSetContext()` and `PCShellSetApply()`
19000b4b7b1cSBarry Smith 
190120f4b53cSBarry Smith   `PCRegister()` may be called multiple times to add several user-defined preconditioners.
190220f4b53cSBarry Smith 
19030b4b7b1cSBarry Smith .seealso: [](ch_ksp), `PC`, `PCType`, `PCRegisterAll()`, `PCSetType()`, `PCShellSetContext()`, `PCShellSetApply()`, `PCSHELL`
19044b9ad928SBarry Smith @*/
1905d71ae5a4SJacob Faibussowitsch PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC))
1906d71ae5a4SJacob Faibussowitsch {
19074b9ad928SBarry Smith   PetscFunctionBegin;
19089566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
19099566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCList, sname, function));
19103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19114b9ad928SBarry Smith }
19124b9ad928SBarry Smith 
1913d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y)
1914d71ae5a4SJacob Faibussowitsch {
1915186a3e20SStefano Zampini   PC pc;
1916186a3e20SStefano Zampini 
1917186a3e20SStefano Zampini   PetscFunctionBegin;
19189566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &pc));
19199566063dSJacob Faibussowitsch   PetscCall(PCApply(pc, X, Y));
19203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1921186a3e20SStefano Zampini }
1922186a3e20SStefano Zampini 
19235d83a8b1SBarry Smith /*@
19240b4b7b1cSBarry Smith   PCComputeOperator - Computes the explicit preconditioned operator as a matrix `Mat`.
19254b9ad928SBarry Smith 
1926c3339decSBarry Smith   Collective
19274b9ad928SBarry Smith 
1928d8d19677SJose E. Roman   Input Parameters:
19290b4b7b1cSBarry Smith + pc      - the `PC` preconditioner object
1930562efe2eSBarry Smith - mattype - the `MatType` to be used for the operator
19314b9ad928SBarry Smith 
19324b9ad928SBarry Smith   Output Parameter:
1933a5b23f4aSJose E. Roman . mat - the explicit preconditioned operator
19344b9ad928SBarry Smith 
193520f4b53cSBarry Smith   Level: advanced
193620f4b53cSBarry Smith 
1937f1580f4eSBarry Smith   Note:
1938186a3e20SStefano Zampini   This computation is done by applying the operators to columns of the identity matrix.
1939186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
1940562efe2eSBarry Smith   Currently, this routine uses a dense matrix format when `mattype` == `NULL`
19414b9ad928SBarry Smith 
19420b4b7b1cSBarry Smith   Developer Note:
19430b4b7b1cSBarry Smith   This should be called `PCCreateExplicitOperator()`
19440b4b7b1cSBarry Smith 
1945562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPComputeOperator()`, `MatType`
19464b9ad928SBarry Smith @*/
1947d71ae5a4SJacob Faibussowitsch PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat)
1948d71ae5a4SJacob Faibussowitsch {
1949186a3e20SStefano Zampini   PetscInt N, M, m, n;
1950186a3e20SStefano Zampini   Mat      A, Apc;
19514b9ad928SBarry Smith 
19524b9ad928SBarry Smith   PetscFunctionBegin;
19530700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
19544f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
19559566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, &A, NULL));
19569566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
19579566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
19589566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc));
19599566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (void (*)(void))MatMult_PC));
19609566063dSJacob Faibussowitsch   PetscCall(MatComputeOperator(Apc, mattype, mat));
19619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Apc));
19623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19634b9ad928SBarry Smith }
19644b9ad928SBarry Smith 
19656c237d78SBarry Smith /*@
19660b4b7b1cSBarry Smith   PCSetCoordinates - sets the coordinates of all the nodes (degrees of freedom in the vector) on the local process
19676c237d78SBarry Smith 
1968c3339decSBarry Smith   Collective
19696c237d78SBarry Smith 
19706c237d78SBarry Smith   Input Parameters:
19710b4b7b1cSBarry Smith + pc     - the `PC` preconditioner context
19726c237d78SBarry Smith . dim    - the dimension of the coordinates 1, 2, or 3
197314893cbeSStefano Zampini . nloc   - the blocked size of the coordinates array
197414893cbeSStefano Zampini - coords - the coordinates array
19756c237d78SBarry Smith 
19766c237d78SBarry Smith   Level: intermediate
19776c237d78SBarry Smith 
19780b4b7b1cSBarry Smith   Notes:
197920f4b53cSBarry Smith   `coords` is an array of the dim coordinates for the nodes on
198020f4b53cSBarry Smith   the local processor, of size `dim`*`nloc`.
19816e415bd2SNuno Nobre   If there are 108 equations (dofs) on a processor
19826e415bd2SNuno Nobre   for a 3d displacement finite element discretization of elasticity (so
198314893cbeSStefano Zampini   that there are nloc = 36 = 108/3 nodes) then the array must have 108
19846c237d78SBarry Smith   double precision values (ie, 3 * 36).  These x y z coordinates
19856c237d78SBarry Smith   should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
19866c237d78SBarry Smith   ... , N-1.z ].
19876c237d78SBarry Smith 
19880b4b7b1cSBarry Smith   The information provided here can be used by some preconditioners, such as `PCGAMG`, to produce a better preconditioner.
19890b4b7b1cSBarry Smith   See also  `MatSetNearNullSpace()`.
19900b4b7b1cSBarry Smith 
1991562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `MatSetNearNullSpace()`
19926c237d78SBarry Smith @*/
1993d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[])
1994d71ae5a4SJacob Faibussowitsch {
19956c237d78SBarry Smith   PetscFunctionBegin;
199632954da3SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
199732954da3SStefano Zampini   PetscValidLogicalCollectiveInt(pc, dim, 2);
199822794d57SStefano Zampini   PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords));
19993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20006c237d78SBarry Smith }
2001fd2dd295SFande Kong 
2002fd2dd295SFande Kong /*@
2003fd2dd295SFande Kong   PCGetInterpolations - Gets interpolation matrices for all levels (except level 0)
2004fd2dd295SFande Kong 
2005c3339decSBarry Smith   Logically Collective
2006fd2dd295SFande Kong 
2007d8d19677SJose E. Roman   Input Parameter:
2008d8d19677SJose E. Roman . pc - the precondition context
2009fd2dd295SFande Kong 
2010d8d19677SJose E. Roman   Output Parameters:
2011d8d19677SJose E. Roman + num_levels     - the number of levels
2012562efe2eSBarry Smith - interpolations - the interpolation matrices (size of `num_levels`-1)
2013fd2dd295SFande Kong 
2014fd2dd295SFande Kong   Level: advanced
2015fd2dd295SFande Kong 
2016562efe2eSBarry Smith   Developer Note:
2017f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
2018fd2dd295SFande Kong 
2019562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()`
2020fd2dd295SFande Kong @*/
2021d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[])
2022d71ae5a4SJacob Faibussowitsch {
2023fd2dd295SFande Kong   PetscFunctionBegin;
2024fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20254f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20264f572ea9SToby Isaac   PetscAssertPointer(interpolations, 3);
2027cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations));
20283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2029fd2dd295SFande Kong }
2030fd2dd295SFande Kong 
2031fd2dd295SFande Kong /*@
2032fd2dd295SFande Kong   PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level)
2033fd2dd295SFande Kong 
2034c3339decSBarry Smith   Logically Collective
2035fd2dd295SFande Kong 
2036d8d19677SJose E. Roman   Input Parameter:
2037d8d19677SJose E. Roman . pc - the precondition context
2038fd2dd295SFande Kong 
2039d8d19677SJose E. Roman   Output Parameters:
2040d8d19677SJose E. Roman + num_levels      - the number of levels
2041562efe2eSBarry Smith - coarseOperators - the coarse operator matrices (size of `num_levels`-1)
2042fd2dd295SFande Kong 
2043fd2dd295SFande Kong   Level: advanced
2044fd2dd295SFande Kong 
2045562efe2eSBarry Smith   Developer Note:
2046f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
2047fd2dd295SFande Kong 
2048562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()`
2049fd2dd295SFande Kong @*/
2050d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
2051d71ae5a4SJacob Faibussowitsch {
2052fd2dd295SFande Kong   PetscFunctionBegin;
2053fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20544f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20554f572ea9SToby Isaac   PetscAssertPointer(coarseOperators, 3);
2056cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators));
20573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2058fd2dd295SFande Kong }
2059