xref: /petsc/include/petscksp.h (revision c2b71004dfc13f0a93945b0cd567a028e08616d3)
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
62c8e378dSBarry Smith #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*/
3019fd82e9SBarry Smith typedef const char* KSPType;
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"
47*c2b71004SHong Zhang #define   KSPFBCGSR     "fbcgsr"
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 *);
6719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,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));
80f550243cSJed Brown PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(const char[]);
8130de9b25SBarry Smith 
8230de9b25SBarry Smith /*MC
8330de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
8430de9b25SBarry Smith 
8530de9b25SBarry Smith    Synopsis:
861890ba74SBarry Smith    PetscErrorCode KSPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(KSP))
8730de9b25SBarry Smith 
8830de9b25SBarry Smith    Not Collective
8930de9b25SBarry Smith 
9030de9b25SBarry Smith    Input Parameters:
9130de9b25SBarry Smith +  name_solver - name of a new user-defined solver
9230de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
9330de9b25SBarry Smith .  name_create - name of routine to create method context
9430de9b25SBarry Smith -  routine_create - routine to create method context
9530de9b25SBarry Smith 
9630de9b25SBarry Smith    Notes:
9730de9b25SBarry Smith    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
9830de9b25SBarry Smith 
9930de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
10030de9b25SBarry Smith    is ignored.
10130de9b25SBarry Smith 
10230de9b25SBarry Smith    Sample usage:
10330de9b25SBarry Smith .vb
10430de9b25SBarry Smith    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
10530de9b25SBarry Smith                "MySolverCreate",MySolverCreate);
10630de9b25SBarry Smith .ve
10730de9b25SBarry Smith 
10830de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
10930de9b25SBarry Smith $     KSPSetType(ksp,"my_solver")
11030de9b25SBarry Smith    or at runtime via the option
11130de9b25SBarry Smith $     -ksp_type my_solver
11230de9b25SBarry Smith 
11330de9b25SBarry Smith    Level: advanced
11430de9b25SBarry Smith 
115ab901514SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
11630de9b25SBarry Smith           and others of the form ${any_environmental_variable} occuring in pathname will be
11730de9b25SBarry Smith           replaced with appropriate values.
11830de9b25SBarry Smith          If your function is not being put into a shared library then use KSPRegister() instead
11930de9b25SBarry Smith 
12030de9b25SBarry Smith .keywords: KSP, register
12130de9b25SBarry Smith 
12230de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy()
12330de9b25SBarry Smith 
12430de9b25SBarry Smith M*/
125aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
126f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
1276df38c32SLois Curfman McInnes #else
128f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
1296df38c32SLois Curfman McInnes #endif
13082bf6240SBarry Smith 
13119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool );
137014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool  *);
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool );
139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *);
140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool );
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool  *);
142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *);
143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool );
144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *);
145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool );
146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace);
151014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*);
152014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1532eac72dbSBarry Smith 
154014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
155014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
156aabeff55SBarry Smith 
157014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
158014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
159014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
160014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
161014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
162014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
1634b0e389bSBarry Smith 
1640e33f6ddSBarry Smith /* not sure where to put this */
165014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
166014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
167014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
168014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
169014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
1702eac72dbSBarry Smith 
171014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *);
1720a71e23eSMatthew Knepley 
173014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
174014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
1752eac72dbSBarry Smith 
176014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
178014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
179014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
180014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP);
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
1844b0e389bSBarry Smith 
185014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
1889f236934SBarry Smith 
189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
193014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1941d73ed98SBarry Smith 
195014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
1971d73ed98SBarry Smith 
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
199014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
200014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
20158ad63f4SBarry Smith 
202b49fd9e1SBarry Smith /*E
203b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
204b49fd9e1SBarry Smith 
205b49fd9e1SBarry Smith    Level: advanced
206b49fd9e1SBarry Smith 
207bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
208e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
209b49fd9e1SBarry Smith 
210b49fd9e1SBarry Smith E*/
21178d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2126a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2131f7e983dSSatish Balay /*MC
2148c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
2158c5b8ba0SBarry Smith 
2168c5b8ba0SBarry Smith    Level: advanced
2178c5b8ba0SBarry Smith 
2188c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2198c5b8ba0SBarry Smith 
220bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
221e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2228c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2238c5b8ba0SBarry Smith M*/
2248c5b8ba0SBarry Smith 
2251f7e983dSSatish Balay /*MC
2268c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2278c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2288c5b8ba0SBarry Smith           poor orthogonality.
2298c5b8ba0SBarry Smith 
2308c5b8ba0SBarry Smith    Level: advanced
2318c5b8ba0SBarry Smith 
2328c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2338c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2348c5b8ba0SBarry Smith 
235bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
236e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2378c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2388c5b8ba0SBarry Smith M*/
2398c5b8ba0SBarry Smith 
2401f7e983dSSatish Balay /*MC
2418c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2428c5b8ba0SBarry Smith 
2438c5b8ba0SBarry Smith    Level: advanced
2448c5b8ba0SBarry Smith 
2458c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2468c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2478c5b8ba0SBarry Smith 
2488c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2498c5b8ba0SBarry Smith 
250bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
251e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2528c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2538c5b8ba0SBarry Smith M*/
2548c5b8ba0SBarry Smith 
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
25708480c60SBarry Smith 
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
261c38d4ed2SBarry Smith 
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
265121fd945SKris Buschelman 
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
269d9492815SBarry Smith 
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2722eac72dbSBarry Smith 
273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
274014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
277390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
278390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
280816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *);
281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMS(KSP,PetscInt,PetscReal,void*);
284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSCreate(KSP,const char*,void**);
285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSDestroy(void**);
286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
28784cb2905SBarry Smith 
288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*);
290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
292c01c455dSBarry Smith 
293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure);
294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
297014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
2997287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
3007287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
3011eb62cbbSBarry Smith 
302014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
303014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
304014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
305014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
3061f7f0c4fSBarry Smith 
307014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
3081eb62cbbSBarry Smith 
309014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
310014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
311db9b2ab1SHong Zhang 
312014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
313014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
31483ab6a24SBarry Smith 
31528ce4d24SBarry Smith /*E
3168a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
3178a4b9c5cSBarry Smith        test routines.
3188a4b9c5cSBarry Smith 
3198a4b9c5cSBarry Smith    Level: advanced
3208a4b9c5cSBarry Smith 
321a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
3221718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
323a3f661c8SBarry Smith 
3248a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
3258a4b9c5cSBarry Smith 
32694b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
3271718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
3288a4b9c5cSBarry Smith E*/
3299e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
3309e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
331014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
332a21b2a99SBarry Smith 
3331f7e983dSSatish Balay /*MC
3349793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
3358c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
3368c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3378c5b8ba0SBarry Smith 
3388c5b8ba0SBarry Smith    Level: advanced
3398c5b8ba0SBarry Smith 
340085a36d4SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
3418c5b8ba0SBarry Smith 
342ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3438c5b8ba0SBarry Smith M*/
3448c5b8ba0SBarry Smith 
3451f7e983dSSatish Balay /*MC
346ce9499c7SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the
3478c5b8ba0SBarry Smith        convergence test routine.
3488c5b8ba0SBarry Smith 
3498c5b8ba0SBarry Smith    Level: advanced
3508c5b8ba0SBarry Smith 
3519793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3528c5b8ba0SBarry Smith M*/
3538c5b8ba0SBarry Smith 
3541f7e983dSSatish Balay /*MC
355ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3568c5b8ba0SBarry Smith        convergence test routine.
3578c5b8ba0SBarry Smith 
3588c5b8ba0SBarry Smith    Level: advanced
3598c5b8ba0SBarry Smith 
3609793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3618c5b8ba0SBarry Smith M*/
3628c5b8ba0SBarry Smith 
3631f7e983dSSatish Balay /*MC
364ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
365a3f661c8SBarry Smith        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS
3668c5b8ba0SBarry Smith 
3678c5b8ba0SBarry Smith    Level: advanced
3688c5b8ba0SBarry Smith 
3699793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3708c5b8ba0SBarry Smith M*/
3718c5b8ba0SBarry Smith 
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
376014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool );
3778a4b9c5cSBarry Smith 
3788a4b9c5cSBarry Smith /*E
37928ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
38028ce4d24SBarry Smith          have converged or diverged
38128ce4d24SBarry Smith 
38228ce4d24SBarry Smith    Level: beginner
38328ce4d24SBarry Smith 
3844d0a8057SBarry Smith    Notes: See KSPGetConvergedReason() for explanation of each value
38528ce4d24SBarry Smith 
3864d0a8057SBarry Smith    Developer notes: this must match finclude/petscksp.h
3874d0a8057SBarry Smith 
3884d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
3894d0a8057SBarry Smith       any of the values here also change them that array of names.
39086c02ca4SBarry Smith 
391c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
39228ce4d24SBarry Smith E*/
393d15094e1SBarry Smith typedef enum {/* converged */
3949ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
3959ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
396d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
397d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
398b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3998031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
4008031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
401329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
402af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
403d15094e1SBarry Smith               /* diverged */
404b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
405d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
406d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
407d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
408b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
409b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
410b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
4116aee1118SBarry Smith               KSP_DIVERGED_NAN                 = -9,
4126aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
413d15094e1SBarry Smith 
414d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
415014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
416d15094e1SBarry Smith 
417c838673bSBarry Smith /*MC
418c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
419c838673bSBarry Smith 
420c838673bSBarry Smith    Level: beginner
421c838673bSBarry Smith 
422c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
423c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
424c838673bSBarry Smith        2-norm of the residual for right preconditioning
425c838673bSBarry Smith 
426c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
427c838673bSBarry Smith 
428c838673bSBarry Smith M*/
429c838673bSBarry Smith 
430c838673bSBarry Smith /*MC
431c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
432c838673bSBarry Smith 
433c838673bSBarry Smith    Level: beginner
434c838673bSBarry Smith 
435c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
436c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
437c838673bSBarry Smith        2-norm of the residual for right preconditioning
438c838673bSBarry Smith 
439c838673bSBarry Smith    Level: beginner
440c838673bSBarry Smith 
441c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
442c838673bSBarry Smith 
443c838673bSBarry Smith M*/
444c838673bSBarry Smith 
445c838673bSBarry Smith /*MC
446c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
447c838673bSBarry Smith 
448c838673bSBarry Smith    Level: beginner
449c838673bSBarry Smith 
450c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
451c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
452c838673bSBarry Smith        2-norm of the residual for right preconditioning
453c838673bSBarry Smith 
454c838673bSBarry Smith    Level: beginner
455c838673bSBarry Smith 
456c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
457c838673bSBarry Smith 
458c838673bSBarry Smith M*/
459c838673bSBarry Smith 
460c838673bSBarry Smith /*MC
461c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
462c838673bSBarry Smith       reached
463c838673bSBarry Smith 
464c838673bSBarry Smith    Level: beginner
465c838673bSBarry Smith 
466c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
467c838673bSBarry Smith 
468c838673bSBarry Smith M*/
469c838673bSBarry Smith 
470c838673bSBarry Smith /*MC
4718c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
4728c9406f6SLisandro Dalcin            the preconditioner is applied. Also used when the KSPSkipConverged() convergence
4734d0a8057SBarry Smith            test routine is set in KSP.
474c838673bSBarry Smith 
475c838673bSBarry Smith 
476c838673bSBarry Smith    Level: beginner
477c838673bSBarry Smith 
478c838673bSBarry Smith 
479c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
480c838673bSBarry Smith 
481c838673bSBarry Smith M*/
482c838673bSBarry Smith 
483c838673bSBarry Smith /*MC
484c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
4853014e516SBarry Smith           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
4863014e516SBarry Smith           preconditioner.
487c838673bSBarry Smith 
488c838673bSBarry Smith    Level: beginner
489c838673bSBarry Smith 
490c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
491c838673bSBarry Smith 
492c838673bSBarry Smith M*/
493c838673bSBarry Smith 
494c838673bSBarry Smith /*MC
495c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
496c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
497c838673bSBarry Smith 
498c838673bSBarry Smith 
499c838673bSBarry Smith    Level: beginner
500c838673bSBarry Smith 
501c838673bSBarry Smith 
502c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
503c838673bSBarry Smith 
504c838673bSBarry Smith M*/
505c838673bSBarry Smith 
506c838673bSBarry Smith /*MC
507c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
508c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
509c838673bSBarry Smith 
510c838673bSBarry Smith    Level: beginner
511c838673bSBarry Smith 
512c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
513c838673bSBarry Smith 
514c838673bSBarry Smith M*/
515c838673bSBarry Smith 
516c838673bSBarry Smith /*MC
517c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
518c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
519c838673bSBarry Smith         be positive definite
520c838673bSBarry Smith 
521c838673bSBarry Smith    Level: beginner
522c838673bSBarry Smith 
5232401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
524c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
525c838673bSBarry Smith 
526c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
527c838673bSBarry Smith 
528c838673bSBarry Smith M*/
529c838673bSBarry Smith 
530c838673bSBarry Smith /*MC
531c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
532c838673bSBarry Smith         while the KSPSolve() is still running.
533c838673bSBarry Smith 
534c838673bSBarry Smith    Level: beginner
535c838673bSBarry Smith 
536c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
537c838673bSBarry Smith 
538c838673bSBarry Smith M*/
539c838673bSBarry Smith 
540014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
542014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
543014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
544014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *);
545014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **);
546014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP);
547014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP);
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
549014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
550abef13c0SSatish Balay 
551014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
552d4fbbf0eSBarry Smith 
55328ce4d24SBarry Smith /*E
55428ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
55528ce4d24SBarry Smith 
55628ce4d24SBarry Smith    Level: beginner
55728ce4d24SBarry Smith 
55828ce4d24SBarry Smith .seealso: KSPCGSetType()
55928ce4d24SBarry Smith E*/
5609dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5616a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
56228ce4d24SBarry Smith 
563014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
564014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
5658031f4adStmunson 
566014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
567014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
568014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
569fcae7a14Stmunson 
570014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
571014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
572014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
573e559a7feSLois Curfman McInnes 
574014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
576014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
577014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
578014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
5798031f4adStmunson 
580014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
5811d6018f0SLisandro Dalcin 
582014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
583014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
5843369ce9aSBarry Smith 
5851d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
5861d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
5871d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*);
588014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
589014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
590014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*);
591014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
5922f2e5d10SKris Buschelman 
593014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
594014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
59503919abeSBarry Smith 
596f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
597ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
598f8a50e2bSBarry Smith 
599014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
600014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
601014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
602014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
603014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
604014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
605f8a50e2bSBarry Smith 
606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
609f8a50e2bSBarry Smith 
610014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
611014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
612d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
6131324cf18SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementSet(Mat,Mat,Mat,Mat,Mat,Mat);
614014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure);
615014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
616014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *);
6173f22127dSBarry Smith 
618014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat);
619fcfd50ebSBarry Smith 
620014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
621014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
622014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
623014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
624014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
625014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*);
626014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
627014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*);
628014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*);
629014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
630014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
6316c699258SBarry Smith 
6322eac72dbSBarry Smith #endif
633