xref: /petsc/include/petscksp.h (revision 0a71e23eb7bd078e49ca454da70af7f66de52c5a)
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
60a835dfdSSatish Balay #include "petscpc.h"
7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
82eac72dbSBarry Smith 
9dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitializePackage(const char[]);
101dbb0a54SBarry Smith 
1128ce4d24SBarry Smith /*S
1228ce4d24SBarry Smith      KSP - Abstract PETSc object that manages all Krylov methods
1328ce4d24SBarry Smith 
1428ce4d24SBarry Smith    Level: beginner
1528ce4d24SBarry Smith 
1628ce4d24SBarry Smith   Concepts: Krylov methods
1728ce4d24SBarry Smith 
1894b7f48cSBarry Smith .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP
1928ce4d24SBarry Smith S*/
20e2a1c21fSSatish Balay typedef struct _p_KSP*     KSP;
212eac72dbSBarry Smith 
2228ce4d24SBarry Smith /*E
2328ce4d24SBarry Smith     KSPType - String with the name of a PETSc Krylov method or the creation function
2428ce4d24SBarry Smith        with an optional dynamic library name, for example
2528ce4d24SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mykspcreate()
2628ce4d24SBarry Smith 
2728ce4d24SBarry Smith    Level: beginner
2828ce4d24SBarry Smith 
2928ce4d24SBarry Smith .seealso: KSPSetType(), KSP
3028ce4d24SBarry Smith E*/
3182bf6240SBarry Smith #define KSPRICHARDSON "richardson"
3282bf6240SBarry Smith #define KSPCHEBYCHEV  "chebychev"
3382bf6240SBarry Smith #define KSPCG         "cg"
34df2a969aSvictorle #define   KSPCGNE       "cgne"
3580e17431SMatthew Knepley #define   KSPSTCG       "stcg"
3682bf6240SBarry Smith #define KSPGMRES      "gmres"
379dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
389dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
3982bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
4082bf6240SBarry Smith #define KSPBCGS       "bcgs"
41f0037002Svictorle #define KSPBCGSL      "bcgsl"
4282bf6240SBarry Smith #define KSPCGS        "cgs"
4382bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
4482bf6240SBarry Smith #define KSPCR         "cr"
4582bf6240SBarry Smith #define KSPLSQR       "lsqr"
4682bf6240SBarry Smith #define KSPPREONLY    "preonly"
4782bf6240SBarry Smith #define KSPQCG        "qcg"
48c9cf9b11SBarry Smith #define KSPBICG       "bicg"
49b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
5001c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
5121356919SSatish Balay #define KSPLCD        "lcd"
52f69a0ea3SMatthew Knepley #define KSPType const char*
532eac72dbSBarry Smith 
548ba1e511SMatthew Knepley /* Logging support */
55dba47a55SKris Buschelman extern PetscCookie PETSCKSP_DLLEXPORT KSP_COOKIE;
568ba1e511SMatthew Knepley 
57dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *);
58f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,KSPType);
59dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP);
60dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP);
61dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec);
62dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec);
63dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP);
642eac72dbSBarry Smith 
65b0a32e0cSBarry Smith extern PetscFList KSPList;
66dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]);
67dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void);
68dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));
6930de9b25SBarry Smith 
7030de9b25SBarry Smith /*MC
7130de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
7230de9b25SBarry Smith 
7330de9b25SBarry Smith    Synopsis:
74d360dc6fSBarry Smith    PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP))
7530de9b25SBarry Smith 
7630de9b25SBarry Smith    Not Collective
7730de9b25SBarry Smith 
7830de9b25SBarry Smith    Input Parameters:
7930de9b25SBarry Smith +  name_solver - name of a new user-defined solver
8030de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
8130de9b25SBarry Smith .  name_create - name of routine to create method context
8230de9b25SBarry Smith -  routine_create - routine to create method context
8330de9b25SBarry Smith 
8430de9b25SBarry Smith    Notes:
8530de9b25SBarry Smith    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
8630de9b25SBarry Smith 
8730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
8830de9b25SBarry Smith    is ignored.
8930de9b25SBarry Smith 
9030de9b25SBarry Smith    Sample usage:
9130de9b25SBarry Smith .vb
9230de9b25SBarry Smith    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
9330de9b25SBarry Smith                "MySolverCreate",MySolverCreate);
9430de9b25SBarry Smith .ve
9530de9b25SBarry Smith 
9630de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
9730de9b25SBarry Smith $     KSPSetType(ksp,"my_solver")
9830de9b25SBarry Smith    or at runtime via the option
9930de9b25SBarry Smith $     -ksp_type my_solver
10030de9b25SBarry Smith 
10130de9b25SBarry Smith    Level: advanced
10230de9b25SBarry Smith 
103ab901514SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
10430de9b25SBarry Smith           and others of the form ${any_environmental_variable} occuring in pathname will be
10530de9b25SBarry Smith           replaced with appropriate values.
10630de9b25SBarry Smith          If your function is not being put into a shared library then use KSPRegister() instead
10730de9b25SBarry Smith 
10830de9b25SBarry Smith .keywords: KSP, register
10930de9b25SBarry Smith 
11030de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy()
11130de9b25SBarry Smith 
11230de9b25SBarry Smith M*/
113aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
114f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
1156df38c32SLois Curfman McInnes #else
116f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
1176df38c32SLois Curfman McInnes #endif
11882bf6240SBarry Smith 
119dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,KSPType *);
120dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPreconditionerSide(KSP,PCSide);
121dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPreconditionerSide(KSP,PCSide*);
122dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
123dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
124dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscTruth);
125dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscTruth *);
126dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscTruth);
127dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscTruth*);
128a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscTruth*);
129dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth);
130a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscTruth*);
131dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth);
132dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *);
133dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *);
134dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*);
135dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*);
136dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace);
137dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*);
1382eac72dbSBarry Smith 
139dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC);
140dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*);
141aabeff55SBarry Smith 
142dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetMonitor(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
143dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPClearMonitor(KSP);
144dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **);
145dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
146dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);
1474b0e389bSBarry Smith 
1480e33f6ddSBarry Smith /* not sure where to put this */
149dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*);
150dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
151dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
152dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
1532eac72dbSBarry Smith 
154*0a71e23eSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC,KSP *);
155*0a71e23eSMatthew Knepley 
156dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *);
157dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *);
1582eac72dbSBarry Smith 
159dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal);
160dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
161dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
162dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
163dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
1644b0e389bSBarry Smith 
165dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt);
166dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal);
1679f236934SBarry Smith 
168dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP);
169dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
170dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
171dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1721d73ed98SBarry Smith 
173dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt);
174dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP);
1751d73ed98SBarry Smith 
176b49fd9e1SBarry Smith /*E
177b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
178b49fd9e1SBarry Smith 
179b49fd9e1SBarry Smith    Level: advanced
180b49fd9e1SBarry Smith 
181b49fd9e1SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1828c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
183b49fd9e1SBarry Smith 
184b49fd9e1SBarry Smith E*/
18578d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
1869dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[];
1871f7e983dSSatish Balay /*MC
1888c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
1898c5b8ba0SBarry Smith 
1908c5b8ba0SBarry Smith    Level: advanced
1918c5b8ba0SBarry Smith 
1928c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
1938c5b8ba0SBarry Smith 
1948c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1958c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
1968c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
1978c5b8ba0SBarry Smith M*/
1988c5b8ba0SBarry Smith 
1991f7e983dSSatish Balay /*MC
2008c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2018c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2028c5b8ba0SBarry Smith           poor orthogonality.
2038c5b8ba0SBarry Smith 
2048c5b8ba0SBarry Smith    Level: advanced
2058c5b8ba0SBarry Smith 
2068c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2078c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2088c5b8ba0SBarry Smith 
2098c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2108c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2118c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2128c5b8ba0SBarry Smith M*/
2138c5b8ba0SBarry Smith 
2141f7e983dSSatish Balay /*MC
2158c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2168c5b8ba0SBarry Smith 
2178c5b8ba0SBarry Smith    Level: advanced
2188c5b8ba0SBarry Smith 
2198c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2208c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2218c5b8ba0SBarry Smith 
2228c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2238c5b8ba0SBarry Smith 
2248c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2258c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2268c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2278c5b8ba0SBarry Smith M*/
2288c5b8ba0SBarry Smith 
229dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
23008480c60SBarry Smith 
231dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
232dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
233dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
234c38d4ed2SBarry Smith 
235dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal);
236dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*);
237dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*);
238121fd945SKris Buschelman 
239d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal);
240d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth);
241d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,int);
242d9492815SBarry Smith 
243dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP);
244dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2452eac72dbSBarry Smith 
246dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSingularValueMonitor(KSP,PetscInt,PetscReal,void *);
247dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultMonitor(KSP,PetscInt,PetscReal,void *);
248dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPTrueMonitor(KSP,PetscInt,PetscReal,void *);
249dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultSMonitor(KSP,PetscInt,PetscReal,void *);
250dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPVecViewMonitor(KSP,PetscInt,PetscReal,void *);
251dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESKrylovMonitor(KSP,PetscInt,PetscReal,void *);
25284cb2905SBarry Smith 
253dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec);
254dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*);
255dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
256dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
257c01c455dSBarry Smith 
258dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure);
259dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
260dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]);
261dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]);
2622dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]);
2631eb62cbbSBarry Smith 
264dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth);
265dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*);
266dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth);
267dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*);
2681f7f0c4fSBarry Smith 
269dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer);
2701eb62cbbSBarry Smith 
27128ce4d24SBarry Smith /*E
2728a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
2738a4b9c5cSBarry Smith        test routines.
2748a4b9c5cSBarry Smith 
2758a4b9c5cSBarry Smith    Level: advanced
2768a4b9c5cSBarry Smith 
2778a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
2788a4b9c5cSBarry Smith 
27994b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
2808a4b9c5cSBarry Smith           KSPSetConvergenceTest()
2818a4b9c5cSBarry Smith E*/
282954b3ec6SBarry Smith typedef enum {KSP_NO_NORM = 0,KSP_PRECONDITIONED_NORM = 1,KSP_UNPRECONDITIONED_NORM = 2,KSP_NATURAL_NORM = 3} KSPNormType;
2839dcbbd2bSBarry Smith extern const char *KSPNormTypes[];
2841f7e983dSSatish Balay /*MC
2858c5b8ba0SBarry Smith     KSP_NO_NORM - Do not compute a norm during the Krylov process. This will
2868c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
2878c5b8ba0SBarry Smith           be based on a norm of a residual etc.
2888c5b8ba0SBarry Smith 
2898c5b8ba0SBarry Smith    Level: advanced
2908c5b8ba0SBarry Smith 
291085a36d4SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
2928c5b8ba0SBarry Smith 
2938c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM
2948c5b8ba0SBarry Smith M*/
2958c5b8ba0SBarry Smith 
2961f7e983dSSatish Balay /*MC
2978c5b8ba0SBarry Smith     KSP_PRECONDITIONED_NORM - Compute the norm of the preconditioned residual and pass that to the
2988c5b8ba0SBarry Smith        convergence test routine.
2998c5b8ba0SBarry Smith 
3008c5b8ba0SBarry Smith    Level: advanced
3018c5b8ba0SBarry Smith 
3028c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
3038c5b8ba0SBarry Smith M*/
3048c5b8ba0SBarry Smith 
3051f7e983dSSatish Balay /*MC
3068c5b8ba0SBarry Smith     KSP_UNPRECONDITIONED_NORM - Compute the norm of the true residual (b - A*x) and pass that to the
3078c5b8ba0SBarry Smith        convergence test routine.
3088c5b8ba0SBarry Smith 
3098c5b8ba0SBarry Smith    Level: advanced
3108c5b8ba0SBarry Smith 
3118c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
3128c5b8ba0SBarry Smith M*/
3138c5b8ba0SBarry Smith 
3141f7e983dSSatish Balay /*MC
3158c5b8ba0SBarry Smith     KSP_NATURAL_NORM - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
3168c5b8ba0SBarry Smith        convergence test routine.
3178c5b8ba0SBarry Smith 
3188c5b8ba0SBarry Smith    Level: advanced
3198c5b8ba0SBarry Smith 
3208c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSPSetConvergenceTest()
3218c5b8ba0SBarry Smith M*/
3228c5b8ba0SBarry Smith 
323dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType);
3248a4b9c5cSBarry Smith 
3258a4b9c5cSBarry Smith /*E
32628ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
32728ce4d24SBarry Smith          have converged or diverged
32828ce4d24SBarry Smith 
32928ce4d24SBarry Smith    Level: beginner
33028ce4d24SBarry Smith 
33128ce4d24SBarry Smith    Notes: this must match finclude/petscksp.h
33228ce4d24SBarry Smith 
33386c02ca4SBarry Smith    Developer note: The string versions of these are in
3344b9ad928SBarry Smith      src/ksp/ksp/interface/itfunc.c called convergedreasons.
33550f1acb2SBarry Smith      If these enums are changed you must change those.
33686c02ca4SBarry Smith 
337c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
33828ce4d24SBarry Smith E*/
339d15094e1SBarry Smith typedef enum {/* converged */
340d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
341d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
342b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
343a4ce1dd3SMatthew Knepley               KSP_CONVERGED_STCG_NEG_CURVE     =  5,
344a4ce1dd3SMatthew Knepley               KSP_CONVERGED_STCG_CONSTRAINED   =  6,
345329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
346d15094e1SBarry Smith               /* diverged */
347b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
348d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
349d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
350d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
351b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
352b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
353b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
3546aee1118SBarry Smith               KSP_DIVERGED_NAN                 = -9,
3556aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
356d15094e1SBarry Smith 
357d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
3589dcbbd2bSBarry Smith extern const char **KSPConvergedReasons;
359d15094e1SBarry Smith 
360c838673bSBarry Smith /*MC
361c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
362c838673bSBarry Smith 
363c838673bSBarry Smith    Level: beginner
364c838673bSBarry Smith 
365c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
366c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
367c838673bSBarry Smith        2-norm of the residual for right preconditioning
368c838673bSBarry Smith 
369c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
370c838673bSBarry Smith 
371c838673bSBarry Smith M*/
372c838673bSBarry Smith 
373c838673bSBarry Smith /*MC
374c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
375c838673bSBarry Smith 
376c838673bSBarry Smith    Level: beginner
377c838673bSBarry Smith 
378c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
379c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
380c838673bSBarry Smith        2-norm of the residual for right preconditioning
381c838673bSBarry Smith 
382c838673bSBarry Smith    Level: beginner
383c838673bSBarry Smith 
384c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
385c838673bSBarry Smith 
386c838673bSBarry Smith M*/
387c838673bSBarry Smith 
388c838673bSBarry Smith /*MC
389c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
390c838673bSBarry Smith 
391c838673bSBarry Smith    Level: beginner
392c838673bSBarry Smith 
393c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
394c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
395c838673bSBarry Smith        2-norm of the residual for right preconditioning
396c838673bSBarry Smith 
397c838673bSBarry Smith    Level: beginner
398c838673bSBarry Smith 
399c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
400c838673bSBarry Smith 
401c838673bSBarry Smith M*/
402c838673bSBarry Smith 
403c838673bSBarry Smith /*MC
404c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
405c838673bSBarry Smith       reached
406c838673bSBarry Smith 
407c838673bSBarry Smith    Level: beginner
408c838673bSBarry Smith 
409c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
410c838673bSBarry Smith 
411c838673bSBarry Smith M*/
412c838673bSBarry Smith 
413c838673bSBarry Smith /*MC
414c838673bSBarry Smith      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of the
415c838673bSBarry Smith            preconditioner is applied.
416c838673bSBarry Smith 
417c838673bSBarry Smith 
418c838673bSBarry Smith    Level: beginner
419c838673bSBarry Smith 
420c838673bSBarry Smith 
421c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
422c838673bSBarry Smith 
423c838673bSBarry Smith M*/
424c838673bSBarry Smith 
425c838673bSBarry Smith /*MC
426c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
427c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
428c838673bSBarry Smith 
429c838673bSBarry Smith    Level: beginner
430c838673bSBarry Smith 
431c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
432c838673bSBarry Smith 
433c838673bSBarry Smith M*/
434c838673bSBarry Smith 
435c838673bSBarry Smith /*MC
436c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
437c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
438c838673bSBarry Smith 
439c838673bSBarry Smith 
440c838673bSBarry Smith    Level: beginner
441c838673bSBarry Smith 
442c838673bSBarry Smith 
443c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
444c838673bSBarry Smith 
445c838673bSBarry Smith M*/
446c838673bSBarry Smith 
447c838673bSBarry Smith /*MC
448c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
449c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
450c838673bSBarry Smith 
451c838673bSBarry Smith    Level: beginner
452c838673bSBarry Smith 
453c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
454c838673bSBarry Smith 
455c838673bSBarry Smith M*/
456c838673bSBarry Smith 
457c838673bSBarry Smith /*MC
458c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
459c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
460c838673bSBarry Smith         be positive definite
461c838673bSBarry Smith 
462c838673bSBarry Smith    Level: beginner
463c838673bSBarry Smith 
4642401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
465c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
466c838673bSBarry Smith 
467c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
468c838673bSBarry Smith 
469c838673bSBarry Smith M*/
470c838673bSBarry Smith 
471c838673bSBarry Smith /*MC
472c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
473c838673bSBarry Smith         while the KSPSolve() is still running.
474c838673bSBarry Smith 
475c838673bSBarry Smith    Level: beginner
476c838673bSBarry Smith 
477c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
478c838673bSBarry Smith 
479c838673bSBarry Smith M*/
480c838673bSBarry Smith 
481dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *);
482dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **);
483dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
48494b02367SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUIRNorm(KSP);
485dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
486dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *);
487abef13c0SSatish Balay 
488dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *);
489d4fbbf0eSBarry Smith 
49028ce4d24SBarry Smith /*E
49128ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
49228ce4d24SBarry Smith 
49328ce4d24SBarry Smith    Level: beginner
49428ce4d24SBarry Smith 
49528ce4d24SBarry Smith .seealso: KSPCGSetType()
49628ce4d24SBarry Smith E*/
4979dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
4989dcbbd2bSBarry Smith extern const char *KSPCGTypes[];
49928ce4d24SBarry Smith 
500dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType);
5010d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal);
5020d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetQuadratic(KSP,PetscReal*);
503e559a7feSLois Curfman McInnes 
504dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP);
505dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP);
5063369ce9aSBarry Smith 
507dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
508dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitor(KSP,PetscInt,PetscReal,void*);
509dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorDestroy(PetscDrawLG);
510dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
511dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitor(KSP,PetscInt,PetscReal,void*);
512dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorDestroy(PetscDrawLG);
5132f2e5d10SKris Buschelman 
5140338f08bSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
51503919abeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
51603919abeSBarry Smith 
517e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
5182eac72dbSBarry Smith #endif
519