xref: /petsc/include/petscksp.h (revision 2f2e5d10202bd0387eb8c302ae412b45563fb849)
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 
9dfbe8321SBarry Smith EXTERN PetscErrorCode 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"
3582bf6240SBarry Smith #define KSPGMRES      "gmres"
3682bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
3782bf6240SBarry Smith #define KSPBCGS       "bcgs"
38f0037002Svictorle #define KSPBCGSL      "bcgsl"
3982bf6240SBarry Smith #define KSPCGS        "cgs"
4082bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
4182bf6240SBarry Smith #define KSPCR         "cr"
4282bf6240SBarry Smith #define KSPLSQR       "lsqr"
4382bf6240SBarry Smith #define KSPPREONLY    "preonly"
4482bf6240SBarry Smith #define KSPQCG        "qcg"
45c9cf9b11SBarry Smith #define KSPBICG       "bicg"
46c38d4ed2SBarry Smith #define KSPFGMRES     "fgmres"
47b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
4801c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
49bae98175SBarry Smith #define KSPLGMRES     "lgmres"
5049773a63SBarry Smith #define KSPType char*
512eac72dbSBarry Smith 
528ba1e511SMatthew Knepley /* Logging support */
536849ba73SBarry Smith extern PetscCookie KSP_COOKIE;
546849ba73SBarry Smith extern PetscEvent    KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve;
558ba1e511SMatthew Knepley 
56dfbe8321SBarry Smith EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
57dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetType(KSP,const KSPType);
58dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetUp(KSP);
59dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
60dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
61dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
62dfbe8321SBarry Smith EXTERN PetscErrorCode KSPDestroy(KSP);
632eac72dbSBarry Smith 
64b0a32e0cSBarry Smith extern PetscFList KSPList;
65dfbe8321SBarry Smith EXTERN PetscErrorCode KSPRegisterAll(const char[]);
66dfbe8321SBarry Smith EXTERN PetscErrorCode KSPRegisterDestroy(void);
672eac72dbSBarry Smith 
686849ba73SBarry Smith EXTERN PetscErrorCode 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 
119dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
120dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetPreconditionerSide(KSP,PCSide);
121dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetPreconditionerSide(KSP,PCSide*);
12213f74950SBarry Smith EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
12313f74950SBarry Smith EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
124dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscTruth);
125dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscTruth *);
126dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscTruth);
127dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscTruth*);
128dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscTruth);
129dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscTruth);
130dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
131dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
132dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
13313f74950SBarry Smith EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
134dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace);
135dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*);
1362eac72dbSBarry Smith 
137dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetPC(KSP,PC);
138dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
139aabeff55SBarry Smith 
14013f74950SBarry Smith EXTERN PetscErrorCode KSPSetMonitor(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
141dfbe8321SBarry Smith EXTERN PetscErrorCode KSPClearMonitor(KSP);
142dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
14313f74950SBarry Smith EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
14413f74950SBarry Smith EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);
1454b0e389bSBarry Smith 
1460e33f6ddSBarry Smith /* not sure where to put this */
147dfbe8321SBarry Smith EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
14813f74950SBarry Smith EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
14913f74950SBarry Smith EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
15069a612a9SBarry Smith EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
1512eac72dbSBarry Smith 
152dfbe8321SBarry Smith EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
153dfbe8321SBarry Smith EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
1542eac72dbSBarry Smith 
155dfbe8321SBarry Smith EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
156dfbe8321SBarry Smith EXTERN PetscErrorCode KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
157dfbe8321SBarry Smith EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
15813f74950SBarry Smith EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
15913f74950SBarry Smith EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
1604b0e389bSBarry Smith 
16113f74950SBarry Smith EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
162dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
1639f236934SBarry Smith 
164dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
16513f74950SBarry Smith EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
16613f74950SBarry Smith EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
16713f74950SBarry Smith EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1681d73ed98SBarry Smith 
16913f74950SBarry Smith EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
170dfbe8321SBarry Smith EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
1711d73ed98SBarry Smith 
172b49fd9e1SBarry Smith /*E
173b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
174b49fd9e1SBarry Smith 
175b49fd9e1SBarry Smith    Level: advanced
176b49fd9e1SBarry Smith 
177b49fd9e1SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1788c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
179b49fd9e1SBarry Smith 
180b49fd9e1SBarry Smith E*/
18178d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
1828c5b8ba0SBarry Smith 
1838c5b8ba0SBarry Smith /*M
1848c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
1858c5b8ba0SBarry Smith 
1868c5b8ba0SBarry Smith    Level: advanced
1878c5b8ba0SBarry Smith 
1888c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
1898c5b8ba0SBarry Smith 
1908c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1918c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
1928c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
1938c5b8ba0SBarry Smith M*/
1948c5b8ba0SBarry Smith 
1958c5b8ba0SBarry Smith /*M
1968c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
1978c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
1988c5b8ba0SBarry Smith           poor orthogonality.
1998c5b8ba0SBarry Smith 
2008c5b8ba0SBarry Smith    Level: advanced
2018c5b8ba0SBarry Smith 
2028c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2038c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2048c5b8ba0SBarry Smith 
2058c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2068c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2078c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2088c5b8ba0SBarry Smith M*/
2098c5b8ba0SBarry Smith 
2108c5b8ba0SBarry Smith /*M
2118c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2128c5b8ba0SBarry Smith 
2138c5b8ba0SBarry Smith    Level: advanced
2148c5b8ba0SBarry Smith 
2158c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2168c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2178c5b8ba0SBarry Smith 
2188c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2198c5b8ba0SBarry Smith 
2208c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2218c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2228c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2238c5b8ba0SBarry Smith M*/
2248c5b8ba0SBarry Smith 
225dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
22608480c60SBarry Smith 
22713f74950SBarry Smith EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
22813f74950SBarry Smith EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
22913f74950SBarry Smith EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
230c38d4ed2SBarry Smith 
231dfbe8321SBarry Smith EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
232dfbe8321SBarry Smith EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
233dfbe8321SBarry Smith EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
234121fd945SKris Buschelman 
235dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetFromOptions(KSP);
2366849ba73SBarry Smith EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2372eac72dbSBarry Smith 
23813f74950SBarry Smith EXTERN PetscErrorCode KSPSingularValueMonitor(KSP,PetscInt,PetscReal,void *);
23913f74950SBarry Smith EXTERN PetscErrorCode KSPDefaultMonitor(KSP,PetscInt,PetscReal,void *);
24013f74950SBarry Smith EXTERN PetscErrorCode KSPTrueMonitor(KSP,PetscInt,PetscReal,void *);
24113f74950SBarry Smith EXTERN PetscErrorCode KSPDefaultSMonitor(KSP,PetscInt,PetscReal,void *);
24213f74950SBarry Smith EXTERN PetscErrorCode KSPVecViewMonitor(KSP,PetscInt,PetscReal,void *);
24313f74950SBarry Smith EXTERN PetscErrorCode KSPGMRESKrylovMonitor(KSP,PetscInt,PetscReal,void *);
24484cb2905SBarry Smith 
245dfbe8321SBarry Smith EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
246dfbe8321SBarry Smith EXTERN PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*);
247dfbe8321SBarry Smith EXTERN PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
248c01c455dSBarry Smith 
249dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure);
25061f53403SBarry Smith EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
251dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
252dfbe8321SBarry Smith EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
253dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,char*[]);
2541eb62cbbSBarry Smith 
255dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscTruth);
256dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscTruth*);
257dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscTruth);
258dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscTruth*);
2591f7f0c4fSBarry Smith 
260dfbe8321SBarry Smith EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
2611eb62cbbSBarry Smith 
26228ce4d24SBarry Smith /*E
2638a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
2648a4b9c5cSBarry Smith        test routines.
2658a4b9c5cSBarry Smith 
2668a4b9c5cSBarry Smith    Level: advanced
2678a4b9c5cSBarry Smith 
2688a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
2698a4b9c5cSBarry Smith 
27094b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
2718a4b9c5cSBarry Smith           KSPSetConvergenceTest()
2728a4b9c5cSBarry Smith E*/
2738a4b9c5cSBarry Smith typedef enum {KSP_NO_NORM               = 0,
2748a4b9c5cSBarry Smith               KSP_PRECONDITIONED_NORM   = 1,
2758a4b9c5cSBarry Smith               KSP_UNPRECONDITIONED_NORM = 2,
2768a4b9c5cSBarry Smith               KSP_NATURAL_NORM          = 3} KSPNormType;
2778c5b8ba0SBarry Smith 
2788c5b8ba0SBarry Smith /*M
2798c5b8ba0SBarry Smith     KSP_NO_NORM - Do not compute a norm during the Krylov process. This will
2808c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
2818c5b8ba0SBarry Smith           be based on a norm of a residual etc.
2828c5b8ba0SBarry Smith 
2838c5b8ba0SBarry Smith    Level: advanced
2848c5b8ba0SBarry Smith 
2858c5b8ba0SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this is ignored
2868c5b8ba0SBarry Smith 
2878c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM
2888c5b8ba0SBarry Smith M*/
2898c5b8ba0SBarry Smith 
2908c5b8ba0SBarry Smith /*M
2918c5b8ba0SBarry Smith     KSP_PRECONDITIONED_NORM - Compute the norm of the preconditioned residual and pass that to the
2928c5b8ba0SBarry Smith        convergence test routine.
2938c5b8ba0SBarry Smith 
2948c5b8ba0SBarry Smith    Level: advanced
2958c5b8ba0SBarry Smith 
2968c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
2978c5b8ba0SBarry Smith M*/
2988c5b8ba0SBarry Smith 
2998c5b8ba0SBarry Smith /*M
3008c5b8ba0SBarry Smith     KSP_UNPRECONDITIONED_NORM - Compute the norm of the true residual (b - A*x) and pass that to the
3018c5b8ba0SBarry Smith        convergence test routine.
3028c5b8ba0SBarry Smith 
3038c5b8ba0SBarry Smith    Level: advanced
3048c5b8ba0SBarry Smith 
3058c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
3068c5b8ba0SBarry Smith M*/
3078c5b8ba0SBarry Smith 
3088c5b8ba0SBarry Smith /*M
3098c5b8ba0SBarry Smith     KSP_NATURAL_NORM - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
3108c5b8ba0SBarry Smith        convergence test routine.
3118c5b8ba0SBarry Smith 
3128c5b8ba0SBarry Smith    Level: advanced
3138c5b8ba0SBarry Smith 
3148c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSPSetConvergenceTest()
3158c5b8ba0SBarry Smith M*/
3168c5b8ba0SBarry Smith 
317dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
3188a4b9c5cSBarry Smith 
3198a4b9c5cSBarry Smith /*E
32028ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
32128ce4d24SBarry Smith          have converged or diverged
32228ce4d24SBarry Smith 
32328ce4d24SBarry Smith    Level: beginner
32428ce4d24SBarry Smith 
32528ce4d24SBarry Smith    Notes: this must match finclude/petscksp.h
32628ce4d24SBarry Smith 
32786c02ca4SBarry Smith    Developer note: The string versions of these are in
3284b9ad928SBarry Smith      src/ksp/ksp/interface/itfunc.c called convergedreasons.
32950f1acb2SBarry Smith      If these enums are changed you must change those.
33086c02ca4SBarry Smith 
331c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
33228ce4d24SBarry Smith E*/
333d15094e1SBarry Smith typedef enum {/* converged */
334d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
335d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
336b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3378d9717b1SSatish Balay               KSP_CONVERGED_QCG_NEG_CURVE      =  5,
338329f5518SBarry Smith               KSP_CONVERGED_QCG_CONSTRAINED    =  6,
339329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
340d15094e1SBarry Smith               /* diverged */
341b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
342d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
343d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
344d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
345b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
346b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
347b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
3486aee1118SBarry Smith               KSP_DIVERGED_NAN                 = -9,
3496aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
350d15094e1SBarry Smith 
351d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
352d15094e1SBarry Smith 
353c838673bSBarry Smith /*MC
354c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
355c838673bSBarry Smith 
356c838673bSBarry Smith    Level: beginner
357c838673bSBarry Smith 
358c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
359c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
360c838673bSBarry Smith        2-norm of the residual for right preconditioning
361c838673bSBarry Smith 
362c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
363c838673bSBarry Smith 
364c838673bSBarry Smith M*/
365c838673bSBarry Smith 
366c838673bSBarry Smith /*MC
367c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
368c838673bSBarry Smith 
369c838673bSBarry Smith    Level: beginner
370c838673bSBarry Smith 
371c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
372c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
373c838673bSBarry Smith        2-norm of the residual for right preconditioning
374c838673bSBarry Smith 
375c838673bSBarry Smith    Level: beginner
376c838673bSBarry Smith 
377c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
378c838673bSBarry Smith 
379c838673bSBarry Smith M*/
380c838673bSBarry Smith 
381c838673bSBarry Smith /*MC
382c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
383c838673bSBarry Smith 
384c838673bSBarry Smith    Level: beginner
385c838673bSBarry Smith 
386c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
387c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
388c838673bSBarry Smith        2-norm of the residual for right preconditioning
389c838673bSBarry Smith 
390c838673bSBarry Smith    Level: beginner
391c838673bSBarry Smith 
392c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
393c838673bSBarry Smith 
394c838673bSBarry Smith M*/
395c838673bSBarry Smith 
396c838673bSBarry Smith /*MC
397c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
398c838673bSBarry Smith       reached
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_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of the
408c838673bSBarry Smith            preconditioner is applied.
409c838673bSBarry Smith 
410c838673bSBarry Smith 
411c838673bSBarry Smith    Level: beginner
412c838673bSBarry Smith 
413c838673bSBarry Smith 
414c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
415c838673bSBarry Smith 
416c838673bSBarry Smith M*/
417c838673bSBarry Smith 
418c838673bSBarry Smith /*MC
419c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
420c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
421c838673bSBarry Smith 
422c838673bSBarry Smith    Level: beginner
423c838673bSBarry Smith 
424c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
425c838673bSBarry Smith 
426c838673bSBarry Smith M*/
427c838673bSBarry Smith 
428c838673bSBarry Smith /*MC
429c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
430c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
431c838673bSBarry Smith 
432c838673bSBarry Smith 
433c838673bSBarry Smith    Level: beginner
434c838673bSBarry Smith 
435c838673bSBarry Smith 
436c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
437c838673bSBarry Smith 
438c838673bSBarry Smith M*/
439c838673bSBarry Smith 
440c838673bSBarry Smith /*MC
441c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
442c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
443c838673bSBarry Smith 
444c838673bSBarry Smith    Level: beginner
445c838673bSBarry Smith 
446c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
447c838673bSBarry Smith 
448c838673bSBarry Smith M*/
449c838673bSBarry Smith 
450c838673bSBarry Smith /*MC
451c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
452c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
453c838673bSBarry Smith         be positive definite
454c838673bSBarry Smith 
455c838673bSBarry Smith    Level: beginner
456c838673bSBarry Smith 
457c838673bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_icc_shift to force
458c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
459c838673bSBarry Smith 
460c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
461c838673bSBarry Smith 
462c838673bSBarry Smith M*/
463c838673bSBarry Smith 
464c838673bSBarry Smith /*MC
465c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
466c838673bSBarry Smith         while the KSPSolve() is still running.
467c838673bSBarry Smith 
468c838673bSBarry Smith    Level: beginner
469c838673bSBarry Smith 
470c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
471c838673bSBarry Smith 
472c838673bSBarry Smith M*/
473c838673bSBarry Smith 
47413f74950SBarry Smith EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *);
475dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
47613f74950SBarry Smith EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
47713f74950SBarry Smith EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
478dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
479abef13c0SSatish Balay 
480dfbe8321SBarry Smith EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
481d4fbbf0eSBarry Smith 
48228ce4d24SBarry Smith /*E
48328ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
48428ce4d24SBarry Smith 
48528ce4d24SBarry Smith    Level: beginner
48628ce4d24SBarry Smith 
48728ce4d24SBarry Smith .seealso: KSPCGSetType()
48828ce4d24SBarry Smith E*/
4896d4a8577SBarry Smith typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType;
49028ce4d24SBarry Smith 
491dfbe8321SBarry Smith EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
492e559a7feSLois Curfman McInnes 
493dfbe8321SBarry Smith EXTERN PetscErrorCode PCPreSolve(PC,KSP);
494dfbe8321SBarry Smith EXTERN PetscErrorCode PCPostSolve(PC,KSP);
4953369ce9aSBarry Smith 
496dfbe8321SBarry Smith EXTERN PetscErrorCode KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
49713f74950SBarry Smith EXTERN PetscErrorCode KSPLGMonitor(KSP,PetscInt,PetscReal,void*);
498dfbe8321SBarry Smith EXTERN PetscErrorCode KSPLGMonitorDestroy(PetscDrawLG);
499dfbe8321SBarry Smith EXTERN PetscErrorCode KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
50013f74950SBarry Smith EXTERN PetscErrorCode KSPLGTrueMonitor(KSP,PetscInt,PetscReal,void*);
501dfbe8321SBarry Smith EXTERN PetscErrorCode KSPLGTrueMonitorDestroy(PetscDrawLG);
5021eb62cbbSBarry Smith 
503*2f2e5d10SKris Buschelman 
504e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
5052eac72dbSBarry Smith #endif
506