xref: /petsc/src/ksp/pc/interface/precon.c (revision 5d83a8b16d06840f96948f1a43aa9c83c769a60a)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith     The PC (preconditioner) interface routines, callable by users.
34b9ad928SBarry Smith */
4af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscksp.h" I*/
51e25c274SJed Brown #include <petscdm.h>
64b9ad928SBarry Smith 
74b9ad928SBarry Smith /* Logging support */
80700a824SBarry Smith PetscClassId  PC_CLASSID;
9c677e75fSPierre Jolivet PetscLogEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_MatApply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
10011ea8aeSBarry Smith PetscLogEvent PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks;
11ab83eea4SMatthew G. Knepley PetscInt      PetscMGLevelId;
120316ec64SBarry Smith PetscLogStage PCMPIStage;
134b9ad928SBarry Smith 
14b4f8a55aSStefano Zampini PETSC_INTERN PetscErrorCode PCGetDefaultType_Private(PC pc, const char *type[])
15d71ae5a4SJacob Faibussowitsch {
162e70eadcSBarry Smith   PetscMPIInt size;
17b4f8a55aSStefano Zampini   PetscBool   hasopblock, hasopsolve, flg1, flg2, set, flg3, isnormal;
18566e8bf2SBarry Smith 
19566e8bf2SBarry Smith   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
21566e8bf2SBarry Smith   if (pc->pmat) {
22b4f8a55aSStefano Zampini     PetscCall(MatHasOperation(pc->pmat, MATOP_GET_DIAGONAL_BLOCK, &hasopblock));
23b4f8a55aSStefano Zampini     PetscCall(MatHasOperation(pc->pmat, MATOP_SOLVE, &hasopsolve));
24566e8bf2SBarry Smith     if (size == 1) {
259566063dSJacob Faibussowitsch       PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ICC, &flg1));
269566063dSJacob Faibussowitsch       PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ILU, &flg2));
279566063dSJacob Faibussowitsch       PetscCall(MatIsSymmetricKnown(pc->pmat, &set, &flg3));
28da74c943SJose E. Roman       PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isnormal, MATNORMAL, MATNORMALHERMITIAN, NULL));
29ba8a7149SBarry Smith       if (flg1 && (!flg2 || (set && flg3))) {
30566e8bf2SBarry Smith         *type = PCICC;
31566e8bf2SBarry Smith       } else if (flg2) {
32566e8bf2SBarry Smith         *type = PCILU;
33da74c943SJose E. Roman       } else if (isnormal) {
34da74c943SJose E. Roman         *type = PCNONE;
35b4f8a55aSStefano Zampini       } else if (hasopblock) { /* likely is a parallel matrix run on one processor */
36566e8bf2SBarry Smith         *type = PCBJACOBI;
37b4f8a55aSStefano Zampini       } else if (hasopsolve) {
38b4f8a55aSStefano Zampini         *type = PCMAT;
39566e8bf2SBarry Smith       } else {
40566e8bf2SBarry Smith         *type = PCNONE;
41566e8bf2SBarry Smith       }
42566e8bf2SBarry Smith     } else {
43b4f8a55aSStefano Zampini       if (hasopblock) {
44566e8bf2SBarry Smith         *type = PCBJACOBI;
45b4f8a55aSStefano Zampini       } else if (hasopsolve) {
46b4f8a55aSStefano Zampini         *type = PCMAT;
47566e8bf2SBarry Smith       } else {
48566e8bf2SBarry Smith         *type = PCNONE;
49566e8bf2SBarry Smith       }
50566e8bf2SBarry Smith     }
51b4f8a55aSStefano Zampini   } else *type = NULL;
523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53566e8bf2SBarry Smith }
544b9ad928SBarry Smith 
5588584be7SBarry Smith /*@
56562efe2eSBarry Smith   PCReset - Resets a `PC` context to the pcsetupcalled = 0 state and removes any allocated `Vec`s and `Mat`s
5788584be7SBarry Smith 
58c3339decSBarry Smith   Collective
5988584be7SBarry Smith 
6088584be7SBarry Smith   Input Parameter:
6188584be7SBarry Smith . pc - the preconditioner context
6288584be7SBarry Smith 
6388584be7SBarry Smith   Level: developer
6488584be7SBarry Smith 
65f1580f4eSBarry Smith   Note:
66562efe2eSBarry Smith   This allows a `PC` to be reused for a different sized linear system but using the same options that have been previously set in `pc`
6788584be7SBarry Smith 
68562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`
6988584be7SBarry Smith @*/
70d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset(PC pc)
71d71ae5a4SJacob Faibussowitsch {
7288584be7SBarry Smith   PetscFunctionBegin;
7388584be7SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
74dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, reset);
759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleright));
769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleleft));
779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->pmat));
789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->mat));
792fa5cd67SKarl Rupp 
800ce6c5a2SBarry Smith   pc->setupcalled = 0;
813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8288584be7SBarry Smith }
8388584be7SBarry Smith 
840764c050SBarry Smith /*@
85f1580f4eSBarry Smith   PCDestroy - Destroys `PC` context that was created with `PCCreate()`.
864b9ad928SBarry Smith 
87c3339decSBarry Smith   Collective
884b9ad928SBarry Smith 
894b9ad928SBarry Smith   Input Parameter:
904b9ad928SBarry Smith . pc - the preconditioner context
914b9ad928SBarry Smith 
924b9ad928SBarry Smith   Level: developer
934b9ad928SBarry Smith 
94562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`
954b9ad928SBarry Smith @*/
96d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy(PC *pc)
97d71ae5a4SJacob Faibussowitsch {
984b9ad928SBarry Smith   PetscFunctionBegin;
993ba16761SJacob Faibussowitsch   if (!*pc) PetscFunctionReturn(PETSC_SUCCESS);
100f4f49eeaSPierre Jolivet   PetscValidHeaderSpecific(*pc, PC_CLASSID, 1);
101f4f49eeaSPierre Jolivet   if (--((PetscObject)*pc)->refct > 0) {
1029371c9d4SSatish Balay     *pc = NULL;
1033ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1049371c9d4SSatish Balay   }
1054b9ad928SBarry Smith 
1069566063dSJacob Faibussowitsch   PetscCall(PCReset(*pc));
107241cb3abSLisandro Dalcin 
108e04113cfSBarry Smith   /* if memory was published with SAWs then destroy it */
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*pc));
110f4f49eeaSPierre Jolivet   PetscTryTypeMethod(*pc, destroy);
1119566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&(*pc)->dm));
1129566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(pc));
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1144b9ad928SBarry Smith }
1154b9ad928SBarry Smith 
116cc4c1da9SBarry Smith /*@
117c5eb9154SBarry Smith   PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right
1184b9ad928SBarry Smith   scaling as needed by certain time-stepping codes.
1194b9ad928SBarry Smith 
120c3339decSBarry Smith   Logically Collective
1214b9ad928SBarry Smith 
1224b9ad928SBarry Smith   Input Parameter:
1234b9ad928SBarry Smith . pc - the preconditioner context
1244b9ad928SBarry Smith 
1254b9ad928SBarry Smith   Output Parameter:
126f1580f4eSBarry Smith . flag - `PETSC_TRUE` if it applies the scaling
1274b9ad928SBarry Smith 
1284b9ad928SBarry Smith   Level: developer
1294b9ad928SBarry Smith 
130f1580f4eSBarry Smith   Note:
131562efe2eSBarry Smith   If this returns `PETSC_TRUE` then the system solved via the Krylov method is, for left and right preconditioning,
1324b9ad928SBarry Smith 
133562efe2eSBarry Smith   $$
134562efe2eSBarry Smith   \begin{align*}
135562efe2eSBarry Smith   D M A D^{-1} y = D M b  \\
136562efe2eSBarry Smith   D A M D^{-1} z = D b.
137562efe2eSBarry Smith   \end{align*}
138562efe2eSBarry Smith   $$
139562efe2eSBarry Smith 
140562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCSetDiagonalScale()`
1414b9ad928SBarry Smith @*/
142d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetDiagonalScale(PC pc, PetscBool *flag)
143d71ae5a4SJacob Faibussowitsch {
1444b9ad928SBarry Smith   PetscFunctionBegin;
1450700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1464f572ea9SToby Isaac   PetscAssertPointer(flag, 2);
1474b9ad928SBarry Smith   *flag = pc->diagonalscale;
1483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1494b9ad928SBarry Smith }
1504b9ad928SBarry Smith 
1514b9ad928SBarry Smith /*@
152c5eb9154SBarry Smith   PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right
1534b9ad928SBarry Smith   scaling as needed by certain time-stepping codes.
1544b9ad928SBarry Smith 
155c3339decSBarry Smith   Logically Collective
1564b9ad928SBarry Smith 
1574b9ad928SBarry Smith   Input Parameters:
1584b9ad928SBarry Smith + pc - the preconditioner context
1594b9ad928SBarry Smith - s  - scaling vector
1604b9ad928SBarry Smith 
1614b9ad928SBarry Smith   Level: intermediate
1624b9ad928SBarry Smith 
16395452b02SPatrick Sanan   Notes:
164562efe2eSBarry Smith   The system solved via the Krylov method is, for left and right preconditioning,
165562efe2eSBarry Smith   $$
166562efe2eSBarry Smith   \begin{align*}
167562efe2eSBarry Smith   D M A D^{-1} y = D M b \\
168562efe2eSBarry Smith   D A M D^{-1} z = D b.
169562efe2eSBarry Smith   \end{align*}
170562efe2eSBarry Smith   $$
1714b9ad928SBarry Smith 
172562efe2eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.
1734b9ad928SBarry Smith 
174562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCGetDiagonalScale()`
1754b9ad928SBarry Smith @*/
176d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetDiagonalScale(PC pc, Vec s)
177d71ae5a4SJacob Faibussowitsch {
1784b9ad928SBarry Smith   PetscFunctionBegin;
1790700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1800700a824SBarry Smith   PetscValidHeaderSpecific(s, VEC_CLASSID, 2);
1814b9ad928SBarry Smith   pc->diagonalscale = PETSC_TRUE;
1822fa5cd67SKarl Rupp 
1839566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)s));
1849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleleft));
1852fa5cd67SKarl Rupp 
1864b9ad928SBarry Smith   pc->diagonalscaleleft = s;
1872fa5cd67SKarl Rupp 
1889566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(s, &pc->diagonalscaleright));
1899566063dSJacob Faibussowitsch   PetscCall(VecCopy(s, pc->diagonalscaleright));
1909566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(pc->diagonalscaleright));
1913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1924b9ad928SBarry Smith }
1934b9ad928SBarry Smith 
194bc08b0f1SBarry Smith /*@
195c5eb9154SBarry Smith   PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes.
1964b9ad928SBarry Smith 
197c3339decSBarry Smith   Logically Collective
1984b9ad928SBarry Smith 
1994b9ad928SBarry Smith   Input Parameters:
2004b9ad928SBarry Smith + pc  - the preconditioner context
2014b9ad928SBarry Smith . in  - input vector
202a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in)
2034b9ad928SBarry Smith 
2044b9ad928SBarry Smith   Level: intermediate
2054b9ad928SBarry Smith 
20695452b02SPatrick Sanan   Notes:
207562efe2eSBarry Smith   The system solved via the Krylov method is, for left and right preconditioning,
2084b9ad928SBarry Smith 
209562efe2eSBarry Smith   $$
210562efe2eSBarry Smith   \begin{align*}
211562efe2eSBarry Smith   D M A D^{-1} y = D M b  \\
212562efe2eSBarry Smith   D A M D^{-1} z = D b.
213562efe2eSBarry Smith   \end{align*}
214562efe2eSBarry Smith   $$
2154b9ad928SBarry Smith 
216562efe2eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.
2174b9ad928SBarry Smith 
218562efe2eSBarry Smith   If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out`
219562efe2eSBarry Smith 
220562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCSetDiagonalScale()`, `PCDiagonalScaleRight()`, `PCDiagonalScale()`
2214b9ad928SBarry Smith @*/
222d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleLeft(PC pc, Vec in, Vec out)
223d71ae5a4SJacob Faibussowitsch {
2244b9ad928SBarry Smith   PetscFunctionBegin;
2250700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2260700a824SBarry Smith   PetscValidHeaderSpecific(in, VEC_CLASSID, 2);
2270700a824SBarry Smith   PetscValidHeaderSpecific(out, VEC_CLASSID, 3);
2284b9ad928SBarry Smith   if (pc->diagonalscale) {
2299566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(out, pc->diagonalscaleleft, in));
2304b9ad928SBarry Smith   } else if (in != out) {
2319566063dSJacob Faibussowitsch     PetscCall(VecCopy(in, out));
2324b9ad928SBarry Smith   }
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2344b9ad928SBarry Smith }
2354b9ad928SBarry Smith 
236bc08b0f1SBarry Smith /*@
2374b9ad928SBarry Smith   PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
2384b9ad928SBarry Smith 
239c3339decSBarry Smith   Logically Collective
2404b9ad928SBarry Smith 
2414b9ad928SBarry Smith   Input Parameters:
2424b9ad928SBarry Smith + pc  - the preconditioner context
2434b9ad928SBarry Smith . in  - input vector
244a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in)
2454b9ad928SBarry Smith 
2464b9ad928SBarry Smith   Level: intermediate
2474b9ad928SBarry Smith 
24895452b02SPatrick Sanan   Notes:
249562efe2eSBarry Smith   The system solved via the Krylov method is, for left and right preconditioning,
2504b9ad928SBarry Smith 
251562efe2eSBarry Smith   $$
252562efe2eSBarry Smith   \begin{align*}
253562efe2eSBarry Smith   D M A D^{-1} y = D M b  \\
254562efe2eSBarry Smith   D A M D^{-1} z = D b.
255a3524536SPierre Jolivet   \end{align*}
256562efe2eSBarry Smith   $$
2574b9ad928SBarry Smith 
258562efe2eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.
2594b9ad928SBarry Smith 
260562efe2eSBarry Smith   If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out`
261562efe2eSBarry Smith 
262562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCSetDiagonalScale()`, `PCDiagonalScale()`
2634b9ad928SBarry Smith @*/
264d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleRight(PC pc, Vec in, Vec out)
265d71ae5a4SJacob Faibussowitsch {
2664b9ad928SBarry Smith   PetscFunctionBegin;
2670700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2680700a824SBarry Smith   PetscValidHeaderSpecific(in, VEC_CLASSID, 2);
2690700a824SBarry Smith   PetscValidHeaderSpecific(out, VEC_CLASSID, 3);
2704b9ad928SBarry Smith   if (pc->diagonalscale) {
2719566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(out, pc->diagonalscaleright, in));
2724b9ad928SBarry Smith   } else if (in != out) {
2739566063dSJacob Faibussowitsch     PetscCall(VecCopy(in, out));
2744b9ad928SBarry Smith   }
2753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2764b9ad928SBarry Smith }
2774b9ad928SBarry Smith 
27849517cdeSBarry Smith /*@
27949517cdeSBarry Smith   PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the
280f1580f4eSBarry Smith   operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
281f1580f4eSBarry Smith   `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.
28249517cdeSBarry Smith 
283c3339decSBarry Smith   Logically Collective
28449517cdeSBarry Smith 
28549517cdeSBarry Smith   Input Parameters:
28649517cdeSBarry Smith + pc  - the preconditioner context
287f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)
28849517cdeSBarry Smith 
28949517cdeSBarry Smith   Options Database Key:
290562efe2eSBarry Smith . -pc_use_amat <true,false> - use the amat argument to `KSPSetOperators()` or `PCSetOperators()` to apply the operator
29149517cdeSBarry Smith 
29220f4b53cSBarry Smith   Level: intermediate
29320f4b53cSBarry Smith 
294f1580f4eSBarry Smith   Note:
29549517cdeSBarry Smith   For the common case in which the linear system matrix and the matrix used to construct the
296562efe2eSBarry Smith   preconditioner are identical, this routine has no affect.
29749517cdeSBarry Smith 
298562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetUseAmat()`, `PCBJACOBI`, `PGMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`,
299562efe2eSBarry Smith           `KSPSetOperators()`, `PCSetOperators()`
30049517cdeSBarry Smith @*/
301d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUseAmat(PC pc, PetscBool flg)
302d71ae5a4SJacob Faibussowitsch {
30349517cdeSBarry Smith   PetscFunctionBegin;
30449517cdeSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
30549517cdeSBarry Smith   pc->useAmat = flg;
3063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30749517cdeSBarry Smith }
30849517cdeSBarry Smith 
309422a814eSBarry Smith /*@
310562efe2eSBarry Smith   PCSetErrorIfFailure - Causes `PC` to generate an error if a floating point exception, for example a zero pivot, is detected.
311422a814eSBarry Smith 
312c3339decSBarry Smith   Logically Collective
313422a814eSBarry Smith 
314422a814eSBarry Smith   Input Parameters:
315562efe2eSBarry Smith + pc  - iterative context obtained from `PCCreate()`
316f1580f4eSBarry Smith - flg - `PETSC_TRUE` indicates you want the error generated
317422a814eSBarry Smith 
318422a814eSBarry Smith   Level: advanced
319422a814eSBarry Smith 
320422a814eSBarry Smith   Notes:
321f1580f4eSBarry Smith   Normally PETSc continues if a linear solver fails due to a failed setup of a preconditioner, you can call `KSPGetConvergedReason()` after a `KSPSolve()`
322422a814eSBarry Smith   to determine if it has converged or failed. Or use -ksp_error_if_not_converged to cause the program to terminate as soon as lack of convergence is
323422a814eSBarry Smith   detected.
324422a814eSBarry Smith 
325562efe2eSBarry Smith   This is propagated into `KSP`s used by this `PC`, which then propagate it into `PC`s used by those `KSP`s
326422a814eSBarry Smith 
327562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPSetErrorIfNotConverged()`, `PCGetInitialGuessNonzero()`, `PCSetInitialGuessKnoll()`, `PCGetInitialGuessKnoll()`
328422a814eSBarry Smith @*/
329d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetErrorIfFailure(PC pc, PetscBool flg)
330d71ae5a4SJacob Faibussowitsch {
331422a814eSBarry Smith   PetscFunctionBegin;
332422a814eSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
333422a814eSBarry Smith   PetscValidLogicalCollectiveBool(pc, flg, 2);
334422a814eSBarry Smith   pc->erroriffailure = flg;
3353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
336422a814eSBarry Smith }
337422a814eSBarry Smith 
33849517cdeSBarry Smith /*@
33949517cdeSBarry Smith   PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the
340f1580f4eSBarry Smith   operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
341f1580f4eSBarry Smith   `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.
34249517cdeSBarry Smith 
343c3339decSBarry Smith   Logically Collective
34449517cdeSBarry Smith 
34549517cdeSBarry Smith   Input Parameter:
34649517cdeSBarry Smith . pc - the preconditioner context
34749517cdeSBarry Smith 
34849517cdeSBarry Smith   Output Parameter:
349f1580f4eSBarry Smith . flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)
35049517cdeSBarry Smith 
35120f4b53cSBarry Smith   Level: intermediate
35220f4b53cSBarry Smith 
353f1580f4eSBarry Smith   Note:
35449517cdeSBarry Smith   For the common case in which the linear system matrix and the matrix used to construct the
35549517cdeSBarry Smith   preconditioner are identical, this routine is does nothing.
35649517cdeSBarry Smith 
357562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUseAmat()`, `PCBJACOBI`, `PGMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`
35849517cdeSBarry Smith @*/
359d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetUseAmat(PC pc, PetscBool *flg)
360d71ae5a4SJacob Faibussowitsch {
36149517cdeSBarry Smith   PetscFunctionBegin;
36249517cdeSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
36349517cdeSBarry Smith   *flg = pc->useAmat;
3643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36549517cdeSBarry Smith }
36649517cdeSBarry Smith 
367f39d8e23SSatish Balay /*@
3683821be0aSBarry Smith   PCSetKSPNestLevel - sets the amount of nesting the `KSP` that contains this `PC` has
3693821be0aSBarry Smith 
3703821be0aSBarry Smith   Collective
3713821be0aSBarry Smith 
3723821be0aSBarry Smith   Input Parameters:
3733821be0aSBarry Smith + pc    - the `PC`
3743821be0aSBarry Smith - level - the nest level
3753821be0aSBarry Smith 
3763821be0aSBarry Smith   Level: developer
3773821be0aSBarry Smith 
3783821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCGetKSPNestLevel()`, `KSPSetNestLevel()`
3793821be0aSBarry Smith @*/
3803821be0aSBarry Smith PetscErrorCode PCSetKSPNestLevel(PC pc, PetscInt level)
3813821be0aSBarry Smith {
3823821be0aSBarry Smith   PetscFunctionBegin;
3837a99bfcaSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
3847a99bfcaSBarry Smith   PetscValidLogicalCollectiveInt(pc, level, 2);
3853821be0aSBarry Smith   pc->kspnestlevel = level;
3863821be0aSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3873821be0aSBarry Smith }
3883821be0aSBarry Smith 
3893821be0aSBarry Smith /*@
3903821be0aSBarry Smith   PCGetKSPNestLevel - gets the amount of nesting the `KSP` that contains this `PC` has
3913821be0aSBarry Smith 
3927a99bfcaSBarry Smith   Not Collective
3933821be0aSBarry Smith 
3943821be0aSBarry Smith   Input Parameter:
3953821be0aSBarry Smith . pc - the `PC`
3963821be0aSBarry Smith 
3973821be0aSBarry Smith   Output Parameter:
3983821be0aSBarry Smith . level - the nest level
3993821be0aSBarry Smith 
4003821be0aSBarry Smith   Level: developer
4013821be0aSBarry Smith 
4023821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `KSPGetNestLevel()`
4033821be0aSBarry Smith @*/
4043821be0aSBarry Smith PetscErrorCode PCGetKSPNestLevel(PC pc, PetscInt *level)
4053821be0aSBarry Smith {
4063821be0aSBarry Smith   PetscFunctionBegin;
4077a99bfcaSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4087a99bfcaSBarry Smith   PetscAssertPointer(level, 2);
4093821be0aSBarry Smith   *level = pc->kspnestlevel;
4103821be0aSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
4113821be0aSBarry Smith }
4123821be0aSBarry Smith 
4133821be0aSBarry Smith /*@
414f1580f4eSBarry Smith   PCCreate - Creates a preconditioner context, `PC`
4154b9ad928SBarry Smith 
416d083f849SBarry Smith   Collective
4174b9ad928SBarry Smith 
4184b9ad928SBarry Smith   Input Parameter:
4194b9ad928SBarry Smith . comm - MPI communicator
4204b9ad928SBarry Smith 
4214b9ad928SBarry Smith   Output Parameter:
4227a99bfcaSBarry Smith . newpc - location to put the preconditioner context
4234b9ad928SBarry Smith 
42420f4b53cSBarry Smith   Level: developer
42520f4b53cSBarry Smith 
426f1580f4eSBarry Smith   Note:
427f1580f4eSBarry Smith   The default preconditioner for sparse matrices is `PCILU` or `PCICC` with 0 fill on one process and block Jacobi (`PCBJACOBI`) with `PCILU` or `PCICC`
428f1580f4eSBarry Smith   in parallel. For dense matrices it is always `PCNONE`.
4294b9ad928SBarry Smith 
430562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCApply()`, `PCDestroy()`
4314b9ad928SBarry Smith @*/
432d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreate(MPI_Comm comm, PC *newpc)
433d71ae5a4SJacob Faibussowitsch {
4344b9ad928SBarry Smith   PC pc;
4354b9ad928SBarry Smith 
4364b9ad928SBarry Smith   PetscFunctionBegin;
4374f572ea9SToby Isaac   PetscAssertPointer(newpc, 2);
4389566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
4394b9ad928SBarry Smith 
4409566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(pc, PC_CLASSID, "PC", "Preconditioner", "PC", comm, PCDestroy, PCView));
4410a545947SLisandro Dalcin   pc->mat                  = NULL;
4420a545947SLisandro Dalcin   pc->pmat                 = NULL;
4434b9ad928SBarry Smith   pc->setupcalled          = 0;
4440c24e6a1SHong Zhang   pc->setfromoptionscalled = 0;
4450a545947SLisandro Dalcin   pc->data                 = NULL;
4464b9ad928SBarry Smith   pc->diagonalscale        = PETSC_FALSE;
4470a545947SLisandro Dalcin   pc->diagonalscaleleft    = NULL;
4480a545947SLisandro Dalcin   pc->diagonalscaleright   = NULL;
4494b9ad928SBarry Smith 
4500a545947SLisandro Dalcin   pc->modifysubmatrices  = NULL;
4510a545947SLisandro Dalcin   pc->modifysubmatricesP = NULL;
4522fa5cd67SKarl Rupp 
453a7d21a52SLisandro Dalcin   *newpc = pc;
4543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4554b9ad928SBarry Smith }
4564b9ad928SBarry Smith 
4574b9ad928SBarry Smith /*@
4584b9ad928SBarry Smith   PCApply - Applies the preconditioner to a vector.
4594b9ad928SBarry Smith 
460c3339decSBarry Smith   Collective
4614b9ad928SBarry Smith 
4624b9ad928SBarry Smith   Input Parameters:
4634b9ad928SBarry Smith + pc - the preconditioner context
4644b9ad928SBarry Smith - x  - input vector
4654b9ad928SBarry Smith 
4664b9ad928SBarry Smith   Output Parameter:
4674b9ad928SBarry Smith . y - output vector
4684b9ad928SBarry Smith 
4694b9ad928SBarry Smith   Level: developer
4704b9ad928SBarry Smith 
471562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`, `PCApplyBAorAB()`
4724b9ad928SBarry Smith @*/
473d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApply(PC pc, Vec x, Vec y)
474d71ae5a4SJacob Faibussowitsch {
475106e20bfSBarry Smith   PetscInt m, n, mv, nv;
4764b9ad928SBarry Smith 
4774b9ad928SBarry Smith   PetscFunctionBegin;
4780700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4790700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
4800700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
4817827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
482e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
483540bdfbdSVaclav Hapla   /* use pmat to check vector sizes since for KSPLSQR the pmat may be of a different size than mat */
4849566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &m, &n));
4859566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x, &mv));
4869566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(y, &nv));
487540bdfbdSVaclav Hapla   /* check pmat * y = x is feasible */
48863a3b9bcSJacob 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);
48963a3b9bcSJacob 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);
4909566063dSJacob Faibussowitsch   PetscCall(VecSetErrorIfLocked(y, 3));
491106e20bfSBarry Smith 
4929566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
4939566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
4949566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
495dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, apply, x, y);
4969566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
497e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
4989566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
4993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5004b9ad928SBarry Smith }
5014b9ad928SBarry Smith 
5024b9ad928SBarry Smith /*@
503562efe2eSBarry Smith   PCMatApply - Applies the preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApply()`, `Y` and `X` must be different matrices.
504c41ea70eSPierre Jolivet 
505c3339decSBarry Smith   Collective
506c677e75fSPierre Jolivet 
507c677e75fSPierre Jolivet   Input Parameters:
508c677e75fSPierre Jolivet + pc - the preconditioner context
509c677e75fSPierre Jolivet - X  - block of input vectors
510c677e75fSPierre Jolivet 
511c677e75fSPierre Jolivet   Output Parameter:
512c677e75fSPierre Jolivet . Y - block of output vectors
513c677e75fSPierre Jolivet 
514c677e75fSPierre Jolivet   Level: developer
515c677e75fSPierre Jolivet 
516562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `KSPMatSolve()`
517c677e75fSPierre Jolivet @*/
518d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMatApply(PC pc, Mat X, Mat Y)
519d71ae5a4SJacob Faibussowitsch {
520c677e75fSPierre Jolivet   Mat       A;
521c677e75fSPierre Jolivet   Vec       cy, cx;
522bd82155bSPierre Jolivet   PetscInt  m1, M1, m2, M2, n1, N1, n2, N2, m3, M3, n3, N3;
523c677e75fSPierre Jolivet   PetscBool match;
524c677e75fSPierre Jolivet 
525c677e75fSPierre Jolivet   PetscFunctionBegin;
526c677e75fSPierre Jolivet   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
527e2b200f6SPierre Jolivet   PetscValidHeaderSpecific(X, MAT_CLASSID, 2);
528e2b200f6SPierre Jolivet   PetscValidHeaderSpecific(Y, MAT_CLASSID, 3);
529e2b200f6SPierre Jolivet   PetscCheckSameComm(pc, 1, X, 2);
530e2b200f6SPierre Jolivet   PetscCheckSameComm(pc, 1, Y, 3);
5317827d75bSBarry Smith   PetscCheck(Y != X, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "Y and X must be different matrices");
5329566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &A));
5339566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m3, &n3));
5349566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m2, &n2));
5359566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &m1, &n1));
5369566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M3, &N3));
5379566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &M2, &N2));
5389566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Y, &M1, &N1));
53963a3b9bcSJacob 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);
54063a3b9bcSJacob 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);
54163a3b9bcSJacob 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);
5429566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &match, MATSEQDENSE, MATMPIDENSE, ""));
54328b400f6SJacob Faibussowitsch   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of output vectors not stored in a dense Mat");
5449566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
54528b400f6SJacob Faibussowitsch   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of input vectors not stored in a dense Mat");
5469566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
547c677e75fSPierre Jolivet   if (pc->ops->matapply) {
5489566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0));
549dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, matapply, X, Y);
5509566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0));
551c677e75fSPierre Jolivet   } else {
5529566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "PC type %s applying column by column\n", ((PetscObject)pc)->type_name));
553bd82155bSPierre Jolivet     for (n1 = 0; n1 < N1; ++n1) {
5549566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecRead(X, n1, &cx));
5559566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecWrite(Y, n1, &cy));
5569566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, cx, cy));
5579566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecWrite(Y, n1, &cy));
5589566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecRead(X, n1, &cx));
559c677e75fSPierre Jolivet     }
560c677e75fSPierre Jolivet   }
5613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
562c677e75fSPierre Jolivet }
563c677e75fSPierre Jolivet 
564c677e75fSPierre Jolivet /*@
5654b9ad928SBarry Smith   PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
5664b9ad928SBarry Smith 
567c3339decSBarry Smith   Collective
5684b9ad928SBarry Smith 
5694b9ad928SBarry Smith   Input Parameters:
5704b9ad928SBarry Smith + pc - the preconditioner context
5714b9ad928SBarry Smith - x  - input vector
5724b9ad928SBarry Smith 
5734b9ad928SBarry Smith   Output Parameter:
5744b9ad928SBarry Smith . y - output vector
5754b9ad928SBarry Smith 
57620f4b53cSBarry Smith   Level: developer
57720f4b53cSBarry Smith 
578f1580f4eSBarry Smith   Note:
579f1580f4eSBarry Smith   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
5804b9ad928SBarry Smith 
581562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricRight()`
5824b9ad928SBarry Smith @*/
583d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricLeft(PC pc, Vec x, Vec y)
584d71ae5a4SJacob Faibussowitsch {
5854b9ad928SBarry Smith   PetscFunctionBegin;
5860700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
5870700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
5880700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
5897827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
590e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
5919566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
5929566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
5939566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplySymmetricLeft, pc, x, y, 0));
594dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applysymmetricleft, x, y);
5959566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplySymmetricLeft, pc, x, y, 0));
5969566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
597e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
5983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5994b9ad928SBarry Smith }
6004b9ad928SBarry Smith 
6014b9ad928SBarry Smith /*@
6024b9ad928SBarry Smith   PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
6034b9ad928SBarry Smith 
604c3339decSBarry Smith   Collective
6054b9ad928SBarry Smith 
6064b9ad928SBarry Smith   Input Parameters:
6074b9ad928SBarry Smith + pc - the preconditioner context
6084b9ad928SBarry Smith - x  - input vector
6094b9ad928SBarry Smith 
6104b9ad928SBarry Smith   Output Parameter:
6114b9ad928SBarry Smith . y - output vector
6124b9ad928SBarry Smith 
6134b9ad928SBarry Smith   Level: developer
6144b9ad928SBarry Smith 
615f1580f4eSBarry Smith   Note:
616f1580f4eSBarry Smith   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
6174b9ad928SBarry Smith 
618562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricLeft()`
6194b9ad928SBarry Smith @*/
620d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricRight(PC pc, Vec x, Vec y)
621d71ae5a4SJacob Faibussowitsch {
6224b9ad928SBarry Smith   PetscFunctionBegin;
6230700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6240700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6250700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6267827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
627e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
6289566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
6299566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
6309566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplySymmetricRight, pc, x, y, 0));
631dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applysymmetricright, x, y);
6329566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplySymmetricRight, pc, x, y, 0));
6339566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
634e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
6353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6364b9ad928SBarry Smith }
6374b9ad928SBarry Smith 
6384b9ad928SBarry Smith /*@
6394b9ad928SBarry Smith   PCApplyTranspose - Applies the transpose of preconditioner to a vector.
6404b9ad928SBarry Smith 
641c3339decSBarry Smith   Collective
6424b9ad928SBarry Smith 
6434b9ad928SBarry Smith   Input Parameters:
6444b9ad928SBarry Smith + pc - the preconditioner context
6454b9ad928SBarry Smith - x  - input vector
6464b9ad928SBarry Smith 
6474b9ad928SBarry Smith   Output Parameter:
6484b9ad928SBarry Smith . y - output vector
6494b9ad928SBarry Smith 
65020f4b53cSBarry Smith   Level: developer
65120f4b53cSBarry Smith 
652f1580f4eSBarry Smith   Note:
65395452b02SPatrick Sanan   For complex numbers this applies the non-Hermitian transpose.
6544c97465dSBarry Smith 
655562efe2eSBarry Smith   Developer Note:
656f1580f4eSBarry Smith   We need to implement a `PCApplyHermitianTranspose()`
6574c97465dSBarry Smith 
658562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyBAorAB()`, `PCApplyBAorABTranspose()`, `PCApplyTransposeExists()`
6594b9ad928SBarry Smith @*/
660d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTranspose(PC pc, Vec x, Vec y)
661d71ae5a4SJacob Faibussowitsch {
6624b9ad928SBarry Smith   PetscFunctionBegin;
6630700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6640700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6650700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6667827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
667e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
6689566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
6699566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
6709566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
671dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applytranspose, x, y);
6729566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
6739566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
674e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
6753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6764b9ad928SBarry Smith }
6774b9ad928SBarry Smith 
6784b9ad928SBarry Smith /*@
6799cc050e5SLisandro Dalcin   PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
6804b9ad928SBarry Smith 
681c3339decSBarry Smith   Collective
6824b9ad928SBarry Smith 
683f1580f4eSBarry Smith   Input Parameter:
6844b9ad928SBarry Smith . pc - the preconditioner context
6854b9ad928SBarry Smith 
6864b9ad928SBarry Smith   Output Parameter:
687f1580f4eSBarry Smith . flg - `PETSC_TRUE` if a transpose operation is defined
6884b9ad928SBarry Smith 
6894b9ad928SBarry Smith   Level: developer
6904b9ad928SBarry Smith 
691562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`
6924b9ad928SBarry Smith @*/
693d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTransposeExists(PC pc, PetscBool *flg)
694d71ae5a4SJacob Faibussowitsch {
6954b9ad928SBarry Smith   PetscFunctionBegin;
6960700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6974f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
6986ce5e81cSLisandro Dalcin   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
6996ce5e81cSLisandro Dalcin   else *flg = PETSC_FALSE;
7003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7014b9ad928SBarry Smith }
7024b9ad928SBarry Smith 
7034b9ad928SBarry Smith /*@
704562efe2eSBarry Smith   PCApplyBAorAB - Applies the preconditioner and operator to a vector. $y = B*A*x $ or $ y = A*B*x$.
7054b9ad928SBarry Smith 
706c3339decSBarry Smith   Collective
7074b9ad928SBarry Smith 
7084b9ad928SBarry Smith   Input Parameters:
7094b9ad928SBarry Smith + pc   - the preconditioner context
710f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
7114b9ad928SBarry Smith . x    - input vector
7124b9ad928SBarry Smith - work - work vector
7134b9ad928SBarry Smith 
7144b9ad928SBarry Smith   Output Parameter:
7154b9ad928SBarry Smith . y - output vector
7164b9ad928SBarry Smith 
7174b9ad928SBarry Smith   Level: developer
7184b9ad928SBarry Smith 
719f1580f4eSBarry Smith   Note:
720562efe2eSBarry 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.
721562efe2eSBarry Smith   The specific `KSPSolve()` method must also be written to handle the post-solve "correction" for the diagonal scaling.
72213e0d083SBarry Smith 
723562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorABTranspose()`
7244b9ad928SBarry Smith @*/
725d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorAB(PC pc, PCSide side, Vec x, Vec y, Vec work)
726d71ae5a4SJacob Faibussowitsch {
7274b9ad928SBarry Smith   PetscFunctionBegin;
7280700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
729a69c7061SStefano Zampini   PetscValidLogicalCollectiveEnum(pc, side, 2);
7300700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
7310700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
7320700a824SBarry Smith   PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
733a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, x, 3);
734a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, y, 4);
735a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, work, 5);
7367827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
7377827d75bSBarry 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");
7387827d75bSBarry Smith   PetscCheck(!pc->diagonalscale || side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot include diagonal scaling with symmetric preconditioner application");
739e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
7404b9ad928SBarry Smith 
7419566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
7424b9ad928SBarry Smith   if (pc->diagonalscale) {
7434b9ad928SBarry Smith     if (pc->ops->applyBA) {
7444b9ad928SBarry Smith       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
7459566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &work2));
7469566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, work2));
747dbbe0bcdSBarry Smith       PetscUseTypeMethod(pc, applyBA, side, work2, y, work);
7489566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7499566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&work2));
7504b9ad928SBarry Smith     } else if (side == PC_RIGHT) {
7519566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, y));
7529566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, y, work));
7539566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7549566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7554b9ad928SBarry Smith     } else if (side == PC_LEFT) {
7569566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, y));
7579566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, y, work));
7589566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, work, y));
7599566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7607827d75bSBarry Smith     } else PetscCheck(side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot provide diagonal scaling with symmetric application of preconditioner");
7614b9ad928SBarry Smith   } else {
7624b9ad928SBarry Smith     if (pc->ops->applyBA) {
763dbbe0bcdSBarry Smith       PetscUseTypeMethod(pc, applyBA, side, x, y, work);
7644b9ad928SBarry Smith     } else if (side == PC_RIGHT) {
7659566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, x, work));
7669566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7674b9ad928SBarry Smith     } else if (side == PC_LEFT) {
7689566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, x, work));
7699566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, work, y));
7704b9ad928SBarry Smith     } else if (side == PC_SYMMETRIC) {
7714b9ad928SBarry Smith       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
7729566063dSJacob Faibussowitsch       PetscCall(PCApplySymmetricRight(pc, x, work));
7739566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7749566063dSJacob Faibussowitsch       PetscCall(VecCopy(y, work));
7759566063dSJacob Faibussowitsch       PetscCall(PCApplySymmetricLeft(pc, work, y));
7764b9ad928SBarry Smith     }
7774b9ad928SBarry Smith   }
778e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
7793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7804b9ad928SBarry Smith }
7814b9ad928SBarry Smith 
7824b9ad928SBarry Smith /*@
7834b9ad928SBarry Smith   PCApplyBAorABTranspose - Applies the transpose of the preconditioner
784562efe2eSBarry Smith   and operator to a vector. That is, applies $B^T * A^T$ with left preconditioning,
785562efe2eSBarry Smith   NOT $(B*A)^T = A^T*B^T$.
7864b9ad928SBarry Smith 
787c3339decSBarry Smith   Collective
7884b9ad928SBarry Smith 
7894b9ad928SBarry Smith   Input Parameters:
7904b9ad928SBarry Smith + pc   - the preconditioner context
791f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
7924b9ad928SBarry Smith . x    - input vector
7934b9ad928SBarry Smith - work - work vector
7944b9ad928SBarry Smith 
7954b9ad928SBarry Smith   Output Parameter:
7964b9ad928SBarry Smith . y - output vector
7974b9ad928SBarry Smith 
79820f4b53cSBarry Smith   Level: developer
79920f4b53cSBarry Smith 
800f1580f4eSBarry Smith   Note:
801562efe2eSBarry 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
802562efe2eSBarry Smith   defined by $B^T$. This is why this has the funny form that it computes $B^T * A^T$
8039b3038f0SBarry Smith 
804562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorAB()`
8054b9ad928SBarry Smith @*/
806d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorABTranspose(PC pc, PCSide side, Vec x, Vec y, Vec work)
807d71ae5a4SJacob Faibussowitsch {
8084b9ad928SBarry Smith   PetscFunctionBegin;
8090700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8100700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
8110700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
8120700a824SBarry Smith   PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
8137827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
814e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
8154b9ad928SBarry Smith   if (pc->ops->applyBAtranspose) {
816dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, applyBAtranspose, side, x, y, work);
817e0f629ddSJacob Faibussowitsch     if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8183ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
8194b9ad928SBarry Smith   }
8207827d75bSBarry Smith   PetscCheck(side == PC_LEFT || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right or left");
8214b9ad928SBarry Smith 
8229566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
8234b9ad928SBarry Smith   if (side == PC_RIGHT) {
8249566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(pc, x, work));
8259566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->mat, work, y));
8264b9ad928SBarry Smith   } else if (side == PC_LEFT) {
8279566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->mat, x, work));
8289566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(pc, work, y));
8294b9ad928SBarry Smith   }
8304b9ad928SBarry Smith   /* add support for PC_SYMMETRIC */
831e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8334b9ad928SBarry Smith }
8344b9ad928SBarry Smith 
8354b9ad928SBarry Smith /*@
8364b9ad928SBarry Smith   PCApplyRichardsonExists - Determines whether a particular preconditioner has a
8374b9ad928SBarry Smith   built-in fast application of Richardson's method.
8384b9ad928SBarry Smith 
8394b9ad928SBarry Smith   Not Collective
8404b9ad928SBarry Smith 
8414b9ad928SBarry Smith   Input Parameter:
8424b9ad928SBarry Smith . pc - the preconditioner
8434b9ad928SBarry Smith 
8444b9ad928SBarry Smith   Output Parameter:
845f1580f4eSBarry Smith . exists - `PETSC_TRUE` or `PETSC_FALSE`
8464b9ad928SBarry Smith 
8474b9ad928SBarry Smith   Level: developer
8484b9ad928SBarry Smith 
84939b1ba1bSPierre Jolivet .seealso: [](ch_ksp), `PC`, `KSPRICHARDSON`, `PCApplyRichardson()`
8504b9ad928SBarry Smith @*/
851d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyRichardsonExists(PC pc, PetscBool *exists)
852d71ae5a4SJacob Faibussowitsch {
8534b9ad928SBarry Smith   PetscFunctionBegin;
8540700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8554f572ea9SToby Isaac   PetscAssertPointer(exists, 2);
8564b9ad928SBarry Smith   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
8574b9ad928SBarry Smith   else *exists = PETSC_FALSE;
8583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8594b9ad928SBarry Smith }
8604b9ad928SBarry Smith 
8614b9ad928SBarry Smith /*@
8624b9ad928SBarry Smith   PCApplyRichardson - Applies several steps of Richardson iteration with
8634b9ad928SBarry Smith   the particular preconditioner. This routine is usually used by the
8644b9ad928SBarry Smith   Krylov solvers and not the application code directly.
8654b9ad928SBarry Smith 
866c3339decSBarry Smith   Collective
8674b9ad928SBarry Smith 
8684b9ad928SBarry Smith   Input Parameters:
8694b9ad928SBarry Smith + pc        - the preconditioner context
870dd8e379bSPierre Jolivet . b         - the right-hand side
8714b9ad928SBarry Smith . w         - one work vector
8724b9ad928SBarry Smith . rtol      - relative decrease in residual norm convergence criteria
87370441072SBarry Smith . abstol    - absolute residual norm convergence criteria
8744b9ad928SBarry Smith . dtol      - divergence residual norm increase criteria
8757319c654SBarry Smith . its       - the number of iterations to apply.
8767319c654SBarry Smith - guesszero - if the input x contains nonzero initial guess
8774b9ad928SBarry Smith 
878d8d19677SJose E. Roman   Output Parameters:
8794d0a8057SBarry Smith + outits - number of iterations actually used (for SOR this always equals its)
8804d0a8057SBarry Smith . reason - the reason the apply terminated
881f1580f4eSBarry Smith - y      - the solution (also contains initial guess if guesszero is `PETSC_FALSE`
8824b9ad928SBarry Smith 
88320f4b53cSBarry Smith   Level: developer
88420f4b53cSBarry Smith 
8854b9ad928SBarry Smith   Notes:
8864b9ad928SBarry Smith   Most preconditioners do not support this function. Use the command
887f1580f4eSBarry Smith   `PCApplyRichardsonExists()` to determine if one does.
8884b9ad928SBarry Smith 
889f1580f4eSBarry Smith   Except for the `PCMG` this routine ignores the convergence tolerances
8904b9ad928SBarry Smith   and always runs for the number of iterations
8914b9ad928SBarry Smith 
892562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyRichardsonExists()`
8934b9ad928SBarry Smith @*/
894d71ae5a4SJacob 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)
895d71ae5a4SJacob Faibussowitsch {
8964b9ad928SBarry Smith   PetscFunctionBegin;
8970700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8980700a824SBarry Smith   PetscValidHeaderSpecific(b, VEC_CLASSID, 2);
8990700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
9000700a824SBarry Smith   PetscValidHeaderSpecific(w, VEC_CLASSID, 4);
9017827d75bSBarry Smith   PetscCheck(b != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "b and y must be different vectors");
9029566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
903dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applyrichardson, b, y, w, rtol, abstol, dtol, its, guesszero, outits, reason);
9043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9054b9ad928SBarry Smith }
9064b9ad928SBarry Smith 
907422a814eSBarry Smith /*@
908f1580f4eSBarry Smith   PCSetFailedReason - Sets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
9091b2b9847SBarry Smith 
910c3339decSBarry Smith   Logically Collective
9111b2b9847SBarry Smith 
912d8d19677SJose E. Roman   Input Parameters:
9131b2b9847SBarry Smith + pc     - the preconditioner context
9141b2b9847SBarry Smith - reason - the reason it failedx
9151b2b9847SBarry Smith 
9161b2b9847SBarry Smith   Level: advanced
9171b2b9847SBarry Smith 
918562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCFailedReason`
9191b2b9847SBarry Smith @*/
920d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFailedReason(PC pc, PCFailedReason reason)
921d71ae5a4SJacob Faibussowitsch {
9221b2b9847SBarry Smith   PetscFunctionBegin;
9236479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9241b2b9847SBarry Smith   pc->failedreason = reason;
9253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9261b2b9847SBarry Smith }
9271b2b9847SBarry Smith 
9281b2b9847SBarry Smith /*@
929f1580f4eSBarry Smith   PCGetFailedReason - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
930422a814eSBarry Smith 
931c3339decSBarry Smith   Logically Collective
932422a814eSBarry Smith 
933422a814eSBarry Smith   Input Parameter:
934422a814eSBarry Smith . pc - the preconditioner context
935422a814eSBarry Smith 
936422a814eSBarry Smith   Output Parameter:
9371b2b9847SBarry Smith . reason - the reason it failed
938422a814eSBarry Smith 
939422a814eSBarry Smith   Level: advanced
940422a814eSBarry Smith 
941f1580f4eSBarry Smith   Note:
942f1580f4eSBarry Smith   This is the maximum over reason over all ranks in the PC communicator. It is only valid after
9436479e835SStefano Zampini   a call `KSPCheckDot()` or  `KSPCheckNorm()` inside a `KSPSolve()` or `PCReduceFailedReason()`.
9446479e835SStefano Zampini   It is not valid immediately after a `PCSetUp()` or `PCApply()`, then use `PCGetFailedReasonRank()`
9451b2b9847SBarry Smith 
946a94f484eSPierre Jolivet .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReasonRank()`, `PCSetFailedReason()`
947422a814eSBarry Smith @*/
948d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason)
949d71ae5a4SJacob Faibussowitsch {
950422a814eSBarry Smith   PetscFunctionBegin;
9516479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9529754764cSHong Zhang   if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled;
9530ae0b32dSHong Zhang   else *reason = pc->failedreason;
9543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
955422a814eSBarry Smith }
956422a814eSBarry Smith 
9571b2b9847SBarry Smith /*@
958f1580f4eSBarry Smith   PCGetFailedReasonRank - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail on this MPI rank
9591b2b9847SBarry Smith 
960f1580f4eSBarry Smith   Not Collective
9611b2b9847SBarry Smith 
9621b2b9847SBarry Smith   Input Parameter:
9631b2b9847SBarry Smith . pc - the preconditioner context
9641b2b9847SBarry Smith 
9651b2b9847SBarry Smith   Output Parameter:
9661b2b9847SBarry Smith . reason - the reason it failed
9671b2b9847SBarry Smith 
96820f4b53cSBarry Smith   Level: advanced
96920f4b53cSBarry Smith 
970f1580f4eSBarry Smith   Note:
971562efe2eSBarry Smith   Different processes may have different reasons or no reason, see `PCGetFailedReason()`
9721b2b9847SBarry Smith 
973562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`, `PCReduceFailedReason()`
9741b2b9847SBarry Smith @*/
975d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReasonRank(PC pc, PCFailedReason *reason)
976d71ae5a4SJacob Faibussowitsch {
9771b2b9847SBarry Smith   PetscFunctionBegin;
9786479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9791b2b9847SBarry Smith   if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled;
9801b2b9847SBarry Smith   else *reason = pc->failedreason;
9813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9821b2b9847SBarry Smith }
983422a814eSBarry Smith 
9846479e835SStefano Zampini /*@
9856479e835SStefano Zampini   PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC`
9866479e835SStefano Zampini 
9876479e835SStefano Zampini   Collective
9886479e835SStefano Zampini 
9896479e835SStefano Zampini   Input Parameter:
9906479e835SStefano Zampini . pc - the preconditioner context
9916479e835SStefano Zampini 
9926479e835SStefano Zampini   Level: advanced
9936479e835SStefano Zampini 
9946479e835SStefano Zampini   Note:
995562efe2eSBarry Smith   Different MPI processes may have different reasons or no reason, see `PCGetFailedReason()`. This routine
9966479e835SStefano Zampini   makes them have a common value (failure if any MPI process had a failure).
9976479e835SStefano Zampini 
998562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`
9996479e835SStefano Zampini @*/
10006479e835SStefano Zampini PetscErrorCode PCReduceFailedReason(PC pc)
10016479e835SStefano Zampini {
10026479e835SStefano Zampini   PetscInt buf;
10036479e835SStefano Zampini 
10046479e835SStefano Zampini   PetscFunctionBegin;
10056479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
10066479e835SStefano Zampini   buf = (PetscInt)pc->failedreason;
10076479e835SStefano Zampini   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
10086479e835SStefano Zampini   pc->failedreason = (PCFailedReason)buf;
10096479e835SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
10106479e835SStefano Zampini }
10116479e835SStefano Zampini 
1012c00cb57fSBarry Smith /*  Next line needed to deactivate KSP_Solve logging */
1013c00cb57fSBarry Smith #include <petsc/private/kspimpl.h>
1014c00cb57fSBarry Smith 
10154b9ad928SBarry Smith /*
10164b9ad928SBarry Smith       a setupcall of 0 indicates never setup,
101723ee1639SBarry Smith                      1 indicates has been previously setup
1018422a814eSBarry Smith                     -1 indicates a PCSetUp() was attempted and failed
10194b9ad928SBarry Smith */
10204b9ad928SBarry Smith /*@
10214b9ad928SBarry Smith   PCSetUp - Prepares for the use of a preconditioner.
10224b9ad928SBarry Smith 
1023c3339decSBarry Smith   Collective
10244b9ad928SBarry Smith 
10254b9ad928SBarry Smith   Input Parameter:
10264b9ad928SBarry Smith . pc - the preconditioner context
10274b9ad928SBarry Smith 
10284b9ad928SBarry Smith   Level: developer
10294b9ad928SBarry Smith 
1030562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`
10314b9ad928SBarry Smith @*/
1032d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp(PC pc)
1033d71ae5a4SJacob Faibussowitsch {
1034566e8bf2SBarry Smith   const char      *def;
1035fc9b51d3SBarry Smith   PetscObjectState matstate, matnonzerostate;
10364b9ad928SBarry Smith 
10374b9ad928SBarry Smith   PetscFunctionBegin;
10380700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
103928b400f6SJacob Faibussowitsch   PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first");
10406ce5e81cSLisandro Dalcin 
104123ee1639SBarry Smith   if (pc->setupcalled && pc->reusepreconditioner) {
10429566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n"));
10433ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10444b9ad928SBarry Smith   }
10454b9ad928SBarry Smith 
10469566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate));
10479566063dSJacob Faibussowitsch   PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate));
104823ee1639SBarry Smith   if (!pc->setupcalled) {
10495e62d202SMark Adams     //PetscCall(PetscInfo(pc, "Setting up PC for first time\n"));
105023ee1639SBarry Smith     pc->flag = DIFFERENT_NONZERO_PATTERN;
10519df67409SStefano Zampini   } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS);
10529df67409SStefano Zampini   else {
10539df67409SStefano Zampini     if (matnonzerostate != pc->matnonzerostate) {
10549566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n"));
105523ee1639SBarry Smith       pc->flag = DIFFERENT_NONZERO_PATTERN;
105623ee1639SBarry Smith     } else {
10575e62d202SMark Adams       //PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n"));
105823ee1639SBarry Smith       pc->flag = SAME_NONZERO_PATTERN;
105923ee1639SBarry Smith     }
106023ee1639SBarry Smith   }
106123ee1639SBarry Smith   pc->matstate        = matstate;
1062fc9b51d3SBarry Smith   pc->matnonzerostate = matnonzerostate;
106323ee1639SBarry Smith 
10647adad957SLisandro Dalcin   if (!((PetscObject)pc)->type_name) {
10659566063dSJacob Faibussowitsch     PetscCall(PCGetDefaultType_Private(pc, &def));
10669566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, def));
1067566e8bf2SBarry Smith   }
10684b9ad928SBarry Smith 
10699566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
10709566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure));
10719566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0));
10724b9ad928SBarry Smith   if (pc->ops->setup) {
1073c00cb57fSBarry Smith     /* do not log solves and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */
10749566063dSJacob Faibussowitsch     PetscCall(KSPInitializePackage());
10759566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePush(KSP_Solve));
10769566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePush(PC_Apply));
1077dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, setup);
10789566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePop(KSP_Solve));
10799566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePop(PC_Apply));
10804b9ad928SBarry Smith   }
10819566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0));
1082422a814eSBarry Smith   if (!pc->setupcalled) pc->setupcalled = 1;
10833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10844b9ad928SBarry Smith }
10854b9ad928SBarry Smith 
10864b9ad928SBarry Smith /*@
10874b9ad928SBarry Smith   PCSetUpOnBlocks - Sets up the preconditioner for each block in
10884b9ad928SBarry Smith   the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
10894b9ad928SBarry Smith   methods.
10904b9ad928SBarry Smith 
1091c3339decSBarry Smith   Collective
10924b9ad928SBarry Smith 
1093f1580f4eSBarry Smith   Input Parameter:
10944b9ad928SBarry Smith . pc - the preconditioner context
10954b9ad928SBarry Smith 
10964b9ad928SBarry Smith   Level: developer
10974b9ad928SBarry Smith 
1098f1580f4eSBarry Smith   Note:
1099f1580f4eSBarry Smith   For nested preconditioners such as `PCBJACOBI` `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is
1100f1580f4eSBarry Smith   called on the outer `PC`, this routine ensures it is called.
1101f1580f4eSBarry Smith 
1102562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()`
11034b9ad928SBarry Smith @*/
1104d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUpOnBlocks(PC pc)
1105d71ae5a4SJacob Faibussowitsch {
11064b9ad928SBarry Smith   PetscFunctionBegin;
11070700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11083ba16761SJacob Faibussowitsch   if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS);
11099566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0));
1110dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, setuponblocks);
11119566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0));
11123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11134b9ad928SBarry Smith }
11144b9ad928SBarry Smith 
11154b9ad928SBarry Smith /*@C
11164b9ad928SBarry Smith   PCSetModifySubMatrices - Sets a user-defined routine for modifying the
111704c3f3b8SBarry Smith   submatrices that arise within certain subdomain-based preconditioners such as `PCASM`
11184b9ad928SBarry Smith 
1119c3339decSBarry Smith   Logically Collective
11204b9ad928SBarry Smith 
11214b9ad928SBarry Smith   Input Parameters:
11224b9ad928SBarry Smith + pc   - the preconditioner context
11234b9ad928SBarry Smith . func - routine for modifying the submatrices
112404c3f3b8SBarry Smith - ctx  - optional user-defined context (may be `NULL`)
11254b9ad928SBarry Smith 
112620f4b53cSBarry Smith   Calling sequence of `func`:
112720f4b53cSBarry Smith + pc     - the preconditioner context
112804c3f3b8SBarry Smith . nsub   - number of index sets
112920f4b53cSBarry Smith . row    - an array of index sets that contain the global row numbers
11304b9ad928SBarry Smith          that comprise each local submatrix
11314b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
11324b9ad928SBarry Smith          that comprise each local submatrix
11334b9ad928SBarry Smith . submat - array of local submatrices
11344b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
113504c3f3b8SBarry Smith          user-defined func routine (may be `NULL`)
11364b9ad928SBarry Smith 
113720f4b53cSBarry Smith   Level: advanced
113820f4b53cSBarry Smith 
11394b9ad928SBarry Smith   Notes:
114004c3f3b8SBarry Smith   The basic submatrices are extracted from the preconditioner matrix as
114104c3f3b8SBarry Smith   usual; the user can then alter these (for example, to set different boundary
114204c3f3b8SBarry Smith   conditions for each submatrix) before they are used for the local solves.
114304c3f3b8SBarry Smith 
1144f1580f4eSBarry Smith   `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and
1145f1580f4eSBarry Smith   `KSPSolve()`.
11464b9ad928SBarry Smith 
1147f1580f4eSBarry Smith   A routine set by `PCSetModifySubMatrices()` is currently called within
1148f1580f4eSBarry Smith   the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`)
11494b9ad928SBarry Smith   preconditioners.  All other preconditioners ignore this routine.
11504b9ad928SBarry Smith 
1151562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()`
11524b9ad928SBarry Smith @*/
115304c3f3b8SBarry Smith PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx), void *ctx)
1154d71ae5a4SJacob Faibussowitsch {
11554b9ad928SBarry Smith   PetscFunctionBegin;
11560700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11574b9ad928SBarry Smith   pc->modifysubmatrices  = func;
11584b9ad928SBarry Smith   pc->modifysubmatricesP = ctx;
11593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11604b9ad928SBarry Smith }
11614b9ad928SBarry Smith 
11624b9ad928SBarry Smith /*@C
11634b9ad928SBarry Smith   PCModifySubMatrices - Calls an optional user-defined routine within
1164f1580f4eSBarry Smith   certain preconditioners if one has been set with `PCSetModifySubMatrices()`.
11654b9ad928SBarry Smith 
1166c3339decSBarry Smith   Collective
11674b9ad928SBarry Smith 
11684b9ad928SBarry Smith   Input Parameters:
11694b9ad928SBarry Smith + pc     - the preconditioner context
11704b9ad928SBarry Smith . nsub   - the number of local submatrices
11714b9ad928SBarry Smith . row    - an array of index sets that contain the global row numbers
11724b9ad928SBarry Smith          that comprise each local submatrix
11734b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
11744b9ad928SBarry Smith          that comprise each local submatrix
11754b9ad928SBarry Smith . submat - array of local submatrices
11764b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
1177562efe2eSBarry Smith          user-defined routine (may be `NULL`)
11784b9ad928SBarry Smith 
11794b9ad928SBarry Smith   Output Parameter:
11804b9ad928SBarry Smith . submat - array of local submatrices (the entries of which may
11814b9ad928SBarry Smith             have been modified)
11824b9ad928SBarry Smith 
118320f4b53cSBarry Smith   Level: developer
118420f4b53cSBarry Smith 
118504c3f3b8SBarry Smith   Note:
11864b9ad928SBarry Smith   The user should NOT generally call this routine, as it will
118704c3f3b8SBarry Smith   automatically be called within certain preconditioners.
11884b9ad928SBarry Smith 
1189562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetModifySubMatrices()`
11904b9ad928SBarry Smith @*/
1191d71ae5a4SJacob Faibussowitsch PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx)
1192d71ae5a4SJacob Faibussowitsch {
11934b9ad928SBarry Smith   PetscFunctionBegin;
11940700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11953ba16761SJacob Faibussowitsch   if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS);
11969566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0));
11979566063dSJacob Faibussowitsch   PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx));
11989566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0));
11993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12004b9ad928SBarry Smith }
12014b9ad928SBarry Smith 
12024b9ad928SBarry Smith /*@
12034b9ad928SBarry Smith   PCSetOperators - Sets the matrix associated with the linear system and
12044b9ad928SBarry Smith   a (possibly) different one associated with the preconditioner.
12054b9ad928SBarry Smith 
1206c3339decSBarry Smith   Logically Collective
12074b9ad928SBarry Smith 
12084b9ad928SBarry Smith   Input Parameters:
12094b9ad928SBarry Smith + pc   - the preconditioner context
1210e5d3d808SBarry Smith . Amat - the matrix that defines the linear system
1211d1e9a80fSBarry Smith - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
12124b9ad928SBarry Smith 
121320f4b53cSBarry Smith   Level: intermediate
1214189c0b8aSBarry Smith 
121520f4b53cSBarry Smith   Notes:
121620f4b53cSBarry Smith   Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used.
121720f4b53cSBarry Smith 
121820f4b53cSBarry Smith   If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then
1219f1580f4eSBarry Smith   first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()`
1220f1580f4eSBarry Smith   on it and then pass it back in in your call to `KSPSetOperators()`.
1221189c0b8aSBarry Smith 
12224b9ad928SBarry Smith   More Notes about Repeated Solution of Linear Systems:
122320f4b53cSBarry Smith   PETSc does NOT reset the matrix entries of either `Amat` or `Pmat`
12244b9ad928SBarry Smith   to zero after a linear solve; the user is completely responsible for
1225f1580f4eSBarry Smith   matrix assembly.  See the routine `MatZeroEntries()` if desiring to
12264b9ad928SBarry Smith   zero all elements of a matrix.
12274b9ad928SBarry Smith 
1228562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`
12294b9ad928SBarry Smith  @*/
1230d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat)
1231d71ae5a4SJacob Faibussowitsch {
12323b3f6333SBarry Smith   PetscInt m1, n1, m2, n2;
12334b9ad928SBarry Smith 
12344b9ad928SBarry Smith   PetscFunctionBegin;
12350700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12360700a824SBarry Smith   if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
12370700a824SBarry Smith   if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3);
12383fc8bf9cSBarry Smith   if (Amat) PetscCheckSameComm(pc, 1, Amat, 2);
12393fc8bf9cSBarry Smith   if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3);
124031641f1aSBarry Smith   if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) {
12419566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Amat, &m1, &n1));
12429566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->mat, &m2, &n2));
124363a3b9bcSJacob 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);
12449566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Pmat, &m1, &n1));
12459566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2));
124663a3b9bcSJacob 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);
12473b3f6333SBarry Smith   }
12484b9ad928SBarry Smith 
124923ee1639SBarry Smith   if (Pmat != pc->pmat) {
125023ee1639SBarry Smith     /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
125123ee1639SBarry Smith     pc->matnonzerostate = -1;
125223ee1639SBarry Smith     pc->matstate        = -1;
125323ee1639SBarry Smith   }
125423ee1639SBarry Smith 
1255906ed7ccSBarry Smith   /* reference first in case the matrices are the same */
12569566063dSJacob Faibussowitsch   if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat));
12579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->mat));
12589566063dSJacob Faibussowitsch   if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat));
12599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->pmat));
12604b9ad928SBarry Smith   pc->mat  = Amat;
12614b9ad928SBarry Smith   pc->pmat = Pmat;
12623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1263e56f5c9eSBarry Smith }
1264e56f5c9eSBarry Smith 
126523ee1639SBarry Smith /*@
126623ee1639SBarry Smith   PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed.
126723ee1639SBarry Smith 
1268c3339decSBarry Smith   Logically Collective
126923ee1639SBarry Smith 
127023ee1639SBarry Smith   Input Parameters:
127123ee1639SBarry Smith + pc   - the preconditioner context
1272f1580f4eSBarry Smith - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
127323ee1639SBarry Smith 
12742b26979fSBarry Smith   Level: intermediate
12752b26979fSBarry Smith 
1276f1580f4eSBarry Smith   Note:
1277f1580f4eSBarry Smith   Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option
1278f1580f4eSBarry Smith   prevents this.
1279f1580f4eSBarry Smith 
1280562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()`
128123ee1639SBarry Smith  @*/
1282d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag)
1283d71ae5a4SJacob Faibussowitsch {
128423ee1639SBarry Smith   PetscFunctionBegin;
128523ee1639SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1286f9177522SStefano Zampini   PetscValidLogicalCollectiveBool(pc, flag, 2);
128723ee1639SBarry Smith   pc->reusepreconditioner = flag;
12883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12894b9ad928SBarry Smith }
12904b9ad928SBarry Smith 
1291c60c7ad4SBarry Smith /*@
1292f1580f4eSBarry Smith   PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed.
1293c60c7ad4SBarry Smith 
1294bba28a21SBarry Smith   Not Collective
1295c60c7ad4SBarry Smith 
1296c60c7ad4SBarry Smith   Input Parameter:
1297c60c7ad4SBarry Smith . pc - the preconditioner context
1298c60c7ad4SBarry Smith 
1299c60c7ad4SBarry Smith   Output Parameter:
1300f1580f4eSBarry Smith . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
1301c60c7ad4SBarry Smith 
1302d0418729SSatish Balay   Level: intermediate
1303d0418729SSatish Balay 
1304562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()`
1305c60c7ad4SBarry Smith  @*/
1306d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag)
1307d71ae5a4SJacob Faibussowitsch {
1308c60c7ad4SBarry Smith   PetscFunctionBegin;
1309c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
13104f572ea9SToby Isaac   PetscAssertPointer(flag, 2);
1311c60c7ad4SBarry Smith   *flag = pc->reusepreconditioner;
13123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1313c60c7ad4SBarry Smith }
1314c60c7ad4SBarry Smith 
1315487a658cSBarry Smith /*@
13164b9ad928SBarry Smith   PCGetOperators - Gets the matrix associated with the linear system and
13174b9ad928SBarry Smith   possibly a different one associated with the preconditioner.
13184b9ad928SBarry Smith 
1319562efe2eSBarry Smith   Not Collective, though parallel `Mat`s are returned if `pc` is parallel
13204b9ad928SBarry Smith 
13214b9ad928SBarry Smith   Input Parameter:
13224b9ad928SBarry Smith . pc - the preconditioner context
13234b9ad928SBarry Smith 
13244b9ad928SBarry Smith   Output Parameters:
1325e5d3d808SBarry Smith + Amat - the matrix defining the linear system
132623ee1639SBarry Smith - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
13274b9ad928SBarry Smith 
13284b9ad928SBarry Smith   Level: intermediate
13294b9ad928SBarry Smith 
1330f1580f4eSBarry Smith   Note:
133195452b02SPatrick Sanan   Does not increase the reference count of the matrices, so you should not destroy them
1332298cc208SBarry Smith 
1333f1580f4eSBarry Smith   Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators
1334f1580f4eSBarry Smith   are created in `PC` and returned to the user. In this case, if both operators
133573950996SBarry Smith   mat and pmat are requested, two DIFFERENT operators will be returned. If
133673950996SBarry Smith   only one is requested both operators in the PC will be the same (i.e. as
1337f1580f4eSBarry Smith   if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats).
133873950996SBarry Smith   The user must set the sizes of the returned matrices and their type etc just
1339f1580f4eSBarry Smith   as if the user created them with `MatCreate()`. For example,
134073950996SBarry Smith 
1341f1580f4eSBarry Smith .vb
1342f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1343f1580f4eSBarry Smith            set size, type, etc of Amat
134473950996SBarry Smith 
1345f1580f4eSBarry Smith          MatCreate(comm,&mat);
1346f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Amat);
1347f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)mat);
1348f1580f4eSBarry Smith            set size, type, etc of Amat
1349f1580f4eSBarry Smith .ve
135073950996SBarry Smith 
135173950996SBarry Smith   and
135273950996SBarry Smith 
1353f1580f4eSBarry Smith .vb
1354f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1355f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
135673950996SBarry Smith 
1357f1580f4eSBarry Smith          MatCreate(comm,&Amat);
1358f1580f4eSBarry Smith          MatCreate(comm,&Pmat);
1359f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1360f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Amat);
1361f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Pmat);
1362f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
1363f1580f4eSBarry Smith .ve
136473950996SBarry Smith 
1365f1580f4eSBarry Smith   The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy
1366b8d709abSRichard Tran Mills   of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely
1367f1580f4eSBarry Smith   managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look
1368f1580f4eSBarry Smith   at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to
1369f1580f4eSBarry Smith   the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP
1370f1580f4eSBarry Smith   you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you).
1371f1580f4eSBarry Smith   Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when
137273950996SBarry Smith   it can be created for you?
137373950996SBarry Smith 
1374562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()`
13754b9ad928SBarry Smith @*/
1376d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat)
1377d71ae5a4SJacob Faibussowitsch {
13784b9ad928SBarry Smith   PetscFunctionBegin;
13790700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1380e5d3d808SBarry Smith   if (Amat) {
138173950996SBarry Smith     if (!pc->mat) {
13829fca8c71SStefano Zampini       if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */
13839a4708feSJed Brown         pc->mat = pc->pmat;
13849566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->mat));
1385e5d3d808SBarry Smith       } else { /* both Amat and Pmat are empty */
13869566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat));
1387e5d3d808SBarry Smith         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
138873950996SBarry Smith           pc->pmat = pc->mat;
13899566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->pmat));
139073950996SBarry Smith         }
139173950996SBarry Smith       }
13929a4708feSJed Brown     }
1393e5d3d808SBarry Smith     *Amat = pc->mat;
139473950996SBarry Smith   }
1395e5d3d808SBarry Smith   if (Pmat) {
139673950996SBarry Smith     if (!pc->pmat) {
1397e5d3d808SBarry Smith       if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */
13989a4708feSJed Brown         pc->pmat = pc->mat;
13999566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->pmat));
14009a4708feSJed Brown       } else {
14019566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat));
1402e5d3d808SBarry Smith         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
140373950996SBarry Smith           pc->mat = pc->pmat;
14049566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->mat));
140573950996SBarry Smith         }
140673950996SBarry Smith       }
14079a4708feSJed Brown     }
1408e5d3d808SBarry Smith     *Pmat = pc->pmat;
140973950996SBarry Smith   }
14103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14114b9ad928SBarry Smith }
14124b9ad928SBarry Smith 
1413*5d83a8b1SBarry Smith /*@
1414906ed7ccSBarry Smith   PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1415f1580f4eSBarry Smith   possibly a different one associated with the preconditioner have been set in the `PC`.
1416906ed7ccSBarry Smith 
141720f4b53cSBarry Smith   Not Collective, though the results on all processes should be the same
1418906ed7ccSBarry Smith 
1419906ed7ccSBarry Smith   Input Parameter:
1420906ed7ccSBarry Smith . pc - the preconditioner context
1421906ed7ccSBarry Smith 
1422906ed7ccSBarry Smith   Output Parameters:
1423906ed7ccSBarry Smith + mat  - the matrix associated with the linear system was set
1424906ed7ccSBarry Smith - pmat - matrix associated with the preconditioner was set, usually the same
1425906ed7ccSBarry Smith 
1426906ed7ccSBarry Smith   Level: intermediate
1427906ed7ccSBarry Smith 
1428562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`
1429906ed7ccSBarry Smith @*/
1430d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat)
1431d71ae5a4SJacob Faibussowitsch {
1432906ed7ccSBarry Smith   PetscFunctionBegin;
14330700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1434906ed7ccSBarry Smith   if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1435906ed7ccSBarry Smith   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
14363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1437906ed7ccSBarry Smith }
1438906ed7ccSBarry Smith 
1439f39d8e23SSatish Balay /*@
1440a4fd02acSBarry Smith   PCFactorGetMatrix - Gets the factored matrix from the
1441f1580f4eSBarry Smith   preconditioner context.  This routine is valid only for the `PCLU`,
1442f1580f4eSBarry Smith   `PCILU`, `PCCHOLESKY`, and `PCICC` methods.
14434b9ad928SBarry Smith 
1444562efe2eSBarry Smith   Not Collective though `mat` is parallel if `pc` is parallel
14454b9ad928SBarry Smith 
1446f1580f4eSBarry Smith   Input Parameter:
14474b9ad928SBarry Smith . pc - the preconditioner context
14484b9ad928SBarry Smith 
1449feefa0e1SJacob Faibussowitsch   Output Parameters:
14504b9ad928SBarry Smith . mat - the factored matrix
14514b9ad928SBarry Smith 
14524b9ad928SBarry Smith   Level: advanced
14534b9ad928SBarry Smith 
1454f1580f4eSBarry Smith   Note:
1455562efe2eSBarry Smith   Does not increase the reference count for `mat` so DO NOT destroy it
14569405f653SBarry Smith 
1457562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC`
14584b9ad928SBarry Smith @*/
1459d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat)
1460d71ae5a4SJacob Faibussowitsch {
14614b9ad928SBarry Smith   PetscFunctionBegin;
14620700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14634f572ea9SToby Isaac   PetscAssertPointer(mat, 2);
14647def7855SStefano Zampini   PetscCall(PCFactorSetUpMatSolverType(pc));
1465dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, getfactoredmatrix, mat);
14663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14674b9ad928SBarry Smith }
14684b9ad928SBarry Smith 
1469cc4c1da9SBarry Smith /*@
14704b9ad928SBarry Smith   PCSetOptionsPrefix - Sets the prefix used for searching for all
1471f1580f4eSBarry Smith   `PC` options in the database.
14724b9ad928SBarry Smith 
1473c3339decSBarry Smith   Logically Collective
14744b9ad928SBarry Smith 
14754b9ad928SBarry Smith   Input Parameters:
14764b9ad928SBarry Smith + pc     - the preconditioner context
1477f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
14784b9ad928SBarry Smith 
1479f1580f4eSBarry Smith   Note:
14804b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
14814b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
14824b9ad928SBarry Smith   hyphen.
14834b9ad928SBarry Smith 
14844b9ad928SBarry Smith   Level: advanced
14854b9ad928SBarry Smith 
1486562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()`
14874b9ad928SBarry Smith @*/
1488d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[])
1489d71ae5a4SJacob Faibussowitsch {
14904b9ad928SBarry Smith   PetscFunctionBegin;
14910700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14929566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix));
14933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14944b9ad928SBarry Smith }
14954b9ad928SBarry Smith 
1496cc4c1da9SBarry Smith /*@
14974b9ad928SBarry Smith   PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1498f1580f4eSBarry Smith   `PC` options in the database.
14994b9ad928SBarry Smith 
1500c3339decSBarry Smith   Logically Collective
15014b9ad928SBarry Smith 
15024b9ad928SBarry Smith   Input Parameters:
15034b9ad928SBarry Smith + pc     - the preconditioner context
1504f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
15054b9ad928SBarry Smith 
1506f1580f4eSBarry Smith   Note:
15074b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
15084b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
15094b9ad928SBarry Smith   hyphen.
15104b9ad928SBarry Smith 
15114b9ad928SBarry Smith   Level: advanced
15124b9ad928SBarry Smith 
1513562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()`
15144b9ad928SBarry Smith @*/
1515d71ae5a4SJacob Faibussowitsch PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[])
1516d71ae5a4SJacob Faibussowitsch {
15174b9ad928SBarry Smith   PetscFunctionBegin;
15180700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15199566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix));
15203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15214b9ad928SBarry Smith }
15224b9ad928SBarry Smith 
1523cc4c1da9SBarry Smith /*@
15244b9ad928SBarry Smith   PCGetOptionsPrefix - Gets the prefix used for searching for all
15254b9ad928SBarry Smith   PC options in the database.
15264b9ad928SBarry Smith 
15274b9ad928SBarry Smith   Not Collective
15284b9ad928SBarry Smith 
1529f1580f4eSBarry Smith   Input Parameter:
15304b9ad928SBarry Smith . pc - the preconditioner context
15314b9ad928SBarry Smith 
1532f1580f4eSBarry Smith   Output Parameter:
15334b9ad928SBarry Smith . prefix - pointer to the prefix string used, is returned
15344b9ad928SBarry Smith 
15354b9ad928SBarry Smith   Level: advanced
15364b9ad928SBarry Smith 
1537562efe2eSBarry Smith   Fortran Note:
1538562efe2eSBarry Smith   The user should pass in a string `prefix` of
1539562efe2eSBarry Smith   sufficient length to hold the prefix.
1540562efe2eSBarry Smith 
1541562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()`
15424b9ad928SBarry Smith @*/
1543d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[])
1544d71ae5a4SJacob Faibussowitsch {
15454b9ad928SBarry Smith   PetscFunctionBegin;
15460700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15474f572ea9SToby Isaac   PetscAssertPointer(prefix, 2);
15489566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix));
15493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15504b9ad928SBarry Smith }
15514b9ad928SBarry Smith 
15528066bbecSBarry Smith /*
1553dd8e379bSPierre Jolivet    Indicates the right-hand side will be changed by KSPSolve(), this occurs for a few
15548066bbecSBarry Smith   preconditioners including BDDC and Eisentat that transform the equations before applying
15558066bbecSBarry Smith   the Krylov methods
15568066bbecSBarry Smith */
1557d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change)
1558d71ae5a4SJacob Faibussowitsch {
1559a06fd7f2SStefano Zampini   PetscFunctionBegin;
1560a06fd7f2SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15614f572ea9SToby Isaac   PetscAssertPointer(change, 2);
1562a06fd7f2SStefano Zampini   *change = PETSC_FALSE;
1563cac4c232SBarry Smith   PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change));
15643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1565a06fd7f2SStefano Zampini }
1566a06fd7f2SStefano Zampini 
15674b9ad928SBarry Smith /*@
15684b9ad928SBarry Smith   PCPreSolve - Optional pre-solve phase, intended for any
15694b9ad928SBarry Smith   preconditioner-specific actions that must be performed before
15704b9ad928SBarry Smith   the iterative solve itself.
15714b9ad928SBarry Smith 
1572c3339decSBarry Smith   Collective
15734b9ad928SBarry Smith 
15744b9ad928SBarry Smith   Input Parameters:
15754b9ad928SBarry Smith + pc  - the preconditioner context
15764b9ad928SBarry Smith - ksp - the Krylov subspace context
15774b9ad928SBarry Smith 
15784b9ad928SBarry Smith   Level: developer
15794b9ad928SBarry Smith 
1580feefa0e1SJacob Faibussowitsch   Example Usage:
15814b9ad928SBarry Smith .vb
15824b9ad928SBarry Smith     PCPreSolve(pc,ksp);
158323ce1328SBarry Smith     KSPSolve(ksp,b,x);
15844b9ad928SBarry Smith     PCPostSolve(pc,ksp);
15854b9ad928SBarry Smith .ve
15864b9ad928SBarry Smith 
15874b9ad928SBarry Smith   Notes:
1588f1580f4eSBarry Smith   The pre-solve phase is distinct from the `PCSetUp()` phase.
15894b9ad928SBarry Smith 
1590f1580f4eSBarry Smith   `KSPSolve()` calls this directly, so is rarely called by the user.
15914b9ad928SBarry Smith 
1592562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCPostSolve()`
15934b9ad928SBarry Smith @*/
1594d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPreSolve(PC pc, KSP ksp)
1595d71ae5a4SJacob Faibussowitsch {
15964b9ad928SBarry Smith   Vec x, rhs;
15974b9ad928SBarry Smith 
15984b9ad928SBarry Smith   PetscFunctionBegin;
15990700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16000700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1601d9a03883SBarry Smith   pc->presolvedone++;
16027827d75bSBarry Smith   PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice");
16039566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
16049566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
16054b9ad928SBarry Smith 
1606dbbe0bcdSBarry Smith   if (pc->ops->presolve) PetscUseTypeMethod(pc, presolve, ksp, rhs, x);
1607f4f49eeaSPierre Jolivet   else if (pc->presolve) PetscCall(pc->presolve(pc, ksp));
16083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16094b9ad928SBarry Smith }
16104b9ad928SBarry Smith 
1611f560b561SHong Zhang /*@C
1612f1580f4eSBarry Smith   PCSetPreSolve - Sets function used by `PCPreSolve()` which is intended for any
1613f560b561SHong Zhang   preconditioner-specific actions that must be performed before
1614f560b561SHong Zhang   the iterative solve itself.
1615f560b561SHong Zhang 
1616c3339decSBarry Smith   Logically Collective
1617f560b561SHong Zhang 
1618f560b561SHong Zhang   Input Parameters:
1619f560b561SHong Zhang + pc       - the preconditioner object
1620f560b561SHong Zhang - presolve - the function to call before the solve
1621f560b561SHong Zhang 
162220f4b53cSBarry Smith   Calling sequence of `presolve`:
162320f4b53cSBarry Smith + pc  - the `PC` context
162420f4b53cSBarry Smith - ksp - the `KSP` context
1625f560b561SHong Zhang 
1626f560b561SHong Zhang   Level: developer
1627f560b561SHong Zhang 
1628562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCPreSolve()`
1629f560b561SHong Zhang @*/
163004c3f3b8SBarry Smith PetscErrorCode PCSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC pc, KSP ksp))
1631d71ae5a4SJacob Faibussowitsch {
1632f560b561SHong Zhang   PetscFunctionBegin;
1633f560b561SHong Zhang   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1634f560b561SHong Zhang   pc->presolve = presolve;
16353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1636f560b561SHong Zhang }
1637f560b561SHong Zhang 
16384b9ad928SBarry Smith /*@
16394b9ad928SBarry Smith   PCPostSolve - Optional post-solve phase, intended for any
16404b9ad928SBarry Smith   preconditioner-specific actions that must be performed after
16414b9ad928SBarry Smith   the iterative solve itself.
16424b9ad928SBarry Smith 
1643c3339decSBarry Smith   Collective
16444b9ad928SBarry Smith 
16454b9ad928SBarry Smith   Input Parameters:
16464b9ad928SBarry Smith + pc  - the preconditioner context
16474b9ad928SBarry Smith - ksp - the Krylov subspace context
16484b9ad928SBarry Smith 
1649feefa0e1SJacob Faibussowitsch   Example Usage:
16504b9ad928SBarry Smith .vb
16514b9ad928SBarry Smith     PCPreSolve(pc,ksp);
165223ce1328SBarry Smith     KSPSolve(ksp,b,x);
16534b9ad928SBarry Smith     PCPostSolve(pc,ksp);
16544b9ad928SBarry Smith .ve
16554b9ad928SBarry Smith 
1656562efe2eSBarry Smith   Level: developer
1657562efe2eSBarry Smith 
16584b9ad928SBarry Smith   Note:
1659f1580f4eSBarry Smith   `KSPSolve()` calls this routine directly, so it is rarely called by the user.
16604b9ad928SBarry Smith 
1661562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetPostSolve()`, `PCSetPresolve()`, `PCPreSolve()`, `KSPSolve()`
16624b9ad928SBarry Smith @*/
1663d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPostSolve(PC pc, KSP ksp)
1664d71ae5a4SJacob Faibussowitsch {
16654b9ad928SBarry Smith   Vec x, rhs;
16664b9ad928SBarry Smith 
16674b9ad928SBarry Smith   PetscFunctionBegin;
16680700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16690700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1670d9a03883SBarry Smith   pc->presolvedone--;
16719566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
16729566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
1673dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, postsolve, ksp, rhs, x);
16743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16754b9ad928SBarry Smith }
16764b9ad928SBarry Smith 
1677ffeef943SBarry Smith /*@
1678f1580f4eSBarry Smith   PCLoad - Loads a `PC` that has been stored in binary  with `PCView()`.
167955849f57SBarry Smith 
1680c3339decSBarry Smith   Collective
168155849f57SBarry Smith 
168255849f57SBarry Smith   Input Parameters:
1683f1580f4eSBarry Smith + newdm  - the newly loaded `PC`, this needs to have been created with `PCCreate()` or
1684f1580f4eSBarry Smith            some related function before a call to `PCLoad()`.
1685f1580f4eSBarry Smith - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()`
168655849f57SBarry Smith 
168755849f57SBarry Smith   Level: intermediate
168855849f57SBarry Smith 
1689f1580f4eSBarry Smith   Note:
1690562efe2eSBarry Smith   The type is determined by the data in the file, any `PCType` set into the `PC` before this call is ignored.
169155849f57SBarry Smith 
1692562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()`
169355849f57SBarry Smith @*/
1694d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLoad(PC newdm, PetscViewer viewer)
1695d71ae5a4SJacob Faibussowitsch {
169655849f57SBarry Smith   PetscBool isbinary;
1697060da220SMatthew G. Knepley   PetscInt  classid;
169855849f57SBarry Smith   char      type[256];
169955849f57SBarry Smith 
170055849f57SBarry Smith   PetscFunctionBegin;
170155849f57SBarry Smith   PetscValidHeaderSpecific(newdm, PC_CLASSID, 1);
170255849f57SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
17039566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
170428b400f6SJacob Faibussowitsch   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
170555849f57SBarry Smith 
17069566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
17077827d75bSBarry Smith   PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file");
17089566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
17099566063dSJacob Faibussowitsch   PetscCall(PCSetType(newdm, type));
1710dbbe0bcdSBarry Smith   PetscTryTypeMethod(newdm, load, viewer);
17113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
171255849f57SBarry Smith }
171355849f57SBarry Smith 
17149804daf3SBarry Smith #include <petscdraw.h>
1715e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1716e04113cfSBarry Smith   #include <petscviewersaws.h>
17170acecf5bSBarry Smith #endif
1718fe2efc57SMark 
1719ffeef943SBarry Smith /*@
1720562efe2eSBarry Smith   PCViewFromOptions - View from the `PC` based on options in the options database
1721fe2efc57SMark 
1722c3339decSBarry Smith   Collective
1723fe2efc57SMark 
1724fe2efc57SMark   Input Parameters:
172520f4b53cSBarry Smith + A    - the `PC` context
1726f1580f4eSBarry Smith . obj  - Optional object that provides the options prefix
1727736c3998SJose E. Roman - name - command line option
1728fe2efc57SMark 
1729fe2efc57SMark   Level: intermediate
1730f1580f4eSBarry Smith 
1731562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()`
1732fe2efc57SMark @*/
1733d71ae5a4SJacob Faibussowitsch PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[])
1734d71ae5a4SJacob Faibussowitsch {
1735fe2efc57SMark   PetscFunctionBegin;
1736fe2efc57SMark   PetscValidHeaderSpecific(A, PC_CLASSID, 1);
17379566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
17383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1739fe2efc57SMark }
1740fe2efc57SMark 
1741ffeef943SBarry Smith /*@
1742f1580f4eSBarry Smith   PCView - Prints information about the `PC`
17434b9ad928SBarry Smith 
1744c3339decSBarry Smith   Collective
17454b9ad928SBarry Smith 
17464b9ad928SBarry Smith   Input Parameters:
1747feefa0e1SJacob Faibussowitsch + pc     - the `PC` context
17484b9ad928SBarry Smith - viewer - optional visualization context
17494b9ad928SBarry Smith 
175020f4b53cSBarry Smith   Level: developer
175120f4b53cSBarry Smith 
1752f1580f4eSBarry Smith   Notes:
17534b9ad928SBarry Smith   The available visualization contexts include
1754f1580f4eSBarry Smith +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
1755f1580f4eSBarry Smith -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
17564b9ad928SBarry Smith   output where only the first processor opens
17574b9ad928SBarry Smith   the file. All other processors send their
17584b9ad928SBarry Smith   data to the first processor to print.
17594b9ad928SBarry Smith 
17604b9ad928SBarry Smith   The user can open an alternative visualization contexts with
1761f1580f4eSBarry Smith   `PetscViewerASCIIOpen()` (output to a specified file).
17624b9ad928SBarry Smith 
1763562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewer`, `KSPView()`, `PetscViewerASCIIOpen()`
17644b9ad928SBarry Smith @*/
1765d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView(PC pc, PetscViewer viewer)
1766d71ae5a4SJacob Faibussowitsch {
176719fd82e9SBarry Smith   PCType            cstr;
17686cd81132SPierre Jolivet   PetscViewerFormat format;
17696cd81132SPierre Jolivet   PetscBool         iascii, isstring, isbinary, isdraw, pop = PETSC_FALSE;
1770e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1771536b137fSBarry Smith   PetscBool issaws;
17720acecf5bSBarry Smith #endif
17734b9ad928SBarry Smith 
17744b9ad928SBarry Smith   PetscFunctionBegin;
17750700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
177648a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
17770700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1778c9780b6fSBarry Smith   PetscCheckSameComm(pc, 1, viewer, 2);
17794b9ad928SBarry Smith 
17809566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
17819566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
17829566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
17839566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1784e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
17859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
17860acecf5bSBarry Smith #endif
1787219b1045SBarry Smith 
178832077d6dSBarry Smith   if (iascii) {
17899566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer));
179048a46eb9SPierre Jolivet     if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, "  PC has not been set up so information may be incomplete\n"));
17919566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1792dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
17939566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
1794834dbeb0SBarry Smith     if (pc->mat) {
17956cd81132SPierre Jolivet       PetscCall(PetscViewerGetFormat(viewer, &format));
17966cd81132SPierre Jolivet       if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
17979566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
17986cd81132SPierre Jolivet         pop = PETSC_TRUE;
17996cd81132SPierre Jolivet       }
18004b9ad928SBarry Smith       if (pc->pmat == pc->mat) {
18019566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix = precond matrix:\n"));
18029566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
18039566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
18049566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
18054b9ad928SBarry Smith       } else {
1806834dbeb0SBarry Smith         if (pc->pmat) {
18079566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix followed by preconditioner matrix:\n"));
18084b9ad928SBarry Smith         } else {
18099566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix:\n"));
18104b9ad928SBarry Smith         }
18119566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
18129566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
18139566063dSJacob Faibussowitsch         if (pc->pmat) PetscCall(MatView(pc->pmat, viewer));
18149566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
18154b9ad928SBarry Smith       }
18166cd81132SPierre Jolivet       if (pop) PetscCall(PetscViewerPopFormat(viewer));
18174b9ad928SBarry Smith     }
18184b9ad928SBarry Smith   } else if (isstring) {
18199566063dSJacob Faibussowitsch     PetscCall(PCGetType(pc, &cstr));
18209566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr));
1821dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18229566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18239566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18245b0b0462SBarry Smith   } else if (isbinary) {
182555849f57SBarry Smith     PetscInt    classid = PC_FILE_CLASSID;
182655849f57SBarry Smith     MPI_Comm    comm;
182755849f57SBarry Smith     PetscMPIInt rank;
182855849f57SBarry Smith     char        type[256];
182955849f57SBarry Smith 
18309566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
18319566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1832dd400576SPatrick Sanan     if (rank == 0) {
18339566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
18349566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256));
18359566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
183655849f57SBarry Smith     }
1837dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
1838219b1045SBarry Smith   } else if (isdraw) {
1839219b1045SBarry Smith     PetscDraw draw;
1840d9884438SBarry Smith     char      str[25];
184189fd9fafSBarry Smith     PetscReal x, y, bottom, h;
1842d9884438SBarry Smith     PetscInt  n;
1843219b1045SBarry Smith 
18449566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
18459566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
18461d840656SPeter Brune     if (pc->mat) {
18479566063dSJacob Faibussowitsch       PetscCall(MatGetSize(pc->mat, &n, NULL));
184863a3b9bcSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n));
18491d840656SPeter Brune     } else {
18509566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name));
18511d840656SPeter Brune     }
18529566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
185389fd9fafSBarry Smith     bottom = y - h;
18549566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
1855dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18569566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
1857e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1858536b137fSBarry Smith   } else if (issaws) {
1859d45a07a7SBarry Smith     PetscMPIInt rank;
1860d45a07a7SBarry Smith 
18619566063dSJacob Faibussowitsch     PetscCall(PetscObjectName((PetscObject)pc));
18629566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
186348a46eb9SPierre Jolivet     if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer));
18649566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18659566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18660acecf5bSBarry Smith #endif
18674b9ad928SBarry Smith   }
18683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18694b9ad928SBarry Smith }
18704b9ad928SBarry Smith 
18714b9ad928SBarry Smith /*@C
1872562efe2eSBarry Smith   PCRegister -  Adds a method (`PCType`) to the preconditioner package.
18731c84c290SBarry Smith 
1874cc4c1da9SBarry Smith   Not collective. No Fortran Support
18751c84c290SBarry Smith 
18761c84c290SBarry Smith   Input Parameters:
187720f4b53cSBarry Smith + sname    - name of a new user-defined solver
187820f4b53cSBarry Smith - function - routine to create method context
18791c84c290SBarry Smith 
1880feefa0e1SJacob Faibussowitsch   Example Usage:
18811c84c290SBarry Smith .vb
1882bdf89e91SBarry Smith    PCRegister("my_solver", MySolverCreate);
18831c84c290SBarry Smith .ve
18841c84c290SBarry Smith 
18851c84c290SBarry Smith   Then, your solver can be chosen with the procedural interface via
18861c84c290SBarry Smith $     PCSetType(pc, "my_solver")
18871c84c290SBarry Smith   or at runtime via the option
18881c84c290SBarry Smith $     -pc_type my_solver
18894b9ad928SBarry Smith 
18904b9ad928SBarry Smith   Level: advanced
18911c84c290SBarry Smith 
189220f4b53cSBarry Smith   Note:
189320f4b53cSBarry Smith   `PCRegister()` may be called multiple times to add several user-defined preconditioners.
189420f4b53cSBarry Smith 
1895562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCType`, `PCRegisterAll()`
18964b9ad928SBarry Smith @*/
1897d71ae5a4SJacob Faibussowitsch PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC))
1898d71ae5a4SJacob Faibussowitsch {
18994b9ad928SBarry Smith   PetscFunctionBegin;
19009566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
19019566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCList, sname, function));
19023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19034b9ad928SBarry Smith }
19044b9ad928SBarry Smith 
1905d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y)
1906d71ae5a4SJacob Faibussowitsch {
1907186a3e20SStefano Zampini   PC pc;
1908186a3e20SStefano Zampini 
1909186a3e20SStefano Zampini   PetscFunctionBegin;
19109566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &pc));
19119566063dSJacob Faibussowitsch   PetscCall(PCApply(pc, X, Y));
19123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1913186a3e20SStefano Zampini }
1914186a3e20SStefano Zampini 
1915*5d83a8b1SBarry Smith /*@
19160bacdadaSStefano Zampini   PCComputeOperator - Computes the explicit preconditioned operator.
19174b9ad928SBarry Smith 
1918c3339decSBarry Smith   Collective
19194b9ad928SBarry Smith 
1920d8d19677SJose E. Roman   Input Parameters:
1921186a3e20SStefano Zampini + pc      - the preconditioner object
1922562efe2eSBarry Smith - mattype - the `MatType` to be used for the operator
19234b9ad928SBarry Smith 
19244b9ad928SBarry Smith   Output Parameter:
1925a5b23f4aSJose E. Roman . mat - the explicit preconditioned operator
19264b9ad928SBarry Smith 
192720f4b53cSBarry Smith   Level: advanced
192820f4b53cSBarry Smith 
1929f1580f4eSBarry Smith   Note:
1930186a3e20SStefano Zampini   This computation is done by applying the operators to columns of the identity matrix.
1931186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
1932562efe2eSBarry Smith   Currently, this routine uses a dense matrix format when `mattype` == `NULL`
19334b9ad928SBarry Smith 
1934562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPComputeOperator()`, `MatType`
19354b9ad928SBarry Smith @*/
1936d71ae5a4SJacob Faibussowitsch PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat)
1937d71ae5a4SJacob Faibussowitsch {
1938186a3e20SStefano Zampini   PetscInt N, M, m, n;
1939186a3e20SStefano Zampini   Mat      A, Apc;
19404b9ad928SBarry Smith 
19414b9ad928SBarry Smith   PetscFunctionBegin;
19420700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
19434f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
19449566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, &A, NULL));
19459566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
19469566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
19479566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc));
19489566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (void (*)(void))MatMult_PC));
19499566063dSJacob Faibussowitsch   PetscCall(MatComputeOperator(Apc, mattype, mat));
19509566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Apc));
19513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19524b9ad928SBarry Smith }
19534b9ad928SBarry Smith 
19546c237d78SBarry Smith /*@
19556c237d78SBarry Smith   PCSetCoordinates - sets the coordinates of all the nodes on the local process
19566c237d78SBarry Smith 
1957c3339decSBarry Smith   Collective
19586c237d78SBarry Smith 
19596c237d78SBarry Smith   Input Parameters:
19606c237d78SBarry Smith + pc     - the solver context
19616c237d78SBarry Smith . dim    - the dimension of the coordinates 1, 2, or 3
196214893cbeSStefano Zampini . nloc   - the blocked size of the coordinates array
196314893cbeSStefano Zampini - coords - the coordinates array
19646c237d78SBarry Smith 
19656c237d78SBarry Smith   Level: intermediate
19666c237d78SBarry Smith 
1967f1580f4eSBarry Smith   Note:
196820f4b53cSBarry Smith   `coords` is an array of the dim coordinates for the nodes on
196920f4b53cSBarry Smith   the local processor, of size `dim`*`nloc`.
197014893cbeSStefano Zampini   If there are 108 equation on a processor
19716c237d78SBarry Smith   for a displacement finite element discretization of elasticity (so
197214893cbeSStefano Zampini   that there are nloc = 36 = 108/3 nodes) then the array must have 108
19736c237d78SBarry Smith   double precision values (ie, 3 * 36).  These x y z coordinates
19746c237d78SBarry Smith   should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
19756c237d78SBarry Smith   ... , N-1.z ].
19766c237d78SBarry Smith 
1977562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `MatSetNearNullSpace()`
19786c237d78SBarry Smith @*/
1979d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[])
1980d71ae5a4SJacob Faibussowitsch {
19816c237d78SBarry Smith   PetscFunctionBegin;
198232954da3SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
198332954da3SStefano Zampini   PetscValidLogicalCollectiveInt(pc, dim, 2);
198422794d57SStefano Zampini   PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords));
19853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19866c237d78SBarry Smith }
1987fd2dd295SFande Kong 
1988fd2dd295SFande Kong /*@
1989fd2dd295SFande Kong   PCGetInterpolations - Gets interpolation matrices for all levels (except level 0)
1990fd2dd295SFande Kong 
1991c3339decSBarry Smith   Logically Collective
1992fd2dd295SFande Kong 
1993d8d19677SJose E. Roman   Input Parameter:
1994d8d19677SJose E. Roman . pc - the precondition context
1995fd2dd295SFande Kong 
1996d8d19677SJose E. Roman   Output Parameters:
1997d8d19677SJose E. Roman + num_levels     - the number of levels
1998562efe2eSBarry Smith - interpolations - the interpolation matrices (size of `num_levels`-1)
1999fd2dd295SFande Kong 
2000fd2dd295SFande Kong   Level: advanced
2001fd2dd295SFande Kong 
2002562efe2eSBarry Smith   Developer Note:
2003f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
2004fd2dd295SFande Kong 
2005562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()`
2006fd2dd295SFande Kong @*/
2007d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[])
2008d71ae5a4SJacob Faibussowitsch {
2009fd2dd295SFande Kong   PetscFunctionBegin;
2010fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20114f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20124f572ea9SToby Isaac   PetscAssertPointer(interpolations, 3);
2013cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations));
20143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2015fd2dd295SFande Kong }
2016fd2dd295SFande Kong 
2017fd2dd295SFande Kong /*@
2018fd2dd295SFande Kong   PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level)
2019fd2dd295SFande Kong 
2020c3339decSBarry Smith   Logically Collective
2021fd2dd295SFande Kong 
2022d8d19677SJose E. Roman   Input Parameter:
2023d8d19677SJose E. Roman . pc - the precondition context
2024fd2dd295SFande Kong 
2025d8d19677SJose E. Roman   Output Parameters:
2026d8d19677SJose E. Roman + num_levels      - the number of levels
2027562efe2eSBarry Smith - coarseOperators - the coarse operator matrices (size of `num_levels`-1)
2028fd2dd295SFande Kong 
2029fd2dd295SFande Kong   Level: advanced
2030fd2dd295SFande Kong 
2031562efe2eSBarry Smith   Developer Note:
2032f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
2033fd2dd295SFande Kong 
2034562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()`
2035fd2dd295SFande Kong @*/
2036d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
2037d71ae5a4SJacob Faibussowitsch {
2038fd2dd295SFande Kong   PetscFunctionBegin;
2039fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20404f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20414f572ea9SToby Isaac   PetscAssertPointer(coarseOperators, 3);
2042cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators));
20433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2044fd2dd295SFande Kong }
2045