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