xref: /petsc/include/petscksp.h (revision b3cc6726fbddf13cd805f591b2b0fc9b163819ab)
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 
91dbb0a54SBarry Smith EXTERN int 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"
3882bf6240SBarry Smith #define KSPCGS        "cgs"
3982bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
4082bf6240SBarry Smith #define KSPCR         "cr"
4182bf6240SBarry Smith #define KSPLSQR       "lsqr"
4282bf6240SBarry Smith #define KSPPREONLY    "preonly"
4382bf6240SBarry Smith #define KSPQCG        "qcg"
44c9cf9b11SBarry Smith #define KSPBICG       "bicg"
45c38d4ed2SBarry Smith #define KSPFGMRES     "fgmres"
46b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
4701c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
48bae98175SBarry Smith #define KSPLGMRES     "lgmres"
4949773a63SBarry Smith #define KSPType char*
502eac72dbSBarry Smith 
518ba1e511SMatthew Knepley /* Logging support */
528ba1e511SMatthew Knepley extern int KSP_COOKIE;
53d5ba7fb7SMatthew Knepley extern int KSP_GMRESOrthogonalization;
541dbb0a54SBarry Smith extern int KSP_SetUp, KSP_Solve;
558ba1e511SMatthew Knepley 
56ca44d042SBarry Smith EXTERN int KSPCreate(MPI_Comm,KSP *);
570e33f6ddSBarry Smith EXTERN int KSPSetType(KSP,const KSPType);
58ca44d042SBarry Smith EXTERN int KSPSetUp(KSP);
59b5d62d44SBarry Smith EXTERN int KSPSetUpOnBlocks(KSP);
6023ce1328SBarry Smith EXTERN int KSPSolve(KSP,Vec,Vec);
6123ce1328SBarry Smith EXTERN int KSPSolveTranspose(KSP,Vec,Vec);
62ca44d042SBarry Smith EXTERN int KSPDestroy(KSP);
632eac72dbSBarry Smith 
64b0a32e0cSBarry Smith extern PetscFList KSPList;
651836bdbcSSatish Balay EXTERN int KSPRegisterAll(const char[]);
66ca44d042SBarry Smith EXTERN int KSPRegisterDestroy(void);
672eac72dbSBarry Smith 
6879cabe49SKris Buschelman EXTERN int KSPRegister(const char[],const char[],const char[],int(*)(KSP));
6930de9b25SBarry Smith 
7030de9b25SBarry Smith /*MC
7130de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
7230de9b25SBarry Smith 
7330de9b25SBarry Smith    Synopsis:
7430de9b25SBarry Smith    int KSPRegisterDynamic(char *name_solver,char *path,char *name_create,int (*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 
119ca44d042SBarry Smith EXTERN int KSPGetType(KSP,KSPType *);
120ca44d042SBarry Smith EXTERN int KSPSetPreconditionerSide(KSP,PCSide);
121ca44d042SBarry Smith EXTERN int KSPGetPreconditionerSide(KSP,PCSide*);
12287828ca2SBarry Smith EXTERN int KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,int*);
12387828ca2SBarry Smith EXTERN int KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,int);
1243a7fca6bSBarry Smith EXTERN int KSPSetInitialGuessNonzero(KSP,PetscTruth);
125ca44d042SBarry Smith EXTERN int KSPGetInitialGuessNonzero(KSP,PetscTruth *);
1266ab422ffSBarry Smith EXTERN int KSPSetInitialGuessKnoll(KSP,PetscTruth);
1276ab422ffSBarry Smith EXTERN int KSPGetInitialGuessKnoll(KSP,PetscTruth*);
1283a7fca6bSBarry Smith EXTERN int KSPSetComputeEigenvalues(KSP,PetscTruth);
1293a7fca6bSBarry Smith EXTERN int KSPSetComputeSingularValues(KSP,PetscTruth);
130ca44d042SBarry Smith EXTERN int KSPGetRhs(KSP,Vec *);
131ca44d042SBarry Smith EXTERN int KSPGetSolution(KSP,Vec *);
13287828ca2SBarry Smith EXTERN int KSPGetResidualNorm(KSP,PetscReal*);
133ca44d042SBarry Smith EXTERN int KSPGetIterationNumber(KSP,int*);
134d8fd42c4SBarry Smith EXTERN int KSPSetNullSpace(KSP,MatNullSpace);
135d8fd42c4SBarry Smith EXTERN int KSPGetNullSpace(KSP,MatNullSpace*);
1362eac72dbSBarry Smith 
137ca44d042SBarry Smith EXTERN int KSPSetPC(KSP,PC);
138ca44d042SBarry Smith EXTERN int KSPGetPC(KSP,PC*);
139aabeff55SBarry Smith 
14087828ca2SBarry Smith EXTERN int KSPSetMonitor(KSP,int (*)(KSP,int,PetscReal,void*),void *,int (*)(void*));
141ca44d042SBarry Smith EXTERN int KSPClearMonitor(KSP);
142ca44d042SBarry Smith EXTERN int KSPGetMonitorContext(KSP,void **);
1431836bdbcSSatish Balay EXTERN int KSPGetResidualHistory(KSP,PetscReal*[],int *);
1441836bdbcSSatish Balay EXTERN int KSPSetResidualHistory(KSP,PetscReal[],int,PetscTruth);
1454b0e389bSBarry Smith 
1460e33f6ddSBarry Smith /* not sure where to put this */
1470e33f6ddSBarry Smith EXTERN int PCKSPGetKSP(PC,KSP*);
1480e33f6ddSBarry Smith EXTERN int PCBJacobiGetSubKSP(PC,int*,int*,KSP*[]);
1490e33f6ddSBarry Smith EXTERN int PCASMGetSubKSP(PC,int*,int*,KSP*[]);
1502eac72dbSBarry Smith 
151ca44d042SBarry Smith EXTERN int KSPBuildSolution(KSP,Vec,Vec *);
152ca44d042SBarry Smith EXTERN int KSPBuildResidual(KSP,Vec,Vec,Vec *);
1532eac72dbSBarry Smith 
15487828ca2SBarry Smith EXTERN int KSPRichardsonSetScale(KSP,PetscReal);
15587828ca2SBarry Smith EXTERN int KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
15687828ca2SBarry Smith EXTERN int KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
15787828ca2SBarry Smith EXTERN int KSPComputeEigenvalues(KSP,int,PetscReal*,PetscReal*,int *);
15887828ca2SBarry Smith EXTERN int KSPComputeEigenvaluesExplicitly(KSP,int,PetscReal*,PetscReal*);
1594b0e389bSBarry Smith 
160a24cce64SKris Buschelman EXTERN int KSPGMRESSetRestart(KSP, int);
161a24cce64SKris Buschelman EXTERN int KSPGMRESSetHapTol(KSP,PetscReal);
1629f236934SBarry Smith 
163ca44d042SBarry Smith EXTERN int KSPGMRESSetPreAllocateVectors(KSP);
164ca44d042SBarry Smith EXTERN int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int));
165ca44d042SBarry Smith EXTERN int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int);
166b49fd9e1SBarry Smith EXTERN int KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,int);
1671d73ed98SBarry Smith 
1681d73ed98SBarry Smith EXTERN int KSPLGMRESSetAugDim(KSP,int);
1691d73ed98SBarry Smith EXTERN int KSPLGMRESSetConstant(KSP);
1701d73ed98SBarry Smith 
171b49fd9e1SBarry Smith /*E
172b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
173b49fd9e1SBarry Smith 
174b49fd9e1SBarry Smith    Level: advanced
175b49fd9e1SBarry Smith 
176b49fd9e1SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1778c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
178b49fd9e1SBarry Smith 
179b49fd9e1SBarry Smith E*/
18078d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
1818c5b8ba0SBarry Smith 
1828c5b8ba0SBarry Smith /*M
1838c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
1848c5b8ba0SBarry Smith 
1858c5b8ba0SBarry Smith    Level: advanced
1868c5b8ba0SBarry Smith 
1878c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
1888c5b8ba0SBarry Smith 
1898c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1908c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
1918c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
1928c5b8ba0SBarry Smith M*/
1938c5b8ba0SBarry Smith 
1948c5b8ba0SBarry Smith /*M
1958c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
1968c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
1978c5b8ba0SBarry Smith           poor orthogonality.
1988c5b8ba0SBarry Smith 
1998c5b8ba0SBarry Smith    Level: advanced
2008c5b8ba0SBarry Smith 
2018c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2028c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2038c5b8ba0SBarry Smith 
2048c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2058c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2068c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2078c5b8ba0SBarry Smith M*/
2088c5b8ba0SBarry Smith 
2098c5b8ba0SBarry Smith /*M
2108c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2118c5b8ba0SBarry Smith 
2128c5b8ba0SBarry Smith    Level: advanced
2138c5b8ba0SBarry Smith 
2148c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2158c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2168c5b8ba0SBarry Smith 
2178c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2188c5b8ba0SBarry Smith 
2198c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2208c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2218c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2228c5b8ba0SBarry Smith M*/
2238c5b8ba0SBarry Smith 
224b49fd9e1SBarry Smith EXTERN int KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
22508480c60SBarry Smith 
22687828ca2SBarry Smith EXTERN int KSPFGMRESModifyPCNoChange(KSP,int,int,PetscReal,void*);
22794b7f48cSBarry Smith EXTERN int KSPFGMRESModifyPCKSP(KSP,int,int,PetscReal,void*);
22887828ca2SBarry Smith EXTERN int KSPFGMRESSetModifyPC(KSP,int (*)(KSP,int,int,PetscReal,void*),void*,int(*)(void*));
229c38d4ed2SBarry Smith 
23087828ca2SBarry Smith EXTERN int KSPQCGSetTrustRegionRadius(KSP,PetscReal);
23187828ca2SBarry Smith EXTERN int KSPQCGGetQuadratic(KSP,PetscReal*);
23287828ca2SBarry Smith EXTERN int KSPQCGGetTrialStepNorm(KSP,PetscReal*);
233121fd945SKris Buschelman 
234ca44d042SBarry Smith EXTERN int KSPSetFromOptions(KSP);
235ca44d042SBarry Smith EXTERN int KSPAddOptionsChecker(int (*)(KSP));
2362eac72dbSBarry Smith 
23787828ca2SBarry Smith EXTERN int KSPSingularValueMonitor(KSP,int,PetscReal,void *);
23887828ca2SBarry Smith EXTERN int KSPDefaultMonitor(KSP,int,PetscReal,void *);
23987828ca2SBarry Smith EXTERN int KSPTrueMonitor(KSP,int,PetscReal,void *);
24087828ca2SBarry Smith EXTERN int KSPDefaultSMonitor(KSP,int,PetscReal,void *);
24187828ca2SBarry Smith EXTERN int KSPVecViewMonitor(KSP,int,PetscReal,void *);
24287828ca2SBarry Smith EXTERN int KSPGMRESKrylovMonitor(KSP,int,PetscReal,void *);
24384cb2905SBarry Smith 
244ca44d042SBarry Smith EXTERN int KSPUnwindPreconditioner(KSP,Vec,Vec);
245ca44d042SBarry Smith EXTERN int KSPDefaultBuildSolution(KSP,Vec,Vec*);
246ca44d042SBarry Smith EXTERN int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
247c01c455dSBarry Smith 
248f0821613SBarry Smith EXTERN int KSPSetOperators(KSP,Mat,Mat,MatStructure);
2491836bdbcSSatish Balay EXTERN int KSPSetOptionsPrefix(KSP,const char[]);
2501836bdbcSSatish Balay EXTERN int KSPAppendOptionsPrefix(KSP,const char[]);
2511836bdbcSSatish Balay EXTERN int KSPGetOptionsPrefix(KSP,char*[]);
2521eb62cbbSBarry Smith 
2531f7f0c4fSBarry Smith EXTERN int KSPSetDiagonalScale(KSP,PetscTruth);
2541f7f0c4fSBarry Smith EXTERN int KSPGetDiagonalScale(KSP,PetscTruth*);
2551f7f0c4fSBarry Smith EXTERN int KSPSetDiagonalScaleFix(KSP,PetscTruth);
2561f7f0c4fSBarry Smith EXTERN int KSPGetDiagonalScaleFix(KSP,PetscTruth*);
2571f7f0c4fSBarry Smith 
258b0a32e0cSBarry Smith EXTERN int KSPView(KSP,PetscViewer);
2591eb62cbbSBarry Smith 
26028ce4d24SBarry Smith /*E
2618a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
2628a4b9c5cSBarry Smith        test routines.
2638a4b9c5cSBarry Smith 
2648a4b9c5cSBarry Smith    Level: advanced
2658a4b9c5cSBarry Smith 
2668a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
2678a4b9c5cSBarry Smith 
26894b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
2698a4b9c5cSBarry Smith           KSPSetConvergenceTest()
2708a4b9c5cSBarry Smith E*/
2718a4b9c5cSBarry Smith typedef enum {KSP_NO_NORM               = 0,
2728a4b9c5cSBarry Smith               KSP_PRECONDITIONED_NORM   = 1,
2738a4b9c5cSBarry Smith               KSP_UNPRECONDITIONED_NORM = 2,
2748a4b9c5cSBarry Smith               KSP_NATURAL_NORM          = 3} KSPNormType;
2758c5b8ba0SBarry Smith 
2768c5b8ba0SBarry Smith /*M
2778c5b8ba0SBarry Smith     KSP_NO_NORM - Do not compute a norm during the Krylov process. This will
2788c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
2798c5b8ba0SBarry Smith           be based on a norm of a residual etc.
2808c5b8ba0SBarry Smith 
2818c5b8ba0SBarry Smith    Level: advanced
2828c5b8ba0SBarry Smith 
2838c5b8ba0SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this is ignored
2848c5b8ba0SBarry Smith 
2858c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM
2868c5b8ba0SBarry Smith M*/
2878c5b8ba0SBarry Smith 
2888c5b8ba0SBarry Smith /*M
2898c5b8ba0SBarry Smith     KSP_PRECONDITIONED_NORM - Compute the norm of the preconditioned residual and pass that to the
2908c5b8ba0SBarry Smith        convergence test routine.
2918c5b8ba0SBarry Smith 
2928c5b8ba0SBarry Smith    Level: advanced
2938c5b8ba0SBarry Smith 
2948c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
2958c5b8ba0SBarry Smith M*/
2968c5b8ba0SBarry Smith 
2978c5b8ba0SBarry Smith /*M
2988c5b8ba0SBarry Smith     KSP_UNPRECONDITIONED_NORM - Compute the norm of the true residual (b - A*x) and pass that to the
2998c5b8ba0SBarry Smith        convergence test routine.
3008c5b8ba0SBarry Smith 
3018c5b8ba0SBarry Smith    Level: advanced
3028c5b8ba0SBarry Smith 
3038c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest()
3048c5b8ba0SBarry Smith M*/
3058c5b8ba0SBarry Smith 
3068c5b8ba0SBarry Smith /*M
3078c5b8ba0SBarry Smith     KSP_NATURAL_NORM - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
3088c5b8ba0SBarry Smith        convergence test routine.
3098c5b8ba0SBarry Smith 
3108c5b8ba0SBarry Smith    Level: advanced
3118c5b8ba0SBarry Smith 
3128c5b8ba0SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSPSetConvergenceTest()
3138c5b8ba0SBarry Smith M*/
3148c5b8ba0SBarry Smith 
3158a4b9c5cSBarry Smith EXTERN int KSPSetNormType(KSP,KSPNormType);
3168a4b9c5cSBarry Smith 
3178a4b9c5cSBarry Smith /*E
31828ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
31928ce4d24SBarry Smith          have converged or diverged
32028ce4d24SBarry Smith 
32128ce4d24SBarry Smith    Level: beginner
32228ce4d24SBarry Smith 
32328ce4d24SBarry Smith    Notes: this must match finclude/petscksp.h
32428ce4d24SBarry Smith 
32586c02ca4SBarry Smith    Developer note: The string versions of these are in
3264b9ad928SBarry Smith      src/ksp/ksp/interface/itfunc.c called convergedreasons.
32786c02ca4SBarry Smith      If these enums are changed you much change those.
32886c02ca4SBarry Smith 
329718d5cbfSvictorle .seealso: KSPSolve(), KSPGetConvergedReason()
33028ce4d24SBarry Smith E*/
331d15094e1SBarry Smith typedef enum {/* converged */
332d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
333d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
334b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3358d9717b1SSatish Balay               KSP_CONVERGED_QCG_NEG_CURVE      =  5,
336329f5518SBarry Smith               KSP_CONVERGED_QCG_CONSTRAINED    =  6,
337329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
338d15094e1SBarry Smith               /* diverged */
339*b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
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