xref: /petsc/include/petscksp.h (revision 1e7308c7d0b3c5375d581985f0b0fcd07f4ab00a)
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*/
31a313700dSBarry Smith #define KSPType char*
3282bf6240SBarry Smith #define KSPRICHARDSON "richardson"
3382bf6240SBarry Smith #define KSPCHEBYCHEV  "chebychev"
3482bf6240SBarry Smith #define KSPCG         "cg"
35df2a969aSvictorle #define   KSPCGNE       "cgne"
36fcae7a14Stmunson #define   KSPNASH       "nash"
3780e17431SMatthew Knepley #define   KSPSTCG       "stcg"
388031f4adStmunson #define   KSPGLTR       "gltr"
3982bf6240SBarry Smith #define KSPGMRES      "gmres"
409dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
419dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
4282bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
4382bf6240SBarry Smith #define KSPBCGS       "bcgs"
44e1c61ce8SBarry Smith #define KSPIBCGS        "ibcgs"
45f0037002Svictorle #define KSPBCGSL        "bcgsl"
4682bf6240SBarry Smith #define KSPCGS        "cgs"
4782bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
4882bf6240SBarry Smith #define KSPCR         "cr"
4982bf6240SBarry Smith #define KSPLSQR       "lsqr"
5082bf6240SBarry Smith #define KSPPREONLY    "preonly"
5182bf6240SBarry Smith #define KSPQCG        "qcg"
52c9cf9b11SBarry Smith #define KSPBICG       "bicg"
53b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
5401c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
5521356919SSatish Balay #define KSPLCD        "lcd"
562eac72dbSBarry Smith 
578ba1e511SMatthew Knepley /* Logging support */
58dba47a55SKris Buschelman extern PetscCookie PETSCKSP_DLLEXPORT KSP_COOKIE;
598ba1e511SMatthew Knepley 
60dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *);
61a313700dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,const KSPType);
62dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP);
63dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP);
64dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec);
65dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec);
66dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP);
672eac72dbSBarry Smith 
68b0a32e0cSBarry Smith extern PetscFList KSPList;
69dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]);
70dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void);
71dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));
7230de9b25SBarry Smith 
7330de9b25SBarry Smith /*MC
7430de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
7530de9b25SBarry Smith 
7630de9b25SBarry Smith    Synopsis:
77d360dc6fSBarry Smith    PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP))
7830de9b25SBarry Smith 
7930de9b25SBarry Smith    Not Collective
8030de9b25SBarry Smith 
8130de9b25SBarry Smith    Input Parameters:
8230de9b25SBarry Smith +  name_solver - name of a new user-defined solver
8330de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
8430de9b25SBarry Smith .  name_create - name of routine to create method context
8530de9b25SBarry Smith -  routine_create - routine to create method context
8630de9b25SBarry Smith 
8730de9b25SBarry Smith    Notes:
8830de9b25SBarry Smith    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
8930de9b25SBarry Smith 
9030de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
9130de9b25SBarry Smith    is ignored.
9230de9b25SBarry Smith 
9330de9b25SBarry Smith    Sample usage:
9430de9b25SBarry Smith .vb
9530de9b25SBarry Smith    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
9630de9b25SBarry Smith                "MySolverCreate",MySolverCreate);
9730de9b25SBarry Smith .ve
9830de9b25SBarry Smith 
9930de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
10030de9b25SBarry Smith $     KSPSetType(ksp,"my_solver")
10130de9b25SBarry Smith    or at runtime via the option
10230de9b25SBarry Smith $     -ksp_type my_solver
10330de9b25SBarry Smith 
10430de9b25SBarry Smith    Level: advanced
10530de9b25SBarry Smith 
106ab901514SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
10730de9b25SBarry Smith           and others of the form ${any_environmental_variable} occuring in pathname will be
10830de9b25SBarry Smith           replaced with appropriate values.
10930de9b25SBarry Smith          If your function is not being put into a shared library then use KSPRegister() instead
11030de9b25SBarry Smith 
11130de9b25SBarry Smith .keywords: KSP, register
11230de9b25SBarry Smith 
11330de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy()
11430de9b25SBarry Smith 
11530de9b25SBarry Smith M*/
116aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
117f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
1186df38c32SLois Curfman McInnes #else
119f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
1206df38c32SLois Curfman McInnes #endif
12182bf6240SBarry Smith 
122a313700dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,const KSPType *);
123dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPreconditionerSide(KSP,PCSide);
124dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPreconditionerSide(KSP,PCSide*);
125dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
126dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
127dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscTruth);
128dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscTruth *);
129dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscTruth);
130dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscTruth*);
131a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscTruth*);
132dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth);
133a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscTruth*);
134dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth);
135dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *);
136dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *);
137dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*);
138dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*);
139dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace);
140dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*);
141906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1422eac72dbSBarry Smith 
143dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC);
144dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*);
145aabeff55SBarry Smith 
146a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
147a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorCancel(KSP);
148dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **);
149dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
150dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);
1514b0e389bSBarry Smith 
1520e33f6ddSBarry Smith /* not sure where to put this */
153dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*);
154dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
155dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
156dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
1572eac72dbSBarry Smith 
1580a71e23eSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC,KSP *);
1590a71e23eSMatthew Knepley 
160dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *);
161dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *);
1622eac72dbSBarry Smith 
163dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal);
164dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
165dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
166dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
167dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
1684b0e389bSBarry Smith 
169dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt);
170dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal);
1719f236934SBarry Smith 
172dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP);
173dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
174dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
175dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1761d73ed98SBarry Smith 
177dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt);
178dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP);
1791d73ed98SBarry Smith 
180b49fd9e1SBarry Smith /*E
181b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
182b49fd9e1SBarry Smith 
183b49fd9e1SBarry Smith    Level: advanced
184b49fd9e1SBarry Smith 
185b49fd9e1SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1868c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
187b49fd9e1SBarry Smith 
188b49fd9e1SBarry Smith E*/
18978d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
1909dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[];
1911f7e983dSSatish Balay /*MC
1928c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
1938c5b8ba0SBarry Smith 
1948c5b8ba0SBarry Smith    Level: advanced
1958c5b8ba0SBarry Smith 
1968c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
1978c5b8ba0SBarry Smith 
1988c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1998c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2008c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2018c5b8ba0SBarry Smith M*/
2028c5b8ba0SBarry Smith 
2031f7e983dSSatish Balay /*MC
2048c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2058c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2068c5b8ba0SBarry Smith           poor orthogonality.
2078c5b8ba0SBarry Smith 
2088c5b8ba0SBarry Smith    Level: advanced
2098c5b8ba0SBarry Smith 
2108c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2118c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2128c5b8ba0SBarry Smith 
2138c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2148c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2158c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2168c5b8ba0SBarry Smith M*/
2178c5b8ba0SBarry Smith 
2181f7e983dSSatish Balay /*MC
2198c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2208c5b8ba0SBarry Smith 
2218c5b8ba0SBarry Smith    Level: advanced
2228c5b8ba0SBarry Smith 
2238c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2248c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2258c5b8ba0SBarry Smith 
2268c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2278c5b8ba0SBarry Smith 
2288c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2298c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2308c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2318c5b8ba0SBarry Smith M*/
2328c5b8ba0SBarry Smith 
233dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
23408480c60SBarry Smith 
235dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
236dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
237dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
238c38d4ed2SBarry Smith 
239dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal);
240dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*);
241dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*);
242121fd945SKris Buschelman 
243d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal);
244d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth);
245d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,int);
246d9492815SBarry Smith 
247dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP);
248dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2492eac72dbSBarry Smith 
250a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
251a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
2521a2f9f99SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
2537517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
254a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
255a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
256a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
25784cb2905SBarry Smith 
258dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec);
259dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*);
260dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
261dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
262c01c455dSBarry Smith 
263dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure);
264dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
265906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperatorsSet(KSP,PetscTruth*,PetscTruth*);
266dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]);
267dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]);
2682dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]);
2691eb62cbbSBarry Smith 
270dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth);
271dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*);
272dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth);
273dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*);
2741f7f0c4fSBarry Smith 
275dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer);
2761eb62cbbSBarry Smith 
277db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRSetStandardErrorVec(KSP,Vec);
278db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRGetStandardErrorVec(KSP,Vec*);
279db9b2ab1SHong Zhang 
28028ce4d24SBarry Smith /*E
2818a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
2828a4b9c5cSBarry Smith        test routines.
2838a4b9c5cSBarry Smith 
2848a4b9c5cSBarry Smith    Level: advanced
2858a4b9c5cSBarry Smith 
286a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
287a3f661c8SBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
288a3f661c8SBarry Smith 
2898a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
2908a4b9c5cSBarry Smith 
29194b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
292a3f661c8SBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
2938a4b9c5cSBarry Smith E*/
294ce9499c7SBarry Smith typedef enum {KSP_NORM_NO = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
2959dcbbd2bSBarry Smith extern const char *KSPNormTypes[];
2961f7e983dSSatish Balay /*MC
297ce9499c7SBarry Smith     KSP_NORM_NO - Do not compute a norm during the Krylov process. This will
2988c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
2998c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3008c5b8ba0SBarry Smith 
3018c5b8ba0SBarry Smith    Level: advanced
3028c5b8ba0SBarry Smith 
303085a36d4SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
3048c5b8ba0SBarry Smith 
305ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3068c5b8ba0SBarry Smith M*/
3078c5b8ba0SBarry Smith 
3081f7e983dSSatish Balay /*MC
309ce9499c7SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the
3108c5b8ba0SBarry Smith        convergence test routine.
3118c5b8ba0SBarry Smith 
3128c5b8ba0SBarry Smith    Level: advanced
3138c5b8ba0SBarry Smith 
314ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3158c5b8ba0SBarry Smith M*/
3168c5b8ba0SBarry Smith 
3171f7e983dSSatish Balay /*MC
318ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3198c5b8ba0SBarry Smith        convergence test routine.
3208c5b8ba0SBarry Smith 
3218c5b8ba0SBarry Smith    Level: advanced
3228c5b8ba0SBarry Smith 
323ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3248c5b8ba0SBarry Smith M*/
3258c5b8ba0SBarry Smith 
3261f7e983dSSatish Balay /*MC
327ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
328a3f661c8SBarry Smith        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS
3298c5b8ba0SBarry Smith 
3308c5b8ba0SBarry Smith    Level: advanced
3318c5b8ba0SBarry Smith 
332ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3338c5b8ba0SBarry Smith M*/
3348c5b8ba0SBarry Smith 
335dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType);
336863a724eSLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNormType(KSP,KSPNormType*);
3373c5daeb9SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetCheckNormIteration(KSP,PetscInt);
338*1e7308c7SRichard Tran Mills EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetLagNorm(KSP,PetscTruth);
3398a4b9c5cSBarry Smith 
3408a4b9c5cSBarry Smith /*E
34128ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
34228ce4d24SBarry Smith          have converged or diverged
34328ce4d24SBarry Smith 
34428ce4d24SBarry Smith    Level: beginner
34528ce4d24SBarry Smith 
3464d0a8057SBarry Smith    Notes: See KSPGetConvergedReason() for explanation of each value
34728ce4d24SBarry Smith 
3484d0a8057SBarry Smith    Developer notes: this must match finclude/petscksp.h
3494d0a8057SBarry Smith 
3504d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
3514d0a8057SBarry Smith       any of the values here also change them that array of names.
35286c02ca4SBarry Smith 
353c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
35428ce4d24SBarry Smith E*/
355d15094e1SBarry Smith typedef enum {/* converged */
356d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
357d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
358b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3598031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
3608031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
361329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
362af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
363d15094e1SBarry Smith               /* diverged */
364b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
365d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
366d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
367d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
368b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
369b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
370b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
3716aee1118SBarry Smith               KSP_DIVERGED_NAN                 = -9,
3726aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
373d15094e1SBarry Smith 
374d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
3759dcbbd2bSBarry Smith extern const char **KSPConvergedReasons;
376d15094e1SBarry Smith 
377c838673bSBarry Smith /*MC
378c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
379c838673bSBarry Smith 
380c838673bSBarry Smith    Level: beginner
381c838673bSBarry Smith 
382c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
383c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
384c838673bSBarry Smith        2-norm of the residual for right preconditioning
385c838673bSBarry Smith 
386c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
387c838673bSBarry Smith 
388c838673bSBarry Smith M*/
389c838673bSBarry Smith 
390c838673bSBarry Smith /*MC
391c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
392c838673bSBarry Smith 
393c838673bSBarry Smith    Level: beginner
394c838673bSBarry Smith 
395c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
396c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
397c838673bSBarry Smith        2-norm of the residual for right preconditioning
398c838673bSBarry Smith 
399c838673bSBarry Smith    Level: beginner
400c838673bSBarry Smith 
401c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
402c838673bSBarry Smith 
403c838673bSBarry Smith M*/
404c838673bSBarry Smith 
405c838673bSBarry Smith /*MC
406c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
407c838673bSBarry Smith 
408c838673bSBarry Smith    Level: beginner
409c838673bSBarry Smith 
410c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
411c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
412c838673bSBarry Smith        2-norm of the residual for right preconditioning
413c838673bSBarry Smith 
414c838673bSBarry Smith    Level: beginner
415c838673bSBarry Smith 
416c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
417c838673bSBarry Smith 
418c838673bSBarry Smith M*/
419c838673bSBarry Smith 
420c838673bSBarry Smith /*MC
421c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
422c838673bSBarry Smith       reached
423c838673bSBarry Smith 
424c838673bSBarry Smith    Level: beginner
425c838673bSBarry Smith 
426c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
427c838673bSBarry Smith 
428c838673bSBarry Smith M*/
429c838673bSBarry Smith 
430c838673bSBarry Smith /*MC
4318c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
4328c9406f6SLisandro Dalcin            the preconditioner is applied. Also used when the KSPSkipConverged() convergence
4334d0a8057SBarry Smith            test routine is set in KSP.
434c838673bSBarry Smith 
435c838673bSBarry Smith 
436c838673bSBarry Smith    Level: beginner
437c838673bSBarry Smith 
438c838673bSBarry Smith 
439c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
440c838673bSBarry Smith 
441c838673bSBarry Smith M*/
442c838673bSBarry Smith 
443c838673bSBarry Smith /*MC
444c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
445c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
446c838673bSBarry Smith 
447c838673bSBarry Smith    Level: beginner
448c838673bSBarry Smith 
449c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
450c838673bSBarry Smith 
451c838673bSBarry Smith M*/
452c838673bSBarry Smith 
453c838673bSBarry Smith /*MC
454c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
455c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
456c838673bSBarry Smith 
457c838673bSBarry Smith 
458c838673bSBarry Smith    Level: beginner
459c838673bSBarry Smith 
460c838673bSBarry Smith 
461c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
462c838673bSBarry Smith 
463c838673bSBarry Smith M*/
464c838673bSBarry Smith 
465c838673bSBarry Smith /*MC
466c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
467c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
468c838673bSBarry Smith 
469c838673bSBarry Smith    Level: beginner
470c838673bSBarry Smith 
471c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
472c838673bSBarry Smith 
473c838673bSBarry Smith M*/
474c838673bSBarry Smith 
475c838673bSBarry Smith /*MC
476c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
477c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
478c838673bSBarry Smith         be positive definite
479c838673bSBarry Smith 
480c838673bSBarry Smith    Level: beginner
481c838673bSBarry Smith 
4822401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
483c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
484c838673bSBarry Smith 
485c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
486c838673bSBarry Smith 
487c838673bSBarry Smith M*/
488c838673bSBarry Smith 
489c838673bSBarry Smith /*MC
490c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
491c838673bSBarry Smith         while the KSPSolve() is still running.
492c838673bSBarry Smith 
493c838673bSBarry Smith    Level: beginner
494c838673bSBarry Smith 
495c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
496c838673bSBarry Smith 
497c838673bSBarry Smith M*/
498c838673bSBarry Smith 
499971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
500dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **);
501dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
502971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedDestroy(void *);
503971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedCreate(void **);
50494b02367SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUIRNorm(KSP);
5055f4984edSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUMIRNorm(KSP);
506dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
507dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *);
508abef13c0SSatish Balay 
509dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *);
510d4fbbf0eSBarry Smith 
51128ce4d24SBarry Smith /*E
51228ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
51328ce4d24SBarry Smith 
51428ce4d24SBarry Smith    Level: beginner
51528ce4d24SBarry Smith 
51628ce4d24SBarry Smith .seealso: KSPCGSetType()
51728ce4d24SBarry Smith E*/
5189dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5199dcbbd2bSBarry Smith extern const char *KSPCGTypes[];
52028ce4d24SBarry Smith 
521dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType);
5228031f4adStmunson 
523fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHSetRadius(KSP,PetscReal);
524fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetNormD(KSP,PetscReal *);
525fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetObjFcn(KSP,PetscReal *);
526fcae7a14Stmunson 
5270d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal);
52811f224b9Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetNormD(KSP,PetscReal *);
529d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetObjFcn(KSP,PetscReal *);
530e559a7feSLois Curfman McInnes 
5318031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRSetRadius(KSP,PetscReal);
5328031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetNormD(KSP,PetscReal *);
533d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetObjFcn(KSP,PetscReal *);
534d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetMinEig(KSP,PetscReal *);
5353c851360Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetLambda(KSP,PetscReal *);
5368031f4adStmunson 
537dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP);
538dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP);
5393369ce9aSBarry Smith 
540a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
541a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLG(KSP,PetscInt,PetscReal,void*);
542a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGDestroy(PetscDrawLG);
5437517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
5447517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
5457517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG);
546b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
547b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
548b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeDestroy(PetscDrawLG);
5492f2e5d10SKris Buschelman 
5500338f08bSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
55103919abeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
55203919abeSBarry Smith 
553e884886fSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatMFFDKSPMonitor(KSP,PetscInt,PetscReal,void *);
554e884886fSBarry Smith 
555f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
556436851aaSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscTruth monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
557f8a50e2bSBarry Smith 
558f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
559f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessDestroy(KSPFischerGuess);
560f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessReset(KSPFischerGuess);
561f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessUpdate(KSPFischerGuess,Vec);
562f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
563436851aaSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessSetFromOptions(KSPFischerGuess);
564f8a50e2bSBarry Smith 
565f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
566f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFischerGuess(KSP,KSPFischerGuess);
567f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetFischerGuess(KSP,KSPFischerGuess*);
568f8a50e2bSBarry Smith 
5693f22127dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatSchurComplementGetKSP(Mat,KSP*);
5700fe3d645SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,MatStructure);
5710fe3d645SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*);
5723f22127dSBarry Smith 
573e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
5742eac72dbSBarry Smith #endif
575