xref: /petsc/include/petscksp.h (revision d4211eb91971e37398c9f1d2162f6d95e03d3b29)
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"
342c8fc646SJed Brown #define KSPGROPPCG    "groppcg"
352c8fc646SJed Brown #define KSPPIPECG     "pipecg"
36df2a969aSvictorle #define   KSPCGNE       "cgne"
37fcae7a14Stmunson #define   KSPNASH       "nash"
3880e17431SMatthew Knepley #define   KSPSTCG       "stcg"
398031f4adStmunson #define   KSPGLTR       "gltr"
4082bf6240SBarry Smith #define KSPGMRES      "gmres"
419dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
429dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
434f8e6cd9SSatish Balay #define   KSPDGMRES     "dgmres"
4461b468f9SJed Brown #define   KSPPGMRES     "pgmres"
4582bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
4682bf6240SBarry Smith #define KSPBCGS       "bcgs"
47e1c61ce8SBarry Smith #define   KSPIBCGS      "ibcgs"
4818ac38e6SHong Zhang #define   KSPFBCGS      "fbcgs"
49c2b71004SHong Zhang #define   KSPFBCGSR     "fbcgsr"
50f0037002Svictorle #define   KSPBCGSL      "bcgsl"
5182bf6240SBarry Smith #define KSPCGS        "cgs"
5282bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
5382bf6240SBarry Smith #define KSPCR         "cr"
542c8fc646SJed Brown #define KSPPIPECR     "pipecr"
5582bf6240SBarry Smith #define KSPLSQR       "lsqr"
5682bf6240SBarry Smith #define KSPPREONLY    "preonly"
5782bf6240SBarry Smith #define KSPQCG        "qcg"
58c9cf9b11SBarry Smith #define KSPBICG       "bicg"
59b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
6001c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
6121356919SSatish Balay #define KSPLCD        "lcd"
621d6018f0SLisandro Dalcin #define KSPPYTHON     "python"
6358ad63f4SBarry Smith #define KSPGCR        "gcr"
64bbce1389SJed Brown #define KSPSPECEST    "specest"
652eac72dbSBarry Smith 
668ba1e511SMatthew Knepley /* Logging support */
67014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID;
685358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID;
698ba1e511SMatthew Knepley 
70014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *);
7119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType);
72014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP);
73014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP);
74014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec);
75014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec);
76014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP);
77014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*);
782eac72dbSBarry Smith 
79140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList;
80014dd563SJed Brown PETSC_EXTERN PetscBool         KSPRegisterAllCalled;
81014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRegisterAll(const char[]);
82014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRegisterDestroy(void);
83014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));
84f550243cSJed Brown PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(const char[]);
8530de9b25SBarry Smith 
8630de9b25SBarry Smith /*MC
8730de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
8830de9b25SBarry Smith 
8930de9b25SBarry Smith    Synopsis:
90f2ba6396SBarry Smith    #include "petscksp.h"
911890ba74SBarry Smith    PetscErrorCode KSPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(KSP))
9230de9b25SBarry Smith 
9330de9b25SBarry Smith    Not Collective
9430de9b25SBarry Smith 
9530de9b25SBarry Smith    Input Parameters:
9630de9b25SBarry Smith +  name_solver - name of a new user-defined solver
9730de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
9830de9b25SBarry Smith .  name_create - name of routine to create method context
9930de9b25SBarry Smith -  routine_create - routine to create method context
10030de9b25SBarry Smith 
10130de9b25SBarry Smith    Notes:
10230de9b25SBarry Smith    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
10330de9b25SBarry Smith 
10430de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
10530de9b25SBarry Smith    is ignored.
10630de9b25SBarry Smith 
10730de9b25SBarry Smith    Sample usage:
10830de9b25SBarry Smith .vb
10930de9b25SBarry Smith    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
11030de9b25SBarry Smith                "MySolverCreate",MySolverCreate);
11130de9b25SBarry Smith .ve
11230de9b25SBarry Smith 
11330de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
11430de9b25SBarry Smith $     KSPSetType(ksp,"my_solver")
11530de9b25SBarry Smith    or at runtime via the option
11630de9b25SBarry Smith $     -ksp_type my_solver
11730de9b25SBarry Smith 
11830de9b25SBarry Smith    Level: advanced
11930de9b25SBarry Smith 
120ab901514SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
12130de9b25SBarry Smith           and others of the form ${any_environmental_variable} occuring in pathname will be
12230de9b25SBarry Smith           replaced with appropriate values.
12330de9b25SBarry Smith          If your function is not being put into a shared library then use KSPRegister() instead
12430de9b25SBarry Smith 
12530de9b25SBarry Smith .keywords: KSP, register
12630de9b25SBarry Smith 
12730de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy()
12830de9b25SBarry Smith 
12930de9b25SBarry Smith M*/
130aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
131f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
1326df38c32SLois Curfman McInnes #else
133f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
1346df38c32SLois Curfman McInnes #endif
13582bf6240SBarry Smith 
13619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *);
137014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide);
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*);
139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool );
142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool  *);
143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool );
144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *);
145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool );
146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool  *);
147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *);
148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool );
149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *);
150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool );
151014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *);
152014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *);
153014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*);
154014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*);
155014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace);
156014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*);
157014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1582eac72dbSBarry Smith 
159*d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
160*d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*);
161*d4211eb9SBarry Smith 
162014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC);
163014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*);
164aabeff55SBarry Smith 
165014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal);
166014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
167014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP);
168014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **);
169014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
170014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool );
1714b0e389bSBarry Smith 
172fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*);
173fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *);
174fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP);
175fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt);
176fa0ddf94SBarry Smith 
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*);
178014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
179014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
180014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
182b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
183b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
184b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
185b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);
186014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *);
1870a71e23eSMatthew Knepley 
188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *);
189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *);
1902eac72dbSBarry Smith 
191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal);
192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool );
193014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal);
194014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal);
195ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom);
196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP);
197014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
199014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
2004b0e389bSBarry Smith 
201014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt);
202014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*);
203014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal);
2049f236934SBarry Smith 
205014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP);
206014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
207014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt));
208014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
209014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
2101d73ed98SBarry Smith 
211014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt);
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP);
2131d73ed98SBarry Smith 
214014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt);
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*);
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
21758ad63f4SBarry Smith 
218b49fd9e1SBarry Smith /*E
219b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
220b49fd9e1SBarry Smith 
221b49fd9e1SBarry Smith    Level: advanced
222b49fd9e1SBarry Smith 
223bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
224e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
225b49fd9e1SBarry Smith 
226b49fd9e1SBarry Smith E*/
22778d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
2286a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[];
2291f7e983dSSatish Balay /*MC
2308c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
2318c5b8ba0SBarry Smith 
2328c5b8ba0SBarry Smith    Level: advanced
2338c5b8ba0SBarry Smith 
2348c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
2358c5b8ba0SBarry Smith 
236bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
237e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2388c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2398c5b8ba0SBarry Smith M*/
2408c5b8ba0SBarry Smith 
2411f7e983dSSatish Balay /*MC
2428c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2438c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2448c5b8ba0SBarry Smith           poor orthogonality.
2458c5b8ba0SBarry Smith 
2468c5b8ba0SBarry Smith    Level: advanced
2478c5b8ba0SBarry Smith 
2488c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2498c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2508c5b8ba0SBarry Smith 
251bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
252e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2538c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2548c5b8ba0SBarry Smith M*/
2558c5b8ba0SBarry Smith 
2561f7e983dSSatish Balay /*MC
2578c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2588c5b8ba0SBarry Smith 
2598c5b8ba0SBarry Smith    Level: advanced
2608c5b8ba0SBarry Smith 
2618c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2628c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2638c5b8ba0SBarry Smith 
2648c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2658c5b8ba0SBarry Smith 
266bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
267e928de20SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2688c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2698c5b8ba0SBarry Smith M*/
2708c5b8ba0SBarry Smith 
271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
272014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*);
27308480c60SBarry Smith 
274014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
277c38d4ed2SBarry Smith 
278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal);
279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*);
280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*);
281121fd945SKris Buschelman 
282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal);
283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool );
284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt);
285e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool);
286d9492815SBarry Smith 
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP);
288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2892eac72dbSBarry Smith 
290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *);
293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *);
294390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy);
295390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy);
296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
297816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *);
298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
299014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMS(KSP,PetscInt,PetscReal,void*);
301014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSCreate(KSP,const char*,void**);
302014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSDestroy(void**);
303014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
30484cb2905SBarry Smith 
305014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec);
306014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
307c01c455dSBarry Smith 
308014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure);
309014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
310014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *);
311014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]);
312014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]);
313014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]);
3147287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt);
3157287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*);
3161eb62cbbSBarry Smith 
317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool );
318014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *);
319014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool );
320014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *);
3211f7f0c4fSBarry Smith 
322014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer);
32355849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer);
32455849f57SBarry Smith 
32555849f57SBarry Smith #define KSP_FILE_CLASSID 1211223
3261eb62cbbSBarry Smith 
327014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec);
328014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*);
329db9b2ab1SHong Zhang 
330014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*);
331014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*);
33283ab6a24SBarry Smith 
33328ce4d24SBarry Smith /*E
3348a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
3358a4b9c5cSBarry Smith        test routines.
3368a4b9c5cSBarry Smith 
3378a4b9c5cSBarry Smith    Level: advanced
3388a4b9c5cSBarry Smith 
339a3f661c8SBarry Smith    Each solver only supports a subset of these and some may support different ones
3401718446fSBarry Smith    depending on left or right preconditioning, see KSPSetPCSide()
341a3f661c8SBarry Smith 
3428a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
3438a4b9c5cSBarry Smith 
34494b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
3451718446fSBarry Smith           KSPSetConvergenceTest(), KSPSetPCSide()
3468a4b9c5cSBarry Smith E*/
3479e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
3489e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1)
349014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes;
350a21b2a99SBarry Smith 
3511f7e983dSSatish Balay /*MC
3529793b452SBarry Smith     KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will
3538c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
3548c5b8ba0SBarry Smith           be based on a norm of a residual etc.
3558c5b8ba0SBarry Smith 
3568c5b8ba0SBarry Smith    Level: advanced
3578c5b8ba0SBarry Smith 
358085a36d4SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
3598c5b8ba0SBarry Smith 
360ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3618c5b8ba0SBarry Smith M*/
3628c5b8ba0SBarry Smith 
3631f7e983dSSatish Balay /*MC
364ce9499c7SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the
3658c5b8ba0SBarry Smith        convergence test routine.
3668c5b8ba0SBarry Smith 
3678c5b8ba0SBarry Smith    Level: advanced
3688c5b8ba0SBarry Smith 
3699793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3708c5b8ba0SBarry Smith M*/
3718c5b8ba0SBarry Smith 
3721f7e983dSSatish Balay /*MC
373ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3748c5b8ba0SBarry Smith        convergence test routine.
3758c5b8ba0SBarry Smith 
3768c5b8ba0SBarry Smith    Level: advanced
3778c5b8ba0SBarry Smith 
3789793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3798c5b8ba0SBarry Smith M*/
3808c5b8ba0SBarry Smith 
3811f7e983dSSatish Balay /*MC
382ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
383a3f661c8SBarry Smith        convergence test routine. This is only supported by  KSPCG, KSPCR, KSPCGNE, KSPCGS
3848c5b8ba0SBarry Smith 
3858c5b8ba0SBarry Smith    Level: advanced
3868c5b8ba0SBarry Smith 
3879793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3888c5b8ba0SBarry Smith M*/
3898c5b8ba0SBarry Smith 
390014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType);
391014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*);
392014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt);
393014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt);
394014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool);
3958a4b9c5cSBarry Smith 
396b30237c6SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
397b30237c6SBarry Smith 
3988a4b9c5cSBarry Smith /*E
39928ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
40028ce4d24SBarry Smith          have converged or diverged
40128ce4d24SBarry Smith 
40228ce4d24SBarry Smith    Level: beginner
40328ce4d24SBarry Smith 
4044d0a8057SBarry Smith    Notes: See KSPGetConvergedReason() for explanation of each value
40528ce4d24SBarry Smith 
4064d0a8057SBarry Smith    Developer notes: this must match finclude/petscksp.h
4074d0a8057SBarry Smith 
4084d0a8057SBarry Smith       The string versions of these are KSPConvergedReasons; if you change
4094d0a8057SBarry Smith       any of the values here also change them that array of names.
41086c02ca4SBarry Smith 
411c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
41228ce4d24SBarry Smith E*/
413d15094e1SBarry Smith typedef enum {/* converged */
4149ecb1286SBarry Smith               KSP_CONVERGED_RTOL_NORMAL        =  1,
4159ecb1286SBarry Smith               KSP_CONVERGED_ATOL_NORMAL        =  9,
416d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
417d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
418b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
4198031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
4208031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
421329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
422af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
423d15094e1SBarry Smith               /* diverged */
424b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
425d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
426d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
427d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
428b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
429b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
430b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
4314d51c080SBarry Smith               KSP_DIVERGED_NANORINF            = -9,
4326aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
433d15094e1SBarry Smith 
434d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
435014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons;
436d15094e1SBarry Smith 
437c838673bSBarry Smith /*MC
438c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
439c838673bSBarry Smith 
440c838673bSBarry Smith    Level: beginner
441c838673bSBarry Smith 
442c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
443c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
444c838673bSBarry Smith        2-norm of the residual for right preconditioning
445c838673bSBarry Smith 
446c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
447c838673bSBarry Smith 
448c838673bSBarry Smith M*/
449c838673bSBarry Smith 
450c838673bSBarry Smith /*MC
451c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
452c838673bSBarry Smith 
453c838673bSBarry Smith    Level: beginner
454c838673bSBarry Smith 
455c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
456c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
457c838673bSBarry Smith        2-norm of the residual for right preconditioning
458c838673bSBarry Smith 
459c838673bSBarry Smith    Level: beginner
460c838673bSBarry Smith 
461c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
462c838673bSBarry Smith 
463c838673bSBarry Smith M*/
464c838673bSBarry Smith 
465c838673bSBarry Smith /*MC
466c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
467c838673bSBarry Smith 
468c838673bSBarry Smith    Level: beginner
469c838673bSBarry Smith 
470c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
471c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
472c838673bSBarry Smith        2-norm of the residual for right preconditioning
473c838673bSBarry Smith 
474c838673bSBarry Smith    Level: beginner
475c838673bSBarry Smith 
476c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
477c838673bSBarry Smith 
478c838673bSBarry Smith M*/
479c838673bSBarry Smith 
480c838673bSBarry Smith /*MC
481c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
482c838673bSBarry Smith       reached
483c838673bSBarry Smith 
484c838673bSBarry Smith    Level: beginner
485c838673bSBarry Smith 
486c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
487c838673bSBarry Smith 
488c838673bSBarry Smith M*/
489c838673bSBarry Smith 
490c838673bSBarry Smith /*MC
4918c9406f6SLisandro Dalcin      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of
4928c9406f6SLisandro Dalcin            the preconditioner is applied. Also used when the KSPSkipConverged() convergence
4934d0a8057SBarry Smith            test routine is set in KSP.
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_BREAKDOWN - A breakdown in the Krylov method was detected so the
5053014e516SBarry Smith           method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or
5063014e516SBarry Smith           preconditioner.
507c838673bSBarry Smith 
508c838673bSBarry Smith    Level: beginner
509c838673bSBarry Smith 
510c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
511c838673bSBarry Smith 
512c838673bSBarry Smith M*/
513c838673bSBarry Smith 
514c838673bSBarry Smith /*MC
515c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
516c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
517c838673bSBarry Smith 
518c838673bSBarry Smith 
519c838673bSBarry Smith    Level: beginner
520c838673bSBarry Smith 
521c838673bSBarry Smith 
522c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
523c838673bSBarry Smith 
524c838673bSBarry Smith M*/
525c838673bSBarry Smith 
526c838673bSBarry Smith /*MC
527c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
528c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
529c838673bSBarry Smith 
530c838673bSBarry Smith    Level: beginner
531c838673bSBarry Smith 
532c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
533c838673bSBarry Smith 
534c838673bSBarry Smith M*/
535c838673bSBarry Smith 
536c838673bSBarry Smith /*MC
537c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
538c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
539c838673bSBarry Smith         be positive definite
540c838673bSBarry Smith 
541c838673bSBarry Smith    Level: beginner
542c838673bSBarry Smith 
5432401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
544c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
545c838673bSBarry Smith 
546c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
547c838673bSBarry Smith 
548c838673bSBarry Smith M*/
549c838673bSBarry Smith 
550c838673bSBarry Smith /*MC
551c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
552c838673bSBarry Smith         while the KSPSolve() is still running.
553c838673bSBarry Smith 
554c838673bSBarry Smith    Level: beginner
555c838673bSBarry Smith 
556c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
557c838673bSBarry Smith 
558c838673bSBarry Smith M*/
559c838673bSBarry Smith 
560014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*));
561014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **);
562014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
563014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
564014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *);
565014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **);
566014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP);
567014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP);
568014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
569014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *);
570abef13c0SSatish Balay 
571014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *);
572d4fbbf0eSBarry Smith 
57328ce4d24SBarry Smith /*E
57428ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
57528ce4d24SBarry Smith 
57628ce4d24SBarry Smith    Level: beginner
57728ce4d24SBarry Smith 
57828ce4d24SBarry Smith .seealso: KSPCGSetType()
57928ce4d24SBarry Smith E*/
5809dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5816a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[];
58228ce4d24SBarry Smith 
583014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType);
584014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool );
5858031f4adStmunson 
586014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal);
587014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *);
588014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *);
589fcae7a14Stmunson 
590014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal);
591014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *);
592014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *);
593e559a7feSLois Curfman McInnes 
594014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal);
595014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *);
596014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *);
597014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *);
598014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *);
5998031f4adStmunson 
600014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]);
6011d6018f0SLisandro Dalcin 
602014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP);
603014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP);
6043369ce9aSBarry Smith 
6059804daf3SBarry Smith #include <petscdrawtypes.h>
6061d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
6071d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*);
6081d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*);
609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
611014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*);
612014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*);
6132f2e5d10SKris Buschelman 
614014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
615014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec));
61603919abeSBarry Smith 
617f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */
618ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool  monitor;Mat mat; KSP ksp;}* KSPFischerGuess;
619f8a50e2bSBarry Smith 
620014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*);
621014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*);
622014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess);
623014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec);
624014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec);
625014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess);
626f8a50e2bSBarry Smith 
627014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt);
628014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess);
629014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*);
630f8a50e2bSBarry Smith 
631014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*);
632014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*);
633d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP);
6341324cf18SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementSet(Mat,Mat,Mat,Mat,Mat,Mat);
635014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure);
636014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*);
637014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *);
6383f22127dSBarry Smith 
639014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM);
640014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool );
641014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*);
642014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*);
643014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*);
644fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *);
645014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*);
646fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*);
647014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*);
648014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*);
649014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
650014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
651fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*);
652fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*);
6536c699258SBarry Smith 
6542eac72dbSBarry Smith #endif
655