mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-05-30 01:53:59 +08:00
Fix: Correct Lapacke bindings for BDCSVD and JacobiSVD to match the updated API
This commit is contained in:
parent
434a2fc4a4
commit
f3e7d64f3d
@ -182,7 +182,9 @@ class BDCSVD : public SVDBase<BDCSVD<MatrixType_, Options_> > {
|
||||
* \deprecated Will be removed in the next major Eigen version. Options should
|
||||
* be specified in the \a Options template parameter.
|
||||
*/
|
||||
EIGEN_DEPRECATED BDCSVD(const MatrixType& matrix, unsigned int computationOptions) : m_algoswap(16), m_numIters(0) {
|
||||
template <typename Derived>
|
||||
EIGEN_DEPRECATED BDCSVD(const MatrixBase<Derived>& matrix, unsigned int computationOptions)
|
||||
: m_algoswap(16), m_numIters(0) {
|
||||
internal::check_svd_options_assertions<MatrixType, Options>(computationOptions, matrix.rows(), matrix.cols());
|
||||
compute_impl(matrix, computationOptions);
|
||||
}
|
||||
|
@ -58,7 +58,8 @@ class BDCSVD_LAPACKE : public BDCSVD<MatrixType_, Options> {
|
||||
// construct this by moving from a parent object
|
||||
BDCSVD_LAPACKE(SVD&& svd) : SVD(std::move(svd)) {}
|
||||
|
||||
void compute_impl_lapacke(const MatrixType& matrix, unsigned int computationOptions) {
|
||||
template <typename Derived>
|
||||
void compute_impl_lapacke(const MatrixBase<Derived>& matrix, unsigned int computationOptions) {
|
||||
SVD::allocate(matrix.rows(), matrix.cols(), computationOptions);
|
||||
|
||||
SVD::m_nonzeroSingularValues = SVD::m_diagSize;
|
||||
@ -120,8 +121,8 @@ class BDCSVD_LAPACKE : public BDCSVD<MatrixType_, Options> {
|
||||
}
|
||||
};
|
||||
|
||||
template <typename MatrixType_, int Options>
|
||||
BDCSVD<MatrixType_, Options>& BDCSVD_wrapper(BDCSVD<MatrixType_, Options>& svd, const MatrixType_& matrix,
|
||||
template <typename MatrixType_, int Options, typename Derived>
|
||||
BDCSVD<MatrixType_, Options>& BDCSVD_wrapper(BDCSVD<MatrixType_, Options>& svd, const MatrixBase<Derived>& matrix,
|
||||
int computationOptions) {
|
||||
// we need to move to the wrapper type and back
|
||||
BDCSVD_LAPACKE<MatrixType_, Options> tmpSvd(std::move(svd));
|
||||
@ -134,12 +135,13 @@ BDCSVD<MatrixType_, Options>& BDCSVD_wrapper(BDCSVD<MatrixType_, Options>& svd,
|
||||
|
||||
} // end namespace internal
|
||||
|
||||
#define EIGEN_LAPACKE_SDD(EIGTYPE, EIGCOLROW, OPTIONS) \
|
||||
template <> \
|
||||
inline BDCSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>& \
|
||||
BDCSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>::compute_impl( \
|
||||
const Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>& matrix, unsigned int computationOptions) { \
|
||||
return internal::lapacke_helpers::BDCSVD_wrapper(*this, matrix, computationOptions); \
|
||||
#define EIGEN_LAPACKE_SDD(EIGTYPE, EIGCOLROW, OPTIONS) \
|
||||
template <> \
|
||||
template <typename Derived> \
|
||||
inline BDCSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>& \
|
||||
BDCSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>::compute_impl( \
|
||||
const MatrixBase<Derived>& matrix, unsigned int computationOptions) { \
|
||||
return internal::lapacke_helpers::BDCSVD_wrapper(*this, matrix, computationOptions); \
|
||||
}
|
||||
|
||||
#define EIGEN_LAPACK_SDD_OPTIONS(OPTIONS) \
|
||||
|
@ -40,65 +40,65 @@ namespace Eigen {
|
||||
|
||||
/** \internal Specialization for the data types supported by LAPACKe */
|
||||
|
||||
#define EIGEN_LAPACKE_SVD(EIGTYPE, LAPACKE_TYPE, LAPACKE_RTYPE, LAPACKE_PREFIX, EIGCOLROW, LAPACKE_COLROW, OPTIONS) \
|
||||
template <> \
|
||||
inline JacobiSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>& \
|
||||
JacobiSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>::compute_impl( \
|
||||
const Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>& matrix, unsigned int computationOptions) { \
|
||||
typedef Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic> MatrixType; \
|
||||
/*typedef MatrixType::Scalar Scalar;*/ \
|
||||
/*typedef MatrixType::RealScalar RealScalar;*/ \
|
||||
allocate(matrix.rows(), matrix.cols(), computationOptions); \
|
||||
\
|
||||
/*const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();*/ \
|
||||
m_nonzeroSingularValues = diagSize(); \
|
||||
\
|
||||
lapack_int lda = internal::convert_index<lapack_int>(matrix.outerStride()), ldu, ldvt; \
|
||||
lapack_int matrix_order = LAPACKE_COLROW; \
|
||||
char jobu, jobvt; \
|
||||
LAPACKE_TYPE *u, *vt, dummy; \
|
||||
jobu = (m_computeFullU) ? 'A' : (m_computeThinU) ? 'S' : 'N'; \
|
||||
jobvt = (m_computeFullV) ? 'A' : (m_computeThinV) ? 'S' : 'N'; \
|
||||
if (computeU()) { \
|
||||
ldu = internal::convert_index<lapack_int>(m_matrixU.outerStride()); \
|
||||
u = (LAPACKE_TYPE*)m_matrixU.data(); \
|
||||
} else { \
|
||||
ldu = 1; \
|
||||
u = &dummy; \
|
||||
} \
|
||||
MatrixType localV; \
|
||||
lapack_int vt_rows = (m_computeFullV) ? internal::convert_index<lapack_int>(cols()) \
|
||||
: (m_computeThinV) ? internal::convert_index<lapack_int>(diagSize()) \
|
||||
: 1; \
|
||||
if (computeV()) { \
|
||||
localV.resize(vt_rows, cols()); \
|
||||
ldvt = internal::convert_index<lapack_int>(localV.outerStride()); \
|
||||
vt = (LAPACKE_TYPE*)localV.data(); \
|
||||
} else { \
|
||||
ldvt = 1; \
|
||||
vt = &dummy; \
|
||||
} \
|
||||
Matrix<LAPACKE_RTYPE, Dynamic, Dynamic> superb; \
|
||||
superb.resize(diagSize(), 1); \
|
||||
MatrixType m_temp; \
|
||||
m_temp = matrix; \
|
||||
lapack_int info = LAPACKE_##LAPACKE_PREFIX##gesvd( \
|
||||
matrix_order, jobu, jobvt, internal::convert_index<lapack_int>(rows()), \
|
||||
internal::convert_index<lapack_int>(cols()), (LAPACKE_TYPE*)m_temp.data(), lda, \
|
||||
(LAPACKE_RTYPE*)m_singularValues.data(), u, ldu, vt, ldvt, superb.data()); \
|
||||
/* Check the result of the LAPACK call */ \
|
||||
if (info < 0 || !m_singularValues.allFinite()) { \
|
||||
m_info = InvalidInput; \
|
||||
} else if (info > 0) { \
|
||||
m_info = NoConvergence; \
|
||||
} else { \
|
||||
m_info = Success; \
|
||||
if (computeV()) m_matrixV = localV.adjoint(); \
|
||||
} \
|
||||
/* for(int i=0;i<diagSize();i++) if (m_singularValues.coeffRef(i) < precision) { m_nonzeroSingularValues--; \
|
||||
* m_singularValues.coeffRef(i)=RealScalar(0);}*/ \
|
||||
m_isInitialized = true; \
|
||||
return *this; \
|
||||
#define EIGEN_LAPACKE_SVD(EIGTYPE, LAPACKE_TYPE, LAPACKE_RTYPE, LAPACKE_PREFIX, EIGCOLROW, LAPACKE_COLROW, OPTIONS) \
|
||||
template <> \
|
||||
template <typename Derived> \
|
||||
inline JacobiSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>& \
|
||||
JacobiSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>::compute_impl( \
|
||||
const MatrixBase<Derived>& matrix, unsigned int computationOptions) { \
|
||||
/*typedef MatrixType::Scalar Scalar;*/ \
|
||||
/*typedef MatrixType::RealScalar RealScalar;*/ \
|
||||
allocate(matrix.rows(), matrix.cols(), computationOptions); \
|
||||
\
|
||||
/*const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();*/ \
|
||||
m_nonzeroSingularValues = diagSize(); \
|
||||
\
|
||||
lapack_int lda = internal::convert_index<lapack_int>(matrix.outerStride()), ldu, ldvt; \
|
||||
lapack_int matrix_order = LAPACKE_COLROW; \
|
||||
char jobu, jobvt; \
|
||||
LAPACKE_TYPE *u, *vt, dummy; \
|
||||
jobu = (m_computeFullU) ? 'A' : (m_computeThinU) ? 'S' : 'N'; \
|
||||
jobvt = (m_computeFullV) ? 'A' : (m_computeThinV) ? 'S' : 'N'; \
|
||||
if (computeU()) { \
|
||||
ldu = internal::convert_index<lapack_int>(m_matrixU.outerStride()); \
|
||||
u = (LAPACKE_TYPE*)m_matrixU.data(); \
|
||||
} else { \
|
||||
ldu = 1; \
|
||||
u = &dummy; \
|
||||
} \
|
||||
MatrixType localV; \
|
||||
lapack_int vt_rows = (m_computeFullV) ? internal::convert_index<lapack_int>(cols()) \
|
||||
: (m_computeThinV) ? internal::convert_index<lapack_int>(diagSize()) \
|
||||
: 1; \
|
||||
if (computeV()) { \
|
||||
localV.resize(vt_rows, cols()); \
|
||||
ldvt = internal::convert_index<lapack_int>(localV.outerStride()); \
|
||||
vt = (LAPACKE_TYPE*)localV.data(); \
|
||||
} else { \
|
||||
ldvt = 1; \
|
||||
vt = &dummy; \
|
||||
} \
|
||||
Matrix<LAPACKE_RTYPE, Dynamic, Dynamic> superb; \
|
||||
superb.resize(diagSize(), 1); \
|
||||
MatrixType m_temp; \
|
||||
m_temp = matrix; \
|
||||
lapack_int info = LAPACKE_##LAPACKE_PREFIX##gesvd( \
|
||||
matrix_order, jobu, jobvt, internal::convert_index<lapack_int>(rows()), \
|
||||
internal::convert_index<lapack_int>(cols()), (LAPACKE_TYPE*)m_temp.data(), lda, \
|
||||
(LAPACKE_RTYPE*)m_singularValues.data(), u, ldu, vt, ldvt, superb.data()); \
|
||||
/* Check the result of the LAPACK call */ \
|
||||
if (info < 0 || !m_singularValues.allFinite()) { \
|
||||
m_info = InvalidInput; \
|
||||
} else if (info > 0) { \
|
||||
m_info = NoConvergence; \
|
||||
} else { \
|
||||
m_info = Success; \
|
||||
if (computeV()) m_matrixV = localV.adjoint(); \
|
||||
} \
|
||||
/* for(int i=0;i<diagSize();i++) if (m_singularValues.coeffRef(i) < precision) { m_nonzeroSingularValues--; \
|
||||
* m_singularValues.coeffRef(i)=RealScalar(0);}*/ \
|
||||
m_isInitialized = true; \
|
||||
return *this; \
|
||||
}
|
||||
|
||||
#define EIGEN_LAPACK_SVD_OPTIONS(OPTIONS) \
|
||||
|
Loading…
x
Reference in New Issue
Block a user