aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/QR
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/QR
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/QR')
-rw-r--r--Eigen/src/QR/ColPivHouseholderQR.h38
-rw-r--r--Eigen/src/QR/FullPivHouseholderQR.h31
-rw-r--r--Eigen/src/QR/HouseholderQR.h13
3 files changed, 55 insertions, 27 deletions
diff --git a/Eigen/src/QR/ColPivHouseholderQR.h b/Eigen/src/QR/ColPivHouseholderQR.h
index 7b3842cbe..172e4a89f 100644
--- a/Eigen/src/QR/ColPivHouseholderQR.h
+++ b/Eigen/src/QR/ColPivHouseholderQR.h
@@ -118,7 +118,8 @@ template<typename _MatrixType> class ColPivHouseholderQR
*
* \sa compute()
*/
- explicit ColPivHouseholderQR(const MatrixType& matrix)
+ template<typename InputType>
+ explicit ColPivHouseholderQR(const EigenBase<InputType>& matrix)
: m_qr(matrix.rows(), matrix.cols()),
m_hCoeffs((std::min)(matrix.rows(),matrix.cols())),
m_colsPermutation(PermIndexType(matrix.cols())),
@@ -128,7 +129,7 @@ template<typename _MatrixType> class ColPivHouseholderQR
m_isInitialized(false),
m_usePrescribedThreshold(false)
{
- compute(matrix);
+ compute(matrix.derived());
}
/** This method finds a solution x to the equation Ax=b, where A is the matrix of which
@@ -185,7 +186,8 @@ template<typename _MatrixType> class ColPivHouseholderQR
return m_qr;
}
- ColPivHouseholderQR& compute(const MatrixType& matrix);
+ template<typename InputType>
+ ColPivHouseholderQR& compute(const EigenBase<InputType>& matrix);
/** \returns a const reference to the column permutation matrix */
const PermutationType& colsPermutation() const
@@ -404,6 +406,8 @@ template<typename _MatrixType> class ColPivHouseholderQR
EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
}
+ void computeInPlace();
+
MatrixType m_qr;
HCoeffsType m_hCoeffs;
PermutationType m_colsPermutation;
@@ -440,24 +444,34 @@ typename MatrixType::RealScalar ColPivHouseholderQR<MatrixType>::logAbsDetermina
* \sa class ColPivHouseholderQR, ColPivHouseholderQR(const MatrixType&)
*/
template<typename MatrixType>
-ColPivHouseholderQR<MatrixType>& ColPivHouseholderQR<MatrixType>::compute(const MatrixType& matrix)
+template<typename InputType>
+ColPivHouseholderQR<MatrixType>& ColPivHouseholderQR<MatrixType>::compute(const EigenBase<InputType>& matrix)
{
check_template_parameters();
- using std::abs;
- Index rows = matrix.rows();
- Index cols = matrix.cols();
- Index size = matrix.diagonalSize();
-
// the column permutation is stored as int indices, so just to be sure:
- eigen_assert(cols<=NumTraits<int>::highest());
+ eigen_assert(matrix.cols()<=NumTraits<int>::highest());
m_qr = matrix;
+
+ computeInPlace();
+
+ return *this;
+}
+
+template<typename MatrixType>
+void ColPivHouseholderQR<MatrixType>::computeInPlace()
+{
+ using std::abs;
+ Index rows = m_qr.rows();
+ Index cols = m_qr.cols();
+ Index size = m_qr.diagonalSize();
+
m_hCoeffs.resize(size);
m_temp.resize(cols);
- m_colsTranspositions.resize(matrix.cols());
+ m_colsTranspositions.resize(m_qr.cols());
Index number_of_transpositions = 0;
m_colSqNorms.resize(cols);
@@ -522,8 +536,6 @@ ColPivHouseholderQR<MatrixType>& ColPivHouseholderQR<MatrixType>::compute(const
m_det_pq = (number_of_transpositions%2) ? -1 : 1;
m_isInitialized = true;
-
- return *this;
}
#ifndef EIGEN_PARSED_BY_DOXYGEN
diff --git a/Eigen/src/QR/FullPivHouseholderQR.h b/Eigen/src/QR/FullPivHouseholderQR.h
index 4c2c958a8..64fe6b7b8 100644
--- a/Eigen/src/QR/FullPivHouseholderQR.h
+++ b/Eigen/src/QR/FullPivHouseholderQR.h
@@ -121,7 +121,8 @@ template<typename _MatrixType> class FullPivHouseholderQR
*
* \sa compute()
*/
- explicit FullPivHouseholderQR(const MatrixType& matrix)
+ template<typename InputType>
+ explicit FullPivHouseholderQR(const EigenBase<InputType>& matrix)
: m_qr(matrix.rows(), matrix.cols()),
m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
m_rows_transpositions((std::min)(matrix.rows(), matrix.cols())),
@@ -131,7 +132,7 @@ template<typename _MatrixType> class FullPivHouseholderQR
m_isInitialized(false),
m_usePrescribedThreshold(false)
{
- compute(matrix);
+ compute(matrix.derived());
}
/** This method finds a solution x to the equation Ax=b, where A is the matrix of which
@@ -172,7 +173,8 @@ template<typename _MatrixType> class FullPivHouseholderQR
return m_qr;
}
- FullPivHouseholderQR& compute(const MatrixType& matrix);
+ template<typename InputType>
+ FullPivHouseholderQR& compute(const EigenBase<InputType>& matrix);
/** \returns a const reference to the column permutation matrix */
const PermutationType& colsPermutation() const
@@ -386,6 +388,8 @@ template<typename _MatrixType> class FullPivHouseholderQR
EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
}
+ void computeInPlace();
+
MatrixType m_qr;
HCoeffsType m_hCoeffs;
IntDiagSizeVectorType m_rows_transpositions;
@@ -423,16 +427,27 @@ typename MatrixType::RealScalar FullPivHouseholderQR<MatrixType>::logAbsDetermin
* \sa class FullPivHouseholderQR, FullPivHouseholderQR(const MatrixType&)
*/
template<typename MatrixType>
-FullPivHouseholderQR<MatrixType>& FullPivHouseholderQR<MatrixType>::compute(const MatrixType& matrix)
+template<typename InputType>
+FullPivHouseholderQR<MatrixType>& FullPivHouseholderQR<MatrixType>::compute(const EigenBase<InputType>& matrix)
{
check_template_parameters();
+ m_qr = matrix.derived();
+
+ computeInPlace();
+
+ return *this;
+}
+
+template<typename MatrixType>
+void FullPivHouseholderQR<MatrixType>::computeInPlace()
+{
using std::abs;
- Index rows = matrix.rows();
- Index cols = matrix.cols();
+ Index rows = m_qr.rows();
+ Index cols = m_qr.cols();
Index size = (std::min)(rows,cols);
- m_qr = matrix;
+
m_hCoeffs.resize(size);
m_temp.resize(cols);
@@ -503,8 +518,6 @@ FullPivHouseholderQR<MatrixType>& FullPivHouseholderQR<MatrixType>::compute(cons
m_det_pq = (number_of_transpositions%2) ? -1 : 1;
m_isInitialized = true;
-
- return *this;
}
#ifndef EIGEN_PARSED_BY_DOXYGEN
diff --git a/Eigen/src/QR/HouseholderQR.h b/Eigen/src/QR/HouseholderQR.h
index 878654be5..1eb861025 100644
--- a/Eigen/src/QR/HouseholderQR.h
+++ b/Eigen/src/QR/HouseholderQR.h
@@ -92,13 +92,14 @@ template<typename _MatrixType> class HouseholderQR
*
* \sa compute()
*/
- explicit HouseholderQR(const MatrixType& matrix)
+ template<typename InputType>
+ explicit HouseholderQR(const EigenBase<InputType>& matrix)
: m_qr(matrix.rows(), matrix.cols()),
m_hCoeffs((std::min)(matrix.rows(),matrix.cols())),
m_temp(matrix.cols()),
m_isInitialized(false)
{
- compute(matrix);
+ compute(matrix.derived());
}
/** This method finds a solution x to the equation Ax=b, where A is the matrix of which
@@ -149,7 +150,8 @@ template<typename _MatrixType> class HouseholderQR
return m_qr;
}
- HouseholderQR& compute(const MatrixType& matrix);
+ template<typename InputType>
+ HouseholderQR& compute(const EigenBase<InputType>& matrix);
/** \returns the absolute value of the determinant of the matrix of which
* *this is the QR decomposition. It has only linear complexity
@@ -352,7 +354,8 @@ void HouseholderQR<_MatrixType>::_solve_impl(const RhsType &rhs, DstType &dst) c
* \sa class HouseholderQR, HouseholderQR(const MatrixType&)
*/
template<typename MatrixType>
-HouseholderQR<MatrixType>& HouseholderQR<MatrixType>::compute(const MatrixType& matrix)
+template<typename InputType>
+HouseholderQR<MatrixType>& HouseholderQR<MatrixType>::compute(const EigenBase<InputType>& matrix)
{
check_template_parameters();
@@ -360,7 +363,7 @@ HouseholderQR<MatrixType>& HouseholderQR<MatrixType>::compute(const MatrixType&
Index cols = matrix.cols();
Index size = (std::min)(rows,cols);
- m_qr = matrix;
+ m_qr = matrix.derived();
m_hCoeffs.resize(size);
m_temp.resize(cols);