xref: /petsc/include/petscksp.h (revision b61692a6f9161375ec08c9563f854e524b42d5c4)
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"
561d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
579f6b3dcaSBarry Smith #define KSPBROYDEN    "broyden"
5858ad63f4SBarry Smith #define KSPGCR        "gcr"
592eac72dbSBarry Smith 
608ba1e511SMatthew Knepley /* Logging support */
610700a824SBarry Smith extern PetscClassId PETSCKSP_DLLEXPORT KSP_CLASSID;
628ba1e511SMatthew Knepley 
63dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *);
64a313700dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,const KSPType);
65dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP);
66dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP);
67dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec);
68dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec);
69dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP);
702eac72dbSBarry Smith 
71b0a32e0cSBarry Smith extern PetscFList KSPList;
72b022a5c1SBarry Smith extern PetscTruth KSPRegisterAllCalled;
73dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]);
74dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void);
75dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));
7630de9b25SBarry Smith 
7730de9b25SBarry Smith /*MC
7830de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
7930de9b25SBarry Smith 
8030de9b25SBarry Smith    Synopsis:
811890ba74SBarry Smith    PetscErrorCode KSPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(KSP))
8230de9b25SBarry Smith 
8330de9b25SBarry Smith    Not Collective
8430de9b25SBarry Smith 
8530de9b25SBarry Smith    Input Parameters:
8630de9b25SBarry Smith +  name_solver - name of a new user-defined solver
8730de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
8830de9b25SBarry Smith .  name_create - name of routine to create method context
8930de9b25SBarry Smith -  routine_create - routine to create method context
9030de9b25SBarry Smith 
9130de9b25SBarry Smith    Notes:
9230de9b25SBarry Smith    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
9330de9b25SBarry Smith 
9430de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
9530de9b25SBarry Smith    is ignored.
9630de9b25SBarry Smith 
9730de9b25SBarry Smith    Sample usage:
9830de9b25SBarry Smith .vb
9930de9b25SBarry Smith    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
10030de9b25SBarry Smith                "MySolverCreate",MySolverCreate);
10130de9b25SBarry Smith .ve
10230de9b25SBarry Smith 
10330de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
10430de9b25SBarry Smith $     KSPSetType(ksp,"my_solver")
10530de9b25SBarry Smith    or at runtime via the option
10630de9b25SBarry Smith $     -ksp_type my_solver
10730de9b25SBarry Smith 
10830de9b25SBarry Smith    Level: advanced
10930de9b25SBarry Smith 
110ab901514SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
11130de9b25SBarry Smith           and others of the form ${any_environmental_variable} occuring in pathname will be
11230de9b25SBarry Smith           replaced with appropriate values.
11330de9b25SBarry Smith          If your function is not being put into a shared library then use KSPRegister() instead
11430de9b25SBarry Smith 
11530de9b25SBarry Smith .keywords: KSP, register
11630de9b25SBarry Smith 
11730de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy()
11830de9b25SBarry Smith 
11930de9b25SBarry Smith M*/
120aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
121f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
1226df38c32SLois Curfman McInnes #else
123f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
1246df38c32SLois Curfman McInnes #endif
12582bf6240SBarry Smith 
126a313700dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,const KSPType *);
1271718446fSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPCSide(KSP,PCSide);
128b037da10SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPCSide(KSP,PCSide*);
129dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
130dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
131dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscTruth);
132dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscTruth *);
133dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscTruth);
134dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscTruth*);
135a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscTruth*);
136dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth);
137a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscTruth*);
138dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth);
139dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *);
140dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *);
141dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*);
142dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*);
143dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace);
144dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*);
145906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1462eac72dbSBarry Smith 
147dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC);
148dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*);
149aabeff55SBarry Smith 
150a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
151a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorCancel(KSP);
152dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **);
153dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
154dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);
1554b0e389bSBarry Smith 
1560e33f6ddSBarry Smith /* not sure where to put this */
157dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*);
158dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
159dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
160dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
1612eac72dbSBarry Smith 
1620a71e23eSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC,KSP *);
1630a71e23eSMatthew Knepley 
164dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *);
165dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *);
1662eac72dbSBarry Smith 
167dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal);
168*b61692a6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetSelfScale(KSP,PetscTruth);
169dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
170dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
171dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
172dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
1734b0e389bSBarry Smith 
174dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt);
175dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal);
1769f236934SBarry Smith 
177dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP);
178dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
179dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
180dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1811d73ed98SBarry Smith 
182dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt);
183dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP);
1841d73ed98SBarry Smith 
18558ad63f4SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRSetRestart(KSP,PetscInt);
1867ae9e5e1SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
18758ad63f4SBarry Smith 
188b49fd9e1SBarry Smith /*E
189b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
190b49fd9e1SBarry Smith 
191b49fd9e1SBarry Smith    Level: advanced
192b49fd9e1SBarry Smith 
193b49fd9e1SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1948c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
195b49fd9e1SBarry Smith 
196b49fd9e1SBarry Smith E*/
19778d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
1989dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[];
1991f7e983dSSatish Balay /*MC
2008c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
2018c5b8ba0SBarry Smith 
2028c5b8ba0SBarry Smith    Level: advanced
2038c5b8ba0SBarry Smith 
2048c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2058c5b8ba0SBarry Smith 
2068c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2078c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2088c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2098c5b8ba0SBarry Smith M*/
2108c5b8ba0SBarry Smith 
2111f7e983dSSatish Balay /*MC
2128c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2138c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2148c5b8ba0SBarry Smith           poor orthogonality.
2158c5b8ba0SBarry Smith 
2168c5b8ba0SBarry Smith    Level: advanced
2178c5b8ba0SBarry Smith 
2188c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2198c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2208c5b8ba0SBarry Smith 
2218c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2228c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2238c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2248c5b8ba0SBarry Smith M*/
2258c5b8ba0SBarry Smith 
2261f7e983dSSatish Balay /*MC
2278c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2288c5b8ba0SBarry Smith 
2298c5b8ba0SBarry Smith    Level: advanced
2308c5b8ba0SBarry Smith 
2318c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2328c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2338c5b8ba0SBarry Smith 
2348c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2358c5b8ba0SBarry Smith 
2368c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2378c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2388c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2398c5b8ba0SBarry Smith M*/
2408c5b8ba0SBarry Smith 
241dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
24208480c60SBarry Smith 
243dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
244dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
245dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
246c38d4ed2SBarry Smith 
247dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal);
248dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*);
249dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*);
250121fd945SKris Buschelman 
251d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal);
252d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth);
253d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,int);
254d9492815SBarry Smith 
255dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP);
256dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2572eac72dbSBarry Smith 
258a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
259a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
26090984118SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultLSQR(KSP,PetscInt,PetscReal,void *);
2611a2f9f99SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
2627517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
263a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
264a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
265a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
26684cb2905SBarry Smith 
267dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec);
268dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*);
269dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
270dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
271c01c455dSBarry Smith 
272dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure);
273dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
274906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperatorsSet(KSP,PetscTruth*,PetscTruth*);
275dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]);
276dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]);
2772dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]);
2781eb62cbbSBarry Smith 
279dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth);
280dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*);
281dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth);
282dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*);
2831f7f0c4fSBarry Smith 
284dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer);
2851eb62cbbSBarry Smith 
286db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRSetStandardErrorVec(KSP,Vec);
287db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRGetStandardErrorVec(KSP,Vec*);
288db9b2ab1SHong Zhang 
28928ce4d24SBarry Smith /*E
2908a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
2918a4b9c5cSBarry Smith        test routines.
2928a4b9c5cSBarry Smith 
2938a4b9c5cSBarry Smith    Level: advanced
2948a4b9c5cSBarry Smith 
295a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
2961718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
297a3f661c8SBarry Smith 
2988a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
2998a4b9c5cSBarry Smith 
30094b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
3011718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
3028a4b9c5cSBarry Smith E*/
303ce9499c7SBarry Smith typedef enum {KSP_NORM_NO = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
3049dcbbd2bSBarry Smith extern const char *KSPNormTypes[];
3051f7e983dSSatish Balay /*MC
306ce9499c7SBarry Smith     KSP_NORM_NO - Do not compute a norm during the Krylov process. This will
3078c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
3088c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3098c5b8ba0SBarry Smith 
3108c5b8ba0SBarry Smith    Level: advanced
3118c5b8ba0SBarry Smith 
312085a36d4SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
3138c5b8ba0SBarry Smith 
314ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3158c5b8ba0SBarry Smith M*/
3168c5b8ba0SBarry Smith 
3171f7e983dSSatish Balay /*MC
318ce9499c7SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual 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_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3248c5b8ba0SBarry Smith M*/
3258c5b8ba0SBarry Smith 
3261f7e983dSSatish Balay /*MC
327ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3288c5b8ba0SBarry Smith        convergence test routine.
3298c5b8ba0SBarry Smith 
3308c5b8ba0SBarry Smith    Level: advanced
3318c5b8ba0SBarry Smith 
332ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3338c5b8ba0SBarry Smith M*/
3348c5b8ba0SBarry Smith 
3351f7e983dSSatish Balay /*MC
336ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
337a3f661c8SBarry Smith        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS
3388c5b8ba0SBarry Smith 
3398c5b8ba0SBarry Smith    Level: advanced
3408c5b8ba0SBarry Smith 
341ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3428c5b8ba0SBarry Smith M*/
3438c5b8ba0SBarry Smith 
344dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType);
345863a724eSLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNormType(KSP,KSPNormType*);
3463c5daeb9SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetCheckNormIteration(KSP,PetscInt);
3471e7308c7SRichard Tran Mills EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetLagNorm(KSP,PetscTruth);
3488a4b9c5cSBarry Smith 
3498a4b9c5cSBarry Smith /*E
35028ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
35128ce4d24SBarry Smith          have converged or diverged
35228ce4d24SBarry Smith 
35328ce4d24SBarry Smith    Level: beginner
35428ce4d24SBarry Smith 
3554d0a8057SBarry Smith    Notes: See KSPGetConvergedReason() for explanation of each value
35628ce4d24SBarry Smith 
3574d0a8057SBarry Smith    Developer notes: this must match finclude/petscksp.h
3584d0a8057SBarry Smith 
3594d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
3604d0a8057SBarry Smith       any of the values here also change them that array of names.
36186c02ca4SBarry Smith 
362c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
36328ce4d24SBarry Smith E*/
364d15094e1SBarry Smith typedef enum {/* converged */
365d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
366d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
367b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3688031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
3698031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
370329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
371af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
372d15094e1SBarry Smith               /* diverged */
373b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
374d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
375d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
376d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
377b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
378b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
379b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
3806aee1118SBarry Smith               KSP_DIVERGED_NAN                 = -9,
3816aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
382d15094e1SBarry Smith 
383d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
3849dcbbd2bSBarry Smith extern const char **KSPConvergedReasons;
385d15094e1SBarry Smith 
386c838673bSBarry Smith /*MC
387c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
388c838673bSBarry Smith 
389c838673bSBarry Smith    Level: beginner
390c838673bSBarry Smith 
391c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
392c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
393c838673bSBarry Smith        2-norm of the residual for right preconditioning
394c838673bSBarry Smith 
395c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
396c838673bSBarry Smith 
397c838673bSBarry Smith M*/
398c838673bSBarry Smith 
399c838673bSBarry Smith /*MC
400c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
401c838673bSBarry Smith 
402c838673bSBarry Smith    Level: beginner
403c838673bSBarry Smith 
404c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
405c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
406c838673bSBarry Smith        2-norm of the residual for right preconditioning
407c838673bSBarry Smith 
408c838673bSBarry Smith    Level: beginner
409c838673bSBarry Smith 
410c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
411c838673bSBarry Smith 
412c838673bSBarry Smith M*/
413c838673bSBarry Smith 
414c838673bSBarry Smith /*MC
415c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
416c838673bSBarry Smith 
417c838673bSBarry Smith    Level: beginner
418c838673bSBarry Smith 
419c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
420c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
421c838673bSBarry Smith        2-norm of the residual for right preconditioning
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_ITS - Ran out of iterations before any convergence criteria was
431c838673bSBarry Smith       reached
432c838673bSBarry Smith 
433c838673bSBarry Smith    Level: beginner
434c838673bSBarry Smith 
435c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
436c838673bSBarry Smith 
437c838673bSBarry Smith M*/
438c838673bSBarry Smith 
439c838673bSBarry Smith /*MC
4408c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
4418c9406f6SLisandro Dalcin            the preconditioner is applied. Also used when the KSPSkipConverged() convergence
4424d0a8057SBarry Smith            test routine is set in KSP.
443c838673bSBarry Smith 
444c838673bSBarry Smith 
445c838673bSBarry Smith    Level: beginner
446c838673bSBarry Smith 
447c838673bSBarry Smith 
448c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
449c838673bSBarry Smith 
450c838673bSBarry Smith M*/
451c838673bSBarry Smith 
452c838673bSBarry Smith /*MC
453c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
454c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
455c838673bSBarry Smith 
456c838673bSBarry Smith    Level: beginner
457c838673bSBarry Smith 
458c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
459c838673bSBarry Smith 
460c838673bSBarry Smith M*/
461c838673bSBarry Smith 
462c838673bSBarry Smith /*MC
463c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
464c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
465c838673bSBarry Smith 
466c838673bSBarry Smith 
467c838673bSBarry Smith    Level: beginner
468c838673bSBarry Smith 
469c838673bSBarry Smith 
470c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
471c838673bSBarry Smith 
472c838673bSBarry Smith M*/
473c838673bSBarry Smith 
474c838673bSBarry Smith /*MC
475c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
476c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
477c838673bSBarry Smith 
478c838673bSBarry Smith    Level: beginner
479c838673bSBarry Smith 
480c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
481c838673bSBarry Smith 
482c838673bSBarry Smith M*/
483c838673bSBarry Smith 
484c838673bSBarry Smith /*MC
485c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
486c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
487c838673bSBarry Smith         be positive definite
488c838673bSBarry Smith 
489c838673bSBarry Smith    Level: beginner
490c838673bSBarry Smith 
4912401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
492c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
493c838673bSBarry Smith 
494c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
495c838673bSBarry Smith 
496c838673bSBarry Smith M*/
497c838673bSBarry Smith 
498c838673bSBarry Smith /*MC
499c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
500c838673bSBarry Smith         while the KSPSolve() is still running.
501c838673bSBarry Smith 
502c838673bSBarry Smith    Level: beginner
503c838673bSBarry Smith 
504c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
505c838673bSBarry Smith 
506c838673bSBarry Smith M*/
507c838673bSBarry Smith 
508971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
509dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **);
510dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
51190984118SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
512971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedDestroy(void *);
513971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedCreate(void **);
51494b02367SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUIRNorm(KSP);
5155f4984edSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUMIRNorm(KSP);
516dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
517dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *);
518abef13c0SSatish Balay 
519dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *);
520d4fbbf0eSBarry Smith 
52128ce4d24SBarry Smith /*E
52228ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
52328ce4d24SBarry Smith 
52428ce4d24SBarry Smith    Level: beginner
52528ce4d24SBarry Smith 
52628ce4d24SBarry Smith .seealso: KSPCGSetType()
52728ce4d24SBarry Smith E*/
5289dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5299dcbbd2bSBarry Smith extern const char *KSPCGTypes[];
53028ce4d24SBarry Smith 
531dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType);
532608df7e2SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGUseSingleReduction(KSP,PetscTruth);
5338031f4adStmunson 
534fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHSetRadius(KSP,PetscReal);
535fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetNormD(KSP,PetscReal *);
536fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetObjFcn(KSP,PetscReal *);
537fcae7a14Stmunson 
5380d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal);
53911f224b9Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetNormD(KSP,PetscReal *);
540d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetObjFcn(KSP,PetscReal *);
541e559a7feSLois Curfman McInnes 
5428031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRSetRadius(KSP,PetscReal);
5438031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetNormD(KSP,PetscReal *);
544d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetObjFcn(KSP,PetscReal *);
545d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetMinEig(KSP,PetscReal *);
5463c851360Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetLambda(KSP,PetscReal *);
5478031f4adStmunson 
5481d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPPythonSetType(KSP,const char[]);
5491d6018f0SLisandro Dalcin 
550dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP);
551dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP);
5523369ce9aSBarry Smith 
553a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
554a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLG(KSP,PetscInt,PetscReal,void*);
555a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGDestroy(PetscDrawLG);
5567517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
5577517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
5587517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG);
559b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
560b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
561b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeDestroy(PetscDrawLG);
5622f2e5d10SKris Buschelman 
5636891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
5646891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
56503919abeSBarry Smith 
566f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
567436851aaSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscTruth monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
568f8a50e2bSBarry Smith 
569f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
570f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessDestroy(KSPFischerGuess);
571f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessReset(KSPFischerGuess);
572f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessUpdate(KSPFischerGuess,Vec);
573f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
574436851aaSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessSetFromOptions(KSPFischerGuess);
575f8a50e2bSBarry Smith 
576f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
577f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFischerGuess(KSP,KSPFischerGuess);
578f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetFischerGuess(KSP,KSPFischerGuess*);
579f8a50e2bSBarry Smith 
580084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
5813f22127dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatSchurComplementGetKSP(Mat,KSP*);
582084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure);
583084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
584dfe085dbSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *);
5853f22127dSBarry Smith 
5866c699258SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDM(KSP,DM);
587f22f69f0SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDMActive(KSP,PetscTruth);
5886c699258SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDM(KSP,DM*);
5896c699258SBarry Smith 
590e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
5912eac72dbSBarry Smith #endif
592