Actual source code: slepcsvd.h
1: /*
2: User interface for the SLEPC singular value solvers.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
26: #include slepcsys.h
27: #include slepceps.h
28: PETSC_EXTERN_CXX_BEGIN
30: extern PetscErrorCode SVDInitializePackage(const char[]);
32: /*S
33: SVD - Abstract SLEPc object that manages all the singular value
34: problem solvers.
36: Level: beginner
38: .seealso: SVDCreate()
39: S*/
40: typedef struct _p_SVD* SVD;
42: /*E
43: SVDType - String with the name of a SLEPc singular value solver
45: Level: beginner
47: .seealso: SVDSetType(), SVD
48: E*/
49: #define SVDType char*
50: #define SVDCROSS "cross"
51: #define SVDCYCLIC "cyclic"
52: #define SVDLAPACK "lapack"
53: #define SVDLANCZOS "lanczos"
54: #define SVDTRLANCZOS "trlanczos"
56: /* Logging support */
57: extern PetscClassId SVD_CLASSID;
59: /*E
60: SVDTransposeMode - determines how to handle the transpose of the matrix
62: Level: advanced
64: .seealso: SVDSetTransposeMode(), SVDGetTransposeMode()
65: E*/
66: typedef enum { SVD_TRANSPOSE_EXPLICIT,
67: SVD_TRANSPOSE_IMPLICIT } SVDTransposeMode;
69: /*E
70: SVDWhich - determines whether largest or smallest singular triplets
71: are to be computed
73: Level: intermediate
75: .seealso: SVDSetWhichSingularTriplets(), SVDGetWhichSingularTriplets()
76: E*/
77: typedef enum { SVD_LARGEST,
78: SVD_SMALLEST } SVDWhich;
80: /*E
81: SVDConvergedReason - reason a singular value solver was said to
82: have converged or diverged
84: Level: beginner
86: .seealso: SVDSolve(), SVDGetConvergedReason(), SVDSetTolerances()
87: E*/
88: typedef enum {/* converged */
89: SVD_CONVERGED_TOL = 2,
90: /* diverged */
91: SVD_DIVERGED_ITS = -3,
92: SVD_DIVERGED_BREAKDOWN = -4,
93: SVD_CONVERGED_ITERATING = 0 } SVDConvergedReason;
95: extern PetscErrorCode SVDCreate(MPI_Comm,SVD*);
96: extern PetscErrorCode SVDSetIP(SVD,IP);
97: extern PetscErrorCode SVDGetIP(SVD,IP*);
98: extern PetscErrorCode SVDSetType(SVD,const SVDType);
99: extern PetscErrorCode SVDGetType(SVD,const SVDType*);
100: extern PetscErrorCode SVDSetOperator(SVD,Mat);
101: extern PetscErrorCode SVDGetOperator(SVD,Mat*);
102: extern PetscErrorCode SVDSetInitialSpace(SVD,PetscInt,Vec*);
103: extern PetscErrorCode SVDSetTransposeMode(SVD,SVDTransposeMode);
104: extern PetscErrorCode SVDGetTransposeMode(SVD,SVDTransposeMode*);
105: extern PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
106: extern PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
107: extern PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
108: extern PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
109: extern PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
110: extern PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
111: extern PetscErrorCode SVDSetFromOptions(SVD);
112: extern PetscErrorCode SVDSetOptionsPrefix(SVD,const char*);
113: extern PetscErrorCode SVDAppendOptionsPrefix(SVD,const char*);
114: extern PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
115: extern PetscErrorCode SVDSetUp(SVD);
116: extern PetscErrorCode SVDSolve(SVD);
117: extern PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
118: extern PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
119: extern PetscErrorCode SVDGetConverged(SVD,PetscInt*);
120: extern PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
121: extern PetscErrorCode SVDComputeResidualNorms(SVD,PetscInt,PetscReal*,PetscReal*);
122: extern PetscErrorCode SVDComputeRelativeError(SVD,PetscInt,PetscReal*);
123: extern PetscErrorCode SVDGetOperationCounters(SVD,PetscInt*,PetscInt*);
124: extern PetscErrorCode SVDView(SVD,PetscViewer);
125: extern PetscErrorCode SVDPrintSolution(SVD,PetscViewer);
126: extern PetscErrorCode SVDDestroy(SVD*);
127: extern PetscErrorCode SVDReset(SVD);
129: extern PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
130: extern PetscErrorCode SVDMonitorCancel(SVD);
131: extern PetscErrorCode SVDGetMonitorContext(SVD,void **);
132: extern PetscErrorCode SVDMonitorAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
133: extern PetscErrorCode SVDMonitorFirst(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
134: extern PetscErrorCode SVDMonitorConverged(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
135: extern PetscErrorCode SVDMonitorLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
136: extern PetscErrorCode SVDMonitorLGAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
138: extern PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
139: extern PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);
141: extern PetscErrorCode SVDDense(PetscInt,PetscInt,PetscScalar*,PetscReal*,PetscScalar*,PetscScalar*);
143: extern PetscErrorCode SVDCrossSetEPS(SVD,EPS);
144: extern PetscErrorCode SVDCrossGetEPS(SVD,EPS*);
146: extern PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
147: extern PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
148: extern PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
149: extern PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);
151: extern PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
152: extern PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);
154: extern PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
155: extern PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
157: extern PetscFList SVDList;
158: extern PetscBool SVDRegisterAllCalled;
159: extern PetscErrorCode SVDRegisterAll(const char[]);
160: extern PetscErrorCode SVDRegisterDestroy(void);
161: extern PetscErrorCode SVDRegister(const char[],const char[],const char[],PetscErrorCode(*)(SVD));
163: /*MC
164: SVDRegisterDynamic - Adds a method to the singular value solver package.
166: Synopsis:
167: PetscErrorCode SVDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(SVD))
169: Not Collective
171: Input Parameters:
172: + name_solver - name of a new user-defined solver
173: . path - path (either absolute or relative) the library containing this solver
174: . name_create - name of routine to create the solver context
175: - routine_create - routine to create the solver context
177: Notes:
178: SVDRegisterDynamic() may be called multiple times to add several user-defined solvers.
180: If dynamic libraries are used, then the fourth input argument (routine_create)
181: is ignored.
183: Sample usage:
184: .vb
185: SVDRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
186: "MySolverCreate",MySolverCreate);
187: .ve
189: Then, your solver can be chosen with the procedural interface via
190: $ SVDSetType(svd,"my_solver")
191: or at runtime via the option
192: $ -svd_type my_solver
194: Level: advanced
196: .seealso: SVDRegisterDestroy(), SVDRegisterAll()
197: M*/
198: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
199: #define SVDRegisterDynamic(a,b,c,d) SVDRegister(a,b,c,0)
200: #else
201: #define SVDRegisterDynamic(a,b,c,d) SVDRegister(a,b,c,d)
202: #endif
204: PETSC_EXTERN_CXX_END
205: #endif