aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Eigenvalues/EigenSolver.h
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2015-09-07 10:42:04 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2015-09-07 10:42:04 +0200
commit7031a851d45a8526474ac1ac972ad12a48e99f1a (patch)
tree8a387115aa6617e8b17862430aab4546a7c24674 /Eigen/src/Eigenvalues/EigenSolver.h
parent1702fcb72e14026af14d8af400b1a5fd4d959d19 (diff)
Generalize matrix ctor and compute() method of dense decomposition to 1) limit temporaries, 2) forward expressions to nested decompositions, 3) fix ambiguous ctor instanciation for square decomposition
Diffstat (limited to 'Eigen/src/Eigenvalues/EigenSolver.h')
-rw-r--r--Eigen/src/Eigenvalues/EigenSolver.h15
1 files changed, 9 insertions, 6 deletions
diff --git a/Eigen/src/Eigenvalues/EigenSolver.h b/Eigen/src/Eigenvalues/EigenSolver.h
index 8bc82df83..532ca7d63 100644
--- a/Eigen/src/Eigenvalues/EigenSolver.h
+++ b/Eigen/src/Eigenvalues/EigenSolver.h
@@ -110,7 +110,7 @@ template<typename _MatrixType> class EigenSolver
*
* \sa compute() for an example.
*/
- EigenSolver() : m_eivec(), m_eivalues(), m_isInitialized(false), m_realSchur(), m_matT(), m_tmp() {}
+ EigenSolver() : m_eivec(), m_eivalues(), m_isInitialized(false), m_realSchur(), m_matT(), m_tmp() {}
/** \brief Default constructor with memory preallocation
*
@@ -143,7 +143,8 @@ template<typename _MatrixType> class EigenSolver
*
* \sa compute()
*/
- explicit EigenSolver(const MatrixType& matrix, bool computeEigenvectors = true)
+ template<typename InputType>
+ explicit EigenSolver(const EigenBase<InputType>& matrix, bool computeEigenvectors = true)
: m_eivec(matrix.rows(), matrix.cols()),
m_eivalues(matrix.cols()),
m_isInitialized(false),
@@ -152,7 +153,7 @@ template<typename _MatrixType> class EigenSolver
m_matT(matrix.rows(), matrix.cols()),
m_tmp(matrix.cols())
{
- compute(matrix, computeEigenvectors);
+ compute(matrix.derived(), computeEigenvectors);
}
/** \brief Returns the eigenvectors of given matrix.
@@ -273,7 +274,8 @@ template<typename _MatrixType> class EigenSolver
* Example: \include EigenSolver_compute.cpp
* Output: \verbinclude EigenSolver_compute.out
*/
- EigenSolver& compute(const MatrixType& matrix, bool computeEigenvectors = true);
+ template<typename InputType>
+ EigenSolver& compute(const EigenBase<InputType>& matrix, bool computeEigenvectors = true);
/** \returns NumericalIssue if the input contains INF or NaN values or overflow occured. Returns Success otherwise. */
ComputationInfo info() const
@@ -370,8 +372,9 @@ typename EigenSolver<MatrixType>::EigenvectorsType EigenSolver<MatrixType>::eige
}
template<typename MatrixType>
+template<typename InputType>
EigenSolver<MatrixType>&
-EigenSolver<MatrixType>::compute(const MatrixType& matrix, bool computeEigenvectors)
+EigenSolver<MatrixType>::compute(const EigenBase<InputType>& matrix, bool computeEigenvectors)
{
check_template_parameters();
@@ -381,7 +384,7 @@ EigenSolver<MatrixType>::compute(const MatrixType& matrix, bool computeEigenvect
eigen_assert(matrix.cols() == matrix.rows());
// Reduce to real Schur form.
- m_realSchur.compute(matrix, computeEigenvectors);
+ m_realSchur.compute(matrix.derived(), computeEigenvectors);
m_info = m_realSchur.info();