xref: /petsc/include/petscksp.h (revision 2dcb1b2a648e6798dfef446393ab7396b79ac556)
1f26ada1bSBarry Smith /*
2f26ada1bSBarry Smith    Defines the interface functions for the Krylov subspace accelerators.
3f26ada1bSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCKSP_H
50a835dfdSSatish Balay #define __PETSCKSP_H
60a835dfdSSatish Balay #include "petscpc.h"
7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
82eac72dbSBarry Smith 
9dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitializePackage(const char[]);
101dbb0a54SBarry Smith 
1128ce4d24SBarry Smith /*S
1228ce4d24SBarry Smith      KSP - Abstract PETSc object that manages all Krylov methods
1328ce4d24SBarry Smith 
1428ce4d24SBarry Smith    Level: beginner
1528ce4d24SBarry Smith 
1628ce4d24SBarry Smith   Concepts: Krylov methods
1728ce4d24SBarry Smith 
1894b7f48cSBarry Smith .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP
1928ce4d24SBarry Smith S*/
20e2a1c21fSSatish Balay typedef struct _p_KSP*     KSP;
212eac72dbSBarry Smith 
2228ce4d24SBarry Smith /*E
2328ce4d24SBarry Smith     KSPType - String with the name of a PETSc Krylov method or the creation function
2428ce4d24SBarry Smith        with an optional dynamic library name, for example
2528ce4d24SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mykspcreate()
2628ce4d24SBarry Smith 
2728ce4d24SBarry Smith    Level: beginner
2828ce4d24SBarry Smith 
2928ce4d24SBarry Smith .seealso: KSPSetType(), KSP
3028ce4d24SBarry Smith E*/
3182bf6240SBarry Smith #define KSPRICHARDSON "richardson"
3282bf6240SBarry Smith #define KSPCHEBYCHEV  "chebychev"
3382bf6240SBarry Smith #define KSPCG         "cg"
34df2a969aSvictorle #define   KSPCGNE       "cgne"
3582bf6240SBarry Smith #define KSPGMRES      "gmres"
369dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
379dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
3882bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
3982bf6240SBarry Smith #define KSPBCGS       "bcgs"
40f0037002Svictorle #define KSPBCGSL      "bcgsl"
4182bf6240SBarry Smith #define KSPCGS        "cgs"
4282bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
4382bf6240SBarry Smith #define KSPCR         "cr"
4482bf6240SBarry Smith #define KSPLSQR       "lsqr"
4582bf6240SBarry Smith #define KSPPREONLY    "preonly"
4682bf6240SBarry Smith #define KSPQCG        "qcg"
47c9cf9b11SBarry Smith #define KSPBICG       "bicg"
48b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
4901c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
5049773a63SBarry Smith #define KSPType char*
512eac72dbSBarry Smith 
528ba1e511SMatthew Knepley /* Logging support */
53dba47a55SKris Buschelman extern PetscCookie PETSCKSP_DLLEXPORT KSP_COOKIE;
546849ba73SBarry Smith extern PetscEvent    KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve;
558ba1e511SMatthew Knepley 
56dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *);
57dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,const KSPType);
58dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP);
59dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP);
60dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec);
61dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec);
62dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP);
632eac72dbSBarry Smith 
64b0a32e0cSBarry Smith extern PetscFList KSPList;
65dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]);
66dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void);
67dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));
6830de9b25SBarry Smith 
6930de9b25SBarry Smith /*MC
7030de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
7130de9b25SBarry Smith 
7230de9b25SBarry Smith    Synopsis:
73d360dc6fSBarry Smith    PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP))
7430de9b25SBarry Smith 
7530de9b25SBarry Smith    Not Collective
7630de9b25SBarry Smith 
7730de9b25SBarry Smith    Input Parameters:
7830de9b25SBarry Smith +  name_solver - name of a new user-defined solver
7930de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
8030de9b25SBarry Smith .  name_create - name of routine to create method context
8130de9b25SBarry Smith -  routine_create - routine to create method context
8230de9b25SBarry Smith 
8330de9b25SBarry Smith    Notes:
8430de9b25SBarry Smith    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
8530de9b25SBarry Smith 
8630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
8730de9b25SBarry Smith    is ignored.
8830de9b25SBarry Smith 
8930de9b25SBarry Smith    Sample usage:
9030de9b25SBarry Smith .vb
9130de9b25SBarry Smith    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
9230de9b25SBarry Smith                "MySolverCreate",MySolverCreate);
9330de9b25SBarry Smith .ve
9430de9b25SBarry Smith 
9530de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
9630de9b25SBarry Smith $     KSPSetType(ksp,"my_solver")
9730de9b25SBarry Smith    or at runtime via the option
9830de9b25SBarry Smith $     -ksp_type my_solver
9930de9b25SBarry Smith 
10030de9b25SBarry Smith    Level: advanced
10130de9b25SBarry Smith 
102ab901514SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
10330de9b25SBarry Smith           and others of the form ${any_environmental_variable} occuring in pathname will be
10430de9b25SBarry Smith           replaced with appropriate values.
10530de9b25SBarry Smith          If your function is not being put into a shared library then use KSPRegister() instead
10630de9b25SBarry Smith 
10730de9b25SBarry Smith .keywords: KSP, register
10830de9b25SBarry Smith 
10930de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy()
11030de9b25SBarry Smith 
11130de9b25SBarry Smith M*/
112aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
113f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
1146df38c32SLois Curfman McInnes #else
115f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
1166df38c32SLois Curfman McInnes #endif
11782bf6240SBarry Smith 
118dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,KSPType *);
119dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPreconditionerSide(KSP,PCSide);
120dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPreconditionerSide(KSP,PCSide*);
121dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
122dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
123dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscTruth);
124dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscTruth *);
125dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscTruth);
126dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscTruth*);
127dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth);
128dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth);
129dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *);
130dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *);
131dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*);
132dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*);
133dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace);
134dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*);
1352eac72dbSBarry Smith 
136dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC);
137dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*);
138aabeff55SBarry Smith 
139dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetMonitor(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
140dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPClearMonitor(KSP);
141dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **);
142dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
143dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);
1444b0e389bSBarry Smith 
1450e33f6ddSBarry Smith /* not sure where to put this */
146dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*);
147dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
148dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
149dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
1502eac72dbSBarry Smith 
151dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *);
152dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *);
1532eac72dbSBarry Smith 
154dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal);
155dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
156dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
157dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
158dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
1594b0e389bSBarry Smith 
160dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt);
161dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal);
1629f236934SBarry Smith 
163dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP);
164dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
165dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
166dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1671d73ed98SBarry Smith 
168dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt);
169dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT 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;
1819dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[];
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 
224dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
22508480c60SBarry Smith 
226dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
227dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
228dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
229c38d4ed2SBarry Smith 
230dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal);
231dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*);
232dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*);
233121fd945SKris Buschelman 
234dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP);
235dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2362eac72dbSBarry Smith 
237dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSingularValueMonitor(KSP,PetscInt,PetscReal,void *);
238dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultMonitor(KSP,PetscInt,PetscReal,void *);
239dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPTrueMonitor(KSP,PetscInt,PetscReal,void *);
240dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultSMonitor(KSP,PetscInt,PetscReal,void *);
241dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPVecViewMonitor(KSP,PetscInt,PetscReal,void *);
242dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESKrylovMonitor(KSP,PetscInt,PetscReal,void *);
24384cb2905SBarry Smith 
244dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec);
245dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*);
246dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
247dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
248c01c455dSBarry Smith 
249dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure);
250dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
251dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]);
252dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]);
253*2dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]);
2541eb62cbbSBarry Smith 
255dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth);
256dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*);
257dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth);
258dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*);
2591f7f0c4fSBarry Smith 
260dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT 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;
2779dcbbd2bSBarry Smith extern const char *KSPNormTypes[];
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 
285085a36d4SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this option 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 
317dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT 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;
3529dcbbd2bSBarry Smith extern const char **KSPConvergedReasons;
353d15094e1SBarry Smith 
354c838673bSBarry Smith /*MC
355c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
356c838673bSBarry Smith 
357c838673bSBarry Smith    Level: beginner
358c838673bSBarry Smith 
359c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
360c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
361c838673bSBarry Smith        2-norm of the residual for right preconditioning
362c838673bSBarry Smith 
363c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
364c838673bSBarry Smith 
365c838673bSBarry Smith M*/
366c838673bSBarry Smith 
367c838673bSBarry Smith /*MC
368c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
369c838673bSBarry Smith 
370c838673bSBarry Smith    Level: beginner
371c838673bSBarry Smith 
372c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
373c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
374c838673bSBarry Smith        2-norm of the residual for right preconditioning
375c838673bSBarry Smith 
376c838673bSBarry Smith    Level: beginner
377c838673bSBarry Smith 
378c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
379c838673bSBarry Smith 
380c838673bSBarry Smith M*/
381c838673bSBarry Smith 
382c838673bSBarry Smith /*MC
383c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
384c838673bSBarry Smith 
385c838673bSBarry Smith    Level: beginner
386c838673bSBarry Smith 
387c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
388c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
389c838673bSBarry Smith        2-norm of the residual for right preconditioning
390c838673bSBarry Smith 
391c838673bSBarry Smith    Level: beginner
392c838673bSBarry Smith 
393c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
394c838673bSBarry Smith 
395c838673bSBarry Smith M*/
396c838673bSBarry Smith 
397c838673bSBarry Smith /*MC
398c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
399c838673bSBarry Smith       reached
400c838673bSBarry Smith 
401c838673bSBarry Smith    Level: beginner
402c838673bSBarry Smith 
403c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
404c838673bSBarry Smith 
405c838673bSBarry Smith M*/
406c838673bSBarry Smith 
407c838673bSBarry Smith /*MC
408c838673bSBarry Smith      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of the
409c838673bSBarry Smith            preconditioner is applied.
410c838673bSBarry Smith 
411c838673bSBarry Smith 
412c838673bSBarry Smith    Level: beginner
413c838673bSBarry Smith 
414c838673bSBarry Smith 
415c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
416c838673bSBarry Smith 
417c838673bSBarry Smith M*/
418c838673bSBarry Smith 
419c838673bSBarry Smith /*MC
420c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
421c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
422c838673bSBarry Smith 
423c838673bSBarry Smith    Level: beginner
424c838673bSBarry Smith 
425c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
426c838673bSBarry Smith 
427c838673bSBarry Smith M*/
428c838673bSBarry Smith 
429c838673bSBarry Smith /*MC
430c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
431c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
432c838673bSBarry Smith 
433c838673bSBarry Smith 
434c838673bSBarry Smith    Level: beginner
435c838673bSBarry Smith 
436c838673bSBarry Smith 
437c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
438c838673bSBarry Smith 
439c838673bSBarry Smith M*/
440c838673bSBarry Smith 
441c838673bSBarry Smith /*MC
442c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
443c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
444c838673bSBarry Smith 
445c838673bSBarry Smith    Level: beginner
446c838673bSBarry Smith 
447c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
448c838673bSBarry Smith 
449c838673bSBarry Smith M*/
450c838673bSBarry Smith 
451c838673bSBarry Smith /*MC
452c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
453c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
454c838673bSBarry Smith         be positive definite
455c838673bSBarry Smith 
456c838673bSBarry Smith    Level: beginner
457c838673bSBarry Smith 
458c838673bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_icc_shift to force
459c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
460c838673bSBarry Smith 
461c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
462c838673bSBarry Smith 
463c838673bSBarry Smith M*/
464c838673bSBarry Smith 
465c838673bSBarry Smith /*MC
466c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
467c838673bSBarry Smith         while the KSPSolve() is still running.
468c838673bSBarry Smith 
469c838673bSBarry Smith    Level: beginner
470c838673bSBarry Smith 
471c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
472c838673bSBarry Smith 
473c838673bSBarry Smith M*/
474c838673bSBarry Smith 
475dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *);
476dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **);
477dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
478dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
479dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *);
480abef13c0SSatish Balay 
481dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *);
482d4fbbf0eSBarry Smith 
48328ce4d24SBarry Smith /*E
48428ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
48528ce4d24SBarry Smith 
48628ce4d24SBarry Smith    Level: beginner
48728ce4d24SBarry Smith 
48828ce4d24SBarry Smith .seealso: KSPCGSetType()
48928ce4d24SBarry Smith E*/
4909dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
4919dcbbd2bSBarry Smith extern const char *KSPCGTypes[];
49228ce4d24SBarry Smith 
493dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType);
494e559a7feSLois Curfman McInnes 
495dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP);
496dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP);
4973369ce9aSBarry Smith 
498dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
499dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitor(KSP,PetscInt,PetscReal,void*);
500dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorDestroy(PetscDrawLG);
501dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
502dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitor(KSP,PetscInt,PetscReal,void*);
503dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorDestroy(PetscDrawLG);
5042f2e5d10SKris Buschelman 
50503919abeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostPreSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
50603919abeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
50703919abeSBarry Smith 
508e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
5092eac72dbSBarry Smith #endif
510