xref: /petsc/include/petscksp.h (revision 6a6fc655e6b1b6c3ff84be1691b1c5064ef5f4c2)
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"
72eac72dbSBarry Smith 
8014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitializePackage(const char[]);
91dbb0a54SBarry Smith 
1028ce4d24SBarry Smith /*S
1128ce4d24SBarry Smith      KSP - Abstract PETSc object that manages all Krylov methods
1228ce4d24SBarry Smith 
1328ce4d24SBarry Smith    Level: beginner
1428ce4d24SBarry Smith 
1528ce4d24SBarry Smith   Concepts: Krylov methods
1628ce4d24SBarry Smith 
1794b7f48cSBarry Smith .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP
1828ce4d24SBarry Smith S*/
19e2a1c21fSSatish Balay typedef struct _p_KSP*     KSP;
202eac72dbSBarry Smith 
2176bdecfbSBarry Smith /*J
2228ce4d24SBarry Smith     KSPType - String with the name of a PETSc Krylov method or the creation function
2328ce4d24SBarry Smith        with an optional dynamic library name, for example
2428ce4d24SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mykspcreate()
2528ce4d24SBarry Smith 
2628ce4d24SBarry Smith    Level: beginner
2728ce4d24SBarry Smith 
2828ce4d24SBarry Smith .seealso: KSPSetType(), KSP
2976bdecfbSBarry Smith J*/
30a313700dSBarry Smith #define KSPType char*
3182bf6240SBarry Smith #define KSPRICHARDSON "richardson"
326c9de887SHong Zhang #define KSPCHEBYSHEV  "chebyshev"
3382bf6240SBarry Smith #define KSPCG         "cg"
34df2a969aSvictorle #define   KSPCGNE       "cgne"
35fcae7a14Stmunson #define   KSPNASH       "nash"
3680e17431SMatthew Knepley #define   KSPSTCG       "stcg"
378031f4adStmunson #define   KSPGLTR       "gltr"
3882bf6240SBarry Smith #define KSPGMRES      "gmres"
399dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
409dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
414f8e6cd9SSatish Balay #define   KSPDGMRES     "dgmres"
4261b468f9SJed Brown #define   KSPPGMRES     "pgmres"
4382bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
4482bf6240SBarry Smith #define KSPBCGS       "bcgs"
45e1c61ce8SBarry Smith #define   KSPIBCGS      "ibcgs"
4618ac38e6SHong Zhang #define   KSPFBCGS      "fbcgs"
4718ac38e6SHong Zhang #define   KSPIFBCGS     "ifbcgs"
48f0037002Svictorle #define   KSPBCGSL      "bcgsl"
4982bf6240SBarry Smith #define KSPCGS        "cgs"
5082bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
5182bf6240SBarry Smith #define KSPCR         "cr"
5282bf6240SBarry Smith #define KSPLSQR       "lsqr"
5382bf6240SBarry Smith #define KSPPREONLY    "preonly"
5482bf6240SBarry Smith #define KSPQCG        "qcg"
55c9cf9b11SBarry Smith #define KSPBICG       "bicg"
56b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
5701c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
5821356919SSatish Balay #define KSPLCD        "lcd"
591d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
6058ad63f4SBarry Smith #define KSPGCR        "gcr"
61bbce1389SJed Brown #define KSPSPECEST    "specest"
622eac72dbSBarry Smith 
638ba1e511SMatthew Knepley /* Logging support */
64014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
658ba1e511SMatthew Knepley 
66014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
67014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetType(KSP,const KSPType);
68014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
69014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
70014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
71014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
72014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
73014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
742eac72dbSBarry Smith 
75014dd563SJed Brown PETSC_EXTERN PetscFList KSPList;
76014dd563SJed Brown PETSC_EXTERN PetscBool KSPRegisterAllCalled;
77014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRegisterAll(const char[]);
78014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRegisterDestroy(void);
79014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));
8030de9b25SBarry Smith 
8130de9b25SBarry Smith /*MC
8230de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
8330de9b25SBarry Smith 
8430de9b25SBarry Smith    Synopsis:
851890ba74SBarry Smith    PetscErrorCode KSPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(KSP))
8630de9b25SBarry Smith 
8730de9b25SBarry Smith    Not Collective
8830de9b25SBarry Smith 
8930de9b25SBarry Smith    Input Parameters:
9030de9b25SBarry Smith +  name_solver - name of a new user-defined solver
9130de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
9230de9b25SBarry Smith .  name_create - name of routine to create method context
9330de9b25SBarry Smith -  routine_create - routine to create method context
9430de9b25SBarry Smith 
9530de9b25SBarry Smith    Notes:
9630de9b25SBarry Smith    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
9730de9b25SBarry Smith 
9830de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
9930de9b25SBarry Smith    is ignored.
10030de9b25SBarry Smith 
10130de9b25SBarry Smith    Sample usage:
10230de9b25SBarry Smith .vb
10330de9b25SBarry Smith    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
10430de9b25SBarry Smith                "MySolverCreate",MySolverCreate);
10530de9b25SBarry Smith .ve
10630de9b25SBarry Smith 
10730de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
10830de9b25SBarry Smith $     KSPSetType(ksp,"my_solver")
10930de9b25SBarry Smith    or at runtime via the option
11030de9b25SBarry Smith $     -ksp_type my_solver
11130de9b25SBarry Smith 
11230de9b25SBarry Smith    Level: advanced
11330de9b25SBarry Smith 
114ab901514SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
11530de9b25SBarry Smith           and others of the form ${any_environmental_variable} occuring in pathname will be
11630de9b25SBarry Smith           replaced with appropriate values.
11730de9b25SBarry Smith          If your function is not being put into a shared library then use KSPRegister() instead
11830de9b25SBarry Smith 
11930de9b25SBarry Smith .keywords: KSP, register
12030de9b25SBarry Smith 
12130de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy()
12230de9b25SBarry Smith 
12330de9b25SBarry Smith M*/
124aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
125f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
1266df38c32SLois Curfman McInnes #else
127f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
1286df38c32SLois Curfman McInnes #endif
12982bf6240SBarry Smith 
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetType(KSP,const KSPType *);
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool );
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool  *);
137014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool );
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *);
139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool );
140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool  *);
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *);
142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool );
143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *);
144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool );
145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace);
150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*);
151014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1522eac72dbSBarry Smith 
153014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
154014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
155aabeff55SBarry Smith 
156014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
157014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
158014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
159014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
160014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
161014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
1624b0e389bSBarry Smith 
1630e33f6ddSBarry Smith /* not sure where to put this */
164014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
165014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
166014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
167014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
168014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
1692eac72dbSBarry Smith 
170014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *);
1710a71e23eSMatthew Knepley 
172014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
173014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
1742eac72dbSBarry Smith 
175014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
176014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
178014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
179014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP);
180014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
1834b0e389bSBarry Smith 
184014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
185014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
1879f236934SBarry Smith 
188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1931d73ed98SBarry Smith 
194014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
195014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
1961d73ed98SBarry Smith 
197014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
199014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
20058ad63f4SBarry Smith 
201b49fd9e1SBarry Smith /*E
202b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
203b49fd9e1SBarry Smith 
204b49fd9e1SBarry Smith    Level: advanced
205b49fd9e1SBarry Smith 
206bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
207e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
208b49fd9e1SBarry Smith 
209b49fd9e1SBarry Smith E*/
21078d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
211*6a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2121f7e983dSSatish Balay /*MC
2138c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
2148c5b8ba0SBarry Smith 
2158c5b8ba0SBarry Smith    Level: advanced
2168c5b8ba0SBarry Smith 
2178c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2188c5b8ba0SBarry Smith 
219bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
220e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2218c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2228c5b8ba0SBarry Smith M*/
2238c5b8ba0SBarry Smith 
2241f7e983dSSatish Balay /*MC
2258c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2268c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2278c5b8ba0SBarry Smith           poor orthogonality.
2288c5b8ba0SBarry Smith 
2298c5b8ba0SBarry Smith    Level: advanced
2308c5b8ba0SBarry Smith 
2318c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2328c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2338c5b8ba0SBarry Smith 
234bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
235e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2368c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2378c5b8ba0SBarry Smith M*/
2388c5b8ba0SBarry Smith 
2391f7e983dSSatish Balay /*MC
2408c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2418c5b8ba0SBarry Smith 
2428c5b8ba0SBarry Smith    Level: advanced
2438c5b8ba0SBarry Smith 
2448c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2458c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2468c5b8ba0SBarry Smith 
2478c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2488c5b8ba0SBarry Smith 
249bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
250e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2518c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2528c5b8ba0SBarry Smith M*/
2538c5b8ba0SBarry Smith 
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
25608480c60SBarry Smith 
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
260c38d4ed2SBarry Smith 
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
264121fd945SKris Buschelman 
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
268d9492815SBarry Smith 
269014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2712eac72dbSBarry Smith 
272014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
274014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMS(KSP,PetscInt,PetscReal,void*);
280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSCreate(KSP,const char*,void**);
281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSDestroy(void**);
282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
28384cb2905SBarry Smith 
284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*);
286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
288c01c455dSBarry Smith 
289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure);
290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
2957287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPIncrementTabLevel(KSP,PetscObject,PetscInt);
2967287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
2977287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
2981eb62cbbSBarry Smith 
299014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
301014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
302014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
3031f7f0c4fSBarry Smith 
304014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
3051eb62cbbSBarry Smith 
306014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
307014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
308db9b2ab1SHong Zhang 
309014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
310014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
31183ab6a24SBarry Smith 
31228ce4d24SBarry Smith /*E
3138a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
3148a4b9c5cSBarry Smith        test routines.
3158a4b9c5cSBarry Smith 
3168a4b9c5cSBarry Smith    Level: advanced
3178a4b9c5cSBarry Smith 
318a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
3191718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
320a3f661c8SBarry Smith 
3218a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
3228a4b9c5cSBarry Smith 
32394b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
3241718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
3258a4b9c5cSBarry Smith E*/
3269e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
3279e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
328014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
329a21b2a99SBarry Smith 
3301f7e983dSSatish Balay /*MC
3319793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
3328c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
3338c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3348c5b8ba0SBarry Smith 
3358c5b8ba0SBarry Smith    Level: advanced
3368c5b8ba0SBarry Smith 
337085a36d4SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
3388c5b8ba0SBarry Smith 
339ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3408c5b8ba0SBarry Smith M*/
3418c5b8ba0SBarry Smith 
3421f7e983dSSatish Balay /*MC
343ce9499c7SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the
3448c5b8ba0SBarry Smith        convergence test routine.
3458c5b8ba0SBarry Smith 
3468c5b8ba0SBarry Smith    Level: advanced
3478c5b8ba0SBarry Smith 
3489793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3498c5b8ba0SBarry Smith M*/
3508c5b8ba0SBarry Smith 
3511f7e983dSSatish Balay /*MC
352ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3538c5b8ba0SBarry Smith        convergence test routine.
3548c5b8ba0SBarry Smith 
3558c5b8ba0SBarry Smith    Level: advanced
3568c5b8ba0SBarry Smith 
3579793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3588c5b8ba0SBarry Smith M*/
3598c5b8ba0SBarry Smith 
3601f7e983dSSatish Balay /*MC
361ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
362a3f661c8SBarry Smith        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS
3638c5b8ba0SBarry Smith 
3648c5b8ba0SBarry Smith    Level: advanced
3658c5b8ba0SBarry Smith 
3669793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3678c5b8ba0SBarry Smith M*/
3688c5b8ba0SBarry Smith 
369014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool );
3748a4b9c5cSBarry Smith 
3758a4b9c5cSBarry Smith /*E
37628ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
37728ce4d24SBarry Smith          have converged or diverged
37828ce4d24SBarry Smith 
37928ce4d24SBarry Smith    Level: beginner
38028ce4d24SBarry Smith 
3814d0a8057SBarry Smith    Notes: See KSPGetConvergedReason() for explanation of each value
38228ce4d24SBarry Smith 
3834d0a8057SBarry Smith    Developer notes: this must match finclude/petscksp.h
3844d0a8057SBarry Smith 
3854d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
3864d0a8057SBarry Smith       any of the values here also change them that array of names.
38786c02ca4SBarry Smith 
388c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
38928ce4d24SBarry Smith E*/
390d15094e1SBarry Smith typedef enum {/* converged */
3919ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
3929ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
393d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
394d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
395b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3968031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
3978031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
398329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
399af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
400d15094e1SBarry Smith               /* diverged */
401b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
402d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
403d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
404d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
405b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
406b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
407b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
4086aee1118SBarry Smith               KSP_DIVERGED_NAN                 = -9,
4096aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
410d15094e1SBarry Smith 
411d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
412014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
413d15094e1SBarry Smith 
414c838673bSBarry Smith /*MC
415c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*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 .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
424c838673bSBarry Smith 
425c838673bSBarry Smith M*/
426c838673bSBarry Smith 
427c838673bSBarry Smith /*MC
428c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
429c838673bSBarry Smith 
430c838673bSBarry Smith    Level: beginner
431c838673bSBarry Smith 
432c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
433c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
434c838673bSBarry Smith        2-norm of the residual for right preconditioning
435c838673bSBarry Smith 
436c838673bSBarry Smith    Level: beginner
437c838673bSBarry Smith 
438c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
439c838673bSBarry Smith 
440c838673bSBarry Smith M*/
441c838673bSBarry Smith 
442c838673bSBarry Smith /*MC
443c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
444c838673bSBarry Smith 
445c838673bSBarry Smith    Level: beginner
446c838673bSBarry Smith 
447c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
448c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
449c838673bSBarry Smith        2-norm of the residual for right preconditioning
450c838673bSBarry Smith 
451c838673bSBarry Smith    Level: beginner
452c838673bSBarry Smith 
453c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
454c838673bSBarry Smith 
455c838673bSBarry Smith M*/
456c838673bSBarry Smith 
457c838673bSBarry Smith /*MC
458c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
459c838673bSBarry Smith       reached
460c838673bSBarry Smith 
461c838673bSBarry Smith    Level: beginner
462c838673bSBarry Smith 
463c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
464c838673bSBarry Smith 
465c838673bSBarry Smith M*/
466c838673bSBarry Smith 
467c838673bSBarry Smith /*MC
4688c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
4698c9406f6SLisandro Dalcin            the preconditioner is applied. Also used when the KSPSkipConverged() convergence
4704d0a8057SBarry Smith            test routine is set in KSP.
471c838673bSBarry Smith 
472c838673bSBarry Smith 
473c838673bSBarry Smith    Level: beginner
474c838673bSBarry Smith 
475c838673bSBarry Smith 
476c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
477c838673bSBarry Smith 
478c838673bSBarry Smith M*/
479c838673bSBarry Smith 
480c838673bSBarry Smith /*MC
481c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
4823014e516SBarry Smith           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
4833014e516SBarry Smith           preconditioner.
484c838673bSBarry Smith 
485c838673bSBarry Smith    Level: beginner
486c838673bSBarry Smith 
487c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
488c838673bSBarry Smith 
489c838673bSBarry Smith M*/
490c838673bSBarry Smith 
491c838673bSBarry Smith /*MC
492c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
493c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
494c838673bSBarry Smith 
495c838673bSBarry Smith 
496c838673bSBarry Smith    Level: beginner
497c838673bSBarry Smith 
498c838673bSBarry Smith 
499c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
500c838673bSBarry Smith 
501c838673bSBarry Smith M*/
502c838673bSBarry Smith 
503c838673bSBarry Smith /*MC
504c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
505c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
506c838673bSBarry Smith 
507c838673bSBarry Smith    Level: beginner
508c838673bSBarry Smith 
509c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
510c838673bSBarry Smith 
511c838673bSBarry Smith M*/
512c838673bSBarry Smith 
513c838673bSBarry Smith /*MC
514c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
515c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
516c838673bSBarry Smith         be positive definite
517c838673bSBarry Smith 
518c838673bSBarry Smith    Level: beginner
519c838673bSBarry Smith 
5202401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
521c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
522c838673bSBarry Smith 
523c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
524c838673bSBarry Smith 
525c838673bSBarry Smith M*/
526c838673bSBarry Smith 
527c838673bSBarry Smith /*MC
528c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
529c838673bSBarry Smith         while the KSPSolve() is still running.
530c838673bSBarry Smith 
531c838673bSBarry Smith    Level: beginner
532c838673bSBarry Smith 
533c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
534c838673bSBarry Smith 
535c838673bSBarry Smith M*/
536c838673bSBarry Smith 
537014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
538014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
540014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *);
542014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **);
543014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP);
544014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP);
545014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
546014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
547abef13c0SSatish Balay 
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
549d4fbbf0eSBarry Smith 
55028ce4d24SBarry Smith /*E
55128ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
55228ce4d24SBarry Smith 
55328ce4d24SBarry Smith    Level: beginner
55428ce4d24SBarry Smith 
55528ce4d24SBarry Smith .seealso: KSPCGSetType()
55628ce4d24SBarry Smith E*/
5579dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
558*6a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
55928ce4d24SBarry Smith 
560014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
561014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
5628031f4adStmunson 
563014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
564014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
565014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
566fcae7a14Stmunson 
567014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
568014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
569014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
570e559a7feSLois Curfman McInnes 
571014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
572014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
573014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
574014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
5768031f4adStmunson 
577014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
5781d6018f0SLisandro Dalcin 
579014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
580014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
5813369ce9aSBarry Smith 
582014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
583014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLG(KSP,PetscInt,PetscReal,void*);
584014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGDestroy(PetscDrawLG*);
585014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
586014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
587014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*);
588014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
589014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
590014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRangeDestroy(PetscDrawLG*);
5912f2e5d10SKris Buschelman 
592014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
593014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
59403919abeSBarry Smith 
595f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
596ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
597f8a50e2bSBarry Smith 
598014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
599014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
600014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
601014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
602014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
603014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
604f8a50e2bSBarry Smith 
605014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
608f8a50e2bSBarry Smith 
609014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
610014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
611014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure);
612014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
613014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *);
6143f22127dSBarry Smith 
615014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat);
616fcfd50ebSBarry Smith 
617014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
618014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
619014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
620014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
621014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
622014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*);
623014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
624014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*);
625014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*);
626014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
627014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
6286c699258SBarry Smith 
6292eac72dbSBarry Smith #endif
630