Actual source code: slepcsvd.h
slepc-3.5.2 2014-10-10
1: /*
2: User interface for SLEPc's singular value solvers.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2014, 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>
27: #include <slepcbv.h>
28: #include <slepcds.h>
30: PETSC_EXTERN PetscErrorCode SVDInitializePackage(void);
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: /*J
43: SVDType - String with the name of a SLEPc singular value solver
45: Level: beginner
47: .seealso: SVDSetType(), SVD
48: J*/
49: typedef const char* SVDType;
50: #define SVDCROSS "cross"
51: #define SVDCYCLIC "cyclic"
52: #define SVDLAPACK "lapack"
53: #define SVDLANCZOS "lanczos"
54: #define SVDTRLANCZOS "trlanczos"
56: /* Logging support */
57: PETSC_EXTERN PetscClassId SVD_CLASSID;
59: /*E
60: SVDWhich - Determines whether largest or smallest singular triplets
61: are to be computed
63: Level: intermediate
65: .seealso: SVDSetWhichSingularTriplets(), SVDGetWhichSingularTriplets()
66: E*/
67: typedef enum { SVD_LARGEST,
68: SVD_SMALLEST } SVDWhich;
70: /*E
71: SVDConvergedReason - Reason a singular value solver was said to
72: have converged or diverged
74: Level: beginner
76: .seealso: SVDSolve(), SVDGetConvergedReason(), SVDSetTolerances()
77: E*/
78: typedef enum {/* converged */
79: SVD_CONVERGED_TOL = 2,
80: /* diverged */
81: SVD_DIVERGED_ITS = -3,
82: SVD_DIVERGED_BREAKDOWN = -4,
83: SVD_CONVERGED_ITERATING = 0 } SVDConvergedReason;
85: PETSC_EXTERN PetscErrorCode SVDCreate(MPI_Comm,SVD*);
86: PETSC_EXTERN PetscErrorCode SVDSetBV(SVD,BV,BV);
87: PETSC_EXTERN PetscErrorCode SVDGetBV(SVD,BV*,BV*);
88: PETSC_EXTERN PetscErrorCode SVDSetDS(SVD,DS);
89: PETSC_EXTERN PetscErrorCode SVDGetDS(SVD,DS*);
90: PETSC_EXTERN PetscErrorCode SVDSetType(SVD,SVDType);
91: PETSC_EXTERN PetscErrorCode SVDGetType(SVD,SVDType*);
92: PETSC_EXTERN PetscErrorCode SVDSetOperator(SVD,Mat);
93: PETSC_EXTERN PetscErrorCode SVDGetOperator(SVD,Mat*);
94: PETSC_EXTERN PetscErrorCode SVDSetInitialSpace(SVD,PetscInt,Vec*);
95: PETSC_EXTERN PetscErrorCode SVDSetInitialSpaceLeft(SVD,PetscInt,Vec*);
96: PETSC_EXTERN PetscErrorCode SVDSetImplicitTranspose(SVD,PetscBool);
97: PETSC_EXTERN PetscErrorCode SVDGetImplicitTranspose(SVD,PetscBool*);
98: PETSC_EXTERN PetscErrorCode SVDSetDimensions(SVD,PetscInt,PetscInt,PetscInt);
99: PETSC_EXTERN PetscErrorCode SVDGetDimensions(SVD,PetscInt*,PetscInt*,PetscInt*);
100: PETSC_EXTERN PetscErrorCode SVDSetTolerances(SVD,PetscReal,PetscInt);
101: PETSC_EXTERN PetscErrorCode SVDGetTolerances(SVD,PetscReal*,PetscInt*);
102: PETSC_EXTERN PetscErrorCode SVDSetWhichSingularTriplets(SVD,SVDWhich);
103: PETSC_EXTERN PetscErrorCode SVDGetWhichSingularTriplets(SVD,SVDWhich*);
104: PETSC_EXTERN PetscErrorCode SVDSetFromOptions(SVD);
105: PETSC_EXTERN PetscErrorCode SVDSetOptionsPrefix(SVD,const char*);
106: PETSC_EXTERN PetscErrorCode SVDAppendOptionsPrefix(SVD,const char*);
107: PETSC_EXTERN PetscErrorCode SVDGetOptionsPrefix(SVD,const char*[]);
108: PETSC_EXTERN PetscErrorCode SVDSetUp(SVD);
109: PETSC_EXTERN PetscErrorCode SVDSolve(SVD);
110: PETSC_EXTERN PetscErrorCode SVDGetIterationNumber(SVD,PetscInt*);
111: PETSC_EXTERN PetscErrorCode SVDGetConvergedReason(SVD,SVDConvergedReason*);
112: PETSC_EXTERN PetscErrorCode SVDGetConverged(SVD,PetscInt*);
113: PETSC_EXTERN PetscErrorCode SVDGetSingularTriplet(SVD,PetscInt,PetscReal*,Vec,Vec);
114: PETSC_EXTERN PetscErrorCode SVDComputeResidualNorms(SVD,PetscInt,PetscReal*,PetscReal*);
115: PETSC_EXTERN PetscErrorCode SVDComputeRelativeError(SVD,PetscInt,PetscReal*);
116: PETSC_EXTERN PetscErrorCode SVDView(SVD,PetscViewer);
117: PETSC_EXTERN PetscErrorCode SVDPrintSolution(SVD,PetscViewer);
118: PETSC_EXTERN PetscErrorCode SVDDestroy(SVD*);
119: PETSC_EXTERN PetscErrorCode SVDReset(SVD);
121: PETSC_EXTERN PetscErrorCode SVDMonitor(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt);
122: PETSC_EXTERN PetscErrorCode SVDMonitorSet(SVD,PetscErrorCode (*)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
123: PETSC_EXTERN PetscErrorCode SVDMonitorCancel(SVD);
124: PETSC_EXTERN PetscErrorCode SVDGetMonitorContext(SVD,void **);
125: PETSC_EXTERN PetscErrorCode SVDMonitorAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
126: PETSC_EXTERN PetscErrorCode SVDMonitorFirst(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
127: PETSC_EXTERN PetscErrorCode SVDMonitorConverged(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
128: PETSC_EXTERN PetscErrorCode SVDMonitorLG(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
129: PETSC_EXTERN PetscErrorCode SVDMonitorLGAll(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
131: PETSC_EXTERN PetscErrorCode SVDSetTrackAll(SVD,PetscBool);
132: PETSC_EXTERN PetscErrorCode SVDGetTrackAll(SVD,PetscBool*);
134: PETSC_EXTERN PetscErrorCode SVDCrossSetEPS(SVD,EPS);
135: PETSC_EXTERN PetscErrorCode SVDCrossGetEPS(SVD,EPS*);
137: PETSC_EXTERN PetscErrorCode SVDCyclicSetExplicitMatrix(SVD,PetscBool);
138: PETSC_EXTERN PetscErrorCode SVDCyclicGetExplicitMatrix(SVD,PetscBool*);
139: PETSC_EXTERN PetscErrorCode SVDCyclicSetEPS(SVD,EPS);
140: PETSC_EXTERN PetscErrorCode SVDCyclicGetEPS(SVD,EPS*);
142: PETSC_EXTERN PetscErrorCode SVDLanczosSetOneSide(SVD,PetscBool);
143: PETSC_EXTERN PetscErrorCode SVDLanczosGetOneSide(SVD,PetscBool*);
145: PETSC_EXTERN PetscErrorCode SVDTRLanczosSetOneSide(SVD,PetscBool);
146: PETSC_EXTERN PetscErrorCode SVDTRLanczosGetOneSide(SVD,PetscBool*);
148: PETSC_EXTERN PetscFunctionList SVDList;
149: PETSC_EXTERN PetscBool SVDRegisterAllCalled;
150: PETSC_EXTERN PetscErrorCode SVDRegisterAll(void);
151: PETSC_EXTERN PetscErrorCode SVDRegister(const char[],PetscErrorCode(*)(SVD));
153: PETSC_EXTERN PetscErrorCode SVDAllocateSolution(SVD,PetscInt);
155: #endif