xref: /petsc/include/petscksp.h (revision 906ed7cc33fecbafab01746eec64dcdcc8a4842f)
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*);
138*906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1392eac72dbSBarry Smith 
140dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC);
141dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*);
142aabeff55SBarry Smith 
143dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetMonitor(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
144dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPClearMonitor(KSP);
145dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **);
146dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
147dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);
1484b0e389bSBarry Smith 
1490e33f6ddSBarry Smith /* not sure where to put this */
150dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*);
151dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
152dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
153dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
1542eac72dbSBarry Smith 
1550a71e23eSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC,KSP *);
1560a71e23eSMatthew Knepley 
157dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *);
158dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *);
1592eac72dbSBarry Smith 
160dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal);
161dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
162dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
163dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
164dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
1654b0e389bSBarry Smith 
166dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt);
167dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal);
1689f236934SBarry Smith 
169dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP);
170dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
171dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
172dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1731d73ed98SBarry Smith 
174dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt);
175dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP);
1761d73ed98SBarry Smith 
177b49fd9e1SBarry Smith /*E
178b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
179b49fd9e1SBarry Smith 
180b49fd9e1SBarry Smith    Level: advanced
181b49fd9e1SBarry Smith 
182b49fd9e1SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1838c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
184b49fd9e1SBarry Smith 
185b49fd9e1SBarry Smith E*/
18678d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
1879dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[];
1881f7e983dSSatish Balay /*MC
1898c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
1908c5b8ba0SBarry Smith 
1918c5b8ba0SBarry Smith    Level: advanced
1928c5b8ba0SBarry Smith 
1938c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
1948c5b8ba0SBarry Smith 
1958c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1968c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
1978c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
1988c5b8ba0SBarry Smith M*/
1998c5b8ba0SBarry Smith 
2001f7e983dSSatish Balay /*MC
2018c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2028c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2038c5b8ba0SBarry Smith           poor orthogonality.
2048c5b8ba0SBarry Smith 
2058c5b8ba0SBarry Smith    Level: advanced
2068c5b8ba0SBarry Smith 
2078c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2088c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2098c5b8ba0SBarry Smith 
2108c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2118c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2128c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2138c5b8ba0SBarry Smith M*/
2148c5b8ba0SBarry Smith 
2151f7e983dSSatish Balay /*MC
2168c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2178c5b8ba0SBarry Smith 
2188c5b8ba0SBarry Smith    Level: advanced
2198c5b8ba0SBarry Smith 
2208c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2218c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2228c5b8ba0SBarry Smith 
2238c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2248c5b8ba0SBarry Smith 
2258c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2268c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2278c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2288c5b8ba0SBarry Smith M*/
2298c5b8ba0SBarry Smith 
230dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
23108480c60SBarry Smith 
232dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
233dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
234dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
235c38d4ed2SBarry Smith 
236dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal);
237dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*);
238dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*);
239121fd945SKris Buschelman 
240d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal);
241d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth);
242d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,int);
243d9492815SBarry Smith 
244dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP);
245dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2462eac72dbSBarry Smith 
247dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSingularValueMonitor(KSP,PetscInt,PetscReal,void *);
248dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultMonitor(KSP,PetscInt,PetscReal,void *);
249dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPTrueMonitor(KSP,PetscInt,PetscReal,void *);
250dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultSMonitor(KSP,PetscInt,PetscReal,void *);
251dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPVecViewMonitor(KSP,PetscInt,PetscReal,void *);
252dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESKrylovMonitor(KSP,PetscInt,PetscReal,void *);
25384cb2905SBarry Smith 
254dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec);
255dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*);
256dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
257dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
258c01c455dSBarry Smith 
259dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure);
260dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
261*906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperatorsSet(KSP,PetscTruth*,PetscTruth*);
262dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]);
263dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]);
2642dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]);
2651eb62cbbSBarry Smith 
266dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth);
267dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*);
268dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth);
269dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*);
2701f7f0c4fSBarry Smith 
271dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer);
2721eb62cbbSBarry Smith 
27328ce4d24SBarry Smith /*E
2748a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
2758a4b9c5cSBarry Smith        test routines.
2768a4b9c5cSBarry Smith 
2778a4b9c5cSBarry Smith    Level: advanced
2788a4b9c5cSBarry Smith 
2798a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
2808a4b9c5cSBarry Smith 
28194b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
2828a4b9c5cSBarry Smith           KSPSetConvergenceTest()
2838a4b9c5cSBarry Smith E*/
284954b3ec6SBarry Smith typedef enum {KSP_NO_NORM = 0,KSP_PRECONDITIONED_NORM = 1,KSP_UNPRECONDITIONED_NORM = 2,KSP_NATURAL_NORM = 3} KSPNormType;
2859dcbbd2bSBarry Smith extern const char *KSPNormTypes[];
2861f7e983dSSatish Balay /*MC
2878c5b8ba0SBarry Smith     KSP_NO_NORM - Do not compute a norm during the Krylov process. This will
2888c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
2898c5b8ba0SBarry Smith           be based on a norm of a residual etc.
2908c5b8ba0SBarry Smith 
2918c5b8ba0SBarry Smith    Level: advanced
2928c5b8ba0SBarry Smith 
293085a36d4SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
2948c5b8ba0SBarry Smith 
2958c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM
2968c5b8ba0SBarry Smith M*/
2978c5b8ba0SBarry Smith 
2981f7e983dSSatish Balay /*MC
2998c5b8ba0SBarry Smith     KSP_PRECONDITIONED_NORM - Compute the norm of the preconditioned residual and pass that to the
3008c5b8ba0SBarry Smith        convergence test routine.
3018c5b8ba0SBarry Smith 
3028c5b8ba0SBarry Smith    Level: advanced
3038c5b8ba0SBarry Smith 
3048c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
3058c5b8ba0SBarry Smith M*/
3068c5b8ba0SBarry Smith 
3071f7e983dSSatish Balay /*MC
3088c5b8ba0SBarry Smith     KSP_UNPRECONDITIONED_NORM - Compute the norm of the true residual (b - A*x) and pass that to the
3098c5b8ba0SBarry Smith        convergence test routine.
3108c5b8ba0SBarry Smith 
3118c5b8ba0SBarry Smith    Level: advanced
3128c5b8ba0SBarry Smith 
3138c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
3148c5b8ba0SBarry Smith M*/
3158c5b8ba0SBarry Smith 
3161f7e983dSSatish Balay /*MC
3178c5b8ba0SBarry Smith     KSP_NATURAL_NORM - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
3188c5b8ba0SBarry Smith        convergence test routine.
3198c5b8ba0SBarry Smith 
3208c5b8ba0SBarry Smith    Level: advanced
3218c5b8ba0SBarry Smith 
3228c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSPSetConvergenceTest()
3238c5b8ba0SBarry Smith M*/
3248c5b8ba0SBarry Smith 
325dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType);
3268a4b9c5cSBarry Smith 
3278a4b9c5cSBarry Smith /*E
32828ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
32928ce4d24SBarry Smith          have converged or diverged
33028ce4d24SBarry Smith 
33128ce4d24SBarry Smith    Level: beginner
33228ce4d24SBarry Smith 
33328ce4d24SBarry Smith    Notes: this must match finclude/petscksp.h
33428ce4d24SBarry Smith 
33586c02ca4SBarry Smith    Developer note: The string versions of these are in
3364b9ad928SBarry Smith      src/ksp/ksp/interface/itfunc.c called convergedreasons.
33750f1acb2SBarry Smith      If these enums are changed you must change those.
33886c02ca4SBarry Smith 
339c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
34028ce4d24SBarry Smith E*/
341d15094e1SBarry Smith typedef enum {/* converged */
342d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
343d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
344b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
345a4ce1dd3SMatthew Knepley               KSP_CONVERGED_STCG_NEG_CURVE     =  5,
346a4ce1dd3SMatthew Knepley               KSP_CONVERGED_STCG_CONSTRAINED   =  6,
347329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
348af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
349d15094e1SBarry Smith               /* diverged */
350b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
351d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
352d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
353d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
354b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
355b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
356b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
3576aee1118SBarry Smith               KSP_DIVERGED_NAN                 = -9,
3586aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
359d15094e1SBarry Smith 
360d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
3619dcbbd2bSBarry Smith extern const char **KSPConvergedReasons;
362d15094e1SBarry Smith 
363c838673bSBarry Smith /*MC
364c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
365c838673bSBarry Smith 
366c838673bSBarry Smith    Level: beginner
367c838673bSBarry Smith 
368c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
369c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
370c838673bSBarry Smith        2-norm of the residual for right preconditioning
371c838673bSBarry Smith 
372c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
373c838673bSBarry Smith 
374c838673bSBarry Smith M*/
375c838673bSBarry Smith 
376c838673bSBarry Smith /*MC
377c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
378c838673bSBarry Smith 
379c838673bSBarry Smith    Level: beginner
380c838673bSBarry Smith 
381c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
382c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
383c838673bSBarry Smith        2-norm of the residual for right preconditioning
384c838673bSBarry Smith 
385c838673bSBarry Smith    Level: beginner
386c838673bSBarry Smith 
387c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
388c838673bSBarry Smith 
389c838673bSBarry Smith M*/
390c838673bSBarry Smith 
391c838673bSBarry Smith /*MC
392c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
393c838673bSBarry Smith 
394c838673bSBarry Smith    Level: beginner
395c838673bSBarry Smith 
396c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
397c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
398c838673bSBarry Smith        2-norm of the residual for right preconditioning
399c838673bSBarry Smith 
400c838673bSBarry Smith    Level: beginner
401c838673bSBarry Smith 
402c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
403c838673bSBarry Smith 
404c838673bSBarry Smith M*/
405c838673bSBarry Smith 
406c838673bSBarry Smith /*MC
407c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
408c838673bSBarry Smith       reached
409c838673bSBarry Smith 
410c838673bSBarry Smith    Level: beginner
411c838673bSBarry Smith 
412c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
413c838673bSBarry Smith 
414c838673bSBarry Smith M*/
415c838673bSBarry Smith 
416c838673bSBarry Smith /*MC
417c838673bSBarry Smith      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of the
418c838673bSBarry Smith            preconditioner is applied.
419c838673bSBarry Smith 
420c838673bSBarry Smith 
421c838673bSBarry Smith    Level: beginner
422c838673bSBarry Smith 
423c838673bSBarry Smith 
424c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
425c838673bSBarry Smith 
426c838673bSBarry Smith M*/
427c838673bSBarry Smith 
428c838673bSBarry Smith /*MC
429c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
430c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
431c838673bSBarry Smith 
432c838673bSBarry Smith    Level: beginner
433c838673bSBarry Smith 
434c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
435c838673bSBarry Smith 
436c838673bSBarry Smith M*/
437c838673bSBarry Smith 
438c838673bSBarry Smith /*MC
439c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
440c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
441c838673bSBarry Smith 
442c838673bSBarry Smith 
443c838673bSBarry Smith    Level: beginner
444c838673bSBarry Smith 
445c838673bSBarry Smith 
446c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
447c838673bSBarry Smith 
448c838673bSBarry Smith M*/
449c838673bSBarry Smith 
450c838673bSBarry Smith /*MC
451c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
452c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
453c838673bSBarry Smith 
454c838673bSBarry Smith    Level: beginner
455c838673bSBarry Smith 
456c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
457c838673bSBarry Smith 
458c838673bSBarry Smith M*/
459c838673bSBarry Smith 
460c838673bSBarry Smith /*MC
461c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
462c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
463c838673bSBarry Smith         be positive definite
464c838673bSBarry Smith 
465c838673bSBarry Smith    Level: beginner
466c838673bSBarry Smith 
4672401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
468c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
469c838673bSBarry Smith 
470c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
471c838673bSBarry Smith 
472c838673bSBarry Smith M*/
473c838673bSBarry Smith 
474c838673bSBarry Smith /*MC
475c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
476c838673bSBarry Smith         while the KSPSolve() is still running.
477c838673bSBarry Smith 
478c838673bSBarry Smith    Level: beginner
479c838673bSBarry Smith 
480c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
481c838673bSBarry Smith 
482c838673bSBarry Smith M*/
483c838673bSBarry Smith 
484dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *);
485dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **);
486dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
48794b02367SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUIRNorm(KSP);
488dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
489dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *);
490abef13c0SSatish Balay 
491dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *);
492d4fbbf0eSBarry Smith 
49328ce4d24SBarry Smith /*E
49428ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
49528ce4d24SBarry Smith 
49628ce4d24SBarry Smith    Level: beginner
49728ce4d24SBarry Smith 
49828ce4d24SBarry Smith .seealso: KSPCGSetType()
49928ce4d24SBarry Smith E*/
5009dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5019dcbbd2bSBarry Smith extern const char *KSPCGTypes[];
50228ce4d24SBarry Smith 
503dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType);
5040d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal);
5050d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetQuadratic(KSP,PetscReal*);
506e559a7feSLois Curfman McInnes 
507dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP);
508dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP);
5093369ce9aSBarry Smith 
510dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
511dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitor(KSP,PetscInt,PetscReal,void*);
512dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorDestroy(PetscDrawLG);
513dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
514dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitor(KSP,PetscInt,PetscReal,void*);
515dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorDestroy(PetscDrawLG);
5162f2e5d10SKris Buschelman 
5170338f08bSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
51803919abeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
51903919abeSBarry Smith 
520e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
5212eac72dbSBarry Smith #endif
522