xref: /petsc/include/petscksp.h (revision 639f9d9dbbc54d6ac4e42e98283c540b41bb2cee)
1 /* $Id: ksp.h,v 1.45 1996/09/14 03:10:32 bsmith Exp bsmith $ */
2 /*
3    Defines the interface functions for the Krylov subspace accelerators.
4 */
5 #ifndef __KSP_PACKAGE
6 #define __KSP_PACKAGE
7 #include "petsc.h"
8 #include "vec.h"
9 #include "mat.h"
10 #include "pc.h"
11 
12 #define KSP_COOKIE  PETSC_COOKIE+8
13 
14 typedef struct _KSP*     KSP;
15 
16 typedef enum { KSPRICHARDSON, KSPCHEBYCHEV, KSPCG, KSPGMRES, KSPTCQMR, KSPBCGS,
17                KSPCGS, KSPTFQMR, KSPCR, KSPLSQR, KSPPREONLY, KSPQCG} KSPType;
18 
19 extern int KSPCreate(MPI_Comm,KSP *);
20 extern int KSPSetType(KSP,KSPType);
21 extern int KSPSetUp(KSP);
22 extern int KSPSolve(KSP,int *);
23 extern int KSPDestroy(KSP);
24 
25 extern int KSPRegisterAll();
26 extern int KSPRegisterDestroy();
27 extern int KSPRegister(KSPType,char *,int (*)(KSP));
28 
29 extern int KSPGetType(KSP, KSPType *,char **);
30 extern int KSPSetPreconditionerSide(KSP,PCSide);
31 extern int KSPGetPreconditionerSide(KSP,PCSide*);
32 extern int KSPGetTolerances(KSP,double*,double*,double*,int*);
33 extern int KSPSetTolerances(KSP,double,double,double,int);
34 extern int KSPSetComputeResidual(KSP,PetscTruth);
35 extern int KSPSetUsePreconditionedResidual(KSP);
36 extern int KSPSetInitialGuessNonzero(KSP);
37 #define KSPSetComputeEigenvalues(a) KSPSetComputeSingularValues(a)
38 extern int KSPSetComputeSingularValues(KSP);
39 extern int KSPSetRhs(KSP,Vec);
40 extern int KSPGetRhs(KSP,Vec *);
41 extern int KSPSetSolution(KSP,Vec);
42 extern int KSPGetSolution(KSP,Vec *);
43 
44 extern int KSPSetPC(KSP,PC);
45 extern int KSPGetPC(KSP,PC*);
46 
47 extern int KSPSetMonitor(KSP,int (*)(KSP,int,double, void*), void *);
48 extern int KSPGetMonitorContext(KSP,void **);
49 extern int KSPSetResidualHistory(KSP, double *,int);
50 
51 extern int KSPSetConvergenceTest(KSP,int (*)(KSP,int,double, void*), void *);
52 extern int KSPGetConvergenceContext(KSP,void **);
53 
54 extern int KSPBuildSolution(KSP, Vec,Vec *);
55 extern int KSPBuildResidual(KSP, Vec, Vec,Vec *);
56 
57 extern int KSPRichardsonSetScale(KSP , double);
58 extern int KSPChebychevSetEigenvalues(KSP , double, double);
59 extern int KSPComputeExtremeSingularValues(KSP, double*,double*);
60 extern int KSPComputeEigenvalues(KSP,int,double*,double*);
61 extern int KSPComputeEigenvaluesExplicitly(KSP,int,double*,double*);
62 
63 extern int KSPGMRESSetRestart(KSP, int);
64 extern int KSPGMRESSetPreAllocateVectors(KSP);
65 extern int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int));
66 extern int KSPGMRESUnmodifiedGramSchmidtOrthogonalization(KSP,int);
67 extern int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int);
68 extern int KSPGMRESIROrthogonalization(KSP,int);
69 
70 extern int KSPSetFromOptions(KSP);
71 extern int KSPAddOptionsChecker(int (*)(KSP) );
72 
73 extern int KSPSingularValueMonitor(KSP,int,double, void * );
74 extern int KSPCGDefaultConverged(KSP,int,double, void *);
75 extern int KSPDefaultMonitor(KSP,int,double, void *);
76 extern int KSPTrueMonitor(KSP,int,double, void *);
77 extern int KSPDefaultSMonitor(KSP,int,double, void *);
78 extern int KSPDefaultConverged(KSP,int,double, void *);
79 
80 extern int KSPResidual(KSP,Vec,Vec,Vec,Vec,Vec,Vec);
81 extern int KSPUnwindPreconditioner(KSP,Vec,Vec);
82 extern int KSPDefaultBuildSolution(KSP,Vec,Vec*);
83 extern int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
84 
85 extern int KSPPrintHelp(KSP);
86 
87 extern int KSPSetOptionsPrefix(KSP,char*);
88 extern int KSPAppendOptionsPrefix(KSP,char*);
89 extern int KSPGetOptionsPrefix(KSP,char**);
90 
91 extern int KSPView(KSP,Viewer);
92 
93 extern int KSPComputeExplicitOperator(KSP,Mat *);
94 
95 typedef enum {KSP_CG_SYMMETRIC=1, KSP_CG_HERMITIAN=2} KSPCGType;
96 extern int KSPCGSetType(KSP,KSPCGType);
97 
98 #if defined(__DRAW_PACKAGE)
99 extern int KSPLGMonitorCreate(char*,char*,int,int,int,int,DrawLG*);
100 extern int KSPLGMonitor(KSP,int,double,void*);
101 extern int KSPLGMonitorDestroy(DrawLG);
102 extern int KSPLGTrueMonitorCreate(MPI_Comm,char*,char*,int,int,int,int,DrawLG*);
103 extern int KSPLGTrueMonitor(KSP,int,double,void*);
104 extern int KSPLGTrueMonitorDestroy(DrawLG);
105 #endif
106 
107 #endif
108 
109 
110