xref: /petsc/include/petscksp.h (revision 0d2a7ff2b865a94f6b70adb47ee5ab4919353b54)
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"
51f69a0ea3SMatthew Knepley #define KSPType const char*
522eac72dbSBarry Smith 
538ba1e511SMatthew Knepley /* Logging support */
54dba47a55SKris Buschelman extern PetscCookie PETSCKSP_DLLEXPORT KSP_COOKIE;
558ba1e511SMatthew Knepley 
56dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *);
57f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,KSPType);
58dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP);
59dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP);
60dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec);
61dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec);
62dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP);
632eac72dbSBarry Smith 
64b0a32e0cSBarry Smith extern PetscFList KSPList;
65dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]);
66dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void);
67dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));
6830de9b25SBarry Smith 
6930de9b25SBarry Smith /*MC
7030de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
7130de9b25SBarry Smith 
7230de9b25SBarry Smith    Synopsis:
73d360dc6fSBarry Smith    PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP))
7430de9b25SBarry Smith 
7530de9b25SBarry Smith    Not Collective
7630de9b25SBarry Smith 
7730de9b25SBarry Smith    Input Parameters:
7830de9b25SBarry Smith +  name_solver - name of a new user-defined solver
7930de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
8030de9b25SBarry Smith .  name_create - name of routine to create method context
8130de9b25SBarry Smith -  routine_create - routine to create method context
8230de9b25SBarry Smith 
8330de9b25SBarry Smith    Notes:
8430de9b25SBarry Smith    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
8530de9b25SBarry Smith 
8630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
8730de9b25SBarry Smith    is ignored.
8830de9b25SBarry Smith 
8930de9b25SBarry Smith    Sample usage:
9030de9b25SBarry Smith .vb
9130de9b25SBarry Smith    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
9230de9b25SBarry Smith                "MySolverCreate",MySolverCreate);
9330de9b25SBarry Smith .ve
9430de9b25SBarry Smith 
9530de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
9630de9b25SBarry Smith $     KSPSetType(ksp,"my_solver")
9730de9b25SBarry Smith    or at runtime via the option
9830de9b25SBarry Smith $     -ksp_type my_solver
9930de9b25SBarry Smith 
10030de9b25SBarry Smith    Level: advanced
10130de9b25SBarry Smith 
102ab901514SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
10330de9b25SBarry Smith           and others of the form ${any_environmental_variable} occuring in pathname will be
10430de9b25SBarry Smith           replaced with appropriate values.
10530de9b25SBarry Smith          If your function is not being put into a shared library then use KSPRegister() instead
10630de9b25SBarry Smith 
10730de9b25SBarry Smith .keywords: KSP, register
10830de9b25SBarry Smith 
10930de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy()
11030de9b25SBarry Smith 
11130de9b25SBarry Smith M*/
112aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
113f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
1146df38c32SLois Curfman McInnes #else
115f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
1166df38c32SLois Curfman McInnes #endif
11782bf6240SBarry Smith 
118dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,KSPType *);
119dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPreconditionerSide(KSP,PCSide);
120dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPreconditionerSide(KSP,PCSide*);
121dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
122dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
123dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscTruth);
124dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscTruth *);
125dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscTruth);
126dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscTruth*);
127a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscTruth*);
128dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth);
129a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscTruth*);
130dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth);
131dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *);
132dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *);
133dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*);
134dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*);
135dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace);
136dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*);
1372eac72dbSBarry Smith 
138dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC);
139dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*);
140aabeff55SBarry Smith 
141dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetMonitor(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
142dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPClearMonitor(KSP);
143dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **);
144dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
145dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);
1464b0e389bSBarry Smith 
1470e33f6ddSBarry Smith /* not sure where to put this */
148dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*);
149dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
150dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
151dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
1522eac72dbSBarry Smith 
153dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *);
154dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *);
1552eac72dbSBarry Smith 
156dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal);
157dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
158dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
159dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
160dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
1614b0e389bSBarry Smith 
162dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt);
163dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal);
1649f236934SBarry Smith 
165dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP);
166dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
167dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
168dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1691d73ed98SBarry Smith 
170dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt);
171dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP);
1721d73ed98SBarry Smith 
173b49fd9e1SBarry Smith /*E
174b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
175b49fd9e1SBarry Smith 
176b49fd9e1SBarry Smith    Level: advanced
177b49fd9e1SBarry Smith 
178b49fd9e1SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1798c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
180b49fd9e1SBarry Smith 
181b49fd9e1SBarry Smith E*/
18278d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
1839dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[];
1841f7e983dSSatish Balay /*MC
1858c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
1868c5b8ba0SBarry Smith 
1878c5b8ba0SBarry Smith    Level: advanced
1888c5b8ba0SBarry Smith 
1898c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
1908c5b8ba0SBarry Smith 
1918c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1928c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
1938c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
1948c5b8ba0SBarry Smith M*/
1958c5b8ba0SBarry Smith 
1961f7e983dSSatish Balay /*MC
1978c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
1988c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
1998c5b8ba0SBarry Smith           poor orthogonality.
2008c5b8ba0SBarry Smith 
2018c5b8ba0SBarry Smith    Level: advanced
2028c5b8ba0SBarry Smith 
2038c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2048c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2058c5b8ba0SBarry Smith 
2068c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2078c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2088c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2098c5b8ba0SBarry Smith M*/
2108c5b8ba0SBarry Smith 
2111f7e983dSSatish Balay /*MC
2128c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2138c5b8ba0SBarry Smith 
2148c5b8ba0SBarry Smith    Level: advanced
2158c5b8ba0SBarry Smith 
2168c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2178c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2188c5b8ba0SBarry Smith 
2198c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2208c5b8ba0SBarry Smith 
2218c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2228c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2238c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2248c5b8ba0SBarry Smith M*/
2258c5b8ba0SBarry Smith 
226dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
22708480c60SBarry Smith 
228dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
229dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
230dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
231c38d4ed2SBarry Smith 
232dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal);
233dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*);
234dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*);
235121fd945SKris Buschelman 
236d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal);
237d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth);
238d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,int);
239d9492815SBarry Smith 
240dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP);
241dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2422eac72dbSBarry Smith 
243dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSingularValueMonitor(KSP,PetscInt,PetscReal,void *);
244dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultMonitor(KSP,PetscInt,PetscReal,void *);
245dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPTrueMonitor(KSP,PetscInt,PetscReal,void *);
246dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultSMonitor(KSP,PetscInt,PetscReal,void *);
247dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPVecViewMonitor(KSP,PetscInt,PetscReal,void *);
248dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESKrylovMonitor(KSP,PetscInt,PetscReal,void *);
24984cb2905SBarry Smith 
250dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec);
251dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*);
252dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
253dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
254c01c455dSBarry Smith 
255dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure);
256dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
257dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]);
258dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]);
2592dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]);
2601eb62cbbSBarry Smith 
261dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth);
262dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*);
263dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth);
264dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*);
2651f7f0c4fSBarry Smith 
266dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer);
2671eb62cbbSBarry Smith 
26828ce4d24SBarry Smith /*E
2698a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
2708a4b9c5cSBarry Smith        test routines.
2718a4b9c5cSBarry Smith 
2728a4b9c5cSBarry Smith    Level: advanced
2738a4b9c5cSBarry Smith 
2748a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
2758a4b9c5cSBarry Smith 
27694b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
2778a4b9c5cSBarry Smith           KSPSetConvergenceTest()
2788a4b9c5cSBarry Smith E*/
279954b3ec6SBarry Smith typedef enum {KSP_NO_NORM = 0,KSP_PRECONDITIONED_NORM = 1,KSP_UNPRECONDITIONED_NORM = 2,KSP_NATURAL_NORM = 3} KSPNormType;
2809dcbbd2bSBarry Smith extern const char *KSPNormTypes[];
2811f7e983dSSatish Balay /*MC
2828c5b8ba0SBarry Smith     KSP_NO_NORM - Do not compute a norm during the Krylov process. This will
2838c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
2848c5b8ba0SBarry Smith           be based on a norm of a residual etc.
2858c5b8ba0SBarry Smith 
2868c5b8ba0SBarry Smith    Level: advanced
2878c5b8ba0SBarry Smith 
288085a36d4SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
2898c5b8ba0SBarry Smith 
2908c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM
2918c5b8ba0SBarry Smith M*/
2928c5b8ba0SBarry Smith 
2931f7e983dSSatish Balay /*MC
2948c5b8ba0SBarry Smith     KSP_PRECONDITIONED_NORM - Compute the norm of the preconditioned residual and pass that to the
2958c5b8ba0SBarry Smith        convergence test routine.
2968c5b8ba0SBarry Smith 
2978c5b8ba0SBarry Smith    Level: advanced
2988c5b8ba0SBarry Smith 
2998c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
3008c5b8ba0SBarry Smith M*/
3018c5b8ba0SBarry Smith 
3021f7e983dSSatish Balay /*MC
3038c5b8ba0SBarry Smith     KSP_UNPRECONDITIONED_NORM - Compute the norm of the true residual (b - A*x) and pass that to the
3048c5b8ba0SBarry Smith        convergence test routine.
3058c5b8ba0SBarry Smith 
3068c5b8ba0SBarry Smith    Level: advanced
3078c5b8ba0SBarry Smith 
3088c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
3098c5b8ba0SBarry Smith M*/
3108c5b8ba0SBarry Smith 
3111f7e983dSSatish Balay /*MC
3128c5b8ba0SBarry Smith     KSP_NATURAL_NORM - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
3138c5b8ba0SBarry Smith        convergence test routine.
3148c5b8ba0SBarry Smith 
3158c5b8ba0SBarry Smith    Level: advanced
3168c5b8ba0SBarry Smith 
3178c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSPSetConvergenceTest()
3188c5b8ba0SBarry Smith M*/
3198c5b8ba0SBarry Smith 
320dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType);
3218a4b9c5cSBarry Smith 
3228a4b9c5cSBarry Smith /*E
32328ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
32428ce4d24SBarry Smith          have converged or diverged
32528ce4d24SBarry Smith 
32628ce4d24SBarry Smith    Level: beginner
32728ce4d24SBarry Smith 
32828ce4d24SBarry Smith    Notes: this must match finclude/petscksp.h
32928ce4d24SBarry Smith 
33086c02ca4SBarry Smith    Developer note: The string versions of these are in
3314b9ad928SBarry Smith      src/ksp/ksp/interface/itfunc.c called convergedreasons.
33250f1acb2SBarry Smith      If these enums are changed you must change those.
33386c02ca4SBarry Smith 
334c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
33528ce4d24SBarry Smith E*/
336d15094e1SBarry Smith typedef enum {/* converged */
337d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
338d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
339b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
340a4ce1dd3SMatthew Knepley               KSP_CONVERGED_STCG_NEG_CURVE     =  5,
341a4ce1dd3SMatthew Knepley               KSP_CONVERGED_STCG_CONSTRAINED   =  6,
342329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
343d15094e1SBarry Smith               /* diverged */
344b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
345d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
346d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
347d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
348b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
349b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
350b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
3516aee1118SBarry Smith               KSP_DIVERGED_NAN                 = -9,
3526aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
353d15094e1SBarry Smith 
354d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
3559dcbbd2bSBarry Smith extern const char **KSPConvergedReasons;
356d15094e1SBarry Smith 
357c838673bSBarry Smith /*MC
358c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
359c838673bSBarry Smith 
360c838673bSBarry Smith    Level: beginner
361c838673bSBarry Smith 
362c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
363c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
364c838673bSBarry Smith        2-norm of the residual for right preconditioning
365c838673bSBarry Smith 
366c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
367c838673bSBarry Smith 
368c838673bSBarry Smith M*/
369c838673bSBarry Smith 
370c838673bSBarry Smith /*MC
371c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
372c838673bSBarry Smith 
373c838673bSBarry Smith    Level: beginner
374c838673bSBarry Smith 
375c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
376c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
377c838673bSBarry Smith        2-norm of the residual for right preconditioning
378c838673bSBarry Smith 
379c838673bSBarry Smith    Level: beginner
380c838673bSBarry Smith 
381c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
382c838673bSBarry Smith 
383c838673bSBarry Smith M*/
384c838673bSBarry Smith 
385c838673bSBarry Smith /*MC
386c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
387c838673bSBarry Smith 
388c838673bSBarry Smith    Level: beginner
389c838673bSBarry Smith 
390c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
391c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
392c838673bSBarry Smith        2-norm of the residual for right preconditioning
393c838673bSBarry Smith 
394c838673bSBarry Smith    Level: beginner
395c838673bSBarry Smith 
396c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
397c838673bSBarry Smith 
398c838673bSBarry Smith M*/
399c838673bSBarry Smith 
400c838673bSBarry Smith /*MC
401c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
402c838673bSBarry Smith       reached
403c838673bSBarry Smith 
404c838673bSBarry Smith    Level: beginner
405c838673bSBarry Smith 
406c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
407c838673bSBarry Smith 
408c838673bSBarry Smith M*/
409c838673bSBarry Smith 
410c838673bSBarry Smith /*MC
411c838673bSBarry Smith      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of the
412c838673bSBarry Smith            preconditioner is applied.
413c838673bSBarry Smith 
414c838673bSBarry Smith 
415c838673bSBarry Smith    Level: beginner
416c838673bSBarry Smith 
417c838673bSBarry Smith 
418c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
419c838673bSBarry Smith 
420c838673bSBarry Smith M*/
421c838673bSBarry Smith 
422c838673bSBarry Smith /*MC
423c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
424c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
425c838673bSBarry Smith 
426c838673bSBarry Smith    Level: beginner
427c838673bSBarry Smith 
428c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
429c838673bSBarry Smith 
430c838673bSBarry Smith M*/
431c838673bSBarry Smith 
432c838673bSBarry Smith /*MC
433c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
434c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
435c838673bSBarry Smith 
436c838673bSBarry Smith 
437c838673bSBarry Smith    Level: beginner
438c838673bSBarry Smith 
439c838673bSBarry Smith 
440c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
441c838673bSBarry Smith 
442c838673bSBarry Smith M*/
443c838673bSBarry Smith 
444c838673bSBarry Smith /*MC
445c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
446c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
447c838673bSBarry Smith 
448c838673bSBarry Smith    Level: beginner
449c838673bSBarry Smith 
450c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
451c838673bSBarry Smith 
452c838673bSBarry Smith M*/
453c838673bSBarry Smith 
454c838673bSBarry Smith /*MC
455c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
456c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
457c838673bSBarry Smith         be positive definite
458c838673bSBarry Smith 
459c838673bSBarry Smith    Level: beginner
460c838673bSBarry Smith 
4612401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
462c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
463c838673bSBarry Smith 
464c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
465c838673bSBarry Smith 
466c838673bSBarry Smith M*/
467c838673bSBarry Smith 
468c838673bSBarry Smith /*MC
469c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
470c838673bSBarry Smith         while the KSPSolve() is still running.
471c838673bSBarry Smith 
472c838673bSBarry Smith    Level: beginner
473c838673bSBarry Smith 
474c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
475c838673bSBarry Smith 
476c838673bSBarry Smith M*/
477c838673bSBarry Smith 
478dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *);
479dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **);
480dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
481dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
482dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *);
483abef13c0SSatish Balay 
484dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *);
485d4fbbf0eSBarry Smith 
48628ce4d24SBarry Smith /*E
48728ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
48828ce4d24SBarry Smith 
48928ce4d24SBarry Smith    Level: beginner
49028ce4d24SBarry Smith 
49128ce4d24SBarry Smith .seealso: KSPCGSetType()
49228ce4d24SBarry Smith E*/
4939dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
4949dcbbd2bSBarry Smith extern const char *KSPCGTypes[];
49528ce4d24SBarry Smith 
496dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType);
497*0d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal);
498*0d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetQuadratic(KSP,PetscReal*);
499e559a7feSLois Curfman McInnes 
500dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP);
501dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP);
5023369ce9aSBarry Smith 
503dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
504dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitor(KSP,PetscInt,PetscReal,void*);
505dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorDestroy(PetscDrawLG);
506dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
507dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitor(KSP,PetscInt,PetscReal,void*);
508dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorDestroy(PetscDrawLG);
5092f2e5d10SKris Buschelman 
5100338f08bSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
51103919abeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
51203919abeSBarry Smith 
513e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
5142eac72dbSBarry Smith #endif
515