xref: /petsc/include/petscksp.h (revision 94b7f48cc472a54ea2ce57edf1fe19e8a254237c)
1 /* $Id: petscksp.h,v 1.107 2001/08/06 21:16:38 bsmith Exp $ */
2 /*
3    Defines the interface functions for the Krylov subspace accelerators.
4 */
5 #ifndef __PETSCKSP_H
6 #define __PETSCKSP_H
7 #include "petscpc.h"
8 PETSC_EXTERN_CXX_BEGIN
9 
10 EXTERN int KSPInitializePackage(const char[]);
11 
12 /*S
13      KSP - Abstract PETSc object that manages all Krylov methods
14 
15    Level: beginner
16 
17   Concepts: Krylov methods
18 
19 .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP
20 S*/
21 typedef struct _p_KSP*     KSP;
22 
23 /*E
24     KSPType - String with the name of a PETSc Krylov method or the creation function
25        with an optional dynamic library name, for example
26        http://www.mcs.anl.gov/petsc/lib.a:mykspcreate()
27 
28    Level: beginner
29 
30 .seealso: KSPSetType(), KSP
31 E*/
32 #define KSPRICHARDSON "richardson"
33 #define KSPCHEBYCHEV  "chebychev"
34 #define KSPCG         "cg"
35 #define KSPCGNE       "cgne"
36 #define KSPGMRES      "gmres"
37 #define KSPTCQMR      "tcqmr"
38 #define KSPBCGS       "bcgs"
39 #define KSPCGS        "cgs"
40 #define KSPTFQMR      "tfqmr"
41 #define KSPCR         "cr"
42 #define KSPLSQR       "lsqr"
43 #define KSPPREONLY    "preonly"
44 #define KSPQCG        "qcg"
45 #define KSPBICG       "bicg"
46 #define KSPFGMRES     "fgmres"
47 #define KSPMINRES     "minres"
48 #define KSPSYMMLQ     "symmlq"
49 #define KSPLGMRES     "lgmres"
50 typedef char * KSPType;
51 
52 /* Logging support */
53 extern int KSP_COOKIE;
54 extern int KSP_GMRESOrthogonalization;
55 extern int KSP_SetUp, KSP_Solve;
56 
57 EXTERN int KSPCreate(MPI_Comm,KSP *);
58 EXTERN int KSPSetType(KSP,KSPType);
59 EXTERN int KSPSetUp(KSP);
60 EXTERN int KSPSetUpOnBlocks(KSP);
61 EXTERN int KSPSolve(KSP);
62 EXTERN int KSPSolveTranspose(KSP);
63 EXTERN int KSPDestroy(KSP);
64 
65 extern PetscFList KSPList;
66 EXTERN int KSPRegisterAll(const char[]);
67 EXTERN int KSPRegisterDestroy(void);
68 
69 EXTERN int KSPRegister(const char[],const char[],const char[],int(*)(KSP));
70 
71 /*MC
72    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
73 
74    Synopsis:
75    int KSPRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(KSP))
76 
77    Not Collective
78 
79    Input Parameters:
80 +  name_solver - name of a new user-defined solver
81 .  path - path (either absolute or relative) the library containing this solver
82 .  name_create - name of routine to create method context
83 -  routine_create - routine to create method context
84 
85    Notes:
86    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
87 
88    If dynamic libraries are used, then the fourth input argument (routine_create)
89    is ignored.
90 
91    Sample usage:
92 .vb
93    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
94                "MySolverCreate",MySolverCreate);
95 .ve
96 
97    Then, your solver can be chosen with the procedural interface via
98 $     KSPSetType(ksp,"my_solver")
99    or at runtime via the option
100 $     -ksp_type my_solver
101 
102    Level: advanced
103 
104    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT},
105           and others of the form ${any_environmental_variable} occuring in pathname will be
106           replaced with appropriate values.
107          If your function is not being put into a shared library then use KSPRegister() instead
108 
109 .keywords: KSP, register
110 
111 .seealso: KSPRegisterAll(), KSPRegisterDestroy()
112 
113 M*/
114 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
115 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
116 #else
117 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
118 #endif
119 
120 EXTERN int KSPGetType(KSP,KSPType *);
121 EXTERN int KSPSetPreconditionerSide(KSP,PCSide);
122 EXTERN int KSPGetPreconditionerSide(KSP,PCSide*);
123 EXTERN int KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,int*);
124 EXTERN int KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,int);
125 EXTERN int KSPSetInitialGuessNonzero(KSP,PetscTruth);
126 EXTERN int KSPGetInitialGuessNonzero(KSP,PetscTruth *);
127 EXTERN int KSPSetInitialGuessKnoll(KSP,PetscTruth);
128 EXTERN int KSPGetInitialGuessKnoll(KSP,PetscTruth*);
129 EXTERN int KSPSetComputeEigenvalues(KSP,PetscTruth);
130 EXTERN int KSPSetComputeSingularValues(KSP,PetscTruth);
131 EXTERN int KSPSetRhs(KSP,Vec);
132 EXTERN int KSPGetRhs(KSP,Vec *);
133 EXTERN int KSPSetSolution(KSP,Vec);
134 EXTERN int KSPGetSolution(KSP,Vec *);
135 EXTERN int KSPGetResidualNorm(KSP,PetscReal*);
136 EXTERN int KSPGetIterationNumber(KSP,int*);
137 
138 EXTERN int KSPSetPC(KSP,PC);
139 EXTERN int KSPGetPC(KSP,PC*);
140 
141 EXTERN int KSPSetMonitor(KSP,int (*)(KSP,int,PetscReal,void*),void *,int (*)(void*));
142 EXTERN int KSPClearMonitor(KSP);
143 EXTERN int KSPGetMonitorContext(KSP,void **);
144 EXTERN int KSPGetResidualHistory(KSP,PetscReal*[],int *);
145 EXTERN int KSPSetResidualHistory(KSP,PetscReal[],int,PetscTruth);
146 
147 
148 EXTERN int KSPBuildSolution(KSP,Vec,Vec *);
149 EXTERN int KSPBuildResidual(KSP,Vec,Vec,Vec *);
150 
151 EXTERN int KSPRichardsonSetScale(KSP,PetscReal);
152 EXTERN int KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
153 EXTERN int KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
154 EXTERN int KSPComputeEigenvalues(KSP,int,PetscReal*,PetscReal*,int *);
155 EXTERN int KSPComputeEigenvaluesExplicitly(KSP,int,PetscReal*,PetscReal*);
156 
157 #define KSPGMRESSetRestart(ksp,r) PetscTryMethod((ksp),KSPGMRESSetRestart_C,(KSP,int),((ksp),(r)))
158 #define KSPGMRESSetHapTol(ksp,tol) PetscTryMethod((ksp),KSPGMRESSetHapTol_C,(KSP,PetscReal),((ksp),(tol)))
159 
160 EXTERN int KSPGMRESSetPreAllocateVectors(KSP);
161 EXTERN int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int));
162 EXTERN int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int);
163 EXTERN int KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,int);
164 /*E
165     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
166 
167    Level: advanced
168 
169 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
170           KSPGMRESSetCGSRefinementType()
171 
172 E*/
173 typedef enum {KSP_GMRES_CGS_REFINEMENT_NONE, KSP_GMRES_CGS_REFINEMENT_IFNEEDED, KSP_GMRES_CGS_REFINEMENT_ALWAYS} KSPGMRESCGSRefinementType;
174 EXTERN int KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
175 
176 EXTERN int KSPFGMRESModifyPCNoChange(KSP,int,int,PetscReal,void*);
177 EXTERN int KSPFGMRESModifyPCKSP(KSP,int,int,PetscReal,void*);
178 EXTERN int KSPFGMRESSetModifyPC(KSP,int (*)(KSP,int,int,PetscReal,void*),void*,int(*)(void*));
179 
180 EXTERN int KSPQCGSetTrustRegionRadius(KSP,PetscReal);
181 EXTERN int KSPQCGGetQuadratic(KSP,PetscReal*);
182 EXTERN int KSPQCGGetTrialStepNorm(KSP,PetscReal*);
183 
184 EXTERN int KSPSetFromOptions(KSP);
185 EXTERN int KSPAddOptionsChecker(int (*)(KSP));
186 
187 EXTERN int KSPSingularValueMonitor(KSP,int,PetscReal,void *);
188 EXTERN int KSPDefaultMonitor(KSP,int,PetscReal,void *);
189 EXTERN int KSPTrueMonitor(KSP,int,PetscReal,void *);
190 EXTERN int KSPDefaultSMonitor(KSP,int,PetscReal,void *);
191 EXTERN int KSPVecViewMonitor(KSP,int,PetscReal,void *);
192 EXTERN int KSPGMRESKrylovMonitor(KSP,int,PetscReal,void *);
193 
194 EXTERN int KSPUnwindPreconditioner(KSP,Vec,Vec);
195 EXTERN int KSPDefaultBuildSolution(KSP,Vec,Vec*);
196 EXTERN int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
197 
198 EXTERN int KSPSetOperators(KSP,Mat,Mat,MatStructure);
199 EXTERN int KSPSetOptionsPrefix(KSP,const char[]);
200 EXTERN int KSPAppendOptionsPrefix(KSP,const char[]);
201 EXTERN int KSPGetOptionsPrefix(KSP,char*[]);
202 
203 EXTERN int KSPSetDiagonalScale(KSP,PetscTruth);
204 EXTERN int KSPGetDiagonalScale(KSP,PetscTruth*);
205 EXTERN int KSPSetDiagonalScaleFix(KSP,PetscTruth);
206 EXTERN int KSPGetDiagonalScaleFix(KSP,PetscTruth*);
207 
208 EXTERN int KSPView(KSP,PetscViewer);
209 
210 /*E
211     KSPNormType - Norm that is passed in the Krylov convergence
212        test routines.
213 
214    Level: advanced
215 
216    Notes: this must match finclude/petscksp.h
217 
218 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
219           KSPSetConvergenceTest()
220 E*/
221 typedef enum {KSP_NO_NORM               = 0,
222               KSP_PRECONDITIONED_NORM   = 1,
223               KSP_UNPRECONDITIONED_NORM = 2,
224               KSP_NATURAL_NORM          = 3} KSPNormType;
225 EXTERN int KSPSetNormType(KSP,KSPNormType);
226 
227 /*E
228     KSPConvergedReason - reason a Krylov method was said to
229          have converged or diverged
230 
231    Level: beginner
232 
233    Notes: this must match finclude/petscksp.h
234 
235    Developer note: The string versions of these are in
236      src/sles/ksp/interface/itfunc.c called convergedreasons.
237      If these enums are changed you much change those.
238 
239 .seealso: SLESSolve(), KSPSolve(), KSPGetConvergedReason()
240 E*/
241 typedef enum {/* converged */
242               KSP_CONVERGED_RTOL               =  2,
243               KSP_CONVERGED_ATOL               =  3,
244               KSP_CONVERGED_ITS                =  4,
245               KSP_CONVERGED_QCG_NEG_CURVE      =  5,
246               KSP_CONVERGED_QCG_CONSTRAINED    =  6,
247               KSP_CONVERGED_STEP_LENGTH        =  7,
248               /* diverged */
249               KSP_DIVERGED_ITS                 = -3,
250               KSP_DIVERGED_DTOL                = -4,
251               KSP_DIVERGED_BREAKDOWN           = -5,
252               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
253               KSP_DIVERGED_NONSYMMETRIC        = -7,
254               KSP_DIVERGED_INDEFINITE_PC       = -8,
255 
256               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
257 
258 EXTERN int KSPSetConvergenceTest(KSP,int (*)(KSP,int,PetscReal,KSPConvergedReason*,void*),void *);
259 EXTERN int KSPGetConvergenceContext(KSP,void **);
260 EXTERN int KSPDefaultConverged(KSP,int,PetscReal,KSPConvergedReason*,void *);
261 EXTERN int KSPSkipConverged(KSP,int,PetscReal,KSPConvergedReason*,void *);
262 EXTERN int KSPGetConvergedReason(KSP,KSPConvergedReason *);
263 
264 EXTERN int KSPComputeExplicitOperator(KSP,Mat *);
265 
266 /*E
267     KSPCGType - Determines what type of CG to use
268 
269    Level: beginner
270 
271 .seealso: KSPCGSetType()
272 E*/
273 typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType;
274 
275 EXTERN int KSPCGSetType(KSP,KSPCGType);
276 
277 EXTERN int PCPreSolve(PC,KSP);
278 EXTERN int PCPostSolve(PC,KSP);
279 
280 EXTERN int KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
281 EXTERN int KSPLGMonitor(KSP,int,PetscReal,void*);
282 EXTERN int KSPLGMonitorDestroy(PetscDrawLG);
283 EXTERN int KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
284 EXTERN int KSPLGTrueMonitor(KSP,int,PetscReal,void*);
285 EXTERN int KSPLGTrueMonitorDestroy(PetscDrawLG);
286 
287 PETSC_EXTERN_CXX_END
288 #endif
289