xref: /petsc/include/petscksp.h (revision 15e8a5b3e3ad67526f50ce4d2e8f9da59aa8d948)
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 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
65 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
66 #else
67 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
68 #endif
69 
70 EXTERN int KSPGetType(KSP,KSPType *);
71 EXTERN int KSPSetPreconditionerSide(KSP,PCSide);
72 EXTERN int KSPGetPreconditionerSide(KSP,PCSide*);
73 EXTERN int KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,int*);
74 EXTERN int KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,int);
75 EXTERN int KSPSetInitialGuessNonzero(KSP,PetscTruth);
76 EXTERN int KSPGetInitialGuessNonzero(KSP,PetscTruth *);
77 EXTERN int KSPSetInitialGuessKnoll(KSP,PetscTruth);
78 EXTERN int KSPGetInitialGuessKnoll(KSP,PetscTruth*);
79 EXTERN int KSPSetComputeEigenvalues(KSP,PetscTruth);
80 EXTERN int KSPSetComputeSingularValues(KSP,PetscTruth);
81 EXTERN int KSPSetRhs(KSP,Vec);
82 EXTERN int KSPGetRhs(KSP,Vec *);
83 EXTERN int KSPSetSolution(KSP,Vec);
84 EXTERN int KSPGetSolution(KSP,Vec *);
85 EXTERN int KSPGetResidualNorm(KSP,PetscReal*);
86 EXTERN int KSPGetIterationNumber(KSP,int*);
87 
88 EXTERN int KSPSetPC(KSP,PC);
89 EXTERN int KSPGetPC(KSP,PC*);
90 
91 EXTERN int KSPSetMonitor(KSP,int (*)(KSP,int,PetscReal,void*),void *,int (*)(void*));
92 EXTERN int KSPClearMonitor(KSP);
93 EXTERN int KSPGetMonitorContext(KSP,void **);
94 EXTERN int KSPGetResidualHistory(KSP,PetscReal **,int *);
95 EXTERN int KSPSetResidualHistory(KSP,PetscReal *,int,PetscTruth);
96 
97 
98 EXTERN int KSPBuildSolution(KSP,Vec,Vec *);
99 EXTERN int KSPBuildResidual(KSP,Vec,Vec,Vec *);
100 
101 EXTERN int KSPRichardsonSetScale(KSP,PetscReal);
102 EXTERN int KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
103 EXTERN int KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
104 EXTERN int KSPComputeEigenvalues(KSP,int,PetscReal*,PetscReal*,int *);
105 EXTERN int KSPComputeEigenvaluesExplicitly(KSP,int,PetscReal*,PetscReal*);
106 
107 #define KSPGMRESSetRestart(ksp,r) PetscTryMethod((ksp),KSPGMRESSetRestart_C,(KSP,int),((ksp),(r)))
108 #define KSPGMRESSetHapTol(ksp,tol) PetscTryMethod((ksp),KSPGMRESSetHapTol_C,(KSP,PetscReal),((ksp),(tol)))
109 
110 EXTERN int KSPGMRESSetPreAllocateVectors(KSP);
111 EXTERN int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int));
112 EXTERN int KSPGMRESUnmodifiedGramSchmidtOrthogonalization(KSP,int);
113 EXTERN int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int);
114 EXTERN int KSPGMRESIROrthogonalization(KSP,int);
115 
116 EXTERN int KSPFGMRESModifyPCNoChange(KSP,int,int,PetscReal,void*);
117 EXTERN int KSPFGMRESModifyPCSLES(KSP,int,int,PetscReal,void*);
118 EXTERN int KSPFGMRESSetModifyPC(KSP,int (*)(KSP,int,int,PetscReal,void*),void*,int(*)(void*));
119 
120 EXTERN int KSPQCGSetTrustRegionRadius(KSP,PetscReal);
121 EXTERN int KSPQCGGetQuadratic(KSP,PetscReal*);
122 EXTERN int KSPQCGGetTrialStepNorm(KSP,PetscReal*);
123 
124 EXTERN int KSPSetFromOptions(KSP);
125 EXTERN int KSPAddOptionsChecker(int (*)(KSP));
126 
127 EXTERN int KSPSingularValueMonitor(KSP,int,PetscReal,void *);
128 EXTERN int KSPDefaultMonitor(KSP,int,PetscReal,void *);
129 EXTERN int KSPTrueMonitor(KSP,int,PetscReal,void *);
130 EXTERN int KSPDefaultSMonitor(KSP,int,PetscReal,void *);
131 EXTERN int KSPVecViewMonitor(KSP,int,PetscReal,void *);
132 EXTERN int KSPGMRESKrylovMonitor(KSP,int,PetscReal,void *);
133 
134 EXTERN int KSPUnwindPreconditioner(KSP,Vec,Vec);
135 EXTERN int KSPDefaultBuildSolution(KSP,Vec,Vec*);
136 EXTERN int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
137 
138 EXTERN int KSPSetOptionsPrefix(KSP,char*);
139 EXTERN int KSPAppendOptionsPrefix(KSP,char*);
140 EXTERN int KSPGetOptionsPrefix(KSP,char**);
141 
142 EXTERN int KSPView(KSP,PetscViewer);
143 
144 /*E
145     KSPNormType - Norm that is passed in the Krylov convergence
146        test routines.
147 
148    Level: advanced
149 
150    Notes: this must match finclude/petscksp.h
151 
152 .seealso: SLESSolve(), KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
153           KSPSetConvergenceTest()
154 E*/
155 typedef enum {KSP_NO_NORM               = 0,
156               KSP_PRECONDITIONED_NORM   = 1,
157               KSP_UNPRECONDITIONED_NORM = 2,
158               KSP_NATURAL_NORM          = 3} KSPNormType;
159 EXTERN int KSPSetNormType(KSP,KSPNormType);
160 
161 /*E
162     KSPConvergedReason - reason a Krylov method was said to
163          have converged or diverged
164 
165    Level: beginner
166 
167    Notes: this must match finclude/petscksp.h
168 
169 .seealso: SLESSolve(), KSPSolve(), KSPGetConvergedReason()
170 E*/
171 typedef enum {/* converged */
172               KSP_CONVERGED_RTOL               =  2,
173               KSP_CONVERGED_ATOL               =  3,
174               KSP_CONVERGED_ITS                =  4,
175               KSP_CONVERGED_QCG_NEG_CURVE      =  5,
176               KSP_CONVERGED_QCG_CONSTRAINED    =  6,
177               KSP_CONVERGED_STEP_LENGTH        =  7,
178               /* diverged */
179               KSP_DIVERGED_ITS                 = -3,
180               KSP_DIVERGED_DTOL                = -4,
181               KSP_DIVERGED_BREAKDOWN           = -5,
182               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
183               KSP_DIVERGED_NONSYMMETRIC        = -7,
184               KSP_DIVERGED_INDEFINITE_PC       = -8,
185 
186               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
187 
188 EXTERN int KSPSetConvergenceTest(KSP,int (*)(KSP,int,PetscReal,KSPConvergedReason*,void*),void *);
189 EXTERN int KSPGetConvergenceContext(KSP,void **);
190 EXTERN int KSPDefaultConverged(KSP,int,PetscReal,KSPConvergedReason*,void *);
191 EXTERN int KSPSkipConverged(KSP,int,PetscReal,KSPConvergedReason*,void *);
192 EXTERN int KSPGetConvergedReason(KSP,KSPConvergedReason *);
193 
194 EXTERN int KSPComputeExplicitOperator(KSP,Mat *);
195 
196 /*E
197     KSPCGType - Determines what type of CG to use
198 
199    Level: beginner
200 
201 .seealso: KSPCGSetType()
202 E*/
203 typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType;
204 
205 EXTERN int KSPCGSetType(KSP,KSPCGType);
206 
207 EXTERN int PCPreSolve(PC,KSP);
208 EXTERN int PCPostSolve(PC,KSP);
209 
210 EXTERN int KSPLGMonitorCreate(char*,char*,int,int,int,int,PetscDrawLG*);
211 EXTERN int KSPLGMonitor(KSP,int,PetscReal,void*);
212 EXTERN int KSPLGMonitorDestroy(PetscDrawLG);
213 EXTERN int KSPLGTrueMonitorCreate(MPI_Comm,char*,char*,int,int,int,int,PetscDrawLG*);
214 EXTERN int KSPLGTrueMonitor(KSP,int,PetscReal,void*);
215 EXTERN int KSPLGTrueMonitorDestroy(PetscDrawLG);
216 
217 #endif
218 
219 
220