xref: /petsc/include/petscksp.h (revision 483d6965044e70e2c4df7306809d831344a09676)
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"
45*483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres"
469dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
479dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
484f8e6cd9SSatish Balay #define   KSPDGMRES     "dgmres"
4961b468f9SJed Brown #define   KSPPGMRES     "pgmres"
5082bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
5182bf6240SBarry Smith #define KSPBCGS       "bcgs"
52e1c61ce8SBarry Smith #define   KSPIBCGS      "ibcgs"
5318ac38e6SHong Zhang #define   KSPFBCGS      "fbcgs"
54c2b71004SHong Zhang #define   KSPFBCGSR     "fbcgsr"
55f0037002Svictorle #define   KSPBCGSL      "bcgsl"
5682bf6240SBarry Smith #define KSPCGS        "cgs"
5782bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
5882bf6240SBarry Smith #define KSPCR         "cr"
592c8fc646SJed Brown #define KSPPIPECR     "pipecr"
6082bf6240SBarry Smith #define KSPLSQR       "lsqr"
6182bf6240SBarry Smith #define KSPPREONLY    "preonly"
6282bf6240SBarry Smith #define KSPQCG        "qcg"
63c9cf9b11SBarry Smith #define KSPBICG       "bicg"
64b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
6501c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
6621356919SSatish Balay #define KSPLCD        "lcd"
671d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
6858ad63f4SBarry Smith #define KSPGCR        "gcr"
69fad47a0aSPatrick Sanan #define KSPPIPEGCR    "pipegcr"
70e4d80e07Szianekhodja #define KSPTSIRM      "tsirm"
71e4d80e07Szianekhodja #define KSPCGLS       "cgls"
722eac72dbSBarry Smith 
738ba1e511SMatthew Knepley /* Logging support */
74014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
755358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
768ba1e511SMatthew Knepley 
77014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
7819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
79c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
80014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
81014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
82014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
83014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
84014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
85014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
8623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool);
8725c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool);
882eac72dbSBarry Smith 
89140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
90bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP));
9130de9b25SBarry Smith 
92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
93014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
94014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
95c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool );
97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool  *);
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool );
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *);
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool );
101014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool  *);
102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool );
1037d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool );
104c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *);
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool );
106c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
110014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
111734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*);
1122a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1132a7a6963SBarry 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);}
1142eac72dbSBarry Smith 
115d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
116d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
117d4211eb9SBarry Smith 
118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
120aabeff55SBarry Smith 
121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
124014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
125014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
126014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
1274b0e389bSBarry Smith 
128fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
129fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *);
130fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
131fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
132fa0ddf94SBarry Smith 
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
137014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
138b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
139b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
140b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
141b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
142014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *);
1430a71e23eSMatthew Knepley 
144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
1462eac72dbSBarry Smith 
147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
15058450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
15182fe98a5SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseRandom(KSP,PetscBool);
152ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom);
15358450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*);
154014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
155d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *);
156d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]);
1577d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt *,Vec[],PetscReal[],PetscReal[]);
1584b0e389bSBarry Smith 
159640f4f53SPatrick Sanan /*E
160640f4f53SPatrick Sanan 
16106137d0aSPatrick Sanan   KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods
162640f4f53SPatrick Sanan 
16306137d0aSPatrick Sanan   KSP_FCD_TRUNCATION uses all (up to mmax) stored directions
16406137d0aSPatrick Sanan   KSP_FCD_TRUNCATION_RESTART uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1..
165640f4f53SPatrick Sanan 
1662b26979fSBarry Smith    Level: intermediate
16706137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()
168640f4f53SPatrick Sanan 
169640f4f53SPatrick Sanan E*/
17006137d0aSPatrick Sanan typedef enum {KSP_FCD_TRUNCATION,KSP_FCD_TRUNCATION_RESTART} KSPFCDTruncationType;
17106137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
17206137d0aSPatrick Sanan 
17306137d0aSPatrick Sanan #define KSPFCDGetNumOldDirections(ctx,i,mi) ({                                                    \
17406137d0aSPatrick Sanan   if(ctx->trunctype == KSP_FCD_TRUNCATION_RESTART){                                               \
17506137d0aSPatrick Sanan     mi = ((i-1) % ctx->mmax)+1;                                                                   \
17606137d0aSPatrick Sanan     if (mi==1 && i!=1) ++(ctx->n_search_space_resets);                                            \
17706137d0aSPatrick Sanan   } else if (ctx->trunctype == KSP_FCD_TRUNCATION)                                                \
17806137d0aSPatrick Sanan     mi = ctx->mmax;                                                                               \
17906137d0aSPatrick Sanan  else {                                                                                           \
18006137d0aSPatrick Sanan    SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unrecognized Truncation Strategy");CHKERRQ(ierr); \
18106137d0aSPatrick Sanan   }                                                                                               \
18206137d0aSPatrick Sanan })
183640f4f53SPatrick Sanan 
184640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt);
185640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*);
186640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt);
187640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*);
18806137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType);
18906137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*);
190640f4f53SPatrick Sanan 
191390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt);
192390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*);
193390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt);
194390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*);
195390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType);
196390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*);
197390d8e47SPatrick Sanan 
198fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
199fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
200fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
201fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
202fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
203fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
204fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
205fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
206fad47a0aSPatrick Sanan 
207014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
208014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
209014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
2109f236934SBarry Smith 
211014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
214014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
2161d73ed98SBarry Smith 
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
2191d73ed98SBarry Smith 
220*483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
221*483d6965SPatrick Sanan 
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
22558ad63f4SBarry Smith 
226b49fd9e1SBarry Smith /*E
227b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
228b49fd9e1SBarry Smith 
229b49fd9e1SBarry Smith    Level: advanced
230b49fd9e1SBarry Smith 
231bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
232e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
233b49fd9e1SBarry Smith 
234b49fd9e1SBarry Smith E*/
23578d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2366a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2371f7e983dSSatish Balay /*MC
2381957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
2398c5b8ba0SBarry Smith 
2408c5b8ba0SBarry Smith    Level: advanced
2418c5b8ba0SBarry Smith 
2428c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2438c5b8ba0SBarry Smith 
244bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
245e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2468c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2478c5b8ba0SBarry Smith M*/
2488c5b8ba0SBarry Smith 
2491f7e983dSSatish Balay /*MC
2508c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2518c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2528c5b8ba0SBarry Smith           poor orthogonality.
2538c5b8ba0SBarry Smith 
2548c5b8ba0SBarry Smith    Level: advanced
2558c5b8ba0SBarry Smith 
2568c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2578c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2588c5b8ba0SBarry Smith 
259bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
260e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2618c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2628c5b8ba0SBarry Smith M*/
2638c5b8ba0SBarry Smith 
2641f7e983dSSatish Balay /*MC
2658c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2668c5b8ba0SBarry Smith 
2678c5b8ba0SBarry Smith    Level: advanced
2688c5b8ba0SBarry Smith 
2698c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2708c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2718c5b8ba0SBarry Smith 
2728c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2738c5b8ba0SBarry Smith 
274bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
275e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2768c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2778c5b8ba0SBarry Smith M*/
2788c5b8ba0SBarry Smith 
279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
28108480c60SBarry Smith 
282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
285c38d4ed2SBarry Smith 
286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
289121fd945SKris Buschelman 
290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
293e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
294d9492815SBarry Smith 
295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2972eac72dbSBarry Smith 
298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
299014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
301014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
302390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
303390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
304014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
305816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *);
306014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
307014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
308e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
309e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
310e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
311014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
31284cb2905SBarry Smith 
313014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
314014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
315c01c455dSBarry Smith 
31623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
31723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
318014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
319014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
320014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
321014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
3227287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
3237287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
3241eb62cbbSBarry Smith 
325014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
326014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
327014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
328014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
3291f7f0c4fSBarry Smith 
330014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
33155849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
332685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
3332a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
3342a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
33555849f57SBarry Smith 
33655849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
3371eb62cbbSBarry Smith 
338014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
339014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
340db9b2ab1SHong Zhang 
341014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
342014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
34368ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
34483ab6a24SBarry Smith 
34528ce4d24SBarry Smith /*E
3468a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
3478a4b9c5cSBarry Smith        test routines.
3488a4b9c5cSBarry Smith 
3498a4b9c5cSBarry Smith    Level: advanced
3508a4b9c5cSBarry Smith 
351a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
3521718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
353a3f661c8SBarry Smith 
354af0996ceSBarry Smith    Notes: this must match petsc/finclude/petscksp.h
3558a4b9c5cSBarry Smith 
35694b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
3571718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
3588a4b9c5cSBarry Smith E*/
3599e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
3609e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
361014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
362a21b2a99SBarry Smith 
3631f7e983dSSatish Balay /*MC
3649793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
3658c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
3668c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3678c5b8ba0SBarry Smith 
3688c5b8ba0SBarry Smith    Level: advanced
3698c5b8ba0SBarry Smith 
3701957e957SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
3718c5b8ba0SBarry Smith 
372ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3738c5b8ba0SBarry Smith M*/
3748c5b8ba0SBarry Smith 
3751f7e983dSSatish Balay /*MC
3761957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
3778c5b8ba0SBarry Smith        convergence test routine.
3788c5b8ba0SBarry Smith 
3798c5b8ba0SBarry Smith    Level: advanced
3808c5b8ba0SBarry Smith 
3819793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3828c5b8ba0SBarry Smith M*/
3838c5b8ba0SBarry Smith 
3841f7e983dSSatish Balay /*MC
385ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3868c5b8ba0SBarry Smith        convergence test routine.
3878c5b8ba0SBarry Smith 
3888c5b8ba0SBarry Smith    Level: advanced
3898c5b8ba0SBarry Smith 
3909793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3918c5b8ba0SBarry Smith M*/
3928c5b8ba0SBarry Smith 
3931f7e983dSSatish Balay /*MC
394ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
395390d8e47SPatrick Sanan        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
3968c5b8ba0SBarry Smith 
3978c5b8ba0SBarry Smith    Level: advanced
3988c5b8ba0SBarry Smith 
3999793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
4008c5b8ba0SBarry Smith M*/
4018c5b8ba0SBarry Smith 
402014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
405014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
406014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
4078a4b9c5cSBarry Smith 
4088a4b9c5cSBarry Smith /*E
4091957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
41028ce4d24SBarry Smith 
41128ce4d24SBarry Smith    Level: beginner
41228ce4d24SBarry Smith 
4134d0a8057SBarry Smith    Notes: See KSPGetConvergedReason() for explanation of each value
41428ce4d24SBarry Smith 
415af0996ceSBarry Smith    Developer notes: this must match petsc/finclude/petscksp.h
4164d0a8057SBarry Smith 
4174d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
4184d0a8057SBarry Smith       any of the values here also change them that array of names.
41986c02ca4SBarry Smith 
420c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
42128ce4d24SBarry Smith E*/
422d15094e1SBarry Smith typedef enum {/* converged */
4239ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
4249ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
425d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
426d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
427b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
4288031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
4298031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
430329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
431af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
432d15094e1SBarry Smith               /* diverged */
433b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
434d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
435d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
436d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
437b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
438b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
439b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
4404d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
4416aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
442422a814eSBarry Smith               KSP_DIVERGED_PCSETUP_FAILED      = -11,
443d15094e1SBarry Smith 
444d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
445014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
446d15094e1SBarry Smith 
447c838673bSBarry Smith /*MC
448c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
449c838673bSBarry Smith 
450c838673bSBarry Smith    Level: beginner
451c838673bSBarry Smith 
452c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
453c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
454c838673bSBarry Smith        2-norm of the residual for right preconditioning
455c838673bSBarry Smith 
456c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
457c838673bSBarry Smith 
458c838673bSBarry Smith M*/
459c838673bSBarry Smith 
460c838673bSBarry Smith /*MC
461c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
462c838673bSBarry Smith 
463c838673bSBarry Smith    Level: beginner
464c838673bSBarry Smith 
465c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
466c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
467c838673bSBarry Smith        2-norm of the residual for right preconditioning
468c838673bSBarry Smith 
469c838673bSBarry Smith    Level: beginner
470c838673bSBarry Smith 
471c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
472c838673bSBarry Smith 
473c838673bSBarry Smith M*/
474c838673bSBarry Smith 
475c838673bSBarry Smith /*MC
476c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
477c838673bSBarry Smith 
478c838673bSBarry Smith    Level: beginner
479c838673bSBarry Smith 
480c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
481c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
482c838673bSBarry Smith        2-norm of the residual for right preconditioning
483c838673bSBarry Smith 
484c838673bSBarry Smith    Level: beginner
485c838673bSBarry Smith 
486c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
487c838673bSBarry Smith 
488c838673bSBarry Smith M*/
489c838673bSBarry Smith 
490c838673bSBarry Smith /*MC
491c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
492c838673bSBarry Smith       reached
493c838673bSBarry Smith 
494c838673bSBarry Smith    Level: beginner
495c838673bSBarry Smith 
496c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
497c838673bSBarry Smith 
498c838673bSBarry Smith M*/
499c838673bSBarry Smith 
500c838673bSBarry Smith /*MC
5018c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
5020059c7bdSJed Brown            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
5034d0a8057SBarry Smith            test routine is set in KSP.
504c838673bSBarry Smith 
505c838673bSBarry Smith    Level: beginner
506c838673bSBarry Smith 
507c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
508c838673bSBarry Smith 
509c838673bSBarry Smith M*/
510c838673bSBarry Smith 
511c838673bSBarry Smith /*MC
512c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
5133014e516SBarry Smith           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
5143014e516SBarry Smith           preconditioner.
515c838673bSBarry Smith 
516c838673bSBarry Smith    Level: beginner
517c838673bSBarry Smith 
518c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
519c838673bSBarry Smith 
520c838673bSBarry Smith M*/
521c838673bSBarry Smith 
522c838673bSBarry Smith /*MC
523c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
524c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
525c838673bSBarry Smith 
526c838673bSBarry Smith    Level: beginner
527c838673bSBarry Smith 
528c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
529c838673bSBarry Smith 
530c838673bSBarry Smith M*/
531c838673bSBarry Smith 
532c838673bSBarry Smith /*MC
533c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
534c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
535c838673bSBarry Smith 
536c838673bSBarry Smith    Level: beginner
537c838673bSBarry Smith 
538c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
539c838673bSBarry Smith 
540c838673bSBarry Smith M*/
541c838673bSBarry Smith 
542c838673bSBarry Smith /*MC
543c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
544c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
545c838673bSBarry Smith         be positive definite
546c838673bSBarry Smith 
547c838673bSBarry Smith    Level: beginner
548c838673bSBarry Smith 
5492401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
550c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
551c838673bSBarry Smith 
552c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
553c838673bSBarry Smith 
554c838673bSBarry Smith M*/
555c838673bSBarry Smith 
556c838673bSBarry Smith /*MC
557c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
558c838673bSBarry Smith         while the KSPSolve() is still running.
559c838673bSBarry Smith 
560c838673bSBarry Smith    Level: beginner
561c838673bSBarry Smith 
562c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
563c838673bSBarry Smith 
564c838673bSBarry Smith M*/
565c838673bSBarry Smith 
566014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
567014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
5688de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
569014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
5708de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
5718de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
5728de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
5738de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
5740059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
576abef13c0SSatish Balay 
5778ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
5788ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
5798ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
5808ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
5818ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
5828ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
5838ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
5848ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
5858ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
5868ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
5878ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
5888ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
5898ea1b3e6SJed Brown 
590014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
591d4fbbf0eSBarry Smith 
59228ce4d24SBarry Smith /*E
59328ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
59428ce4d24SBarry Smith 
59528ce4d24SBarry Smith    Level: beginner
59628ce4d24SBarry Smith 
59728ce4d24SBarry Smith .seealso: KSPCGSetType()
59828ce4d24SBarry Smith E*/
5999dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
6006a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
60128ce4d24SBarry Smith 
602014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
603014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
6048031f4adStmunson 
605014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
608fcae7a14Stmunson 
609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
611014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
612e559a7feSLois Curfman McInnes 
613014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
614014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
615014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
616014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
617014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
6188031f4adStmunson 
619014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
6201d6018f0SLisandro Dalcin 
621014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
622014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
6233369ce9aSBarry Smith 
6249804daf3SBarry Smith #include <petscdrawtypes.h>
625d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
626d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
627d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
628d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
629014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
6302f2e5d10SKris Buschelman 
631014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
632014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
63303919abeSBarry Smith 
634f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
635ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
636f8a50e2bSBarry Smith 
637014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
638014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
639014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
640014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
641014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
642014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
643f8a50e2bSBarry Smith 
644014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
645014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
646014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
647f8a50e2bSBarry Smith 
648470b340bSDmitry Karpeev /*E
649470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
650470b340bSDmitry Karpeev 
651470b340bSDmitry Karpeev     Level: intermediate
652470b340bSDmitry Karpeev 
653470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
654470b340bSDmitry Karpeev E*/
655470b340bSDmitry Karpeev typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType;
656470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
657470b340bSDmitry Karpeev 
658014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
659014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
660d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
661bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
662aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
663bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
664470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
665470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
6665bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
6675a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
668470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *);
669470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
6703f22127dSBarry Smith 
671014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
672014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
673014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
674014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
675014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
676fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
67723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
678fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
67923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
68023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
681014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
682014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
683fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
684fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
6856c699258SBarry Smith 
68602b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
68702b02e71SToby 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);
6882eac72dbSBarry Smith #endif
689