xref: /petsc/src/ksp/pc/interface/precon.c (revision b4f8a55a69b2805e0b0c023bc6ec977a7834d796)
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 
14*b4f8a55aSStefano Zampini PETSC_INTERN PetscErrorCode PCGetDefaultType_Private(PC pc, const char *type[])
15d71ae5a4SJacob Faibussowitsch {
162e70eadcSBarry Smith   PetscMPIInt size;
17*b4f8a55aSStefano Zampini   PetscBool   hasopblock, hasopsolve, flg1, flg2, set, flg3, isnormal;
18566e8bf2SBarry Smith 
19566e8bf2SBarry Smith   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
21566e8bf2SBarry Smith   if (pc->pmat) {
22*b4f8a55aSStefano Zampini     PetscCall(MatHasOperation(pc->pmat, MATOP_GET_DIAGONAL_BLOCK, &hasopblock));
23*b4f8a55aSStefano Zampini     PetscCall(MatHasOperation(pc->pmat, MATOP_SOLVE, &hasopsolve));
24566e8bf2SBarry Smith     if (size == 1) {
259566063dSJacob Faibussowitsch       PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ICC, &flg1));
269566063dSJacob Faibussowitsch       PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ILU, &flg2));
279566063dSJacob Faibussowitsch       PetscCall(MatIsSymmetricKnown(pc->pmat, &set, &flg3));
28da74c943SJose E. Roman       PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isnormal, MATNORMAL, MATNORMALHERMITIAN, NULL));
29ba8a7149SBarry Smith       if (flg1 && (!flg2 || (set && flg3))) {
30566e8bf2SBarry Smith         *type = PCICC;
31566e8bf2SBarry Smith       } else if (flg2) {
32566e8bf2SBarry Smith         *type = PCILU;
33da74c943SJose E. Roman       } else if (isnormal) {
34da74c943SJose E. Roman         *type = PCNONE;
35*b4f8a55aSStefano Zampini       } else if (hasopblock) { /* likely is a parallel matrix run on one processor */
36566e8bf2SBarry Smith         *type = PCBJACOBI;
37*b4f8a55aSStefano Zampini       } else if (hasopsolve) {
38*b4f8a55aSStefano Zampini         *type = PCMAT;
39566e8bf2SBarry Smith       } else {
40566e8bf2SBarry Smith         *type = PCNONE;
41566e8bf2SBarry Smith       }
42566e8bf2SBarry Smith     } else {
43*b4f8a55aSStefano Zampini       if (hasopblock) {
44566e8bf2SBarry Smith         *type = PCBJACOBI;
45*b4f8a55aSStefano Zampini       } else if (hasopsolve) {
46*b4f8a55aSStefano Zampini         *type = PCMAT;
47566e8bf2SBarry Smith       } else {
48566e8bf2SBarry Smith         *type = PCNONE;
49566e8bf2SBarry Smith       }
50566e8bf2SBarry Smith     }
51*b4f8a55aSStefano Zampini   } else *type = NULL;
523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53566e8bf2SBarry Smith }
544b9ad928SBarry Smith 
5588584be7SBarry Smith /*@
56f1580f4eSBarry Smith   PCReset - Resets a PC context to the pcsetupcalled = 0 state and removes any allocated `Vec`s and `Mat`s
5788584be7SBarry Smith 
58c3339decSBarry Smith   Collective
5988584be7SBarry Smith 
6088584be7SBarry Smith   Input Parameter:
6188584be7SBarry Smith . pc - the preconditioner context
6288584be7SBarry Smith 
6388584be7SBarry Smith   Level: developer
6488584be7SBarry Smith 
65f1580f4eSBarry Smith   Note:
66f1580f4eSBarry 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 the PC
6788584be7SBarry Smith 
68f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCSetUp()`
6988584be7SBarry Smith @*/
70d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset(PC pc)
71d71ae5a4SJacob Faibussowitsch {
7288584be7SBarry Smith   PetscFunctionBegin;
7388584be7SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
74dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, reset);
759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleright));
769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleleft));
779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->pmat));
789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->mat));
792fa5cd67SKarl Rupp 
800ce6c5a2SBarry Smith   pc->setupcalled = 0;
813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8288584be7SBarry Smith }
8388584be7SBarry Smith 
841fb7b255SJunchao Zhang /*@C
85f1580f4eSBarry Smith   PCDestroy - Destroys `PC` context that was created with `PCCreate()`.
864b9ad928SBarry Smith 
87c3339decSBarry Smith   Collective
884b9ad928SBarry Smith 
894b9ad928SBarry Smith   Input Parameter:
904b9ad928SBarry Smith . pc - the preconditioner context
914b9ad928SBarry Smith 
924b9ad928SBarry Smith   Level: developer
934b9ad928SBarry Smith 
94f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCSetUp()`
954b9ad928SBarry Smith @*/
96d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy(PC *pc)
97d71ae5a4SJacob Faibussowitsch {
984b9ad928SBarry Smith   PetscFunctionBegin;
993ba16761SJacob Faibussowitsch   if (!*pc) PetscFunctionReturn(PETSC_SUCCESS);
1006bf464f9SBarry Smith   PetscValidHeaderSpecific((*pc), PC_CLASSID, 1);
1019371c9d4SSatish Balay   if (--((PetscObject)(*pc))->refct > 0) {
1029371c9d4SSatish Balay     *pc = NULL;
1033ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1049371c9d4SSatish Balay   }
1054b9ad928SBarry Smith 
1069566063dSJacob Faibussowitsch   PetscCall(PCReset(*pc));
107241cb3abSLisandro Dalcin 
108e04113cfSBarry Smith   /* if memory was published with SAWs then destroy it */
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*pc));
110dbbe0bcdSBarry Smith   PetscTryTypeMethod((*pc), destroy);
1119566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&(*pc)->dm));
1129566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(pc));
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1144b9ad928SBarry Smith }
1154b9ad928SBarry Smith 
1164b9ad928SBarry Smith /*@C
117c5eb9154SBarry Smith   PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right
1184b9ad928SBarry Smith   scaling as needed by certain time-stepping codes.
1194b9ad928SBarry Smith 
120c3339decSBarry Smith   Logically Collective
1214b9ad928SBarry Smith 
1224b9ad928SBarry Smith   Input Parameter:
1234b9ad928SBarry Smith . pc - the preconditioner context
1244b9ad928SBarry Smith 
1254b9ad928SBarry Smith   Output Parameter:
126f1580f4eSBarry Smith . flag - `PETSC_TRUE` if it applies the scaling
1274b9ad928SBarry Smith 
1284b9ad928SBarry Smith   Level: developer
1294b9ad928SBarry Smith 
130f1580f4eSBarry Smith   Note:
131f1580f4eSBarry Smith   If this returns `PETSC_TRUE` then the system solved via the Krylov method is
132f1580f4eSBarry Smith .vb
133f1580f4eSBarry Smith       D M A D^{-1} y = D M b  for left preconditioning or
134f1580f4eSBarry Smith       D A M D^{-1} z = D b for right preconditioning
135f1580f4eSBarry Smith .ve
1364b9ad928SBarry Smith 
137f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCSetDiagonalScale()`
1384b9ad928SBarry Smith @*/
139d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetDiagonalScale(PC pc, PetscBool *flag)
140d71ae5a4SJacob Faibussowitsch {
1414b9ad928SBarry Smith   PetscFunctionBegin;
1420700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1434f572ea9SToby Isaac   PetscAssertPointer(flag, 2);
1444b9ad928SBarry Smith   *flag = pc->diagonalscale;
1453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1464b9ad928SBarry Smith }
1474b9ad928SBarry Smith 
1484b9ad928SBarry Smith /*@
149c5eb9154SBarry Smith   PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right
1504b9ad928SBarry Smith   scaling as needed by certain time-stepping codes.
1514b9ad928SBarry Smith 
152c3339decSBarry Smith   Logically Collective
1534b9ad928SBarry Smith 
1544b9ad928SBarry Smith   Input Parameters:
1554b9ad928SBarry Smith + pc - the preconditioner context
1564b9ad928SBarry Smith - s  - scaling vector
1574b9ad928SBarry Smith 
1584b9ad928SBarry Smith   Level: intermediate
1594b9ad928SBarry Smith 
16095452b02SPatrick Sanan   Notes:
16195452b02SPatrick Sanan   The system solved via the Krylov method is
162f1580f4eSBarry Smith .vb
163f1580f4eSBarry Smith            D M A D^{-1} y = D M b  for left preconditioning or
164f1580f4eSBarry Smith            D A M D^{-1} z = D b for right preconditioning
165f1580f4eSBarry Smith .ve
1664b9ad928SBarry Smith 
167f1580f4eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by D. `PCDiagonalScaleRight()` scales a vector by D^{-1}.
1684b9ad928SBarry Smith 
169db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCGetDiagonalScale()`
1704b9ad928SBarry Smith @*/
171d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetDiagonalScale(PC pc, Vec s)
172d71ae5a4SJacob Faibussowitsch {
1734b9ad928SBarry Smith   PetscFunctionBegin;
1740700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1750700a824SBarry Smith   PetscValidHeaderSpecific(s, VEC_CLASSID, 2);
1764b9ad928SBarry Smith   pc->diagonalscale = PETSC_TRUE;
1772fa5cd67SKarl Rupp 
1789566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)s));
1799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleleft));
1802fa5cd67SKarl Rupp 
1814b9ad928SBarry Smith   pc->diagonalscaleleft = s;
1822fa5cd67SKarl Rupp 
1839566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(s, &pc->diagonalscaleright));
1849566063dSJacob Faibussowitsch   PetscCall(VecCopy(s, pc->diagonalscaleright));
1859566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(pc->diagonalscaleright));
1863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1874b9ad928SBarry Smith }
1884b9ad928SBarry Smith 
189bc08b0f1SBarry Smith /*@
190c5eb9154SBarry Smith   PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes.
1914b9ad928SBarry Smith 
192c3339decSBarry Smith   Logically Collective
1934b9ad928SBarry Smith 
1944b9ad928SBarry Smith   Input Parameters:
1954b9ad928SBarry Smith + pc  - the preconditioner context
1964b9ad928SBarry Smith . in  - input vector
197a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in)
1984b9ad928SBarry Smith 
1994b9ad928SBarry Smith   Level: intermediate
2004b9ad928SBarry Smith 
20195452b02SPatrick Sanan   Notes:
20295452b02SPatrick Sanan   The system solved via the Krylov method is
203f1580f4eSBarry Smith .vb
204f1580f4eSBarry Smith         D M A D^{-1} y = D M b  for left preconditioning or
205f1580f4eSBarry Smith         D A M D^{-1} z = D b for right preconditioning
206f1580f4eSBarry Smith .ve
2074b9ad928SBarry Smith 
208f1580f4eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by D. `PCDiagonalScaleRight()` scales a vector by D^{-1}.
2094b9ad928SBarry Smith 
2104b9ad928SBarry Smith   If diagonal scaling is turned off and in is not out then in is copied to out
2114b9ad928SBarry Smith 
212db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleSet()`, `PCDiagonalScaleRight()`, `PCDiagonalScale()`
2134b9ad928SBarry Smith @*/
214d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleLeft(PC pc, Vec in, Vec out)
215d71ae5a4SJacob Faibussowitsch {
2164b9ad928SBarry Smith   PetscFunctionBegin;
2170700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2180700a824SBarry Smith   PetscValidHeaderSpecific(in, VEC_CLASSID, 2);
2190700a824SBarry Smith   PetscValidHeaderSpecific(out, VEC_CLASSID, 3);
2204b9ad928SBarry Smith   if (pc->diagonalscale) {
2219566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(out, pc->diagonalscaleleft, in));
2224b9ad928SBarry Smith   } else if (in != out) {
2239566063dSJacob Faibussowitsch     PetscCall(VecCopy(in, out));
2244b9ad928SBarry Smith   }
2253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2264b9ad928SBarry Smith }
2274b9ad928SBarry Smith 
228bc08b0f1SBarry Smith /*@
2294b9ad928SBarry Smith   PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
2304b9ad928SBarry Smith 
231c3339decSBarry Smith   Logically Collective
2324b9ad928SBarry Smith 
2334b9ad928SBarry Smith   Input Parameters:
2344b9ad928SBarry Smith + pc  - the preconditioner context
2354b9ad928SBarry Smith . in  - input vector
236a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in)
2374b9ad928SBarry Smith 
2384b9ad928SBarry Smith   Level: intermediate
2394b9ad928SBarry Smith 
24095452b02SPatrick Sanan   Notes:
24195452b02SPatrick Sanan   The system solved via the Krylov method is
242f1580f4eSBarry Smith .vb
243f1580f4eSBarry Smith         D M A D^{-1} y = D M b  for left preconditioning or
244f1580f4eSBarry Smith         D A M D^{-1} z = D b for right preconditioning
245f1580f4eSBarry Smith .ve
2464b9ad928SBarry Smith 
247f1580f4eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by D. `PCDiagonalScaleRight()` scales a vector by D^{-1}.
2484b9ad928SBarry Smith 
2494b9ad928SBarry Smith   If diagonal scaling is turned off and in is not out then in is copied to out
2504b9ad928SBarry Smith 
251db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleSet()`, `PCDiagonalScale()`
2524b9ad928SBarry Smith @*/
253d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleRight(PC pc, Vec in, Vec out)
254d71ae5a4SJacob Faibussowitsch {
2554b9ad928SBarry Smith   PetscFunctionBegin;
2560700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2570700a824SBarry Smith   PetscValidHeaderSpecific(in, VEC_CLASSID, 2);
2580700a824SBarry Smith   PetscValidHeaderSpecific(out, VEC_CLASSID, 3);
2594b9ad928SBarry Smith   if (pc->diagonalscale) {
2609566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(out, pc->diagonalscaleright, in));
2614b9ad928SBarry Smith   } else if (in != out) {
2629566063dSJacob Faibussowitsch     PetscCall(VecCopy(in, out));
2634b9ad928SBarry Smith   }
2643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2654b9ad928SBarry Smith }
2664b9ad928SBarry Smith 
26749517cdeSBarry Smith /*@
26849517cdeSBarry Smith   PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the
269f1580f4eSBarry Smith   operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
270f1580f4eSBarry Smith   `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.
27149517cdeSBarry Smith 
272c3339decSBarry Smith   Logically Collective
27349517cdeSBarry Smith 
27449517cdeSBarry Smith   Input Parameters:
27549517cdeSBarry Smith + pc  - the preconditioner context
276f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)
27749517cdeSBarry Smith 
27849517cdeSBarry Smith   Options Database Key:
279147403d9SBarry Smith . -pc_use_amat <true,false> - use the amat to apply the operator
28049517cdeSBarry Smith 
28120f4b53cSBarry Smith   Level: intermediate
28220f4b53cSBarry Smith 
283f1580f4eSBarry Smith   Note:
28449517cdeSBarry Smith   For the common case in which the linear system matrix and the matrix used to construct the
28549517cdeSBarry Smith   preconditioner are identical, this routine is does nothing.
28649517cdeSBarry Smith 
287f1580f4eSBarry Smith .seealso: `PC`, `PCGetUseAmat()`, `PCBJACOBI`, `PGMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`
28849517cdeSBarry Smith @*/
289d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUseAmat(PC pc, PetscBool flg)
290d71ae5a4SJacob Faibussowitsch {
29149517cdeSBarry Smith   PetscFunctionBegin;
29249517cdeSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
29349517cdeSBarry Smith   pc->useAmat = flg;
2943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29549517cdeSBarry Smith }
29649517cdeSBarry Smith 
297422a814eSBarry Smith /*@
298f1580f4eSBarry Smith   PCSetErrorIfFailure - Causes `PC` to generate an error if a FPE, for example a zero pivot, is detected.
299422a814eSBarry Smith 
300c3339decSBarry Smith   Logically Collective
301422a814eSBarry Smith 
302422a814eSBarry Smith   Input Parameters:
303422a814eSBarry Smith + pc  - iterative context obtained from PCCreate()
304f1580f4eSBarry Smith - flg - `PETSC_TRUE` indicates you want the error generated
305422a814eSBarry Smith 
306422a814eSBarry Smith   Level: advanced
307422a814eSBarry Smith 
308422a814eSBarry Smith   Notes:
309f1580f4eSBarry Smith   Normally PETSc continues if a linear solver fails due to a failed setup of a preconditioner, you can call `KSPGetConvergedReason()` after a `KSPSolve()`
310422a814eSBarry 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
311422a814eSBarry Smith   detected.
312422a814eSBarry Smith 
313422a814eSBarry Smith   This is propagated into KSPs used by this PC, which then propagate it into PCs used by those KSPs
314422a814eSBarry Smith 
315f1580f4eSBarry Smith .seealso: `PC`, `KSPSetErrorIfNotConverged()`, `PCGetInitialGuessNonzero()`, `PCSetInitialGuessKnoll()`, `PCGetInitialGuessKnoll()`
316422a814eSBarry Smith @*/
317d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetErrorIfFailure(PC pc, PetscBool flg)
318d71ae5a4SJacob Faibussowitsch {
319422a814eSBarry Smith   PetscFunctionBegin;
320422a814eSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
321422a814eSBarry Smith   PetscValidLogicalCollectiveBool(pc, flg, 2);
322422a814eSBarry Smith   pc->erroriffailure = flg;
3233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
324422a814eSBarry Smith }
325422a814eSBarry Smith 
32649517cdeSBarry Smith /*@
32749517cdeSBarry Smith   PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the
328f1580f4eSBarry Smith   operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
329f1580f4eSBarry Smith   `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.
33049517cdeSBarry Smith 
331c3339decSBarry Smith   Logically Collective
33249517cdeSBarry Smith 
33349517cdeSBarry Smith   Input Parameter:
33449517cdeSBarry Smith . pc - the preconditioner context
33549517cdeSBarry Smith 
33649517cdeSBarry Smith   Output Parameter:
337f1580f4eSBarry Smith . flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)
33849517cdeSBarry Smith 
33920f4b53cSBarry Smith   Level: intermediate
34020f4b53cSBarry Smith 
341f1580f4eSBarry Smith   Note:
34249517cdeSBarry Smith   For the common case in which the linear system matrix and the matrix used to construct the
34349517cdeSBarry Smith   preconditioner are identical, this routine is does nothing.
34449517cdeSBarry Smith 
345f1580f4eSBarry Smith .seealso: `PC`, `PCSetUseAmat()`, `PCBJACOBI`, `PGMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`
34649517cdeSBarry Smith @*/
347d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetUseAmat(PC pc, PetscBool *flg)
348d71ae5a4SJacob Faibussowitsch {
34949517cdeSBarry Smith   PetscFunctionBegin;
35049517cdeSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
35149517cdeSBarry Smith   *flg = pc->useAmat;
3523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35349517cdeSBarry Smith }
35449517cdeSBarry Smith 
355f39d8e23SSatish Balay /*@
3563821be0aSBarry Smith   PCSetKSPNestLevel - sets the amount of nesting the `KSP` that contains this `PC` has
3573821be0aSBarry Smith 
3583821be0aSBarry Smith   Collective
3593821be0aSBarry Smith 
3603821be0aSBarry Smith   Input Parameters:
3613821be0aSBarry Smith + pc    - the `PC`
3623821be0aSBarry Smith - level - the nest level
3633821be0aSBarry Smith 
3643821be0aSBarry Smith   Level: developer
3653821be0aSBarry Smith 
3663821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCGetKSPNestLevel()`, `KSPSetNestLevel()`
3673821be0aSBarry Smith @*/
3683821be0aSBarry Smith PetscErrorCode PCSetKSPNestLevel(PC pc, PetscInt level)
3693821be0aSBarry Smith {
3703821be0aSBarry Smith   PetscFunctionBegin;
3717a99bfcaSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
3727a99bfcaSBarry Smith   PetscValidLogicalCollectiveInt(pc, level, 2);
3733821be0aSBarry Smith   pc->kspnestlevel = level;
3743821be0aSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3753821be0aSBarry Smith }
3763821be0aSBarry Smith 
3773821be0aSBarry Smith /*@
3783821be0aSBarry Smith   PCGetKSPNestLevel - gets the amount of nesting the `KSP` that contains this `PC` has
3793821be0aSBarry Smith 
3807a99bfcaSBarry Smith   Not Collective
3813821be0aSBarry Smith 
3823821be0aSBarry Smith   Input Parameter:
3833821be0aSBarry Smith . pc - the `PC`
3843821be0aSBarry Smith 
3853821be0aSBarry Smith   Output Parameter:
3863821be0aSBarry Smith . level - the nest level
3873821be0aSBarry Smith 
3883821be0aSBarry Smith   Level: developer
3893821be0aSBarry Smith 
3903821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `KSPGetNestLevel()`
3913821be0aSBarry Smith @*/
3923821be0aSBarry Smith PetscErrorCode PCGetKSPNestLevel(PC pc, PetscInt *level)
3933821be0aSBarry Smith {
3943821be0aSBarry Smith   PetscFunctionBegin;
3957a99bfcaSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
3967a99bfcaSBarry Smith   PetscAssertPointer(level, 2);
3973821be0aSBarry Smith   *level = pc->kspnestlevel;
3983821be0aSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3993821be0aSBarry Smith }
4003821be0aSBarry Smith 
4013821be0aSBarry Smith /*@
402f1580f4eSBarry Smith   PCCreate - Creates a preconditioner context, `PC`
4034b9ad928SBarry Smith 
404d083f849SBarry Smith   Collective
4054b9ad928SBarry Smith 
4064b9ad928SBarry Smith   Input Parameter:
4074b9ad928SBarry Smith . comm - MPI communicator
4084b9ad928SBarry Smith 
4094b9ad928SBarry Smith   Output Parameter:
4107a99bfcaSBarry Smith . newpc - location to put the preconditioner context
4114b9ad928SBarry Smith 
41220f4b53cSBarry Smith   Level: developer
41320f4b53cSBarry Smith 
414f1580f4eSBarry Smith   Note:
415f1580f4eSBarry 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`
416f1580f4eSBarry Smith   in parallel. For dense matrices it is always `PCNONE`.
4174b9ad928SBarry Smith 
418f1580f4eSBarry Smith .seealso: `PC`, `PCSetUp()`, `PCApply()`, `PCDestroy()`
4194b9ad928SBarry Smith @*/
420d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreate(MPI_Comm comm, PC *newpc)
421d71ae5a4SJacob Faibussowitsch {
4224b9ad928SBarry Smith   PC pc;
4234b9ad928SBarry Smith 
4244b9ad928SBarry Smith   PetscFunctionBegin;
4254f572ea9SToby Isaac   PetscAssertPointer(newpc, 2);
4260a545947SLisandro Dalcin   *newpc = NULL;
4279566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
4284b9ad928SBarry Smith 
4299566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(pc, PC_CLASSID, "PC", "Preconditioner", "PC", comm, PCDestroy, PCView));
430a7d21a52SLisandro Dalcin 
4310a545947SLisandro Dalcin   pc->mat                  = NULL;
4320a545947SLisandro Dalcin   pc->pmat                 = NULL;
4334b9ad928SBarry Smith   pc->setupcalled          = 0;
4340c24e6a1SHong Zhang   pc->setfromoptionscalled = 0;
4350a545947SLisandro Dalcin   pc->data                 = NULL;
4364b9ad928SBarry Smith   pc->diagonalscale        = PETSC_FALSE;
4370a545947SLisandro Dalcin   pc->diagonalscaleleft    = NULL;
4380a545947SLisandro Dalcin   pc->diagonalscaleright   = NULL;
4394b9ad928SBarry Smith 
4400a545947SLisandro Dalcin   pc->modifysubmatrices  = NULL;
4410a545947SLisandro Dalcin   pc->modifysubmatricesP = NULL;
4422fa5cd67SKarl Rupp 
443a7d21a52SLisandro Dalcin   *newpc = pc;
4443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4454b9ad928SBarry Smith }
4464b9ad928SBarry Smith 
4474b9ad928SBarry Smith /*@
4484b9ad928SBarry Smith   PCApply - Applies the preconditioner to a vector.
4494b9ad928SBarry Smith 
450c3339decSBarry Smith   Collective
4514b9ad928SBarry Smith 
4524b9ad928SBarry Smith   Input Parameters:
4534b9ad928SBarry Smith + pc - the preconditioner context
4544b9ad928SBarry Smith - x  - input vector
4554b9ad928SBarry Smith 
4564b9ad928SBarry Smith   Output Parameter:
4574b9ad928SBarry Smith . y - output vector
4584b9ad928SBarry Smith 
4594b9ad928SBarry Smith   Level: developer
4604b9ad928SBarry Smith 
461f1580f4eSBarry Smith .seealso: `PC`, `PCApplyTranspose()`, `PCApplyBAorAB()`
4624b9ad928SBarry Smith @*/
463d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApply(PC pc, Vec x, Vec y)
464d71ae5a4SJacob Faibussowitsch {
465106e20bfSBarry Smith   PetscInt m, n, mv, nv;
4664b9ad928SBarry Smith 
4674b9ad928SBarry Smith   PetscFunctionBegin;
4680700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4690700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
4700700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
4717827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
472e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
473540bdfbdSVaclav Hapla   /* use pmat to check vector sizes since for KSPLSQR the pmat may be of a different size than mat */
4749566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &m, &n));
4759566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x, &mv));
4769566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(y, &nv));
477540bdfbdSVaclav Hapla   /* check pmat * y = x is feasible */
47863a3b9bcSJacob 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);
47963a3b9bcSJacob 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);
4809566063dSJacob Faibussowitsch   PetscCall(VecSetErrorIfLocked(y, 3));
481106e20bfSBarry Smith 
4829566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
4839566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
4849566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
485dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, apply, x, y);
4869566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
487e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
4889566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
4893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4904b9ad928SBarry Smith }
4914b9ad928SBarry Smith 
4924b9ad928SBarry Smith /*@
493f1580f4eSBarry Smith   PCMatApply - Applies the preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApply()`, Y and X must be different matrices.
494c41ea70eSPierre Jolivet 
495c3339decSBarry Smith   Collective
496c677e75fSPierre Jolivet 
497c677e75fSPierre Jolivet   Input Parameters:
498c677e75fSPierre Jolivet + pc - the preconditioner context
499c677e75fSPierre Jolivet - X  - block of input vectors
500c677e75fSPierre Jolivet 
501c677e75fSPierre Jolivet   Output Parameter:
502c677e75fSPierre Jolivet . Y - block of output vectors
503c677e75fSPierre Jolivet 
504c677e75fSPierre Jolivet   Level: developer
505c677e75fSPierre Jolivet 
506f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `KSPMatSolve()`
507c677e75fSPierre Jolivet @*/
508d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMatApply(PC pc, Mat X, Mat Y)
509d71ae5a4SJacob Faibussowitsch {
510c677e75fSPierre Jolivet   Mat       A;
511c677e75fSPierre Jolivet   Vec       cy, cx;
512bd82155bSPierre Jolivet   PetscInt  m1, M1, m2, M2, n1, N1, n2, N2, m3, M3, n3, N3;
513c677e75fSPierre Jolivet   PetscBool match;
514c677e75fSPierre Jolivet 
515c677e75fSPierre Jolivet   PetscFunctionBegin;
516c677e75fSPierre Jolivet   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
517e2b200f6SPierre Jolivet   PetscValidHeaderSpecific(X, MAT_CLASSID, 2);
518e2b200f6SPierre Jolivet   PetscValidHeaderSpecific(Y, MAT_CLASSID, 3);
519e2b200f6SPierre Jolivet   PetscCheckSameComm(pc, 1, X, 2);
520e2b200f6SPierre Jolivet   PetscCheckSameComm(pc, 1, Y, 3);
5217827d75bSBarry Smith   PetscCheck(Y != X, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "Y and X must be different matrices");
5229566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &A));
5239566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m3, &n3));
5249566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m2, &n2));
5259566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &m1, &n1));
5269566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M3, &N3));
5279566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &M2, &N2));
5289566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Y, &M1, &N1));
52963a3b9bcSJacob 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);
53063a3b9bcSJacob 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);
53163a3b9bcSJacob 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);
5329566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &match, MATSEQDENSE, MATMPIDENSE, ""));
53328b400f6SJacob Faibussowitsch   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of output vectors not stored in a dense Mat");
5349566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
53528b400f6SJacob Faibussowitsch   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of input vectors not stored in a dense Mat");
5369566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
537c677e75fSPierre Jolivet   if (pc->ops->matapply) {
5389566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0));
539dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, matapply, X, Y);
5409566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0));
541c677e75fSPierre Jolivet   } else {
5429566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "PC type %s applying column by column\n", ((PetscObject)pc)->type_name));
543bd82155bSPierre Jolivet     for (n1 = 0; n1 < N1; ++n1) {
5449566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecRead(X, n1, &cx));
5459566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecWrite(Y, n1, &cy));
5469566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, cx, cy));
5479566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecWrite(Y, n1, &cy));
5489566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecRead(X, n1, &cx));
549c677e75fSPierre Jolivet     }
550c677e75fSPierre Jolivet   }
5513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
552c677e75fSPierre Jolivet }
553c677e75fSPierre Jolivet 
554c677e75fSPierre Jolivet /*@
5554b9ad928SBarry Smith   PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
5564b9ad928SBarry Smith 
557c3339decSBarry Smith   Collective
5584b9ad928SBarry Smith 
5594b9ad928SBarry Smith   Input Parameters:
5604b9ad928SBarry Smith + pc - the preconditioner context
5614b9ad928SBarry Smith - x  - input vector
5624b9ad928SBarry Smith 
5634b9ad928SBarry Smith   Output Parameter:
5644b9ad928SBarry Smith . y - output vector
5654b9ad928SBarry Smith 
56620f4b53cSBarry Smith   Level: developer
56720f4b53cSBarry Smith 
568f1580f4eSBarry Smith   Note:
569f1580f4eSBarry Smith   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
5704b9ad928SBarry Smith 
571f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplySymmetricRight()`
5724b9ad928SBarry Smith @*/
573d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricLeft(PC pc, Vec x, Vec y)
574d71ae5a4SJacob Faibussowitsch {
5754b9ad928SBarry Smith   PetscFunctionBegin;
5760700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
5770700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
5780700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
5797827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
580e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
5819566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
5829566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
5839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplySymmetricLeft, pc, x, y, 0));
584dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applysymmetricleft, x, y);
5859566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplySymmetricLeft, pc, x, y, 0));
5869566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
587e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
5883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5894b9ad928SBarry Smith }
5904b9ad928SBarry Smith 
5914b9ad928SBarry Smith /*@
5924b9ad928SBarry Smith   PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
5934b9ad928SBarry Smith 
594c3339decSBarry Smith   Collective
5954b9ad928SBarry Smith 
5964b9ad928SBarry Smith   Input Parameters:
5974b9ad928SBarry Smith + pc - the preconditioner context
5984b9ad928SBarry Smith - x  - input vector
5994b9ad928SBarry Smith 
6004b9ad928SBarry Smith   Output Parameter:
6014b9ad928SBarry Smith . y - output vector
6024b9ad928SBarry Smith 
6034b9ad928SBarry Smith   Level: developer
6044b9ad928SBarry Smith 
605f1580f4eSBarry Smith   Note:
606f1580f4eSBarry Smith   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
6074b9ad928SBarry Smith 
608f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplySymmetricLeft()`
6094b9ad928SBarry Smith @*/
610d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricRight(PC pc, Vec x, Vec y)
611d71ae5a4SJacob Faibussowitsch {
6124b9ad928SBarry Smith   PetscFunctionBegin;
6130700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6140700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6150700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6167827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
617e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
6189566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
6199566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
6209566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplySymmetricRight, pc, x, y, 0));
621dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applysymmetricright, x, y);
6229566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplySymmetricRight, pc, x, y, 0));
6239566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
624e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
6253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6264b9ad928SBarry Smith }
6274b9ad928SBarry Smith 
6284b9ad928SBarry Smith /*@
6294b9ad928SBarry Smith   PCApplyTranspose - Applies the transpose of preconditioner to a vector.
6304b9ad928SBarry Smith 
631c3339decSBarry Smith   Collective
6324b9ad928SBarry Smith 
6334b9ad928SBarry Smith   Input Parameters:
6344b9ad928SBarry Smith + pc - the preconditioner context
6354b9ad928SBarry Smith - x  - input vector
6364b9ad928SBarry Smith 
6374b9ad928SBarry Smith   Output Parameter:
6384b9ad928SBarry Smith . y - output vector
6394b9ad928SBarry Smith 
64020f4b53cSBarry Smith   Level: developer
64120f4b53cSBarry Smith 
642f1580f4eSBarry Smith   Note:
64395452b02SPatrick Sanan   For complex numbers this applies the non-Hermitian transpose.
6444c97465dSBarry Smith 
645feefa0e1SJacob Faibussowitsch   Developer Notes:
646f1580f4eSBarry Smith   We need to implement a `PCApplyHermitianTranspose()`
6474c97465dSBarry Smith 
648f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplyBAorAB()`, `PCApplyBAorABTranspose()`, `PCApplyTransposeExists()`
6494b9ad928SBarry Smith @*/
650d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTranspose(PC pc, Vec x, Vec y)
651d71ae5a4SJacob Faibussowitsch {
6524b9ad928SBarry Smith   PetscFunctionBegin;
6530700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6540700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6550700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6567827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
657e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
6589566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
6599566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
6609566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
661dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applytranspose, x, y);
6629566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
6639566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
664e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
6653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6664b9ad928SBarry Smith }
6674b9ad928SBarry Smith 
6684b9ad928SBarry Smith /*@
6699cc050e5SLisandro Dalcin   PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
6704b9ad928SBarry Smith 
671c3339decSBarry Smith   Collective
6724b9ad928SBarry Smith 
673f1580f4eSBarry Smith   Input Parameter:
6744b9ad928SBarry Smith . pc - the preconditioner context
6754b9ad928SBarry Smith 
6764b9ad928SBarry Smith   Output Parameter:
677f1580f4eSBarry Smith . flg - `PETSC_TRUE` if a transpose operation is defined
6784b9ad928SBarry Smith 
6794b9ad928SBarry Smith   Level: developer
6804b9ad928SBarry Smith 
681f1580f4eSBarry Smith .seealso: `PC`, `PCApplyTranspose()`
6824b9ad928SBarry Smith @*/
683d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTransposeExists(PC pc, PetscBool *flg)
684d71ae5a4SJacob Faibussowitsch {
6854b9ad928SBarry Smith   PetscFunctionBegin;
6860700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6874f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
6886ce5e81cSLisandro Dalcin   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
6896ce5e81cSLisandro Dalcin   else *flg = PETSC_FALSE;
6903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6914b9ad928SBarry Smith }
6924b9ad928SBarry Smith 
6934b9ad928SBarry Smith /*@
6941a2f9f99SBarry Smith   PCApplyBAorAB - Applies the preconditioner and operator to a vector. y = B*A*x or y = A*B*x.
6954b9ad928SBarry Smith 
696c3339decSBarry Smith   Collective
6974b9ad928SBarry Smith 
6984b9ad928SBarry Smith   Input Parameters:
6994b9ad928SBarry Smith + pc   - the preconditioner context
700f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
7014b9ad928SBarry Smith . x    - input vector
7024b9ad928SBarry Smith - work - work vector
7034b9ad928SBarry Smith 
7044b9ad928SBarry Smith   Output Parameter:
7054b9ad928SBarry Smith . y - output vector
7064b9ad928SBarry Smith 
7074b9ad928SBarry Smith   Level: developer
7084b9ad928SBarry Smith 
709f1580f4eSBarry Smith   Note:
710f1580f4eSBarry 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. Note that the
711f1580f4eSBarry Smith   specific `KSPSolve()` method must also be written to handle the post-solve "correction" for the diagonal scaling.
71213e0d083SBarry Smith 
713f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorABTranspose()`
7144b9ad928SBarry Smith @*/
715d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorAB(PC pc, PCSide side, Vec x, Vec y, Vec work)
716d71ae5a4SJacob Faibussowitsch {
7174b9ad928SBarry Smith   PetscFunctionBegin;
7180700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
719a69c7061SStefano Zampini   PetscValidLogicalCollectiveEnum(pc, side, 2);
7200700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
7210700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
7220700a824SBarry Smith   PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
723a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, x, 3);
724a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, y, 4);
725a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, work, 5);
7267827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
7277827d75bSBarry 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");
7287827d75bSBarry Smith   PetscCheck(!pc->diagonalscale || side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot include diagonal scaling with symmetric preconditioner application");
729e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
7304b9ad928SBarry Smith 
7319566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
7324b9ad928SBarry Smith   if (pc->diagonalscale) {
7334b9ad928SBarry Smith     if (pc->ops->applyBA) {
7344b9ad928SBarry Smith       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
7359566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &work2));
7369566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, work2));
737dbbe0bcdSBarry Smith       PetscUseTypeMethod(pc, applyBA, side, work2, y, work);
7389566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7399566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&work2));
7404b9ad928SBarry Smith     } else if (side == PC_RIGHT) {
7419566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, y));
7429566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, y, work));
7439566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7449566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7454b9ad928SBarry Smith     } else if (side == PC_LEFT) {
7469566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, y));
7479566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, y, work));
7489566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, work, y));
7499566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7507827d75bSBarry Smith     } else PetscCheck(side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot provide diagonal scaling with symmetric application of preconditioner");
7514b9ad928SBarry Smith   } else {
7524b9ad928SBarry Smith     if (pc->ops->applyBA) {
753dbbe0bcdSBarry Smith       PetscUseTypeMethod(pc, applyBA, side, x, y, work);
7544b9ad928SBarry Smith     } else if (side == PC_RIGHT) {
7559566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, x, work));
7569566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7574b9ad928SBarry Smith     } else if (side == PC_LEFT) {
7589566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, x, work));
7599566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, work, y));
7604b9ad928SBarry Smith     } else if (side == PC_SYMMETRIC) {
7614b9ad928SBarry Smith       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
7629566063dSJacob Faibussowitsch       PetscCall(PCApplySymmetricRight(pc, x, work));
7639566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7649566063dSJacob Faibussowitsch       PetscCall(VecCopy(y, work));
7659566063dSJacob Faibussowitsch       PetscCall(PCApplySymmetricLeft(pc, work, y));
7664b9ad928SBarry Smith     }
7674b9ad928SBarry Smith   }
768e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
7693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7704b9ad928SBarry Smith }
7714b9ad928SBarry Smith 
7724b9ad928SBarry Smith /*@
7734b9ad928SBarry Smith   PCApplyBAorABTranspose - Applies the transpose of the preconditioner
7744b9ad928SBarry Smith   and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
7759b3038f0SBarry Smith   NOT tr(B*A) = tr(A)*tr(B).
7764b9ad928SBarry Smith 
777c3339decSBarry Smith   Collective
7784b9ad928SBarry Smith 
7794b9ad928SBarry Smith   Input Parameters:
7804b9ad928SBarry Smith + pc   - the preconditioner context
781f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
7824b9ad928SBarry Smith . x    - input vector
7834b9ad928SBarry Smith - work - work vector
7844b9ad928SBarry Smith 
7854b9ad928SBarry Smith   Output Parameter:
7864b9ad928SBarry Smith . y - output vector
7874b9ad928SBarry Smith 
78820f4b53cSBarry Smith   Level: developer
78920f4b53cSBarry Smith 
790f1580f4eSBarry Smith   Note:
791f1580f4eSBarry Smith   This routine is used internally so that the same Krylov code can be used to solve A x = b and A' x = b, with a preconditioner
7929b3038f0SBarry Smith   defined by B'. This is why this has the funny form that it computes tr(B) * tr(A)
7939b3038f0SBarry Smith 
794f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorAB()`
7954b9ad928SBarry Smith @*/
796d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorABTranspose(PC pc, PCSide side, Vec x, Vec y, Vec work)
797d71ae5a4SJacob Faibussowitsch {
7984b9ad928SBarry Smith   PetscFunctionBegin;
7990700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8000700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
8010700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
8020700a824SBarry Smith   PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
8037827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
804e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
8054b9ad928SBarry Smith   if (pc->ops->applyBAtranspose) {
806dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, applyBAtranspose, side, x, y, work);
807e0f629ddSJacob Faibussowitsch     if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8083ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
8094b9ad928SBarry Smith   }
8107827d75bSBarry Smith   PetscCheck(side == PC_LEFT || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right or left");
8114b9ad928SBarry Smith 
8129566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
8134b9ad928SBarry Smith   if (side == PC_RIGHT) {
8149566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(pc, x, work));
8159566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->mat, work, y));
8164b9ad928SBarry Smith   } else if (side == PC_LEFT) {
8179566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->mat, x, work));
8189566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(pc, work, y));
8194b9ad928SBarry Smith   }
8204b9ad928SBarry Smith   /* add support for PC_SYMMETRIC */
821e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8234b9ad928SBarry Smith }
8244b9ad928SBarry Smith 
8254b9ad928SBarry Smith /*@
8264b9ad928SBarry Smith   PCApplyRichardsonExists - Determines whether a particular preconditioner has a
8274b9ad928SBarry Smith   built-in fast application of Richardson's method.
8284b9ad928SBarry Smith 
8294b9ad928SBarry Smith   Not Collective
8304b9ad928SBarry Smith 
8314b9ad928SBarry Smith   Input Parameter:
8324b9ad928SBarry Smith . pc - the preconditioner
8334b9ad928SBarry Smith 
8344b9ad928SBarry Smith   Output Parameter:
835f1580f4eSBarry Smith . exists - `PETSC_TRUE` or `PETSC_FALSE`
8364b9ad928SBarry Smith 
8374b9ad928SBarry Smith   Level: developer
8384b9ad928SBarry Smith 
839f1580f4eSBarry Smith .seealso: `PC`, `PCRICHARDSON`, `PCApplyRichardson()`
8404b9ad928SBarry Smith @*/
841d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyRichardsonExists(PC pc, PetscBool *exists)
842d71ae5a4SJacob Faibussowitsch {
8434b9ad928SBarry Smith   PetscFunctionBegin;
8440700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8454f572ea9SToby Isaac   PetscAssertPointer(exists, 2);
8464b9ad928SBarry Smith   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
8474b9ad928SBarry Smith   else *exists = PETSC_FALSE;
8483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8494b9ad928SBarry Smith }
8504b9ad928SBarry Smith 
8514b9ad928SBarry Smith /*@
8524b9ad928SBarry Smith   PCApplyRichardson - Applies several steps of Richardson iteration with
8534b9ad928SBarry Smith   the particular preconditioner. This routine is usually used by the
8544b9ad928SBarry Smith   Krylov solvers and not the application code directly.
8554b9ad928SBarry Smith 
856c3339decSBarry Smith   Collective
8574b9ad928SBarry Smith 
8584b9ad928SBarry Smith   Input Parameters:
8594b9ad928SBarry Smith + pc        - the preconditioner context
8607319c654SBarry Smith . b         - the right hand side
8614b9ad928SBarry Smith . w         - one work vector
8624b9ad928SBarry Smith . rtol      - relative decrease in residual norm convergence criteria
86370441072SBarry Smith . abstol    - absolute residual norm convergence criteria
8644b9ad928SBarry Smith . dtol      - divergence residual norm increase criteria
8657319c654SBarry Smith . its       - the number of iterations to apply.
8667319c654SBarry Smith - guesszero - if the input x contains nonzero initial guess
8674b9ad928SBarry Smith 
868d8d19677SJose E. Roman   Output Parameters:
8694d0a8057SBarry Smith + outits - number of iterations actually used (for SOR this always equals its)
8704d0a8057SBarry Smith . reason - the reason the apply terminated
871f1580f4eSBarry Smith - y      - the solution (also contains initial guess if guesszero is `PETSC_FALSE`
8724b9ad928SBarry Smith 
87320f4b53cSBarry Smith   Level: developer
87420f4b53cSBarry Smith 
8754b9ad928SBarry Smith   Notes:
8764b9ad928SBarry Smith   Most preconditioners do not support this function. Use the command
877f1580f4eSBarry Smith   `PCApplyRichardsonExists()` to determine if one does.
8784b9ad928SBarry Smith 
879f1580f4eSBarry Smith   Except for the `PCMG` this routine ignores the convergence tolerances
8804b9ad928SBarry Smith   and always runs for the number of iterations
8814b9ad928SBarry Smith 
882f1580f4eSBarry Smith .seealso: `PC`, `PCApplyRichardsonExists()`
8834b9ad928SBarry Smith @*/
884d71ae5a4SJacob 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)
885d71ae5a4SJacob Faibussowitsch {
8864b9ad928SBarry Smith   PetscFunctionBegin;
8870700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8880700a824SBarry Smith   PetscValidHeaderSpecific(b, VEC_CLASSID, 2);
8890700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
8900700a824SBarry Smith   PetscValidHeaderSpecific(w, VEC_CLASSID, 4);
8917827d75bSBarry Smith   PetscCheck(b != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "b and y must be different vectors");
8929566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
893dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applyrichardson, b, y, w, rtol, abstol, dtol, its, guesszero, outits, reason);
8943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8954b9ad928SBarry Smith }
8964b9ad928SBarry Smith 
897422a814eSBarry Smith /*@
898f1580f4eSBarry Smith   PCSetFailedReason - Sets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
8991b2b9847SBarry Smith 
900c3339decSBarry Smith   Logically Collective
9011b2b9847SBarry Smith 
902d8d19677SJose E. Roman   Input Parameters:
9031b2b9847SBarry Smith + pc     - the preconditioner context
9041b2b9847SBarry Smith - reason - the reason it failedx
9051b2b9847SBarry Smith 
9061b2b9847SBarry Smith   Level: advanced
9071b2b9847SBarry Smith 
908f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCFailedReason`
9091b2b9847SBarry Smith @*/
910d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFailedReason(PC pc, PCFailedReason reason)
911d71ae5a4SJacob Faibussowitsch {
9121b2b9847SBarry Smith   PetscFunctionBegin;
9136479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9141b2b9847SBarry Smith   pc->failedreason = reason;
9153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9161b2b9847SBarry Smith }
9171b2b9847SBarry Smith 
9181b2b9847SBarry Smith /*@
919f1580f4eSBarry Smith   PCGetFailedReason - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
920422a814eSBarry Smith 
921c3339decSBarry Smith   Logically Collective
922422a814eSBarry Smith 
923422a814eSBarry Smith   Input Parameter:
924422a814eSBarry Smith . pc - the preconditioner context
925422a814eSBarry Smith 
926422a814eSBarry Smith   Output Parameter:
9271b2b9847SBarry Smith . reason - the reason it failed
928422a814eSBarry Smith 
929422a814eSBarry Smith   Level: advanced
930422a814eSBarry Smith 
931f1580f4eSBarry Smith   Note:
932f1580f4eSBarry Smith   This is the maximum over reason over all ranks in the PC communicator. It is only valid after
9336479e835SStefano Zampini   a call `KSPCheckDot()` or  `KSPCheckNorm()` inside a `KSPSolve()` or `PCReduceFailedReason()`.
9346479e835SStefano Zampini   It is not valid immediately after a `PCSetUp()` or `PCApply()`, then use `PCGetFailedReasonRank()`
9351b2b9847SBarry Smith 
936feefa0e1SJacob Faibussowitsch .seealso: `PC`, ``PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReasonRank()`, `PCSetFailedReason()`
937422a814eSBarry Smith @*/
938d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason)
939d71ae5a4SJacob Faibussowitsch {
940422a814eSBarry Smith   PetscFunctionBegin;
9416479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9429754764cSHong Zhang   if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled;
9430ae0b32dSHong Zhang   else *reason = pc->failedreason;
9443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
945422a814eSBarry Smith }
946422a814eSBarry Smith 
9471b2b9847SBarry Smith /*@
948f1580f4eSBarry Smith   PCGetFailedReasonRank - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail on this MPI rank
9491b2b9847SBarry Smith 
950f1580f4eSBarry Smith   Not Collective
9511b2b9847SBarry Smith 
9521b2b9847SBarry Smith   Input Parameter:
9531b2b9847SBarry Smith . pc - the preconditioner context
9541b2b9847SBarry Smith 
9551b2b9847SBarry Smith   Output Parameter:
9561b2b9847SBarry Smith . reason - the reason it failed
9571b2b9847SBarry Smith 
95820f4b53cSBarry Smith   Level: advanced
95920f4b53cSBarry Smith 
960f1580f4eSBarry Smith   Note:
961f1580f4eSBarry Smith   Different ranks may have different reasons or no reason, see `PCGetFailedReason()`
9621b2b9847SBarry Smith 
9636479e835SStefano Zampini .seealso: `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`, `PCReduceFailedReason()`
9641b2b9847SBarry Smith @*/
965d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReasonRank(PC pc, PCFailedReason *reason)
966d71ae5a4SJacob Faibussowitsch {
9671b2b9847SBarry Smith   PetscFunctionBegin;
9686479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9691b2b9847SBarry Smith   if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled;
9701b2b9847SBarry Smith   else *reason = pc->failedreason;
9713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9721b2b9847SBarry Smith }
973422a814eSBarry Smith 
9746479e835SStefano Zampini /*@
9756479e835SStefano Zampini   PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC`
9766479e835SStefano Zampini 
9776479e835SStefano Zampini   Collective
9786479e835SStefano Zampini 
9796479e835SStefano Zampini   Input Parameter:
9806479e835SStefano Zampini . pc - the preconditioner context
9816479e835SStefano Zampini 
9826479e835SStefano Zampini   Level: advanced
9836479e835SStefano Zampini 
9846479e835SStefano Zampini   Note:
9856479e835SStefano Zampini   Different MPI ranks may have different reasons or no reason, see `PCGetFailedReason()`. This routine
9866479e835SStefano Zampini   makes them have a common value (failure if any MPI process had a failure).
9876479e835SStefano Zampini 
9886479e835SStefano Zampini .seealso: `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`
9896479e835SStefano Zampini @*/
9906479e835SStefano Zampini PetscErrorCode PCReduceFailedReason(PC pc)
9916479e835SStefano Zampini {
9926479e835SStefano Zampini   PetscInt buf;
9936479e835SStefano Zampini 
9946479e835SStefano Zampini   PetscFunctionBegin;
9956479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9966479e835SStefano Zampini   buf = (PetscInt)pc->failedreason;
9976479e835SStefano Zampini   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
9986479e835SStefano Zampini   pc->failedreason = (PCFailedReason)buf;
9996479e835SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
10006479e835SStefano Zampini }
10016479e835SStefano Zampini 
1002c00cb57fSBarry Smith /*  Next line needed to deactivate KSP_Solve logging */
1003c00cb57fSBarry Smith #include <petsc/private/kspimpl.h>
1004c00cb57fSBarry Smith 
10054b9ad928SBarry Smith /*
10064b9ad928SBarry Smith       a setupcall of 0 indicates never setup,
100723ee1639SBarry Smith                      1 indicates has been previously setup
1008422a814eSBarry Smith                     -1 indicates a PCSetUp() was attempted and failed
10094b9ad928SBarry Smith */
10104b9ad928SBarry Smith /*@
10114b9ad928SBarry Smith   PCSetUp - Prepares for the use of a preconditioner.
10124b9ad928SBarry Smith 
1013c3339decSBarry Smith   Collective
10144b9ad928SBarry Smith 
10154b9ad928SBarry Smith   Input Parameter:
10164b9ad928SBarry Smith . pc - the preconditioner context
10174b9ad928SBarry Smith 
10184b9ad928SBarry Smith   Level: developer
10194b9ad928SBarry Smith 
1020f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`
10214b9ad928SBarry Smith @*/
1022d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp(PC pc)
1023d71ae5a4SJacob Faibussowitsch {
1024566e8bf2SBarry Smith   const char      *def;
1025fc9b51d3SBarry Smith   PetscObjectState matstate, matnonzerostate;
10264b9ad928SBarry Smith 
10274b9ad928SBarry Smith   PetscFunctionBegin;
10280700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
102928b400f6SJacob Faibussowitsch   PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first");
10306ce5e81cSLisandro Dalcin 
103123ee1639SBarry Smith   if (pc->setupcalled && pc->reusepreconditioner) {
10329566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n"));
10333ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10344b9ad928SBarry Smith   }
10354b9ad928SBarry Smith 
10369566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate));
10379566063dSJacob Faibussowitsch   PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate));
103823ee1639SBarry Smith   if (!pc->setupcalled) {
10395e62d202SMark Adams     //PetscCall(PetscInfo(pc, "Setting up PC for first time\n"));
104023ee1639SBarry Smith     pc->flag = DIFFERENT_NONZERO_PATTERN;
10419df67409SStefano Zampini   } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS);
10429df67409SStefano Zampini   else {
10439df67409SStefano Zampini     if (matnonzerostate != pc->matnonzerostate) {
10449566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n"));
104523ee1639SBarry Smith       pc->flag = DIFFERENT_NONZERO_PATTERN;
104623ee1639SBarry Smith     } else {
10475e62d202SMark Adams       //PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n"));
104823ee1639SBarry Smith       pc->flag = SAME_NONZERO_PATTERN;
104923ee1639SBarry Smith     }
105023ee1639SBarry Smith   }
105123ee1639SBarry Smith   pc->matstate        = matstate;
1052fc9b51d3SBarry Smith   pc->matnonzerostate = matnonzerostate;
105323ee1639SBarry Smith 
10547adad957SLisandro Dalcin   if (!((PetscObject)pc)->type_name) {
10559566063dSJacob Faibussowitsch     PetscCall(PCGetDefaultType_Private(pc, &def));
10569566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, def));
1057566e8bf2SBarry Smith   }
10584b9ad928SBarry Smith 
10599566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
10609566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure));
10619566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0));
10624b9ad928SBarry Smith   if (pc->ops->setup) {
1063c00cb57fSBarry Smith     /* do not log solves and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */
10649566063dSJacob Faibussowitsch     PetscCall(KSPInitializePackage());
10659566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePush(KSP_Solve));
10669566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePush(PC_Apply));
1067dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, setup);
10689566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePop(KSP_Solve));
10699566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePop(PC_Apply));
10704b9ad928SBarry Smith   }
10719566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0));
1072422a814eSBarry Smith   if (!pc->setupcalled) pc->setupcalled = 1;
10733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10744b9ad928SBarry Smith }
10754b9ad928SBarry Smith 
10764b9ad928SBarry Smith /*@
10774b9ad928SBarry Smith   PCSetUpOnBlocks - Sets up the preconditioner for each block in
10784b9ad928SBarry Smith   the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
10794b9ad928SBarry Smith   methods.
10804b9ad928SBarry Smith 
1081c3339decSBarry Smith   Collective
10824b9ad928SBarry Smith 
1083f1580f4eSBarry Smith   Input Parameter:
10844b9ad928SBarry Smith . pc - the preconditioner context
10854b9ad928SBarry Smith 
10864b9ad928SBarry Smith   Level: developer
10874b9ad928SBarry Smith 
1088f1580f4eSBarry Smith   Note:
1089f1580f4eSBarry Smith   For nested preconditioners such as `PCBJACOBI` `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is
1090f1580f4eSBarry Smith   called on the outer `PC`, this routine ensures it is called.
1091f1580f4eSBarry Smith 
1092feefa0e1SJacob Faibussowitsch .seealso: `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()`
10934b9ad928SBarry Smith @*/
1094d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUpOnBlocks(PC pc)
1095d71ae5a4SJacob Faibussowitsch {
10964b9ad928SBarry Smith   PetscFunctionBegin;
10970700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
10983ba16761SJacob Faibussowitsch   if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS);
10999566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0));
1100dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, setuponblocks);
11019566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0));
11023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11034b9ad928SBarry Smith }
11044b9ad928SBarry Smith 
11054b9ad928SBarry Smith /*@C
11064b9ad928SBarry Smith   PCSetModifySubMatrices - Sets a user-defined routine for modifying the
110704c3f3b8SBarry Smith   submatrices that arise within certain subdomain-based preconditioners such as `PCASM`
11084b9ad928SBarry Smith 
1109c3339decSBarry Smith   Logically Collective
11104b9ad928SBarry Smith 
11114b9ad928SBarry Smith   Input Parameters:
11124b9ad928SBarry Smith + pc   - the preconditioner context
11134b9ad928SBarry Smith . func - routine for modifying the submatrices
111404c3f3b8SBarry Smith - ctx  - optional user-defined context (may be `NULL`)
11154b9ad928SBarry Smith 
111620f4b53cSBarry Smith   Calling sequence of `func`:
111720f4b53cSBarry Smith + pc     - the preconditioner context
111804c3f3b8SBarry Smith . nsub   - number of index sets
111920f4b53cSBarry Smith . row    - an array of index sets that contain the global row numbers
11204b9ad928SBarry Smith          that comprise each local submatrix
11214b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
11224b9ad928SBarry Smith          that comprise each local submatrix
11234b9ad928SBarry Smith . submat - array of local submatrices
11244b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
112504c3f3b8SBarry Smith          user-defined func routine (may be `NULL`)
11264b9ad928SBarry Smith 
112720f4b53cSBarry Smith   Level: advanced
112820f4b53cSBarry Smith 
11294b9ad928SBarry Smith   Notes:
113004c3f3b8SBarry Smith   The basic submatrices are extracted from the preconditioner matrix as
113104c3f3b8SBarry Smith   usual; the user can then alter these (for example, to set different boundary
113204c3f3b8SBarry Smith   conditions for each submatrix) before they are used for the local solves.
113304c3f3b8SBarry Smith 
1134f1580f4eSBarry Smith   `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and
1135f1580f4eSBarry Smith   `KSPSolve()`.
11364b9ad928SBarry Smith 
1137f1580f4eSBarry Smith   A routine set by `PCSetModifySubMatrices()` is currently called within
1138f1580f4eSBarry Smith   the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`)
11394b9ad928SBarry Smith   preconditioners.  All other preconditioners ignore this routine.
11404b9ad928SBarry Smith 
1141f1580f4eSBarry Smith .seealso: `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()`
11424b9ad928SBarry Smith @*/
114304c3f3b8SBarry Smith PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx), void *ctx)
1144d71ae5a4SJacob Faibussowitsch {
11454b9ad928SBarry Smith   PetscFunctionBegin;
11460700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11474b9ad928SBarry Smith   pc->modifysubmatrices  = func;
11484b9ad928SBarry Smith   pc->modifysubmatricesP = ctx;
11493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11504b9ad928SBarry Smith }
11514b9ad928SBarry Smith 
11524b9ad928SBarry Smith /*@C
11534b9ad928SBarry Smith   PCModifySubMatrices - Calls an optional user-defined routine within
1154f1580f4eSBarry Smith   certain preconditioners if one has been set with `PCSetModifySubMatrices()`.
11554b9ad928SBarry Smith 
1156c3339decSBarry Smith   Collective
11574b9ad928SBarry Smith 
11584b9ad928SBarry Smith   Input Parameters:
11594b9ad928SBarry Smith + pc     - the preconditioner context
11604b9ad928SBarry Smith . nsub   - the number of local submatrices
11614b9ad928SBarry Smith . row    - an array of index sets that contain the global row numbers
11624b9ad928SBarry Smith          that comprise each local submatrix
11634b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
11644b9ad928SBarry Smith          that comprise each local submatrix
11654b9ad928SBarry Smith . submat - array of local submatrices
11664b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
11674b9ad928SBarry Smith          user-defined routine (may be null)
11684b9ad928SBarry Smith 
11694b9ad928SBarry Smith   Output Parameter:
11704b9ad928SBarry Smith . submat - array of local submatrices (the entries of which may
11714b9ad928SBarry Smith             have been modified)
11724b9ad928SBarry Smith 
117320f4b53cSBarry Smith   Level: developer
117420f4b53cSBarry Smith 
117504c3f3b8SBarry Smith   Note:
11764b9ad928SBarry Smith   The user should NOT generally call this routine, as it will
117704c3f3b8SBarry Smith   automatically be called within certain preconditioners.
11784b9ad928SBarry Smith 
1179f1580f4eSBarry Smith .seealso: `PC`, `PCSetModifySubMatrices()`
11804b9ad928SBarry Smith @*/
1181d71ae5a4SJacob Faibussowitsch PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx)
1182d71ae5a4SJacob Faibussowitsch {
11834b9ad928SBarry Smith   PetscFunctionBegin;
11840700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11853ba16761SJacob Faibussowitsch   if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS);
11869566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0));
11879566063dSJacob Faibussowitsch   PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx));
11889566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0));
11893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11904b9ad928SBarry Smith }
11914b9ad928SBarry Smith 
11924b9ad928SBarry Smith /*@
11934b9ad928SBarry Smith   PCSetOperators - Sets the matrix associated with the linear system and
11944b9ad928SBarry Smith   a (possibly) different one associated with the preconditioner.
11954b9ad928SBarry Smith 
1196c3339decSBarry Smith   Logically Collective
11974b9ad928SBarry Smith 
11984b9ad928SBarry Smith   Input Parameters:
11994b9ad928SBarry Smith + pc   - the preconditioner context
1200e5d3d808SBarry Smith . Amat - the matrix that defines the linear system
1201d1e9a80fSBarry Smith - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
12024b9ad928SBarry Smith 
120320f4b53cSBarry Smith   Level: intermediate
1204189c0b8aSBarry Smith 
120520f4b53cSBarry Smith   Notes:
120620f4b53cSBarry Smith   Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used.
120720f4b53cSBarry Smith 
120820f4b53cSBarry Smith   If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then
1209f1580f4eSBarry Smith   first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()`
1210f1580f4eSBarry Smith   on it and then pass it back in in your call to `KSPSetOperators()`.
1211189c0b8aSBarry Smith 
12124b9ad928SBarry Smith   More Notes about Repeated Solution of Linear Systems:
121320f4b53cSBarry Smith   PETSc does NOT reset the matrix entries of either `Amat` or `Pmat`
12144b9ad928SBarry Smith   to zero after a linear solve; the user is completely responsible for
1215f1580f4eSBarry Smith   matrix assembly.  See the routine `MatZeroEntries()` if desiring to
12164b9ad928SBarry Smith   zero all elements of a matrix.
12174b9ad928SBarry Smith 
1218f1580f4eSBarry Smith .seealso: `PC`, `PCGetOperators()`, `MatZeroEntries()`
12194b9ad928SBarry Smith  @*/
1220d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat)
1221d71ae5a4SJacob Faibussowitsch {
12223b3f6333SBarry Smith   PetscInt m1, n1, m2, n2;
12234b9ad928SBarry Smith 
12244b9ad928SBarry Smith   PetscFunctionBegin;
12250700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12260700a824SBarry Smith   if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
12270700a824SBarry Smith   if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3);
12283fc8bf9cSBarry Smith   if (Amat) PetscCheckSameComm(pc, 1, Amat, 2);
12293fc8bf9cSBarry Smith   if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3);
123031641f1aSBarry Smith   if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) {
12319566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Amat, &m1, &n1));
12329566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->mat, &m2, &n2));
123363a3b9bcSJacob 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);
12349566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Pmat, &m1, &n1));
12359566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2));
123663a3b9bcSJacob 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);
12373b3f6333SBarry Smith   }
12384b9ad928SBarry Smith 
123923ee1639SBarry Smith   if (Pmat != pc->pmat) {
124023ee1639SBarry Smith     /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
124123ee1639SBarry Smith     pc->matnonzerostate = -1;
124223ee1639SBarry Smith     pc->matstate        = -1;
124323ee1639SBarry Smith   }
124423ee1639SBarry Smith 
1245906ed7ccSBarry Smith   /* reference first in case the matrices are the same */
12469566063dSJacob Faibussowitsch   if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat));
12479566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->mat));
12489566063dSJacob Faibussowitsch   if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat));
12499566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->pmat));
12504b9ad928SBarry Smith   pc->mat  = Amat;
12514b9ad928SBarry Smith   pc->pmat = Pmat;
12523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1253e56f5c9eSBarry Smith }
1254e56f5c9eSBarry Smith 
125523ee1639SBarry Smith /*@
125623ee1639SBarry Smith   PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed.
125723ee1639SBarry Smith 
1258c3339decSBarry Smith   Logically Collective
125923ee1639SBarry Smith 
126023ee1639SBarry Smith   Input Parameters:
126123ee1639SBarry Smith + pc   - the preconditioner context
1262f1580f4eSBarry Smith - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
126323ee1639SBarry Smith 
12642b26979fSBarry Smith   Level: intermediate
12652b26979fSBarry Smith 
1266f1580f4eSBarry Smith   Note:
1267f1580f4eSBarry Smith   Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option
1268f1580f4eSBarry Smith   prevents this.
1269f1580f4eSBarry Smith 
1270f1580f4eSBarry Smith .seealso: `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()`
127123ee1639SBarry Smith  @*/
1272d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag)
1273d71ae5a4SJacob Faibussowitsch {
127423ee1639SBarry Smith   PetscFunctionBegin;
127523ee1639SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1276f9177522SStefano Zampini   PetscValidLogicalCollectiveBool(pc, flag, 2);
127723ee1639SBarry Smith   pc->reusepreconditioner = flag;
12783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12794b9ad928SBarry Smith }
12804b9ad928SBarry Smith 
1281c60c7ad4SBarry Smith /*@
1282f1580f4eSBarry Smith   PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed.
1283c60c7ad4SBarry Smith 
1284bba28a21SBarry Smith   Not Collective
1285c60c7ad4SBarry Smith 
1286c60c7ad4SBarry Smith   Input Parameter:
1287c60c7ad4SBarry Smith . pc - the preconditioner context
1288c60c7ad4SBarry Smith 
1289c60c7ad4SBarry Smith   Output Parameter:
1290f1580f4eSBarry Smith . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
1291c60c7ad4SBarry Smith 
1292d0418729SSatish Balay   Level: intermediate
1293d0418729SSatish Balay 
1294f1580f4eSBarry Smith .seealso: `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()`
1295c60c7ad4SBarry Smith  @*/
1296d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag)
1297d71ae5a4SJacob Faibussowitsch {
1298c60c7ad4SBarry Smith   PetscFunctionBegin;
1299c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
13004f572ea9SToby Isaac   PetscAssertPointer(flag, 2);
1301c60c7ad4SBarry Smith   *flag = pc->reusepreconditioner;
13023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1303c60c7ad4SBarry Smith }
1304c60c7ad4SBarry Smith 
1305487a658cSBarry Smith /*@
13064b9ad928SBarry Smith   PCGetOperators - Gets the matrix associated with the linear system and
13074b9ad928SBarry Smith   possibly a different one associated with the preconditioner.
13084b9ad928SBarry Smith 
130920f4b53cSBarry Smith   Not Collective, though parallel mats are returned if pc is parallel
13104b9ad928SBarry Smith 
13114b9ad928SBarry Smith   Input Parameter:
13124b9ad928SBarry Smith . pc - the preconditioner context
13134b9ad928SBarry Smith 
13144b9ad928SBarry Smith   Output Parameters:
1315e5d3d808SBarry Smith + Amat - the matrix defining the linear system
131623ee1639SBarry Smith - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
13174b9ad928SBarry Smith 
13184b9ad928SBarry Smith   Level: intermediate
13194b9ad928SBarry Smith 
1320f1580f4eSBarry Smith   Note:
132195452b02SPatrick Sanan   Does not increase the reference count of the matrices, so you should not destroy them
1322298cc208SBarry Smith 
1323f1580f4eSBarry Smith   Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators
1324f1580f4eSBarry Smith   are created in `PC` and returned to the user. In this case, if both operators
132573950996SBarry Smith   mat and pmat are requested, two DIFFERENT operators will be returned. If
132673950996SBarry Smith   only one is requested both operators in the PC will be the same (i.e. as
1327f1580f4eSBarry Smith   if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats).
132873950996SBarry Smith   The user must set the sizes of the returned matrices and their type etc just
1329f1580f4eSBarry Smith   as if the user created them with `MatCreate()`. For example,
133073950996SBarry Smith 
1331f1580f4eSBarry Smith .vb
1332f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1333f1580f4eSBarry Smith            set size, type, etc of Amat
133473950996SBarry Smith 
1335f1580f4eSBarry Smith          MatCreate(comm,&mat);
1336f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Amat);
1337f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)mat);
1338f1580f4eSBarry Smith            set size, type, etc of Amat
1339f1580f4eSBarry Smith .ve
134073950996SBarry Smith 
134173950996SBarry Smith   and
134273950996SBarry Smith 
1343f1580f4eSBarry Smith .vb
1344f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1345f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
134673950996SBarry Smith 
1347f1580f4eSBarry Smith          MatCreate(comm,&Amat);
1348f1580f4eSBarry Smith          MatCreate(comm,&Pmat);
1349f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1350f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Amat);
1351f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Pmat);
1352f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
1353f1580f4eSBarry Smith .ve
135473950996SBarry Smith 
1355f1580f4eSBarry Smith   The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy
1356b8d709abSRichard Tran Mills   of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely
1357f1580f4eSBarry Smith   managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look
1358f1580f4eSBarry Smith   at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to
1359f1580f4eSBarry Smith   the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP
1360f1580f4eSBarry Smith   you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you).
1361f1580f4eSBarry Smith   Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when
136273950996SBarry Smith   it can be created for you?
136373950996SBarry Smith 
1364f1580f4eSBarry Smith .seealso: `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()`
13654b9ad928SBarry Smith @*/
1366d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat)
1367d71ae5a4SJacob Faibussowitsch {
13684b9ad928SBarry Smith   PetscFunctionBegin;
13690700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1370e5d3d808SBarry Smith   if (Amat) {
137173950996SBarry Smith     if (!pc->mat) {
13729fca8c71SStefano Zampini       if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */
13739a4708feSJed Brown         pc->mat = pc->pmat;
13749566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->mat));
1375e5d3d808SBarry Smith       } else { /* both Amat and Pmat are empty */
13769566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat));
1377e5d3d808SBarry Smith         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
137873950996SBarry Smith           pc->pmat = pc->mat;
13799566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->pmat));
138073950996SBarry Smith         }
138173950996SBarry Smith       }
13829a4708feSJed Brown     }
1383e5d3d808SBarry Smith     *Amat = pc->mat;
138473950996SBarry Smith   }
1385e5d3d808SBarry Smith   if (Pmat) {
138673950996SBarry Smith     if (!pc->pmat) {
1387e5d3d808SBarry Smith       if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */
13889a4708feSJed Brown         pc->pmat = pc->mat;
13899566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->pmat));
13909a4708feSJed Brown       } else {
13919566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat));
1392e5d3d808SBarry Smith         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
139373950996SBarry Smith           pc->mat = pc->pmat;
13949566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->mat));
139573950996SBarry Smith         }
139673950996SBarry Smith       }
13979a4708feSJed Brown     }
1398e5d3d808SBarry Smith     *Pmat = pc->pmat;
139973950996SBarry Smith   }
14003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14014b9ad928SBarry Smith }
14024b9ad928SBarry Smith 
1403906ed7ccSBarry Smith /*@C
1404906ed7ccSBarry Smith   PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1405f1580f4eSBarry Smith   possibly a different one associated with the preconditioner have been set in the `PC`.
1406906ed7ccSBarry Smith 
140720f4b53cSBarry Smith   Not Collective, though the results on all processes should be the same
1408906ed7ccSBarry Smith 
1409906ed7ccSBarry Smith   Input Parameter:
1410906ed7ccSBarry Smith . pc - the preconditioner context
1411906ed7ccSBarry Smith 
1412906ed7ccSBarry Smith   Output Parameters:
1413906ed7ccSBarry Smith + mat  - the matrix associated with the linear system was set
1414906ed7ccSBarry Smith - pmat - matrix associated with the preconditioner was set, usually the same
1415906ed7ccSBarry Smith 
1416906ed7ccSBarry Smith   Level: intermediate
1417906ed7ccSBarry Smith 
1418f1580f4eSBarry Smith .seealso: `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`
1419906ed7ccSBarry Smith @*/
1420d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat)
1421d71ae5a4SJacob Faibussowitsch {
1422906ed7ccSBarry Smith   PetscFunctionBegin;
14230700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1424906ed7ccSBarry Smith   if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1425906ed7ccSBarry Smith   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
14263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1427906ed7ccSBarry Smith }
1428906ed7ccSBarry Smith 
1429f39d8e23SSatish Balay /*@
1430a4fd02acSBarry Smith   PCFactorGetMatrix - Gets the factored matrix from the
1431f1580f4eSBarry Smith   preconditioner context.  This routine is valid only for the `PCLU`,
1432f1580f4eSBarry Smith   `PCILU`, `PCCHOLESKY`, and `PCICC` methods.
14334b9ad928SBarry Smith 
1434c3339decSBarry Smith   Not Collective though mat is parallel if pc is parallel
14354b9ad928SBarry Smith 
1436f1580f4eSBarry Smith   Input Parameter:
14374b9ad928SBarry Smith . pc - the preconditioner context
14384b9ad928SBarry Smith 
1439feefa0e1SJacob Faibussowitsch   Output Parameters:
14404b9ad928SBarry Smith . mat - the factored matrix
14414b9ad928SBarry Smith 
14424b9ad928SBarry Smith   Level: advanced
14434b9ad928SBarry Smith 
1444f1580f4eSBarry Smith   Note:
144595452b02SPatrick Sanan   Does not increase the reference count for the matrix so DO NOT destroy it
14469405f653SBarry Smith 
1447f1580f4eSBarry Smith .seealso: `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC`
14484b9ad928SBarry Smith @*/
1449d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat)
1450d71ae5a4SJacob Faibussowitsch {
14514b9ad928SBarry Smith   PetscFunctionBegin;
14520700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14534f572ea9SToby Isaac   PetscAssertPointer(mat, 2);
14547def7855SStefano Zampini   PetscCall(PCFactorSetUpMatSolverType(pc));
1455dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, getfactoredmatrix, mat);
14563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14574b9ad928SBarry Smith }
14584b9ad928SBarry Smith 
14594b9ad928SBarry Smith /*@C
14604b9ad928SBarry Smith   PCSetOptionsPrefix - Sets the prefix used for searching for all
1461f1580f4eSBarry Smith   `PC` options in the database.
14624b9ad928SBarry Smith 
1463c3339decSBarry Smith   Logically Collective
14644b9ad928SBarry Smith 
14654b9ad928SBarry Smith   Input Parameters:
14664b9ad928SBarry Smith + pc     - the preconditioner context
1467f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
14684b9ad928SBarry Smith 
1469f1580f4eSBarry Smith   Note:
14704b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
14714b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
14724b9ad928SBarry Smith   hyphen.
14734b9ad928SBarry Smith 
14744b9ad928SBarry Smith   Level: advanced
14754b9ad928SBarry Smith 
1476f1580f4eSBarry Smith .seealso: `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()`
14774b9ad928SBarry Smith @*/
1478d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[])
1479d71ae5a4SJacob Faibussowitsch {
14804b9ad928SBarry Smith   PetscFunctionBegin;
14810700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14829566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix));
14833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14844b9ad928SBarry Smith }
14854b9ad928SBarry Smith 
14864b9ad928SBarry Smith /*@C
14874b9ad928SBarry Smith   PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1488f1580f4eSBarry Smith   `PC` options in the database.
14894b9ad928SBarry Smith 
1490c3339decSBarry Smith   Logically Collective
14914b9ad928SBarry Smith 
14924b9ad928SBarry Smith   Input Parameters:
14934b9ad928SBarry Smith + pc     - the preconditioner context
1494f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
14954b9ad928SBarry Smith 
1496f1580f4eSBarry Smith   Note:
14974b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
14984b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
14994b9ad928SBarry Smith   hyphen.
15004b9ad928SBarry Smith 
15014b9ad928SBarry Smith   Level: advanced
15024b9ad928SBarry Smith 
1503f1580f4eSBarry Smith .seealso: `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()`
15044b9ad928SBarry Smith @*/
1505d71ae5a4SJacob Faibussowitsch PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[])
1506d71ae5a4SJacob Faibussowitsch {
15074b9ad928SBarry Smith   PetscFunctionBegin;
15080700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15099566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix));
15103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15114b9ad928SBarry Smith }
15124b9ad928SBarry Smith 
15134b9ad928SBarry Smith /*@C
15144b9ad928SBarry Smith   PCGetOptionsPrefix - Gets the prefix used for searching for all
15154b9ad928SBarry Smith   PC options in the database.
15164b9ad928SBarry Smith 
15174b9ad928SBarry Smith   Not Collective
15184b9ad928SBarry Smith 
1519f1580f4eSBarry Smith   Input Parameter:
15204b9ad928SBarry Smith . pc - the preconditioner context
15214b9ad928SBarry Smith 
1522f1580f4eSBarry Smith   Output Parameter:
15234b9ad928SBarry Smith . prefix - pointer to the prefix string used, is returned
15244b9ad928SBarry Smith 
1525feefa0e1SJacob Faibussowitsch   Fortran Notes:
1526f1580f4eSBarry Smith   The user should pass in a string 'prefix' of
15274b9ad928SBarry Smith   sufficient length to hold the prefix.
15284b9ad928SBarry Smith 
15294b9ad928SBarry Smith   Level: advanced
15304b9ad928SBarry Smith 
1531f1580f4eSBarry Smith .seealso: `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()`
15324b9ad928SBarry Smith @*/
1533d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[])
1534d71ae5a4SJacob Faibussowitsch {
15354b9ad928SBarry Smith   PetscFunctionBegin;
15360700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15374f572ea9SToby Isaac   PetscAssertPointer(prefix, 2);
15389566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix));
15393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15404b9ad928SBarry Smith }
15414b9ad928SBarry Smith 
15428066bbecSBarry Smith /*
15438066bbecSBarry Smith    Indicates the right hand side will be changed by KSPSolve(), this occurs for a few
15448066bbecSBarry Smith   preconditioners including BDDC and Eisentat that transform the equations before applying
15458066bbecSBarry Smith   the Krylov methods
15468066bbecSBarry Smith */
1547d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change)
1548d71ae5a4SJacob Faibussowitsch {
1549a06fd7f2SStefano Zampini   PetscFunctionBegin;
1550a06fd7f2SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15514f572ea9SToby Isaac   PetscAssertPointer(change, 2);
1552a06fd7f2SStefano Zampini   *change = PETSC_FALSE;
1553cac4c232SBarry Smith   PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change));
15543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1555a06fd7f2SStefano Zampini }
1556a06fd7f2SStefano Zampini 
15574b9ad928SBarry Smith /*@
15584b9ad928SBarry Smith   PCPreSolve - Optional pre-solve phase, intended for any
15594b9ad928SBarry Smith   preconditioner-specific actions that must be performed before
15604b9ad928SBarry Smith   the iterative solve itself.
15614b9ad928SBarry Smith 
1562c3339decSBarry Smith   Collective
15634b9ad928SBarry Smith 
15644b9ad928SBarry Smith   Input Parameters:
15654b9ad928SBarry Smith + pc  - the preconditioner context
15664b9ad928SBarry Smith - ksp - the Krylov subspace context
15674b9ad928SBarry Smith 
15684b9ad928SBarry Smith   Level: developer
15694b9ad928SBarry Smith 
1570feefa0e1SJacob Faibussowitsch   Example Usage:
15714b9ad928SBarry Smith .vb
15724b9ad928SBarry Smith     PCPreSolve(pc,ksp);
157323ce1328SBarry Smith     KSPSolve(ksp,b,x);
15744b9ad928SBarry Smith     PCPostSolve(pc,ksp);
15754b9ad928SBarry Smith .ve
15764b9ad928SBarry Smith 
15774b9ad928SBarry Smith   Notes:
1578f1580f4eSBarry Smith   The pre-solve phase is distinct from the `PCSetUp()` phase.
15794b9ad928SBarry Smith 
1580f1580f4eSBarry Smith   `KSPSolve()` calls this directly, so is rarely called by the user.
15814b9ad928SBarry Smith 
1582f1580f4eSBarry Smith .seealso: `PC`, `PCPostSolve()`
15834b9ad928SBarry Smith @*/
1584d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPreSolve(PC pc, KSP ksp)
1585d71ae5a4SJacob Faibussowitsch {
15864b9ad928SBarry Smith   Vec x, rhs;
15874b9ad928SBarry Smith 
15884b9ad928SBarry Smith   PetscFunctionBegin;
15890700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15900700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1591d9a03883SBarry Smith   pc->presolvedone++;
15927827d75bSBarry Smith   PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice");
15939566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
15949566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
15954b9ad928SBarry Smith 
1596dbbe0bcdSBarry Smith   if (pc->ops->presolve) PetscUseTypeMethod(pc, presolve, ksp, rhs, x);
15971baa6e33SBarry Smith   else if (pc->presolve) PetscCall((pc->presolve)(pc, ksp));
15983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15994b9ad928SBarry Smith }
16004b9ad928SBarry Smith 
1601f560b561SHong Zhang /*@C
1602f1580f4eSBarry Smith   PCSetPreSolve - Sets function used by `PCPreSolve()` which is intended for any
1603f560b561SHong Zhang   preconditioner-specific actions that must be performed before
1604f560b561SHong Zhang   the iterative solve itself.
1605f560b561SHong Zhang 
1606c3339decSBarry Smith   Logically Collective
1607f560b561SHong Zhang 
1608f560b561SHong Zhang   Input Parameters:
1609f560b561SHong Zhang + pc       - the preconditioner object
1610f560b561SHong Zhang - presolve - the function to call before the solve
1611f560b561SHong Zhang 
161220f4b53cSBarry Smith   Calling sequence of `presolve`:
161320f4b53cSBarry Smith + pc  - the `PC` context
161420f4b53cSBarry Smith - ksp - the `KSP` context
1615f560b561SHong Zhang 
1616f560b561SHong Zhang   Level: developer
1617f560b561SHong Zhang 
1618db781477SPatrick Sanan .seealso: `PC`, `PCSetUp()`, `PCPreSolve()`
1619f560b561SHong Zhang @*/
162004c3f3b8SBarry Smith PetscErrorCode PCSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC pc, KSP ksp))
1621d71ae5a4SJacob Faibussowitsch {
1622f560b561SHong Zhang   PetscFunctionBegin;
1623f560b561SHong Zhang   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1624f560b561SHong Zhang   pc->presolve = presolve;
16253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1626f560b561SHong Zhang }
1627f560b561SHong Zhang 
16284b9ad928SBarry Smith /*@
16294b9ad928SBarry Smith   PCPostSolve - Optional post-solve phase, intended for any
16304b9ad928SBarry Smith   preconditioner-specific actions that must be performed after
16314b9ad928SBarry Smith   the iterative solve itself.
16324b9ad928SBarry Smith 
1633c3339decSBarry Smith   Collective
16344b9ad928SBarry Smith 
16354b9ad928SBarry Smith   Input Parameters:
16364b9ad928SBarry Smith + pc  - the preconditioner context
16374b9ad928SBarry Smith - ksp - the Krylov subspace context
16384b9ad928SBarry Smith 
1639feefa0e1SJacob Faibussowitsch   Example Usage:
16404b9ad928SBarry Smith .vb
16414b9ad928SBarry Smith     PCPreSolve(pc,ksp);
164223ce1328SBarry Smith     KSPSolve(ksp,b,x);
16434b9ad928SBarry Smith     PCPostSolve(pc,ksp);
16444b9ad928SBarry Smith .ve
16454b9ad928SBarry Smith 
16464b9ad928SBarry Smith   Note:
1647f1580f4eSBarry Smith   `KSPSolve()` calls this routine directly, so it is rarely called by the user.
16484b9ad928SBarry Smith 
16494b9ad928SBarry Smith   Level: developer
16504b9ad928SBarry Smith 
1651f1580f4eSBarry Smith .seealso: `PC`, `PCSetPostSolve()`, `PCSetPresolve()`, `PCPreSolve()`, `KSPSolve()`
16524b9ad928SBarry Smith @*/
1653d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPostSolve(PC pc, KSP ksp)
1654d71ae5a4SJacob Faibussowitsch {
16554b9ad928SBarry Smith   Vec x, rhs;
16564b9ad928SBarry Smith 
16574b9ad928SBarry Smith   PetscFunctionBegin;
16580700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16590700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1660d9a03883SBarry Smith   pc->presolvedone--;
16619566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
16629566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
1663dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, postsolve, ksp, rhs, x);
16643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16654b9ad928SBarry Smith }
16664b9ad928SBarry Smith 
166755849f57SBarry Smith /*@C
1668f1580f4eSBarry Smith   PCLoad - Loads a `PC` that has been stored in binary  with `PCView()`.
166955849f57SBarry Smith 
1670c3339decSBarry Smith   Collective
167155849f57SBarry Smith 
167255849f57SBarry Smith   Input Parameters:
1673f1580f4eSBarry Smith + newdm  - the newly loaded `PC`, this needs to have been created with `PCCreate()` or
1674f1580f4eSBarry Smith            some related function before a call to `PCLoad()`.
1675f1580f4eSBarry Smith - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()`
167655849f57SBarry Smith 
167755849f57SBarry Smith   Level: intermediate
167855849f57SBarry Smith 
1679f1580f4eSBarry Smith   Note:
168055849f57SBarry Smith   The type is determined by the data in the file, any type set into the PC before this call is ignored.
168155849f57SBarry Smith 
1682f1580f4eSBarry Smith .seealso: `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()`
168355849f57SBarry Smith @*/
1684d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLoad(PC newdm, PetscViewer viewer)
1685d71ae5a4SJacob Faibussowitsch {
168655849f57SBarry Smith   PetscBool isbinary;
1687060da220SMatthew G. Knepley   PetscInt  classid;
168855849f57SBarry Smith   char      type[256];
168955849f57SBarry Smith 
169055849f57SBarry Smith   PetscFunctionBegin;
169155849f57SBarry Smith   PetscValidHeaderSpecific(newdm, PC_CLASSID, 1);
169255849f57SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
16939566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
169428b400f6SJacob Faibussowitsch   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
169555849f57SBarry Smith 
16969566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
16977827d75bSBarry Smith   PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file");
16989566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
16999566063dSJacob Faibussowitsch   PetscCall(PCSetType(newdm, type));
1700dbbe0bcdSBarry Smith   PetscTryTypeMethod(newdm, load, viewer);
17013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170255849f57SBarry Smith }
170355849f57SBarry Smith 
17049804daf3SBarry Smith #include <petscdraw.h>
1705e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1706e04113cfSBarry Smith   #include <petscviewersaws.h>
17070acecf5bSBarry Smith #endif
1708fe2efc57SMark 
1709fe2efc57SMark /*@C
1710f1580f4eSBarry Smith   PCViewFromOptions - View from the `PC` based on options in the database
1711fe2efc57SMark 
1712c3339decSBarry Smith   Collective
1713fe2efc57SMark 
1714fe2efc57SMark   Input Parameters:
171520f4b53cSBarry Smith + A    - the `PC` context
1716f1580f4eSBarry Smith . obj  - Optional object that provides the options prefix
1717736c3998SJose E. Roman - name - command line option
1718fe2efc57SMark 
1719fe2efc57SMark   Level: intermediate
1720f1580f4eSBarry Smith 
1721db781477SPatrick Sanan .seealso: `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()`
1722fe2efc57SMark @*/
1723d71ae5a4SJacob Faibussowitsch PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[])
1724d71ae5a4SJacob Faibussowitsch {
1725fe2efc57SMark   PetscFunctionBegin;
1726fe2efc57SMark   PetscValidHeaderSpecific(A, PC_CLASSID, 1);
17279566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
17283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1729fe2efc57SMark }
1730fe2efc57SMark 
17314b9ad928SBarry Smith /*@C
1732f1580f4eSBarry Smith   PCView - Prints information about the `PC`
17334b9ad928SBarry Smith 
1734c3339decSBarry Smith   Collective
17354b9ad928SBarry Smith 
17364b9ad928SBarry Smith   Input Parameters:
1737feefa0e1SJacob Faibussowitsch + pc     - the `PC` context
17384b9ad928SBarry Smith - viewer - optional visualization context
17394b9ad928SBarry Smith 
174020f4b53cSBarry Smith   Level: developer
174120f4b53cSBarry Smith 
1742f1580f4eSBarry Smith   Notes:
17434b9ad928SBarry Smith   The available visualization contexts include
1744f1580f4eSBarry Smith +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
1745f1580f4eSBarry Smith -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
17464b9ad928SBarry Smith   output where only the first processor opens
17474b9ad928SBarry Smith   the file. All other processors send their
17484b9ad928SBarry Smith   data to the first processor to print.
17494b9ad928SBarry Smith 
17504b9ad928SBarry Smith   The user can open an alternative visualization contexts with
1751f1580f4eSBarry Smith   `PetscViewerASCIIOpen()` (output to a specified file).
17524b9ad928SBarry Smith 
1753f1580f4eSBarry Smith .seealso: `PC`, `PetscViewer`, `KSPView()`, `PetscViewerASCIIOpen()`
17544b9ad928SBarry Smith @*/
1755d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView(PC pc, PetscViewer viewer)
1756d71ae5a4SJacob Faibussowitsch {
175719fd82e9SBarry Smith   PCType            cstr;
17586cd81132SPierre Jolivet   PetscViewerFormat format;
17596cd81132SPierre Jolivet   PetscBool         iascii, isstring, isbinary, isdraw, pop = PETSC_FALSE;
1760e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1761536b137fSBarry Smith   PetscBool issaws;
17620acecf5bSBarry Smith #endif
17634b9ad928SBarry Smith 
17644b9ad928SBarry Smith   PetscFunctionBegin;
17650700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
176648a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
17670700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1768c9780b6fSBarry Smith   PetscCheckSameComm(pc, 1, viewer, 2);
17694b9ad928SBarry Smith 
17709566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
17719566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
17729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
17739566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1774e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
17759566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
17760acecf5bSBarry Smith #endif
1777219b1045SBarry Smith 
177832077d6dSBarry Smith   if (iascii) {
17799566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer));
178048a46eb9SPierre Jolivet     if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, "  PC has not been set up so information may be incomplete\n"));
17819566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1782dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
17839566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
1784834dbeb0SBarry Smith     if (pc->mat) {
17856cd81132SPierre Jolivet       PetscCall(PetscViewerGetFormat(viewer, &format));
17866cd81132SPierre Jolivet       if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
17879566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
17886cd81132SPierre Jolivet         pop = PETSC_TRUE;
17896cd81132SPierre Jolivet       }
17904b9ad928SBarry Smith       if (pc->pmat == pc->mat) {
17919566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix = precond matrix:\n"));
17929566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
17939566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
17949566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
17954b9ad928SBarry Smith       } else {
1796834dbeb0SBarry Smith         if (pc->pmat) {
17979566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix followed by preconditioner matrix:\n"));
17984b9ad928SBarry Smith         } else {
17999566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix:\n"));
18004b9ad928SBarry Smith         }
18019566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
18029566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
18039566063dSJacob Faibussowitsch         if (pc->pmat) PetscCall(MatView(pc->pmat, viewer));
18049566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
18054b9ad928SBarry Smith       }
18066cd81132SPierre Jolivet       if (pop) PetscCall(PetscViewerPopFormat(viewer));
18074b9ad928SBarry Smith     }
18084b9ad928SBarry Smith   } else if (isstring) {
18099566063dSJacob Faibussowitsch     PetscCall(PCGetType(pc, &cstr));
18109566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr));
1811dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18129566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18139566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18145b0b0462SBarry Smith   } else if (isbinary) {
181555849f57SBarry Smith     PetscInt    classid = PC_FILE_CLASSID;
181655849f57SBarry Smith     MPI_Comm    comm;
181755849f57SBarry Smith     PetscMPIInt rank;
181855849f57SBarry Smith     char        type[256];
181955849f57SBarry Smith 
18209566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
18219566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1822dd400576SPatrick Sanan     if (rank == 0) {
18239566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
18249566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256));
18259566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
182655849f57SBarry Smith     }
1827dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
1828219b1045SBarry Smith   } else if (isdraw) {
1829219b1045SBarry Smith     PetscDraw draw;
1830d9884438SBarry Smith     char      str[25];
183189fd9fafSBarry Smith     PetscReal x, y, bottom, h;
1832d9884438SBarry Smith     PetscInt  n;
1833219b1045SBarry Smith 
18349566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
18359566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
18361d840656SPeter Brune     if (pc->mat) {
18379566063dSJacob Faibussowitsch       PetscCall(MatGetSize(pc->mat, &n, NULL));
183863a3b9bcSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n));
18391d840656SPeter Brune     } else {
18409566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name));
18411d840656SPeter Brune     }
18429566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
184389fd9fafSBarry Smith     bottom = y - h;
18449566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
1845dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18469566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
1847e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1848536b137fSBarry Smith   } else if (issaws) {
1849d45a07a7SBarry Smith     PetscMPIInt rank;
1850d45a07a7SBarry Smith 
18519566063dSJacob Faibussowitsch     PetscCall(PetscObjectName((PetscObject)pc));
18529566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
185348a46eb9SPierre Jolivet     if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer));
18549566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18559566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18560acecf5bSBarry Smith #endif
18574b9ad928SBarry Smith   }
18583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18594b9ad928SBarry Smith }
18604b9ad928SBarry Smith 
18614b9ad928SBarry Smith /*@C
18621c84c290SBarry Smith   PCRegister -  Adds a method to the preconditioner package.
18631c84c290SBarry Smith 
18641c84c290SBarry Smith   Not collective
18651c84c290SBarry Smith 
18661c84c290SBarry Smith   Input Parameters:
186720f4b53cSBarry Smith + sname    - name of a new user-defined solver
186820f4b53cSBarry Smith - function - routine to create method context
18691c84c290SBarry Smith 
1870feefa0e1SJacob Faibussowitsch   Example Usage:
18711c84c290SBarry Smith .vb
1872bdf89e91SBarry Smith    PCRegister("my_solver", MySolverCreate);
18731c84c290SBarry Smith .ve
18741c84c290SBarry Smith 
18751c84c290SBarry Smith   Then, your solver can be chosen with the procedural interface via
18761c84c290SBarry Smith $     PCSetType(pc, "my_solver")
18771c84c290SBarry Smith   or at runtime via the option
18781c84c290SBarry Smith $     -pc_type my_solver
18794b9ad928SBarry Smith 
18804b9ad928SBarry Smith   Level: advanced
18811c84c290SBarry Smith 
188220f4b53cSBarry Smith   Note:
188320f4b53cSBarry Smith   `PCRegister()` may be called multiple times to add several user-defined preconditioners.
188420f4b53cSBarry Smith 
1885f1580f4eSBarry Smith .seealso: `PC`, `PCType`, `PCRegisterAll()`
18864b9ad928SBarry Smith @*/
1887d71ae5a4SJacob Faibussowitsch PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC))
1888d71ae5a4SJacob Faibussowitsch {
18894b9ad928SBarry Smith   PetscFunctionBegin;
18909566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
18919566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCList, sname, function));
18923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18934b9ad928SBarry Smith }
18944b9ad928SBarry Smith 
1895d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y)
1896d71ae5a4SJacob Faibussowitsch {
1897186a3e20SStefano Zampini   PC pc;
1898186a3e20SStefano Zampini 
1899186a3e20SStefano Zampini   PetscFunctionBegin;
19009566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &pc));
19019566063dSJacob Faibussowitsch   PetscCall(PCApply(pc, X, Y));
19023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1903186a3e20SStefano Zampini }
1904186a3e20SStefano Zampini 
1905feefa0e1SJacob Faibussowitsch /*@C
19060bacdadaSStefano Zampini   PCComputeOperator - Computes the explicit preconditioned operator.
19074b9ad928SBarry Smith 
1908c3339decSBarry Smith   Collective
19094b9ad928SBarry Smith 
1910d8d19677SJose E. Roman   Input Parameters:
1911186a3e20SStefano Zampini + pc      - the preconditioner object
1912186a3e20SStefano Zampini - mattype - the matrix type to be used for the operator
19134b9ad928SBarry Smith 
19144b9ad928SBarry Smith   Output Parameter:
1915a5b23f4aSJose E. Roman . mat - the explicit preconditioned operator
19164b9ad928SBarry Smith 
191720f4b53cSBarry Smith   Level: advanced
191820f4b53cSBarry Smith 
1919f1580f4eSBarry Smith   Note:
1920186a3e20SStefano Zampini   This computation is done by applying the operators to columns of the identity matrix.
1921186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
1922186a3e20SStefano Zampini   Currently, this routine uses a dense matrix format when mattype == NULL
19234b9ad928SBarry Smith 
1924f1580f4eSBarry Smith .seealso: `PC`, `KSPComputeOperator()`, `MatType`
19254b9ad928SBarry Smith @*/
1926d71ae5a4SJacob Faibussowitsch PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat)
1927d71ae5a4SJacob Faibussowitsch {
1928186a3e20SStefano Zampini   PetscInt N, M, m, n;
1929186a3e20SStefano Zampini   Mat      A, Apc;
19304b9ad928SBarry Smith 
19314b9ad928SBarry Smith   PetscFunctionBegin;
19320700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
19334f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
19349566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, &A, NULL));
19359566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
19369566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
19379566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc));
19389566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (void (*)(void))MatMult_PC));
19399566063dSJacob Faibussowitsch   PetscCall(MatComputeOperator(Apc, mattype, mat));
19409566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Apc));
19413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19424b9ad928SBarry Smith }
19434b9ad928SBarry Smith 
19446c237d78SBarry Smith /*@
19456c237d78SBarry Smith   PCSetCoordinates - sets the coordinates of all the nodes on the local process
19466c237d78SBarry Smith 
1947c3339decSBarry Smith   Collective
19486c237d78SBarry Smith 
19496c237d78SBarry Smith   Input Parameters:
19506c237d78SBarry Smith + pc     - the solver context
19516c237d78SBarry Smith . dim    - the dimension of the coordinates 1, 2, or 3
195214893cbeSStefano Zampini . nloc   - the blocked size of the coordinates array
195314893cbeSStefano Zampini - coords - the coordinates array
19546c237d78SBarry Smith 
19556c237d78SBarry Smith   Level: intermediate
19566c237d78SBarry Smith 
1957f1580f4eSBarry Smith   Note:
195820f4b53cSBarry Smith   `coords` is an array of the dim coordinates for the nodes on
195920f4b53cSBarry Smith   the local processor, of size `dim`*`nloc`.
196014893cbeSStefano Zampini   If there are 108 equation on a processor
19616c237d78SBarry Smith   for a displacement finite element discretization of elasticity (so
196214893cbeSStefano Zampini   that there are nloc = 36 = 108/3 nodes) then the array must have 108
19636c237d78SBarry Smith   double precision values (ie, 3 * 36).  These x y z coordinates
19646c237d78SBarry Smith   should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
19656c237d78SBarry Smith   ... , N-1.z ].
19666c237d78SBarry Smith 
1967f1580f4eSBarry Smith .seealso: `PC`, `MatSetNearNullSpace()`
19686c237d78SBarry Smith @*/
1969d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[])
1970d71ae5a4SJacob Faibussowitsch {
19716c237d78SBarry Smith   PetscFunctionBegin;
197232954da3SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
197332954da3SStefano Zampini   PetscValidLogicalCollectiveInt(pc, dim, 2);
197422794d57SStefano Zampini   PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords));
19753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19766c237d78SBarry Smith }
1977fd2dd295SFande Kong 
1978fd2dd295SFande Kong /*@
1979fd2dd295SFande Kong   PCGetInterpolations - Gets interpolation matrices for all levels (except level 0)
1980fd2dd295SFande Kong 
1981c3339decSBarry Smith   Logically Collective
1982fd2dd295SFande Kong 
1983d8d19677SJose E. Roman   Input Parameter:
1984d8d19677SJose E. Roman . pc - the precondition context
1985fd2dd295SFande Kong 
1986d8d19677SJose E. Roman   Output Parameters:
1987d8d19677SJose E. Roman + num_levels     - the number of levels
1988d8d19677SJose E. Roman - interpolations - the interpolation matrices (size of num_levels-1)
1989fd2dd295SFande Kong 
1990fd2dd295SFande Kong   Level: advanced
1991fd2dd295SFande Kong 
1992feefa0e1SJacob Faibussowitsch   Developer Notes:
1993f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
1994fd2dd295SFande Kong 
1995f1580f4eSBarry Smith .seealso: `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()`
1996fd2dd295SFande Kong @*/
1997d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[])
1998d71ae5a4SJacob Faibussowitsch {
1999fd2dd295SFande Kong   PetscFunctionBegin;
2000fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20014f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20024f572ea9SToby Isaac   PetscAssertPointer(interpolations, 3);
2003cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations));
20043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2005fd2dd295SFande Kong }
2006fd2dd295SFande Kong 
2007fd2dd295SFande Kong /*@
2008fd2dd295SFande Kong   PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level)
2009fd2dd295SFande Kong 
2010c3339decSBarry Smith   Logically Collective
2011fd2dd295SFande Kong 
2012d8d19677SJose E. Roman   Input Parameter:
2013d8d19677SJose E. Roman . pc - the precondition context
2014fd2dd295SFande Kong 
2015d8d19677SJose E. Roman   Output Parameters:
2016d8d19677SJose E. Roman + num_levels      - the number of levels
2017d8d19677SJose E. Roman - coarseOperators - the coarse operator matrices (size of num_levels-1)
2018fd2dd295SFande Kong 
2019fd2dd295SFande Kong   Level: advanced
2020fd2dd295SFande Kong 
2021feefa0e1SJacob Faibussowitsch   Developer Notes:
2022f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
2023fd2dd295SFande Kong 
2024f1580f4eSBarry Smith .seealso: `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()`
2025fd2dd295SFande Kong @*/
2026d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
2027d71ae5a4SJacob Faibussowitsch {
2028fd2dd295SFande Kong   PetscFunctionBegin;
2029fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20304f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20314f572ea9SToby Isaac   PetscAssertPointer(coarseOperators, 3);
2032cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators));
20333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2034fd2dd295SFande Kong }
2035