xref: /petsc/include/petscksp.h (revision 69a612a9be7fbcc0ab27d96ce41c4e3795ac7f2a)
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 
10330de9b25SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT},
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*);
122dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,int*);
123dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,int);
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*);
133dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetIterationNumber(KSP,int*);
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 
1406849ba73SBarry Smith EXTERN PetscErrorCode KSPSetMonitor(KSP,PetscErrorCode (*)(KSP,int,PetscReal,void*),void *,PetscErrorCode (*)(void*));
141dfbe8321SBarry Smith EXTERN PetscErrorCode KSPClearMonitor(KSP);
142dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
143dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],int *);
144dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],int,PetscTruth);
1454b0e389bSBarry Smith 
1460e33f6ddSBarry Smith /* not sure where to put this */
147dfbe8321SBarry Smith EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
148dfbe8321SBarry Smith EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,int*,int*,KSP*[]);
149dfbe8321SBarry Smith EXTERN PetscErrorCode PCASMGetSubKSP(PC,int*,int*,KSP*[]);
150*69a612a9SBarry 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*);
158dfbe8321SBarry Smith EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,int,PetscReal*,PetscReal*,int *);
159dfbe8321SBarry Smith EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,int,PetscReal*,PetscReal*);
1604b0e389bSBarry Smith 
161dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, int);
162dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
1639f236934SBarry Smith 
164dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
1656849ba73SBarry Smith EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,int));
166dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int);
167dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,int);
1681d73ed98SBarry Smith 
169dfbe8321SBarry Smith EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,int);
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 
227dfbe8321SBarry Smith EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,int,int,PetscReal,void*);
228dfbe8321SBarry Smith EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,int,int,PetscReal,void*);
2296849ba73SBarry Smith EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,int,int,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 
238dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSingularValueMonitor(KSP,int,PetscReal,void *);
239dfbe8321SBarry Smith EXTERN PetscErrorCode KSPDefaultMonitor(KSP,int,PetscReal,void *);
240dfbe8321SBarry Smith EXTERN PetscErrorCode KSPTrueMonitor(KSP,int,PetscReal,void *);
241dfbe8321SBarry Smith EXTERN PetscErrorCode KSPDefaultSMonitor(KSP,int,PetscReal,void *);
242dfbe8321SBarry Smith EXTERN PetscErrorCode KSPVecViewMonitor(KSP,int,PetscReal,void *);
243dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGMRESKrylovMonitor(KSP,int,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.
32986c02ca4SBarry Smith      If these enums are changed you much change those.
33086c02ca4SBarry Smith 
331718d5cbfSvictorle .seealso: KSPSolve(), KSPGetConvergedReason()
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,
348d15094e1SBarry Smith 
349d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
350d15094e1SBarry Smith 
3516849ba73SBarry Smith EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,int,PetscReal,KSPConvergedReason*,void*),void *);
352dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
353dfbe8321SBarry Smith EXTERN PetscErrorCode KSPDefaultConverged(KSP,int,PetscReal,KSPConvergedReason*,void *);
354dfbe8321SBarry Smith EXTERN PetscErrorCode KSPSkipConverged(KSP,int,PetscReal,KSPConvergedReason*,void *);
355dfbe8321SBarry Smith EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
356abef13c0SSatish Balay 
357dfbe8321SBarry Smith EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
358d4fbbf0eSBarry Smith 
35928ce4d24SBarry Smith /*E
36028ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
36128ce4d24SBarry Smith 
36228ce4d24SBarry Smith    Level: beginner
36328ce4d24SBarry Smith 
36428ce4d24SBarry Smith .seealso: KSPCGSetType()
36528ce4d24SBarry Smith E*/
3666d4a8577SBarry Smith typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType;
36728ce4d24SBarry Smith 
368dfbe8321SBarry Smith EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
369e559a7feSLois Curfman McInnes 
370dfbe8321SBarry Smith EXTERN PetscErrorCode PCPreSolve(PC,KSP);
371dfbe8321SBarry Smith EXTERN PetscErrorCode PCPostSolve(PC,KSP);
3723369ce9aSBarry Smith 
373dfbe8321SBarry Smith EXTERN PetscErrorCode KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
374dfbe8321SBarry Smith EXTERN PetscErrorCode KSPLGMonitor(KSP,int,PetscReal,void*);
375dfbe8321SBarry Smith EXTERN PetscErrorCode KSPLGMonitorDestroy(PetscDrawLG);
376dfbe8321SBarry Smith EXTERN PetscErrorCode KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
377dfbe8321SBarry Smith EXTERN PetscErrorCode KSPLGTrueMonitor(KSP,int,PetscReal,void*);
378dfbe8321SBarry Smith EXTERN PetscErrorCode KSPLGTrueMonitorDestroy(PetscDrawLG);
3791eb62cbbSBarry Smith 
380e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
3812eac72dbSBarry Smith #endif
382