xref: /petsc/include/petscksp.h (revision 8c5b8ba0f41b759ad5a83c0b9a7f2092f2c8d324)
173f4d377SMatthew Knepley /* $Id: petscksp.h,v 1.107 2001/08/06 21:16:38 bsmith Exp $ */
2f26ada1bSBarry Smith /*
3f26ada1bSBarry Smith    Defines the interface functions for the Krylov subspace accelerators.
4f26ada1bSBarry Smith */
50a835dfdSSatish Balay #ifndef __PETSCKSP_H
60a835dfdSSatish Balay #define __PETSCKSP_H
70a835dfdSSatish Balay #include "petscpc.h"
8e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
92eac72dbSBarry Smith 
101dbb0a54SBarry Smith EXTERN int KSPInitializePackage(const char[]);
111dbb0a54SBarry Smith 
1228ce4d24SBarry Smith /*S
1328ce4d24SBarry Smith      KSP - Abstract PETSc object that manages all Krylov methods
1428ce4d24SBarry Smith 
1528ce4d24SBarry Smith    Level: beginner
1628ce4d24SBarry Smith 
1728ce4d24SBarry Smith   Concepts: Krylov methods
1828ce4d24SBarry Smith 
1994b7f48cSBarry Smith .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP
2028ce4d24SBarry Smith S*/
21e2a1c21fSSatish Balay typedef struct _p_KSP*     KSP;
222eac72dbSBarry Smith 
2328ce4d24SBarry Smith /*E
2428ce4d24SBarry Smith     KSPType - String with the name of a PETSc Krylov method or the creation function
2528ce4d24SBarry Smith        with an optional dynamic library name, for example
2628ce4d24SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mykspcreate()
2728ce4d24SBarry Smith 
2828ce4d24SBarry Smith    Level: beginner
2928ce4d24SBarry Smith 
3028ce4d24SBarry Smith .seealso: KSPSetType(), KSP
3128ce4d24SBarry Smith E*/
3282bf6240SBarry Smith #define KSPRICHARDSON "richardson"
3382bf6240SBarry Smith #define KSPCHEBYCHEV  "chebychev"
3482bf6240SBarry Smith #define KSPCG         "cg"
35df2a969aSvictorle #define KSPCGNE       "cgne"
3682bf6240SBarry Smith #define KSPGMRES      "gmres"
3782bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
3882bf6240SBarry Smith #define KSPBCGS       "bcgs"
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 */
538ba1e511SMatthew Knepley extern int KSP_COOKIE;
54d5ba7fb7SMatthew Knepley extern int KSP_GMRESOrthogonalization;
551dbb0a54SBarry Smith extern int KSP_SetUp, KSP_Solve;
568ba1e511SMatthew Knepley 
57ca44d042SBarry Smith EXTERN int KSPCreate(MPI_Comm,KSP *);
580e33f6ddSBarry Smith EXTERN int KSPSetType(KSP,const KSPType);
59ca44d042SBarry Smith EXTERN int KSPSetUp(KSP);
60b5d62d44SBarry Smith EXTERN int KSPSetUpOnBlocks(KSP);
61c293cc10SBarry Smith EXTERN int KSPSolve(KSP);
62c293cc10SBarry Smith EXTERN int KSPSolveTranspose(KSP);
63ca44d042SBarry Smith EXTERN int KSPDestroy(KSP);
642eac72dbSBarry Smith 
65b0a32e0cSBarry Smith extern PetscFList KSPList;
661836bdbcSSatish Balay EXTERN int KSPRegisterAll(const char[]);
67ca44d042SBarry Smith EXTERN int KSPRegisterDestroy(void);
682eac72dbSBarry Smith 
6979cabe49SKris Buschelman EXTERN int KSPRegister(const char[],const char[],const char[],int(*)(KSP));
7030de9b25SBarry Smith 
7130de9b25SBarry Smith /*MC
7230de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
7330de9b25SBarry Smith 
7430de9b25SBarry Smith    Synopsis:
7530de9b25SBarry Smith    int KSPRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(KSP))
7630de9b25SBarry Smith 
7730de9b25SBarry Smith    Not Collective
7830de9b25SBarry Smith 
7930de9b25SBarry Smith    Input Parameters:
8030de9b25SBarry Smith +  name_solver - name of a new user-defined solver
8130de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
8230de9b25SBarry Smith .  name_create - name of routine to create method context
8330de9b25SBarry Smith -  routine_create - routine to create method context
8430de9b25SBarry Smith 
8530de9b25SBarry Smith    Notes:
8630de9b25SBarry Smith    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
8730de9b25SBarry Smith 
8830de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
8930de9b25SBarry Smith    is ignored.
9030de9b25SBarry Smith 
9130de9b25SBarry Smith    Sample usage:
9230de9b25SBarry Smith .vb
9330de9b25SBarry Smith    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
9430de9b25SBarry Smith                "MySolverCreate",MySolverCreate);
9530de9b25SBarry Smith .ve
9630de9b25SBarry Smith 
9730de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
9830de9b25SBarry Smith $     KSPSetType(ksp,"my_solver")
9930de9b25SBarry Smith    or at runtime via the option
10030de9b25SBarry Smith $     -ksp_type my_solver
10130de9b25SBarry Smith 
10230de9b25SBarry Smith    Level: advanced
10330de9b25SBarry Smith 
10430de9b25SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT},
10530de9b25SBarry Smith           and others of the form ${any_environmental_variable} occuring in pathname will be
10630de9b25SBarry Smith           replaced with appropriate values.
10730de9b25SBarry Smith          If your function is not being put into a shared library then use KSPRegister() instead
10830de9b25SBarry Smith 
10930de9b25SBarry Smith .keywords: KSP, register
11030de9b25SBarry Smith 
11130de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy()
11230de9b25SBarry Smith 
11330de9b25SBarry Smith M*/
114aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
115f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
1166df38c32SLois Curfman McInnes #else
117f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
1186df38c32SLois Curfman McInnes #endif
11982bf6240SBarry Smith 
120ca44d042SBarry Smith EXTERN int KSPGetType(KSP,KSPType *);
121ca44d042SBarry Smith EXTERN int KSPSetPreconditionerSide(KSP,PCSide);
122ca44d042SBarry Smith EXTERN int KSPGetPreconditionerSide(KSP,PCSide*);
12387828ca2SBarry Smith EXTERN int KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,int*);
12487828ca2SBarry Smith EXTERN int KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,int);
1253a7fca6bSBarry Smith EXTERN int KSPSetInitialGuessNonzero(KSP,PetscTruth);
126ca44d042SBarry Smith EXTERN int KSPGetInitialGuessNonzero(KSP,PetscTruth *);
1276ab422ffSBarry Smith EXTERN int KSPSetInitialGuessKnoll(KSP,PetscTruth);
1286ab422ffSBarry Smith EXTERN int KSPGetInitialGuessKnoll(KSP,PetscTruth*);
1293a7fca6bSBarry Smith EXTERN int KSPSetComputeEigenvalues(KSP,PetscTruth);
1303a7fca6bSBarry Smith EXTERN int KSPSetComputeSingularValues(KSP,PetscTruth);
131ca44d042SBarry Smith EXTERN int KSPSetRhs(KSP,Vec);
132ca44d042SBarry Smith EXTERN int KSPGetRhs(KSP,Vec *);
133ca44d042SBarry Smith EXTERN int KSPSetSolution(KSP,Vec);
134ca44d042SBarry Smith EXTERN int KSPGetSolution(KSP,Vec *);
13587828ca2SBarry Smith EXTERN int KSPGetResidualNorm(KSP,PetscReal*);
136ca44d042SBarry Smith EXTERN int KSPGetIterationNumber(KSP,int*);
1372eac72dbSBarry Smith 
138ca44d042SBarry Smith EXTERN int KSPSetPC(KSP,PC);
139ca44d042SBarry Smith EXTERN int KSPGetPC(KSP,PC*);
140aabeff55SBarry Smith 
14187828ca2SBarry Smith EXTERN int KSPSetMonitor(KSP,int (*)(KSP,int,PetscReal,void*),void *,int (*)(void*));
142ca44d042SBarry Smith EXTERN int KSPClearMonitor(KSP);
143ca44d042SBarry Smith EXTERN int KSPGetMonitorContext(KSP,void **);
1441836bdbcSSatish Balay EXTERN int KSPGetResidualHistory(KSP,PetscReal*[],int *);
1451836bdbcSSatish Balay EXTERN int KSPSetResidualHistory(KSP,PetscReal[],int,PetscTruth);
1464b0e389bSBarry Smith 
1470e33f6ddSBarry Smith /* not sure where to put this */
1480e33f6ddSBarry Smith EXTERN int PCKSPGetKSP(PC,KSP*);
1490e33f6ddSBarry Smith EXTERN int PCBJacobiGetSubKSP(PC,int*,int*,KSP*[]);
1500e33f6ddSBarry Smith EXTERN int PCASMGetSubKSP(PC,int*,int*,KSP*[]);
1512eac72dbSBarry Smith 
152ca44d042SBarry Smith EXTERN int KSPBuildSolution(KSP,Vec,Vec *);
153ca44d042SBarry Smith EXTERN int KSPBuildResidual(KSP,Vec,Vec,Vec *);
1542eac72dbSBarry Smith 
15587828ca2SBarry Smith EXTERN int KSPRichardsonSetScale(KSP,PetscReal);
15687828ca2SBarry Smith EXTERN int KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
15787828ca2SBarry Smith EXTERN int KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
15887828ca2SBarry Smith EXTERN int KSPComputeEigenvalues(KSP,int,PetscReal*,PetscReal*,int *);
15987828ca2SBarry Smith EXTERN int KSPComputeEigenvaluesExplicitly(KSP,int,PetscReal*,PetscReal*);
1604b0e389bSBarry Smith 
161a24cce64SKris Buschelman EXTERN int KSPGMRESSetRestart(KSP, int);
162a24cce64SKris Buschelman EXTERN int KSPGMRESSetHapTol(KSP,PetscReal);
1639f236934SBarry Smith 
164ca44d042SBarry Smith EXTERN int KSPGMRESSetPreAllocateVectors(KSP);
165ca44d042SBarry Smith EXTERN int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int));
166ca44d042SBarry Smith EXTERN int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int);
167b49fd9e1SBarry Smith EXTERN int KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,int);
1681d73ed98SBarry Smith 
1691d73ed98SBarry Smith EXTERN int KSPLGMRESSetAugDim(KSP,int);
1701d73ed98SBarry Smith EXTERN int 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(),
178*8c5b8ba0SBarry 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;
182*8c5b8ba0SBarry Smith 
183*8c5b8ba0SBarry Smith /*M
184*8c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
185*8c5b8ba0SBarry Smith 
186*8c5b8ba0SBarry Smith    Level: advanced
187*8c5b8ba0SBarry Smith 
188*8c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
189*8c5b8ba0SBarry Smith 
190*8c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
191*8c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
192*8c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
193*8c5b8ba0SBarry Smith M*/
194*8c5b8ba0SBarry Smith 
195*8c5b8ba0SBarry Smith /*M
196*8c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
197*8c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
198*8c5b8ba0SBarry Smith           poor orthogonality.
199*8c5b8ba0SBarry Smith 
200*8c5b8ba0SBarry Smith    Level: advanced
201*8c5b8ba0SBarry Smith 
202*8c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
203*8c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
204*8c5b8ba0SBarry Smith 
205*8c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
206*8c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
207*8c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
208*8c5b8ba0SBarry Smith M*/
209*8c5b8ba0SBarry Smith 
210*8c5b8ba0SBarry Smith /*M
211*8c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
212*8c5b8ba0SBarry Smith 
213*8c5b8ba0SBarry Smith    Level: advanced
214*8c5b8ba0SBarry Smith 
215*8c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
216*8c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
217*8c5b8ba0SBarry Smith 
218*8c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
219*8c5b8ba0SBarry Smith 
220*8c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
221*8c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
222*8c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
223*8c5b8ba0SBarry Smith M*/
224*8c5b8ba0SBarry Smith 
225b49fd9e1SBarry Smith EXTERN int KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
22608480c60SBarry Smith 
22787828ca2SBarry Smith EXTERN int KSPFGMRESModifyPCNoChange(KSP,int,int,PetscReal,void*);
22894b7f48cSBarry Smith EXTERN int KSPFGMRESModifyPCKSP(KSP,int,int,PetscReal,void*);
22987828ca2SBarry Smith EXTERN int KSPFGMRESSetModifyPC(KSP,int (*)(KSP,int,int,PetscReal,void*),void*,int(*)(void*));
230c38d4ed2SBarry Smith 
23187828ca2SBarry Smith EXTERN int KSPQCGSetTrustRegionRadius(KSP,PetscReal);
23287828ca2SBarry Smith EXTERN int KSPQCGGetQuadratic(KSP,PetscReal*);
23387828ca2SBarry Smith EXTERN int KSPQCGGetTrialStepNorm(KSP,PetscReal*);
234121fd945SKris Buschelman 
235ca44d042SBarry Smith EXTERN int KSPSetFromOptions(KSP);
236ca44d042SBarry Smith EXTERN int KSPAddOptionsChecker(int (*)(KSP));
2372eac72dbSBarry Smith 
23887828ca2SBarry Smith EXTERN int KSPSingularValueMonitor(KSP,int,PetscReal,void *);
23987828ca2SBarry Smith EXTERN int KSPDefaultMonitor(KSP,int,PetscReal,void *);
24087828ca2SBarry Smith EXTERN int KSPTrueMonitor(KSP,int,PetscReal,void *);
24187828ca2SBarry Smith EXTERN int KSPDefaultSMonitor(KSP,int,PetscReal,void *);
24287828ca2SBarry Smith EXTERN int KSPVecViewMonitor(KSP,int,PetscReal,void *);
24387828ca2SBarry Smith EXTERN int KSPGMRESKrylovMonitor(KSP,int,PetscReal,void *);
24484cb2905SBarry Smith 
245ca44d042SBarry Smith EXTERN int KSPUnwindPreconditioner(KSP,Vec,Vec);
246ca44d042SBarry Smith EXTERN int KSPDefaultBuildSolution(KSP,Vec,Vec*);
247ca44d042SBarry Smith EXTERN int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
248c01c455dSBarry Smith 
249f0821613SBarry Smith EXTERN int KSPSetOperators(KSP,Mat,Mat,MatStructure);
2501836bdbcSSatish Balay EXTERN int KSPSetOptionsPrefix(KSP,const char[]);
2511836bdbcSSatish Balay EXTERN int KSPAppendOptionsPrefix(KSP,const char[]);
2521836bdbcSSatish Balay EXTERN int KSPGetOptionsPrefix(KSP,char*[]);
2531eb62cbbSBarry Smith 
2541f7f0c4fSBarry Smith EXTERN int KSPSetDiagonalScale(KSP,PetscTruth);
2551f7f0c4fSBarry Smith EXTERN int KSPGetDiagonalScale(KSP,PetscTruth*);
2561f7f0c4fSBarry Smith EXTERN int KSPSetDiagonalScaleFix(KSP,PetscTruth);
2571f7f0c4fSBarry Smith EXTERN int KSPGetDiagonalScaleFix(KSP,PetscTruth*);
2581f7f0c4fSBarry Smith 
259b0a32e0cSBarry Smith EXTERN int KSPView(KSP,PetscViewer);
2601eb62cbbSBarry Smith 
26128ce4d24SBarry Smith /*E
2628a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
2638a4b9c5cSBarry Smith        test routines.
2648a4b9c5cSBarry Smith 
2658a4b9c5cSBarry Smith    Level: advanced
2668a4b9c5cSBarry Smith 
2678a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
2688a4b9c5cSBarry Smith 
26994b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
2708a4b9c5cSBarry Smith           KSPSetConvergenceTest()
2718a4b9c5cSBarry Smith E*/
2728a4b9c5cSBarry Smith typedef enum {KSP_NO_NORM               = 0,
2738a4b9c5cSBarry Smith               KSP_PRECONDITIONED_NORM   = 1,
2748a4b9c5cSBarry Smith               KSP_UNPRECONDITIONED_NORM = 2,
2758a4b9c5cSBarry Smith               KSP_NATURAL_NORM          = 3} KSPNormType;
276*8c5b8ba0SBarry Smith 
277*8c5b8ba0SBarry Smith /*M
278*8c5b8ba0SBarry Smith     KSP_NO_NORM - Do not compute a norm during the Krylov process. This will
279*8c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
280*8c5b8ba0SBarry Smith           be based on a norm of a residual etc.
281*8c5b8ba0SBarry Smith 
282*8c5b8ba0SBarry Smith    Level: advanced
283*8c5b8ba0SBarry Smith 
284*8c5b8ba0SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this is ignored
285*8c5b8ba0SBarry Smith 
286*8c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM
287*8c5b8ba0SBarry Smith M*/
288*8c5b8ba0SBarry Smith 
289*8c5b8ba0SBarry Smith /*M
290*8c5b8ba0SBarry Smith     KSP_PRECONDITIONED_NORM - Compute the norm of the preconditioned residual and pass that to the
291*8c5b8ba0SBarry Smith        convergence test routine.
292*8c5b8ba0SBarry Smith 
293*8c5b8ba0SBarry Smith    Level: advanced
294*8c5b8ba0SBarry Smith 
295*8c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
296*8c5b8ba0SBarry Smith M*/
297*8c5b8ba0SBarry Smith 
298*8c5b8ba0SBarry Smith /*M
299*8c5b8ba0SBarry Smith     KSP_UNPRECONDITIONED_NORM - Compute the norm of the true residual (b - A*x) and pass that to the
300*8c5b8ba0SBarry Smith        convergence test routine.
301*8c5b8ba0SBarry Smith 
302*8c5b8ba0SBarry Smith    Level: advanced
303*8c5b8ba0SBarry Smith 
304*8c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
305*8c5b8ba0SBarry Smith M*/
306*8c5b8ba0SBarry Smith 
307*8c5b8ba0SBarry Smith /*M
308*8c5b8ba0SBarry Smith     KSP_NATURAL_NORM - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
309*8c5b8ba0SBarry Smith        convergence test routine.
310*8c5b8ba0SBarry Smith 
311*8c5b8ba0SBarry Smith    Level: advanced
312*8c5b8ba0SBarry Smith 
313*8c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSPSetConvergenceTest()
314*8c5b8ba0SBarry Smith M*/
315*8c5b8ba0SBarry Smith 
3168a4b9c5cSBarry Smith EXTERN int KSPSetNormType(KSP,KSPNormType);
3178a4b9c5cSBarry Smith 
3188a4b9c5cSBarry Smith /*E
31928ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
32028ce4d24SBarry Smith          have converged or diverged
32128ce4d24SBarry Smith 
32228ce4d24SBarry Smith    Level: beginner
32328ce4d24SBarry Smith 
32428ce4d24SBarry Smith    Notes: this must match finclude/petscksp.h
32528ce4d24SBarry Smith 
32686c02ca4SBarry Smith    Developer note: The string versions of these are in
3274b9ad928SBarry Smith      src/ksp/ksp/interface/itfunc.c called convergedreasons.
32886c02ca4SBarry Smith      If these enums are changed you much change those.
32986c02ca4SBarry Smith 
330718d5cbfSvictorle .seealso: KSPSolve(), KSPGetConvergedReason()
33128ce4d24SBarry Smith E*/
332d15094e1SBarry Smith typedef enum {/* converged */
333d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
334d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
335b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3368d9717b1SSatish Balay               KSP_CONVERGED_QCG_NEG_CURVE      =  5,
337329f5518SBarry Smith               KSP_CONVERGED_QCG_CONSTRAINED    =  6,
338329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
339d15094e1SBarry Smith               /* diverged */
340d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
341d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
342d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
343b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
344b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
345b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
346d15094e1SBarry Smith 
347d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
348d15094e1SBarry Smith 
34987828ca2SBarry Smith EXTERN int KSPSetConvergenceTest(KSP,int (*)(KSP,int,PetscReal,KSPConvergedReason*,void*),void *);
350ca44d042SBarry Smith EXTERN int KSPGetConvergenceContext(KSP,void **);
35187828ca2SBarry Smith EXTERN int KSPDefaultConverged(KSP,int,PetscReal,KSPConvergedReason*,void *);
35287828ca2SBarry Smith EXTERN int KSPSkipConverged(KSP,int,PetscReal,KSPConvergedReason*,void *);
353ca44d042SBarry Smith EXTERN int KSPGetConvergedReason(KSP,KSPConvergedReason *);
354abef13c0SSatish Balay 
355ca44d042SBarry Smith EXTERN int KSPComputeExplicitOperator(KSP,Mat *);
356d4fbbf0eSBarry Smith 
35728ce4d24SBarry Smith /*E
35828ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
35928ce4d24SBarry Smith 
36028ce4d24SBarry Smith    Level: beginner
36128ce4d24SBarry Smith 
36228ce4d24SBarry Smith .seealso: KSPCGSetType()
36328ce4d24SBarry Smith E*/
3646d4a8577SBarry Smith typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType;
36528ce4d24SBarry Smith 
366ca44d042SBarry Smith EXTERN int KSPCGSetType(KSP,KSPCGType);
367e559a7feSLois Curfman McInnes 
368ca44d042SBarry Smith EXTERN int PCPreSolve(PC,KSP);
369ca44d042SBarry Smith EXTERN int PCPostSolve(PC,KSP);
3703369ce9aSBarry Smith 
3711836bdbcSSatish Balay EXTERN int KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
37287828ca2SBarry Smith EXTERN int KSPLGMonitor(KSP,int,PetscReal,void*);
373b0a32e0cSBarry Smith EXTERN int KSPLGMonitorDestroy(PetscDrawLG);
3741836bdbcSSatish Balay EXTERN int KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
37587828ca2SBarry Smith EXTERN int KSPLGTrueMonitor(KSP,int,PetscReal,void*);
376b0a32e0cSBarry Smith EXTERN int KSPLGTrueMonitorDestroy(PetscDrawLG);
3771eb62cbbSBarry Smith 
378e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
3792eac72dbSBarry Smith #endif
380