38#ifndef EIGEN_BDCSVD_LAPACKE_H
39#define EIGEN_BDCSVD_LAPACKE_H
45namespace lapacke_helpers {
50template <
typename MatrixType_,
int Options>
51class BDCSVD_LAPACKE :
public BDCSVD<MatrixType_, Options> {
52 typedef BDCSVD<MatrixType_, Options> SVD;
53 typedef typename SVD::MatrixType MatrixType;
54 typedef typename SVD::Scalar Scalar;
55 typedef typename SVD::RealScalar RealScalar;
59 BDCSVD_LAPACKE(SVD&& svd) : SVD(std::move(svd)) {}
61 void compute_impl_lapacke(
const MatrixType& matrix,
unsigned int computationOptions) {
62 SVD::allocate(matrix.rows(), matrix.cols(), computationOptions);
64 SVD::m_nonzeroSingularValues = SVD::m_diagSize;
67 const lapack_int matrix_order = lapack_storage_of(matrix);
68 const char jobz = (SVD::m_computeFullU || SVD::m_computeFullV) ?
'A'
69 : (SVD::m_computeThinU || SVD::m_computeThinV) ?
'S'
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;
74 Scalar *u, *vt, dummy;
76 if (
SVD::computeU() && !(SVD::m_computeThinU && SVD::m_computeFullV)) {
77 ldu = to_lapack(SVD::m_matrixU.outerStride());
78 u = SVD::m_matrixU.data();
80 localU.resize(SVD::rows(), u_cols);
81 ldu = to_lapack(localU.outerStride());
89 localV.resize(vt_rows, SVD::cols());
90 ldvt = to_lapack(localV.outerStride());
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);
105 if (
info < 0 || !SVD::m_singularValues.allFinite()) {
108 }
else if (
info > 0) {
112 if (SVD::m_computeThinU && SVD::m_computeFullV) {
113 SVD::m_matrixU = localU.leftCols(SVD::m_matrixU.cols());
116 SVD::m_matrixV = localV.adjoint().leftCols(SVD::m_matrixV.cols());
119 SVD::m_isInitialized =
true;
123template <
typename MatrixType_,
int Options>
124BDCSVD<MatrixType_, Options>& BDCSVD_wrapper(BDCSVD<MatrixType_, Options>& svd,
const MatrixType_& matrix,
125 int computationOptions) {
127 BDCSVD_LAPACKE<MatrixType_, Options> tmpSvd(std::move(svd));
128 tmpSvd.compute_impl_lapacke(matrix, computationOptions);
129 svd = std::move(tmpSvd);
137#define EIGEN_LAPACKE_SDD(EIGTYPE, EIGCOLROW, OPTIONS) \
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); \
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) \
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)
156EIGEN_LAPACK_SDD_OPTIONS(0)
157EIGEN_LAPACK_SDD_OPTIONS(ComputeThinU)
158EIGEN_LAPACK_SDD_OPTIONS(ComputeThinV)
159EIGEN_LAPACK_SDD_OPTIONS(ComputeFullU)
160EIGEN_LAPACK_SDD_OPTIONS(ComputeFullV)
161EIGEN_LAPACK_SDD_OPTIONS(ComputeThinU | ComputeThinV)
162EIGEN_LAPACK_SDD_OPTIONS(ComputeFullU | ComputeFullV)
163EIGEN_LAPACK_SDD_OPTIONS(ComputeThinU | ComputeFullV)
164EIGEN_LAPACK_SDD_OPTIONS(ComputeFullU | ComputeThinV)
166#undef EIGEN_LAPACK_SDD_OPTIONS
168#undef EIGEN_LAPACKE_SDD
bool computeV() const
Definition SVDBase.h:275
bool computeU() const
Definition SVDBase.h:273
ComputationInfo info() const
Definition SVDBase.h:300
@ InvalidInput
Definition Constants.h:447
@ Success
Definition Constants.h:440
@ NoConvergence
Definition Constants.h:444
Namespace containing all symbols from the Eigen library.
Definition Core:137