xref: /petsc/include/petscksp.h (revision 93f1e87bd04011f4c6aae5303e5e0ec8df4fe3f0)
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"
45483d6965SPatrick 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 
163*93f1e87bSPatrick Sanan   KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
164*93f1e87bSPatrick Sanan   KSP_FCD_TRUNC_TYPE_NOTAY 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*/
170*93f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType;
17106137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[];
17206137d0aSPatrick Sanan 
17306137d0aSPatrick Sanan #define KSPFCDGetNumOldDirections(ctx,i,mi) ({                                                    \
174*93f1e87bSPatrick Sanan   if(ctx->truncstrat == KSP_FCD_TRUNC_TYPE_NOTAY){                                                \
17506137d0aSPatrick Sanan     mi = ((i-1) % ctx->mmax)+1;                                                                   \
176*93f1e87bSPatrick Sanan   } else if (ctx->truncstrat == KSP_FCD_TRUNC_TYPE_STANDARD)                                      \
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 
197fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt);
198fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*);
199fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt);
200fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*);
201fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType);
202fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*);
203fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool);
204fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*);
205fad47a0aSPatrick 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 
219483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar);
220483d6965SPatrick Sanan 
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
22458ad63f4SBarry Smith 
225b49fd9e1SBarry Smith /*E
226b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
227b49fd9e1SBarry Smith 
228b49fd9e1SBarry Smith    Level: advanced
229b49fd9e1SBarry Smith 
230bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
231e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
232b49fd9e1SBarry Smith 
233b49fd9e1SBarry Smith E*/
23478d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2356a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2361f7e983dSSatish Balay /*MC
2371957e957SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process
2388c5b8ba0SBarry Smith 
2398c5b8ba0SBarry Smith    Level: advanced
2408c5b8ba0SBarry Smith 
2418c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2428c5b8ba0SBarry Smith 
243bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
244e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2458c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2468c5b8ba0SBarry Smith M*/
2478c5b8ba0SBarry Smith 
2481f7e983dSSatish Balay /*MC
2498c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2508c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2518c5b8ba0SBarry Smith           poor orthogonality.
2528c5b8ba0SBarry Smith 
2538c5b8ba0SBarry Smith    Level: advanced
2548c5b8ba0SBarry Smith 
2558c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2568c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2578c5b8ba0SBarry Smith 
258bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
259e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2608c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2618c5b8ba0SBarry Smith M*/
2628c5b8ba0SBarry Smith 
2631f7e983dSSatish Balay /*MC
2648c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2658c5b8ba0SBarry Smith 
2668c5b8ba0SBarry Smith    Level: advanced
2678c5b8ba0SBarry Smith 
2688c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2698c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2708c5b8ba0SBarry Smith 
2718c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2728c5b8ba0SBarry Smith 
273bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
274e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2758c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2768c5b8ba0SBarry Smith M*/
2778c5b8ba0SBarry Smith 
278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
28008480c60SBarry Smith 
281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
284c38d4ed2SBarry Smith 
285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
288121fd945SKris Buschelman 
289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
292e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
293d9492815SBarry Smith 
294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2962eac72dbSBarry Smith 
297014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
299014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
301390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
302390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
303014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
304816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *);
305014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
306014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
307e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*);
308e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**);
309e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**);
310014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
31184cb2905SBarry Smith 
312014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
313014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
314c01c455dSBarry Smith 
31523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat);
31623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*);
317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
318014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
319014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
320014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
3217287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
3227287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
3231eb62cbbSBarry Smith 
324014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
325014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
326014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
327014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
3281f7f0c4fSBarry Smith 
329014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
33055849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
331685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
3322a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer);
3332a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP);
33455849f57SBarry Smith 
33555849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
3361eb62cbbSBarry Smith 
337014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
338014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
339db9b2ab1SHong Zhang 
340014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
341014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
34268ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*);
34383ab6a24SBarry Smith 
34428ce4d24SBarry Smith /*E
3458a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
3468a4b9c5cSBarry Smith        test routines.
3478a4b9c5cSBarry Smith 
3488a4b9c5cSBarry Smith    Level: advanced
3498a4b9c5cSBarry Smith 
350a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
3511718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
352a3f661c8SBarry Smith 
353af0996ceSBarry Smith    Notes: this must match petsc/finclude/petscksp.h
3548a4b9c5cSBarry Smith 
35594b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
3561718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
3578a4b9c5cSBarry Smith E*/
3589e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
3599e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
360014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
361a21b2a99SBarry Smith 
3621f7e983dSSatish Balay /*MC
3639793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
3648c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
3658c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3668c5b8ba0SBarry Smith 
3678c5b8ba0SBarry Smith    Level: advanced
3688c5b8ba0SBarry Smith 
3691957e957SBarry Smith     Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored
3708c5b8ba0SBarry Smith 
371ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3728c5b8ba0SBarry Smith M*/
3738c5b8ba0SBarry Smith 
3741f7e983dSSatish Balay /*MC
3751957e957SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the
3768c5b8ba0SBarry Smith        convergence test routine.
3778c5b8ba0SBarry Smith 
3788c5b8ba0SBarry Smith    Level: advanced
3798c5b8ba0SBarry Smith 
3809793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3818c5b8ba0SBarry Smith M*/
3828c5b8ba0SBarry Smith 
3831f7e983dSSatish Balay /*MC
384ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3858c5b8ba0SBarry Smith        convergence test routine.
3868c5b8ba0SBarry Smith 
3878c5b8ba0SBarry Smith    Level: advanced
3888c5b8ba0SBarry Smith 
3899793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3908c5b8ba0SBarry Smith M*/
3918c5b8ba0SBarry Smith 
3921f7e983dSSatish Balay /*MC
393ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
394390d8e47SPatrick Sanan        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR
3958c5b8ba0SBarry Smith 
3968c5b8ba0SBarry Smith    Level: advanced
3978c5b8ba0SBarry Smith 
3989793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3998c5b8ba0SBarry Smith M*/
4008c5b8ba0SBarry Smith 
401014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
402014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
405014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
4068a4b9c5cSBarry Smith 
4078a4b9c5cSBarry Smith /*E
4081957e957SBarry Smith     KSPConvergedReason - reason a Krylov method was said to have converged or diverged
40928ce4d24SBarry Smith 
41028ce4d24SBarry Smith    Level: beginner
41128ce4d24SBarry Smith 
4124d0a8057SBarry Smith    Notes: See KSPGetConvergedReason() for explanation of each value
41328ce4d24SBarry Smith 
414af0996ceSBarry Smith    Developer notes: this must match petsc/finclude/petscksp.h
4154d0a8057SBarry Smith 
4164d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
4174d0a8057SBarry Smith       any of the values here also change them that array of names.
41886c02ca4SBarry Smith 
419c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
42028ce4d24SBarry Smith E*/
421d15094e1SBarry Smith typedef enum {/* converged */
4229ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
4239ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
424d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
425d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
426b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
4278031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
4288031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
429329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
430af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
431d15094e1SBarry Smith               /* diverged */
432b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
433d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
434d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
435d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
436b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
437b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
438b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
4394d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
4406aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
441422a814eSBarry Smith               KSP_DIVERGED_PCSETUP_FAILED      = -11,
442d15094e1SBarry Smith 
443d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
444014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
445d15094e1SBarry Smith 
446c838673bSBarry Smith /*MC
447c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
448c838673bSBarry Smith 
449c838673bSBarry Smith    Level: beginner
450c838673bSBarry Smith 
451c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
452c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
453c838673bSBarry Smith        2-norm of the residual for right preconditioning
454c838673bSBarry Smith 
455c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
456c838673bSBarry Smith 
457c838673bSBarry Smith M*/
458c838673bSBarry Smith 
459c838673bSBarry Smith /*MC
460c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
461c838673bSBarry Smith 
462c838673bSBarry Smith    Level: beginner
463c838673bSBarry Smith 
464c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
465c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
466c838673bSBarry Smith        2-norm of the residual for right preconditioning
467c838673bSBarry Smith 
468c838673bSBarry Smith    Level: beginner
469c838673bSBarry Smith 
470c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
471c838673bSBarry Smith 
472c838673bSBarry Smith M*/
473c838673bSBarry Smith 
474c838673bSBarry Smith /*MC
475c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
476c838673bSBarry Smith 
477c838673bSBarry Smith    Level: beginner
478c838673bSBarry Smith 
479c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
480c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
481c838673bSBarry Smith        2-norm of the residual for right preconditioning
482c838673bSBarry Smith 
483c838673bSBarry Smith    Level: beginner
484c838673bSBarry Smith 
485c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
486c838673bSBarry Smith 
487c838673bSBarry Smith M*/
488c838673bSBarry Smith 
489c838673bSBarry Smith /*MC
490c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
491c838673bSBarry Smith       reached
492c838673bSBarry Smith 
493c838673bSBarry Smith    Level: beginner
494c838673bSBarry Smith 
495c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
496c838673bSBarry Smith 
497c838673bSBarry Smith M*/
498c838673bSBarry Smith 
499c838673bSBarry Smith /*MC
5008c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
5010059c7bdSJed Brown            the preconditioner is applied. Also used when the KSPConvergedSkip() convergence
5024d0a8057SBarry Smith            test routine is set in KSP.
503c838673bSBarry Smith 
504c838673bSBarry Smith    Level: beginner
505c838673bSBarry Smith 
506c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
507c838673bSBarry Smith 
508c838673bSBarry Smith M*/
509c838673bSBarry Smith 
510c838673bSBarry Smith /*MC
511c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
5123014e516SBarry Smith           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
5133014e516SBarry Smith           preconditioner.
514c838673bSBarry Smith 
515c838673bSBarry Smith    Level: beginner
516c838673bSBarry Smith 
517c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
518c838673bSBarry Smith 
519c838673bSBarry Smith M*/
520c838673bSBarry Smith 
521c838673bSBarry Smith /*MC
522c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
523c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
524c838673bSBarry Smith 
525c838673bSBarry Smith    Level: beginner
526c838673bSBarry Smith 
527c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
528c838673bSBarry Smith 
529c838673bSBarry Smith M*/
530c838673bSBarry Smith 
531c838673bSBarry Smith /*MC
532c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
533c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
534c838673bSBarry Smith 
535c838673bSBarry Smith    Level: beginner
536c838673bSBarry Smith 
537c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
538c838673bSBarry Smith 
539c838673bSBarry Smith M*/
540c838673bSBarry Smith 
541c838673bSBarry Smith /*MC
542c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
543c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
544c838673bSBarry Smith         be positive definite
545c838673bSBarry Smith 
546c838673bSBarry Smith    Level: beginner
547c838673bSBarry Smith 
5482401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
549c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
550c838673bSBarry Smith 
551c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
552c838673bSBarry Smith 
553c838673bSBarry Smith M*/
554c838673bSBarry Smith 
555c838673bSBarry Smith /*MC
556c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
557c838673bSBarry Smith         while the KSPSolve() is still running.
558c838673bSBarry Smith 
559c838673bSBarry Smith    Level: beginner
560c838673bSBarry Smith 
561c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
562c838673bSBarry Smith 
563c838673bSBarry Smith M*/
564c838673bSBarry Smith 
565014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
566014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
5678de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
568014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
5698de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *);
5708de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **);
5718de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP);
5728de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP);
5730059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
574014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
575abef13c0SSatish Balay 
5768ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ }
5778ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault)
5788ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ }
5798ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy)
5808ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ }
5818ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate)
5828ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ }
5838ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm)
5848ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ }
5858ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm)
5868ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ }
5878ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip)
5888ea1b3e6SJed Brown 
589014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
590d4fbbf0eSBarry Smith 
59128ce4d24SBarry Smith /*E
59228ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
59328ce4d24SBarry Smith 
59428ce4d24SBarry Smith    Level: beginner
59528ce4d24SBarry Smith 
59628ce4d24SBarry Smith .seealso: KSPCGSetType()
59728ce4d24SBarry Smith E*/
5989dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5996a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
60028ce4d24SBarry Smith 
601014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
602014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
6038031f4adStmunson 
604014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
605014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
607fcae7a14Stmunson 
608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
611e559a7feSLois Curfman McInnes 
612014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
613014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
614014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
615014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
616014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
6178031f4adStmunson 
618014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
6191d6018f0SLisandro Dalcin 
620014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
621014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
6223369ce9aSBarry Smith 
6239804daf3SBarry Smith #include <petscdrawtypes.h>
624d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
625d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
626d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
627d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
628014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
6292f2e5d10SKris Buschelman 
630014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
631014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
63203919abeSBarry Smith 
633f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
634ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
635f8a50e2bSBarry Smith 
636014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
637014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
638014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
639014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
640014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
641014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
642f8a50e2bSBarry Smith 
643014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
644014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
645014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
646f8a50e2bSBarry Smith 
647470b340bSDmitry Karpeev /*E
648470b340bSDmitry Karpeev     MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines
649470b340bSDmitry Karpeev 
650470b340bSDmitry Karpeev     Level: intermediate
651470b340bSDmitry Karpeev 
652470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat()
653470b340bSDmitry Karpeev E*/
654470b340bSDmitry Karpeev typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType;
655470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[];
656470b340bSDmitry Karpeev 
657014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
658014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
659d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
660bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
661aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat);
662bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
663470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType);
664470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*);
6655bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*);
6665a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*);
667470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *);
668470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*);
6693f22127dSBarry Smith 
670014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
671014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
672014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
673014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
674014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
675fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
67623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
677fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
67823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*);
67923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*);
680014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
681014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
682fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
683fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
6846c699258SBarry Smith 
68502b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec);
68602b02e71SToby 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);
6872eac72dbSBarry Smith #endif
688