xref: /petsc/include/petscksp.h (revision 7087cfbefd1a42b179f217f9994fb6cb0d0c1824)
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 
9*7087cfbeSBarry Smith extern PetscErrorCode  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"
597945016fSBarry Smith #define KSPNGMRES     "ngmres"
60bbce1389SJed Brown #define KSPSPECEST    "specest"
612eac72dbSBarry Smith 
628ba1e511SMatthew Knepley /* Logging support */
63*7087cfbeSBarry Smith extern PetscClassId  KSP_CLASSID;
648ba1e511SMatthew Knepley 
65*7087cfbeSBarry Smith extern PetscErrorCode  KSPCreate(MPI_Comm,KSP *);
66*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetType(KSP,const KSPType);
67*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetUp(KSP);
68*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetUpOnBlocks(KSP);
69*7087cfbeSBarry Smith extern PetscErrorCode  KSPSolve(KSP,Vec,Vec);
70*7087cfbeSBarry Smith extern PetscErrorCode  KSPSolveTranspose(KSP,Vec,Vec);
71*7087cfbeSBarry Smith extern PetscErrorCode  KSPDestroy(KSP);
722eac72dbSBarry Smith 
73b0a32e0cSBarry Smith extern PetscFList KSPList;
74ace3abfcSBarry Smith extern PetscBool  KSPRegisterAllCalled;
75*7087cfbeSBarry Smith extern PetscErrorCode  KSPRegisterAll(const char[]);
76*7087cfbeSBarry Smith extern PetscErrorCode  KSPRegisterDestroy(void);
77*7087cfbeSBarry Smith extern PetscErrorCode  KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));
7830de9b25SBarry Smith 
7930de9b25SBarry Smith /*MC
8030de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
8130de9b25SBarry Smith 
8230de9b25SBarry Smith    Synopsis:
831890ba74SBarry Smith    PetscErrorCode KSPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(KSP))
8430de9b25SBarry Smith 
8530de9b25SBarry Smith    Not Collective
8630de9b25SBarry Smith 
8730de9b25SBarry Smith    Input Parameters:
8830de9b25SBarry Smith +  name_solver - name of a new user-defined solver
8930de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
9030de9b25SBarry Smith .  name_create - name of routine to create method context
9130de9b25SBarry Smith -  routine_create - routine to create method context
9230de9b25SBarry Smith 
9330de9b25SBarry Smith    Notes:
9430de9b25SBarry Smith    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
9530de9b25SBarry Smith 
9630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
9730de9b25SBarry Smith    is ignored.
9830de9b25SBarry Smith 
9930de9b25SBarry Smith    Sample usage:
10030de9b25SBarry Smith .vb
10130de9b25SBarry Smith    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
10230de9b25SBarry Smith                "MySolverCreate",MySolverCreate);
10330de9b25SBarry Smith .ve
10430de9b25SBarry Smith 
10530de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
10630de9b25SBarry Smith $     KSPSetType(ksp,"my_solver")
10730de9b25SBarry Smith    or at runtime via the option
10830de9b25SBarry Smith $     -ksp_type my_solver
10930de9b25SBarry Smith 
11030de9b25SBarry Smith    Level: advanced
11130de9b25SBarry Smith 
112ab901514SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
11330de9b25SBarry Smith           and others of the form ${any_environmental_variable} occuring in pathname will be
11430de9b25SBarry Smith           replaced with appropriate values.
11530de9b25SBarry Smith          If your function is not being put into a shared library then use KSPRegister() instead
11630de9b25SBarry Smith 
11730de9b25SBarry Smith .keywords: KSP, register
11830de9b25SBarry Smith 
11930de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy()
12030de9b25SBarry Smith 
12130de9b25SBarry Smith M*/
122aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
123f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
1246df38c32SLois Curfman McInnes #else
125f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
1266df38c32SLois Curfman McInnes #endif
12782bf6240SBarry Smith 
128*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetType(KSP,const KSPType *);
129*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetPCSide(KSP,PCSide);
130*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetPCSide(KSP,PCSide*);
131*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
132*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
133*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetInitialGuessNonzero(KSP,PetscBool );
134*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetInitialGuessNonzero(KSP,PetscBool  *);
135*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetInitialGuessKnoll(KSP,PetscBool );
136*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetInitialGuessKnoll(KSP,PetscBool *);
137*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetErrorIfNotConverged(KSP,PetscBool );
138*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetErrorIfNotConverged(KSP,PetscBool  *);
139*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetComputeEigenvalues(KSP,PetscBool *);
140*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetComputeEigenvalues(KSP,PetscBool );
141*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetComputeSingularValues(KSP,PetscBool *);
142*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetComputeSingularValues(KSP,PetscBool );
143*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetRhs(KSP,Vec *);
144*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetSolution(KSP,Vec *);
145*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetResidualNorm(KSP,PetscReal*);
146*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetIterationNumber(KSP,PetscInt*);
147*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetNullSpace(KSP,MatNullSpace);
148*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetNullSpace(KSP,MatNullSpace*);
149*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1502eac72dbSBarry Smith 
151*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetPC(KSP,PC);
152*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetPC(KSP,PC*);
153aabeff55SBarry Smith 
154*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
155*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorCancel(KSP);
156*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetMonitorContext(KSP,void **);
157*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
158*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
1594b0e389bSBarry Smith 
1600e33f6ddSBarry Smith /* not sure where to put this */
161*7087cfbeSBarry Smith extern PetscErrorCode  PCKSPGetKSP(PC,KSP*);
162*7087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
163*7087cfbeSBarry Smith extern PetscErrorCode  PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
164*7087cfbeSBarry Smith extern PetscErrorCode  PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
165*7087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
1662eac72dbSBarry Smith 
167*7087cfbeSBarry Smith extern PetscErrorCode  PCGalerkinGetKSP(PC,KSP *);
1680a71e23eSMatthew Knepley 
169*7087cfbeSBarry Smith extern PetscErrorCode  KSPBuildSolution(KSP,Vec,Vec *);
170*7087cfbeSBarry Smith extern PetscErrorCode  KSPBuildResidual(KSP,Vec,Vec,Vec *);
1712eac72dbSBarry Smith 
172*7087cfbeSBarry Smith extern PetscErrorCode  KSPRichardsonSetScale(KSP,PetscReal);
173*7087cfbeSBarry Smith extern PetscErrorCode  KSPRichardsonSetSelfScale(KSP,PetscBool );
174*7087cfbeSBarry Smith extern PetscErrorCode  KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
175*7087cfbeSBarry Smith extern PetscErrorCode  KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
176*7087cfbeSBarry Smith extern PetscErrorCode  KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
177*7087cfbeSBarry Smith extern PetscErrorCode  KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
1784b0e389bSBarry Smith 
179*7087cfbeSBarry Smith extern PetscErrorCode  KSPGMRESSetRestart(KSP, PetscInt);
180*7087cfbeSBarry Smith extern PetscErrorCode  KSPGMRESGetRestart(KSP, PetscInt*);
181*7087cfbeSBarry Smith extern PetscErrorCode  KSPGMRESSetHapTol(KSP,PetscReal);
1829f236934SBarry Smith 
183*7087cfbeSBarry Smith extern PetscErrorCode  KSPGMRESSetPreAllocateVectors(KSP);
184*7087cfbeSBarry Smith extern PetscErrorCode  KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
185*7087cfbeSBarry Smith extern PetscErrorCode  KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
186*7087cfbeSBarry Smith extern PetscErrorCode  KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
187*7087cfbeSBarry Smith extern PetscErrorCode  KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1881d73ed98SBarry Smith 
189*7087cfbeSBarry Smith extern PetscErrorCode  KSPLGMRESSetAugDim(KSP,PetscInt);
190*7087cfbeSBarry Smith extern PetscErrorCode  KSPLGMRESSetConstant(KSP);
1911d73ed98SBarry Smith 
192*7087cfbeSBarry Smith extern PetscErrorCode  KSPGCRSetRestart(KSP,PetscInt);
193*7087cfbeSBarry Smith extern PetscErrorCode  KSPGCRGetRestart(KSP,PetscInt*);
194*7087cfbeSBarry Smith extern PetscErrorCode  KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
19558ad63f4SBarry Smith 
196b49fd9e1SBarry Smith /*E
197b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
198b49fd9e1SBarry Smith 
199b49fd9e1SBarry Smith    Level: advanced
200b49fd9e1SBarry Smith 
201bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
202e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
203b49fd9e1SBarry Smith 
204b49fd9e1SBarry Smith E*/
20578d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2069dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[];
2071f7e983dSSatish Balay /*MC
2088c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
2098c5b8ba0SBarry Smith 
2108c5b8ba0SBarry Smith    Level: advanced
2118c5b8ba0SBarry Smith 
2128c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2138c5b8ba0SBarry Smith 
214bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
215e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2168c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2178c5b8ba0SBarry Smith M*/
2188c5b8ba0SBarry Smith 
2191f7e983dSSatish Balay /*MC
2208c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2218c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2228c5b8ba0SBarry Smith           poor orthogonality.
2238c5b8ba0SBarry Smith 
2248c5b8ba0SBarry Smith    Level: advanced
2258c5b8ba0SBarry Smith 
2268c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2278c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2288c5b8ba0SBarry Smith 
229bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
230e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2318c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2328c5b8ba0SBarry Smith M*/
2338c5b8ba0SBarry Smith 
2341f7e983dSSatish Balay /*MC
2358c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2368c5b8ba0SBarry Smith 
2378c5b8ba0SBarry Smith    Level: advanced
2388c5b8ba0SBarry Smith 
2398c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2408c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2418c5b8ba0SBarry Smith 
2428c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2438c5b8ba0SBarry Smith 
244bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
245e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2468c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2478c5b8ba0SBarry Smith M*/
2488c5b8ba0SBarry Smith 
249*7087cfbeSBarry Smith extern PetscErrorCode  KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
250*7087cfbeSBarry Smith extern PetscErrorCode  KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
25108480c60SBarry Smith 
252*7087cfbeSBarry Smith extern PetscErrorCode  KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
253*7087cfbeSBarry Smith extern PetscErrorCode  KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
254*7087cfbeSBarry Smith extern PetscErrorCode  KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
255c38d4ed2SBarry Smith 
256*7087cfbeSBarry Smith extern PetscErrorCode  KSPQCGSetTrustRegionRadius(KSP,PetscReal);
257*7087cfbeSBarry Smith extern PetscErrorCode  KSPQCGGetQuadratic(KSP,PetscReal*);
258*7087cfbeSBarry Smith extern PetscErrorCode  KSPQCGGetTrialStepNorm(KSP,PetscReal*);
259121fd945SKris Buschelman 
260*7087cfbeSBarry Smith extern PetscErrorCode  KSPBCGSLSetXRes(KSP,PetscReal);
261*7087cfbeSBarry Smith extern PetscErrorCode  KSPBCGSLSetPol(KSP,PetscBool );
262*7087cfbeSBarry Smith extern PetscErrorCode  KSPBCGSLSetEll(KSP,PetscInt);
263d9492815SBarry Smith 
264*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetFromOptions(KSP);
265*7087cfbeSBarry Smith extern PetscErrorCode  KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2662eac72dbSBarry Smith 
267*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
268*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
269*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorDefaultLSQR(KSP,PetscInt,PetscReal,void *);
270*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
271*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
272*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
273*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
274*7087cfbeSBarry Smith extern PetscErrorCode  KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
27584cb2905SBarry Smith 
276*7087cfbeSBarry Smith extern PetscErrorCode  KSPUnwindPreconditioner(KSP,Vec,Vec);
277*7087cfbeSBarry Smith extern PetscErrorCode  KSPDefaultBuildSolution(KSP,Vec,Vec*);
278*7087cfbeSBarry Smith extern PetscErrorCode  KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
279*7087cfbeSBarry Smith extern PetscErrorCode  KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
280c01c455dSBarry Smith 
281*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetOperators(KSP,Mat,Mat,MatStructure);
282*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
283*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
284*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetOptionsPrefix(KSP,const char[]);
285*7087cfbeSBarry Smith extern PetscErrorCode  KSPAppendOptionsPrefix(KSP,const char[]);
286*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetOptionsPrefix(KSP,const char*[]);
2871eb62cbbSBarry Smith 
288*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetDiagonalScale(KSP,PetscBool );
289*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetDiagonalScale(KSP,PetscBool *);
290*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetDiagonalScaleFix(KSP,PetscBool );
291*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetDiagonalScaleFix(KSP,PetscBool *);
2921f7f0c4fSBarry Smith 
293*7087cfbeSBarry Smith extern PetscErrorCode  KSPView(KSP,PetscViewer);
2941eb62cbbSBarry Smith 
295*7087cfbeSBarry Smith extern PetscErrorCode  KSPLSQRSetStandardErrorVec(KSP,Vec);
296*7087cfbeSBarry Smith extern PetscErrorCode  KSPLSQRGetStandardErrorVec(KSP,Vec*);
297db9b2ab1SHong Zhang 
29828ce4d24SBarry Smith /*E
2998a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
3008a4b9c5cSBarry Smith        test routines.
3018a4b9c5cSBarry Smith 
3028a4b9c5cSBarry Smith    Level: advanced
3038a4b9c5cSBarry Smith 
304a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
3051718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
306a3f661c8SBarry Smith 
3078a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
3088a4b9c5cSBarry Smith 
30994b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
3101718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
3118a4b9c5cSBarry Smith E*/
3129793b452SBarry Smith typedef enum {KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
3139dcbbd2bSBarry Smith extern const char *KSPNormTypes[];
3141f7e983dSSatish Balay /*MC
3159793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
3168c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
3178c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3188c5b8ba0SBarry Smith 
3198c5b8ba0SBarry Smith    Level: advanced
3208c5b8ba0SBarry Smith 
321085a36d4SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
3228c5b8ba0SBarry Smith 
323ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3248c5b8ba0SBarry Smith M*/
3258c5b8ba0SBarry Smith 
3261f7e983dSSatish Balay /*MC
327ce9499c7SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the
3288c5b8ba0SBarry Smith        convergence test routine.
3298c5b8ba0SBarry Smith 
3308c5b8ba0SBarry Smith    Level: advanced
3318c5b8ba0SBarry Smith 
3329793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3338c5b8ba0SBarry Smith M*/
3348c5b8ba0SBarry Smith 
3351f7e983dSSatish Balay /*MC
336ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3378c5b8ba0SBarry Smith        convergence test routine.
3388c5b8ba0SBarry Smith 
3398c5b8ba0SBarry Smith    Level: advanced
3408c5b8ba0SBarry Smith 
3419793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3428c5b8ba0SBarry Smith M*/
3438c5b8ba0SBarry Smith 
3441f7e983dSSatish Balay /*MC
345ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
346a3f661c8SBarry Smith        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS
3478c5b8ba0SBarry Smith 
3488c5b8ba0SBarry Smith    Level: advanced
3498c5b8ba0SBarry Smith 
3509793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3518c5b8ba0SBarry Smith M*/
3528c5b8ba0SBarry Smith 
353*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetNormType(KSP,KSPNormType);
354*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetNormType(KSP,KSPNormType*);
355*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetCheckNormIteration(KSP,PetscInt);
356*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetLagNorm(KSP,PetscBool );
3578a4b9c5cSBarry Smith 
3588a4b9c5cSBarry Smith /*E
35928ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
36028ce4d24SBarry Smith          have converged or diverged
36128ce4d24SBarry Smith 
36228ce4d24SBarry Smith    Level: beginner
36328ce4d24SBarry Smith 
3644d0a8057SBarry Smith    Notes: See KSPGetConvergedReason() for explanation of each value
36528ce4d24SBarry Smith 
3664d0a8057SBarry Smith    Developer notes: this must match finclude/petscksp.h
3674d0a8057SBarry Smith 
3684d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
3694d0a8057SBarry Smith       any of the values here also change them that array of names.
37086c02ca4SBarry Smith 
371c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
37228ce4d24SBarry Smith E*/
373d15094e1SBarry Smith typedef enum {/* converged */
3749ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
3759ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
376d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
377d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
378b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3798031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
3808031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
381329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
382af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
383d15094e1SBarry Smith               /* diverged */
384b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
385d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
386d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
387d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
388b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
389b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
390b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
3916aee1118SBarry Smith               KSP_DIVERGED_NAN                 = -9,
3926aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
393d15094e1SBarry Smith 
394d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
3959dcbbd2bSBarry Smith extern const char **KSPConvergedReasons;
396d15094e1SBarry Smith 
397c838673bSBarry Smith /*MC
398c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
399c838673bSBarry Smith 
400c838673bSBarry Smith    Level: beginner
401c838673bSBarry Smith 
402c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
403c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
404c838673bSBarry Smith        2-norm of the residual for right preconditioning
405c838673bSBarry Smith 
406c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
407c838673bSBarry Smith 
408c838673bSBarry Smith M*/
409c838673bSBarry Smith 
410c838673bSBarry Smith /*MC
411c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
412c838673bSBarry Smith 
413c838673bSBarry Smith    Level: beginner
414c838673bSBarry Smith 
415c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
416c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
417c838673bSBarry Smith        2-norm of the residual for right preconditioning
418c838673bSBarry Smith 
419c838673bSBarry Smith    Level: beginner
420c838673bSBarry Smith 
421c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
422c838673bSBarry Smith 
423c838673bSBarry Smith M*/
424c838673bSBarry Smith 
425c838673bSBarry Smith /*MC
426c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
427c838673bSBarry Smith 
428c838673bSBarry Smith    Level: beginner
429c838673bSBarry Smith 
430c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
431c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
432c838673bSBarry Smith        2-norm of the residual for right preconditioning
433c838673bSBarry Smith 
434c838673bSBarry Smith    Level: beginner
435c838673bSBarry Smith 
436c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
437c838673bSBarry Smith 
438c838673bSBarry Smith M*/
439c838673bSBarry Smith 
440c838673bSBarry Smith /*MC
441c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
442c838673bSBarry Smith       reached
443c838673bSBarry Smith 
444c838673bSBarry Smith    Level: beginner
445c838673bSBarry Smith 
446c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
447c838673bSBarry Smith 
448c838673bSBarry Smith M*/
449c838673bSBarry Smith 
450c838673bSBarry Smith /*MC
4518c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
4528c9406f6SLisandro Dalcin            the preconditioner is applied. Also used when the KSPSkipConverged() convergence
4534d0a8057SBarry Smith            test routine is set in KSP.
454c838673bSBarry Smith 
455c838673bSBarry Smith 
456c838673bSBarry Smith    Level: beginner
457c838673bSBarry Smith 
458c838673bSBarry Smith 
459c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
460c838673bSBarry Smith 
461c838673bSBarry Smith M*/
462c838673bSBarry Smith 
463c838673bSBarry Smith /*MC
464c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
4653014e516SBarry Smith           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
4663014e516SBarry Smith           preconditioner.
467c838673bSBarry Smith 
468c838673bSBarry Smith    Level: beginner
469c838673bSBarry Smith 
470c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
471c838673bSBarry Smith 
472c838673bSBarry Smith M*/
473c838673bSBarry Smith 
474c838673bSBarry Smith /*MC
475c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
476c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
477c838673bSBarry Smith 
478c838673bSBarry Smith 
479c838673bSBarry Smith    Level: beginner
480c838673bSBarry Smith 
481c838673bSBarry Smith 
482c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
483c838673bSBarry Smith 
484c838673bSBarry Smith M*/
485c838673bSBarry Smith 
486c838673bSBarry Smith /*MC
487c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
488c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
489c838673bSBarry Smith 
490c838673bSBarry Smith    Level: beginner
491c838673bSBarry Smith 
492c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
493c838673bSBarry Smith 
494c838673bSBarry Smith M*/
495c838673bSBarry Smith 
496c838673bSBarry Smith /*MC
497c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
498c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
499c838673bSBarry Smith         be positive definite
500c838673bSBarry Smith 
501c838673bSBarry Smith    Level: beginner
502c838673bSBarry Smith 
5032401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
504c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
505c838673bSBarry Smith 
506c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
507c838673bSBarry Smith 
508c838673bSBarry Smith M*/
509c838673bSBarry Smith 
510c838673bSBarry Smith /*MC
511c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
512c838673bSBarry Smith         while the KSPSolve() is still running.
513c838673bSBarry Smith 
514c838673bSBarry Smith    Level: beginner
515c838673bSBarry Smith 
516c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
517c838673bSBarry Smith 
518c838673bSBarry Smith M*/
519c838673bSBarry Smith 
520*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
521*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetConvergenceContext(KSP,void **);
522*7087cfbeSBarry Smith extern PetscErrorCode  KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
523*7087cfbeSBarry Smith extern PetscErrorCode  KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
524*7087cfbeSBarry Smith extern PetscErrorCode  KSPDefaultConvergedDestroy(void *);
525*7087cfbeSBarry Smith extern PetscErrorCode  KSPDefaultConvergedCreate(void **);
526*7087cfbeSBarry Smith extern PetscErrorCode  KSPDefaultConvergedSetUIRNorm(KSP);
527*7087cfbeSBarry Smith extern PetscErrorCode  KSPDefaultConvergedSetUMIRNorm(KSP);
528*7087cfbeSBarry Smith extern PetscErrorCode  KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
529*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetConvergedReason(KSP,KSPConvergedReason *);
530abef13c0SSatish Balay 
531*7087cfbeSBarry Smith extern PetscErrorCode  KSPComputeExplicitOperator(KSP,Mat *);
532d4fbbf0eSBarry Smith 
53328ce4d24SBarry Smith /*E
53428ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
53528ce4d24SBarry Smith 
53628ce4d24SBarry Smith    Level: beginner
53728ce4d24SBarry Smith 
53828ce4d24SBarry Smith .seealso: KSPCGSetType()
53928ce4d24SBarry Smith E*/
5409dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5419dcbbd2bSBarry Smith extern const char *KSPCGTypes[];
54228ce4d24SBarry Smith 
543*7087cfbeSBarry Smith extern PetscErrorCode  KSPCGSetType(KSP,KSPCGType);
544*7087cfbeSBarry Smith extern PetscErrorCode  KSPCGUseSingleReduction(KSP,PetscBool );
5458031f4adStmunson 
546*7087cfbeSBarry Smith extern PetscErrorCode  KSPNASHSetRadius(KSP,PetscReal);
547*7087cfbeSBarry Smith extern PetscErrorCode  KSPNASHGetNormD(KSP,PetscReal *);
548*7087cfbeSBarry Smith extern PetscErrorCode  KSPNASHGetObjFcn(KSP,PetscReal *);
549fcae7a14Stmunson 
550*7087cfbeSBarry Smith extern PetscErrorCode  KSPSTCGSetRadius(KSP,PetscReal);
551*7087cfbeSBarry Smith extern PetscErrorCode  KSPSTCGGetNormD(KSP,PetscReal *);
552*7087cfbeSBarry Smith extern PetscErrorCode  KSPSTCGGetObjFcn(KSP,PetscReal *);
553e559a7feSLois Curfman McInnes 
554*7087cfbeSBarry Smith extern PetscErrorCode  KSPGLTRSetRadius(KSP,PetscReal);
555*7087cfbeSBarry Smith extern PetscErrorCode  KSPGLTRGetNormD(KSP,PetscReal *);
556*7087cfbeSBarry Smith extern PetscErrorCode  KSPGLTRGetObjFcn(KSP,PetscReal *);
557*7087cfbeSBarry Smith extern PetscErrorCode  KSPGLTRGetMinEig(KSP,PetscReal *);
558*7087cfbeSBarry Smith extern PetscErrorCode  KSPGLTRGetLambda(KSP,PetscReal *);
5598031f4adStmunson 
560*7087cfbeSBarry Smith extern PetscErrorCode  KSPPythonSetType(KSP,const char[]);
5611d6018f0SLisandro Dalcin 
562*7087cfbeSBarry Smith extern PetscErrorCode  PCPreSolve(PC,KSP);
563*7087cfbeSBarry Smith extern PetscErrorCode  PCPostSolve(PC,KSP);
5643369ce9aSBarry Smith 
565*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
566*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorLG(KSP,PetscInt,PetscReal,void*);
567*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorLGDestroy(PetscDrawLG);
568*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
569*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
570*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG);
571*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
572*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
573*7087cfbeSBarry Smith extern PetscErrorCode  KSPMonitorLGRangeDestroy(PetscDrawLG);
5742f2e5d10SKris Buschelman 
575*7087cfbeSBarry Smith extern PetscErrorCode  PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
576*7087cfbeSBarry Smith extern PetscErrorCode  PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
57703919abeSBarry Smith 
578f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
579ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
580f8a50e2bSBarry Smith 
581*7087cfbeSBarry Smith extern PetscErrorCode  KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
582*7087cfbeSBarry Smith extern PetscErrorCode  KSPFischerGuessDestroy(KSPFischerGuess);
583*7087cfbeSBarry Smith extern PetscErrorCode  KSPFischerGuessReset(KSPFischerGuess);
584*7087cfbeSBarry Smith extern PetscErrorCode  KSPFischerGuessUpdate(KSPFischerGuess,Vec);
585*7087cfbeSBarry Smith extern PetscErrorCode  KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
586*7087cfbeSBarry Smith extern PetscErrorCode  KSPFischerGuessSetFromOptions(KSPFischerGuess);
587f8a50e2bSBarry Smith 
588*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
589*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetFischerGuess(KSP,KSPFischerGuess);
590*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetFischerGuess(KSP,KSPFischerGuess*);
591f8a50e2bSBarry Smith 
592*7087cfbeSBarry Smith extern PetscErrorCode  MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
593*7087cfbeSBarry Smith extern PetscErrorCode  MatSchurComplementGetKSP(Mat,KSP*);
594*7087cfbeSBarry Smith extern PetscErrorCode  MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure);
595*7087cfbeSBarry Smith extern PetscErrorCode  MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
596*7087cfbeSBarry Smith extern PetscErrorCode  MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *);
5973f22127dSBarry Smith 
598*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetDM(KSP,DM);
599*7087cfbeSBarry Smith extern PetscErrorCode  KSPSetDMActive(KSP,PetscBool );
600*7087cfbeSBarry Smith extern PetscErrorCode  KSPGetDM(KSP,DM*);
6016c699258SBarry Smith 
602e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
6032eac72dbSBarry Smith #endif
604