PartialSVDSolver.h
Go to the documentation of this file.
1 // Copyright (C) 2018 Yixuan Qiu <yixuan.qiu@cos.name>
2 //
3 // This Source Code Form is subject to the terms of the Mozilla
4 // Public License v. 2.0. If a copy of the MPL was not distributed
5 // with this file, You can obtain one at https://mozilla.org/MPL/2.0/.
6 
7 #ifndef PARTIAL_SVD_SOLVER_H
8 #define PARTIAL_SVD_SOLVER_H
9 
10 #include <Eigen/Core>
11 #include "../SymEigsSolver.h"
12 
13 namespace Spectra {
14 
15 // Abstract class for matrix operation
16 template <typename Scalar>
17 class SVDMatOp
18 {
19 public:
20  virtual int rows() const = 0;
21  virtual int cols() const = 0;
22 
23  // y_out = A' * A * x_in or y_out = A * A' * x_in
24  virtual void perform_op(const Scalar* x_in, Scalar* y_out) = 0;
25 
26  virtual ~SVDMatOp() {}
27 };
28 
29 // Operation of a tall matrix in SVD
30 // We compute the eigenvalues of A' * A
31 // MatrixType is either Eigen::Matrix<Scalar, ...> or Eigen::SparseMatrix<Scalar, ...>
32 template <typename Scalar, typename MatrixType>
33 class SVDTallMatOp : public SVDMatOp<Scalar>
34 {
35 private:
40 
42  const int m_dim;
44 
45 public:
46  // Constructor
48  m_mat(mat),
49  m_dim(std::min(mat.rows(), mat.cols())),
50  m_cache(mat.rows())
51  {}
52 
53  // These are the rows and columns of A' * A
54  int rows() const { return m_dim; }
55  int cols() const { return m_dim; }
56 
57  // y_out = A' * A * x_in
58  void perform_op(const Scalar* x_in, Scalar* y_out)
59  {
60  MapConstVec x(x_in, m_mat.cols());
61  MapVec y(y_out, m_mat.cols());
62  m_cache.noalias() = m_mat * x;
63  y.noalias() = m_mat.transpose() * m_cache;
64  }
65 };
66 
67 // Operation of a wide matrix in SVD
68 // We compute the eigenvalues of A * A'
69 // MatrixType is either Eigen::Matrix<Scalar, ...> or Eigen::SparseMatrix<Scalar, ...>
70 template <typename Scalar, typename MatrixType>
71 class SVDWideMatOp : public SVDMatOp<Scalar>
72 {
73 private:
78 
80  const int m_dim;
82 
83 public:
84  // Constructor
86  m_mat(mat),
87  m_dim(std::min(mat.rows(), mat.cols())),
88  m_cache(mat.cols())
89  {}
90 
91  // These are the rows and columns of A * A'
92  int rows() const { return m_dim; }
93  int cols() const { return m_dim; }
94 
95  // y_out = A * A' * x_in
96  void perform_op(const Scalar* x_in, Scalar* y_out)
97  {
98  MapConstVec x(x_in, m_mat.rows());
99  MapVec y(y_out, m_mat.rows());
100  m_cache.noalias() = m_mat.transpose() * x;
101  y.noalias() = m_mat * m_cache;
102  }
103 };
104 
105 // Partial SVD solver
106 // MatrixType is either Eigen::Matrix<Scalar, ...> or Eigen::SparseMatrix<Scalar, ...>
107 template <typename Scalar = double,
110 {
111 private:
115 
117  const int m_m;
118  const int m_n;
121  int m_nconv;
123 
124 public:
125  // Constructor
126  PartialSVDSolver(ConstGenericMatrix& mat, int ncomp, int ncv) :
127  m_mat(mat), m_m(mat.rows()), m_n(mat.cols()), m_evecs(0, 0)
128  {
129  // Determine the matrix type, tall or wide
130  if (m_m > m_n)
131  {
133  }
134  else
135  {
137  }
138 
139  // Solver object
141  }
142 
143  // Destructor
145  {
146  delete m_eigs;
147  delete m_op;
148  }
149 
150  // Computation
151  int compute(int maxit = 1000, Scalar tol = 1e-10)
152  {
153  m_eigs->init();
154  m_nconv = m_eigs->compute(maxit, tol);
155 
156  return m_nconv;
157  }
158 
159  // The converged singular values
161  {
162  Vector svals = m_eigs->eigenvalues().cwiseSqrt();
163 
164  return svals;
165  }
166 
167  // The converged left singular vectors
168  Matrix matrix_U(int nu)
169  {
170  if (m_evecs.cols() < 1)
171  {
172  m_evecs = m_eigs->eigenvectors();
173  }
174  nu = std::min(nu, m_nconv);
175  if (m_m <= m_n)
176  {
177  return m_evecs.leftCols(nu);
178  }
179 
180  return m_mat * (m_evecs.leftCols(nu).array().rowwise() / m_eigs->eigenvalues().head(nu).transpose().array().sqrt()).matrix();
181  }
182 
183  // The converged right singular vectors
184  Matrix matrix_V(int nv)
185  {
186  if (m_evecs.cols() < 1)
187  {
188  m_evecs = m_eigs->eigenvectors();
189  }
190  nv = std::min(nv, m_nconv);
191  if (m_m > m_n)
192  {
193  return m_evecs.leftCols(nv);
194  }
195 
196  return m_mat.transpose() * (m_evecs.leftCols(nv).array().rowwise() / m_eigs->eigenvalues().head(nv).transpose().array().sqrt()).matrix();
197  }
198 };
199 
200 } // namespace Spectra
201 
202 #endif // PARTIAL_SVD_SOLVER_H
Spectra::PartialSVDSolver::matrix_V
Matrix matrix_V(int nv)
Definition: PartialSVDSolver.h:184
Spectra::SVDWideMatOp::m_cache
Vector m_cache
Definition: PartialSVDSolver.h:81
Spectra::SVDMatOp::~SVDMatOp
virtual ~SVDMatOp()
Definition: PartialSVDSolver.h:26
Spectra::SVDWideMatOp::cols
int cols() const
Definition: PartialSVDSolver.h:93
e
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Spectra::SVDTallMatOp::m_mat
ConstGenericMatrix m_mat
Definition: PartialSVDSolver.h:41
MatrixType
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
Spectra::SVDTallMatOp::rows
int rows() const
Definition: PartialSVDSolver.h:54
Spectra::PartialSVDSolver::~PartialSVDSolver
virtual ~PartialSVDSolver()
Definition: PartialSVDSolver.h:144
Spectra::SVDWideMatOp::rows
int rows() const
Definition: PartialSVDSolver.h:92
Spectra::PartialSVDSolver::m_evecs
Matrix m_evecs
Definition: PartialSVDSolver.h:122
Spectra::SVDTallMatOp::Vector
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
Definition: PartialSVDSolver.h:36
x
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition: gnuplot_common_settings.hh:12
Spectra::PartialSVDSolver::m_nconv
int m_nconv
Definition: PartialSVDSolver.h:121
Spectra::SymEigsSolver
Definition: SymEigsSolver.h:141
Spectra::SVDMatOp::cols
virtual int cols() const =0
mat
MatrixXf mat
Definition: Tutorial_AdvancedInitialization_CommaTemporary.cpp:1
Spectra::SVDTallMatOp::ConstGenericMatrix
const typedef Eigen::Ref< const MatrixType > ConstGenericMatrix
Definition: PartialSVDSolver.h:39
Spectra::SVDTallMatOp::SVDTallMatOp
SVDTallMatOp(ConstGenericMatrix &mat)
Definition: PartialSVDSolver.h:47
rows
int rows
Definition: Tutorial_commainit_02.cpp:1
Spectra::PartialSVDSolver::Matrix
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > Matrix
Definition: PartialSVDSolver.h:112
Spectra::PartialSVDSolver::m_n
const int m_n
Definition: PartialSVDSolver.h:118
Spectra::SVDWideMatOp::m_mat
ConstGenericMatrix m_mat
Definition: PartialSVDSolver.h:79
Spectra::PartialSVDSolver
Definition: PartialSVDSolver.h:109
Spectra::SVDWideMatOp
Definition: PartialSVDSolver.h:71
Spectra::PartialSVDSolver::compute
int compute(int maxit=1000, Scalar tol=1e-10)
Definition: PartialSVDSolver.h:151
Spectra::PartialSVDSolver::matrix_U
Matrix matrix_U(int nu)
Definition: PartialSVDSolver.h:168
Spectra::SVDTallMatOp::MapConstVec
Eigen::Map< const Vector > MapConstVec
Definition: PartialSVDSolver.h:37
Spectra::PartialSVDSolver::m_eigs
SymEigsSolver< Scalar, LARGEST_ALGE, SVDMatOp< Scalar > > * m_eigs
Definition: PartialSVDSolver.h:120
Spectra::SVDMatOp::perform_op
virtual void perform_op(const Scalar *x_in, Scalar *y_out)=0
Spectra::PartialSVDSolver::singular_values
Vector singular_values() const
Definition: PartialSVDSolver.h:160
Spectra::SVDWideMatOp::m_dim
const int m_dim
Definition: PartialSVDSolver.h:80
Spectra::SVDTallMatOp::m_cache
Vector m_cache
Definition: PartialSVDSolver.h:43
Spectra::SVDTallMatOp::MapVec
Eigen::Map< Vector > MapVec
Definition: PartialSVDSolver.h:38
Spectra::SVDTallMatOp
Definition: PartialSVDSolver.h:33
Spectra::PartialSVDSolver::Vector
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
Definition: PartialSVDSolver.h:113
Spectra::SVDWideMatOp::perform_op
void perform_op(const Scalar *x_in, Scalar *y_out)
Definition: PartialSVDSolver.h:96
Spectra::PartialSVDSolver::ConstGenericMatrix
const typedef Eigen::Ref< const MatrixType > ConstGenericMatrix
Definition: PartialSVDSolver.h:114
Spectra::SVDTallMatOp::cols
int cols() const
Definition: PartialSVDSolver.h:55
Spectra::SVDTallMatOp::m_dim
const int m_dim
Definition: PartialSVDSolver.h:42
Eigen::Map
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
Spectra::SVDWideMatOp::MapVec
Eigen::Map< Vector > MapVec
Definition: PartialSVDSolver.h:76
Spectra::SVDWideMatOp::MapConstVec
Eigen::Map< const Vector > MapConstVec
Definition: PartialSVDSolver.h:75
y
Scalar * y
Definition: level1_cplx_impl.h:124
matrix
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: gtsam/3rdparty/Eigen/blas/common.h:110
Spectra::SVDMatOp::rows
virtual int rows() const =0
Spectra::PartialSVDSolver::m_m
const int m_m
Definition: PartialSVDSolver.h:117
Spectra::SVDMatOp
Definition: PartialSVDSolver.h:17
Eigen::Ref< const MatrixType >
Spectra::PartialSVDSolver::PartialSVDSolver
PartialSVDSolver(ConstGenericMatrix &mat, int ncomp, int ncv)
Definition: PartialSVDSolver.h:126
Eigen::PlainObjectBase::cols
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:145
Spectra::PartialSVDSolver::m_op
SVDMatOp< Scalar > * m_op
Definition: PartialSVDSolver.h:119
Spectra::SVDWideMatOp::Vector
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
Definition: PartialSVDSolver.h:74
std
Definition: BFloat16.h:88
Spectra
Definition: LOBPCGSolver.h:19
min
#define min(a, b)
Definition: datatypes.h:19
gtsam::tol
const G double tol
Definition: Group.h:79
Spectra::SVDTallMatOp::perform_op
void perform_op(const Scalar *x_in, Scalar *y_out)
Definition: PartialSVDSolver.h:58
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 >
Spectra::PartialSVDSolver::m_mat
ConstGenericMatrix m_mat
Definition: PartialSVDSolver.h:116
cols
int cols
Definition: Tutorial_commainit_02.cpp:1
Spectra::SVDWideMatOp::SVDWideMatOp
SVDWideMatOp(ConstGenericMatrix &mat)
Definition: PartialSVDSolver.h:85
Scalar
SCALAR Scalar
Definition: bench_gemm.cpp:46
Spectra::SVDWideMatOp::ConstGenericMatrix
const typedef Eigen::Ref< const MatrixType > ConstGenericMatrix
Definition: PartialSVDSolver.h:77


gtsam
Author(s):
autogenerated on Sun Dec 22 2024 04:12:41