diff options
Diffstat (limited to 'Eigen/src/Eigenvalues')
-rw-r--r-- | Eigen/src/Eigenvalues/ComplexEigenSolver.h | 4 | ||||
-rw-r--r-- | Eigen/src/Eigenvalues/ComplexSchur.h | 4 | ||||
-rw-r--r-- | Eigen/src/Eigenvalues/ComplexSchur_MKL.h | 1 | ||||
-rw-r--r-- | Eigen/src/Eigenvalues/EigenSolver.h | 10 | ||||
-rw-r--r-- | Eigen/src/Eigenvalues/GeneralizedEigenSolver.h | 4 | ||||
-rw-r--r-- | Eigen/src/Eigenvalues/GeneralizedSelfAdjointEigenSolver.h | 2 | ||||
-rw-r--r-- | Eigen/src/Eigenvalues/HessenbergDecomposition.h | 4 | ||||
-rw-r--r-- | Eigen/src/Eigenvalues/RealQZ.h | 6 | ||||
-rw-r--r-- | Eigen/src/Eigenvalues/RealSchur.h | 16 | ||||
-rw-r--r-- | Eigen/src/Eigenvalues/RealSchur_MKL.h | 4 | ||||
-rw-r--r-- | Eigen/src/Eigenvalues/SelfAdjointEigenSolver.h | 7 | ||||
-rw-r--r-- | Eigen/src/Eigenvalues/Tridiagonalization.h | 19 |
12 files changed, 34 insertions, 47 deletions
diff --git a/Eigen/src/Eigenvalues/ComplexEigenSolver.h b/Eigen/src/Eigenvalues/ComplexEigenSolver.h index af434bc9b..25082546e 100644 --- a/Eigen/src/Eigenvalues/ComplexEigenSolver.h +++ b/Eigen/src/Eigenvalues/ComplexEigenSolver.h @@ -104,7 +104,7 @@ template<typename _MatrixType> class ComplexEigenSolver * according to the specified problem \a size. * \sa ComplexEigenSolver() */ - ComplexEigenSolver(Index size) + explicit ComplexEigenSolver(Index size) : m_eivec(size, size), m_eivalues(size), m_schur(size), @@ -122,7 +122,7 @@ template<typename _MatrixType> class ComplexEigenSolver * * This constructor calls compute() to compute the eigendecomposition. */ - ComplexEigenSolver(const MatrixType& matrix, bool computeEigenvectors = true) + explicit ComplexEigenSolver(const MatrixType& matrix, bool computeEigenvectors = true) : m_eivec(matrix.rows(),matrix.cols()), m_eivalues(matrix.cols()), m_schur(matrix.rows()), diff --git a/Eigen/src/Eigenvalues/ComplexSchur.h b/Eigen/src/Eigenvalues/ComplexSchur.h index 89e6cade3..a3a5a4649 100644 --- a/Eigen/src/Eigenvalues/ComplexSchur.h +++ b/Eigen/src/Eigenvalues/ComplexSchur.h @@ -91,7 +91,7 @@ template<typename _MatrixType> class ComplexSchur * * \sa compute() for an example. */ - ComplexSchur(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime) + explicit ComplexSchur(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime) : m_matT(size,size), m_matU(size,size), m_hess(size), @@ -109,7 +109,7 @@ template<typename _MatrixType> class ComplexSchur * * \sa matrixT() and matrixU() for examples. */ - ComplexSchur(const MatrixType& matrix, bool computeU = true) + explicit ComplexSchur(const MatrixType& matrix, bool computeU = true) : m_matT(matrix.rows(),matrix.cols()), m_matU(matrix.rows(),matrix.cols()), m_hess(matrix.rows()), diff --git a/Eigen/src/Eigenvalues/ComplexSchur_MKL.h b/Eigen/src/Eigenvalues/ComplexSchur_MKL.h index 91496ae5b..27aed923c 100644 --- a/Eigen/src/Eigenvalues/ComplexSchur_MKL.h +++ b/Eigen/src/Eigenvalues/ComplexSchur_MKL.h @@ -45,7 +45,6 @@ ComplexSchur<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW> >& \ ComplexSchur<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW> >::compute(const Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW>& matrix, bool computeU) \ { \ typedef Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW> MatrixType; \ - typedef MatrixType::Scalar Scalar; \ typedef MatrixType::RealScalar RealScalar; \ typedef std::complex<RealScalar> ComplexScalar; \ \ diff --git a/Eigen/src/Eigenvalues/EigenSolver.h b/Eigen/src/Eigenvalues/EigenSolver.h index d2563d470..9372021ff 100644 --- a/Eigen/src/Eigenvalues/EigenSolver.h +++ b/Eigen/src/Eigenvalues/EigenSolver.h @@ -118,7 +118,7 @@ template<typename _MatrixType> class EigenSolver * according to the specified problem \a size. * \sa EigenSolver() */ - EigenSolver(Index size) + explicit EigenSolver(Index size) : m_eivec(size, size), m_eivalues(size), m_isInitialized(false), @@ -143,7 +143,7 @@ template<typename _MatrixType> class EigenSolver * * \sa compute() */ - EigenSolver(const MatrixType& matrix, bool computeEigenvectors = true) + explicit EigenSolver(const MatrixType& matrix, bool computeEigenvectors = true) : m_eivec(matrix.rows(), matrix.cols()), m_eivalues(matrix.cols()), m_isInitialized(false), @@ -368,7 +368,6 @@ EigenSolver<MatrixType>::compute(const MatrixType& matrix, bool computeEigenvect { using std::sqrt; using std::abs; - using std::max; using numext::isfinite; eigen_assert(matrix.cols() == matrix.rows()); @@ -409,7 +408,7 @@ EigenSolver<MatrixType>::compute(const MatrixType& matrix, bool computeEigenvect { Scalar t0 = m_matT.coeff(i+1, i); Scalar t1 = m_matT.coeff(i, i+1); - Scalar maxval = (max)(abs(p),(max)(abs(t0),abs(t1))); + Scalar maxval = numext::maxi(abs(p),numext::maxi(abs(t0),abs(t1))); t0 /= maxval; t1 /= maxval; Scalar p0 = p/maxval; @@ -600,8 +599,7 @@ void EigenSolver<MatrixType>::doComputeEigenvectors() } // Overflow control - EIGEN_USING_STD_MATH(max); - Scalar t = (max)(abs(m_matT.coeff(i,n-1)),abs(m_matT.coeff(i,n))); + Scalar t = numext::maxi(abs(m_matT.coeff(i,n-1)),abs(m_matT.coeff(i,n))); if ((eps * t) * t > Scalar(1)) m_matT.block(i, n-1, size-i, 2) /= t; diff --git a/Eigen/src/Eigenvalues/GeneralizedEigenSolver.h b/Eigen/src/Eigenvalues/GeneralizedEigenSolver.h index dc240e13e..c20ea03e6 100644 --- a/Eigen/src/Eigenvalues/GeneralizedEigenSolver.h +++ b/Eigen/src/Eigenvalues/GeneralizedEigenSolver.h @@ -122,7 +122,7 @@ template<typename _MatrixType> class GeneralizedEigenSolver * according to the specified problem \a size. * \sa GeneralizedEigenSolver() */ - GeneralizedEigenSolver(Index size) + explicit GeneralizedEigenSolver(Index size) : m_eivec(size, size), m_alphas(size), m_betas(size), @@ -145,7 +145,7 @@ template<typename _MatrixType> class GeneralizedEigenSolver * * \sa compute() */ - GeneralizedEigenSolver(const MatrixType& A, const MatrixType& B, bool computeEigenvectors = true) + explicit GeneralizedEigenSolver(const MatrixType& A, const MatrixType& B, bool computeEigenvectors = true) : m_eivec(A.rows(), A.cols()), m_alphas(A.cols()), m_betas(A.cols()), diff --git a/Eigen/src/Eigenvalues/GeneralizedSelfAdjointEigenSolver.h b/Eigen/src/Eigenvalues/GeneralizedSelfAdjointEigenSolver.h index 07bf1ea09..1ce1f5f58 100644 --- a/Eigen/src/Eigenvalues/GeneralizedSelfAdjointEigenSolver.h +++ b/Eigen/src/Eigenvalues/GeneralizedSelfAdjointEigenSolver.h @@ -74,7 +74,7 @@ class GeneralizedSelfAdjointEigenSolver : public SelfAdjointEigenSolver<_MatrixT * * \sa compute() for an example */ - GeneralizedSelfAdjointEigenSolver(Index size) + explicit GeneralizedSelfAdjointEigenSolver(Index size) : Base(size) {} diff --git a/Eigen/src/Eigenvalues/HessenbergDecomposition.h b/Eigen/src/Eigenvalues/HessenbergDecomposition.h index 3db0c0106..2615a9f23 100644 --- a/Eigen/src/Eigenvalues/HessenbergDecomposition.h +++ b/Eigen/src/Eigenvalues/HessenbergDecomposition.h @@ -97,7 +97,7 @@ template<typename _MatrixType> class HessenbergDecomposition * * \sa compute() for an example. */ - HessenbergDecomposition(Index size = Size==Dynamic ? 2 : Size) + explicit HessenbergDecomposition(Index size = Size==Dynamic ? 2 : Size) : m_matrix(size,size), m_temp(size), m_isInitialized(false) @@ -115,7 +115,7 @@ template<typename _MatrixType> class HessenbergDecomposition * * \sa matrixH() for an example. */ - HessenbergDecomposition(const MatrixType& matrix) + explicit HessenbergDecomposition(const MatrixType& matrix) : m_matrix(matrix), m_temp(matrix.rows()), m_isInitialized(false) diff --git a/Eigen/src/Eigenvalues/RealQZ.h b/Eigen/src/Eigenvalues/RealQZ.h index 5706eeebe..128ef9028 100644 --- a/Eigen/src/Eigenvalues/RealQZ.h +++ b/Eigen/src/Eigenvalues/RealQZ.h @@ -83,7 +83,7 @@ namespace Eigen { * * \sa compute() for an example. */ - RealQZ(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime) : + explicit RealQZ(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime) : m_S(size, size), m_T(size, size), m_Q(size, size), @@ -101,7 +101,7 @@ namespace Eigen { * * This constructor calls compute() to compute the QZ decomposition. */ - RealQZ(const MatrixType& A, const MatrixType& B, bool computeQZ = true) : + explicit RealQZ(const MatrixType& A, const MatrixType& B, bool computeQZ = true) : m_S(A.rows(),A.cols()), m_T(A.rows(),A.cols()), m_Q(A.rows(),A.cols()), @@ -313,7 +313,7 @@ namespace Eigen { using std::abs; using std::sqrt; const Index dim=m_S.cols(); - if (abs(m_S.coeff(i+1,i)==Scalar(0))) + if (abs(m_S.coeff(i+1,i))==Scalar(0)) return; Index z = findSmallDiagEntry(i,i+1); if (z==i-1) diff --git a/Eigen/src/Eigenvalues/RealSchur.h b/Eigen/src/Eigenvalues/RealSchur.h index 64d136341..51e61ba38 100644 --- a/Eigen/src/Eigenvalues/RealSchur.h +++ b/Eigen/src/Eigenvalues/RealSchur.h @@ -80,7 +80,7 @@ template<typename _MatrixType> class RealSchur * * \sa compute() for an example. */ - RealSchur(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime) + explicit RealSchur(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime) : m_matT(size, size), m_matU(size, size), m_workspaceVector(size), @@ -100,7 +100,7 @@ template<typename _MatrixType> class RealSchur * Example: \include RealSchur_RealSchur_MatrixType.cpp * Output: \verbinclude RealSchur_RealSchur_MatrixType.out */ - RealSchur(const MatrixType& matrix, bool computeU = true) + explicit RealSchur(const MatrixType& matrix, bool computeU = true) : m_matT(matrix.rows(),matrix.cols()), m_matU(matrix.rows(),matrix.cols()), m_workspaceVector(matrix.rows()), @@ -234,7 +234,7 @@ template<typename _MatrixType> class RealSchur typedef Matrix<Scalar,3,1> Vector3s; Scalar computeNormOfT(); - Index findSmallSubdiagEntry(Index iu, const Scalar& norm); + Index findSmallSubdiagEntry(Index iu); void splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift); void computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo); void initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector); @@ -286,7 +286,7 @@ RealSchur<MatrixType>& RealSchur<MatrixType>::computeFromHessenberg(const HessMa { while (iu >= 0) { - Index il = findSmallSubdiagEntry(iu, norm); + Index il = findSmallSubdiagEntry(iu); // Check for convergence if (il == iu) // One root found @@ -343,16 +343,14 @@ inline typename MatrixType::Scalar RealSchur<MatrixType>::computeNormOfT() /** \internal Look for single small sub-diagonal element and returns its index */ template<typename MatrixType> -inline typename MatrixType::Index RealSchur<MatrixType>::findSmallSubdiagEntry(Index iu, const Scalar& norm) +inline typename MatrixType::Index RealSchur<MatrixType>::findSmallSubdiagEntry(Index iu) { using std::abs; Index res = iu; while (res > 0) { Scalar s = abs(m_matT.coeff(res-1,res-1)) + abs(m_matT.coeff(res,res)); - if (s == 0.0) - s = norm; - if (abs(m_matT.coeff(res,res-1)) < NumTraits<Scalar>::epsilon() * s) + if (abs(m_matT.coeff(res,res-1)) <= NumTraits<Scalar>::epsilon() * s) break; res--; } @@ -457,9 +455,7 @@ inline void RealSchur<MatrixType>::initFrancisQRStep(Index il, Index iu, const V const Scalar lhs = m_matT.coeff(im,im-1) * (abs(v.coeff(1)) + abs(v.coeff(2))); const Scalar rhs = v.coeff(0) * (abs(m_matT.coeff(im-1,im-1)) + abs(Tmm) + abs(m_matT.coeff(im+1,im+1))); if (abs(lhs) < NumTraits<Scalar>::epsilon() * rhs) - { break; - } } } diff --git a/Eigen/src/Eigenvalues/RealSchur_MKL.h b/Eigen/src/Eigenvalues/RealSchur_MKL.h index ad9736460..c3089b468 100644 --- a/Eigen/src/Eigenvalues/RealSchur_MKL.h +++ b/Eigen/src/Eigenvalues/RealSchur_MKL.h @@ -44,10 +44,6 @@ template<> inline \ RealSchur<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW> >& \ RealSchur<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW> >::compute(const Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW>& matrix, bool computeU) \ { \ - typedef Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW> MatrixType; \ - typedef MatrixType::Scalar Scalar; \ - typedef MatrixType::RealScalar RealScalar; \ -\ eigen_assert(matrix.cols() == matrix.rows()); \ \ lapack_int n = matrix.cols(), sdim, info; \ diff --git a/Eigen/src/Eigenvalues/SelfAdjointEigenSolver.h b/Eigen/src/Eigenvalues/SelfAdjointEigenSolver.h index a6bbdac6b..54f60b197 100644 --- a/Eigen/src/Eigenvalues/SelfAdjointEigenSolver.h +++ b/Eigen/src/Eigenvalues/SelfAdjointEigenSolver.h @@ -133,7 +133,7 @@ template<typename _MatrixType> class SelfAdjointEigenSolver * \sa compute() for an example */ EIGEN_DEVICE_FUNC - SelfAdjointEigenSolver(Index size) + explicit SelfAdjointEigenSolver(Index size) : m_eivec(size, size), m_eivalues(size), m_subdiag(size > 1 ? size - 1 : 1), @@ -156,7 +156,7 @@ template<typename _MatrixType> class SelfAdjointEigenSolver * \sa compute(const MatrixType&, int) */ EIGEN_DEVICE_FUNC - SelfAdjointEigenSolver(const MatrixType& matrix, int options = ComputeEigenvectors) + explicit SelfAdjointEigenSolver(const MatrixType& matrix, int options = ComputeEigenvectors) : m_eivec(matrix.rows(), matrix.cols()), m_eivalues(matrix.cols()), m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1), @@ -732,7 +732,6 @@ struct direct_selfadjoint_eigenvalues<SolverType,2,false> EIGEN_DEVICE_FUNC static inline void run(SolverType& solver, const MatrixType& mat, int options) { - EIGEN_USING_STD_MATH(max) EIGEN_USING_STD_MATH(sqrt); eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows()); @@ -746,7 +745,7 @@ struct direct_selfadjoint_eigenvalues<SolverType,2,false> // map the matrix coefficients to [-1:1] to avoid over- and underflow. Scalar scale = mat.cwiseAbs().maxCoeff(); - scale = (max)(scale,Scalar(1)); + scale = numext::maxi(scale,Scalar(1)); MatrixType scaledMat = mat / scale; // Compute the eigenvalues diff --git a/Eigen/src/Eigenvalues/Tridiagonalization.h b/Eigen/src/Eigenvalues/Tridiagonalization.h index 192278d68..bedd1cb34 100644 --- a/Eigen/src/Eigenvalues/Tridiagonalization.h +++ b/Eigen/src/Eigenvalues/Tridiagonalization.h @@ -18,8 +18,10 @@ namespace internal { template<typename MatrixType> struct TridiagonalizationMatrixTReturnType; template<typename MatrixType> struct traits<TridiagonalizationMatrixTReturnType<MatrixType> > + : public traits<typename MatrixType::PlainObject> { - typedef typename MatrixType::PlainObject ReturnType; + typedef typename MatrixType::PlainObject ReturnType; // FIXME shall it be a BandMatrix? + enum { Flags = 0 }; }; template<typename MatrixType, typename CoeffVectorType> @@ -89,10 +91,8 @@ template<typename _MatrixType> class Tridiagonalization >::type DiagonalReturnType; typedef typename internal::conditional<NumTraits<Scalar>::IsComplex, - typename internal::add_const_on_value_type<typename Diagonal< - Block<const MatrixType,SizeMinusOne,SizeMinusOne> >::RealReturnType>::type, - const Diagonal< - Block<const MatrixType,SizeMinusOne,SizeMinusOne> > + typename internal::add_const_on_value_type<typename Diagonal<const MatrixType, -1>::RealReturnType>::type, + const Diagonal<const MatrixType, -1> >::type SubDiagonalReturnType; /** \brief Return type of matrixQ() */ @@ -110,7 +110,7 @@ template<typename _MatrixType> class Tridiagonalization * * \sa compute() for an example. */ - Tridiagonalization(Index size = Size==Dynamic ? 2 : Size) + explicit Tridiagonalization(Index size = Size==Dynamic ? 2 : Size) : m_matrix(size,size), m_hCoeffs(size > 1 ? size-1 : 1), m_isInitialized(false) @@ -126,7 +126,7 @@ template<typename _MatrixType> class Tridiagonalization * Example: \include Tridiagonalization_Tridiagonalization_MatrixType.cpp * Output: \verbinclude Tridiagonalization_Tridiagonalization_MatrixType.out */ - Tridiagonalization(const MatrixType& matrix) + explicit Tridiagonalization(const MatrixType& matrix) : m_matrix(matrix), m_hCoeffs(matrix.cols() > 1 ? matrix.cols()-1 : 1), m_isInitialized(false) @@ -305,7 +305,7 @@ typename Tridiagonalization<MatrixType>::DiagonalReturnType Tridiagonalization<MatrixType>::diagonal() const { eigen_assert(m_isInitialized && "Tridiagonalization is not initialized."); - return m_matrix.diagonal(); + return m_matrix.diagonal().real(); } template<typename MatrixType> @@ -313,8 +313,7 @@ typename Tridiagonalization<MatrixType>::SubDiagonalReturnType Tridiagonalization<MatrixType>::subDiagonal() const { eigen_assert(m_isInitialized && "Tridiagonalization is not initialized."); - Index n = m_matrix.rows(); - return Block<const MatrixType,SizeMinusOne,SizeMinusOne>(m_matrix, 1, 0, n-1,n-1).diagonal(); + return m_matrix.template diagonal<-1>().real(); } namespace internal { |