xref: /petsc/include/petscksp.h (revision 436851aaa7e03ea0d9fc31aefde2b6b3f89a5be9)
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*/
31e5a9bf91SBarry Smith #define KSPType const 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"
44f0037002Svictorle #define KSPBCGSL      "bcgsl"
4582bf6240SBarry Smith #define KSPCGS        "cgs"
4682bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
4782bf6240SBarry Smith #define KSPCR         "cr"
4882bf6240SBarry Smith #define KSPLSQR       "lsqr"
4982bf6240SBarry Smith #define KSPPREONLY    "preonly"
5082bf6240SBarry Smith #define KSPQCG        "qcg"
51c9cf9b11SBarry Smith #define KSPBICG       "bicg"
52b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
5301c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
5421356919SSatish Balay #define KSPLCD        "lcd"
552eac72dbSBarry Smith 
568ba1e511SMatthew Knepley /* Logging support */
57dba47a55SKris Buschelman extern PetscCookie PETSCKSP_DLLEXPORT KSP_COOKIE;
588ba1e511SMatthew Knepley 
59dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *);
60f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,KSPType);
61dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP);
62dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP);
63dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec);
64dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec);
65dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP);
662eac72dbSBarry Smith 
67b0a32e0cSBarry Smith extern PetscFList KSPList;
68dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]);
69dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void);
70dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));
7130de9b25SBarry Smith 
7230de9b25SBarry Smith /*MC
7330de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
7430de9b25SBarry Smith 
7530de9b25SBarry Smith    Synopsis:
76d360dc6fSBarry Smith    PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP))
7730de9b25SBarry Smith 
7830de9b25SBarry Smith    Not Collective
7930de9b25SBarry Smith 
8030de9b25SBarry Smith    Input Parameters:
8130de9b25SBarry Smith +  name_solver - name of a new user-defined solver
8230de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
8330de9b25SBarry Smith .  name_create - name of routine to create method context
8430de9b25SBarry Smith -  routine_create - routine to create method context
8530de9b25SBarry Smith 
8630de9b25SBarry Smith    Notes:
8730de9b25SBarry Smith    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
8830de9b25SBarry Smith 
8930de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
9030de9b25SBarry Smith    is ignored.
9130de9b25SBarry Smith 
9230de9b25SBarry Smith    Sample usage:
9330de9b25SBarry Smith .vb
9430de9b25SBarry Smith    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
9530de9b25SBarry Smith                "MySolverCreate",MySolverCreate);
9630de9b25SBarry Smith .ve
9730de9b25SBarry Smith 
9830de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
9930de9b25SBarry Smith $     KSPSetType(ksp,"my_solver")
10030de9b25SBarry Smith    or at runtime via the option
10130de9b25SBarry Smith $     -ksp_type my_solver
10230de9b25SBarry Smith 
10330de9b25SBarry Smith    Level: advanced
10430de9b25SBarry Smith 
105ab901514SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
10630de9b25SBarry Smith           and others of the form ${any_environmental_variable} occuring in pathname will be
10730de9b25SBarry Smith           replaced with appropriate values.
10830de9b25SBarry Smith          If your function is not being put into a shared library then use KSPRegister() instead
10930de9b25SBarry Smith 
11030de9b25SBarry Smith .keywords: KSP, register
11130de9b25SBarry Smith 
11230de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy()
11330de9b25SBarry Smith 
11430de9b25SBarry Smith M*/
115aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
116f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
1176df38c32SLois Curfman McInnes #else
118f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
1196df38c32SLois Curfman McInnes #endif
12082bf6240SBarry Smith 
121dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,KSPType *);
122dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPreconditionerSide(KSP,PCSide);
123dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPreconditionerSide(KSP,PCSide*);
124dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
125dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
126dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscTruth);
127dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscTruth *);
128dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscTruth);
129dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscTruth*);
130a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscTruth*);
131dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth);
132a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscTruth*);
133dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth);
134dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *);
135dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *);
136dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*);
137dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*);
138dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace);
139dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*);
140906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1412eac72dbSBarry Smith 
142dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC);
143dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*);
144aabeff55SBarry Smith 
145a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
146a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorCancel(KSP);
147dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **);
148dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
149dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);
1504b0e389bSBarry Smith 
1510e33f6ddSBarry Smith /* not sure where to put this */
152dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*);
153dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
154dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
155dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
1562eac72dbSBarry Smith 
1570a71e23eSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC,KSP *);
1580a71e23eSMatthew Knepley 
159dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *);
160dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *);
1612eac72dbSBarry Smith 
162dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal);
163dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
164dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
165dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
166dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
1674b0e389bSBarry Smith 
168dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt);
169dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal);
1709f236934SBarry Smith 
171dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP);
172dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
173dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
174dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1751d73ed98SBarry Smith 
176dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt);
177dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP);
1781d73ed98SBarry Smith 
179b49fd9e1SBarry Smith /*E
180b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
181b49fd9e1SBarry Smith 
182b49fd9e1SBarry Smith    Level: advanced
183b49fd9e1SBarry Smith 
184b49fd9e1SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1858c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
186b49fd9e1SBarry Smith 
187b49fd9e1SBarry Smith E*/
18878d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
1899dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[];
1901f7e983dSSatish Balay /*MC
1918c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
1928c5b8ba0SBarry Smith 
1938c5b8ba0SBarry Smith    Level: advanced
1948c5b8ba0SBarry Smith 
1958c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
1968c5b8ba0SBarry Smith 
1978c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1988c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
1998c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2008c5b8ba0SBarry Smith M*/
2018c5b8ba0SBarry Smith 
2021f7e983dSSatish Balay /*MC
2038c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2048c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2058c5b8ba0SBarry Smith           poor orthogonality.
2068c5b8ba0SBarry Smith 
2078c5b8ba0SBarry Smith    Level: advanced
2088c5b8ba0SBarry Smith 
2098c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2108c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2118c5b8ba0SBarry Smith 
2128c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2138c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2148c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2158c5b8ba0SBarry Smith M*/
2168c5b8ba0SBarry Smith 
2171f7e983dSSatish Balay /*MC
2188c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2198c5b8ba0SBarry Smith 
2208c5b8ba0SBarry Smith    Level: advanced
2218c5b8ba0SBarry Smith 
2228c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2238c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2248c5b8ba0SBarry Smith 
2258c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2268c5b8ba0SBarry Smith 
2278c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2288c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2298c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2308c5b8ba0SBarry Smith M*/
2318c5b8ba0SBarry Smith 
232dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
23308480c60SBarry Smith 
234dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
235dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
236dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
237c38d4ed2SBarry Smith 
238dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal);
239dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*);
240dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*);
241121fd945SKris Buschelman 
242d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal);
243d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth);
244d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,int);
245d9492815SBarry Smith 
246dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP);
247dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2482eac72dbSBarry Smith 
249a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
250a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
2517517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
252a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
253a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
254a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
25584cb2905SBarry Smith 
256dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec);
257dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*);
258dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
259dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
260c01c455dSBarry Smith 
261dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure);
262dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
263906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperatorsSet(KSP,PetscTruth*,PetscTruth*);
264dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]);
265dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]);
2662dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]);
2671eb62cbbSBarry Smith 
268dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth);
269dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*);
270dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth);
271dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*);
2721f7f0c4fSBarry Smith 
273dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer);
2741eb62cbbSBarry Smith 
275db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRSetStandardErrorVec(KSP,Vec);
276db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRGetStandardErrorVec(KSP,Vec*);
277db9b2ab1SHong Zhang 
27828ce4d24SBarry Smith /*E
2798a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
2808a4b9c5cSBarry Smith        test routines.
2818a4b9c5cSBarry Smith 
2828a4b9c5cSBarry Smith    Level: advanced
2838a4b9c5cSBarry Smith 
2848a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
2858a4b9c5cSBarry Smith 
28694b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
2878a4b9c5cSBarry Smith           KSPSetConvergenceTest()
2888a4b9c5cSBarry Smith E*/
289ce9499c7SBarry Smith typedef enum {KSP_NORM_NO = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
2909dcbbd2bSBarry Smith extern const char *KSPNormTypes[];
2911f7e983dSSatish Balay /*MC
292ce9499c7SBarry Smith     KSP_NORM_NO - Do not compute a norm during the Krylov process. This will
2938c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
2948c5b8ba0SBarry Smith           be based on a norm of a residual etc.
2958c5b8ba0SBarry Smith 
2968c5b8ba0SBarry Smith    Level: advanced
2978c5b8ba0SBarry Smith 
298085a36d4SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
2998c5b8ba0SBarry Smith 
300ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3018c5b8ba0SBarry Smith M*/
3028c5b8ba0SBarry Smith 
3031f7e983dSSatish Balay /*MC
304ce9499c7SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the
3058c5b8ba0SBarry Smith        convergence test routine.
3068c5b8ba0SBarry Smith 
3078c5b8ba0SBarry Smith    Level: advanced
3088c5b8ba0SBarry Smith 
309ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3108c5b8ba0SBarry Smith M*/
3118c5b8ba0SBarry Smith 
3121f7e983dSSatish Balay /*MC
313ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3148c5b8ba0SBarry Smith        convergence test routine.
3158c5b8ba0SBarry Smith 
3168c5b8ba0SBarry Smith    Level: advanced
3178c5b8ba0SBarry Smith 
318ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3198c5b8ba0SBarry Smith M*/
3208c5b8ba0SBarry Smith 
3211f7e983dSSatish Balay /*MC
322ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
3238c5b8ba0SBarry Smith        convergence test routine.
3248c5b8ba0SBarry Smith 
3258c5b8ba0SBarry Smith    Level: advanced
3268c5b8ba0SBarry Smith 
327ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3288c5b8ba0SBarry Smith M*/
3298c5b8ba0SBarry Smith 
330dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType);
331863a724eSLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNormType(KSP,KSPNormType*);
3323c5daeb9SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetCheckNormIteration(KSP,PetscInt);
3338a4b9c5cSBarry Smith 
3348a4b9c5cSBarry Smith /*E
33528ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
33628ce4d24SBarry Smith          have converged or diverged
33728ce4d24SBarry Smith 
33828ce4d24SBarry Smith    Level: beginner
33928ce4d24SBarry Smith 
34028ce4d24SBarry Smith    Notes: this must match finclude/petscksp.h
34128ce4d24SBarry Smith 
34286c02ca4SBarry Smith    Developer note: The string versions of these are in
3434b9ad928SBarry Smith      src/ksp/ksp/interface/itfunc.c called convergedreasons.
34450f1acb2SBarry Smith      If these enums are changed you must change those.
34586c02ca4SBarry Smith 
346c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
34728ce4d24SBarry Smith E*/
348d15094e1SBarry Smith typedef enum {/* converged */
349d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
350d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
351b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3528031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
3538031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
354329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
355af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
356d15094e1SBarry Smith               /* diverged */
357b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
358d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
359d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
360d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
361b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
362b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
363b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
3646aee1118SBarry Smith               KSP_DIVERGED_NAN                 = -9,
3656aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
366d15094e1SBarry Smith 
367d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
3689dcbbd2bSBarry Smith extern const char **KSPConvergedReasons;
369d15094e1SBarry Smith 
370c838673bSBarry Smith /*MC
371c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
372c838673bSBarry Smith 
373c838673bSBarry Smith    Level: beginner
374c838673bSBarry Smith 
375c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
376c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
377c838673bSBarry Smith        2-norm of the residual for right preconditioning
378c838673bSBarry Smith 
379c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
380c838673bSBarry Smith 
381c838673bSBarry Smith M*/
382c838673bSBarry Smith 
383c838673bSBarry Smith /*MC
384c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
385c838673bSBarry Smith 
386c838673bSBarry Smith    Level: beginner
387c838673bSBarry Smith 
388c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
389c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
390c838673bSBarry Smith        2-norm of the residual for right preconditioning
391c838673bSBarry Smith 
392c838673bSBarry Smith    Level: beginner
393c838673bSBarry Smith 
394c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
395c838673bSBarry Smith 
396c838673bSBarry Smith M*/
397c838673bSBarry Smith 
398c838673bSBarry Smith /*MC
399c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
400c838673bSBarry Smith 
401c838673bSBarry Smith    Level: beginner
402c838673bSBarry Smith 
403c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
404c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
405c838673bSBarry Smith        2-norm of the residual for right preconditioning
406c838673bSBarry Smith 
407c838673bSBarry Smith    Level: beginner
408c838673bSBarry Smith 
409c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
410c838673bSBarry Smith 
411c838673bSBarry Smith M*/
412c838673bSBarry Smith 
413c838673bSBarry Smith /*MC
414c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
415c838673bSBarry Smith       reached
416c838673bSBarry Smith 
417c838673bSBarry Smith    Level: beginner
418c838673bSBarry Smith 
419c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
420c838673bSBarry Smith 
421c838673bSBarry Smith M*/
422c838673bSBarry Smith 
423c838673bSBarry Smith /*MC
4248c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
4258c9406f6SLisandro Dalcin            the preconditioner is applied. Also used when the KSPSkipConverged() convergence
4268c9406f6SLisandro Dalcin            test rutine is set in KSP.
427c838673bSBarry Smith 
428c838673bSBarry Smith 
429c838673bSBarry Smith    Level: beginner
430c838673bSBarry Smith 
431c838673bSBarry Smith 
432c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
433c838673bSBarry Smith 
434c838673bSBarry Smith M*/
435c838673bSBarry Smith 
436c838673bSBarry Smith /*MC
437c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
438c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
439c838673bSBarry Smith 
440c838673bSBarry Smith    Level: beginner
441c838673bSBarry Smith 
442c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
443c838673bSBarry Smith 
444c838673bSBarry Smith M*/
445c838673bSBarry Smith 
446c838673bSBarry Smith /*MC
447c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
448c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
449c838673bSBarry Smith 
450c838673bSBarry Smith 
451c838673bSBarry Smith    Level: beginner
452c838673bSBarry Smith 
453c838673bSBarry Smith 
454c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
455c838673bSBarry Smith 
456c838673bSBarry Smith M*/
457c838673bSBarry Smith 
458c838673bSBarry Smith /*MC
459c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
460c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
461c838673bSBarry Smith 
462c838673bSBarry Smith    Level: beginner
463c838673bSBarry Smith 
464c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
465c838673bSBarry Smith 
466c838673bSBarry Smith M*/
467c838673bSBarry Smith 
468c838673bSBarry Smith /*MC
469c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
470c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
471c838673bSBarry Smith         be positive definite
472c838673bSBarry Smith 
473c838673bSBarry Smith    Level: beginner
474c838673bSBarry Smith 
4752401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
476c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
477c838673bSBarry Smith 
478c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
479c838673bSBarry Smith 
480c838673bSBarry Smith M*/
481c838673bSBarry Smith 
482c838673bSBarry Smith /*MC
483c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
484c838673bSBarry Smith         while the KSPSolve() is still running.
485c838673bSBarry Smith 
486c838673bSBarry Smith    Level: beginner
487c838673bSBarry Smith 
488c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
489c838673bSBarry Smith 
490c838673bSBarry Smith M*/
491c838673bSBarry Smith 
492971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
493dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **);
494dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
495971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedDestroy(void *);
496971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedCreate(void **);
49794b02367SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUIRNorm(KSP);
4985f4984edSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUMIRNorm(KSP);
499dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
500dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *);
501abef13c0SSatish Balay 
502dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *);
503d4fbbf0eSBarry Smith 
50428ce4d24SBarry Smith /*E
50528ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
50628ce4d24SBarry Smith 
50728ce4d24SBarry Smith    Level: beginner
50828ce4d24SBarry Smith 
50928ce4d24SBarry Smith .seealso: KSPCGSetType()
51028ce4d24SBarry Smith E*/
5119dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5129dcbbd2bSBarry Smith extern const char *KSPCGTypes[];
51328ce4d24SBarry Smith 
514dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType);
5158031f4adStmunson 
516fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHSetRadius(KSP,PetscReal);
517fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetNormD(KSP,PetscReal *);
518fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetObjFcn(KSP,PetscReal *);
519fcae7a14Stmunson 
5200d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal);
52111f224b9Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetNormD(KSP,PetscReal *);
522d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetObjFcn(KSP,PetscReal *);
523e559a7feSLois Curfman McInnes 
5248031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRSetRadius(KSP,PetscReal);
5258031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetNormD(KSP,PetscReal *);
526d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetObjFcn(KSP,PetscReal *);
527d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetMinEig(KSP,PetscReal *);
5283c851360Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetLambda(KSP,PetscReal *);
5298031f4adStmunson 
530dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP);
531dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP);
5323369ce9aSBarry Smith 
533a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
534a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLG(KSP,PetscInt,PetscReal,void*);
535a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGDestroy(PetscDrawLG);
5367517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
5377517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
5387517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG);
5392f2e5d10SKris Buschelman 
5400338f08bSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
54103919abeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
54203919abeSBarry Smith 
543e884886fSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatMFFDKSPMonitor(KSP,PetscInt,PetscReal,void *);
544e884886fSBarry Smith 
545f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
546*436851aaSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscTruth monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
547f8a50e2bSBarry Smith 
548f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
549f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessDestroy(KSPFischerGuess);
550f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessReset(KSPFischerGuess);
551f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessUpdate(KSPFischerGuess,Vec);
552f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
553*436851aaSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessSetFromOptions(KSPFischerGuess);
554f8a50e2bSBarry Smith 
555f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
556f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFischerGuess(KSP,KSPFischerGuess);
557f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetFischerGuess(KSP,KSPFischerGuess*);
558f8a50e2bSBarry Smith 
559e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
5602eac72dbSBarry Smith #endif
561