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