xref: /petsc/src/ksp/pc/interface/precon.c (revision 562efe2ef922487c6beae96ba39e19afd4eefbe6)
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 
14d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetDefaultType_Private(PC pc, const char *type[])
15d71ae5a4SJacob Faibussowitsch {
162e70eadcSBarry Smith   PetscMPIInt size;
17da74c943SJose E. Roman   PetscBool   hasop, 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) {
229566063dSJacob Faibussowitsch     PetscCall(MatHasOperation(pc->pmat, MATOP_GET_DIAGONAL_BLOCK, &hasop));
23566e8bf2SBarry Smith     if (size == 1) {
249566063dSJacob Faibussowitsch       PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ICC, &flg1));
259566063dSJacob Faibussowitsch       PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ILU, &flg2));
269566063dSJacob Faibussowitsch       PetscCall(MatIsSymmetricKnown(pc->pmat, &set, &flg3));
27da74c943SJose E. Roman       PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isnormal, MATNORMAL, MATNORMALHERMITIAN, NULL));
28ba8a7149SBarry Smith       if (flg1 && (!flg2 || (set && flg3))) {
29566e8bf2SBarry Smith         *type = PCICC;
30566e8bf2SBarry Smith       } else if (flg2) {
31566e8bf2SBarry Smith         *type = PCILU;
32da74c943SJose E. Roman       } else if (isnormal) {
33da74c943SJose E. Roman         *type = PCNONE;
34976e8c5aSLisandro Dalcin       } else if (hasop) { /* likely is a parallel matrix run on one processor */
35566e8bf2SBarry Smith         *type = PCBJACOBI;
36566e8bf2SBarry Smith       } else {
37566e8bf2SBarry Smith         *type = PCNONE;
38566e8bf2SBarry Smith       }
39566e8bf2SBarry Smith     } else {
40976e8c5aSLisandro Dalcin       if (hasop) {
41566e8bf2SBarry Smith         *type = PCBJACOBI;
42566e8bf2SBarry Smith       } else {
43566e8bf2SBarry Smith         *type = PCNONE;
44566e8bf2SBarry Smith       }
45566e8bf2SBarry Smith     }
46566e8bf2SBarry Smith   } else {
47566e8bf2SBarry Smith     if (size == 1) {
48566e8bf2SBarry Smith       *type = PCILU;
49566e8bf2SBarry Smith     } else {
50566e8bf2SBarry Smith       *type = PCBJACOBI;
51566e8bf2SBarry Smith     }
52566e8bf2SBarry Smith   }
533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54566e8bf2SBarry Smith }
554b9ad928SBarry Smith 
5688584be7SBarry Smith /*@
57*562efe2eSBarry Smith   PCReset - Resets a `PC` context to the pcsetupcalled = 0 state and removes any allocated `Vec`s and `Mat`s
5888584be7SBarry Smith 
59c3339decSBarry Smith   Collective
6088584be7SBarry Smith 
6188584be7SBarry Smith   Input Parameter:
6288584be7SBarry Smith . pc - the preconditioner context
6388584be7SBarry Smith 
6488584be7SBarry Smith   Level: developer
6588584be7SBarry Smith 
66f1580f4eSBarry Smith   Note:
67*562efe2eSBarry 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`
6888584be7SBarry Smith 
69*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`
7088584be7SBarry Smith @*/
71d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset(PC pc)
72d71ae5a4SJacob Faibussowitsch {
7388584be7SBarry Smith   PetscFunctionBegin;
7488584be7SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
75dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, reset);
769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleright));
779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleleft));
789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->pmat));
799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->mat));
802fa5cd67SKarl Rupp 
810ce6c5a2SBarry Smith   pc->setupcalled = 0;
823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8388584be7SBarry Smith }
8488584be7SBarry Smith 
851fb7b255SJunchao Zhang /*@C
86f1580f4eSBarry Smith   PCDestroy - Destroys `PC` context that was created with `PCCreate()`.
874b9ad928SBarry Smith 
88c3339decSBarry Smith   Collective
894b9ad928SBarry Smith 
904b9ad928SBarry Smith   Input Parameter:
914b9ad928SBarry Smith . pc - the preconditioner context
924b9ad928SBarry Smith 
934b9ad928SBarry Smith   Level: developer
944b9ad928SBarry Smith 
95*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`
964b9ad928SBarry Smith @*/
97d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy(PC *pc)
98d71ae5a4SJacob Faibussowitsch {
994b9ad928SBarry Smith   PetscFunctionBegin;
1003ba16761SJacob Faibussowitsch   if (!*pc) PetscFunctionReturn(PETSC_SUCCESS);
1016bf464f9SBarry Smith   PetscValidHeaderSpecific((*pc), PC_CLASSID, 1);
1029371c9d4SSatish Balay   if (--((PetscObject)(*pc))->refct > 0) {
1039371c9d4SSatish Balay     *pc = NULL;
1043ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1059371c9d4SSatish Balay   }
1064b9ad928SBarry Smith 
1079566063dSJacob Faibussowitsch   PetscCall(PCReset(*pc));
108241cb3abSLisandro Dalcin 
109e04113cfSBarry Smith   /* if memory was published with SAWs then destroy it */
1109566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*pc));
111dbbe0bcdSBarry Smith   PetscTryTypeMethod((*pc), destroy);
1129566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&(*pc)->dm));
1139566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(pc));
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1154b9ad928SBarry Smith }
1164b9ad928SBarry Smith 
1174b9ad928SBarry Smith /*@C
118c5eb9154SBarry Smith   PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right
1194b9ad928SBarry Smith   scaling as needed by certain time-stepping codes.
1204b9ad928SBarry Smith 
121c3339decSBarry Smith   Logically Collective
1224b9ad928SBarry Smith 
1234b9ad928SBarry Smith   Input Parameter:
1244b9ad928SBarry Smith . pc - the preconditioner context
1254b9ad928SBarry Smith 
1264b9ad928SBarry Smith   Output Parameter:
127f1580f4eSBarry Smith . flag - `PETSC_TRUE` if it applies the scaling
1284b9ad928SBarry Smith 
1294b9ad928SBarry Smith   Level: developer
1304b9ad928SBarry Smith 
131f1580f4eSBarry Smith   Note:
132*562efe2eSBarry Smith   If this returns `PETSC_TRUE` then the system solved via the Krylov method is, for left and right preconditioning,
1334b9ad928SBarry Smith 
134*562efe2eSBarry Smith   $$
135*562efe2eSBarry Smith   \begin{align*}
136*562efe2eSBarry Smith   D M A D^{-1} y = D M b  \\
137*562efe2eSBarry Smith   D A M D^{-1} z = D b.
138*562efe2eSBarry Smith   \end{align*}
139*562efe2eSBarry Smith   $$
140*562efe2eSBarry Smith 
141*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCSetDiagonalScale()`
1424b9ad928SBarry Smith @*/
143d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetDiagonalScale(PC pc, PetscBool *flag)
144d71ae5a4SJacob Faibussowitsch {
1454b9ad928SBarry Smith   PetscFunctionBegin;
1460700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1474f572ea9SToby Isaac   PetscAssertPointer(flag, 2);
1484b9ad928SBarry Smith   *flag = pc->diagonalscale;
1493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1504b9ad928SBarry Smith }
1514b9ad928SBarry Smith 
1524b9ad928SBarry Smith /*@
153c5eb9154SBarry Smith   PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right
1544b9ad928SBarry Smith   scaling as needed by certain time-stepping codes.
1554b9ad928SBarry Smith 
156c3339decSBarry Smith   Logically Collective
1574b9ad928SBarry Smith 
1584b9ad928SBarry Smith   Input Parameters:
1594b9ad928SBarry Smith + pc - the preconditioner context
1604b9ad928SBarry Smith - s  - scaling vector
1614b9ad928SBarry Smith 
1624b9ad928SBarry Smith   Level: intermediate
1634b9ad928SBarry Smith 
16495452b02SPatrick Sanan   Notes:
165*562efe2eSBarry Smith   The system solved via the Krylov method is, for left and right preconditioning,
166*562efe2eSBarry Smith   $$
167*562efe2eSBarry Smith   \begin{align*}
168*562efe2eSBarry Smith   D M A D^{-1} y = D M b \\
169*562efe2eSBarry Smith   D A M D^{-1} z = D b.
170*562efe2eSBarry Smith   \end{align*}
171*562efe2eSBarry Smith   $$
1724b9ad928SBarry Smith 
173*562efe2eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.
1744b9ad928SBarry Smith 
175*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCGetDiagonalScale()`
1764b9ad928SBarry Smith @*/
177d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetDiagonalScale(PC pc, Vec s)
178d71ae5a4SJacob Faibussowitsch {
1794b9ad928SBarry Smith   PetscFunctionBegin;
1800700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1810700a824SBarry Smith   PetscValidHeaderSpecific(s, VEC_CLASSID, 2);
1824b9ad928SBarry Smith   pc->diagonalscale = PETSC_TRUE;
1832fa5cd67SKarl Rupp 
1849566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)s));
1859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleleft));
1862fa5cd67SKarl Rupp 
1874b9ad928SBarry Smith   pc->diagonalscaleleft = s;
1882fa5cd67SKarl Rupp 
1899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(s, &pc->diagonalscaleright));
1909566063dSJacob Faibussowitsch   PetscCall(VecCopy(s, pc->diagonalscaleright));
1919566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(pc->diagonalscaleright));
1923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1934b9ad928SBarry Smith }
1944b9ad928SBarry Smith 
195bc08b0f1SBarry Smith /*@
196c5eb9154SBarry Smith   PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes.
1974b9ad928SBarry Smith 
198c3339decSBarry Smith   Logically Collective
1994b9ad928SBarry Smith 
2004b9ad928SBarry Smith   Input Parameters:
2014b9ad928SBarry Smith + pc  - the preconditioner context
2024b9ad928SBarry Smith . in  - input vector
203a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in)
2044b9ad928SBarry Smith 
2054b9ad928SBarry Smith   Level: intermediate
2064b9ad928SBarry Smith 
20795452b02SPatrick Sanan   Notes:
208*562efe2eSBarry Smith   The system solved via the Krylov method is, for left and right preconditioning,
2094b9ad928SBarry Smith 
210*562efe2eSBarry Smith   $$
211*562efe2eSBarry Smith   \begin{align*}
212*562efe2eSBarry Smith   D M A D^{-1} y = D M b  \\
213*562efe2eSBarry Smith   D A M D^{-1} z = D b.
214*562efe2eSBarry Smith   \end{align*}
215*562efe2eSBarry Smith   $$
2164b9ad928SBarry Smith 
217*562efe2eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.
2184b9ad928SBarry Smith 
219*562efe2eSBarry Smith   If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out`
220*562efe2eSBarry Smith 
221*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCSetDiagonalScale()`, `PCDiagonalScaleRight()`, `PCDiagonalScale()`
2224b9ad928SBarry Smith @*/
223d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleLeft(PC pc, Vec in, Vec out)
224d71ae5a4SJacob Faibussowitsch {
2254b9ad928SBarry Smith   PetscFunctionBegin;
2260700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2270700a824SBarry Smith   PetscValidHeaderSpecific(in, VEC_CLASSID, 2);
2280700a824SBarry Smith   PetscValidHeaderSpecific(out, VEC_CLASSID, 3);
2294b9ad928SBarry Smith   if (pc->diagonalscale) {
2309566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(out, pc->diagonalscaleleft, in));
2314b9ad928SBarry Smith   } else if (in != out) {
2329566063dSJacob Faibussowitsch     PetscCall(VecCopy(in, out));
2334b9ad928SBarry Smith   }
2343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2354b9ad928SBarry Smith }
2364b9ad928SBarry Smith 
237bc08b0f1SBarry Smith /*@
2384b9ad928SBarry Smith   PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
2394b9ad928SBarry Smith 
240c3339decSBarry Smith   Logically Collective
2414b9ad928SBarry Smith 
2424b9ad928SBarry Smith   Input Parameters:
2434b9ad928SBarry Smith + pc  - the preconditioner context
2444b9ad928SBarry Smith . in  - input vector
245a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in)
2464b9ad928SBarry Smith 
2474b9ad928SBarry Smith   Level: intermediate
2484b9ad928SBarry Smith 
24995452b02SPatrick Sanan   Notes:
250*562efe2eSBarry Smith   The system solved via the Krylov method is, for left and right preconditioning,
2514b9ad928SBarry Smith 
252*562efe2eSBarry Smith   $$
253*562efe2eSBarry Smith   \begin{align*}
254*562efe2eSBarry Smith   D M A D^{-1} y = D M b  \\
255*562efe2eSBarry Smith   D A M D^{-1} z = D b.
256*562efe2eSBarry Smith   \end{aligne*}
257*562efe2eSBarry Smith   $$
2584b9ad928SBarry Smith 
259*562efe2eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.
2604b9ad928SBarry Smith 
261*562efe2eSBarry Smith   If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out`
262*562efe2eSBarry Smith 
263*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCSetDiagonalScale()`, `PCDiagonalScale()`
2644b9ad928SBarry Smith @*/
265d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleRight(PC pc, Vec in, Vec out)
266d71ae5a4SJacob Faibussowitsch {
2674b9ad928SBarry Smith   PetscFunctionBegin;
2680700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2690700a824SBarry Smith   PetscValidHeaderSpecific(in, VEC_CLASSID, 2);
2700700a824SBarry Smith   PetscValidHeaderSpecific(out, VEC_CLASSID, 3);
2714b9ad928SBarry Smith   if (pc->diagonalscale) {
2729566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(out, pc->diagonalscaleright, in));
2734b9ad928SBarry Smith   } else if (in != out) {
2749566063dSJacob Faibussowitsch     PetscCall(VecCopy(in, out));
2754b9ad928SBarry Smith   }
2763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2774b9ad928SBarry Smith }
2784b9ad928SBarry Smith 
27949517cdeSBarry Smith /*@
28049517cdeSBarry Smith   PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the
281f1580f4eSBarry Smith   operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
282f1580f4eSBarry Smith   `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.
28349517cdeSBarry Smith 
284c3339decSBarry Smith   Logically Collective
28549517cdeSBarry Smith 
28649517cdeSBarry Smith   Input Parameters:
28749517cdeSBarry Smith + pc  - the preconditioner context
288f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)
28949517cdeSBarry Smith 
29049517cdeSBarry Smith   Options Database Key:
291*562efe2eSBarry Smith . -pc_use_amat <true,false> - use the amat argument to `KSPSetOperators()` or `PCSetOperators()` to apply the operator
29249517cdeSBarry Smith 
29320f4b53cSBarry Smith   Level: intermediate
29420f4b53cSBarry Smith 
295f1580f4eSBarry Smith   Note:
29649517cdeSBarry Smith   For the common case in which the linear system matrix and the matrix used to construct the
297*562efe2eSBarry Smith   preconditioner are identical, this routine has no affect.
29849517cdeSBarry Smith 
299*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetUseAmat()`, `PCBJACOBI`, `PGMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`,
300*562efe2eSBarry Smith           `KSPSetOperators()`, `PCSetOperators()`
30149517cdeSBarry Smith @*/
302d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUseAmat(PC pc, PetscBool flg)
303d71ae5a4SJacob Faibussowitsch {
30449517cdeSBarry Smith   PetscFunctionBegin;
30549517cdeSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
30649517cdeSBarry Smith   pc->useAmat = flg;
3073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30849517cdeSBarry Smith }
30949517cdeSBarry Smith 
310422a814eSBarry Smith /*@
311*562efe2eSBarry Smith   PCSetErrorIfFailure - Causes `PC` to generate an error if a floating point exception, for example a zero pivot, is detected.
312422a814eSBarry Smith 
313c3339decSBarry Smith   Logically Collective
314422a814eSBarry Smith 
315422a814eSBarry Smith   Input Parameters:
316*562efe2eSBarry Smith + pc  - iterative context obtained from `PCCreate()`
317f1580f4eSBarry Smith - flg - `PETSC_TRUE` indicates you want the error generated
318422a814eSBarry Smith 
319422a814eSBarry Smith   Level: advanced
320422a814eSBarry Smith 
321422a814eSBarry Smith   Notes:
322f1580f4eSBarry Smith   Normally PETSc continues if a linear solver fails due to a failed setup of a preconditioner, you can call `KSPGetConvergedReason()` after a `KSPSolve()`
323422a814eSBarry 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
324422a814eSBarry Smith   detected.
325422a814eSBarry Smith 
326*562efe2eSBarry Smith   This is propagated into `KSP`s used by this `PC`, which then propagate it into `PC`s used by those `KSP`s
327422a814eSBarry Smith 
328*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPSetErrorIfNotConverged()`, `PCGetInitialGuessNonzero()`, `PCSetInitialGuessKnoll()`, `PCGetInitialGuessKnoll()`
329422a814eSBarry Smith @*/
330d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetErrorIfFailure(PC pc, PetscBool flg)
331d71ae5a4SJacob Faibussowitsch {
332422a814eSBarry Smith   PetscFunctionBegin;
333422a814eSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
334422a814eSBarry Smith   PetscValidLogicalCollectiveBool(pc, flg, 2);
335422a814eSBarry Smith   pc->erroriffailure = flg;
3363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
337422a814eSBarry Smith }
338422a814eSBarry Smith 
33949517cdeSBarry Smith /*@
34049517cdeSBarry Smith   PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the
341f1580f4eSBarry Smith   operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
342f1580f4eSBarry Smith   `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.
34349517cdeSBarry Smith 
344c3339decSBarry Smith   Logically Collective
34549517cdeSBarry Smith 
34649517cdeSBarry Smith   Input Parameter:
34749517cdeSBarry Smith . pc - the preconditioner context
34849517cdeSBarry Smith 
34949517cdeSBarry Smith   Output Parameter:
350f1580f4eSBarry Smith . flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)
35149517cdeSBarry Smith 
35220f4b53cSBarry Smith   Level: intermediate
35320f4b53cSBarry Smith 
354f1580f4eSBarry Smith   Note:
35549517cdeSBarry Smith   For the common case in which the linear system matrix and the matrix used to construct the
35649517cdeSBarry Smith   preconditioner are identical, this routine is does nothing.
35749517cdeSBarry Smith 
358*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUseAmat()`, `PCBJACOBI`, `PGMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`
35949517cdeSBarry Smith @*/
360d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetUseAmat(PC pc, PetscBool *flg)
361d71ae5a4SJacob Faibussowitsch {
36249517cdeSBarry Smith   PetscFunctionBegin;
36349517cdeSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
36449517cdeSBarry Smith   *flg = pc->useAmat;
3653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36649517cdeSBarry Smith }
36749517cdeSBarry Smith 
368f39d8e23SSatish Balay /*@
3693821be0aSBarry Smith   PCSetKSPNestLevel - sets the amount of nesting the `KSP` that contains this `PC` has
3703821be0aSBarry Smith 
3713821be0aSBarry Smith   Collective
3723821be0aSBarry Smith 
3733821be0aSBarry Smith   Input Parameters:
3743821be0aSBarry Smith + pc    - the `PC`
3753821be0aSBarry Smith - level - the nest level
3763821be0aSBarry Smith 
3773821be0aSBarry Smith   Level: developer
3783821be0aSBarry Smith 
3793821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCGetKSPNestLevel()`, `KSPSetNestLevel()`
3803821be0aSBarry Smith @*/
3813821be0aSBarry Smith PetscErrorCode PCSetKSPNestLevel(PC pc, PetscInt level)
3823821be0aSBarry Smith {
3833821be0aSBarry Smith   PetscFunctionBegin;
3847a99bfcaSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
3857a99bfcaSBarry Smith   PetscValidLogicalCollectiveInt(pc, level, 2);
3863821be0aSBarry Smith   pc->kspnestlevel = level;
3873821be0aSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3883821be0aSBarry Smith }
3893821be0aSBarry Smith 
3903821be0aSBarry Smith /*@
3913821be0aSBarry Smith   PCGetKSPNestLevel - gets the amount of nesting the `KSP` that contains this `PC` has
3923821be0aSBarry Smith 
3937a99bfcaSBarry Smith   Not Collective
3943821be0aSBarry Smith 
3953821be0aSBarry Smith   Input Parameter:
3963821be0aSBarry Smith . pc - the `PC`
3973821be0aSBarry Smith 
3983821be0aSBarry Smith   Output Parameter:
3993821be0aSBarry Smith . level - the nest level
4003821be0aSBarry Smith 
4013821be0aSBarry Smith   Level: developer
4023821be0aSBarry Smith 
4033821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `KSPGetNestLevel()`
4043821be0aSBarry Smith @*/
4053821be0aSBarry Smith PetscErrorCode PCGetKSPNestLevel(PC pc, PetscInt *level)
4063821be0aSBarry Smith {
4073821be0aSBarry Smith   PetscFunctionBegin;
4087a99bfcaSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4097a99bfcaSBarry Smith   PetscAssertPointer(level, 2);
4103821be0aSBarry Smith   *level = pc->kspnestlevel;
4113821be0aSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
4123821be0aSBarry Smith }
4133821be0aSBarry Smith 
4143821be0aSBarry Smith /*@
415f1580f4eSBarry Smith   PCCreate - Creates a preconditioner context, `PC`
4164b9ad928SBarry Smith 
417d083f849SBarry Smith   Collective
4184b9ad928SBarry Smith 
4194b9ad928SBarry Smith   Input Parameter:
4204b9ad928SBarry Smith . comm - MPI communicator
4214b9ad928SBarry Smith 
4224b9ad928SBarry Smith   Output Parameter:
4237a99bfcaSBarry Smith . newpc - location to put the preconditioner context
4244b9ad928SBarry Smith 
42520f4b53cSBarry Smith   Level: developer
42620f4b53cSBarry Smith 
427f1580f4eSBarry Smith   Note:
428f1580f4eSBarry 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`
429f1580f4eSBarry Smith   in parallel. For dense matrices it is always `PCNONE`.
4304b9ad928SBarry Smith 
431*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCApply()`, `PCDestroy()`
4324b9ad928SBarry Smith @*/
433d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreate(MPI_Comm comm, PC *newpc)
434d71ae5a4SJacob Faibussowitsch {
4354b9ad928SBarry Smith   PC pc;
4364b9ad928SBarry Smith 
4374b9ad928SBarry Smith   PetscFunctionBegin;
4384f572ea9SToby Isaac   PetscAssertPointer(newpc, 2);
4390a545947SLisandro Dalcin   *newpc = NULL;
4409566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
4414b9ad928SBarry Smith 
4429566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(pc, PC_CLASSID, "PC", "Preconditioner", "PC", comm, PCDestroy, PCView));
443a7d21a52SLisandro Dalcin 
4440a545947SLisandro Dalcin   pc->mat                  = NULL;
4450a545947SLisandro Dalcin   pc->pmat                 = NULL;
4464b9ad928SBarry Smith   pc->setupcalled          = 0;
4470c24e6a1SHong Zhang   pc->setfromoptionscalled = 0;
4480a545947SLisandro Dalcin   pc->data                 = NULL;
4494b9ad928SBarry Smith   pc->diagonalscale        = PETSC_FALSE;
4500a545947SLisandro Dalcin   pc->diagonalscaleleft    = NULL;
4510a545947SLisandro Dalcin   pc->diagonalscaleright   = NULL;
4524b9ad928SBarry Smith 
4530a545947SLisandro Dalcin   pc->modifysubmatrices  = NULL;
4540a545947SLisandro Dalcin   pc->modifysubmatricesP = NULL;
4552fa5cd67SKarl Rupp 
456a7d21a52SLisandro Dalcin   *newpc = pc;
4573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4584b9ad928SBarry Smith }
4594b9ad928SBarry Smith 
4604b9ad928SBarry Smith /*@
4614b9ad928SBarry Smith   PCApply - Applies the preconditioner to a vector.
4624b9ad928SBarry Smith 
463c3339decSBarry Smith   Collective
4644b9ad928SBarry Smith 
4654b9ad928SBarry Smith   Input Parameters:
4664b9ad928SBarry Smith + pc - the preconditioner context
4674b9ad928SBarry Smith - x  - input vector
4684b9ad928SBarry Smith 
4694b9ad928SBarry Smith   Output Parameter:
4704b9ad928SBarry Smith . y - output vector
4714b9ad928SBarry Smith 
4724b9ad928SBarry Smith   Level: developer
4734b9ad928SBarry Smith 
474*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`, `PCApplyBAorAB()`
4754b9ad928SBarry Smith @*/
476d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApply(PC pc, Vec x, Vec y)
477d71ae5a4SJacob Faibussowitsch {
478106e20bfSBarry Smith   PetscInt m, n, mv, nv;
4794b9ad928SBarry Smith 
4804b9ad928SBarry Smith   PetscFunctionBegin;
4810700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4820700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
4830700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
4847827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
485e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
486540bdfbdSVaclav Hapla   /* use pmat to check vector sizes since for KSPLSQR the pmat may be of a different size than mat */
4879566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &m, &n));
4889566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x, &mv));
4899566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(y, &nv));
490540bdfbdSVaclav Hapla   /* check pmat * y = x is feasible */
49163a3b9bcSJacob 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);
49263a3b9bcSJacob 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);
4939566063dSJacob Faibussowitsch   PetscCall(VecSetErrorIfLocked(y, 3));
494106e20bfSBarry Smith 
4959566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
4969566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
4979566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
498dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, apply, x, y);
4999566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
500e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
5019566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
5023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5034b9ad928SBarry Smith }
5044b9ad928SBarry Smith 
5054b9ad928SBarry Smith /*@
506*562efe2eSBarry Smith   PCMatApply - Applies the preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApply()`, `Y` and `X` must be different matrices.
507c41ea70eSPierre Jolivet 
508c3339decSBarry Smith   Collective
509c677e75fSPierre Jolivet 
510c677e75fSPierre Jolivet   Input Parameters:
511c677e75fSPierre Jolivet + pc - the preconditioner context
512c677e75fSPierre Jolivet - X  - block of input vectors
513c677e75fSPierre Jolivet 
514c677e75fSPierre Jolivet   Output Parameter:
515c677e75fSPierre Jolivet . Y - block of output vectors
516c677e75fSPierre Jolivet 
517c677e75fSPierre Jolivet   Level: developer
518c677e75fSPierre Jolivet 
519*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `KSPMatSolve()`
520c677e75fSPierre Jolivet @*/
521d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMatApply(PC pc, Mat X, Mat Y)
522d71ae5a4SJacob Faibussowitsch {
523c677e75fSPierre Jolivet   Mat       A;
524c677e75fSPierre Jolivet   Vec       cy, cx;
525bd82155bSPierre Jolivet   PetscInt  m1, M1, m2, M2, n1, N1, n2, N2, m3, M3, n3, N3;
526c677e75fSPierre Jolivet   PetscBool match;
527c677e75fSPierre Jolivet 
528c677e75fSPierre Jolivet   PetscFunctionBegin;
529c677e75fSPierre Jolivet   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
530e2b200f6SPierre Jolivet   PetscValidHeaderSpecific(X, MAT_CLASSID, 2);
531e2b200f6SPierre Jolivet   PetscValidHeaderSpecific(Y, MAT_CLASSID, 3);
532e2b200f6SPierre Jolivet   PetscCheckSameComm(pc, 1, X, 2);
533e2b200f6SPierre Jolivet   PetscCheckSameComm(pc, 1, Y, 3);
5347827d75bSBarry Smith   PetscCheck(Y != X, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "Y and X must be different matrices");
5359566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &A));
5369566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m3, &n3));
5379566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m2, &n2));
5389566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &m1, &n1));
5399566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M3, &N3));
5409566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &M2, &N2));
5419566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Y, &M1, &N1));
54263a3b9bcSJacob 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);
54363a3b9bcSJacob 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);
54463a3b9bcSJacob 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);
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &match, MATSEQDENSE, MATMPIDENSE, ""));
54628b400f6SJacob Faibussowitsch   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of output vectors not stored in a dense Mat");
5479566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
54828b400f6SJacob Faibussowitsch   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of input vectors not stored in a dense Mat");
5499566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
550c677e75fSPierre Jolivet   if (pc->ops->matapply) {
5519566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0));
552dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, matapply, X, Y);
5539566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0));
554c677e75fSPierre Jolivet   } else {
5559566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "PC type %s applying column by column\n", ((PetscObject)pc)->type_name));
556bd82155bSPierre Jolivet     for (n1 = 0; n1 < N1; ++n1) {
5579566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecRead(X, n1, &cx));
5589566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecWrite(Y, n1, &cy));
5599566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, cx, cy));
5609566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecWrite(Y, n1, &cy));
5619566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecRead(X, n1, &cx));
562c677e75fSPierre Jolivet     }
563c677e75fSPierre Jolivet   }
5643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
565c677e75fSPierre Jolivet }
566c677e75fSPierre Jolivet 
567c677e75fSPierre Jolivet /*@
5684b9ad928SBarry Smith   PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
5694b9ad928SBarry Smith 
570c3339decSBarry Smith   Collective
5714b9ad928SBarry Smith 
5724b9ad928SBarry Smith   Input Parameters:
5734b9ad928SBarry Smith + pc - the preconditioner context
5744b9ad928SBarry Smith - x  - input vector
5754b9ad928SBarry Smith 
5764b9ad928SBarry Smith   Output Parameter:
5774b9ad928SBarry Smith . y - output vector
5784b9ad928SBarry Smith 
57920f4b53cSBarry Smith   Level: developer
58020f4b53cSBarry Smith 
581f1580f4eSBarry Smith   Note:
582f1580f4eSBarry Smith   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
5834b9ad928SBarry Smith 
584*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricRight()`
5854b9ad928SBarry Smith @*/
586d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricLeft(PC pc, Vec x, Vec y)
587d71ae5a4SJacob Faibussowitsch {
5884b9ad928SBarry Smith   PetscFunctionBegin;
5890700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
5900700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
5910700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
5927827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
593e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
5949566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
5959566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
5969566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplySymmetricLeft, pc, x, y, 0));
597dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applysymmetricleft, x, y);
5989566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplySymmetricLeft, pc, x, y, 0));
5999566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
600e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
6013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6024b9ad928SBarry Smith }
6034b9ad928SBarry Smith 
6044b9ad928SBarry Smith /*@
6054b9ad928SBarry Smith   PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
6064b9ad928SBarry Smith 
607c3339decSBarry Smith   Collective
6084b9ad928SBarry Smith 
6094b9ad928SBarry Smith   Input Parameters:
6104b9ad928SBarry Smith + pc - the preconditioner context
6114b9ad928SBarry Smith - x  - input vector
6124b9ad928SBarry Smith 
6134b9ad928SBarry Smith   Output Parameter:
6144b9ad928SBarry Smith . y - output vector
6154b9ad928SBarry Smith 
6164b9ad928SBarry Smith   Level: developer
6174b9ad928SBarry Smith 
618f1580f4eSBarry Smith   Note:
619f1580f4eSBarry Smith   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
6204b9ad928SBarry Smith 
621*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricLeft()`
6224b9ad928SBarry Smith @*/
623d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricRight(PC pc, Vec x, Vec y)
624d71ae5a4SJacob Faibussowitsch {
6254b9ad928SBarry Smith   PetscFunctionBegin;
6260700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6270700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6280700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6297827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
630e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
6319566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
6329566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
6339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplySymmetricRight, pc, x, y, 0));
634dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applysymmetricright, x, y);
6359566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplySymmetricRight, pc, x, y, 0));
6369566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
637e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
6383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6394b9ad928SBarry Smith }
6404b9ad928SBarry Smith 
6414b9ad928SBarry Smith /*@
6424b9ad928SBarry Smith   PCApplyTranspose - Applies the transpose of preconditioner to a vector.
6434b9ad928SBarry Smith 
644c3339decSBarry Smith   Collective
6454b9ad928SBarry Smith 
6464b9ad928SBarry Smith   Input Parameters:
6474b9ad928SBarry Smith + pc - the preconditioner context
6484b9ad928SBarry Smith - x  - input vector
6494b9ad928SBarry Smith 
6504b9ad928SBarry Smith   Output Parameter:
6514b9ad928SBarry Smith . y - output vector
6524b9ad928SBarry Smith 
65320f4b53cSBarry Smith   Level: developer
65420f4b53cSBarry Smith 
655f1580f4eSBarry Smith   Note:
65695452b02SPatrick Sanan   For complex numbers this applies the non-Hermitian transpose.
6574c97465dSBarry Smith 
658*562efe2eSBarry Smith   Developer Note:
659f1580f4eSBarry Smith   We need to implement a `PCApplyHermitianTranspose()`
6604c97465dSBarry Smith 
661*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyBAorAB()`, `PCApplyBAorABTranspose()`, `PCApplyTransposeExists()`
6624b9ad928SBarry Smith @*/
663d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTranspose(PC pc, Vec x, Vec y)
664d71ae5a4SJacob Faibussowitsch {
6654b9ad928SBarry Smith   PetscFunctionBegin;
6660700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6670700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6680700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6697827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
670e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
6719566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
6729566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
6739566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
674dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applytranspose, x, y);
6759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
6769566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
677e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
6783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6794b9ad928SBarry Smith }
6804b9ad928SBarry Smith 
6814b9ad928SBarry Smith /*@
6829cc050e5SLisandro Dalcin   PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
6834b9ad928SBarry Smith 
684c3339decSBarry Smith   Collective
6854b9ad928SBarry Smith 
686f1580f4eSBarry Smith   Input Parameter:
6874b9ad928SBarry Smith . pc - the preconditioner context
6884b9ad928SBarry Smith 
6894b9ad928SBarry Smith   Output Parameter:
690f1580f4eSBarry Smith . flg - `PETSC_TRUE` if a transpose operation is defined
6914b9ad928SBarry Smith 
6924b9ad928SBarry Smith   Level: developer
6934b9ad928SBarry Smith 
694*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`
6954b9ad928SBarry Smith @*/
696d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTransposeExists(PC pc, PetscBool *flg)
697d71ae5a4SJacob Faibussowitsch {
6984b9ad928SBarry Smith   PetscFunctionBegin;
6990700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
7004f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
7016ce5e81cSLisandro Dalcin   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
7026ce5e81cSLisandro Dalcin   else *flg = PETSC_FALSE;
7033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7044b9ad928SBarry Smith }
7054b9ad928SBarry Smith 
7064b9ad928SBarry Smith /*@
707*562efe2eSBarry Smith   PCApplyBAorAB - Applies the preconditioner and operator to a vector. $y = B*A*x $ or $ y = A*B*x$.
7084b9ad928SBarry Smith 
709c3339decSBarry Smith   Collective
7104b9ad928SBarry Smith 
7114b9ad928SBarry Smith   Input Parameters:
7124b9ad928SBarry Smith + pc   - the preconditioner context
713f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
7144b9ad928SBarry Smith . x    - input vector
7154b9ad928SBarry Smith - work - work vector
7164b9ad928SBarry Smith 
7174b9ad928SBarry Smith   Output Parameter:
7184b9ad928SBarry Smith . y - output vector
7194b9ad928SBarry Smith 
7204b9ad928SBarry Smith   Level: developer
7214b9ad928SBarry Smith 
722f1580f4eSBarry Smith   Note:
723*562efe2eSBarry 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.
724*562efe2eSBarry Smith   The specific `KSPSolve()` method must also be written to handle the post-solve "correction" for the diagonal scaling.
72513e0d083SBarry Smith 
726*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorABTranspose()`
7274b9ad928SBarry Smith @*/
728d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorAB(PC pc, PCSide side, Vec x, Vec y, Vec work)
729d71ae5a4SJacob Faibussowitsch {
7304b9ad928SBarry Smith   PetscFunctionBegin;
7310700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
732a69c7061SStefano Zampini   PetscValidLogicalCollectiveEnum(pc, side, 2);
7330700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
7340700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
7350700a824SBarry Smith   PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
736a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, x, 3);
737a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, y, 4);
738a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, work, 5);
7397827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
7407827d75bSBarry 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");
7417827d75bSBarry Smith   PetscCheck(!pc->diagonalscale || side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot include diagonal scaling with symmetric preconditioner application");
742e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
7434b9ad928SBarry Smith 
7449566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
7454b9ad928SBarry Smith   if (pc->diagonalscale) {
7464b9ad928SBarry Smith     if (pc->ops->applyBA) {
7474b9ad928SBarry Smith       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
7489566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &work2));
7499566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, work2));
750dbbe0bcdSBarry Smith       PetscUseTypeMethod(pc, applyBA, side, work2, y, work);
7519566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7529566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&work2));
7534b9ad928SBarry Smith     } else if (side == PC_RIGHT) {
7549566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, y));
7559566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, y, work));
7569566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7579566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7584b9ad928SBarry Smith     } else if (side == PC_LEFT) {
7599566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, y));
7609566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, y, work));
7619566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, work, y));
7629566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7637827d75bSBarry Smith     } else PetscCheck(side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot provide diagonal scaling with symmetric application of preconditioner");
7644b9ad928SBarry Smith   } else {
7654b9ad928SBarry Smith     if (pc->ops->applyBA) {
766dbbe0bcdSBarry Smith       PetscUseTypeMethod(pc, applyBA, side, x, y, work);
7674b9ad928SBarry Smith     } else if (side == PC_RIGHT) {
7689566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, x, work));
7699566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7704b9ad928SBarry Smith     } else if (side == PC_LEFT) {
7719566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, x, work));
7729566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, work, y));
7734b9ad928SBarry Smith     } else if (side == PC_SYMMETRIC) {
7744b9ad928SBarry Smith       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
7759566063dSJacob Faibussowitsch       PetscCall(PCApplySymmetricRight(pc, x, work));
7769566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7779566063dSJacob Faibussowitsch       PetscCall(VecCopy(y, work));
7789566063dSJacob Faibussowitsch       PetscCall(PCApplySymmetricLeft(pc, work, y));
7794b9ad928SBarry Smith     }
7804b9ad928SBarry Smith   }
781e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
7823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7834b9ad928SBarry Smith }
7844b9ad928SBarry Smith 
7854b9ad928SBarry Smith /*@
7864b9ad928SBarry Smith   PCApplyBAorABTranspose - Applies the transpose of the preconditioner
787*562efe2eSBarry Smith   and operator to a vector. That is, applies $B^T * A^T$ with left preconditioning,
788*562efe2eSBarry Smith   NOT $(B*A)^T = A^T*B^T$.
7894b9ad928SBarry Smith 
790c3339decSBarry Smith   Collective
7914b9ad928SBarry Smith 
7924b9ad928SBarry Smith   Input Parameters:
7934b9ad928SBarry Smith + pc   - the preconditioner context
794f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
7954b9ad928SBarry Smith . x    - input vector
7964b9ad928SBarry Smith - work - work vector
7974b9ad928SBarry Smith 
7984b9ad928SBarry Smith   Output Parameter:
7994b9ad928SBarry Smith . y - output vector
8004b9ad928SBarry Smith 
80120f4b53cSBarry Smith   Level: developer
80220f4b53cSBarry Smith 
803f1580f4eSBarry Smith   Note:
804*562efe2eSBarry 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
805*562efe2eSBarry Smith   defined by $B^T$. This is why this has the funny form that it computes $B^T * A^T$
8069b3038f0SBarry Smith 
807*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorAB()`
8084b9ad928SBarry Smith @*/
809d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorABTranspose(PC pc, PCSide side, Vec x, Vec y, Vec work)
810d71ae5a4SJacob Faibussowitsch {
8114b9ad928SBarry Smith   PetscFunctionBegin;
8120700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8130700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
8140700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
8150700a824SBarry Smith   PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
8167827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
817e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
8184b9ad928SBarry Smith   if (pc->ops->applyBAtranspose) {
819dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, applyBAtranspose, side, x, y, work);
820e0f629ddSJacob Faibussowitsch     if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8213ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
8224b9ad928SBarry Smith   }
8237827d75bSBarry Smith   PetscCheck(side == PC_LEFT || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right or left");
8244b9ad928SBarry Smith 
8259566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
8264b9ad928SBarry Smith   if (side == PC_RIGHT) {
8279566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(pc, x, work));
8289566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->mat, work, y));
8294b9ad928SBarry Smith   } else if (side == PC_LEFT) {
8309566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->mat, x, work));
8319566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(pc, work, y));
8324b9ad928SBarry Smith   }
8334b9ad928SBarry Smith   /* add support for PC_SYMMETRIC */
834e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8364b9ad928SBarry Smith }
8374b9ad928SBarry Smith 
8384b9ad928SBarry Smith /*@
8394b9ad928SBarry Smith   PCApplyRichardsonExists - Determines whether a particular preconditioner has a
8404b9ad928SBarry Smith   built-in fast application of Richardson's method.
8414b9ad928SBarry Smith 
8424b9ad928SBarry Smith   Not Collective
8434b9ad928SBarry Smith 
8444b9ad928SBarry Smith   Input Parameter:
8454b9ad928SBarry Smith . pc - the preconditioner
8464b9ad928SBarry Smith 
8474b9ad928SBarry Smith   Output Parameter:
848f1580f4eSBarry Smith . exists - `PETSC_TRUE` or `PETSC_FALSE`
8494b9ad928SBarry Smith 
8504b9ad928SBarry Smith   Level: developer
8514b9ad928SBarry Smith 
852*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCRICHARDSON`, `PCApplyRichardson()`
8534b9ad928SBarry Smith @*/
854d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyRichardsonExists(PC pc, PetscBool *exists)
855d71ae5a4SJacob Faibussowitsch {
8564b9ad928SBarry Smith   PetscFunctionBegin;
8570700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8584f572ea9SToby Isaac   PetscAssertPointer(exists, 2);
8594b9ad928SBarry Smith   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
8604b9ad928SBarry Smith   else *exists = PETSC_FALSE;
8613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8624b9ad928SBarry Smith }
8634b9ad928SBarry Smith 
8644b9ad928SBarry Smith /*@
8654b9ad928SBarry Smith   PCApplyRichardson - Applies several steps of Richardson iteration with
8664b9ad928SBarry Smith   the particular preconditioner. This routine is usually used by the
8674b9ad928SBarry Smith   Krylov solvers and not the application code directly.
8684b9ad928SBarry Smith 
869c3339decSBarry Smith   Collective
8704b9ad928SBarry Smith 
8714b9ad928SBarry Smith   Input Parameters:
8724b9ad928SBarry Smith + pc        - the preconditioner context
8737319c654SBarry Smith . b         - the right hand side
8744b9ad928SBarry Smith . w         - one work vector
8754b9ad928SBarry Smith . rtol      - relative decrease in residual norm convergence criteria
87670441072SBarry Smith . abstol    - absolute residual norm convergence criteria
8774b9ad928SBarry Smith . dtol      - divergence residual norm increase criteria
8787319c654SBarry Smith . its       - the number of iterations to apply.
8797319c654SBarry Smith - guesszero - if the input x contains nonzero initial guess
8804b9ad928SBarry Smith 
881d8d19677SJose E. Roman   Output Parameters:
8824d0a8057SBarry Smith + outits - number of iterations actually used (for SOR this always equals its)
8834d0a8057SBarry Smith . reason - the reason the apply terminated
884f1580f4eSBarry Smith - y      - the solution (also contains initial guess if guesszero is `PETSC_FALSE`
8854b9ad928SBarry Smith 
88620f4b53cSBarry Smith   Level: developer
88720f4b53cSBarry Smith 
8884b9ad928SBarry Smith   Notes:
8894b9ad928SBarry Smith   Most preconditioners do not support this function. Use the command
890f1580f4eSBarry Smith   `PCApplyRichardsonExists()` to determine if one does.
8914b9ad928SBarry Smith 
892f1580f4eSBarry Smith   Except for the `PCMG` this routine ignores the convergence tolerances
8934b9ad928SBarry Smith   and always runs for the number of iterations
8944b9ad928SBarry Smith 
895*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCApplyRichardsonExists()`
8964b9ad928SBarry Smith @*/
897d71ae5a4SJacob 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)
898d71ae5a4SJacob Faibussowitsch {
8994b9ad928SBarry Smith   PetscFunctionBegin;
9000700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9010700a824SBarry Smith   PetscValidHeaderSpecific(b, VEC_CLASSID, 2);
9020700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
9030700a824SBarry Smith   PetscValidHeaderSpecific(w, VEC_CLASSID, 4);
9047827d75bSBarry Smith   PetscCheck(b != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "b and y must be different vectors");
9059566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
906dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applyrichardson, b, y, w, rtol, abstol, dtol, its, guesszero, outits, reason);
9073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9084b9ad928SBarry Smith }
9094b9ad928SBarry Smith 
910422a814eSBarry Smith /*@
911f1580f4eSBarry Smith   PCSetFailedReason - Sets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
9121b2b9847SBarry Smith 
913c3339decSBarry Smith   Logically Collective
9141b2b9847SBarry Smith 
915d8d19677SJose E. Roman   Input Parameters:
9161b2b9847SBarry Smith + pc     - the preconditioner context
9171b2b9847SBarry Smith - reason - the reason it failedx
9181b2b9847SBarry Smith 
9191b2b9847SBarry Smith   Level: advanced
9201b2b9847SBarry Smith 
921*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCFailedReason`
9221b2b9847SBarry Smith @*/
923d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFailedReason(PC pc, PCFailedReason reason)
924d71ae5a4SJacob Faibussowitsch {
9251b2b9847SBarry Smith   PetscFunctionBegin;
9266479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9271b2b9847SBarry Smith   pc->failedreason = reason;
9283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9291b2b9847SBarry Smith }
9301b2b9847SBarry Smith 
9311b2b9847SBarry Smith /*@
932f1580f4eSBarry Smith   PCGetFailedReason - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
933422a814eSBarry Smith 
934c3339decSBarry Smith   Logically Collective
935422a814eSBarry Smith 
936422a814eSBarry Smith   Input Parameter:
937422a814eSBarry Smith . pc - the preconditioner context
938422a814eSBarry Smith 
939422a814eSBarry Smith   Output Parameter:
9401b2b9847SBarry Smith . reason - the reason it failed
941422a814eSBarry Smith 
942422a814eSBarry Smith   Level: advanced
943422a814eSBarry Smith 
944f1580f4eSBarry Smith   Note:
945f1580f4eSBarry Smith   This is the maximum over reason over all ranks in the PC communicator. It is only valid after
9466479e835SStefano Zampini   a call `KSPCheckDot()` or  `KSPCheckNorm()` inside a `KSPSolve()` or `PCReduceFailedReason()`.
9476479e835SStefano Zampini   It is not valid immediately after a `PCSetUp()` or `PCApply()`, then use `PCGetFailedReasonRank()`
9481b2b9847SBarry Smith 
949*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, ``PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReasonRank()`, `PCSetFailedReason()`
950422a814eSBarry Smith @*/
951d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason)
952d71ae5a4SJacob Faibussowitsch {
953422a814eSBarry Smith   PetscFunctionBegin;
9546479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9559754764cSHong Zhang   if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled;
9560ae0b32dSHong Zhang   else *reason = pc->failedreason;
9573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
958422a814eSBarry Smith }
959422a814eSBarry Smith 
9601b2b9847SBarry Smith /*@
961f1580f4eSBarry Smith   PCGetFailedReasonRank - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail on this MPI rank
9621b2b9847SBarry Smith 
963f1580f4eSBarry Smith   Not Collective
9641b2b9847SBarry Smith 
9651b2b9847SBarry Smith   Input Parameter:
9661b2b9847SBarry Smith . pc - the preconditioner context
9671b2b9847SBarry Smith 
9681b2b9847SBarry Smith   Output Parameter:
9691b2b9847SBarry Smith . reason - the reason it failed
9701b2b9847SBarry Smith 
97120f4b53cSBarry Smith   Level: advanced
97220f4b53cSBarry Smith 
973f1580f4eSBarry Smith   Note:
974*562efe2eSBarry Smith   Different processes may have different reasons or no reason, see `PCGetFailedReason()`
9751b2b9847SBarry Smith 
976*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`, `PCReduceFailedReason()`
9771b2b9847SBarry Smith @*/
978d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReasonRank(PC pc, PCFailedReason *reason)
979d71ae5a4SJacob Faibussowitsch {
9801b2b9847SBarry Smith   PetscFunctionBegin;
9816479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9821b2b9847SBarry Smith   if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled;
9831b2b9847SBarry Smith   else *reason = pc->failedreason;
9843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9851b2b9847SBarry Smith }
986422a814eSBarry Smith 
9876479e835SStefano Zampini /*@
9886479e835SStefano Zampini   PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC`
9896479e835SStefano Zampini 
9906479e835SStefano Zampini   Collective
9916479e835SStefano Zampini 
9926479e835SStefano Zampini   Input Parameter:
9936479e835SStefano Zampini . pc - the preconditioner context
9946479e835SStefano Zampini 
9956479e835SStefano Zampini   Level: advanced
9966479e835SStefano Zampini 
9976479e835SStefano Zampini   Note:
998*562efe2eSBarry Smith   Different MPI processes may have different reasons or no reason, see `PCGetFailedReason()`. This routine
9996479e835SStefano Zampini   makes them have a common value (failure if any MPI process had a failure).
10006479e835SStefano Zampini 
1001*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`
10026479e835SStefano Zampini @*/
10036479e835SStefano Zampini PetscErrorCode PCReduceFailedReason(PC pc)
10046479e835SStefano Zampini {
10056479e835SStefano Zampini   PetscInt buf;
10066479e835SStefano Zampini 
10076479e835SStefano Zampini   PetscFunctionBegin;
10086479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
10096479e835SStefano Zampini   buf = (PetscInt)pc->failedreason;
10106479e835SStefano Zampini   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
10116479e835SStefano Zampini   pc->failedreason = (PCFailedReason)buf;
10126479e835SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
10136479e835SStefano Zampini }
10146479e835SStefano Zampini 
1015c00cb57fSBarry Smith /*  Next line needed to deactivate KSP_Solve logging */
1016c00cb57fSBarry Smith #include <petsc/private/kspimpl.h>
1017c00cb57fSBarry Smith 
10184b9ad928SBarry Smith /*
10194b9ad928SBarry Smith       a setupcall of 0 indicates never setup,
102023ee1639SBarry Smith                      1 indicates has been previously setup
1021422a814eSBarry Smith                     -1 indicates a PCSetUp() was attempted and failed
10224b9ad928SBarry Smith */
10234b9ad928SBarry Smith /*@
10244b9ad928SBarry Smith   PCSetUp - Prepares for the use of a preconditioner.
10254b9ad928SBarry Smith 
1026c3339decSBarry Smith   Collective
10274b9ad928SBarry Smith 
10284b9ad928SBarry Smith   Input Parameter:
10294b9ad928SBarry Smith . pc - the preconditioner context
10304b9ad928SBarry Smith 
10314b9ad928SBarry Smith   Level: developer
10324b9ad928SBarry Smith 
1033*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`
10344b9ad928SBarry Smith @*/
1035d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp(PC pc)
1036d71ae5a4SJacob Faibussowitsch {
1037566e8bf2SBarry Smith   const char      *def;
1038fc9b51d3SBarry Smith   PetscObjectState matstate, matnonzerostate;
10394b9ad928SBarry Smith 
10404b9ad928SBarry Smith   PetscFunctionBegin;
10410700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
104228b400f6SJacob Faibussowitsch   PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first");
10436ce5e81cSLisandro Dalcin 
104423ee1639SBarry Smith   if (pc->setupcalled && pc->reusepreconditioner) {
10459566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n"));
10463ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10474b9ad928SBarry Smith   }
10484b9ad928SBarry Smith 
10499566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate));
10509566063dSJacob Faibussowitsch   PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate));
105123ee1639SBarry Smith   if (!pc->setupcalled) {
10529566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Setting up PC for first time\n"));
105323ee1639SBarry Smith     pc->flag = DIFFERENT_NONZERO_PATTERN;
10549df67409SStefano Zampini   } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS);
10559df67409SStefano Zampini   else {
10569df67409SStefano Zampini     if (matnonzerostate != pc->matnonzerostate) {
10579566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n"));
105823ee1639SBarry Smith       pc->flag = DIFFERENT_NONZERO_PATTERN;
105923ee1639SBarry Smith     } else {
10609566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n"));
106123ee1639SBarry Smith       pc->flag = SAME_NONZERO_PATTERN;
106223ee1639SBarry Smith     }
106323ee1639SBarry Smith   }
106423ee1639SBarry Smith   pc->matstate        = matstate;
1065fc9b51d3SBarry Smith   pc->matnonzerostate = matnonzerostate;
106623ee1639SBarry Smith 
10677adad957SLisandro Dalcin   if (!((PetscObject)pc)->type_name) {
10689566063dSJacob Faibussowitsch     PetscCall(PCGetDefaultType_Private(pc, &def));
10699566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, def));
1070566e8bf2SBarry Smith   }
10714b9ad928SBarry Smith 
10729566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
10739566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure));
10749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0));
10754b9ad928SBarry Smith   if (pc->ops->setup) {
1076c00cb57fSBarry Smith     /* do not log solves and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */
10779566063dSJacob Faibussowitsch     PetscCall(KSPInitializePackage());
10789566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePush(KSP_Solve));
10799566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePush(PC_Apply));
1080dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, setup);
10819566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePop(KSP_Solve));
10829566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePop(PC_Apply));
10834b9ad928SBarry Smith   }
10849566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0));
1085422a814eSBarry Smith   if (!pc->setupcalled) pc->setupcalled = 1;
10863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10874b9ad928SBarry Smith }
10884b9ad928SBarry Smith 
10894b9ad928SBarry Smith /*@
10904b9ad928SBarry Smith   PCSetUpOnBlocks - Sets up the preconditioner for each block in
10914b9ad928SBarry Smith   the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
10924b9ad928SBarry Smith   methods.
10934b9ad928SBarry Smith 
1094c3339decSBarry Smith   Collective
10954b9ad928SBarry Smith 
1096f1580f4eSBarry Smith   Input Parameter:
10974b9ad928SBarry Smith . pc - the preconditioner context
10984b9ad928SBarry Smith 
10994b9ad928SBarry Smith   Level: developer
11004b9ad928SBarry Smith 
1101f1580f4eSBarry Smith   Note:
1102f1580f4eSBarry Smith   For nested preconditioners such as `PCBJACOBI` `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is
1103f1580f4eSBarry Smith   called on the outer `PC`, this routine ensures it is called.
1104f1580f4eSBarry Smith 
1105*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()`
11064b9ad928SBarry Smith @*/
1107d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUpOnBlocks(PC pc)
1108d71ae5a4SJacob Faibussowitsch {
11094b9ad928SBarry Smith   PetscFunctionBegin;
11100700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11113ba16761SJacob Faibussowitsch   if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS);
11129566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0));
1113dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, setuponblocks);
11149566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0));
11153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11164b9ad928SBarry Smith }
11174b9ad928SBarry Smith 
11184b9ad928SBarry Smith /*@C
11194b9ad928SBarry Smith   PCSetModifySubMatrices - Sets a user-defined routine for modifying the
112004c3f3b8SBarry Smith   submatrices that arise within certain subdomain-based preconditioners such as `PCASM`
11214b9ad928SBarry Smith 
1122c3339decSBarry Smith   Logically Collective
11234b9ad928SBarry Smith 
11244b9ad928SBarry Smith   Input Parameters:
11254b9ad928SBarry Smith + pc   - the preconditioner context
11264b9ad928SBarry Smith . func - routine for modifying the submatrices
112704c3f3b8SBarry Smith - ctx  - optional user-defined context (may be `NULL`)
11284b9ad928SBarry Smith 
112920f4b53cSBarry Smith   Calling sequence of `func`:
113020f4b53cSBarry Smith + pc     - the preconditioner context
113104c3f3b8SBarry Smith . nsub   - number of index sets
113220f4b53cSBarry Smith . row    - an array of index sets that contain the global row numbers
11334b9ad928SBarry Smith          that comprise each local submatrix
11344b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
11354b9ad928SBarry Smith          that comprise each local submatrix
11364b9ad928SBarry Smith . submat - array of local submatrices
11374b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
113804c3f3b8SBarry Smith          user-defined func routine (may be `NULL`)
11394b9ad928SBarry Smith 
114020f4b53cSBarry Smith   Level: advanced
114120f4b53cSBarry Smith 
11424b9ad928SBarry Smith   Notes:
114304c3f3b8SBarry Smith   The basic submatrices are extracted from the preconditioner matrix as
114404c3f3b8SBarry Smith   usual; the user can then alter these (for example, to set different boundary
114504c3f3b8SBarry Smith   conditions for each submatrix) before they are used for the local solves.
114604c3f3b8SBarry Smith 
1147f1580f4eSBarry Smith   `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and
1148f1580f4eSBarry Smith   `KSPSolve()`.
11494b9ad928SBarry Smith 
1150f1580f4eSBarry Smith   A routine set by `PCSetModifySubMatrices()` is currently called within
1151f1580f4eSBarry Smith   the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`)
11524b9ad928SBarry Smith   preconditioners.  All other preconditioners ignore this routine.
11534b9ad928SBarry Smith 
1154*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()`
11554b9ad928SBarry Smith @*/
115604c3f3b8SBarry Smith PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx), void *ctx)
1157d71ae5a4SJacob Faibussowitsch {
11584b9ad928SBarry Smith   PetscFunctionBegin;
11590700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11604b9ad928SBarry Smith   pc->modifysubmatrices  = func;
11614b9ad928SBarry Smith   pc->modifysubmatricesP = ctx;
11623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11634b9ad928SBarry Smith }
11644b9ad928SBarry Smith 
11654b9ad928SBarry Smith /*@C
11664b9ad928SBarry Smith   PCModifySubMatrices - Calls an optional user-defined routine within
1167f1580f4eSBarry Smith   certain preconditioners if one has been set with `PCSetModifySubMatrices()`.
11684b9ad928SBarry Smith 
1169c3339decSBarry Smith   Collective
11704b9ad928SBarry Smith 
11714b9ad928SBarry Smith   Input Parameters:
11724b9ad928SBarry Smith + pc     - the preconditioner context
11734b9ad928SBarry Smith . nsub   - the number of local submatrices
11744b9ad928SBarry Smith . row    - an array of index sets that contain the global row numbers
11754b9ad928SBarry Smith          that comprise each local submatrix
11764b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
11774b9ad928SBarry Smith          that comprise each local submatrix
11784b9ad928SBarry Smith . submat - array of local submatrices
11794b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
1180*562efe2eSBarry Smith          user-defined routine (may be `NULL`)
11814b9ad928SBarry Smith 
11824b9ad928SBarry Smith   Output Parameter:
11834b9ad928SBarry Smith . submat - array of local submatrices (the entries of which may
11844b9ad928SBarry Smith             have been modified)
11854b9ad928SBarry Smith 
118620f4b53cSBarry Smith   Level: developer
118720f4b53cSBarry Smith 
118804c3f3b8SBarry Smith   Note:
11894b9ad928SBarry Smith   The user should NOT generally call this routine, as it will
119004c3f3b8SBarry Smith   automatically be called within certain preconditioners.
11914b9ad928SBarry Smith 
1192*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetModifySubMatrices()`
11934b9ad928SBarry Smith @*/
1194d71ae5a4SJacob Faibussowitsch PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx)
1195d71ae5a4SJacob Faibussowitsch {
11964b9ad928SBarry Smith   PetscFunctionBegin;
11970700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11983ba16761SJacob Faibussowitsch   if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS);
11999566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0));
12009566063dSJacob Faibussowitsch   PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx));
12019566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0));
12023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12034b9ad928SBarry Smith }
12044b9ad928SBarry Smith 
12054b9ad928SBarry Smith /*@
12064b9ad928SBarry Smith   PCSetOperators - Sets the matrix associated with the linear system and
12074b9ad928SBarry Smith   a (possibly) different one associated with the preconditioner.
12084b9ad928SBarry Smith 
1209c3339decSBarry Smith   Logically Collective
12104b9ad928SBarry Smith 
12114b9ad928SBarry Smith   Input Parameters:
12124b9ad928SBarry Smith + pc   - the preconditioner context
1213e5d3d808SBarry Smith . Amat - the matrix that defines the linear system
1214d1e9a80fSBarry Smith - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
12154b9ad928SBarry Smith 
121620f4b53cSBarry Smith   Level: intermediate
1217189c0b8aSBarry Smith 
121820f4b53cSBarry Smith   Notes:
121920f4b53cSBarry Smith   Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used.
122020f4b53cSBarry Smith 
122120f4b53cSBarry Smith   If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then
1222f1580f4eSBarry Smith   first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()`
1223f1580f4eSBarry Smith   on it and then pass it back in in your call to `KSPSetOperators()`.
1224189c0b8aSBarry Smith 
12254b9ad928SBarry Smith   More Notes about Repeated Solution of Linear Systems:
122620f4b53cSBarry Smith   PETSc does NOT reset the matrix entries of either `Amat` or `Pmat`
12274b9ad928SBarry Smith   to zero after a linear solve; the user is completely responsible for
1228f1580f4eSBarry Smith   matrix assembly.  See the routine `MatZeroEntries()` if desiring to
12294b9ad928SBarry Smith   zero all elements of a matrix.
12304b9ad928SBarry Smith 
1231*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`
12324b9ad928SBarry Smith  @*/
1233d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat)
1234d71ae5a4SJacob Faibussowitsch {
12353b3f6333SBarry Smith   PetscInt m1, n1, m2, n2;
12364b9ad928SBarry Smith 
12374b9ad928SBarry Smith   PetscFunctionBegin;
12380700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12390700a824SBarry Smith   if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
12400700a824SBarry Smith   if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3);
12413fc8bf9cSBarry Smith   if (Amat) PetscCheckSameComm(pc, 1, Amat, 2);
12423fc8bf9cSBarry Smith   if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3);
124331641f1aSBarry Smith   if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) {
12449566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Amat, &m1, &n1));
12459566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->mat, &m2, &n2));
124663a3b9bcSJacob 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);
12479566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Pmat, &m1, &n1));
12489566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2));
124963a3b9bcSJacob 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);
12503b3f6333SBarry Smith   }
12514b9ad928SBarry Smith 
125223ee1639SBarry Smith   if (Pmat != pc->pmat) {
125323ee1639SBarry Smith     /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
125423ee1639SBarry Smith     pc->matnonzerostate = -1;
125523ee1639SBarry Smith     pc->matstate        = -1;
125623ee1639SBarry Smith   }
125723ee1639SBarry Smith 
1258906ed7ccSBarry Smith   /* reference first in case the matrices are the same */
12599566063dSJacob Faibussowitsch   if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat));
12609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->mat));
12619566063dSJacob Faibussowitsch   if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat));
12629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->pmat));
12634b9ad928SBarry Smith   pc->mat  = Amat;
12644b9ad928SBarry Smith   pc->pmat = Pmat;
12653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1266e56f5c9eSBarry Smith }
1267e56f5c9eSBarry Smith 
126823ee1639SBarry Smith /*@
126923ee1639SBarry Smith   PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed.
127023ee1639SBarry Smith 
1271c3339decSBarry Smith   Logically Collective
127223ee1639SBarry Smith 
127323ee1639SBarry Smith   Input Parameters:
127423ee1639SBarry Smith + pc   - the preconditioner context
1275f1580f4eSBarry Smith - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
127623ee1639SBarry Smith 
12772b26979fSBarry Smith   Level: intermediate
12782b26979fSBarry Smith 
1279f1580f4eSBarry Smith   Note:
1280f1580f4eSBarry Smith   Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option
1281f1580f4eSBarry Smith   prevents this.
1282f1580f4eSBarry Smith 
1283*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()`
128423ee1639SBarry Smith  @*/
1285d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag)
1286d71ae5a4SJacob Faibussowitsch {
128723ee1639SBarry Smith   PetscFunctionBegin;
128823ee1639SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1289f9177522SStefano Zampini   PetscValidLogicalCollectiveBool(pc, flag, 2);
129023ee1639SBarry Smith   pc->reusepreconditioner = flag;
12913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12924b9ad928SBarry Smith }
12934b9ad928SBarry Smith 
1294c60c7ad4SBarry Smith /*@
1295f1580f4eSBarry Smith   PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed.
1296c60c7ad4SBarry Smith 
1297bba28a21SBarry Smith   Not Collective
1298c60c7ad4SBarry Smith 
1299c60c7ad4SBarry Smith   Input Parameter:
1300c60c7ad4SBarry Smith . pc - the preconditioner context
1301c60c7ad4SBarry Smith 
1302c60c7ad4SBarry Smith   Output Parameter:
1303f1580f4eSBarry Smith . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
1304c60c7ad4SBarry Smith 
1305d0418729SSatish Balay   Level: intermediate
1306d0418729SSatish Balay 
1307*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()`
1308c60c7ad4SBarry Smith  @*/
1309d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag)
1310d71ae5a4SJacob Faibussowitsch {
1311c60c7ad4SBarry Smith   PetscFunctionBegin;
1312c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
13134f572ea9SToby Isaac   PetscAssertPointer(flag, 2);
1314c60c7ad4SBarry Smith   *flag = pc->reusepreconditioner;
13153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1316c60c7ad4SBarry Smith }
1317c60c7ad4SBarry Smith 
1318487a658cSBarry Smith /*@
13194b9ad928SBarry Smith   PCGetOperators - Gets the matrix associated with the linear system and
13204b9ad928SBarry Smith   possibly a different one associated with the preconditioner.
13214b9ad928SBarry Smith 
1322*562efe2eSBarry Smith   Not Collective, though parallel `Mat`s are returned if `pc` is parallel
13234b9ad928SBarry Smith 
13244b9ad928SBarry Smith   Input Parameter:
13254b9ad928SBarry Smith . pc - the preconditioner context
13264b9ad928SBarry Smith 
13274b9ad928SBarry Smith   Output Parameters:
1328e5d3d808SBarry Smith + Amat - the matrix defining the linear system
132923ee1639SBarry Smith - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
13304b9ad928SBarry Smith 
13314b9ad928SBarry Smith   Level: intermediate
13324b9ad928SBarry Smith 
1333f1580f4eSBarry Smith   Note:
133495452b02SPatrick Sanan   Does not increase the reference count of the matrices, so you should not destroy them
1335298cc208SBarry Smith 
1336f1580f4eSBarry Smith   Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators
1337f1580f4eSBarry Smith   are created in `PC` and returned to the user. In this case, if both operators
133873950996SBarry Smith   mat and pmat are requested, two DIFFERENT operators will be returned. If
133973950996SBarry Smith   only one is requested both operators in the PC will be the same (i.e. as
1340f1580f4eSBarry Smith   if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats).
134173950996SBarry Smith   The user must set the sizes of the returned matrices and their type etc just
1342f1580f4eSBarry Smith   as if the user created them with `MatCreate()`. For example,
134373950996SBarry Smith 
1344f1580f4eSBarry Smith .vb
1345f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1346f1580f4eSBarry Smith            set size, type, etc of Amat
134773950996SBarry Smith 
1348f1580f4eSBarry Smith          MatCreate(comm,&mat);
1349f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Amat);
1350f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)mat);
1351f1580f4eSBarry Smith            set size, type, etc of Amat
1352f1580f4eSBarry Smith .ve
135373950996SBarry Smith 
135473950996SBarry Smith   and
135573950996SBarry Smith 
1356f1580f4eSBarry Smith .vb
1357f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1358f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
135973950996SBarry Smith 
1360f1580f4eSBarry Smith          MatCreate(comm,&Amat);
1361f1580f4eSBarry Smith          MatCreate(comm,&Pmat);
1362f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1363f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Amat);
1364f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Pmat);
1365f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
1366f1580f4eSBarry Smith .ve
136773950996SBarry Smith 
1368f1580f4eSBarry Smith   The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy
1369b8d709abSRichard Tran Mills   of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely
1370f1580f4eSBarry Smith   managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look
1371f1580f4eSBarry Smith   at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to
1372f1580f4eSBarry Smith   the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP
1373f1580f4eSBarry Smith   you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you).
1374f1580f4eSBarry Smith   Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when
137573950996SBarry Smith   it can be created for you?
137673950996SBarry Smith 
1377*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()`
13784b9ad928SBarry Smith @*/
1379d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat)
1380d71ae5a4SJacob Faibussowitsch {
13814b9ad928SBarry Smith   PetscFunctionBegin;
13820700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1383e5d3d808SBarry Smith   if (Amat) {
138473950996SBarry Smith     if (!pc->mat) {
13859fca8c71SStefano Zampini       if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */
13869a4708feSJed Brown         pc->mat = pc->pmat;
13879566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->mat));
1388e5d3d808SBarry Smith       } else { /* both Amat and Pmat are empty */
13899566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat));
1390e5d3d808SBarry Smith         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
139173950996SBarry Smith           pc->pmat = pc->mat;
13929566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->pmat));
139373950996SBarry Smith         }
139473950996SBarry Smith       }
13959a4708feSJed Brown     }
1396e5d3d808SBarry Smith     *Amat = pc->mat;
139773950996SBarry Smith   }
1398e5d3d808SBarry Smith   if (Pmat) {
139973950996SBarry Smith     if (!pc->pmat) {
1400e5d3d808SBarry Smith       if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */
14019a4708feSJed Brown         pc->pmat = pc->mat;
14029566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->pmat));
14039a4708feSJed Brown       } else {
14049566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat));
1405e5d3d808SBarry Smith         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
140673950996SBarry Smith           pc->mat = pc->pmat;
14079566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->mat));
140873950996SBarry Smith         }
140973950996SBarry Smith       }
14109a4708feSJed Brown     }
1411e5d3d808SBarry Smith     *Pmat = pc->pmat;
141273950996SBarry Smith   }
14133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14144b9ad928SBarry Smith }
14154b9ad928SBarry Smith 
1416906ed7ccSBarry Smith /*@C
1417906ed7ccSBarry Smith   PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1418f1580f4eSBarry Smith   possibly a different one associated with the preconditioner have been set in the `PC`.
1419906ed7ccSBarry Smith 
142020f4b53cSBarry Smith   Not Collective, though the results on all processes should be the same
1421906ed7ccSBarry Smith 
1422906ed7ccSBarry Smith   Input Parameter:
1423906ed7ccSBarry Smith . pc - the preconditioner context
1424906ed7ccSBarry Smith 
1425906ed7ccSBarry Smith   Output Parameters:
1426906ed7ccSBarry Smith + mat  - the matrix associated with the linear system was set
1427906ed7ccSBarry Smith - pmat - matrix associated with the preconditioner was set, usually the same
1428906ed7ccSBarry Smith 
1429906ed7ccSBarry Smith   Level: intermediate
1430906ed7ccSBarry Smith 
1431*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`
1432906ed7ccSBarry Smith @*/
1433d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat)
1434d71ae5a4SJacob Faibussowitsch {
1435906ed7ccSBarry Smith   PetscFunctionBegin;
14360700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1437906ed7ccSBarry Smith   if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1438906ed7ccSBarry Smith   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
14393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1440906ed7ccSBarry Smith }
1441906ed7ccSBarry Smith 
1442f39d8e23SSatish Balay /*@
1443a4fd02acSBarry Smith   PCFactorGetMatrix - Gets the factored matrix from the
1444f1580f4eSBarry Smith   preconditioner context.  This routine is valid only for the `PCLU`,
1445f1580f4eSBarry Smith   `PCILU`, `PCCHOLESKY`, and `PCICC` methods.
14464b9ad928SBarry Smith 
1447*562efe2eSBarry Smith   Not Collective though `mat` is parallel if `pc` is parallel
14484b9ad928SBarry Smith 
1449f1580f4eSBarry Smith   Input Parameter:
14504b9ad928SBarry Smith . pc - the preconditioner context
14514b9ad928SBarry Smith 
1452feefa0e1SJacob Faibussowitsch   Output Parameters:
14534b9ad928SBarry Smith . mat - the factored matrix
14544b9ad928SBarry Smith 
14554b9ad928SBarry Smith   Level: advanced
14564b9ad928SBarry Smith 
1457f1580f4eSBarry Smith   Note:
1458*562efe2eSBarry Smith   Does not increase the reference count for `mat` so DO NOT destroy it
14599405f653SBarry Smith 
1460*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC`
14614b9ad928SBarry Smith @*/
1462d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat)
1463d71ae5a4SJacob Faibussowitsch {
14644b9ad928SBarry Smith   PetscFunctionBegin;
14650700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14664f572ea9SToby Isaac   PetscAssertPointer(mat, 2);
14677def7855SStefano Zampini   PetscCall(PCFactorSetUpMatSolverType(pc));
1468dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, getfactoredmatrix, mat);
14693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14704b9ad928SBarry Smith }
14714b9ad928SBarry Smith 
14724b9ad928SBarry Smith /*@C
14734b9ad928SBarry Smith   PCSetOptionsPrefix - Sets the prefix used for searching for all
1474f1580f4eSBarry Smith   `PC` options in the database.
14754b9ad928SBarry Smith 
1476c3339decSBarry Smith   Logically Collective
14774b9ad928SBarry Smith 
14784b9ad928SBarry Smith   Input Parameters:
14794b9ad928SBarry Smith + pc     - the preconditioner context
1480f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
14814b9ad928SBarry Smith 
1482f1580f4eSBarry Smith   Note:
14834b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
14844b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
14854b9ad928SBarry Smith   hyphen.
14864b9ad928SBarry Smith 
14874b9ad928SBarry Smith   Level: advanced
14884b9ad928SBarry Smith 
1489*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()`
14904b9ad928SBarry Smith @*/
1491d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[])
1492d71ae5a4SJacob Faibussowitsch {
14934b9ad928SBarry Smith   PetscFunctionBegin;
14940700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14959566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix));
14963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14974b9ad928SBarry Smith }
14984b9ad928SBarry Smith 
14994b9ad928SBarry Smith /*@C
15004b9ad928SBarry Smith   PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1501f1580f4eSBarry Smith   `PC` options in the database.
15024b9ad928SBarry Smith 
1503c3339decSBarry Smith   Logically Collective
15044b9ad928SBarry Smith 
15054b9ad928SBarry Smith   Input Parameters:
15064b9ad928SBarry Smith + pc     - the preconditioner context
1507f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
15084b9ad928SBarry Smith 
1509f1580f4eSBarry Smith   Note:
15104b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
15114b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
15124b9ad928SBarry Smith   hyphen.
15134b9ad928SBarry Smith 
15144b9ad928SBarry Smith   Level: advanced
15154b9ad928SBarry Smith 
1516*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()`
15174b9ad928SBarry Smith @*/
1518d71ae5a4SJacob Faibussowitsch PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[])
1519d71ae5a4SJacob Faibussowitsch {
15204b9ad928SBarry Smith   PetscFunctionBegin;
15210700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15229566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix));
15233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15244b9ad928SBarry Smith }
15254b9ad928SBarry Smith 
15264b9ad928SBarry Smith /*@C
15274b9ad928SBarry Smith   PCGetOptionsPrefix - Gets the prefix used for searching for all
15284b9ad928SBarry Smith   PC options in the database.
15294b9ad928SBarry Smith 
15304b9ad928SBarry Smith   Not Collective
15314b9ad928SBarry Smith 
1532f1580f4eSBarry Smith   Input Parameter:
15334b9ad928SBarry Smith . pc - the preconditioner context
15344b9ad928SBarry Smith 
1535f1580f4eSBarry Smith   Output Parameter:
15364b9ad928SBarry Smith . prefix - pointer to the prefix string used, is returned
15374b9ad928SBarry Smith 
15384b9ad928SBarry Smith   Level: advanced
15394b9ad928SBarry Smith 
1540*562efe2eSBarry Smith   Fortran Note:
1541*562efe2eSBarry Smith   The user should pass in a string `prefix` of
1542*562efe2eSBarry Smith   sufficient length to hold the prefix.
1543*562efe2eSBarry Smith 
1544*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()`
15454b9ad928SBarry Smith @*/
1546d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[])
1547d71ae5a4SJacob Faibussowitsch {
15484b9ad928SBarry Smith   PetscFunctionBegin;
15490700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15504f572ea9SToby Isaac   PetscAssertPointer(prefix, 2);
15519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix));
15523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15534b9ad928SBarry Smith }
15544b9ad928SBarry Smith 
15558066bbecSBarry Smith /*
15568066bbecSBarry Smith    Indicates the right hand side will be changed by KSPSolve(), this occurs for a few
15578066bbecSBarry Smith   preconditioners including BDDC and Eisentat that transform the equations before applying
15588066bbecSBarry Smith   the Krylov methods
15598066bbecSBarry Smith */
1560d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change)
1561d71ae5a4SJacob Faibussowitsch {
1562a06fd7f2SStefano Zampini   PetscFunctionBegin;
1563a06fd7f2SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15644f572ea9SToby Isaac   PetscAssertPointer(change, 2);
1565a06fd7f2SStefano Zampini   *change = PETSC_FALSE;
1566cac4c232SBarry Smith   PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change));
15673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1568a06fd7f2SStefano Zampini }
1569a06fd7f2SStefano Zampini 
15704b9ad928SBarry Smith /*@
15714b9ad928SBarry Smith   PCPreSolve - Optional pre-solve phase, intended for any
15724b9ad928SBarry Smith   preconditioner-specific actions that must be performed before
15734b9ad928SBarry Smith   the iterative solve itself.
15744b9ad928SBarry Smith 
1575c3339decSBarry Smith   Collective
15764b9ad928SBarry Smith 
15774b9ad928SBarry Smith   Input Parameters:
15784b9ad928SBarry Smith + pc  - the preconditioner context
15794b9ad928SBarry Smith - ksp - the Krylov subspace context
15804b9ad928SBarry Smith 
15814b9ad928SBarry Smith   Level: developer
15824b9ad928SBarry Smith 
1583feefa0e1SJacob Faibussowitsch   Example Usage:
15844b9ad928SBarry Smith .vb
15854b9ad928SBarry Smith     PCPreSolve(pc,ksp);
158623ce1328SBarry Smith     KSPSolve(ksp,b,x);
15874b9ad928SBarry Smith     PCPostSolve(pc,ksp);
15884b9ad928SBarry Smith .ve
15894b9ad928SBarry Smith 
15904b9ad928SBarry Smith   Notes:
1591f1580f4eSBarry Smith   The pre-solve phase is distinct from the `PCSetUp()` phase.
15924b9ad928SBarry Smith 
1593f1580f4eSBarry Smith   `KSPSolve()` calls this directly, so is rarely called by the user.
15944b9ad928SBarry Smith 
1595*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCPostSolve()`
15964b9ad928SBarry Smith @*/
1597d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPreSolve(PC pc, KSP ksp)
1598d71ae5a4SJacob Faibussowitsch {
15994b9ad928SBarry Smith   Vec x, rhs;
16004b9ad928SBarry Smith 
16014b9ad928SBarry Smith   PetscFunctionBegin;
16020700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16030700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1604d9a03883SBarry Smith   pc->presolvedone++;
16057827d75bSBarry Smith   PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice");
16069566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
16079566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
16084b9ad928SBarry Smith 
1609dbbe0bcdSBarry Smith   if (pc->ops->presolve) PetscUseTypeMethod(pc, presolve, ksp, rhs, x);
16101baa6e33SBarry Smith   else if (pc->presolve) PetscCall((pc->presolve)(pc, ksp));
16113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16124b9ad928SBarry Smith }
16134b9ad928SBarry Smith 
1614f560b561SHong Zhang /*@C
1615f1580f4eSBarry Smith   PCSetPreSolve - Sets function used by `PCPreSolve()` which is intended for any
1616f560b561SHong Zhang   preconditioner-specific actions that must be performed before
1617f560b561SHong Zhang   the iterative solve itself.
1618f560b561SHong Zhang 
1619c3339decSBarry Smith   Logically Collective
1620f560b561SHong Zhang 
1621f560b561SHong Zhang   Input Parameters:
1622f560b561SHong Zhang + pc       - the preconditioner object
1623f560b561SHong Zhang - presolve - the function to call before the solve
1624f560b561SHong Zhang 
162520f4b53cSBarry Smith   Calling sequence of `presolve`:
162620f4b53cSBarry Smith + pc  - the `PC` context
162720f4b53cSBarry Smith - ksp - the `KSP` context
1628f560b561SHong Zhang 
1629f560b561SHong Zhang   Level: developer
1630f560b561SHong Zhang 
1631*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCPreSolve()`
1632f560b561SHong Zhang @*/
163304c3f3b8SBarry Smith PetscErrorCode PCSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC pc, KSP ksp))
1634d71ae5a4SJacob Faibussowitsch {
1635f560b561SHong Zhang   PetscFunctionBegin;
1636f560b561SHong Zhang   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1637f560b561SHong Zhang   pc->presolve = presolve;
16383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1639f560b561SHong Zhang }
1640f560b561SHong Zhang 
16414b9ad928SBarry Smith /*@
16424b9ad928SBarry Smith   PCPostSolve - Optional post-solve phase, intended for any
16434b9ad928SBarry Smith   preconditioner-specific actions that must be performed after
16444b9ad928SBarry Smith   the iterative solve itself.
16454b9ad928SBarry Smith 
1646c3339decSBarry Smith   Collective
16474b9ad928SBarry Smith 
16484b9ad928SBarry Smith   Input Parameters:
16494b9ad928SBarry Smith + pc  - the preconditioner context
16504b9ad928SBarry Smith - ksp - the Krylov subspace context
16514b9ad928SBarry Smith 
1652feefa0e1SJacob Faibussowitsch   Example Usage:
16534b9ad928SBarry Smith .vb
16544b9ad928SBarry Smith     PCPreSolve(pc,ksp);
165523ce1328SBarry Smith     KSPSolve(ksp,b,x);
16564b9ad928SBarry Smith     PCPostSolve(pc,ksp);
16574b9ad928SBarry Smith .ve
16584b9ad928SBarry Smith 
1659*562efe2eSBarry Smith   Level: developer
1660*562efe2eSBarry Smith 
16614b9ad928SBarry Smith   Note:
1662f1580f4eSBarry Smith   `KSPSolve()` calls this routine directly, so it is rarely called by the user.
16634b9ad928SBarry Smith 
1664*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCSetPostSolve()`, `PCSetPresolve()`, `PCPreSolve()`, `KSPSolve()`
16654b9ad928SBarry Smith @*/
1666d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPostSolve(PC pc, KSP ksp)
1667d71ae5a4SJacob Faibussowitsch {
16684b9ad928SBarry Smith   Vec x, rhs;
16694b9ad928SBarry Smith 
16704b9ad928SBarry Smith   PetscFunctionBegin;
16710700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16720700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1673d9a03883SBarry Smith   pc->presolvedone--;
16749566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
16759566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
1676dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, postsolve, ksp, rhs, x);
16773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16784b9ad928SBarry Smith }
16794b9ad928SBarry Smith 
168055849f57SBarry Smith /*@C
1681f1580f4eSBarry Smith   PCLoad - Loads a `PC` that has been stored in binary  with `PCView()`.
168255849f57SBarry Smith 
1683c3339decSBarry Smith   Collective
168455849f57SBarry Smith 
168555849f57SBarry Smith   Input Parameters:
1686f1580f4eSBarry Smith + newdm  - the newly loaded `PC`, this needs to have been created with `PCCreate()` or
1687f1580f4eSBarry Smith            some related function before a call to `PCLoad()`.
1688f1580f4eSBarry Smith - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()`
168955849f57SBarry Smith 
169055849f57SBarry Smith   Level: intermediate
169155849f57SBarry Smith 
1692f1580f4eSBarry Smith   Note:
1693*562efe2eSBarry Smith   The type is determined by the data in the file, any `PCType` set into the `PC` before this call is ignored.
169455849f57SBarry Smith 
1695*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()`
169655849f57SBarry Smith @*/
1697d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLoad(PC newdm, PetscViewer viewer)
1698d71ae5a4SJacob Faibussowitsch {
169955849f57SBarry Smith   PetscBool isbinary;
1700060da220SMatthew G. Knepley   PetscInt  classid;
170155849f57SBarry Smith   char      type[256];
170255849f57SBarry Smith 
170355849f57SBarry Smith   PetscFunctionBegin;
170455849f57SBarry Smith   PetscValidHeaderSpecific(newdm, PC_CLASSID, 1);
170555849f57SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
17069566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
170728b400f6SJacob Faibussowitsch   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
170855849f57SBarry Smith 
17099566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
17107827d75bSBarry Smith   PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file");
17119566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
17129566063dSJacob Faibussowitsch   PetscCall(PCSetType(newdm, type));
1713dbbe0bcdSBarry Smith   PetscTryTypeMethod(newdm, load, viewer);
17143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
171555849f57SBarry Smith }
171655849f57SBarry Smith 
17179804daf3SBarry Smith #include <petscdraw.h>
1718e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1719e04113cfSBarry Smith   #include <petscviewersaws.h>
17200acecf5bSBarry Smith #endif
1721fe2efc57SMark 
1722fe2efc57SMark /*@C
1723*562efe2eSBarry Smith   PCViewFromOptions - View from the `PC` based on options in the options database
1724fe2efc57SMark 
1725c3339decSBarry Smith   Collective
1726fe2efc57SMark 
1727fe2efc57SMark   Input Parameters:
172820f4b53cSBarry Smith + A    - the `PC` context
1729f1580f4eSBarry Smith . obj  - Optional object that provides the options prefix
1730736c3998SJose E. Roman - name - command line option
1731fe2efc57SMark 
1732fe2efc57SMark   Level: intermediate
1733f1580f4eSBarry Smith 
1734*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()`
1735fe2efc57SMark @*/
1736d71ae5a4SJacob Faibussowitsch PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[])
1737d71ae5a4SJacob Faibussowitsch {
1738fe2efc57SMark   PetscFunctionBegin;
1739fe2efc57SMark   PetscValidHeaderSpecific(A, PC_CLASSID, 1);
17409566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
17413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1742fe2efc57SMark }
1743fe2efc57SMark 
17444b9ad928SBarry Smith /*@C
1745f1580f4eSBarry Smith   PCView - Prints information about the `PC`
17464b9ad928SBarry Smith 
1747c3339decSBarry Smith   Collective
17484b9ad928SBarry Smith 
17494b9ad928SBarry Smith   Input Parameters:
1750feefa0e1SJacob Faibussowitsch + pc     - the `PC` context
17514b9ad928SBarry Smith - viewer - optional visualization context
17524b9ad928SBarry Smith 
175320f4b53cSBarry Smith   Level: developer
175420f4b53cSBarry Smith 
1755f1580f4eSBarry Smith   Notes:
17564b9ad928SBarry Smith   The available visualization contexts include
1757f1580f4eSBarry Smith +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
1758f1580f4eSBarry Smith -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
17594b9ad928SBarry Smith   output where only the first processor opens
17604b9ad928SBarry Smith   the file. All other processors send their
17614b9ad928SBarry Smith   data to the first processor to print.
17624b9ad928SBarry Smith 
17634b9ad928SBarry Smith   The user can open an alternative visualization contexts with
1764f1580f4eSBarry Smith   `PetscViewerASCIIOpen()` (output to a specified file).
17654b9ad928SBarry Smith 
1766*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PetscViewer`, `KSPView()`, `PetscViewerASCIIOpen()`
17674b9ad928SBarry Smith @*/
1768d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView(PC pc, PetscViewer viewer)
1769d71ae5a4SJacob Faibussowitsch {
177019fd82e9SBarry Smith   PCType            cstr;
17716cd81132SPierre Jolivet   PetscViewerFormat format;
17726cd81132SPierre Jolivet   PetscBool         iascii, isstring, isbinary, isdraw, pop = PETSC_FALSE;
1773e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1774536b137fSBarry Smith   PetscBool issaws;
17750acecf5bSBarry Smith #endif
17764b9ad928SBarry Smith 
17774b9ad928SBarry Smith   PetscFunctionBegin;
17780700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
177948a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
17800700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1781c9780b6fSBarry Smith   PetscCheckSameComm(pc, 1, viewer, 2);
17824b9ad928SBarry Smith 
17839566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
17849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
17859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
17869566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1787e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
17889566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
17890acecf5bSBarry Smith #endif
1790219b1045SBarry Smith 
179132077d6dSBarry Smith   if (iascii) {
17929566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer));
179348a46eb9SPierre Jolivet     if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, "  PC has not been set up so information may be incomplete\n"));
17949566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1795dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
17969566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
1797834dbeb0SBarry Smith     if (pc->mat) {
17986cd81132SPierre Jolivet       PetscCall(PetscViewerGetFormat(viewer, &format));
17996cd81132SPierre Jolivet       if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
18009566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
18016cd81132SPierre Jolivet         pop = PETSC_TRUE;
18026cd81132SPierre Jolivet       }
18034b9ad928SBarry Smith       if (pc->pmat == pc->mat) {
18049566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix = precond matrix:\n"));
18059566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
18069566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
18079566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
18084b9ad928SBarry Smith       } else {
1809834dbeb0SBarry Smith         if (pc->pmat) {
18109566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix followed by preconditioner matrix:\n"));
18114b9ad928SBarry Smith         } else {
18129566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix:\n"));
18134b9ad928SBarry Smith         }
18149566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
18159566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
18169566063dSJacob Faibussowitsch         if (pc->pmat) PetscCall(MatView(pc->pmat, viewer));
18179566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
18184b9ad928SBarry Smith       }
18196cd81132SPierre Jolivet       if (pop) PetscCall(PetscViewerPopFormat(viewer));
18204b9ad928SBarry Smith     }
18214b9ad928SBarry Smith   } else if (isstring) {
18229566063dSJacob Faibussowitsch     PetscCall(PCGetType(pc, &cstr));
18239566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr));
1824dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18259566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18269566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18275b0b0462SBarry Smith   } else if (isbinary) {
182855849f57SBarry Smith     PetscInt    classid = PC_FILE_CLASSID;
182955849f57SBarry Smith     MPI_Comm    comm;
183055849f57SBarry Smith     PetscMPIInt rank;
183155849f57SBarry Smith     char        type[256];
183255849f57SBarry Smith 
18339566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
18349566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1835dd400576SPatrick Sanan     if (rank == 0) {
18369566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
18379566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256));
18389566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
183955849f57SBarry Smith     }
1840dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
1841219b1045SBarry Smith   } else if (isdraw) {
1842219b1045SBarry Smith     PetscDraw draw;
1843d9884438SBarry Smith     char      str[25];
184489fd9fafSBarry Smith     PetscReal x, y, bottom, h;
1845d9884438SBarry Smith     PetscInt  n;
1846219b1045SBarry Smith 
18479566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
18489566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
18491d840656SPeter Brune     if (pc->mat) {
18509566063dSJacob Faibussowitsch       PetscCall(MatGetSize(pc->mat, &n, NULL));
185163a3b9bcSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n));
18521d840656SPeter Brune     } else {
18539566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name));
18541d840656SPeter Brune     }
18559566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
185689fd9fafSBarry Smith     bottom = y - h;
18579566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
1858dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18599566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
1860e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1861536b137fSBarry Smith   } else if (issaws) {
1862d45a07a7SBarry Smith     PetscMPIInt rank;
1863d45a07a7SBarry Smith 
18649566063dSJacob Faibussowitsch     PetscCall(PetscObjectName((PetscObject)pc));
18659566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
186648a46eb9SPierre Jolivet     if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer));
18679566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18689566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18690acecf5bSBarry Smith #endif
18704b9ad928SBarry Smith   }
18713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18724b9ad928SBarry Smith }
18734b9ad928SBarry Smith 
18744b9ad928SBarry Smith /*@C
1875*562efe2eSBarry Smith   PCRegister -  Adds a method (`PCType`) to the preconditioner package.
18761c84c290SBarry Smith 
18771c84c290SBarry Smith   Not collective
18781c84c290SBarry Smith 
18791c84c290SBarry Smith   Input Parameters:
188020f4b53cSBarry Smith + sname    - name of a new user-defined solver
188120f4b53cSBarry Smith - function - routine to create method context
18821c84c290SBarry Smith 
1883feefa0e1SJacob Faibussowitsch   Example Usage:
18841c84c290SBarry Smith .vb
1885bdf89e91SBarry Smith    PCRegister("my_solver", MySolverCreate);
18861c84c290SBarry Smith .ve
18871c84c290SBarry Smith 
18881c84c290SBarry Smith   Then, your solver can be chosen with the procedural interface via
18891c84c290SBarry Smith $     PCSetType(pc, "my_solver")
18901c84c290SBarry Smith   or at runtime via the option
18911c84c290SBarry Smith $     -pc_type my_solver
18924b9ad928SBarry Smith 
18934b9ad928SBarry Smith   Level: advanced
18941c84c290SBarry Smith 
189520f4b53cSBarry Smith   Note:
189620f4b53cSBarry Smith   `PCRegister()` may be called multiple times to add several user-defined preconditioners.
189720f4b53cSBarry Smith 
1898*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCType`, `PCRegisterAll()`
18994b9ad928SBarry Smith @*/
1900d71ae5a4SJacob Faibussowitsch PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC))
1901d71ae5a4SJacob Faibussowitsch {
19024b9ad928SBarry Smith   PetscFunctionBegin;
19039566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
19049566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCList, sname, function));
19053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19064b9ad928SBarry Smith }
19074b9ad928SBarry Smith 
1908d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y)
1909d71ae5a4SJacob Faibussowitsch {
1910186a3e20SStefano Zampini   PC pc;
1911186a3e20SStefano Zampini 
1912186a3e20SStefano Zampini   PetscFunctionBegin;
19139566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &pc));
19149566063dSJacob Faibussowitsch   PetscCall(PCApply(pc, X, Y));
19153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1916186a3e20SStefano Zampini }
1917186a3e20SStefano Zampini 
1918feefa0e1SJacob Faibussowitsch /*@C
19190bacdadaSStefano Zampini   PCComputeOperator - Computes the explicit preconditioned operator.
19204b9ad928SBarry Smith 
1921c3339decSBarry Smith   Collective
19224b9ad928SBarry Smith 
1923d8d19677SJose E. Roman   Input Parameters:
1924186a3e20SStefano Zampini + pc      - the preconditioner object
1925*562efe2eSBarry Smith - mattype - the `MatType` to be used for the operator
19264b9ad928SBarry Smith 
19274b9ad928SBarry Smith   Output Parameter:
1928a5b23f4aSJose E. Roman . mat - the explicit preconditioned operator
19294b9ad928SBarry Smith 
193020f4b53cSBarry Smith   Level: advanced
193120f4b53cSBarry Smith 
1932f1580f4eSBarry Smith   Note:
1933186a3e20SStefano Zampini   This computation is done by applying the operators to columns of the identity matrix.
1934186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
1935*562efe2eSBarry Smith   Currently, this routine uses a dense matrix format when `mattype` == `NULL`
19364b9ad928SBarry Smith 
1937*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSPComputeOperator()`, `MatType`
19384b9ad928SBarry Smith @*/
1939d71ae5a4SJacob Faibussowitsch PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat)
1940d71ae5a4SJacob Faibussowitsch {
1941186a3e20SStefano Zampini   PetscInt N, M, m, n;
1942186a3e20SStefano Zampini   Mat      A, Apc;
19434b9ad928SBarry Smith 
19444b9ad928SBarry Smith   PetscFunctionBegin;
19450700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
19464f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
19479566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, &A, NULL));
19489566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
19499566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
19509566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc));
19519566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (void (*)(void))MatMult_PC));
19529566063dSJacob Faibussowitsch   PetscCall(MatComputeOperator(Apc, mattype, mat));
19539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Apc));
19543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19554b9ad928SBarry Smith }
19564b9ad928SBarry Smith 
19576c237d78SBarry Smith /*@
19586c237d78SBarry Smith   PCSetCoordinates - sets the coordinates of all the nodes on the local process
19596c237d78SBarry Smith 
1960c3339decSBarry Smith   Collective
19616c237d78SBarry Smith 
19626c237d78SBarry Smith   Input Parameters:
19636c237d78SBarry Smith + pc     - the solver context
19646c237d78SBarry Smith . dim    - the dimension of the coordinates 1, 2, or 3
196514893cbeSStefano Zampini . nloc   - the blocked size of the coordinates array
196614893cbeSStefano Zampini - coords - the coordinates array
19676c237d78SBarry Smith 
19686c237d78SBarry Smith   Level: intermediate
19696c237d78SBarry Smith 
1970f1580f4eSBarry Smith   Note:
197120f4b53cSBarry Smith   `coords` is an array of the dim coordinates for the nodes on
197220f4b53cSBarry Smith   the local processor, of size `dim`*`nloc`.
197314893cbeSStefano Zampini   If there are 108 equation on a processor
19746c237d78SBarry Smith   for a displacement finite element discretization of elasticity (so
197514893cbeSStefano Zampini   that there are nloc = 36 = 108/3 nodes) then the array must have 108
19766c237d78SBarry Smith   double precision values (ie, 3 * 36).  These x y z coordinates
19776c237d78SBarry Smith   should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
19786c237d78SBarry Smith   ... , N-1.z ].
19796c237d78SBarry Smith 
1980*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `MatSetNearNullSpace()`
19816c237d78SBarry Smith @*/
1982d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[])
1983d71ae5a4SJacob Faibussowitsch {
19846c237d78SBarry Smith   PetscFunctionBegin;
198532954da3SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
198632954da3SStefano Zampini   PetscValidLogicalCollectiveInt(pc, dim, 2);
198722794d57SStefano Zampini   PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords));
19883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19896c237d78SBarry Smith }
1990fd2dd295SFande Kong 
1991fd2dd295SFande Kong /*@
1992fd2dd295SFande Kong   PCGetInterpolations - Gets interpolation matrices for all levels (except level 0)
1993fd2dd295SFande Kong 
1994c3339decSBarry Smith   Logically Collective
1995fd2dd295SFande Kong 
1996d8d19677SJose E. Roman   Input Parameter:
1997d8d19677SJose E. Roman . pc - the precondition context
1998fd2dd295SFande Kong 
1999d8d19677SJose E. Roman   Output Parameters:
2000d8d19677SJose E. Roman + num_levels     - the number of levels
2001*562efe2eSBarry Smith - interpolations - the interpolation matrices (size of `num_levels`-1)
2002fd2dd295SFande Kong 
2003fd2dd295SFande Kong   Level: advanced
2004fd2dd295SFande Kong 
2005*562efe2eSBarry Smith   Developer Note:
2006f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
2007fd2dd295SFande Kong 
2008*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()`
2009fd2dd295SFande Kong @*/
2010d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[])
2011d71ae5a4SJacob Faibussowitsch {
2012fd2dd295SFande Kong   PetscFunctionBegin;
2013fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20144f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20154f572ea9SToby Isaac   PetscAssertPointer(interpolations, 3);
2016cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations));
20173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2018fd2dd295SFande Kong }
2019fd2dd295SFande Kong 
2020fd2dd295SFande Kong /*@
2021fd2dd295SFande Kong   PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level)
2022fd2dd295SFande Kong 
2023c3339decSBarry Smith   Logically Collective
2024fd2dd295SFande Kong 
2025d8d19677SJose E. Roman   Input Parameter:
2026d8d19677SJose E. Roman . pc - the precondition context
2027fd2dd295SFande Kong 
2028d8d19677SJose E. Roman   Output Parameters:
2029d8d19677SJose E. Roman + num_levels      - the number of levels
2030*562efe2eSBarry Smith - coarseOperators - the coarse operator matrices (size of `num_levels`-1)
2031fd2dd295SFande Kong 
2032fd2dd295SFande Kong   Level: advanced
2033fd2dd295SFande Kong 
2034*562efe2eSBarry Smith   Developer Note:
2035f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
2036fd2dd295SFande Kong 
2037*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()`
2038fd2dd295SFande Kong @*/
2039d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
2040d71ae5a4SJacob Faibussowitsch {
2041fd2dd295SFande Kong   PetscFunctionBegin;
2042fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20434f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20444f572ea9SToby Isaac   PetscAssertPointer(coarseOperators, 3);
2045cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators));
20463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2047fd2dd295SFande Kong }
2048