38 #ifndef EIGEN_BDCSVD_LAPACKE_H 39 #define EIGEN_BDCSVD_LAPACKE_H 45 namespace lapacke_helpers {
50 template <
typename MatrixType_,
int Options>
51 class 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 template <
typename Derived>
62 void compute_impl_lapacke(
const MatrixBase<Derived>& matrix,
unsigned int computationOptions) {
63 SVD::allocate(matrix.rows(), matrix.cols(), computationOptions);
65 SVD::m_nonzeroSingularValues = SVD::m_diagSize;
68 const lapack_int matrix_order = lapack_storage_of(matrix);
69 const char jobz = (SVD::m_computeFullU || SVD::m_computeFullV) ?
'A' 70 : (SVD::m_computeThinU || SVD::m_computeThinV) ?
'S' 72 const lapack_int u_cols = (jobz ==
'A') ? to_lapack(SVD::rows()) : (jobz ==
'S') ? to_lapack(SVD::diagSize()) : 1;
73 const lapack_int vt_rows = (jobz ==
'A') ? to_lapack(SVD::cols()) : (jobz ==
'S') ? to_lapack(SVD::diagSize()) : 1;
75 Scalar *u, *vt, dummy;
77 if (
SVD::computeU() && !(SVD::m_computeThinU && SVD::m_computeFullV)) {
78 ldu = to_lapack(SVD::m_matrixU.outerStride());
79 u = SVD::m_matrixU.data();
81 localU.resize(SVD::rows(), u_cols);
82 ldu = to_lapack(localU.outerStride());
90 localV.resize(vt_rows, SVD::cols());
91 ldvt = to_lapack(localV.outerStride());
101 lapack_int
info = gesdd(matrix_order, jobz, to_lapack(SVD::rows()), to_lapack(SVD::cols()), to_lapack(temp.data()),
102 to_lapack(temp.outerStride()), (RealScalar*)SVD::m_singularValues.data(), to_lapack(u), ldu,
103 to_lapack(vt), ldvt);
106 if (info < 0 || !SVD::m_singularValues.allFinite()) {
109 }
else if (info > 0) {
113 if (SVD::m_computeThinU && SVD::m_computeFullV) {
114 SVD::m_matrixU = localU.leftCols(SVD::m_matrixU.cols());
117 SVD::m_matrixV = localV.adjoint().leftCols(SVD::m_matrixV.cols());
120 SVD::m_isInitialized =
true;
124 template <
typename MatrixType_,
int Options,
typename Derived>
125 BDCSVD<MatrixType_, Options>& BDCSVD_wrapper(BDCSVD<MatrixType_, Options>& svd,
const MatrixBase<Derived>& matrix,
126 int computationOptions) {
128 BDCSVD_LAPACKE<MatrixType_, Options> tmpSvd(std::move(svd));
129 tmpSvd.compute_impl_lapacke(matrix, computationOptions);
130 svd = std::move(tmpSvd);
138 #define EIGEN_LAPACKE_SDD(EIGTYPE, EIGCOLROW, OPTIONS) \ 140 template <typename Derived> \ 141 inline BDCSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>& \ 142 BDCSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>::compute_impl( \ 143 const MatrixBase<Derived>& matrix, unsigned int computationOptions) { \ 144 return internal::lapacke_helpers::BDCSVD_wrapper(*this, matrix, computationOptions); \ 147 #define EIGEN_LAPACK_SDD_OPTIONS(OPTIONS) \ 148 EIGEN_LAPACKE_SDD(double, ColMajor, OPTIONS) \ 149 EIGEN_LAPACKE_SDD(float, ColMajor, OPTIONS) \ 150 EIGEN_LAPACKE_SDD(dcomplex, ColMajor, OPTIONS) \ 151 EIGEN_LAPACKE_SDD(scomplex, ColMajor, OPTIONS) \ 153 EIGEN_LAPACKE_SDD(double, RowMajor, OPTIONS) \ 154 EIGEN_LAPACKE_SDD(float, RowMajor, OPTIONS) \ 155 EIGEN_LAPACKE_SDD(dcomplex, RowMajor, OPTIONS) \ 156 EIGEN_LAPACKE_SDD(scomplex, RowMajor, OPTIONS) 158 EIGEN_LAPACK_SDD_OPTIONS(0)
168 #undef EIGEN_LAPACK_SDD_OPTIONS 170 #undef EIGEN_LAPACKE_SDD 174 #endif // EIGEN_BDCSVD_LAPACKE_H Definition: Constants.h:389
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SVDBase.h:300
Definition: Constants.h:395
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
Definition: BFloat16.h:231
bool computeV() const
Definition: SVDBase.h:275
Definition: Constants.h:447
Definition: Constants.h:440
bool computeU() const
Definition: SVDBase.h:273
Definition: Constants.h:393
Definition: Constants.h:391
Definition: Constants.h:444