xref: /petsc/include/petscksp.h (revision fad47a0a50989ce39c782fb8acaab0e4028b3de9)
1f26ada1bSBarry Smith /*
2f26ada1bSBarry Smith    Defines the interface functions for the Krylov subspace accelerators.
3f26ada1bSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCKSP_H
50a835dfdSSatish Balay #define __PETSCKSP_H
62c8e378dSBarry Smith #include <petscpc.h>
72eac72dbSBarry Smith 
8607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void);
91dbb0a54SBarry Smith 
1028ce4d24SBarry Smith /*S
118f6c3df8SBarry Smith      KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the
128f6c3df8SBarry Smith          linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators).
1328ce4d24SBarry Smith 
1428ce4d24SBarry Smith    Level: beginner
1528ce4d24SBarry Smith 
1628ce4d24SBarry Smith   Concepts: Krylov methods
1728ce4d24SBarry Smith 
188f6c3df8SBarry Smith         Notes: When a direct solver is used but no Krylov solver is used the KSP object is still used by with a
198f6c3df8SBarry Smith        KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver).
208f6c3df8SBarry Smith 
218f6c3df8SBarry Smith .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy()
2228ce4d24SBarry Smith S*/
23e2a1c21fSSatish Balay typedef struct _p_KSP*     KSP;
242eac72dbSBarry Smith 
2576bdecfbSBarry Smith /*J
268f6c3df8SBarry Smith     KSPType - String with the name of a PETSc Krylov method.
2728ce4d24SBarry Smith 
2828ce4d24SBarry Smith    Level: beginner
2928ce4d24SBarry Smith 
308f6c3df8SBarry Smith .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions()
3176bdecfbSBarry Smith J*/
3219fd82e9SBarry Smith typedef const char* KSPType;
3382bf6240SBarry Smith #define KSPRICHARDSON "richardson"
346c9de887SHong Zhang #define KSPCHEBYSHEV  "chebyshev"
3582bf6240SBarry Smith #define KSPCG         "cg"
362c8fc646SJed Brown #define KSPGROPPCG    "groppcg"
372c8fc646SJed Brown #define KSPPIPECG     "pipecg"
38df2a969aSvictorle #define   KSPCGNE       "cgne"
39fcae7a14Stmunson #define   KSPNASH       "nash"
4080e17431SMatthew Knepley #define   KSPSTCG       "stcg"
418031f4adStmunson #define   KSPGLTR       "gltr"
42640f4f53SPatrick Sanan #define KSPFCG        "fcg"
43390d8e47SPatrick Sanan #define KSPPIPEFCG    "pipefcg"
4482bf6240SBarry Smith #define KSPGMRES      "gmres"
459dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
469dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
474f8e6cd9SSatish Balay #define   KSPDGMRES     "dgmres"
4861b468f9SJed Brown #define   KSPPGMRES     "pgmres"
4982bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
5082bf6240SBarry Smith #define KSPBCGS       "bcgs"
51e1c61ce8SBarry Smith #define   KSPIBCGS      "ibcgs"
5218ac38e6SHong Zhang #define   KSPFBCGS      "fbcgs"
53c2b71004SHong Zhang #define   KSPFBCGSR     "fbcgsr"
54f0037002Svictorle #define   KSPBCGSL      "bcgsl"
5582bf6240SBarry Smith #define KSPCGS        "cgs"
5682bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
5782bf6240SBarry Smith #define KSPCR         "cr"
582c8fc646SJed Brown #define KSPPIPECR     "pipecr"
5982bf6240SBarry Smith #define KSPLSQR       "lsqr"
6082bf6240SBarry Smith #define KSPPREONLY    "preonly"
6182bf6240SBarry Smith #define KSPQCG        "qcg"
62c9cf9b11SBarry Smith #define KSPBICG       "bicg"
63b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
6401c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
6521356919SSatish Balay #define KSPLCD        "lcd"
661d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
6758ad63f4SBarry Smith #define KSPGCR        "gcr"
68*fad47a0aSPatrick Sanan #define KSPPIPEGCR    "pipegcr"
69e4d80e07Szianekhodja #define KSPTSIRM      "tsirm"
70e4d80e07Szianekhodja #define KSPCGLS       "cgls"
712eac72dbSBarry Smith 
728ba1e511SMatthew Knepley /* Logging support */
73014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
745358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
758ba1e511SMatthew Knepley 
76014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
7719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
78c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
79014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
80014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
81014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
82014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
83014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
84014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
8523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
8625c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
872eac72dbSBarry Smith 
88140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
89bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
9030de9b25SBarry Smith 
91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
93014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
94c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool );
96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool  *);
97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool );
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *);
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool );
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool  *);
101014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool );
1027d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool );
103c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *);
104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool );
105c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *);
106014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
110734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
1112a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1122a7a6963SBarry Smith PETSC_DEPRECATED("Use KSPCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);}
1132eac72dbSBarry Smith 
114d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
115d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
116d4211eb9SBarry Smith 
117014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
119aabeff55SBarry Smith 
120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
124014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
125014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
1264b0e389bSBarry Smith 
127fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
128fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *);
129fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
130fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
131fa0ddf94SBarry Smith 
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
137b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
138b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
139b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
140b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *);
1420a71e23eSMatthew Knepley 
143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
1452eac72dbSBarry Smith 
146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
14958450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
15082fe98a5SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseRandom(KSP,PetscBool);
151ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom);
15258450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
153014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
154d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *);
155d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
1567d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt *,Vec[],PetscReal[],PetscReal[]);
1574b0e389bSBarry Smith 
158640f4f53SPatrick Sanan /*E
159640f4f53SPatrick Sanan 
16006137d0aSPatrick Sanan   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
161640f4f53SPatrick Sanan 
16206137d0aSPatrick Sanan   KSP_FCD_TRUNCATION uses all (up to mmax) stored directions
16306137d0aSPatrick Sanan   KSP_FCD_TRUNCATION_RESTART uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
164640f4f53SPatrick Sanan 
1652b26979fSBarry Smith    Level: intermediate
16606137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
167640f4f53SPatrick Sanan 
168640f4f53SPatrick Sanan E*/
16906137d0aSPatrick Sanan typedef enum {KSP_FCD_TRUNCATION,KSP_FCD_TRUNCATION_RESTART} KSPFCDTruncationType;
17006137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
17106137d0aSPatrick Sanan 
17206137d0aSPatrick Sanan #define KSPFCDGetNumOldDirections(ctx,i,mi) ({                                                    \
17306137d0aSPatrick Sanan   if(ctx->trunctype == KSP_FCD_TRUNCATION_RESTART){                                               \
17406137d0aSPatrick Sanan     mi = ((i-1) % ctx->mmax)+1;                                                                   \
17506137d0aSPatrick Sanan     if (mi==1 && i!=1) ++(ctx->n_search_space_resets);                                            \
17606137d0aSPatrick Sanan   } else if (ctx->trunctype == KSP_FCD_TRUNCATION)                                                \
17706137d0aSPatrick Sanan     mi = ctx->mmax;                                                                               \
17806137d0aSPatrick Sanan  else {                                                                                           \
17906137d0aSPatrick Sanan    SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unrecognized Truncation Strategy");CHKERRQ(ierr); \
18006137d0aSPatrick Sanan   }                                                                                               \
18106137d0aSPatrick Sanan })
182640f4f53SPatrick Sanan 
183640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
184640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
185640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
186640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
18706137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
18806137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
189640f4f53SPatrick Sanan 
190390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
191390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
192390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
193390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
194390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
195390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
196390d8e47SPatrick Sanan 
197*fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
198*fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
199*fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
200*fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
201*fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
202*fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
203*fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
204*fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
205*fad47a0aSPatrick Sanan 
206014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
207014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
208014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
2099f236934SBarry Smith 
210014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
211014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
214014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
2151d73ed98SBarry Smith 
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
2181d73ed98SBarry Smith 
219014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
220014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
22258ad63f4SBarry Smith 
223b49fd9e1SBarry Smith /*E
224b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
225b49fd9e1SBarry Smith 
226b49fd9e1SBarry Smith    Level: advanced
227b49fd9e1SBarry Smith 
228bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
229e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
230b49fd9e1SBarry Smith 
231b49fd9e1SBarry Smith E*/
23278d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2336a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2341f7e983dSSatish Balay /*MC
2351957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
2368c5b8ba0SBarry Smith 
2378c5b8ba0SBarry Smith    Level: advanced
2388c5b8ba0SBarry Smith 
2398c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2408c5b8ba0SBarry Smith 
241bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
242e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2438c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2448c5b8ba0SBarry Smith M*/
2458c5b8ba0SBarry Smith 
2461f7e983dSSatish Balay /*MC
2478c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2488c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2498c5b8ba0SBarry Smith           poor orthogonality.
2508c5b8ba0SBarry Smith 
2518c5b8ba0SBarry Smith    Level: advanced
2528c5b8ba0SBarry Smith 
2538c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2548c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2558c5b8ba0SBarry Smith 
256bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
257e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2588c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2598c5b8ba0SBarry Smith M*/
2608c5b8ba0SBarry Smith 
2611f7e983dSSatish Balay /*MC
2628c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2638c5b8ba0SBarry Smith 
2648c5b8ba0SBarry Smith    Level: advanced
2658c5b8ba0SBarry Smith 
2668c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2678c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2688c5b8ba0SBarry Smith 
2698c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2708c5b8ba0SBarry Smith 
271bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
272e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2738c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2748c5b8ba0SBarry Smith M*/
2758c5b8ba0SBarry Smith 
276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
27808480c60SBarry Smith 
279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
282c38d4ed2SBarry Smith 
283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
286121fd945SKris Buschelman 
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
290e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
291d9492815SBarry Smith 
292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2942eac72dbSBarry Smith 
295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
297014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
299390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
300390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
301014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
302816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *);
303014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
304014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
305e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
306e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
307e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
308014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
30984cb2905SBarry Smith 
310014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
311014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
312c01c455dSBarry Smith 
31323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
31423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
315014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
316014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
318014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
3197287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
3207287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
3211eb62cbbSBarry Smith 
322014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
323014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
324014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
325014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
3261f7f0c4fSBarry Smith 
327014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
32855849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
329685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
3302a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
3312a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
33255849f57SBarry Smith 
33355849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
3341eb62cbbSBarry Smith 
335014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
336014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
337db9b2ab1SHong Zhang 
338014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
339014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
34068ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
34183ab6a24SBarry Smith 
34228ce4d24SBarry Smith /*E
3438a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
3448a4b9c5cSBarry Smith        test routines.
3458a4b9c5cSBarry Smith 
3468a4b9c5cSBarry Smith    Level: advanced
3478a4b9c5cSBarry Smith 
348a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
3491718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
350a3f661c8SBarry Smith 
351af0996ceSBarry Smith    Notes: this must match petsc/finclude/petscksp.h
3528a4b9c5cSBarry Smith 
35394b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
3541718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
3558a4b9c5cSBarry Smith E*/
3569e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
3579e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
358014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
359a21b2a99SBarry Smith 
3601f7e983dSSatish Balay /*MC
3619793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
3628c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
3638c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3648c5b8ba0SBarry Smith 
3658c5b8ba0SBarry Smith    Level: advanced
3668c5b8ba0SBarry Smith 
3671957e957SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
3688c5b8ba0SBarry Smith 
369ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3708c5b8ba0SBarry Smith M*/
3718c5b8ba0SBarry Smith 
3721f7e983dSSatish Balay /*MC
3731957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
3748c5b8ba0SBarry Smith        convergence test routine.
3758c5b8ba0SBarry Smith 
3768c5b8ba0SBarry Smith    Level: advanced
3778c5b8ba0SBarry Smith 
3789793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3798c5b8ba0SBarry Smith M*/
3808c5b8ba0SBarry Smith 
3811f7e983dSSatish Balay /*MC
382ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3838c5b8ba0SBarry Smith        convergence test routine.
3848c5b8ba0SBarry Smith 
3858c5b8ba0SBarry Smith    Level: advanced
3868c5b8ba0SBarry Smith 
3879793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3888c5b8ba0SBarry Smith M*/
3898c5b8ba0SBarry Smith 
3901f7e983dSSatish Balay /*MC
391ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
392390d8e47SPatrick Sanan        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
3938c5b8ba0SBarry Smith 
3948c5b8ba0SBarry Smith    Level: advanced
3958c5b8ba0SBarry Smith 
3969793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3978c5b8ba0SBarry Smith M*/
3988c5b8ba0SBarry Smith 
399014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
400014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
401014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
402014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
4048a4b9c5cSBarry Smith 
4058a4b9c5cSBarry Smith /*E
4061957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
40728ce4d24SBarry Smith 
40828ce4d24SBarry Smith    Level: beginner
40928ce4d24SBarry Smith 
4104d0a8057SBarry Smith    Notes: See KSPGetConvergedReason() for explanation of each value
41128ce4d24SBarry Smith 
412af0996ceSBarry Smith    Developer notes: this must match petsc/finclude/petscksp.h
4134d0a8057SBarry Smith 
4144d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
4154d0a8057SBarry Smith       any of the values here also change them that array of names.
41686c02ca4SBarry Smith 
417c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
41828ce4d24SBarry Smith E*/
419d15094e1SBarry Smith typedef enum {/* converged */
4209ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
4219ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
422d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
423d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
424b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
4258031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
4268031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
427329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
428af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
429d15094e1SBarry Smith               /* diverged */
430b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
431d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
432d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
433d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
434b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
435b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
436b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
4374d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
4386aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
439422a814eSBarry Smith               KSP_DIVERGED_PCSETUP_FAILED      = -11,
440d15094e1SBarry Smith 
441d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
442014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
443d15094e1SBarry Smith 
444c838673bSBarry Smith /*MC
445c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
446c838673bSBarry Smith 
447c838673bSBarry Smith    Level: beginner
448c838673bSBarry Smith 
449c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
450c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
451c838673bSBarry Smith        2-norm of the residual for right preconditioning
452c838673bSBarry Smith 
453c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
454c838673bSBarry Smith 
455c838673bSBarry Smith M*/
456c838673bSBarry Smith 
457c838673bSBarry Smith /*MC
458c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
459c838673bSBarry Smith 
460c838673bSBarry Smith    Level: beginner
461c838673bSBarry Smith 
462c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
463c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
464c838673bSBarry Smith        2-norm of the residual for right preconditioning
465c838673bSBarry Smith 
466c838673bSBarry Smith    Level: beginner
467c838673bSBarry Smith 
468c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
469c838673bSBarry Smith 
470c838673bSBarry Smith M*/
471c838673bSBarry Smith 
472c838673bSBarry Smith /*MC
473c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
474c838673bSBarry Smith 
475c838673bSBarry Smith    Level: beginner
476c838673bSBarry Smith 
477c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
478c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
479c838673bSBarry Smith        2-norm of the residual for right preconditioning
480c838673bSBarry Smith 
481c838673bSBarry Smith    Level: beginner
482c838673bSBarry Smith 
483c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
484c838673bSBarry Smith 
485c838673bSBarry Smith M*/
486c838673bSBarry Smith 
487c838673bSBarry Smith /*MC
488c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
489c838673bSBarry Smith       reached
490c838673bSBarry Smith 
491c838673bSBarry Smith    Level: beginner
492c838673bSBarry Smith 
493c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
494c838673bSBarry Smith 
495c838673bSBarry Smith M*/
496c838673bSBarry Smith 
497c838673bSBarry Smith /*MC
4988c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
4990059c7bdSJed Brown            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
5004d0a8057SBarry Smith            test routine is set in KSP.
501c838673bSBarry Smith 
502c838673bSBarry Smith    Level: beginner
503c838673bSBarry Smith 
504c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
505c838673bSBarry Smith 
506c838673bSBarry Smith M*/
507c838673bSBarry Smith 
508c838673bSBarry Smith /*MC
509c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
5103014e516SBarry Smith           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
5113014e516SBarry Smith           preconditioner.
512c838673bSBarry Smith 
513c838673bSBarry Smith    Level: beginner
514c838673bSBarry Smith 
515c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
516c838673bSBarry Smith 
517c838673bSBarry Smith M*/
518c838673bSBarry Smith 
519c838673bSBarry Smith /*MC
520c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
521c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
522c838673bSBarry Smith 
523c838673bSBarry Smith    Level: beginner
524c838673bSBarry Smith 
525c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
526c838673bSBarry Smith 
527c838673bSBarry Smith M*/
528c838673bSBarry Smith 
529c838673bSBarry Smith /*MC
530c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
531c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
532c838673bSBarry Smith 
533c838673bSBarry Smith    Level: beginner
534c838673bSBarry Smith 
535c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
536c838673bSBarry Smith 
537c838673bSBarry Smith M*/
538c838673bSBarry Smith 
539c838673bSBarry Smith /*MC
540c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
541c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
542c838673bSBarry Smith         be positive definite
543c838673bSBarry Smith 
544c838673bSBarry Smith    Level: beginner
545c838673bSBarry Smith 
5462401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
547c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
548c838673bSBarry Smith 
549c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
550c838673bSBarry Smith 
551c838673bSBarry Smith M*/
552c838673bSBarry Smith 
553c838673bSBarry Smith /*MC
554c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
555c838673bSBarry Smith         while the KSPSolve() is still running.
556c838673bSBarry Smith 
557c838673bSBarry Smith    Level: beginner
558c838673bSBarry Smith 
559c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
560c838673bSBarry Smith 
561c838673bSBarry Smith M*/
562c838673bSBarry Smith 
563014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
564014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
5658de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
566014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
5678de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
5688de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
5698de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
5708de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
5710059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
572014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
573abef13c0SSatish Balay 
5748ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
5758ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
5768ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
5778ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
5788ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
5798ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
5808ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
5818ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
5828ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
5838ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
5848ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
5858ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
5868ea1b3e6SJed Brown 
587014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
588d4fbbf0eSBarry Smith 
58928ce4d24SBarry Smith /*E
59028ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
59128ce4d24SBarry Smith 
59228ce4d24SBarry Smith    Level: beginner
59328ce4d24SBarry Smith 
59428ce4d24SBarry Smith .seealso: KSPCGSetType()
59528ce4d24SBarry Smith E*/
5969dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5976a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
59828ce4d24SBarry Smith 
599014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
600014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
6018031f4adStmunson 
602014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
603014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
604014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
605fcae7a14Stmunson 
606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
609e559a7feSLois Curfman McInnes 
610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
611014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
612014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
613014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
614014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
6158031f4adStmunson 
616014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
6171d6018f0SLisandro Dalcin 
618014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
619014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
6203369ce9aSBarry Smith 
6219804daf3SBarry Smith #include <petscdrawtypes.h>
622d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
623d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
624d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
625d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
626014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
6272f2e5d10SKris Buschelman 
628014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
629014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
63003919abeSBarry Smith 
631f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
632ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
633f8a50e2bSBarry Smith 
634014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
635014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
636014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
637014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
638014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
639014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
640f8a50e2bSBarry Smith 
641014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
642014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
643014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
644f8a50e2bSBarry Smith 
645470b340bSDmitry Karpeev /*E
646470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
647470b340bSDmitry Karpeev 
648470b340bSDmitry Karpeev     Level: intermediate
649470b340bSDmitry Karpeev 
650470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
651470b340bSDmitry Karpeev E*/
652470b340bSDmitry Karpeev typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType;
653470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
654470b340bSDmitry Karpeev 
655014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
656014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
657d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
658bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
659aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
660bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
661470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
662470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
6635bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
6645a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
665470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *);
666470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
6673f22127dSBarry Smith 
668014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
669014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
670014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
671014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
672014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
673fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
67423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
675fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
67623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
67723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
678014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
679014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
680fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
681fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
6826c699258SBarry Smith 
68302b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
68402b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexProjectField(DM, Vec, void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode, Vec);
6852eac72dbSBarry Smith #endif
686