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