Actual source code: petscksp.h
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
7: #include petscpc.h
8: PETSC_EXTERN_CXX_BEGIN
10: EXTERN int KSPInitializePackage(const char[]);
12: /*S
13: KSP - Abstract PETSc object that manages all Krylov methods
15: Level: beginner
17: Concepts: Krylov methods
19: .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP
20: S*/
21: typedef struct _p_KSP* KSP;
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()
28: Level: beginner
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: #define KSPType char*
52: /* Logging support */
53: extern int KSP_COOKIE;
54: extern int KSP_GMRESOrthogonalization;
55: extern int KSP_SetUp, KSP_Solve;
57: EXTERN int KSPCreate(MPI_Comm,KSP *);
58: EXTERN int KSPSetType(KSP,const 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);
65: extern PetscFList KSPList;
66: EXTERN int KSPRegisterAll(const char[]);
67: EXTERN int KSPRegisterDestroy(void);
69: EXTERN int KSPRegister(const char[],const char[],const char[],int(*)(KSP));
71: /*MC
72: KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
74: Synopsis:
75: int KSPRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(KSP))
77: Not Collective
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
85: Notes:
86: KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
88: If dynamic libraries are used, then the fourth input argument (routine_create)
89: is ignored.
91: Sample usage:
92: .vb
93: KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
94: "MySolverCreate",MySolverCreate);
95: .ve
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
102: Level: advanced
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
109: .keywords: KSP, register
111: .seealso: KSPRegisterAll(), KSPRegisterDestroy()
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
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*);
138: EXTERN int KSPSetPC(KSP,PC);
139: EXTERN int KSPGetPC(KSP,PC*);
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);
147: /* not sure where to put this */
148: EXTERN int PCKSPGetKSP(PC,KSP*);
149: EXTERN int PCBJacobiGetSubKSP(PC,int*,int*,KSP*[]);
150: EXTERN int PCASMGetSubKSP(PC,int*,int*,KSP*[]);
152: EXTERN int KSPBuildSolution(KSP,Vec,Vec *);
153: EXTERN int KSPBuildResidual(KSP,Vec,Vec,Vec *);
155: EXTERN int KSPRichardsonSetScale(KSP,PetscReal);
156: EXTERN int KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
157: EXTERN int KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
158: EXTERN int KSPComputeEigenvalues(KSP,int,PetscReal*,PetscReal*,int *);
159: EXTERN int KSPComputeEigenvaluesExplicitly(KSP,int,PetscReal*,PetscReal*);
161: EXTERN int KSPGMRESSetRestart(KSP, int);
162: EXTERN int KSPGMRESSetHapTol(KSP,PetscReal);
164: EXTERN int KSPGMRESSetPreAllocateVectors(KSP);
165: EXTERN int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int));
166: EXTERN int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int);
167: EXTERN int KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,int);
169: EXTERN int KSPLGMRESSetAugDim(KSP,int);
170: EXTERN int KSPLGMRESSetConstant(KSP);
172: /*E
173: KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
175: Level: advanced
177: .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
178: KSPGMRESSetCGSRefinementType()
180: E*/
181: typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
182: EXTERN int KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
184: EXTERN int KSPFGMRESModifyPCNoChange(KSP,int,int,PetscReal,void*);
185: EXTERN int KSPFGMRESModifyPCKSP(KSP,int,int,PetscReal,void*);
186: EXTERN int KSPFGMRESSetModifyPC(KSP,int (*)(KSP,int,int,PetscReal,void*),void*,int(*)(void*));
188: EXTERN int KSPQCGSetTrustRegionRadius(KSP,PetscReal);
189: EXTERN int KSPQCGGetQuadratic(KSP,PetscReal*);
190: EXTERN int KSPQCGGetTrialStepNorm(KSP,PetscReal*);
192: EXTERN int KSPSetFromOptions(KSP);
193: EXTERN int KSPAddOptionsChecker(int (*)(KSP));
195: EXTERN int KSPSingularValueMonitor(KSP,int,PetscReal,void *);
196: EXTERN int KSPDefaultMonitor(KSP,int,PetscReal,void *);
197: EXTERN int KSPTrueMonitor(KSP,int,PetscReal,void *);
198: EXTERN int KSPDefaultSMonitor(KSP,int,PetscReal,void *);
199: EXTERN int KSPVecViewMonitor(KSP,int,PetscReal,void *);
200: EXTERN int KSPGMRESKrylovMonitor(KSP,int,PetscReal,void *);
202: EXTERN int KSPUnwindPreconditioner(KSP,Vec,Vec);
203: EXTERN int KSPDefaultBuildSolution(KSP,Vec,Vec*);
204: EXTERN int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
206: EXTERN int KSPSetOperators(KSP,Mat,Mat,MatStructure);
207: EXTERN int KSPSetOptionsPrefix(KSP,const char[]);
208: EXTERN int KSPAppendOptionsPrefix(KSP,const char[]);
209: EXTERN int KSPGetOptionsPrefix(KSP,char*[]);
211: EXTERN int KSPSetDiagonalScale(KSP,PetscTruth);
212: EXTERN int KSPGetDiagonalScale(KSP,PetscTruth*);
213: EXTERN int KSPSetDiagonalScaleFix(KSP,PetscTruth);
214: EXTERN int KSPGetDiagonalScaleFix(KSP,PetscTruth*);
216: EXTERN int KSPView(KSP,PetscViewer);
218: /*E
219: KSPNormType - Norm that is passed in the Krylov convergence
220: test routines.
222: Level: advanced
224: Notes: this must match finclude/petscksp.h
226: .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
227: KSPSetConvergenceTest()
228: E*/
229: typedef enum {KSP_NO_NORM = 0,
230: KSP_PRECONDITIONED_NORM = 1,
231: KSP_UNPRECONDITIONED_NORM = 2,
232: KSP_NATURAL_NORM = 3} KSPNormType;
233: EXTERN int KSPSetNormType(KSP,KSPNormType);
235: /*E
236: KSPConvergedReason - reason a Krylov method was said to
237: have converged or diverged
239: Level: beginner
241: Notes: this must match finclude/petscksp.h
243: Developer note: The string versions of these are in
244: src/ksp/ksp/interface/itfunc.c called convergedreasons.
245: If these enums are changed you much change those.
247: .seealso: KSPSolve(), KSPGetConvergedReason()
248: E*/
249: typedef enum {/* converged */
250: KSP_CONVERGED_RTOL = 2,
251: KSP_CONVERGED_ATOL = 3,
252: KSP_CONVERGED_ITS = 4,
253: KSP_CONVERGED_QCG_NEG_CURVE = 5,
254: KSP_CONVERGED_QCG_CONSTRAINED = 6,
255: KSP_CONVERGED_STEP_LENGTH = 7,
256: /* diverged */
257: KSP_DIVERGED_ITS = -3,
258: KSP_DIVERGED_DTOL = -4,
259: KSP_DIVERGED_BREAKDOWN = -5,
260: KSP_DIVERGED_BREAKDOWN_BICG = -6,
261: KSP_DIVERGED_NONSYMMETRIC = -7,
262: KSP_DIVERGED_INDEFINITE_PC = -8,
263:
264: KSP_CONVERGED_ITERATING = 0} KSPConvergedReason;
266: EXTERN int KSPSetConvergenceTest(KSP,int (*)(KSP,int,PetscReal,KSPConvergedReason*,void*),void *);
267: EXTERN int KSPGetConvergenceContext(KSP,void **);
268: EXTERN int KSPDefaultConverged(KSP,int,PetscReal,KSPConvergedReason*,void *);
269: EXTERN int KSPSkipConverged(KSP,int,PetscReal,KSPConvergedReason*,void *);
270: EXTERN int KSPGetConvergedReason(KSP,KSPConvergedReason *);
272: EXTERN int KSPComputeExplicitOperator(KSP,Mat *);
274: /*E
275: KSPCGType - Determines what type of CG to use
277: Level: beginner
279: .seealso: KSPCGSetType()
280: E*/
281: typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType;
283: EXTERN int KSPCGSetType(KSP,KSPCGType);
285: EXTERN int PCPreSolve(PC,KSP);
286: EXTERN int PCPostSolve(PC,KSP);
288: EXTERN int KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
289: EXTERN int KSPLGMonitor(KSP,int,PetscReal,void*);
290: EXTERN int KSPLGMonitorDestroy(PetscDrawLG);
291: EXTERN int KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
292: EXTERN int KSPLGTrueMonitor(KSP,int,PetscReal,void*);
293: EXTERN int KSPLGTrueMonitorDestroy(PetscDrawLG);
295: PETSC_EXTERN_CXX_END
296: #endif