xref: /petsc/src/ksp/pc/interface/precon.c (revision 6cd81132b6b98a06238cc51bfa4557dba4eb398f)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     The PC (preconditioner) interface routines, callable by users.
44b9ad928SBarry Smith */
5af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscksp.h" I*/
61e25c274SJed Brown #include <petscdm.h>
74b9ad928SBarry Smith 
84b9ad928SBarry Smith /* Logging support */
90700a824SBarry Smith PetscClassId  PC_CLASSID;
10c677e75fSPierre Jolivet PetscLogEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_MatApply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
11011ea8aeSBarry Smith PetscLogEvent PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks;
12ab83eea4SMatthew G. Knepley PetscInt      PetscMGLevelId;
130316ec64SBarry Smith PetscLogStage PCMPIStage;
144b9ad928SBarry Smith 
15d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetDefaultType_Private(PC pc, const char *type[])
16d71ae5a4SJacob Faibussowitsch {
172e70eadcSBarry Smith   PetscMPIInt size;
18da74c943SJose E. Roman   PetscBool   hasop, flg1, flg2, set, flg3, isnormal;
19566e8bf2SBarry Smith 
20566e8bf2SBarry Smith   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
22566e8bf2SBarry Smith   if (pc->pmat) {
239566063dSJacob Faibussowitsch     PetscCall(MatHasOperation(pc->pmat, MATOP_GET_DIAGONAL_BLOCK, &hasop));
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;
35976e8c5aSLisandro Dalcin       } else if (hasop) { /* likely is a parallel matrix run on one processor */
36566e8bf2SBarry Smith         *type = PCBJACOBI;
37566e8bf2SBarry Smith       } else {
38566e8bf2SBarry Smith         *type = PCNONE;
39566e8bf2SBarry Smith       }
40566e8bf2SBarry Smith     } else {
41976e8c5aSLisandro Dalcin       if (hasop) {
42566e8bf2SBarry Smith         *type = PCBJACOBI;
43566e8bf2SBarry Smith       } else {
44566e8bf2SBarry Smith         *type = PCNONE;
45566e8bf2SBarry Smith       }
46566e8bf2SBarry Smith     }
47566e8bf2SBarry Smith   } else {
48566e8bf2SBarry Smith     if (size == 1) {
49566e8bf2SBarry Smith       *type = PCILU;
50566e8bf2SBarry Smith     } else {
51566e8bf2SBarry Smith       *type = PCBJACOBI;
52566e8bf2SBarry Smith     }
53566e8bf2SBarry Smith   }
543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
55566e8bf2SBarry Smith }
564b9ad928SBarry Smith 
5788584be7SBarry Smith /*@
58f1580f4eSBarry Smith   PCReset - Resets a PC context to the pcsetupcalled = 0 state and removes any allocated `Vec`s and `Mat`s
5988584be7SBarry Smith 
60c3339decSBarry Smith   Collective
6188584be7SBarry Smith 
6288584be7SBarry Smith   Input Parameter:
6388584be7SBarry Smith . pc - the preconditioner context
6488584be7SBarry Smith 
6588584be7SBarry Smith   Level: developer
6688584be7SBarry Smith 
67f1580f4eSBarry Smith   Note:
68f1580f4eSBarry 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
6988584be7SBarry Smith 
70f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCSetUp()`
7188584be7SBarry Smith @*/
72d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset(PC pc)
73d71ae5a4SJacob Faibussowitsch {
7488584be7SBarry Smith   PetscFunctionBegin;
7588584be7SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
76dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, reset);
779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleright));
789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleleft));
799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->pmat));
809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->mat));
812fa5cd67SKarl Rupp 
820ce6c5a2SBarry Smith   pc->setupcalled = 0;
833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8488584be7SBarry Smith }
8588584be7SBarry Smith 
861fb7b255SJunchao Zhang /*@C
87f1580f4eSBarry Smith   PCDestroy - Destroys `PC` context that was created with `PCCreate()`.
884b9ad928SBarry Smith 
89c3339decSBarry Smith   Collective
904b9ad928SBarry Smith 
914b9ad928SBarry Smith   Input Parameter:
924b9ad928SBarry Smith . pc - the preconditioner context
934b9ad928SBarry Smith 
944b9ad928SBarry Smith   Level: developer
954b9ad928SBarry Smith 
96f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCSetUp()`
974b9ad928SBarry Smith @*/
98d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy(PC *pc)
99d71ae5a4SJacob Faibussowitsch {
1004b9ad928SBarry Smith   PetscFunctionBegin;
1013ba16761SJacob Faibussowitsch   if (!*pc) PetscFunctionReturn(PETSC_SUCCESS);
1026bf464f9SBarry Smith   PetscValidHeaderSpecific((*pc), PC_CLASSID, 1);
1039371c9d4SSatish Balay   if (--((PetscObject)(*pc))->refct > 0) {
1049371c9d4SSatish Balay     *pc = NULL;
1053ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1069371c9d4SSatish Balay   }
1074b9ad928SBarry Smith 
1089566063dSJacob Faibussowitsch   PetscCall(PCReset(*pc));
109241cb3abSLisandro Dalcin 
110e04113cfSBarry Smith   /* if memory was published with SAWs then destroy it */
1119566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsViewOff((PetscObject)*pc));
112dbbe0bcdSBarry Smith   PetscTryTypeMethod((*pc), destroy);
1139566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&(*pc)->dm));
1149566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(pc));
1153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1164b9ad928SBarry Smith }
1174b9ad928SBarry Smith 
1184b9ad928SBarry Smith /*@C
119c5eb9154SBarry Smith   PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right
1204b9ad928SBarry Smith   scaling as needed by certain time-stepping codes.
1214b9ad928SBarry Smith 
122c3339decSBarry Smith   Logically Collective
1234b9ad928SBarry Smith 
1244b9ad928SBarry Smith   Input Parameter:
1254b9ad928SBarry Smith . pc - the preconditioner context
1264b9ad928SBarry Smith 
1274b9ad928SBarry Smith   Output Parameter:
128f1580f4eSBarry Smith . flag - `PETSC_TRUE` if it applies the scaling
1294b9ad928SBarry Smith 
1304b9ad928SBarry Smith   Level: developer
1314b9ad928SBarry Smith 
132f1580f4eSBarry Smith   Note:
133f1580f4eSBarry Smith   If this returns `PETSC_TRUE` then the system solved via the Krylov method is
134f1580f4eSBarry Smith .vb
135f1580f4eSBarry Smith       D M A D^{-1} y = D M b  for left preconditioning or
136f1580f4eSBarry Smith       D A M D^{-1} z = D b for right preconditioning
137f1580f4eSBarry Smith .ve
1384b9ad928SBarry Smith 
139f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCSetDiagonalScale()`
1404b9ad928SBarry Smith @*/
141d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetDiagonalScale(PC pc, PetscBool *flag)
142d71ae5a4SJacob Faibussowitsch {
1434b9ad928SBarry Smith   PetscFunctionBegin;
1440700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1454f572ea9SToby Isaac   PetscAssertPointer(flag, 2);
1464b9ad928SBarry Smith   *flag = pc->diagonalscale;
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1484b9ad928SBarry Smith }
1494b9ad928SBarry Smith 
1504b9ad928SBarry Smith /*@
151c5eb9154SBarry Smith   PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right
1524b9ad928SBarry Smith   scaling as needed by certain time-stepping codes.
1534b9ad928SBarry Smith 
154c3339decSBarry Smith   Logically Collective
1554b9ad928SBarry Smith 
1564b9ad928SBarry Smith   Input Parameters:
1574b9ad928SBarry Smith + pc - the preconditioner context
1584b9ad928SBarry Smith - s  - scaling vector
1594b9ad928SBarry Smith 
1604b9ad928SBarry Smith   Level: intermediate
1614b9ad928SBarry Smith 
16295452b02SPatrick Sanan   Notes:
16395452b02SPatrick Sanan   The system solved via the Krylov method is
164f1580f4eSBarry Smith .vb
165f1580f4eSBarry Smith            D M A D^{-1} y = D M b  for left preconditioning or
166f1580f4eSBarry Smith            D A M D^{-1} z = D b for right preconditioning
167f1580f4eSBarry Smith .ve
1684b9ad928SBarry Smith 
169f1580f4eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by D. `PCDiagonalScaleRight()` scales a vector by D^{-1}.
1704b9ad928SBarry Smith 
171db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCGetDiagonalScale()`
1724b9ad928SBarry Smith @*/
173d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetDiagonalScale(PC pc, Vec s)
174d71ae5a4SJacob Faibussowitsch {
1754b9ad928SBarry Smith   PetscFunctionBegin;
1760700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1770700a824SBarry Smith   PetscValidHeaderSpecific(s, VEC_CLASSID, 2);
1784b9ad928SBarry Smith   pc->diagonalscale = PETSC_TRUE;
1792fa5cd67SKarl Rupp 
1809566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)s));
1819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc->diagonalscaleleft));
1822fa5cd67SKarl Rupp 
1834b9ad928SBarry Smith   pc->diagonalscaleleft = s;
1842fa5cd67SKarl Rupp 
1859566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(s, &pc->diagonalscaleright));
1869566063dSJacob Faibussowitsch   PetscCall(VecCopy(s, pc->diagonalscaleright));
1879566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(pc->diagonalscaleright));
1883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1894b9ad928SBarry Smith }
1904b9ad928SBarry Smith 
191bc08b0f1SBarry Smith /*@
192c5eb9154SBarry Smith   PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes.
1934b9ad928SBarry Smith 
194c3339decSBarry Smith   Logically Collective
1954b9ad928SBarry Smith 
1964b9ad928SBarry Smith   Input Parameters:
1974b9ad928SBarry Smith + pc  - the preconditioner context
1984b9ad928SBarry Smith . in  - input vector
199a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in)
2004b9ad928SBarry Smith 
2014b9ad928SBarry Smith   Level: intermediate
2024b9ad928SBarry Smith 
20395452b02SPatrick Sanan   Notes:
20495452b02SPatrick Sanan   The system solved via the Krylov method is
205f1580f4eSBarry Smith .vb
206f1580f4eSBarry Smith         D M A D^{-1} y = D M b  for left preconditioning or
207f1580f4eSBarry Smith         D A M D^{-1} z = D b for right preconditioning
208f1580f4eSBarry Smith .ve
2094b9ad928SBarry Smith 
210f1580f4eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by D. `PCDiagonalScaleRight()` scales a vector by D^{-1}.
2114b9ad928SBarry Smith 
2124b9ad928SBarry Smith   If diagonal scaling is turned off and in is not out then in is copied to out
2134b9ad928SBarry Smith 
214db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleSet()`, `PCDiagonalScaleRight()`, `PCDiagonalScale()`
2154b9ad928SBarry Smith @*/
216d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleLeft(PC pc, Vec in, Vec out)
217d71ae5a4SJacob Faibussowitsch {
2184b9ad928SBarry Smith   PetscFunctionBegin;
2190700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2200700a824SBarry Smith   PetscValidHeaderSpecific(in, VEC_CLASSID, 2);
2210700a824SBarry Smith   PetscValidHeaderSpecific(out, VEC_CLASSID, 3);
2224b9ad928SBarry Smith   if (pc->diagonalscale) {
2239566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(out, pc->diagonalscaleleft, in));
2244b9ad928SBarry Smith   } else if (in != out) {
2259566063dSJacob Faibussowitsch     PetscCall(VecCopy(in, out));
2264b9ad928SBarry Smith   }
2273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2284b9ad928SBarry Smith }
2294b9ad928SBarry Smith 
230bc08b0f1SBarry Smith /*@
2314b9ad928SBarry Smith   PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
2324b9ad928SBarry Smith 
233c3339decSBarry Smith   Logically Collective
2344b9ad928SBarry Smith 
2354b9ad928SBarry Smith   Input Parameters:
2364b9ad928SBarry Smith + pc  - the preconditioner context
2374b9ad928SBarry Smith . in  - input vector
238a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in)
2394b9ad928SBarry Smith 
2404b9ad928SBarry Smith   Level: intermediate
2414b9ad928SBarry Smith 
24295452b02SPatrick Sanan   Notes:
24395452b02SPatrick Sanan   The system solved via the Krylov method is
244f1580f4eSBarry Smith .vb
245f1580f4eSBarry Smith         D M A D^{-1} y = D M b  for left preconditioning or
246f1580f4eSBarry Smith         D A M D^{-1} z = D b for right preconditioning
247f1580f4eSBarry Smith .ve
2484b9ad928SBarry Smith 
249f1580f4eSBarry Smith   `PCDiagonalScaleLeft()` scales a vector by D. `PCDiagonalScaleRight()` scales a vector by D^{-1}.
2504b9ad928SBarry Smith 
2514b9ad928SBarry Smith   If diagonal scaling is turned off and in is not out then in is copied to out
2524b9ad928SBarry Smith 
253db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleSet()`, `PCDiagonalScale()`
2544b9ad928SBarry Smith @*/
255d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleRight(PC pc, Vec in, Vec out)
256d71ae5a4SJacob Faibussowitsch {
2574b9ad928SBarry Smith   PetscFunctionBegin;
2580700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2590700a824SBarry Smith   PetscValidHeaderSpecific(in, VEC_CLASSID, 2);
2600700a824SBarry Smith   PetscValidHeaderSpecific(out, VEC_CLASSID, 3);
2614b9ad928SBarry Smith   if (pc->diagonalscale) {
2629566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(out, pc->diagonalscaleright, in));
2634b9ad928SBarry Smith   } else if (in != out) {
2649566063dSJacob Faibussowitsch     PetscCall(VecCopy(in, out));
2654b9ad928SBarry Smith   }
2663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2674b9ad928SBarry Smith }
2684b9ad928SBarry Smith 
26949517cdeSBarry Smith /*@
27049517cdeSBarry Smith   PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the
271f1580f4eSBarry Smith   operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
272f1580f4eSBarry Smith   `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.
27349517cdeSBarry Smith 
274c3339decSBarry Smith   Logically Collective
27549517cdeSBarry Smith 
27649517cdeSBarry Smith   Input Parameters:
27749517cdeSBarry Smith + pc  - the preconditioner context
278f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)
27949517cdeSBarry Smith 
28049517cdeSBarry Smith   Options Database Key:
281147403d9SBarry Smith . -pc_use_amat <true,false> - use the amat to apply the operator
28249517cdeSBarry Smith 
28320f4b53cSBarry Smith   Level: intermediate
28420f4b53cSBarry Smith 
285f1580f4eSBarry Smith   Note:
28649517cdeSBarry Smith   For the common case in which the linear system matrix and the matrix used to construct the
28749517cdeSBarry Smith   preconditioner are identical, this routine is does nothing.
28849517cdeSBarry Smith 
289f1580f4eSBarry Smith .seealso: `PC`, `PCGetUseAmat()`, `PCBJACOBI`, `PGMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`
29049517cdeSBarry Smith @*/
291d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUseAmat(PC pc, PetscBool flg)
292d71ae5a4SJacob Faibussowitsch {
29349517cdeSBarry Smith   PetscFunctionBegin;
29449517cdeSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
29549517cdeSBarry Smith   pc->useAmat = flg;
2963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29749517cdeSBarry Smith }
29849517cdeSBarry Smith 
299422a814eSBarry Smith /*@
300f1580f4eSBarry Smith   PCSetErrorIfFailure - Causes `PC` to generate an error if a FPE, for example a zero pivot, is detected.
301422a814eSBarry Smith 
302c3339decSBarry Smith   Logically Collective
303422a814eSBarry Smith 
304422a814eSBarry Smith   Input Parameters:
305422a814eSBarry Smith + pc  - iterative context obtained from PCCreate()
306f1580f4eSBarry Smith - flg - `PETSC_TRUE` indicates you want the error generated
307422a814eSBarry Smith 
308422a814eSBarry Smith   Level: advanced
309422a814eSBarry Smith 
310422a814eSBarry Smith   Notes:
311f1580f4eSBarry Smith   Normally PETSc continues if a linear solver fails due to a failed setup of a preconditioner, you can call `KSPGetConvergedReason()` after a `KSPSolve()`
312422a814eSBarry 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
313422a814eSBarry Smith   detected.
314422a814eSBarry Smith 
315422a814eSBarry Smith   This is propagated into KSPs used by this PC, which then propagate it into PCs used by those KSPs
316422a814eSBarry Smith 
317f1580f4eSBarry Smith .seealso: `PC`, `KSPSetErrorIfNotConverged()`, `PCGetInitialGuessNonzero()`, `PCSetInitialGuessKnoll()`, `PCGetInitialGuessKnoll()`
318422a814eSBarry Smith @*/
319d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetErrorIfFailure(PC pc, PetscBool flg)
320d71ae5a4SJacob Faibussowitsch {
321422a814eSBarry Smith   PetscFunctionBegin;
322422a814eSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
323422a814eSBarry Smith   PetscValidLogicalCollectiveBool(pc, flg, 2);
324422a814eSBarry Smith   pc->erroriffailure = flg;
3253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
326422a814eSBarry Smith }
327422a814eSBarry Smith 
32849517cdeSBarry Smith /*@
32949517cdeSBarry Smith   PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the
330f1580f4eSBarry Smith   operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
331f1580f4eSBarry Smith   `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.
33249517cdeSBarry Smith 
333c3339decSBarry Smith   Logically Collective
33449517cdeSBarry Smith 
33549517cdeSBarry Smith   Input Parameter:
33649517cdeSBarry Smith . pc - the preconditioner context
33749517cdeSBarry Smith 
33849517cdeSBarry Smith   Output Parameter:
339f1580f4eSBarry Smith . flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)
34049517cdeSBarry Smith 
34120f4b53cSBarry Smith   Level: intermediate
34220f4b53cSBarry Smith 
343f1580f4eSBarry Smith   Note:
34449517cdeSBarry Smith   For the common case in which the linear system matrix and the matrix used to construct the
34549517cdeSBarry Smith   preconditioner are identical, this routine is does nothing.
34649517cdeSBarry Smith 
347f1580f4eSBarry Smith .seealso: `PC`, `PCSetUseAmat()`, `PCBJACOBI`, `PGMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`
34849517cdeSBarry Smith @*/
349d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetUseAmat(PC pc, PetscBool *flg)
350d71ae5a4SJacob Faibussowitsch {
35149517cdeSBarry Smith   PetscFunctionBegin;
35249517cdeSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
35349517cdeSBarry Smith   *flg = pc->useAmat;
3543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35549517cdeSBarry Smith }
35649517cdeSBarry Smith 
357f39d8e23SSatish Balay /*@
3583821be0aSBarry Smith   PCSetKSPNestLevel - sets the amount of nesting the `KSP` that contains this `PC` has
3593821be0aSBarry Smith 
3603821be0aSBarry Smith   Collective
3613821be0aSBarry Smith 
3623821be0aSBarry Smith   Input Parameters:
3633821be0aSBarry Smith + pc    - the `PC`
3643821be0aSBarry Smith - level - the nest level
3653821be0aSBarry Smith 
3663821be0aSBarry Smith   Level: developer
3673821be0aSBarry Smith 
3683821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCGetKSPNestLevel()`, `KSPSetNestLevel()`
3693821be0aSBarry Smith @*/
3703821be0aSBarry Smith PetscErrorCode PCSetKSPNestLevel(PC pc, PetscInt level)
3713821be0aSBarry Smith {
3723821be0aSBarry Smith   PetscFunctionBegin;
3737a99bfcaSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
3747a99bfcaSBarry Smith   PetscValidLogicalCollectiveInt(pc, level, 2);
3753821be0aSBarry Smith   pc->kspnestlevel = level;
3763821be0aSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3773821be0aSBarry Smith }
3783821be0aSBarry Smith 
3793821be0aSBarry Smith /*@
3803821be0aSBarry Smith   PCGetKSPNestLevel - gets the amount of nesting the `KSP` that contains this `PC` has
3813821be0aSBarry Smith 
3827a99bfcaSBarry Smith   Not Collective
3833821be0aSBarry Smith 
3843821be0aSBarry Smith   Input Parameter:
3853821be0aSBarry Smith . pc - the `PC`
3863821be0aSBarry Smith 
3873821be0aSBarry Smith   Output Parameter:
3883821be0aSBarry Smith . level - the nest level
3893821be0aSBarry Smith 
3903821be0aSBarry Smith   Level: developer
3913821be0aSBarry Smith 
3923821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `KSPGetNestLevel()`
3933821be0aSBarry Smith @*/
3943821be0aSBarry Smith PetscErrorCode PCGetKSPNestLevel(PC pc, PetscInt *level)
3953821be0aSBarry Smith {
3963821be0aSBarry Smith   PetscFunctionBegin;
3977a99bfcaSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
3987a99bfcaSBarry Smith   PetscAssertPointer(level, 2);
3993821be0aSBarry Smith   *level = pc->kspnestlevel;
4003821be0aSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
4013821be0aSBarry Smith }
4023821be0aSBarry Smith 
4033821be0aSBarry Smith /*@
404f1580f4eSBarry Smith   PCCreate - Creates a preconditioner context, `PC`
4054b9ad928SBarry Smith 
406d083f849SBarry Smith   Collective
4074b9ad928SBarry Smith 
4084b9ad928SBarry Smith   Input Parameter:
4094b9ad928SBarry Smith . comm - MPI communicator
4104b9ad928SBarry Smith 
4114b9ad928SBarry Smith   Output Parameter:
4127a99bfcaSBarry Smith . newpc - location to put the preconditioner context
4134b9ad928SBarry Smith 
41420f4b53cSBarry Smith   Level: developer
41520f4b53cSBarry Smith 
416f1580f4eSBarry Smith   Note:
417f1580f4eSBarry 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`
418f1580f4eSBarry Smith   in parallel. For dense matrices it is always `PCNONE`.
4194b9ad928SBarry Smith 
420f1580f4eSBarry Smith .seealso: `PC`, `PCSetUp()`, `PCApply()`, `PCDestroy()`
4214b9ad928SBarry Smith @*/
422d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreate(MPI_Comm comm, PC *newpc)
423d71ae5a4SJacob Faibussowitsch {
4244b9ad928SBarry Smith   PC pc;
4254b9ad928SBarry Smith 
4264b9ad928SBarry Smith   PetscFunctionBegin;
4274f572ea9SToby Isaac   PetscAssertPointer(newpc, 2);
4280a545947SLisandro Dalcin   *newpc = NULL;
4299566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
4304b9ad928SBarry Smith 
4319566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(pc, PC_CLASSID, "PC", "Preconditioner", "PC", comm, PCDestroy, PCView));
432a7d21a52SLisandro Dalcin 
4330a545947SLisandro Dalcin   pc->mat                  = NULL;
4340a545947SLisandro Dalcin   pc->pmat                 = NULL;
4354b9ad928SBarry Smith   pc->setupcalled          = 0;
4360c24e6a1SHong Zhang   pc->setfromoptionscalled = 0;
4370a545947SLisandro Dalcin   pc->data                 = NULL;
4384b9ad928SBarry Smith   pc->diagonalscale        = PETSC_FALSE;
4390a545947SLisandro Dalcin   pc->diagonalscaleleft    = NULL;
4400a545947SLisandro Dalcin   pc->diagonalscaleright   = NULL;
4414b9ad928SBarry Smith 
4420a545947SLisandro Dalcin   pc->modifysubmatrices  = NULL;
4430a545947SLisandro Dalcin   pc->modifysubmatricesP = NULL;
4442fa5cd67SKarl Rupp 
445a7d21a52SLisandro Dalcin   *newpc = pc;
4463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4474b9ad928SBarry Smith }
4484b9ad928SBarry Smith 
4494b9ad928SBarry Smith /*@
4504b9ad928SBarry Smith   PCApply - Applies the preconditioner to a vector.
4514b9ad928SBarry Smith 
452c3339decSBarry Smith   Collective
4534b9ad928SBarry Smith 
4544b9ad928SBarry Smith   Input Parameters:
4554b9ad928SBarry Smith + pc - the preconditioner context
4564b9ad928SBarry Smith - x  - input vector
4574b9ad928SBarry Smith 
4584b9ad928SBarry Smith   Output Parameter:
4594b9ad928SBarry Smith . y - output vector
4604b9ad928SBarry Smith 
4614b9ad928SBarry Smith   Level: developer
4624b9ad928SBarry Smith 
463f1580f4eSBarry Smith .seealso: `PC`, `PCApplyTranspose()`, `PCApplyBAorAB()`
4644b9ad928SBarry Smith @*/
465d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApply(PC pc, Vec x, Vec y)
466d71ae5a4SJacob Faibussowitsch {
467106e20bfSBarry Smith   PetscInt m, n, mv, nv;
4684b9ad928SBarry Smith 
4694b9ad928SBarry Smith   PetscFunctionBegin;
4700700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4710700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
4720700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
4737827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
474e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
475540bdfbdSVaclav Hapla   /* use pmat to check vector sizes since for KSPLSQR the pmat may be of a different size than mat */
4769566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &m, &n));
4779566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x, &mv));
4789566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(y, &nv));
479540bdfbdSVaclav Hapla   /* check pmat * y = x is feasible */
48063a3b9bcSJacob 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);
48163a3b9bcSJacob 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);
4829566063dSJacob Faibussowitsch   PetscCall(VecSetErrorIfLocked(y, 3));
483106e20bfSBarry Smith 
4849566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
4859566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
4869566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
487dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, apply, x, y);
4889566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
489e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
4909566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
4913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4924b9ad928SBarry Smith }
4934b9ad928SBarry Smith 
4944b9ad928SBarry Smith /*@
495f1580f4eSBarry Smith   PCMatApply - Applies the preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApply()`, Y and X must be different matrices.
496c41ea70eSPierre Jolivet 
497c3339decSBarry Smith   Collective
498c677e75fSPierre Jolivet 
499c677e75fSPierre Jolivet   Input Parameters:
500c677e75fSPierre Jolivet + pc - the preconditioner context
501c677e75fSPierre Jolivet - X  - block of input vectors
502c677e75fSPierre Jolivet 
503c677e75fSPierre Jolivet   Output Parameter:
504c677e75fSPierre Jolivet . Y - block of output vectors
505c677e75fSPierre Jolivet 
506c677e75fSPierre Jolivet   Level: developer
507c677e75fSPierre Jolivet 
508f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `KSPMatSolve()`
509c677e75fSPierre Jolivet @*/
510d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMatApply(PC pc, Mat X, Mat Y)
511d71ae5a4SJacob Faibussowitsch {
512c677e75fSPierre Jolivet   Mat       A;
513c677e75fSPierre Jolivet   Vec       cy, cx;
514bd82155bSPierre Jolivet   PetscInt  m1, M1, m2, M2, n1, N1, n2, N2, m3, M3, n3, N3;
515c677e75fSPierre Jolivet   PetscBool match;
516c677e75fSPierre Jolivet 
517c677e75fSPierre Jolivet   PetscFunctionBegin;
518c677e75fSPierre Jolivet   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
519e2b200f6SPierre Jolivet   PetscValidHeaderSpecific(X, MAT_CLASSID, 2);
520e2b200f6SPierre Jolivet   PetscValidHeaderSpecific(Y, MAT_CLASSID, 3);
521e2b200f6SPierre Jolivet   PetscCheckSameComm(pc, 1, X, 2);
522e2b200f6SPierre Jolivet   PetscCheckSameComm(pc, 1, Y, 3);
5237827d75bSBarry Smith   PetscCheck(Y != X, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "Y and X must be different matrices");
5249566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &A));
5259566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m3, &n3));
5269566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m2, &n2));
5279566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Y, &m1, &n1));
5289566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M3, &N3));
5299566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, &M2, &N2));
5309566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Y, &M1, &N1));
53163a3b9bcSJacob 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);
53263a3b9bcSJacob 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);
53363a3b9bcSJacob 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);
5349566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &match, MATSEQDENSE, MATMPIDENSE, ""));
53528b400f6SJacob Faibussowitsch   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of output vectors not stored in a dense Mat");
5369566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
53728b400f6SJacob Faibussowitsch   PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of input vectors not stored in a dense Mat");
5389566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
539c677e75fSPierre Jolivet   if (pc->ops->matapply) {
5409566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0));
541dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, matapply, X, Y);
5429566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0));
543c677e75fSPierre Jolivet   } else {
5449566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "PC type %s applying column by column\n", ((PetscObject)pc)->type_name));
545bd82155bSPierre Jolivet     for (n1 = 0; n1 < N1; ++n1) {
5469566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecRead(X, n1, &cx));
5479566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecWrite(Y, n1, &cy));
5489566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, cx, cy));
5499566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecWrite(Y, n1, &cy));
5509566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecRead(X, n1, &cx));
551c677e75fSPierre Jolivet     }
552c677e75fSPierre Jolivet   }
5533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
554c677e75fSPierre Jolivet }
555c677e75fSPierre Jolivet 
556c677e75fSPierre Jolivet /*@
5574b9ad928SBarry Smith   PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
5584b9ad928SBarry Smith 
559c3339decSBarry Smith   Collective
5604b9ad928SBarry Smith 
5614b9ad928SBarry Smith   Input Parameters:
5624b9ad928SBarry Smith + pc - the preconditioner context
5634b9ad928SBarry Smith - x  - input vector
5644b9ad928SBarry Smith 
5654b9ad928SBarry Smith   Output Parameter:
5664b9ad928SBarry Smith . y - output vector
5674b9ad928SBarry Smith 
56820f4b53cSBarry Smith   Level: developer
56920f4b53cSBarry Smith 
570f1580f4eSBarry Smith   Note:
571f1580f4eSBarry Smith   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
5724b9ad928SBarry Smith 
573f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplySymmetricRight()`
5744b9ad928SBarry Smith @*/
575d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricLeft(PC pc, Vec x, Vec y)
576d71ae5a4SJacob Faibussowitsch {
5774b9ad928SBarry Smith   PetscFunctionBegin;
5780700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
5790700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
5800700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
5817827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
582e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
5839566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
5849566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
5859566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplySymmetricLeft, pc, x, y, 0));
586dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applysymmetricleft, x, y);
5879566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplySymmetricLeft, pc, x, y, 0));
5889566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
589e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
5903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5914b9ad928SBarry Smith }
5924b9ad928SBarry Smith 
5934b9ad928SBarry Smith /*@
5944b9ad928SBarry Smith   PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
5954b9ad928SBarry Smith 
596c3339decSBarry Smith   Collective
5974b9ad928SBarry Smith 
5984b9ad928SBarry Smith   Input Parameters:
5994b9ad928SBarry Smith + pc - the preconditioner context
6004b9ad928SBarry Smith - x  - input vector
6014b9ad928SBarry Smith 
6024b9ad928SBarry Smith   Output Parameter:
6034b9ad928SBarry Smith . y - output vector
6044b9ad928SBarry Smith 
6054b9ad928SBarry Smith   Level: developer
6064b9ad928SBarry Smith 
607f1580f4eSBarry Smith   Note:
608f1580f4eSBarry Smith   Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
6094b9ad928SBarry Smith 
610f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplySymmetricLeft()`
6114b9ad928SBarry Smith @*/
612d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricRight(PC pc, Vec x, Vec y)
613d71ae5a4SJacob Faibussowitsch {
6144b9ad928SBarry Smith   PetscFunctionBegin;
6150700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6160700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6170700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6187827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
619e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
6209566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
6219566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
6229566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplySymmetricRight, pc, x, y, 0));
623dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applysymmetricright, x, y);
6249566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplySymmetricRight, pc, x, y, 0));
6259566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
626e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
6273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6284b9ad928SBarry Smith }
6294b9ad928SBarry Smith 
6304b9ad928SBarry Smith /*@
6314b9ad928SBarry Smith   PCApplyTranspose - Applies the transpose of preconditioner to a vector.
6324b9ad928SBarry Smith 
633c3339decSBarry Smith   Collective
6344b9ad928SBarry Smith 
6354b9ad928SBarry Smith   Input Parameters:
6364b9ad928SBarry Smith + pc - the preconditioner context
6374b9ad928SBarry Smith - x  - input vector
6384b9ad928SBarry Smith 
6394b9ad928SBarry Smith   Output Parameter:
6404b9ad928SBarry Smith . y - output vector
6414b9ad928SBarry Smith 
64220f4b53cSBarry Smith   Level: developer
64320f4b53cSBarry Smith 
644f1580f4eSBarry Smith   Note:
64595452b02SPatrick Sanan   For complex numbers this applies the non-Hermitian transpose.
6464c97465dSBarry Smith 
647feefa0e1SJacob Faibussowitsch   Developer Notes:
648f1580f4eSBarry Smith   We need to implement a `PCApplyHermitianTranspose()`
6494c97465dSBarry Smith 
650f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplyBAorAB()`, `PCApplyBAorABTranspose()`, `PCApplyTransposeExists()`
6514b9ad928SBarry Smith @*/
652d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTranspose(PC pc, Vec x, Vec y)
653d71ae5a4SJacob Faibussowitsch {
6544b9ad928SBarry Smith   PetscFunctionBegin;
6550700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6560700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
6570700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
6587827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
659e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
6609566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
6619566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(x));
6629566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
663dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applytranspose, x, y);
6649566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
6659566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(x));
666e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
6673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6684b9ad928SBarry Smith }
6694b9ad928SBarry Smith 
6704b9ad928SBarry Smith /*@
6719cc050e5SLisandro Dalcin   PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
6724b9ad928SBarry Smith 
673c3339decSBarry Smith   Collective
6744b9ad928SBarry Smith 
675f1580f4eSBarry Smith   Input Parameter:
6764b9ad928SBarry Smith . pc - the preconditioner context
6774b9ad928SBarry Smith 
6784b9ad928SBarry Smith   Output Parameter:
679f1580f4eSBarry Smith . flg - `PETSC_TRUE` if a transpose operation is defined
6804b9ad928SBarry Smith 
6814b9ad928SBarry Smith   Level: developer
6824b9ad928SBarry Smith 
683f1580f4eSBarry Smith .seealso: `PC`, `PCApplyTranspose()`
6844b9ad928SBarry Smith @*/
685d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTransposeExists(PC pc, PetscBool *flg)
686d71ae5a4SJacob Faibussowitsch {
6874b9ad928SBarry Smith   PetscFunctionBegin;
6880700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6894f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
6906ce5e81cSLisandro Dalcin   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
6916ce5e81cSLisandro Dalcin   else *flg = PETSC_FALSE;
6923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6934b9ad928SBarry Smith }
6944b9ad928SBarry Smith 
6954b9ad928SBarry Smith /*@
6961a2f9f99SBarry Smith   PCApplyBAorAB - Applies the preconditioner and operator to a vector. y = B*A*x or y = A*B*x.
6974b9ad928SBarry Smith 
698c3339decSBarry Smith   Collective
6994b9ad928SBarry Smith 
7004b9ad928SBarry Smith   Input Parameters:
7014b9ad928SBarry Smith + pc   - the preconditioner context
702f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
7034b9ad928SBarry Smith . x    - input vector
7044b9ad928SBarry Smith - work - work vector
7054b9ad928SBarry Smith 
7064b9ad928SBarry Smith   Output Parameter:
7074b9ad928SBarry Smith . y - output vector
7084b9ad928SBarry Smith 
7094b9ad928SBarry Smith   Level: developer
7104b9ad928SBarry Smith 
711f1580f4eSBarry Smith   Note:
712f1580f4eSBarry 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
713f1580f4eSBarry Smith   specific `KSPSolve()` method must also be written to handle the post-solve "correction" for the diagonal scaling.
71413e0d083SBarry Smith 
715f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorABTranspose()`
7164b9ad928SBarry Smith @*/
717d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorAB(PC pc, PCSide side, Vec x, Vec y, Vec work)
718d71ae5a4SJacob Faibussowitsch {
7194b9ad928SBarry Smith   PetscFunctionBegin;
7200700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
721a69c7061SStefano Zampini   PetscValidLogicalCollectiveEnum(pc, side, 2);
7220700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
7230700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
7240700a824SBarry Smith   PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
725a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, x, 3);
726a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, y, 4);
727a69c7061SStefano Zampini   PetscCheckSameComm(pc, 1, work, 5);
7287827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
7297827d75bSBarry 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");
7307827d75bSBarry Smith   PetscCheck(!pc->diagonalscale || side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot include diagonal scaling with symmetric preconditioner application");
731e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
7324b9ad928SBarry Smith 
7339566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
7344b9ad928SBarry Smith   if (pc->diagonalscale) {
7354b9ad928SBarry Smith     if (pc->ops->applyBA) {
7364b9ad928SBarry Smith       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
7379566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &work2));
7389566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, work2));
739dbbe0bcdSBarry Smith       PetscUseTypeMethod(pc, applyBA, side, work2, y, work);
7409566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7419566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&work2));
7424b9ad928SBarry Smith     } else if (side == PC_RIGHT) {
7439566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, y));
7449566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, y, work));
7459566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7469566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7474b9ad928SBarry Smith     } else if (side == PC_LEFT) {
7489566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleRight(pc, x, y));
7499566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, y, work));
7509566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, work, y));
7519566063dSJacob Faibussowitsch       PetscCall(PCDiagonalScaleLeft(pc, y, y));
7527827d75bSBarry Smith     } else PetscCheck(side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot provide diagonal scaling with symmetric application of preconditioner");
7534b9ad928SBarry Smith   } else {
7544b9ad928SBarry Smith     if (pc->ops->applyBA) {
755dbbe0bcdSBarry Smith       PetscUseTypeMethod(pc, applyBA, side, x, y, work);
7564b9ad928SBarry Smith     } else if (side == PC_RIGHT) {
7579566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, x, work));
7589566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7594b9ad928SBarry Smith     } else if (side == PC_LEFT) {
7609566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, x, work));
7619566063dSJacob Faibussowitsch       PetscCall(PCApply(pc, work, y));
7624b9ad928SBarry Smith     } else if (side == PC_SYMMETRIC) {
7634b9ad928SBarry Smith       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
7649566063dSJacob Faibussowitsch       PetscCall(PCApplySymmetricRight(pc, x, work));
7659566063dSJacob Faibussowitsch       PetscCall(MatMult(pc->mat, work, y));
7669566063dSJacob Faibussowitsch       PetscCall(VecCopy(y, work));
7679566063dSJacob Faibussowitsch       PetscCall(PCApplySymmetricLeft(pc, work, y));
7684b9ad928SBarry Smith     }
7694b9ad928SBarry Smith   }
770e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
7713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7724b9ad928SBarry Smith }
7734b9ad928SBarry Smith 
7744b9ad928SBarry Smith /*@
7754b9ad928SBarry Smith   PCApplyBAorABTranspose - Applies the transpose of the preconditioner
7764b9ad928SBarry Smith   and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
7779b3038f0SBarry Smith   NOT tr(B*A) = tr(A)*tr(B).
7784b9ad928SBarry Smith 
779c3339decSBarry Smith   Collective
7804b9ad928SBarry Smith 
7814b9ad928SBarry Smith   Input Parameters:
7824b9ad928SBarry Smith + pc   - the preconditioner context
783f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
7844b9ad928SBarry Smith . x    - input vector
7854b9ad928SBarry Smith - work - work vector
7864b9ad928SBarry Smith 
7874b9ad928SBarry Smith   Output Parameter:
7884b9ad928SBarry Smith . y - output vector
7894b9ad928SBarry Smith 
79020f4b53cSBarry Smith   Level: developer
79120f4b53cSBarry Smith 
792f1580f4eSBarry Smith   Note:
793f1580f4eSBarry 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
7949b3038f0SBarry Smith   defined by B'. This is why this has the funny form that it computes tr(B) * tr(A)
7959b3038f0SBarry Smith 
796f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorAB()`
7974b9ad928SBarry Smith @*/
798d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorABTranspose(PC pc, PCSide side, Vec x, Vec y, Vec work)
799d71ae5a4SJacob Faibussowitsch {
8004b9ad928SBarry Smith   PetscFunctionBegin;
8010700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8020700a824SBarry Smith   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
8030700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
8040700a824SBarry Smith   PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
8057827d75bSBarry Smith   PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
806e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
8074b9ad928SBarry Smith   if (pc->ops->applyBAtranspose) {
808dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, applyBAtranspose, side, x, y, work);
809e0f629ddSJacob Faibussowitsch     if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8103ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
8114b9ad928SBarry Smith   }
8127827d75bSBarry Smith   PetscCheck(side == PC_LEFT || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right or left");
8134b9ad928SBarry Smith 
8149566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
8154b9ad928SBarry Smith   if (side == PC_RIGHT) {
8169566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(pc, x, work));
8179566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->mat, work, y));
8184b9ad928SBarry Smith   } else if (side == PC_LEFT) {
8199566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->mat, x, work));
8209566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(pc, work, y));
8214b9ad928SBarry Smith   }
8224b9ad928SBarry Smith   /* add support for PC_SYMMETRIC */
823e0f629ddSJacob Faibussowitsch   if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
8243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8254b9ad928SBarry Smith }
8264b9ad928SBarry Smith 
8274b9ad928SBarry Smith /*@
8284b9ad928SBarry Smith   PCApplyRichardsonExists - Determines whether a particular preconditioner has a
8294b9ad928SBarry Smith   built-in fast application of Richardson's method.
8304b9ad928SBarry Smith 
8314b9ad928SBarry Smith   Not Collective
8324b9ad928SBarry Smith 
8334b9ad928SBarry Smith   Input Parameter:
8344b9ad928SBarry Smith . pc - the preconditioner
8354b9ad928SBarry Smith 
8364b9ad928SBarry Smith   Output Parameter:
837f1580f4eSBarry Smith . exists - `PETSC_TRUE` or `PETSC_FALSE`
8384b9ad928SBarry Smith 
8394b9ad928SBarry Smith   Level: developer
8404b9ad928SBarry Smith 
841f1580f4eSBarry Smith .seealso: `PC`, `PCRICHARDSON`, `PCApplyRichardson()`
8424b9ad928SBarry Smith @*/
843d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyRichardsonExists(PC pc, PetscBool *exists)
844d71ae5a4SJacob Faibussowitsch {
8454b9ad928SBarry Smith   PetscFunctionBegin;
8460700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8474f572ea9SToby Isaac   PetscAssertPointer(exists, 2);
8484b9ad928SBarry Smith   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
8494b9ad928SBarry Smith   else *exists = PETSC_FALSE;
8503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8514b9ad928SBarry Smith }
8524b9ad928SBarry Smith 
8534b9ad928SBarry Smith /*@
8544b9ad928SBarry Smith   PCApplyRichardson - Applies several steps of Richardson iteration with
8554b9ad928SBarry Smith   the particular preconditioner. This routine is usually used by the
8564b9ad928SBarry Smith   Krylov solvers and not the application code directly.
8574b9ad928SBarry Smith 
858c3339decSBarry Smith   Collective
8594b9ad928SBarry Smith 
8604b9ad928SBarry Smith   Input Parameters:
8614b9ad928SBarry Smith + pc        - the preconditioner context
8627319c654SBarry Smith . b         - the right hand side
8634b9ad928SBarry Smith . w         - one work vector
8644b9ad928SBarry Smith . rtol      - relative decrease in residual norm convergence criteria
86570441072SBarry Smith . abstol    - absolute residual norm convergence criteria
8664b9ad928SBarry Smith . dtol      - divergence residual norm increase criteria
8677319c654SBarry Smith . its       - the number of iterations to apply.
8687319c654SBarry Smith - guesszero - if the input x contains nonzero initial guess
8694b9ad928SBarry Smith 
870d8d19677SJose E. Roman   Output Parameters:
8714d0a8057SBarry Smith + outits - number of iterations actually used (for SOR this always equals its)
8724d0a8057SBarry Smith . reason - the reason the apply terminated
873f1580f4eSBarry Smith - y      - the solution (also contains initial guess if guesszero is `PETSC_FALSE`
8744b9ad928SBarry Smith 
87520f4b53cSBarry Smith   Level: developer
87620f4b53cSBarry Smith 
8774b9ad928SBarry Smith   Notes:
8784b9ad928SBarry Smith   Most preconditioners do not support this function. Use the command
879f1580f4eSBarry Smith   `PCApplyRichardsonExists()` to determine if one does.
8804b9ad928SBarry Smith 
881f1580f4eSBarry Smith   Except for the `PCMG` this routine ignores the convergence tolerances
8824b9ad928SBarry Smith   and always runs for the number of iterations
8834b9ad928SBarry Smith 
884f1580f4eSBarry Smith .seealso: `PC`, `PCApplyRichardsonExists()`
8854b9ad928SBarry Smith @*/
886d71ae5a4SJacob 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)
887d71ae5a4SJacob Faibussowitsch {
8884b9ad928SBarry Smith   PetscFunctionBegin;
8890700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8900700a824SBarry Smith   PetscValidHeaderSpecific(b, VEC_CLASSID, 2);
8910700a824SBarry Smith   PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
8920700a824SBarry Smith   PetscValidHeaderSpecific(w, VEC_CLASSID, 4);
8937827d75bSBarry Smith   PetscCheck(b != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "b and y must be different vectors");
8949566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
895dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, applyrichardson, b, y, w, rtol, abstol, dtol, its, guesszero, outits, reason);
8963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8974b9ad928SBarry Smith }
8984b9ad928SBarry Smith 
899422a814eSBarry Smith /*@
900f1580f4eSBarry Smith   PCSetFailedReason - Sets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
9011b2b9847SBarry Smith 
902c3339decSBarry Smith   Logically Collective
9031b2b9847SBarry Smith 
904d8d19677SJose E. Roman   Input Parameters:
9051b2b9847SBarry Smith + pc     - the preconditioner context
9061b2b9847SBarry Smith - reason - the reason it failedx
9071b2b9847SBarry Smith 
9081b2b9847SBarry Smith   Level: advanced
9091b2b9847SBarry Smith 
910f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCFailedReason`
9111b2b9847SBarry Smith @*/
912d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFailedReason(PC pc, PCFailedReason reason)
913d71ae5a4SJacob Faibussowitsch {
9141b2b9847SBarry Smith   PetscFunctionBegin;
9156479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9161b2b9847SBarry Smith   pc->failedreason = reason;
9173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9181b2b9847SBarry Smith }
9191b2b9847SBarry Smith 
9201b2b9847SBarry Smith /*@
921f1580f4eSBarry Smith   PCGetFailedReason - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
922422a814eSBarry Smith 
923c3339decSBarry Smith   Logically Collective
924422a814eSBarry Smith 
925422a814eSBarry Smith   Input Parameter:
926422a814eSBarry Smith . pc - the preconditioner context
927422a814eSBarry Smith 
928422a814eSBarry Smith   Output Parameter:
9291b2b9847SBarry Smith . reason - the reason it failed
930422a814eSBarry Smith 
931422a814eSBarry Smith   Level: advanced
932422a814eSBarry Smith 
933f1580f4eSBarry Smith   Note:
934f1580f4eSBarry Smith   This is the maximum over reason over all ranks in the PC communicator. It is only valid after
9356479e835SStefano Zampini   a call `KSPCheckDot()` or  `KSPCheckNorm()` inside a `KSPSolve()` or `PCReduceFailedReason()`.
9366479e835SStefano Zampini   It is not valid immediately after a `PCSetUp()` or `PCApply()`, then use `PCGetFailedReasonRank()`
9371b2b9847SBarry Smith 
938feefa0e1SJacob Faibussowitsch .seealso: `PC`, ``PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReasonRank()`, `PCSetFailedReason()`
939422a814eSBarry Smith @*/
940d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason)
941d71ae5a4SJacob Faibussowitsch {
942422a814eSBarry Smith   PetscFunctionBegin;
9436479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9449754764cSHong Zhang   if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled;
9450ae0b32dSHong Zhang   else *reason = pc->failedreason;
9463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
947422a814eSBarry Smith }
948422a814eSBarry Smith 
9491b2b9847SBarry Smith /*@
950f1580f4eSBarry Smith   PCGetFailedReasonRank - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail on this MPI rank
9511b2b9847SBarry Smith 
952f1580f4eSBarry Smith   Not Collective
9531b2b9847SBarry Smith 
9541b2b9847SBarry Smith   Input Parameter:
9551b2b9847SBarry Smith . pc - the preconditioner context
9561b2b9847SBarry Smith 
9571b2b9847SBarry Smith   Output Parameter:
9581b2b9847SBarry Smith . reason - the reason it failed
9591b2b9847SBarry Smith 
96020f4b53cSBarry Smith   Level: advanced
96120f4b53cSBarry Smith 
962f1580f4eSBarry Smith   Note:
963f1580f4eSBarry Smith   Different ranks may have different reasons or no reason, see `PCGetFailedReason()`
9641b2b9847SBarry Smith 
9656479e835SStefano Zampini .seealso: `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`, `PCReduceFailedReason()`
9661b2b9847SBarry Smith @*/
967d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReasonRank(PC pc, PCFailedReason *reason)
968d71ae5a4SJacob Faibussowitsch {
9691b2b9847SBarry Smith   PetscFunctionBegin;
9706479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9711b2b9847SBarry Smith   if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled;
9721b2b9847SBarry Smith   else *reason = pc->failedreason;
9733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9741b2b9847SBarry Smith }
975422a814eSBarry Smith 
9766479e835SStefano Zampini /*@
9776479e835SStefano Zampini   PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC`
9786479e835SStefano Zampini 
9796479e835SStefano Zampini   Collective
9806479e835SStefano Zampini 
9816479e835SStefano Zampini   Input Parameter:
9826479e835SStefano Zampini . pc - the preconditioner context
9836479e835SStefano Zampini 
9846479e835SStefano Zampini   Level: advanced
9856479e835SStefano Zampini 
9866479e835SStefano Zampini   Note:
9876479e835SStefano Zampini   Different MPI ranks may have different reasons or no reason, see `PCGetFailedReason()`. This routine
9886479e835SStefano Zampini   makes them have a common value (failure if any MPI process had a failure).
9896479e835SStefano Zampini 
9906479e835SStefano Zampini .seealso: `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`
9916479e835SStefano Zampini @*/
9926479e835SStefano Zampini PetscErrorCode PCReduceFailedReason(PC pc)
9936479e835SStefano Zampini {
9946479e835SStefano Zampini   PetscInt buf;
9956479e835SStefano Zampini 
9966479e835SStefano Zampini   PetscFunctionBegin;
9976479e835SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9986479e835SStefano Zampini   buf = (PetscInt)pc->failedreason;
9996479e835SStefano Zampini   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
10006479e835SStefano Zampini   pc->failedreason = (PCFailedReason)buf;
10016479e835SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
10026479e835SStefano Zampini }
10036479e835SStefano Zampini 
1004c00cb57fSBarry Smith /*  Next line needed to deactivate KSP_Solve logging */
1005c00cb57fSBarry Smith #include <petsc/private/kspimpl.h>
1006c00cb57fSBarry Smith 
10074b9ad928SBarry Smith /*
10084b9ad928SBarry Smith       a setupcall of 0 indicates never setup,
100923ee1639SBarry Smith                      1 indicates has been previously setup
1010422a814eSBarry Smith                     -1 indicates a PCSetUp() was attempted and failed
10114b9ad928SBarry Smith */
10124b9ad928SBarry Smith /*@
10134b9ad928SBarry Smith   PCSetUp - Prepares for the use of a preconditioner.
10144b9ad928SBarry Smith 
1015c3339decSBarry Smith   Collective
10164b9ad928SBarry Smith 
10174b9ad928SBarry Smith   Input Parameter:
10184b9ad928SBarry Smith . pc - the preconditioner context
10194b9ad928SBarry Smith 
10204b9ad928SBarry Smith   Level: developer
10214b9ad928SBarry Smith 
1022f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`
10234b9ad928SBarry Smith @*/
1024d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp(PC pc)
1025d71ae5a4SJacob Faibussowitsch {
1026566e8bf2SBarry Smith   const char      *def;
1027fc9b51d3SBarry Smith   PetscObjectState matstate, matnonzerostate;
10284b9ad928SBarry Smith 
10294b9ad928SBarry Smith   PetscFunctionBegin;
10300700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
103128b400f6SJacob Faibussowitsch   PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first");
10326ce5e81cSLisandro Dalcin 
103323ee1639SBarry Smith   if (pc->setupcalled && pc->reusepreconditioner) {
10349566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n"));
10353ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10364b9ad928SBarry Smith   }
10374b9ad928SBarry Smith 
10389566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate));
10399566063dSJacob Faibussowitsch   PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate));
104023ee1639SBarry Smith   if (!pc->setupcalled) {
10419566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Setting up PC for first time\n"));
104223ee1639SBarry Smith     pc->flag = DIFFERENT_NONZERO_PATTERN;
10439df67409SStefano Zampini   } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS);
10449df67409SStefano Zampini   else {
10459df67409SStefano Zampini     if (matnonzerostate != pc->matnonzerostate) {
10469566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n"));
104723ee1639SBarry Smith       pc->flag = DIFFERENT_NONZERO_PATTERN;
104823ee1639SBarry Smith     } else {
10499566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n"));
105023ee1639SBarry Smith       pc->flag = SAME_NONZERO_PATTERN;
105123ee1639SBarry Smith     }
105223ee1639SBarry Smith   }
105323ee1639SBarry Smith   pc->matstate        = matstate;
1054fc9b51d3SBarry Smith   pc->matnonzerostate = matnonzerostate;
105523ee1639SBarry Smith 
10567adad957SLisandro Dalcin   if (!((PetscObject)pc)->type_name) {
10579566063dSJacob Faibussowitsch     PetscCall(PCGetDefaultType_Private(pc, &def));
10589566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, def));
1059566e8bf2SBarry Smith   }
10604b9ad928SBarry Smith 
10619566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
10629566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure));
10639566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0));
10644b9ad928SBarry Smith   if (pc->ops->setup) {
1065c00cb57fSBarry Smith     /* do not log solves and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */
10669566063dSJacob Faibussowitsch     PetscCall(KSPInitializePackage());
10679566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePush(KSP_Solve));
10689566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePush(PC_Apply));
1069dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, setup);
10709566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePop(KSP_Solve));
10719566063dSJacob Faibussowitsch     PetscCall(PetscLogEventDeactivatePop(PC_Apply));
10724b9ad928SBarry Smith   }
10739566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0));
1074422a814eSBarry Smith   if (!pc->setupcalled) pc->setupcalled = 1;
10753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10764b9ad928SBarry Smith }
10774b9ad928SBarry Smith 
10784b9ad928SBarry Smith /*@
10794b9ad928SBarry Smith   PCSetUpOnBlocks - Sets up the preconditioner for each block in
10804b9ad928SBarry Smith   the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
10814b9ad928SBarry Smith   methods.
10824b9ad928SBarry Smith 
1083c3339decSBarry Smith   Collective
10844b9ad928SBarry Smith 
1085f1580f4eSBarry Smith   Input Parameter:
10864b9ad928SBarry Smith . pc - the preconditioner context
10874b9ad928SBarry Smith 
10884b9ad928SBarry Smith   Level: developer
10894b9ad928SBarry Smith 
1090f1580f4eSBarry Smith   Note:
1091f1580f4eSBarry Smith   For nested preconditioners such as `PCBJACOBI` `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is
1092f1580f4eSBarry Smith   called on the outer `PC`, this routine ensures it is called.
1093f1580f4eSBarry Smith 
1094feefa0e1SJacob Faibussowitsch .seealso: `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()`
10954b9ad928SBarry Smith @*/
1096d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUpOnBlocks(PC pc)
1097d71ae5a4SJacob Faibussowitsch {
10984b9ad928SBarry Smith   PetscFunctionBegin;
10990700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11003ba16761SJacob Faibussowitsch   if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS);
11019566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0));
1102dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, setuponblocks);
11039566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0));
11043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11054b9ad928SBarry Smith }
11064b9ad928SBarry Smith 
11074b9ad928SBarry Smith /*@C
11084b9ad928SBarry Smith   PCSetModifySubMatrices - Sets a user-defined routine for modifying the
110904c3f3b8SBarry Smith   submatrices that arise within certain subdomain-based preconditioners such as `PCASM`
11104b9ad928SBarry Smith 
1111c3339decSBarry Smith   Logically Collective
11124b9ad928SBarry Smith 
11134b9ad928SBarry Smith   Input Parameters:
11144b9ad928SBarry Smith + pc   - the preconditioner context
11154b9ad928SBarry Smith . func - routine for modifying the submatrices
111604c3f3b8SBarry Smith - ctx  - optional user-defined context (may be `NULL`)
11174b9ad928SBarry Smith 
111820f4b53cSBarry Smith   Calling sequence of `func`:
111920f4b53cSBarry Smith + pc     - the preconditioner context
112004c3f3b8SBarry Smith . nsub   - number of index sets
112120f4b53cSBarry Smith . row    - an array of index sets that contain the global row numbers
11224b9ad928SBarry Smith          that comprise each local submatrix
11234b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
11244b9ad928SBarry Smith          that comprise each local submatrix
11254b9ad928SBarry Smith . submat - array of local submatrices
11264b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
112704c3f3b8SBarry Smith          user-defined func routine (may be `NULL`)
11284b9ad928SBarry Smith 
112920f4b53cSBarry Smith   Level: advanced
113020f4b53cSBarry Smith 
11314b9ad928SBarry Smith   Notes:
113204c3f3b8SBarry Smith   The basic submatrices are extracted from the preconditioner matrix as
113304c3f3b8SBarry Smith   usual; the user can then alter these (for example, to set different boundary
113404c3f3b8SBarry Smith   conditions for each submatrix) before they are used for the local solves.
113504c3f3b8SBarry Smith 
1136f1580f4eSBarry Smith   `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and
1137f1580f4eSBarry Smith   `KSPSolve()`.
11384b9ad928SBarry Smith 
1139f1580f4eSBarry Smith   A routine set by `PCSetModifySubMatrices()` is currently called within
1140f1580f4eSBarry Smith   the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`)
11414b9ad928SBarry Smith   preconditioners.  All other preconditioners ignore this routine.
11424b9ad928SBarry Smith 
1143f1580f4eSBarry Smith .seealso: `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()`
11444b9ad928SBarry Smith @*/
114504c3f3b8SBarry Smith PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx), void *ctx)
1146d71ae5a4SJacob Faibussowitsch {
11474b9ad928SBarry Smith   PetscFunctionBegin;
11480700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11494b9ad928SBarry Smith   pc->modifysubmatrices  = func;
11504b9ad928SBarry Smith   pc->modifysubmatricesP = ctx;
11513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11524b9ad928SBarry Smith }
11534b9ad928SBarry Smith 
11544b9ad928SBarry Smith /*@C
11554b9ad928SBarry Smith   PCModifySubMatrices - Calls an optional user-defined routine within
1156f1580f4eSBarry Smith   certain preconditioners if one has been set with `PCSetModifySubMatrices()`.
11574b9ad928SBarry Smith 
1158c3339decSBarry Smith   Collective
11594b9ad928SBarry Smith 
11604b9ad928SBarry Smith   Input Parameters:
11614b9ad928SBarry Smith + pc     - the preconditioner context
11624b9ad928SBarry Smith . nsub   - the number of local submatrices
11634b9ad928SBarry Smith . row    - an array of index sets that contain the global row numbers
11644b9ad928SBarry Smith          that comprise each local submatrix
11654b9ad928SBarry Smith . col    - an array of index sets that contain the global column numbers
11664b9ad928SBarry Smith          that comprise each local submatrix
11674b9ad928SBarry Smith . submat - array of local submatrices
11684b9ad928SBarry Smith - ctx    - optional user-defined context for private data for the
11694b9ad928SBarry Smith          user-defined routine (may be null)
11704b9ad928SBarry Smith 
11714b9ad928SBarry Smith   Output Parameter:
11724b9ad928SBarry Smith . submat - array of local submatrices (the entries of which may
11734b9ad928SBarry Smith             have been modified)
11744b9ad928SBarry Smith 
117520f4b53cSBarry Smith   Level: developer
117620f4b53cSBarry Smith 
117704c3f3b8SBarry Smith   Note:
11784b9ad928SBarry Smith   The user should NOT generally call this routine, as it will
117904c3f3b8SBarry Smith   automatically be called within certain preconditioners.
11804b9ad928SBarry Smith 
1181f1580f4eSBarry Smith .seealso: `PC`, `PCSetModifySubMatrices()`
11824b9ad928SBarry Smith @*/
1183d71ae5a4SJacob Faibussowitsch PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx)
1184d71ae5a4SJacob Faibussowitsch {
11854b9ad928SBarry Smith   PetscFunctionBegin;
11860700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11873ba16761SJacob Faibussowitsch   if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS);
11889566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0));
11899566063dSJacob Faibussowitsch   PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx));
11909566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0));
11913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11924b9ad928SBarry Smith }
11934b9ad928SBarry Smith 
11944b9ad928SBarry Smith /*@
11954b9ad928SBarry Smith   PCSetOperators - Sets the matrix associated with the linear system and
11964b9ad928SBarry Smith   a (possibly) different one associated with the preconditioner.
11974b9ad928SBarry Smith 
1198c3339decSBarry Smith   Logically Collective
11994b9ad928SBarry Smith 
12004b9ad928SBarry Smith   Input Parameters:
12014b9ad928SBarry Smith + pc   - the preconditioner context
1202e5d3d808SBarry Smith . Amat - the matrix that defines the linear system
1203d1e9a80fSBarry Smith - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
12044b9ad928SBarry Smith 
120520f4b53cSBarry Smith   Level: intermediate
1206189c0b8aSBarry Smith 
120720f4b53cSBarry Smith   Notes:
120820f4b53cSBarry Smith   Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used.
120920f4b53cSBarry Smith 
121020f4b53cSBarry Smith   If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then
1211f1580f4eSBarry Smith   first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()`
1212f1580f4eSBarry Smith   on it and then pass it back in in your call to `KSPSetOperators()`.
1213189c0b8aSBarry Smith 
12144b9ad928SBarry Smith   More Notes about Repeated Solution of Linear Systems:
121520f4b53cSBarry Smith   PETSc does NOT reset the matrix entries of either `Amat` or `Pmat`
12164b9ad928SBarry Smith   to zero after a linear solve; the user is completely responsible for
1217f1580f4eSBarry Smith   matrix assembly.  See the routine `MatZeroEntries()` if desiring to
12184b9ad928SBarry Smith   zero all elements of a matrix.
12194b9ad928SBarry Smith 
1220f1580f4eSBarry Smith .seealso: `PC`, `PCGetOperators()`, `MatZeroEntries()`
12214b9ad928SBarry Smith  @*/
1222d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat)
1223d71ae5a4SJacob Faibussowitsch {
12243b3f6333SBarry Smith   PetscInt m1, n1, m2, n2;
12254b9ad928SBarry Smith 
12264b9ad928SBarry Smith   PetscFunctionBegin;
12270700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12280700a824SBarry Smith   if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
12290700a824SBarry Smith   if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3);
12303fc8bf9cSBarry Smith   if (Amat) PetscCheckSameComm(pc, 1, Amat, 2);
12313fc8bf9cSBarry Smith   if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3);
123231641f1aSBarry Smith   if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) {
12339566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Amat, &m1, &n1));
12349566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->mat, &m2, &n2));
123563a3b9bcSJacob 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);
12369566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Pmat, &m1, &n1));
12379566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2));
123863a3b9bcSJacob 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);
12393b3f6333SBarry Smith   }
12404b9ad928SBarry Smith 
124123ee1639SBarry Smith   if (Pmat != pc->pmat) {
124223ee1639SBarry Smith     /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
124323ee1639SBarry Smith     pc->matnonzerostate = -1;
124423ee1639SBarry Smith     pc->matstate        = -1;
124523ee1639SBarry Smith   }
124623ee1639SBarry Smith 
1247906ed7ccSBarry Smith   /* reference first in case the matrices are the same */
12489566063dSJacob Faibussowitsch   if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat));
12499566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->mat));
12509566063dSJacob Faibussowitsch   if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat));
12519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc->pmat));
12524b9ad928SBarry Smith   pc->mat  = Amat;
12534b9ad928SBarry Smith   pc->pmat = Pmat;
12543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1255e56f5c9eSBarry Smith }
1256e56f5c9eSBarry Smith 
125723ee1639SBarry Smith /*@
125823ee1639SBarry Smith   PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed.
125923ee1639SBarry Smith 
1260c3339decSBarry Smith   Logically Collective
126123ee1639SBarry Smith 
126223ee1639SBarry Smith   Input Parameters:
126323ee1639SBarry Smith + pc   - the preconditioner context
1264f1580f4eSBarry Smith - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
126523ee1639SBarry Smith 
12662b26979fSBarry Smith   Level: intermediate
12672b26979fSBarry Smith 
1268f1580f4eSBarry Smith   Note:
1269f1580f4eSBarry Smith   Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option
1270f1580f4eSBarry Smith   prevents this.
1271f1580f4eSBarry Smith 
1272f1580f4eSBarry Smith .seealso: `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()`
127323ee1639SBarry Smith  @*/
1274d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag)
1275d71ae5a4SJacob Faibussowitsch {
127623ee1639SBarry Smith   PetscFunctionBegin;
127723ee1639SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1278f9177522SStefano Zampini   PetscValidLogicalCollectiveBool(pc, flag, 2);
127923ee1639SBarry Smith   pc->reusepreconditioner = flag;
12803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12814b9ad928SBarry Smith }
12824b9ad928SBarry Smith 
1283c60c7ad4SBarry Smith /*@
1284f1580f4eSBarry Smith   PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed.
1285c60c7ad4SBarry Smith 
1286bba28a21SBarry Smith   Not Collective
1287c60c7ad4SBarry Smith 
1288c60c7ad4SBarry Smith   Input Parameter:
1289c60c7ad4SBarry Smith . pc - the preconditioner context
1290c60c7ad4SBarry Smith 
1291c60c7ad4SBarry Smith   Output Parameter:
1292f1580f4eSBarry Smith . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
1293c60c7ad4SBarry Smith 
1294d0418729SSatish Balay   Level: intermediate
1295d0418729SSatish Balay 
1296f1580f4eSBarry Smith .seealso: `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()`
1297c60c7ad4SBarry Smith  @*/
1298d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag)
1299d71ae5a4SJacob Faibussowitsch {
1300c60c7ad4SBarry Smith   PetscFunctionBegin;
1301c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
13024f572ea9SToby Isaac   PetscAssertPointer(flag, 2);
1303c60c7ad4SBarry Smith   *flag = pc->reusepreconditioner;
13043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1305c60c7ad4SBarry Smith }
1306c60c7ad4SBarry Smith 
1307487a658cSBarry Smith /*@
13084b9ad928SBarry Smith   PCGetOperators - Gets the matrix associated with the linear system and
13094b9ad928SBarry Smith   possibly a different one associated with the preconditioner.
13104b9ad928SBarry Smith 
131120f4b53cSBarry Smith   Not Collective, though parallel mats are returned if pc is parallel
13124b9ad928SBarry Smith 
13134b9ad928SBarry Smith   Input Parameter:
13144b9ad928SBarry Smith . pc - the preconditioner context
13154b9ad928SBarry Smith 
13164b9ad928SBarry Smith   Output Parameters:
1317e5d3d808SBarry Smith + Amat - the matrix defining the linear system
131823ee1639SBarry Smith - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
13194b9ad928SBarry Smith 
13204b9ad928SBarry Smith   Level: intermediate
13214b9ad928SBarry Smith 
1322f1580f4eSBarry Smith   Note:
132395452b02SPatrick Sanan   Does not increase the reference count of the matrices, so you should not destroy them
1324298cc208SBarry Smith 
1325f1580f4eSBarry Smith   Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators
1326f1580f4eSBarry Smith   are created in `PC` and returned to the user. In this case, if both operators
132773950996SBarry Smith   mat and pmat are requested, two DIFFERENT operators will be returned. If
132873950996SBarry Smith   only one is requested both operators in the PC will be the same (i.e. as
1329f1580f4eSBarry Smith   if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats).
133073950996SBarry Smith   The user must set the sizes of the returned matrices and their type etc just
1331f1580f4eSBarry Smith   as if the user created them with `MatCreate()`. For example,
133273950996SBarry Smith 
1333f1580f4eSBarry Smith .vb
1334f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1335f1580f4eSBarry Smith            set size, type, etc of Amat
133673950996SBarry Smith 
1337f1580f4eSBarry Smith          MatCreate(comm,&mat);
1338f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Amat);
1339f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)mat);
1340f1580f4eSBarry Smith            set size, type, etc of Amat
1341f1580f4eSBarry Smith .ve
134273950996SBarry Smith 
134373950996SBarry Smith   and
134473950996SBarry Smith 
1345f1580f4eSBarry Smith .vb
1346f1580f4eSBarry Smith          KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1347f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
134873950996SBarry Smith 
1349f1580f4eSBarry Smith          MatCreate(comm,&Amat);
1350f1580f4eSBarry Smith          MatCreate(comm,&Pmat);
1351f1580f4eSBarry Smith          KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1352f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Amat);
1353f1580f4eSBarry Smith          PetscObjectDereference((PetscObject)Pmat);
1354f1580f4eSBarry Smith            set size, type, etc of Amat and Pmat
1355f1580f4eSBarry Smith .ve
135673950996SBarry Smith 
1357f1580f4eSBarry Smith   The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy
1358b8d709abSRichard Tran Mills   of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely
1359f1580f4eSBarry Smith   managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look
1360f1580f4eSBarry Smith   at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to
1361f1580f4eSBarry Smith   the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP
1362f1580f4eSBarry Smith   you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you).
1363f1580f4eSBarry Smith   Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when
136473950996SBarry Smith   it can be created for you?
136573950996SBarry Smith 
1366f1580f4eSBarry Smith .seealso: `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()`
13674b9ad928SBarry Smith @*/
1368d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat)
1369d71ae5a4SJacob Faibussowitsch {
13704b9ad928SBarry Smith   PetscFunctionBegin;
13710700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1372e5d3d808SBarry Smith   if (Amat) {
137373950996SBarry Smith     if (!pc->mat) {
13749fca8c71SStefano Zampini       if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */
13759a4708feSJed Brown         pc->mat = pc->pmat;
13769566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->mat));
1377e5d3d808SBarry Smith       } else { /* both Amat and Pmat are empty */
13789566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat));
1379e5d3d808SBarry Smith         if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
138073950996SBarry Smith           pc->pmat = pc->mat;
13819566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->pmat));
138273950996SBarry Smith         }
138373950996SBarry Smith       }
13849a4708feSJed Brown     }
1385e5d3d808SBarry Smith     *Amat = pc->mat;
138673950996SBarry Smith   }
1387e5d3d808SBarry Smith   if (Pmat) {
138873950996SBarry Smith     if (!pc->pmat) {
1389e5d3d808SBarry Smith       if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */
13909a4708feSJed Brown         pc->pmat = pc->mat;
13919566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pc->pmat));
13929a4708feSJed Brown       } else {
13939566063dSJacob Faibussowitsch         PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat));
1394e5d3d808SBarry Smith         if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
139573950996SBarry Smith           pc->mat = pc->pmat;
13969566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)pc->mat));
139773950996SBarry Smith         }
139873950996SBarry Smith       }
13999a4708feSJed Brown     }
1400e5d3d808SBarry Smith     *Pmat = pc->pmat;
140173950996SBarry Smith   }
14023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14034b9ad928SBarry Smith }
14044b9ad928SBarry Smith 
1405906ed7ccSBarry Smith /*@C
1406906ed7ccSBarry Smith   PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1407f1580f4eSBarry Smith   possibly a different one associated with the preconditioner have been set in the `PC`.
1408906ed7ccSBarry Smith 
140920f4b53cSBarry Smith   Not Collective, though the results on all processes should be the same
1410906ed7ccSBarry Smith 
1411906ed7ccSBarry Smith   Input Parameter:
1412906ed7ccSBarry Smith . pc - the preconditioner context
1413906ed7ccSBarry Smith 
1414906ed7ccSBarry Smith   Output Parameters:
1415906ed7ccSBarry Smith + mat  - the matrix associated with the linear system was set
1416906ed7ccSBarry Smith - pmat - matrix associated with the preconditioner was set, usually the same
1417906ed7ccSBarry Smith 
1418906ed7ccSBarry Smith   Level: intermediate
1419906ed7ccSBarry Smith 
1420f1580f4eSBarry Smith .seealso: `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`
1421906ed7ccSBarry Smith @*/
1422d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat)
1423d71ae5a4SJacob Faibussowitsch {
1424906ed7ccSBarry Smith   PetscFunctionBegin;
14250700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1426906ed7ccSBarry Smith   if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1427906ed7ccSBarry Smith   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
14283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1429906ed7ccSBarry Smith }
1430906ed7ccSBarry Smith 
1431f39d8e23SSatish Balay /*@
1432a4fd02acSBarry Smith   PCFactorGetMatrix - Gets the factored matrix from the
1433f1580f4eSBarry Smith   preconditioner context.  This routine is valid only for the `PCLU`,
1434f1580f4eSBarry Smith   `PCILU`, `PCCHOLESKY`, and `PCICC` methods.
14354b9ad928SBarry Smith 
1436c3339decSBarry Smith   Not Collective though mat is parallel if pc is parallel
14374b9ad928SBarry Smith 
1438f1580f4eSBarry Smith   Input Parameter:
14394b9ad928SBarry Smith . pc - the preconditioner context
14404b9ad928SBarry Smith 
1441feefa0e1SJacob Faibussowitsch   Output Parameters:
14424b9ad928SBarry Smith . mat - the factored matrix
14434b9ad928SBarry Smith 
14444b9ad928SBarry Smith   Level: advanced
14454b9ad928SBarry Smith 
1446f1580f4eSBarry Smith   Note:
144795452b02SPatrick Sanan   Does not increase the reference count for the matrix so DO NOT destroy it
14489405f653SBarry Smith 
1449f1580f4eSBarry Smith .seealso: `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC`
14504b9ad928SBarry Smith @*/
1451d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat)
1452d71ae5a4SJacob Faibussowitsch {
14534b9ad928SBarry Smith   PetscFunctionBegin;
14540700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14554f572ea9SToby Isaac   PetscAssertPointer(mat, 2);
14567def7855SStefano Zampini   PetscCall(PCFactorSetUpMatSolverType(pc));
1457dbbe0bcdSBarry Smith   PetscUseTypeMethod(pc, getfactoredmatrix, mat);
14583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14594b9ad928SBarry Smith }
14604b9ad928SBarry Smith 
14614b9ad928SBarry Smith /*@C
14624b9ad928SBarry Smith   PCSetOptionsPrefix - Sets the prefix used for searching for all
1463f1580f4eSBarry Smith   `PC` options in the database.
14644b9ad928SBarry Smith 
1465c3339decSBarry Smith   Logically Collective
14664b9ad928SBarry Smith 
14674b9ad928SBarry Smith   Input Parameters:
14684b9ad928SBarry Smith + pc     - the preconditioner context
1469f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
14704b9ad928SBarry Smith 
1471f1580f4eSBarry Smith   Note:
14724b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
14734b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
14744b9ad928SBarry Smith   hyphen.
14754b9ad928SBarry Smith 
14764b9ad928SBarry Smith   Level: advanced
14774b9ad928SBarry Smith 
1478f1580f4eSBarry Smith .seealso: `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()`
14794b9ad928SBarry Smith @*/
1480d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[])
1481d71ae5a4SJacob Faibussowitsch {
14824b9ad928SBarry Smith   PetscFunctionBegin;
14830700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14849566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix));
14853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14864b9ad928SBarry Smith }
14874b9ad928SBarry Smith 
14884b9ad928SBarry Smith /*@C
14894b9ad928SBarry Smith   PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1490f1580f4eSBarry Smith   `PC` options in the database.
14914b9ad928SBarry Smith 
1492c3339decSBarry Smith   Logically Collective
14934b9ad928SBarry Smith 
14944b9ad928SBarry Smith   Input Parameters:
14954b9ad928SBarry Smith + pc     - the preconditioner context
1496f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests
14974b9ad928SBarry Smith 
1498f1580f4eSBarry Smith   Note:
14994b9ad928SBarry Smith   A hyphen (-) must NOT be given at the beginning of the prefix name.
15004b9ad928SBarry Smith   The first character of all runtime options is AUTOMATICALLY the
15014b9ad928SBarry Smith   hyphen.
15024b9ad928SBarry Smith 
15034b9ad928SBarry Smith   Level: advanced
15044b9ad928SBarry Smith 
1505f1580f4eSBarry Smith .seealso: `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()`
15064b9ad928SBarry Smith @*/
1507d71ae5a4SJacob Faibussowitsch PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[])
1508d71ae5a4SJacob Faibussowitsch {
15094b9ad928SBarry Smith   PetscFunctionBegin;
15100700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15119566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix));
15123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15134b9ad928SBarry Smith }
15144b9ad928SBarry Smith 
15154b9ad928SBarry Smith /*@C
15164b9ad928SBarry Smith   PCGetOptionsPrefix - Gets the prefix used for searching for all
15174b9ad928SBarry Smith   PC options in the database.
15184b9ad928SBarry Smith 
15194b9ad928SBarry Smith   Not Collective
15204b9ad928SBarry Smith 
1521f1580f4eSBarry Smith   Input Parameter:
15224b9ad928SBarry Smith . pc - the preconditioner context
15234b9ad928SBarry Smith 
1524f1580f4eSBarry Smith   Output Parameter:
15254b9ad928SBarry Smith . prefix - pointer to the prefix string used, is returned
15264b9ad928SBarry Smith 
1527feefa0e1SJacob Faibussowitsch   Fortran Notes:
1528f1580f4eSBarry Smith   The user should pass in a string 'prefix' of
15294b9ad928SBarry Smith   sufficient length to hold the prefix.
15304b9ad928SBarry Smith 
15314b9ad928SBarry Smith   Level: advanced
15324b9ad928SBarry Smith 
1533f1580f4eSBarry Smith .seealso: `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()`
15344b9ad928SBarry Smith @*/
1535d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[])
1536d71ae5a4SJacob Faibussowitsch {
15374b9ad928SBarry Smith   PetscFunctionBegin;
15380700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15394f572ea9SToby Isaac   PetscAssertPointer(prefix, 2);
15409566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix));
15413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15424b9ad928SBarry Smith }
15434b9ad928SBarry Smith 
15448066bbecSBarry Smith /*
15458066bbecSBarry Smith    Indicates the right hand side will be changed by KSPSolve(), this occurs for a few
15468066bbecSBarry Smith   preconditioners including BDDC and Eisentat that transform the equations before applying
15478066bbecSBarry Smith   the Krylov methods
15488066bbecSBarry Smith */
1549d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change)
1550d71ae5a4SJacob Faibussowitsch {
1551a06fd7f2SStefano Zampini   PetscFunctionBegin;
1552a06fd7f2SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15534f572ea9SToby Isaac   PetscAssertPointer(change, 2);
1554a06fd7f2SStefano Zampini   *change = PETSC_FALSE;
1555cac4c232SBarry Smith   PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change));
15563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1557a06fd7f2SStefano Zampini }
1558a06fd7f2SStefano Zampini 
15594b9ad928SBarry Smith /*@
15604b9ad928SBarry Smith   PCPreSolve - Optional pre-solve phase, intended for any
15614b9ad928SBarry Smith   preconditioner-specific actions that must be performed before
15624b9ad928SBarry Smith   the iterative solve itself.
15634b9ad928SBarry Smith 
1564c3339decSBarry Smith   Collective
15654b9ad928SBarry Smith 
15664b9ad928SBarry Smith   Input Parameters:
15674b9ad928SBarry Smith + pc  - the preconditioner context
15684b9ad928SBarry Smith - ksp - the Krylov subspace context
15694b9ad928SBarry Smith 
15704b9ad928SBarry Smith   Level: developer
15714b9ad928SBarry Smith 
1572feefa0e1SJacob Faibussowitsch   Example Usage:
15734b9ad928SBarry Smith .vb
15744b9ad928SBarry Smith     PCPreSolve(pc,ksp);
157523ce1328SBarry Smith     KSPSolve(ksp,b,x);
15764b9ad928SBarry Smith     PCPostSolve(pc,ksp);
15774b9ad928SBarry Smith .ve
15784b9ad928SBarry Smith 
15794b9ad928SBarry Smith   Notes:
1580f1580f4eSBarry Smith   The pre-solve phase is distinct from the `PCSetUp()` phase.
15814b9ad928SBarry Smith 
1582f1580f4eSBarry Smith   `KSPSolve()` calls this directly, so is rarely called by the user.
15834b9ad928SBarry Smith 
1584f1580f4eSBarry Smith .seealso: `PC`, `PCPostSolve()`
15854b9ad928SBarry Smith @*/
1586d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPreSolve(PC pc, KSP ksp)
1587d71ae5a4SJacob Faibussowitsch {
15884b9ad928SBarry Smith   Vec x, rhs;
15894b9ad928SBarry Smith 
15904b9ad928SBarry Smith   PetscFunctionBegin;
15910700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15920700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1593d9a03883SBarry Smith   pc->presolvedone++;
15947827d75bSBarry Smith   PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice");
15959566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
15969566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
15974b9ad928SBarry Smith 
1598dbbe0bcdSBarry Smith   if (pc->ops->presolve) PetscUseTypeMethod(pc, presolve, ksp, rhs, x);
15991baa6e33SBarry Smith   else if (pc->presolve) PetscCall((pc->presolve)(pc, ksp));
16003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16014b9ad928SBarry Smith }
16024b9ad928SBarry Smith 
1603f560b561SHong Zhang /*@C
1604f1580f4eSBarry Smith   PCSetPreSolve - Sets function used by `PCPreSolve()` which is intended for any
1605f560b561SHong Zhang   preconditioner-specific actions that must be performed before
1606f560b561SHong Zhang   the iterative solve itself.
1607f560b561SHong Zhang 
1608c3339decSBarry Smith   Logically Collective
1609f560b561SHong Zhang 
1610f560b561SHong Zhang   Input Parameters:
1611f560b561SHong Zhang + pc       - the preconditioner object
1612f560b561SHong Zhang - presolve - the function to call before the solve
1613f560b561SHong Zhang 
161420f4b53cSBarry Smith   Calling sequence of `presolve`:
161520f4b53cSBarry Smith + pc  - the `PC` context
161620f4b53cSBarry Smith - ksp - the `KSP` context
1617f560b561SHong Zhang 
1618f560b561SHong Zhang   Level: developer
1619f560b561SHong Zhang 
1620db781477SPatrick Sanan .seealso: `PC`, `PCSetUp()`, `PCPreSolve()`
1621f560b561SHong Zhang @*/
162204c3f3b8SBarry Smith PetscErrorCode PCSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC pc, KSP ksp))
1623d71ae5a4SJacob Faibussowitsch {
1624f560b561SHong Zhang   PetscFunctionBegin;
1625f560b561SHong Zhang   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1626f560b561SHong Zhang   pc->presolve = presolve;
16273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1628f560b561SHong Zhang }
1629f560b561SHong Zhang 
16304b9ad928SBarry Smith /*@
16314b9ad928SBarry Smith   PCPostSolve - Optional post-solve phase, intended for any
16324b9ad928SBarry Smith   preconditioner-specific actions that must be performed after
16334b9ad928SBarry Smith   the iterative solve itself.
16344b9ad928SBarry Smith 
1635c3339decSBarry Smith   Collective
16364b9ad928SBarry Smith 
16374b9ad928SBarry Smith   Input Parameters:
16384b9ad928SBarry Smith + pc  - the preconditioner context
16394b9ad928SBarry Smith - ksp - the Krylov subspace context
16404b9ad928SBarry Smith 
1641feefa0e1SJacob Faibussowitsch   Example Usage:
16424b9ad928SBarry Smith .vb
16434b9ad928SBarry Smith     PCPreSolve(pc,ksp);
164423ce1328SBarry Smith     KSPSolve(ksp,b,x);
16454b9ad928SBarry Smith     PCPostSolve(pc,ksp);
16464b9ad928SBarry Smith .ve
16474b9ad928SBarry Smith 
16484b9ad928SBarry Smith   Note:
1649f1580f4eSBarry Smith   `KSPSolve()` calls this routine directly, so it is rarely called by the user.
16504b9ad928SBarry Smith 
16514b9ad928SBarry Smith   Level: developer
16524b9ad928SBarry Smith 
1653f1580f4eSBarry Smith .seealso: `PC`, `PCSetPostSolve()`, `PCSetPresolve()`, `PCPreSolve()`, `KSPSolve()`
16544b9ad928SBarry Smith @*/
1655d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPostSolve(PC pc, KSP ksp)
1656d71ae5a4SJacob Faibussowitsch {
16574b9ad928SBarry Smith   Vec x, rhs;
16584b9ad928SBarry Smith 
16594b9ad928SBarry Smith   PetscFunctionBegin;
16600700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16610700a824SBarry Smith   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1662d9a03883SBarry Smith   pc->presolvedone--;
16639566063dSJacob Faibussowitsch   PetscCall(KSPGetSolution(ksp, &x));
16649566063dSJacob Faibussowitsch   PetscCall(KSPGetRhs(ksp, &rhs));
1665dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, postsolve, ksp, rhs, x);
16663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16674b9ad928SBarry Smith }
16684b9ad928SBarry Smith 
166955849f57SBarry Smith /*@C
1670f1580f4eSBarry Smith   PCLoad - Loads a `PC` that has been stored in binary  with `PCView()`.
167155849f57SBarry Smith 
1672c3339decSBarry Smith   Collective
167355849f57SBarry Smith 
167455849f57SBarry Smith   Input Parameters:
1675f1580f4eSBarry Smith + newdm  - the newly loaded `PC`, this needs to have been created with `PCCreate()` or
1676f1580f4eSBarry Smith            some related function before a call to `PCLoad()`.
1677f1580f4eSBarry Smith - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()`
167855849f57SBarry Smith 
167955849f57SBarry Smith   Level: intermediate
168055849f57SBarry Smith 
1681f1580f4eSBarry Smith   Note:
168255849f57SBarry Smith   The type is determined by the data in the file, any type set into the PC before this call is ignored.
168355849f57SBarry Smith 
1684f1580f4eSBarry Smith .seealso: `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()`
168555849f57SBarry Smith @*/
1686d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLoad(PC newdm, PetscViewer viewer)
1687d71ae5a4SJacob Faibussowitsch {
168855849f57SBarry Smith   PetscBool isbinary;
1689060da220SMatthew G. Knepley   PetscInt  classid;
169055849f57SBarry Smith   char      type[256];
169155849f57SBarry Smith 
169255849f57SBarry Smith   PetscFunctionBegin;
169355849f57SBarry Smith   PetscValidHeaderSpecific(newdm, PC_CLASSID, 1);
169455849f57SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
16959566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
169628b400f6SJacob Faibussowitsch   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
169755849f57SBarry Smith 
16989566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
16997827d75bSBarry Smith   PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file");
17009566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
17019566063dSJacob Faibussowitsch   PetscCall(PCSetType(newdm, type));
1702dbbe0bcdSBarry Smith   PetscTryTypeMethod(newdm, load, viewer);
17033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170455849f57SBarry Smith }
170555849f57SBarry Smith 
17069804daf3SBarry Smith #include <petscdraw.h>
1707e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1708e04113cfSBarry Smith   #include <petscviewersaws.h>
17090acecf5bSBarry Smith #endif
1710fe2efc57SMark 
1711fe2efc57SMark /*@C
1712f1580f4eSBarry Smith   PCViewFromOptions - View from the `PC` based on options in the database
1713fe2efc57SMark 
1714c3339decSBarry Smith   Collective
1715fe2efc57SMark 
1716fe2efc57SMark   Input Parameters:
171720f4b53cSBarry Smith + A    - the `PC` context
1718f1580f4eSBarry Smith . obj  - Optional object that provides the options prefix
1719736c3998SJose E. Roman - name - command line option
1720fe2efc57SMark 
1721fe2efc57SMark   Level: intermediate
1722f1580f4eSBarry Smith 
1723db781477SPatrick Sanan .seealso: `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()`
1724fe2efc57SMark @*/
1725d71ae5a4SJacob Faibussowitsch PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[])
1726d71ae5a4SJacob Faibussowitsch {
1727fe2efc57SMark   PetscFunctionBegin;
1728fe2efc57SMark   PetscValidHeaderSpecific(A, PC_CLASSID, 1);
17299566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
17303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1731fe2efc57SMark }
1732fe2efc57SMark 
17334b9ad928SBarry Smith /*@C
1734f1580f4eSBarry Smith   PCView - Prints information about the `PC`
17354b9ad928SBarry Smith 
1736c3339decSBarry Smith   Collective
17374b9ad928SBarry Smith 
17384b9ad928SBarry Smith   Input Parameters:
1739feefa0e1SJacob Faibussowitsch + pc     - the `PC` context
17404b9ad928SBarry Smith - viewer - optional visualization context
17414b9ad928SBarry Smith 
174220f4b53cSBarry Smith   Level: developer
174320f4b53cSBarry Smith 
1744f1580f4eSBarry Smith   Notes:
17454b9ad928SBarry Smith   The available visualization contexts include
1746f1580f4eSBarry Smith +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
1747f1580f4eSBarry Smith -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
17484b9ad928SBarry Smith   output where only the first processor opens
17494b9ad928SBarry Smith   the file.  All other processors send their
17504b9ad928SBarry Smith   data to the first processor to print.
17514b9ad928SBarry Smith 
17524b9ad928SBarry Smith   The user can open an alternative visualization contexts with
1753f1580f4eSBarry Smith   `PetscViewerASCIIOpen()` (output to a specified file).
17544b9ad928SBarry Smith 
1755f1580f4eSBarry Smith .seealso: `PC`, `PetscViewer`, `KSPView()`, `PetscViewerASCIIOpen()`
17564b9ad928SBarry Smith @*/
1757d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView(PC pc, PetscViewer viewer)
1758d71ae5a4SJacob Faibussowitsch {
175919fd82e9SBarry Smith   PCType            cstr;
1760*6cd81132SPierre Jolivet   PetscViewerFormat format;
1761*6cd81132SPierre Jolivet   PetscBool         iascii, isstring, isbinary, isdraw, pop = PETSC_FALSE;
1762e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1763536b137fSBarry Smith   PetscBool issaws;
17640acecf5bSBarry Smith #endif
17654b9ad928SBarry Smith 
17664b9ad928SBarry Smith   PetscFunctionBegin;
17670700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
176848a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
17690700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1770c9780b6fSBarry Smith   PetscCheckSameComm(pc, 1, viewer, 2);
17714b9ad928SBarry Smith 
17729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
17739566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
17749566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
17759566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1776e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
17779566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
17780acecf5bSBarry Smith #endif
1779219b1045SBarry Smith 
178032077d6dSBarry Smith   if (iascii) {
17819566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer));
178248a46eb9SPierre Jolivet     if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, "  PC has not been set up so information may be incomplete\n"));
17839566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1784dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
17859566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
1786834dbeb0SBarry Smith     if (pc->mat) {
1787*6cd81132SPierre Jolivet       PetscCall(PetscViewerGetFormat(viewer, &format));
1788*6cd81132SPierre Jolivet       if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
17899566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
1790*6cd81132SPierre Jolivet         pop = PETSC_TRUE;
1791*6cd81132SPierre Jolivet       }
17924b9ad928SBarry Smith       if (pc->pmat == pc->mat) {
17939566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix = precond matrix:\n"));
17949566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
17959566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
17969566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
17974b9ad928SBarry Smith       } else {
1798834dbeb0SBarry Smith         if (pc->pmat) {
17999566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix followed by preconditioner matrix:\n"));
18004b9ad928SBarry Smith         } else {
18019566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  linear system matrix:\n"));
18024b9ad928SBarry Smith         }
18039566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
18049566063dSJacob Faibussowitsch         PetscCall(MatView(pc->mat, viewer));
18059566063dSJacob Faibussowitsch         if (pc->pmat) PetscCall(MatView(pc->pmat, viewer));
18069566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
18074b9ad928SBarry Smith       }
1808*6cd81132SPierre Jolivet       if (pop) PetscCall(PetscViewerPopFormat(viewer));
18094b9ad928SBarry Smith     }
18104b9ad928SBarry Smith   } else if (isstring) {
18119566063dSJacob Faibussowitsch     PetscCall(PCGetType(pc, &cstr));
18129566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr));
1813dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18149566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18159566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18165b0b0462SBarry Smith   } else if (isbinary) {
181755849f57SBarry Smith     PetscInt    classid = PC_FILE_CLASSID;
181855849f57SBarry Smith     MPI_Comm    comm;
181955849f57SBarry Smith     PetscMPIInt rank;
182055849f57SBarry Smith     char        type[256];
182155849f57SBarry Smith 
18229566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
18239566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1824dd400576SPatrick Sanan     if (rank == 0) {
18259566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
18269566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256));
18279566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
182855849f57SBarry Smith     }
1829dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
1830219b1045SBarry Smith   } else if (isdraw) {
1831219b1045SBarry Smith     PetscDraw draw;
1832d9884438SBarry Smith     char      str[25];
183389fd9fafSBarry Smith     PetscReal x, y, bottom, h;
1834d9884438SBarry Smith     PetscInt  n;
1835219b1045SBarry Smith 
18369566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
18379566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
18381d840656SPeter Brune     if (pc->mat) {
18399566063dSJacob Faibussowitsch       PetscCall(MatGetSize(pc->mat, &n, NULL));
184063a3b9bcSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n));
18411d840656SPeter Brune     } else {
18429566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name));
18431d840656SPeter Brune     }
18449566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
184589fd9fafSBarry Smith     bottom = y - h;
18469566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
1847dbbe0bcdSBarry Smith     PetscTryTypeMethod(pc, view, viewer);
18489566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
1849e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS)
1850536b137fSBarry Smith   } else if (issaws) {
1851d45a07a7SBarry Smith     PetscMPIInt rank;
1852d45a07a7SBarry Smith 
18539566063dSJacob Faibussowitsch     PetscCall(PetscObjectName((PetscObject)pc));
18549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
185548a46eb9SPierre Jolivet     if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer));
18569566063dSJacob Faibussowitsch     if (pc->mat) PetscCall(MatView(pc->mat, viewer));
18579566063dSJacob Faibussowitsch     if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
18580acecf5bSBarry Smith #endif
18594b9ad928SBarry Smith   }
18603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18614b9ad928SBarry Smith }
18624b9ad928SBarry Smith 
18634b9ad928SBarry Smith /*@C
18641c84c290SBarry Smith   PCRegister -  Adds a method to the preconditioner package.
18651c84c290SBarry Smith 
18661c84c290SBarry Smith   Not collective
18671c84c290SBarry Smith 
18681c84c290SBarry Smith   Input Parameters:
186920f4b53cSBarry Smith + sname    - name of a new user-defined solver
187020f4b53cSBarry Smith - function - routine to create method context
18711c84c290SBarry Smith 
1872feefa0e1SJacob Faibussowitsch   Example Usage:
18731c84c290SBarry Smith .vb
1874bdf89e91SBarry Smith    PCRegister("my_solver", MySolverCreate);
18751c84c290SBarry Smith .ve
18761c84c290SBarry Smith 
18771c84c290SBarry Smith   Then, your solver can be chosen with the procedural interface via
18781c84c290SBarry Smith $     PCSetType(pc, "my_solver")
18791c84c290SBarry Smith   or at runtime via the option
18801c84c290SBarry Smith $     -pc_type my_solver
18814b9ad928SBarry Smith 
18824b9ad928SBarry Smith   Level: advanced
18831c84c290SBarry Smith 
188420f4b53cSBarry Smith   Note:
188520f4b53cSBarry Smith   `PCRegister()` may be called multiple times to add several user-defined preconditioners.
188620f4b53cSBarry Smith 
1887f1580f4eSBarry Smith .seealso: `PC`, `PCType`, `PCRegisterAll()`
18884b9ad928SBarry Smith @*/
1889d71ae5a4SJacob Faibussowitsch PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC))
1890d71ae5a4SJacob Faibussowitsch {
18914b9ad928SBarry Smith   PetscFunctionBegin;
18929566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
18939566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCList, sname, function));
18943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18954b9ad928SBarry Smith }
18964b9ad928SBarry Smith 
1897d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y)
1898d71ae5a4SJacob Faibussowitsch {
1899186a3e20SStefano Zampini   PC pc;
1900186a3e20SStefano Zampini 
1901186a3e20SStefano Zampini   PetscFunctionBegin;
19029566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &pc));
19039566063dSJacob Faibussowitsch   PetscCall(PCApply(pc, X, Y));
19043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1905186a3e20SStefano Zampini }
1906186a3e20SStefano Zampini 
1907feefa0e1SJacob Faibussowitsch /*@C
19080bacdadaSStefano Zampini   PCComputeOperator - Computes the explicit preconditioned operator.
19094b9ad928SBarry Smith 
1910c3339decSBarry Smith   Collective
19114b9ad928SBarry Smith 
1912d8d19677SJose E. Roman   Input Parameters:
1913186a3e20SStefano Zampini + pc      - the preconditioner object
1914186a3e20SStefano Zampini - mattype - the matrix type to be used for the operator
19154b9ad928SBarry Smith 
19164b9ad928SBarry Smith   Output Parameter:
1917a5b23f4aSJose E. Roman . mat - the explicit preconditioned operator
19184b9ad928SBarry Smith 
191920f4b53cSBarry Smith   Level: advanced
192020f4b53cSBarry Smith 
1921f1580f4eSBarry Smith   Note:
1922186a3e20SStefano Zampini   This computation is done by applying the operators to columns of the identity matrix.
1923186a3e20SStefano Zampini   This routine is costly in general, and is recommended for use only with relatively small systems.
1924186a3e20SStefano Zampini   Currently, this routine uses a dense matrix format when mattype == NULL
19254b9ad928SBarry Smith 
1926f1580f4eSBarry Smith .seealso: `PC`, `KSPComputeOperator()`, `MatType`
19274b9ad928SBarry Smith @*/
1928d71ae5a4SJacob Faibussowitsch PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat)
1929d71ae5a4SJacob Faibussowitsch {
1930186a3e20SStefano Zampini   PetscInt N, M, m, n;
1931186a3e20SStefano Zampini   Mat      A, Apc;
19324b9ad928SBarry Smith 
19334b9ad928SBarry Smith   PetscFunctionBegin;
19340700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
19354f572ea9SToby Isaac   PetscAssertPointer(mat, 3);
19369566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, &A, NULL));
19379566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
19389566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
19399566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc));
19409566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (void (*)(void))MatMult_PC));
19419566063dSJacob Faibussowitsch   PetscCall(MatComputeOperator(Apc, mattype, mat));
19429566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Apc));
19433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19444b9ad928SBarry Smith }
19454b9ad928SBarry Smith 
19466c237d78SBarry Smith /*@
19476c237d78SBarry Smith   PCSetCoordinates - sets the coordinates of all the nodes on the local process
19486c237d78SBarry Smith 
1949c3339decSBarry Smith   Collective
19506c237d78SBarry Smith 
19516c237d78SBarry Smith   Input Parameters:
19526c237d78SBarry Smith + pc     - the solver context
19536c237d78SBarry Smith . dim    - the dimension of the coordinates 1, 2, or 3
195414893cbeSStefano Zampini . nloc   - the blocked size of the coordinates array
195514893cbeSStefano Zampini - coords - the coordinates array
19566c237d78SBarry Smith 
19576c237d78SBarry Smith   Level: intermediate
19586c237d78SBarry Smith 
1959f1580f4eSBarry Smith   Note:
196020f4b53cSBarry Smith   `coords` is an array of the dim coordinates for the nodes on
196120f4b53cSBarry Smith   the local processor, of size `dim`*`nloc`.
196214893cbeSStefano Zampini   If there are 108 equation on a processor
19636c237d78SBarry Smith   for a displacement finite element discretization of elasticity (so
196414893cbeSStefano Zampini   that there are nloc = 36 = 108/3 nodes) then the array must have 108
19656c237d78SBarry Smith   double precision values (ie, 3 * 36).  These x y z coordinates
19666c237d78SBarry Smith   should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
19676c237d78SBarry Smith   ... , N-1.z ].
19686c237d78SBarry Smith 
1969f1580f4eSBarry Smith .seealso: `PC`, `MatSetNearNullSpace()`
19706c237d78SBarry Smith @*/
1971d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[])
1972d71ae5a4SJacob Faibussowitsch {
19736c237d78SBarry Smith   PetscFunctionBegin;
197432954da3SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
197532954da3SStefano Zampini   PetscValidLogicalCollectiveInt(pc, dim, 2);
197622794d57SStefano Zampini   PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords));
19773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19786c237d78SBarry Smith }
1979fd2dd295SFande Kong 
1980fd2dd295SFande Kong /*@
1981fd2dd295SFande Kong   PCGetInterpolations - Gets interpolation matrices for all levels (except level 0)
1982fd2dd295SFande Kong 
1983c3339decSBarry Smith   Logically Collective
1984fd2dd295SFande Kong 
1985d8d19677SJose E. Roman   Input Parameter:
1986d8d19677SJose E. Roman . pc - the precondition context
1987fd2dd295SFande Kong 
1988d8d19677SJose E. Roman   Output Parameters:
1989d8d19677SJose E. Roman + num_levels     - the number of levels
1990d8d19677SJose E. Roman - interpolations - the interpolation matrices (size of num_levels-1)
1991fd2dd295SFande Kong 
1992fd2dd295SFande Kong   Level: advanced
1993fd2dd295SFande Kong 
1994feefa0e1SJacob Faibussowitsch   Developer Notes:
1995f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
1996fd2dd295SFande Kong 
1997f1580f4eSBarry Smith .seealso: `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()`
1998fd2dd295SFande Kong @*/
1999d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[])
2000d71ae5a4SJacob Faibussowitsch {
2001fd2dd295SFande Kong   PetscFunctionBegin;
2002fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20034f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20044f572ea9SToby Isaac   PetscAssertPointer(interpolations, 3);
2005cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations));
20063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2007fd2dd295SFande Kong }
2008fd2dd295SFande Kong 
2009fd2dd295SFande Kong /*@
2010fd2dd295SFande Kong   PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level)
2011fd2dd295SFande Kong 
2012c3339decSBarry Smith   Logically Collective
2013fd2dd295SFande Kong 
2014d8d19677SJose E. Roman   Input Parameter:
2015d8d19677SJose E. Roman . pc - the precondition context
2016fd2dd295SFande Kong 
2017d8d19677SJose E. Roman   Output Parameters:
2018d8d19677SJose E. Roman + num_levels      - the number of levels
2019d8d19677SJose E. Roman - coarseOperators - the coarse operator matrices (size of num_levels-1)
2020fd2dd295SFande Kong 
2021fd2dd295SFande Kong   Level: advanced
2022fd2dd295SFande Kong 
2023feefa0e1SJacob Faibussowitsch   Developer Notes:
2024f1580f4eSBarry Smith   Why is this here instead of in `PCMG` etc?
2025fd2dd295SFande Kong 
2026f1580f4eSBarry Smith .seealso: `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()`
2027fd2dd295SFande Kong @*/
2028d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
2029d71ae5a4SJacob Faibussowitsch {
2030fd2dd295SFande Kong   PetscFunctionBegin;
2031fd2dd295SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
20324f572ea9SToby Isaac   PetscAssertPointer(num_levels, 2);
20334f572ea9SToby Isaac   PetscAssertPointer(coarseOperators, 3);
2034cac4c232SBarry Smith   PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators));
20353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2036fd2dd295SFande Kong }
2037