Actual source code: slepcsvd.h
1: /*
2: User interface for SLEPc's singular value solvers.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
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 <slepceps.h>
28: PETSC_EXTERN PetscErrorCode SVDInitializePackage(void);
30: /*S
31: SVD - Abstract SLEPc object that manages all the singular value
32: problem solvers.
34: Level: beginner
36: .seealso: SVDCreate()
37: S*/
38: typedef struct _p_SVD* SVD;
40: /*J
41: SVDType - String with the name of a SLEPc singular value solver
43: Level: beginner
45: .seealso: SVDSetType(), SVD
46: J*/
47: typedef const char* SVDType;
48: #define SVDCROSS "cross"
49: #define SVDCYCLIC "cyclic"
50: #define SVDLAPACK "lapack"
51: #define SVDLANCZOS "lanczos"
52: #define SVDTRLANCZOS "trlanczos"
54: /* Logging support */
55: PETSC_EXTERN PetscClassId SVD_CLASSID;
57: /*E
58: SVDTransposeMode - Determines how to handle the transpose of the matrix
60: Level: advanced
62: .seealso: SVDSetTransposeMode(), SVDGetTransposeMode()
63: E*/
64: typedef enum { SVD_TRANSPOSE_EXPLICIT,
65: SVD_TRANSPOSE_IMPLICIT } SVDTransposeMode;
67: /*E
68: SVDWhich - Determines whether largest or smallest singular triplets
69: are to be computed
71: Level: intermediate
73: .seealso: SVDSetWhichSingularTriplets(), SVDGetWhichSingularTriplets()
74: E*/
75: typedef enum { SVD_LARGEST,
76: SVD_SMALLEST } SVDWhich;
78: /*E
79: SVDConvergedReason - Reason a singular value solver was said to
80: have converged or diverged
82: Level: beginner
84: .seealso: SVDSolve(), SVDGetConvergedReason(), SVDSetTolerances()
85: E*/
86: typedef enum {/* converged */
87: SVD_CONVERGED_TOL = 2,
88: /* diverged */
89: SVD_DIVERGED_ITS = -3,
90: SVD_DIVERGED_BREAKDOWN = -4,
91: SVD_CONVERGED_ITERATING = 0 } SVDConvergedReason;
93: PETSC_EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
94: PETSC_EXTERN PetscErrorCode SVDSetIP(SVD,IP);
95: PETSC_EXTERN PetscErrorCode SVDGetIP(SVD,IP*);
96: PETSC_EXTERN PetscErrorCode SVDSetDS(SVD,DS);
97: PETSC_EXTERN PetscErrorCode SVDGetDS(SVD,DS*);
98: PETSC_EXTERN PetscErrorCode SVDSetType(SVD,SVDType);
99: PETSC_EXTERN PetscErrorCode SVDGetType(SVD,SVDType*);
100: PETSC_EXTERN PetscErrorCode SVDSetOperator(SVD,Mat);
101: PETSC_EXTERN PetscErrorCode SVDGetOperator(SVD,Mat*);
102: PETSC_EXTERN PetscErrorCode SVDSetInitialSpace(SVD,PetscInt,Vec*);
103: PETSC_EXTERN PetscErrorCode SVDSetInitialSpaceLeft(SVD,PetscInt,Vec*);
104: PETSC_EXTERN PetscErrorCode SVDSetTransposeMode(SVD,SVDTransposeMode);
105: PETSC_EXTERN PetscErrorCode SVDGetTransposeMode(SVD,SVDTransposeMode*);
106: PETSC_EXTERN PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
107: PETSC_EXTERN PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
108: PETSC_EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
109: PETSC_EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
110: PETSC_EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
111: PETSC_EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
112: PETSC_EXTERN PetscErrorCode SVDSetFromOptions(SVD);
113: PETSC_EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char*);
114: PETSC_EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char*);
115: PETSC_EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
116: PETSC_EXTERN PetscErrorCode SVDSetUp(SVD);
117: PETSC_EXTERN PetscErrorCode SVDSolve(SVD);
118: PETSC_EXTERN PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
119: PETSC_EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
120: PETSC_EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
121: PETSC_EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
122: PETSC_EXTERN PetscErrorCode SVDComputeResidualNorms(SVD,PetscInt,PetscReal*,PetscReal*);
123: PETSC_EXTERN PetscErrorCode SVDComputeRelativeError(SVD,PetscInt,PetscReal*);
124: PETSC_EXTERN PetscErrorCode SVDGetOperationCounters(SVD,PetscInt*,PetscInt*);
125: PETSC_EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
126: PETSC_EXTERN PetscErrorCode SVDPrintSolution(SVD,PetscViewer);
127: PETSC_EXTERN PetscErrorCode SVDDestroy(SVD*);
128: PETSC_EXTERN PetscErrorCode SVDReset(SVD);
130: PETSC_EXTERN PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt);
131: PETSC_EXTERN PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
132: PETSC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
133: PETSC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void **);
134: PETSC_EXTERN PetscErrorCode SVDMonitorAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
135: PETSC_EXTERN PetscErrorCode SVDMonitorFirst(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
136: PETSC_EXTERN PetscErrorCode SVDMonitorConverged(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
137: PETSC_EXTERN PetscErrorCode SVDMonitorLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
138: PETSC_EXTERN PetscErrorCode SVDMonitorLGAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
140: PETSC_EXTERN PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
141: PETSC_EXTERN PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);
143: PETSC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
144: PETSC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);
146: PETSC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
147: PETSC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
148: PETSC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
149: PETSC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);
151: PETSC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
152: PETSC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);
154: PETSC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
155: PETSC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
157: PETSC_EXTERN PetscFunctionList SVDList;
158: PETSC_EXTERN PetscBool SVDRegisterAllCalled;
159: PETSC_EXTERN PetscErrorCode SVDRegisterAll(void);
160: PETSC_EXTERN PetscErrorCode SVDRegister(const char[],PetscErrorCode(*)(SVD));
162: #endif