BDCSVD_LAPACKE.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2022 Melven Roehrig-Zoellner <Melven.Roehrig-Zoellner@DLR.de>
5 // Copyright (c) 2011, Intel Corporation. All rights reserved.
6 //
7 // This file is based on the JacobiSVD_LAPACKE.h originally from Intel -
8 // see license notice below:
9 /*
10  Redistribution and use in source and binary forms, with or without modification,
11  are permitted provided that the following conditions are met:
12 
13  * Redistributions of source code must retain the above copyright notice, this
14  list of conditions and the following disclaimer.
15  * Redistributions in binary form must reproduce the above copyright notice,
16  this list of conditions and the following disclaimer in the documentation
17  and/or other materials provided with the distribution.
18  * Neither the name of Intel Corporation nor the names of its contributors may
19  be used to endorse or promote products derived from this software without
20  specific prior written permission.
21 
22  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
23  ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
24  WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
25  DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR
26  ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
27  (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
28  LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
29  ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
30  (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
31  SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
32 
33  ********************************************************************************
34  * Content : Eigen bindings to LAPACKe
35  * Singular Value Decomposition - SVD (divide and conquer variant)
36  ********************************************************************************
37 */
38 #ifndef EIGEN_BDCSVD_LAPACKE_H
39 #define EIGEN_BDCSVD_LAPACKE_H
40 
41 namespace Eigen {
42 
43 namespace internal {
44 
45 namespace lapacke_helpers {
46 
49 // defining a derived class to allow access to protected members
50 template <typename MatrixType_, int Options>
51 class BDCSVD_LAPACKE : public BDCSVD<MatrixType_, Options> {
53  typedef typename SVD::MatrixType MatrixType;
54  typedef typename SVD::Scalar Scalar;
55  typedef typename SVD::RealScalar RealScalar;
56 
57  public:
58  // construct this by moving from a parent object
59  BDCSVD_LAPACKE(SVD&& svd) : SVD(std::move(svd)) {}
60 
61  void compute_impl_lapacke(const MatrixType& matrix, unsigned int computationOptions) {
62  SVD::allocate(matrix.rows(), matrix.cols(), computationOptions);
63 
65 
66  // prepare arguments to ?gesdd
67  const lapack_int matrix_order = lapack_storage_of(matrix);
68  const char jobz = (SVD::m_computeFullU || SVD::m_computeFullV) ? 'A'
70  : 'N';
71  const lapack_int u_cols = (jobz == 'A') ? to_lapack(SVD::rows()) : (jobz == 'S') ? to_lapack(SVD::diagSize()) : 1;
72  const lapack_int vt_rows = (jobz == 'A') ? to_lapack(SVD::cols()) : (jobz == 'S') ? to_lapack(SVD::diagSize()) : 1;
73  lapack_int ldu, ldvt;
74  Scalar *u, *vt, dummy;
75  MatrixType localU;
77  ldu = to_lapack(SVD::m_matrixU.outerStride());
78  u = SVD::m_matrixU.data();
79  } else if (SVD::computeV()) {
80  localU.resize(SVD::rows(), u_cols);
81  ldu = to_lapack(localU.outerStride());
82  u = localU.data();
83  } else {
84  ldu = 1;
85  u = &dummy;
86  }
87  MatrixType localV;
88  if (SVD::computeU() || SVD::computeV()) {
89  localV.resize(vt_rows, SVD::cols());
90  ldvt = to_lapack(localV.outerStride());
91  vt = localV.data();
92  } else {
93  ldvt = 1;
94  vt = &dummy;
95  }
96  MatrixType temp;
97  temp = matrix;
98 
99  // actual call to ?gesdd
100  lapack_int info = gesdd(matrix_order, jobz, to_lapack(SVD::rows()), to_lapack(SVD::cols()), to_lapack(temp.data()),
101  to_lapack(temp.outerStride()), (RealScalar*)SVD::m_singularValues.data(), to_lapack(u), ldu,
102  to_lapack(vt), ldvt);
103 
104  // Check the result of the LAPACK call
105  if (info < 0 || !SVD::m_singularValues.allFinite()) {
106  // this includes info == -4 => NaN entry in A
108  } else if (info > 0) {
110  } else {
113  SVD::m_matrixU = localU.leftCols(SVD::m_matrixU.cols());
114  }
115  if (SVD::computeV()) {
116  SVD::m_matrixV = localV.adjoint().leftCols(SVD::m_matrixV.cols());
117  }
118  }
119  SVD::m_isInitialized = true;
120  }
121 };
122 
123 template <typename MatrixType_, int Options>
125  int computationOptions) {
126  // we need to move to the wrapper type and back
127  BDCSVD_LAPACKE<MatrixType_, Options> tmpSvd(std::move(svd));
128  tmpSvd.compute_impl_lapacke(matrix, computationOptions);
129  svd = std::move(tmpSvd);
130  return svd;
131 }
132 
133 } // end namespace lapacke_helpers
134 
135 } // end namespace internal
136 
137 #define EIGEN_LAPACKE_SDD(EIGTYPE, EIGCOLROW, OPTIONS) \
138  template <> \
139  inline BDCSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>& \
140  BDCSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>::compute_impl( \
141  const Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>& matrix, unsigned int computationOptions) { \
142  return internal::lapacke_helpers::BDCSVD_wrapper(*this, matrix, computationOptions); \
143  }
144 
145 #define EIGEN_LAPACK_SDD_OPTIONS(OPTIONS) \
146  EIGEN_LAPACKE_SDD(double, ColMajor, OPTIONS) \
147  EIGEN_LAPACKE_SDD(float, ColMajor, OPTIONS) \
148  EIGEN_LAPACKE_SDD(dcomplex, ColMajor, OPTIONS) \
149  EIGEN_LAPACKE_SDD(scomplex, ColMajor, OPTIONS) \
150  \
151  EIGEN_LAPACKE_SDD(double, RowMajor, OPTIONS) \
152  EIGEN_LAPACKE_SDD(float, RowMajor, OPTIONS) \
153  EIGEN_LAPACKE_SDD(dcomplex, RowMajor, OPTIONS) \
154  EIGEN_LAPACKE_SDD(scomplex, RowMajor, OPTIONS)
155 
165 
166 #undef EIGEN_LAPACK_SDD_OPTIONS
167 
168 #undef EIGEN_LAPACKE_SDD
169 
170 } // end namespace Eigen
171 
172 #endif // EIGEN_BDCSVD_LAPACKE_H
#define EIGEN_LAPACK_SDD_OPTIONS(OPTIONS)
Definition: BDCSVD_LAPACKE.h:145
cout<< "Here is the matrix m:"<< endl<< m<< endl;JacobiSVD< MatrixXf, ComputeThinU|ComputeThinV > svd(m)
class Bidiagonal Divide and Conquer SVD
Definition: BDCSVD.h:85
void allocate(Index rows, Index cols, unsigned int computationOptions)
Definition: BDCSVD.h:268
Base::RealScalar RealScalar
Definition: BDCSVD.h:97
MatrixType_ MatrixType
Definition: BDCSVD.h:95
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: EigenBase.h:61
Base::Scalar Scalar
Definition: BDCSVD.h:96
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:59
constexpr EIGEN_DEVICE_FUNC const Scalar * data() const
Definition: PlainObjectBase.h:273
bool m_computeFullV
Definition: SVDBase.h:337
internal::variable_if_dynamic< Index, DiagSizeAtCompileTime > m_diagSize
Definition: SVDBase.h:342
ComputationInfo m_info
Definition: SVDBase.h:334
bool m_computeThinU
Definition: SVDBase.h:336
EIGEN_DEVICE_FUNC ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SVDBase.h:300
bool computeV() const
Definition: SVDBase.h:275
bool m_isInitialized
Definition: SVDBase.h:335
MatrixVType m_matrixV
Definition: SVDBase.h:332
bool computeU() const
Definition: SVDBase.h:273
Index m_nonzeroSingularValues
Definition: SVDBase.h:339
SingularValuesType m_singularValues
Definition: SVDBase.h:333
bool m_computeThinV
Definition: SVDBase.h:337
bool m_computeFullU
Definition: SVDBase.h:336
Index diagSize() const
Definition: SVDBase.h:279
MatrixUType m_matrixU
Definition: SVDBase.h:331
internal::traits< Derived >::Scalar Scalar
Definition: SolverBase.h:75
SVD::RealScalar RealScalar
Definition: BDCSVD_LAPACKE.h:55
BDCSVD< MatrixType_, Options > SVD
Definition: BDCSVD_LAPACKE.h:52
BDCSVD_LAPACKE(SVD &&svd)
Definition: BDCSVD_LAPACKE.h:59
SVD::Scalar Scalar
Definition: BDCSVD_LAPACKE.h:54
SVD::MatrixType MatrixType
Definition: BDCSVD_LAPACKE.h:53
void compute_impl_lapacke(const MatrixType &matrix, unsigned int computationOptions)
Definition: BDCSVD_LAPACKE.h:61
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85
@ InvalidInput
Definition: Constants.h:447
@ Success
Definition: Constants.h:440
@ NoConvergence
Definition: Constants.h:444
@ ComputeFullV
Definition: Constants.h:393
@ ComputeThinV
Definition: Constants.h:395
@ ComputeFullU
Definition: Constants.h:389
@ ComputeThinU
Definition: Constants.h:391
#define lapack_int
Definition: lapacke.h:52
EIGEN_ALWAYS_INLINE auto to_lapack(Source value)
Definition: lapacke_helpers.h:61
EIGEN_ALWAYS_INLINE EIGEN_CONSTEXPR lapack_int lapack_storage_of(const EigenBase< Derived > &)
translates storage order of the given Eigen object to the corresponding lapack constant
Definition: lapacke_helpers.h:78
BDCSVD< MatrixType_, Options > & BDCSVD_wrapper(BDCSVD< MatrixType_, Options > &svd, const MatrixType_ &matrix, int computationOptions)
Definition: BDCSVD_LAPACKE.h:124
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
Definition: Eigen_Colamd.h:49