xref: /petsc/include/petscksp.h (revision c38d4ed214221df9ea04de46f7761bef149d00ff)
1 /* $Id: ksp.h,v 1.82 1999/11/05 14:48:27 bsmith Exp bsmith $ */
2 /*
3    Defines the interface functions for the Krylov subspace accelerators.
4 */
5 #ifndef __KSP_H
6 #define __KSP_H
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 _p_KSP*     KSP;
15 
16 #define KSPRICHARDSON "richardson"
17 #define KSPCHEBYCHEV  "chebychev"
18 #define KSPCG         "cg"
19 #define KSPGMRES      "gmres"
20 #define KSPTCQMR      "tcqmr"
21 #define KSPBCGS       "bcgs"
22 #define KSPCGS        "cgs"
23 #define KSPTFQMR      "tfqmr"
24 #define KSPCR         "cr"
25 #define KSPLSQR       "lsqr"
26 #define KSPPREONLY    "preonly"
27 #define KSPQCG        "qcg"
28 #define KSPBICG       "bicg"
29 #define KSPFGMRES     "fgmres"
30 typedef char * KSPType;
31 
32 extern int KSPCreate(MPI_Comm,KSP *);
33 extern int KSPSetType(KSP,KSPType);
34 extern int KSPSetUp(KSP);
35 extern int KSPSolve(KSP,int *);
36 extern int KSPSolveTrans(KSP,int *);
37 extern int KSPDestroy(KSP);
38 
39 extern FList KSPList;
40 extern int KSPRegisterAll(char *);
41 extern int KSPRegisterDestroy(void);
42 
43 extern int KSPRegister(char*,char*,char*,int(*)(KSP));
44 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
45 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
46 #else
47 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
48 #endif
49 
50 extern int KSPGetType(KSP, KSPType *);
51 extern int KSPSetPreconditionerSide(KSP,PCSide);
52 extern int KSPGetPreconditionerSide(KSP,PCSide*);
53 extern int KSPGetTolerances(KSP,double*,double*,double*,int*);
54 extern int KSPSetTolerances(KSP,double,double,double,int);
55 extern int KSPSetComputeResidual(KSP,PetscTruth);
56 extern int KSPSetUsePreconditionedResidual(KSP);
57 extern int KSPSetInitialGuessNonzero(KSP);
58 extern int KSPGetInitialGuessNonzero(KSP,PetscTruth *);
59 extern int KSPSetComputeEigenvalues(KSP);
60 extern int KSPSetComputeSingularValues(KSP);
61 extern int KSPSetRhs(KSP,Vec);
62 extern int KSPGetRhs(KSP,Vec *);
63 extern int KSPSetSolution(KSP,Vec);
64 extern int KSPGetSolution(KSP,Vec *);
65 extern int KSPGetResidualNorm(KSP,double*);
66 extern int KSPGetIterationNumber(KSP,int*);
67 
68 extern int KSPSetPC(KSP,PC);
69 extern int KSPGetPC(KSP,PC*);
70 
71 extern int KSPSetAvoidNorms(KSP);
72 
73 extern int KSPSetMonitor(KSP,int (*)(KSP,int,double, void*), void *,int (*)(void*));
74 extern int KSPClearMonitor(KSP);
75 extern int KSPGetMonitorContext(KSP,void **);
76 extern int KSPGetResidualHistory(KSP, double **, int *);
77 extern int KSPSetResidualHistory(KSP, double *,int,PetscTruth);
78 
79 extern int KSPSetConvergenceTest(KSP,int (*)(KSP,int,double, void*), void *);
80 extern int KSPGetConvergenceContext(KSP,void **);
81 
82 extern int KSPBuildSolution(KSP, Vec,Vec *);
83 extern int KSPBuildResidual(KSP, Vec, Vec,Vec *);
84 
85 extern int KSPRichardsonSetScale(KSP , double);
86 extern int KSPChebychevSetEigenvalues(KSP , double, double);
87 extern int KSPComputeExtremeSingularValues(KSP, double*,double*);
88 extern int KSPComputeEigenvalues(KSP,int,double*,double*,int *);
89 extern int KSPComputeEigenvaluesExplicitly(KSP,int,double*,double*);
90 
91 extern int KSPGMRESSetRestart(KSP, int);
92 extern int KSPGMRESSetPreAllocateVectors(KSP);
93 extern int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int));
94 extern int KSPGMRESUnmodifiedGramSchmidtOrthogonalization(KSP,int);
95 extern int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int);
96 extern int KSPGMRESIROrthogonalization(KSP,int);
97 extern int KSPGMRESPrestartSet(KSP,int);
98 
99 
100 extern int KSPFGMRESSetRestart( KSP, int);
101 extern int KSPFGMRESSetPreAllocateVectors( KSP );
102 extern int KSPFGMRESSetOrthogonalization( KSP, int (*)( KSP, int ) );
103 extern int KSPFGMRESUnmodifiedGramSchmidtOrthogonalization( KSP, int );
104 extern int KSPFGMRESModifiedGramSchmidtOrthogonalization( KSP, int );
105 extern int KSPFGMRESIROrthogonalization( KSP, int );
106 extern int KSPFGMRESPrestartSet( KSP, int );
107 
108 extern int KSPFGMRESModifyPCNoChange( KSP, int, int, int, int, double );
109 extern int KSPFGMRESModifyPCGMRESVariableEx( KSP, int, int, int, int, double );
110 extern int KSPFGMRESModifyPCEx( KSP, int, int, int, int, double );
111 extern int KSPFGMRESSetModifyPC(KSP,int (*)(KSP,int,int,int,int,double));
112 
113 extern int KSPSetFromOptions(KSP);
114 extern int KSPSetTypeFromOptions(KSP);
115 extern int KSPAddOptionsChecker(int (*)(KSP));
116 
117 extern int KSPSingularValueMonitor(KSP,int,double, void * );
118 extern int KSPDefaultMonitor(KSP,int,double, void *);
119 extern int KSPTrueMonitor(KSP,int,double, void *);
120 extern int KSPDefaultSMonitor(KSP,int,double, void *);
121 extern int KSPVecViewMonitor(KSP,int,double,void *);
122 extern int KSPGMRESKrylovMonitor(KSP,int,double,void *);
123 
124 extern int KSPDefaultConverged(KSP,int,double, void *);
125 extern int KSPSkipConverged(KSP,int,double, void *);
126 
127 extern int KSPResidual(KSP,Vec,Vec,Vec,Vec,Vec,Vec);
128 extern int KSPUnwindPreconditioner(KSP,Vec,Vec);
129 extern int KSPDefaultBuildSolution(KSP,Vec,Vec*);
130 extern int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
131 
132 extern int KSPPrintHelp(KSP);
133 
134 extern int KSPSetOptionsPrefix(KSP,char*);
135 extern int KSPAppendOptionsPrefix(KSP,char*);
136 extern int KSPGetOptionsPrefix(KSP,char**);
137 
138 extern int KSPView(KSP,Viewer);
139 
140 extern int KSPComputeExplicitOperator(KSP,Mat *);
141 
142 typedef enum {KSP_CG_SYMMETRIC=1, KSP_CG_HERMITIAN=2} KSPCGType;
143 extern int KSPCGSetType(KSP,KSPCGType);
144 
145 extern int PCPreSolve(PC,KSP);
146 extern int PCPostSolve(PC,KSP);
147 
148 extern int KSPLGMonitorCreate(char*,char*,int,int,int,int,DrawLG*);
149 extern int KSPLGMonitor(KSP,int,double,void*);
150 extern int KSPLGMonitorDestroy(DrawLG);
151 extern int KSPLGTrueMonitorCreate(MPI_Comm,char*,char*,int,int,int,int,DrawLG*);
152 extern int KSPLGTrueMonitor(KSP,int,double,void*);
153 extern int KSPLGTrueMonitorDestroy(DrawLG);
154 
155 #endif
156 
157 
158