mirror of
https://gitlab.com/libeigen/eigen.git
synced 2024-11-21 03:11:25 +08:00
add ColPivotingHouseholderQR
rename RRQR to fullPivotingHouseholderQR
This commit is contained in:
parent
a848ed02ad
commit
97bc1af1f1
3
Eigen/QR
3
Eigen/QR
@ -36,7 +36,8 @@ namespace Eigen {
|
||||
*/
|
||||
|
||||
#include "src/QR/QR.h"
|
||||
#include "src/QR/RRQR.h"
|
||||
#include "src/QR/FullPivotingHouseholderQR.h"
|
||||
#include "src/QR/ColPivotingHouseholderQR.h"
|
||||
#include "src/QR/Tridiagonalization.h"
|
||||
#include "src/QR/EigenSolver.h"
|
||||
#include "src/QR/SelfAdjointEigenSolver.h"
|
||||
|
@ -745,6 +745,8 @@ template<typename Derived> class MatrixBase
|
||||
/////////// QR module ///////////
|
||||
|
||||
const HouseholderQR<PlainMatrixType> householderQr() const;
|
||||
const ColPivotingHouseholderQR<PlainMatrixType> colPivotingHouseholderQr() const;
|
||||
const FullPivotingHouseholderQR<PlainMatrixType> fullPivotingHouseholderQr() const;
|
||||
|
||||
EigenvaluesReturnType eigenvalues() const;
|
||||
RealScalar operatorNorm() const;
|
||||
|
@ -117,6 +117,8 @@ template<typename MatrixType, int Direction = BothDirections> class Reverse;
|
||||
template<typename MatrixType> class LU;
|
||||
template<typename MatrixType> class PartialLU;
|
||||
template<typename MatrixType> class HouseholderQR;
|
||||
template<typename MatrixType> class ColPivotingHouseholderQR;
|
||||
template<typename MatrixType> class FullPivotingHouseholderQR;
|
||||
template<typename MatrixType> class SVD;
|
||||
template<typename MatrixType, int UpLo = LowerTriangular> class LLT;
|
||||
template<typename MatrixType> class LDLT;
|
||||
|
266
Eigen/src/QR/ColPivotingHouseholderQR.h
Normal file
266
Eigen/src/QR/ColPivotingHouseholderQR.h
Normal file
@ -0,0 +1,266 @@
|
||||
// This file is part of Eigen, a lightweight C++ template library
|
||||
// for linear algebra.
|
||||
//
|
||||
// Copyright (C) 2008-2009 Gael Guennebaud <g.gael@free.fr>
|
||||
// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
|
||||
//
|
||||
// Eigen is free software; you can redistribute it and/or
|
||||
// modify it under the terms of the GNU Lesser General Public
|
||||
// License as published by the Free Software Foundation; either
|
||||
// version 3 of the License, or (at your option) any later version.
|
||||
//
|
||||
// Alternatively, you can redistribute it and/or
|
||||
// modify it under the terms of the GNU General Public License as
|
||||
// published by the Free Software Foundation; either version 2 of
|
||||
// the License, or (at your option) any later version.
|
||||
//
|
||||
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
|
||||
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
|
||||
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
|
||||
// GNU General Public License for more details.
|
||||
//
|
||||
// You should have received a copy of the GNU Lesser General Public
|
||||
// License and a copy of the GNU General Public License along with
|
||||
// Eigen. If not, see <http://www.gnu.org/licenses/>.
|
||||
|
||||
#ifndef EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
|
||||
#define EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
|
||||
|
||||
/** \ingroup QR_Module
|
||||
* \nonstableyet
|
||||
*
|
||||
* \class ColPivotingHouseholderQR
|
||||
*
|
||||
* \brief Householder rank-revealing QR decomposition of a matrix
|
||||
*
|
||||
* \param MatrixType the type of the matrix of which we are computing the QR decomposition
|
||||
*
|
||||
* This class performs a rank-revealing QR decomposition using Householder transformations.
|
||||
*
|
||||
* This decomposition performs full-pivoting in order to be rank-revealing and achieve optimal
|
||||
* numerical stability.
|
||||
*
|
||||
* \sa MatrixBase::colPivotingHouseholderQr()
|
||||
*/
|
||||
template<typename MatrixType> class ColPivotingHouseholderQR
|
||||
{
|
||||
public:
|
||||
|
||||
enum {
|
||||
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
|
||||
Options = MatrixType::Options,
|
||||
DiagSizeAtCompileTime = EIGEN_ENUM_MIN(ColsAtCompileTime,RowsAtCompileTime)
|
||||
};
|
||||
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef typename MatrixType::RealScalar RealScalar;
|
||||
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime> MatrixQType;
|
||||
typedef Matrix<Scalar, DiagSizeAtCompileTime, 1> HCoeffsType;
|
||||
typedef Matrix<int, 1, ColsAtCompileTime> IntRowVectorType;
|
||||
typedef Matrix<int, RowsAtCompileTime, 1> IntColVectorType;
|
||||
typedef Matrix<Scalar, 1, ColsAtCompileTime> RowVectorType;
|
||||
typedef Matrix<Scalar, RowsAtCompileTime, 1> ColVectorType;
|
||||
typedef Matrix<RealScalar, 1, ColsAtCompileTime> RealRowVectorType;
|
||||
|
||||
/**
|
||||
* \brief Default Constructor.
|
||||
*
|
||||
* The default constructor is useful in cases in which the user intends to
|
||||
* perform decompositions via ColPivotingHouseholderQR::compute(const MatrixType&).
|
||||
*/
|
||||
ColPivotingHouseholderQR() : m_qr(), m_hCoeffs(), m_isInitialized(false) {}
|
||||
|
||||
ColPivotingHouseholderQR(const MatrixType& matrix)
|
||||
: m_qr(matrix.rows(), matrix.cols()),
|
||||
m_hCoeffs(std::min(matrix.rows(),matrix.cols())),
|
||||
m_isInitialized(false)
|
||||
{
|
||||
compute(matrix);
|
||||
}
|
||||
|
||||
/** This method finds a solution x to the equation Ax=b, where A is the matrix of which
|
||||
* *this is the QR decomposition, if any exists.
|
||||
*
|
||||
* \param b the right-hand-side of the equation to solve.
|
||||
*
|
||||
* \param result a pointer to the vector/matrix in which to store the solution, if any exists.
|
||||
* Resized if necessary, so that result->rows()==A.cols() and result->cols()==b.cols().
|
||||
* If no solution exists, *result is left with undefined coefficients.
|
||||
*
|
||||
* \note The case where b is a matrix is not yet implemented. Also, this
|
||||
* code is space inefficient.
|
||||
*
|
||||
* Example: \include ColPivotingHouseholderQR_solve.cpp
|
||||
* Output: \verbinclude ColPivotingHouseholderQR_solve.out
|
||||
*/
|
||||
template<typename OtherDerived, typename ResultType>
|
||||
void solve(const MatrixBase<OtherDerived>& b, ResultType *result) const;
|
||||
|
||||
MatrixType matrixQ(void) const;
|
||||
|
||||
/** \returns a reference to the matrix where the Householder QR decomposition is stored
|
||||
*/
|
||||
const MatrixType& matrixQR() const { return m_qr; }
|
||||
|
||||
ColPivotingHouseholderQR& compute(const MatrixType& matrix);
|
||||
|
||||
const IntRowVectorType& colsPermutation() const
|
||||
{
|
||||
ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
|
||||
return m_cols_permutation;
|
||||
}
|
||||
|
||||
inline int rank() const
|
||||
{
|
||||
ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
|
||||
return m_rank;
|
||||
}
|
||||
|
||||
protected:
|
||||
MatrixType m_qr;
|
||||
HCoeffsType m_hCoeffs;
|
||||
IntRowVectorType m_cols_permutation;
|
||||
bool m_isInitialized;
|
||||
RealScalar m_precision;
|
||||
int m_rank;
|
||||
int m_det_pq;
|
||||
};
|
||||
|
||||
#ifndef EIGEN_HIDE_HEAVY_CODE
|
||||
|
||||
template<typename MatrixType>
|
||||
ColPivotingHouseholderQR<MatrixType>& ColPivotingHouseholderQR<MatrixType>::compute(const MatrixType& matrix)
|
||||
{
|
||||
int rows = matrix.rows();
|
||||
int cols = matrix.cols();
|
||||
int size = std::min(rows,cols);
|
||||
m_rank = size;
|
||||
|
||||
m_qr = matrix;
|
||||
m_hCoeffs.resize(size);
|
||||
|
||||
RowVectorType temp(cols);
|
||||
|
||||
m_precision = epsilon<Scalar>() * size;
|
||||
|
||||
IntRowVectorType cols_transpositions(matrix.cols());
|
||||
m_cols_permutation.resize(matrix.cols());
|
||||
int number_of_transpositions = 0;
|
||||
|
||||
RealRowVectorType colSqNorms(cols);
|
||||
for(int k = 0; k < cols; ++k)
|
||||
colSqNorms.coeffRef(k) = m_qr.col(k).squaredNorm();
|
||||
RealScalar biggestColSqNorm = colSqNorms.maxCoeff();
|
||||
|
||||
for (int k = 0; k < size; ++k)
|
||||
{
|
||||
int biggest_col_in_corner;
|
||||
RealScalar biggestColSqNormInCorner = colSqNorms.end(cols-k).maxCoeff(&biggest_col_in_corner);
|
||||
biggest_col_in_corner += k;
|
||||
|
||||
// if the corner is negligible, then we have less than full rank, and we can finish early
|
||||
if(ei_isMuchSmallerThan(biggestColSqNormInCorner, biggestColSqNorm, m_precision))
|
||||
{
|
||||
m_rank = k;
|
||||
for(int i = k; i < size; i++)
|
||||
{
|
||||
cols_transpositions.coeffRef(i) = i;
|
||||
m_hCoeffs.coeffRef(i) = Scalar(0);
|
||||
}
|
||||
break;
|
||||
}
|
||||
|
||||
cols_transpositions.coeffRef(k) = biggest_col_in_corner;
|
||||
if(k != biggest_col_in_corner) {
|
||||
m_qr.col(k).swap(m_qr.col(biggest_col_in_corner));
|
||||
++number_of_transpositions;
|
||||
}
|
||||
|
||||
RealScalar beta;
|
||||
m_qr.col(k).end(rows-k).makeHouseholderInPlace(&m_hCoeffs.coeffRef(k), &beta);
|
||||
m_qr.coeffRef(k,k) = beta;
|
||||
|
||||
m_qr.corner(BottomRight, rows-k, cols-k-1)
|
||||
.applyHouseholderOnTheLeft(m_qr.col(k).end(rows-k-1), m_hCoeffs.coeffRef(k), &temp.coeffRef(k+1));
|
||||
|
||||
colSqNorms.end(cols-k-1) -= m_qr.row(k).end(cols-k-1).cwise().abs2();
|
||||
}
|
||||
|
||||
for(int k = 0; k < matrix.cols(); ++k) m_cols_permutation.coeffRef(k) = k;
|
||||
for(int k = 0; k < size; ++k)
|
||||
std::swap(m_cols_permutation.coeffRef(k), m_cols_permutation.coeffRef(cols_transpositions.coeff(k)));
|
||||
|
||||
m_det_pq = (number_of_transpositions%2) ? -1 : 1;
|
||||
m_isInitialized = true;
|
||||
|
||||
return *this;
|
||||
}
|
||||
|
||||
template<typename MatrixType>
|
||||
template<typename OtherDerived, typename ResultType>
|
||||
void ColPivotingHouseholderQR<MatrixType>::solve(
|
||||
const MatrixBase<OtherDerived>& b,
|
||||
ResultType *result
|
||||
) const
|
||||
{
|
||||
ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
|
||||
const int rows = m_qr.rows();
|
||||
const int cols = b.cols();
|
||||
ei_assert(b.rows() == rows);
|
||||
|
||||
typename OtherDerived::PlainMatrixType c(b);
|
||||
|
||||
Matrix<Scalar,1,MatrixType::ColsAtCompileTime> temp(cols);
|
||||
for (int k = 0; k < m_rank; ++k)
|
||||
{
|
||||
int remainingSize = rows-k;
|
||||
c.corner(BottomRight, remainingSize, cols)
|
||||
.applyHouseholderOnTheLeft(m_qr.col(k).end(remainingSize-1), m_hCoeffs.coeff(k), &temp.coeffRef(0));
|
||||
}
|
||||
|
||||
m_qr.corner(TopLeft, m_rank, m_rank)
|
||||
.template triangularView<UpperTriangular>()
|
||||
.solveInPlace(c.corner(TopLeft, m_rank, c.cols()));
|
||||
|
||||
result->resize(m_qr.cols(), b.cols());
|
||||
for(int i = 0; i < m_rank; ++i) result->row(m_cols_permutation.coeff(i)) = c.row(i);
|
||||
for(int i = m_rank; i < m_qr.cols(); ++i) result->row(m_cols_permutation.coeff(i)).setZero();
|
||||
}
|
||||
|
||||
/** \returns the matrix Q */
|
||||
template<typename MatrixType>
|
||||
MatrixType ColPivotingHouseholderQR<MatrixType>::matrixQ() const
|
||||
{
|
||||
ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
|
||||
// compute the product H'_0 H'_1 ... H'_n-1,
|
||||
// where H_k is the k-th Householder transformation I - h_k v_k v_k'
|
||||
// and v_k is the k-th Householder vector [1,m_qr(k+1,k), m_qr(k+2,k), ...]
|
||||
int rows = m_qr.rows();
|
||||
int cols = m_qr.cols();
|
||||
int size = std::min(rows,cols);
|
||||
MatrixType res = MatrixType::Identity(rows, rows);
|
||||
Matrix<Scalar,1,MatrixType::RowsAtCompileTime> temp(rows);
|
||||
for (int k = size-1; k >= 0; k--)
|
||||
{
|
||||
res.block(k, k, rows-k, rows-k)
|
||||
.applyHouseholderOnTheLeft(m_qr.col(k).end(rows-k-1), ei_conj(m_hCoeffs.coeff(k)), &temp.coeffRef(k));
|
||||
}
|
||||
return res;
|
||||
}
|
||||
|
||||
#endif // EIGEN_HIDE_HEAVY_CODE
|
||||
|
||||
/** \return the column-pivoting Householder QR decomposition of \c *this.
|
||||
*
|
||||
* \sa class ColPivotingHouseholderQR
|
||||
*/
|
||||
template<typename Derived>
|
||||
const ColPivotingHouseholderQR<typename MatrixBase<Derived>::PlainMatrixType>
|
||||
MatrixBase<Derived>::colPivotingHouseholderQr() const
|
||||
{
|
||||
return ColPivotingHouseholderQR<PlainMatrixType>(eval());
|
||||
}
|
||||
|
||||
|
||||
#endif // EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
|
@ -23,13 +23,13 @@
|
||||
// License and a copy of the GNU General Public License along with
|
||||
// Eigen. If not, see <http://www.gnu.org/licenses/>.
|
||||
|
||||
#ifndef EIGEN_RRQR_H
|
||||
#define EIGEN_RRQR_H
|
||||
#ifndef EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
|
||||
#define EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
|
||||
|
||||
/** \ingroup QR_Module
|
||||
* \nonstableyet
|
||||
*
|
||||
* \class HouseholderRRQR
|
||||
* \class FullPivotingHouseholderQR
|
||||
*
|
||||
* \brief Householder rank-revealing QR decomposition of a matrix
|
||||
*
|
||||
@ -42,7 +42,7 @@
|
||||
*
|
||||
* \sa MatrixBase::householderRrqr()
|
||||
*/
|
||||
template<typename MatrixType> class HouseholderRRQR
|
||||
template<typename MatrixType> class FullPivotingHouseholderQR
|
||||
{
|
||||
public:
|
||||
|
||||
@ -66,11 +66,11 @@ template<typename MatrixType> class HouseholderRRQR
|
||||
* \brief Default Constructor.
|
||||
*
|
||||
* The default constructor is useful in cases in which the user intends to
|
||||
* perform decompositions via HouseholderRRQR::compute(const MatrixType&).
|
||||
* perform decompositions via FullPivotingHouseholderQR::compute(const MatrixType&).
|
||||
*/
|
||||
HouseholderRRQR() : m_qr(), m_hCoeffs(), m_isInitialized(false) {}
|
||||
FullPivotingHouseholderQR() : m_qr(), m_hCoeffs(), m_isInitialized(false) {}
|
||||
|
||||
HouseholderRRQR(const MatrixType& matrix)
|
||||
FullPivotingHouseholderQR(const MatrixType& matrix)
|
||||
: m_qr(matrix.rows(), matrix.cols()),
|
||||
m_hCoeffs(std::min(matrix.rows(),matrix.cols())),
|
||||
m_isInitialized(false)
|
||||
@ -90,8 +90,8 @@ template<typename MatrixType> class HouseholderRRQR
|
||||
* \note The case where b is a matrix is not yet implemented. Also, this
|
||||
* code is space inefficient.
|
||||
*
|
||||
* Example: \include HouseholderRRQR_solve.cpp
|
||||
* Output: \verbinclude HouseholderRRQR_solve.out
|
||||
* Example: \include FullPivotingHouseholderQR_solve.cpp
|
||||
* Output: \verbinclude FullPivotingHouseholderQR_solve.out
|
||||
*/
|
||||
template<typename OtherDerived, typename ResultType>
|
||||
void solve(const MatrixBase<OtherDerived>& b, ResultType *result) const;
|
||||
@ -102,23 +102,23 @@ template<typename MatrixType> class HouseholderRRQR
|
||||
*/
|
||||
const MatrixType& matrixQR() const { return m_qr; }
|
||||
|
||||
HouseholderRRQR& compute(const MatrixType& matrix);
|
||||
FullPivotingHouseholderQR& compute(const MatrixType& matrix);
|
||||
|
||||
const IntRowVectorType& colsPermutation() const
|
||||
{
|
||||
ei_assert(m_isInitialized && "RRQR is not initialized.");
|
||||
ei_assert(m_isInitialized && "FULLPIVOTINGHOUSEHOLDERQR is not initialized.");
|
||||
return m_cols_permutation;
|
||||
}
|
||||
|
||||
const IntColVectorType& rowsTranspositions() const
|
||||
{
|
||||
ei_assert(m_isInitialized && "RRQR is not initialized.");
|
||||
ei_assert(m_isInitialized && "FULLPIVOTINGHOUSEHOLDERQR is not initialized.");
|
||||
return m_rows_transpositions;
|
||||
}
|
||||
|
||||
inline int rank() const
|
||||
{
|
||||
ei_assert(m_isInitialized && "RRQR is not initialized.");
|
||||
ei_assert(m_isInitialized && "FULLPIVOTINGHOUSEHOLDERQR is not initialized.");
|
||||
return m_rank;
|
||||
}
|
||||
|
||||
@ -136,7 +136,7 @@ template<typename MatrixType> class HouseholderRRQR
|
||||
#ifndef EIGEN_HIDE_HEAVY_CODE
|
||||
|
||||
template<typename MatrixType>
|
||||
HouseholderRRQR<MatrixType>& HouseholderRRQR<MatrixType>::compute(const MatrixType& matrix)
|
||||
FullPivotingHouseholderQR<MatrixType>& FullPivotingHouseholderQR<MatrixType>::compute(const MatrixType& matrix)
|
||||
{
|
||||
int rows = matrix.rows();
|
||||
int cols = matrix.cols();
|
||||
@ -148,7 +148,6 @@ HouseholderRRQR<MatrixType>& HouseholderRRQR<MatrixType>::compute(const MatrixTy
|
||||
|
||||
RowVectorType temp(cols);
|
||||
|
||||
// TODO: experiment to see the best formula
|
||||
m_precision = epsilon<Scalar>() * size;
|
||||
|
||||
m_rows_transpositions.resize(matrix.rows());
|
||||
@ -198,7 +197,6 @@ HouseholderRRQR<MatrixType>& HouseholderRRQR<MatrixType>::compute(const MatrixTy
|
||||
m_qr.col(k).end(rows-k).makeHouseholderInPlace(&m_hCoeffs.coeffRef(k), &beta);
|
||||
m_qr.coeffRef(k,k) = beta;
|
||||
|
||||
// apply H to remaining part of m_qr from the left
|
||||
m_qr.corner(BottomRight, rows-k, cols-k-1)
|
||||
.applyHouseholderOnTheLeft(m_qr.col(k).end(rows-k-1), m_hCoeffs.coeffRef(k), &temp.coeffRef(k+1));
|
||||
}
|
||||
@ -215,12 +213,12 @@ HouseholderRRQR<MatrixType>& HouseholderRRQR<MatrixType>::compute(const MatrixTy
|
||||
|
||||
template<typename MatrixType>
|
||||
template<typename OtherDerived, typename ResultType>
|
||||
void HouseholderRRQR<MatrixType>::solve(
|
||||
void FullPivotingHouseholderQR<MatrixType>::solve(
|
||||
const MatrixBase<OtherDerived>& b,
|
||||
ResultType *result
|
||||
) const
|
||||
{
|
||||
ei_assert(m_isInitialized && "HouseholderRRQR is not initialized.");
|
||||
ei_assert(m_isInitialized && "FullPivotingHouseholderQR is not initialized.");
|
||||
const int rows = m_qr.rows();
|
||||
const int cols = b.cols();
|
||||
ei_assert(b.rows() == rows);
|
||||
@ -247,9 +245,9 @@ void HouseholderRRQR<MatrixType>::solve(
|
||||
|
||||
/** \returns the matrix Q */
|
||||
template<typename MatrixType>
|
||||
MatrixType HouseholderRRQR<MatrixType>::matrixQ() const
|
||||
MatrixType FullPivotingHouseholderQR<MatrixType>::matrixQ() const
|
||||
{
|
||||
ei_assert(m_isInitialized && "HouseholderRRQR is not initialized.");
|
||||
ei_assert(m_isInitialized && "FullPivotingHouseholderQR is not initialized.");
|
||||
// compute the product H'_0 H'_1 ... H'_n-1,
|
||||
// where H_k is the k-th Householder transformation I - h_k v_k v_k'
|
||||
// and v_k is the k-th Householder vector [1,m_qr(k+1,k), m_qr(k+2,k), ...]
|
||||
@ -269,18 +267,15 @@ MatrixType HouseholderRRQR<MatrixType>::matrixQ() const
|
||||
|
||||
#endif // EIGEN_HIDE_HEAVY_CODE
|
||||
|
||||
#if 0
|
||||
/** \return the Householder QR decomposition of \c *this.
|
||||
/** \return the full-pivoting Householder QR decomposition of \c *this.
|
||||
*
|
||||
* \sa class HouseholderRRQR
|
||||
* \sa class FullPivotingHouseholderQR
|
||||
*/
|
||||
template<typename Derived>
|
||||
const HouseholderRRQR<typename MatrixBase<Derived>::PlainMatrixType>
|
||||
MatrixBase<Derived>::householderQr() const
|
||||
const FullPivotingHouseholderQR<typename MatrixBase<Derived>::PlainMatrixType>
|
||||
MatrixBase<Derived>::fullPivotingHouseholderQr() const
|
||||
{
|
||||
return HouseholderRRQR<PlainMatrixType>(eval());
|
||||
return FullPivotingHouseholderQR<PlainMatrixType>(eval());
|
||||
}
|
||||
#endif
|
||||
|
||||
|
||||
#endif // EIGEN_QR_H
|
||||
#endif // EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
|
@ -121,7 +121,8 @@ ei_add_test(lu ${EI_OFLAG})
|
||||
ei_add_test(determinant)
|
||||
ei_add_test(inverse ${EI_OFLAG})
|
||||
ei_add_test(qr)
|
||||
ei_add_test(rrqr)
|
||||
ei_add_test(qr_colpivoting)
|
||||
ei_add_test(qr_fullpivoting)
|
||||
ei_add_test(eigensolver_selfadjoint " " "${GSL_LIBRARIES}")
|
||||
ei_add_test(eigensolver_generic " " "${GSL_LIBRARIES}")
|
||||
ei_add_test(svd)
|
||||
|
@ -37,7 +37,7 @@ template<typename MatrixType> void qr()
|
||||
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
|
||||
MatrixType m1;
|
||||
createRandomMatrixOfRank(rank,rows,cols,m1);
|
||||
HouseholderRRQR<MatrixType> qr(m1);
|
||||
ColPivotingHouseholderQR<MatrixType> qr(m1);
|
||||
VERIFY_IS_APPROX(rank, qr.rank());
|
||||
|
||||
MatrixType r = qr.matrixQR();
|
||||
@ -74,7 +74,7 @@ template<typename MatrixType> void qr_invertible()
|
||||
m1 += a * a.adjoint();
|
||||
}
|
||||
|
||||
HouseholderRRQR<MatrixType> qr(m1);
|
||||
ColPivotingHouseholderQR<MatrixType> qr(m1);
|
||||
m3 = MatrixType::Random(size,size);
|
||||
qr.solve(m3, &m2);
|
||||
VERIFY_IS_APPROX(m3, m1*m2);
|
||||
@ -84,13 +84,13 @@ template<typename MatrixType> void qr_verify_assert()
|
||||
{
|
||||
MatrixType tmp;
|
||||
|
||||
HouseholderRRQR<MatrixType> qr;
|
||||
ColPivotingHouseholderQR<MatrixType> qr;
|
||||
VERIFY_RAISES_ASSERT(qr.matrixR())
|
||||
VERIFY_RAISES_ASSERT(qr.solve(tmp,&tmp))
|
||||
VERIFY_RAISES_ASSERT(qr.matrixQ())
|
||||
}
|
||||
|
||||
void test_rrqr()
|
||||
void test_qr_colpivoting()
|
||||
{
|
||||
for(int i = 0; i < 1; i++) {
|
||||
// FIXME : very weird bug here
|
115
test/qr_fullpivoting.cpp
Normal file
115
test/qr_fullpivoting.cpp
Normal file
@ -0,0 +1,115 @@
|
||||
// This file is part of Eigen, a lightweight C++ template library
|
||||
// for linear algebra.
|
||||
//
|
||||
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
|
||||
// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
|
||||
//
|
||||
// Eigen is free software; you can redistribute it and/or
|
||||
// modify it under the terms of the GNU Lesser General Public
|
||||
// License as published by the Free Software Foundation; either
|
||||
// version 3 of the License, or (at your option) any later version.
|
||||
//
|
||||
// Alternatively, you can redistribute it and/or
|
||||
// modify it under the terms of the GNU General Public License as
|
||||
// published by the Free Software Foundation; either version 2 of
|
||||
// the License, or (at your option) any later version.
|
||||
//
|
||||
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
|
||||
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
|
||||
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
|
||||
// GNU General Public License for more details.
|
||||
//
|
||||
// You should have received a copy of the GNU Lesser General Public
|
||||
// License and a copy of the GNU General Public License along with
|
||||
// Eigen. If not, see <http://www.gnu.org/licenses/>.
|
||||
|
||||
#include "main.h"
|
||||
#include <Eigen/QR>
|
||||
|
||||
template<typename MatrixType> void qr()
|
||||
{
|
||||
/* this test covers the following files: QR.h */
|
||||
int rows = ei_random<int>(20,200), cols = ei_random<int>(20,200), cols2 = ei_random<int>(20,200);
|
||||
int rank = ei_random<int>(1, std::min(rows, cols)-1);
|
||||
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> SquareMatrixType;
|
||||
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
|
||||
MatrixType m1;
|
||||
createRandomMatrixOfRank(rank,rows,cols,m1);
|
||||
FullPivotingHouseholderQR<MatrixType> qr(m1);
|
||||
VERIFY_IS_APPROX(rank, qr.rank());
|
||||
|
||||
MatrixType r = qr.matrixQR();
|
||||
// FIXME need better way to construct trapezoid
|
||||
for(int i = 0; i < rows; i++) for(int j = 0; j < cols; j++) if(i>j) r(i,j) = Scalar(0);
|
||||
|
||||
MatrixType b = qr.matrixQ() * r;
|
||||
|
||||
MatrixType c = MatrixType::Zero(rows,cols);
|
||||
|
||||
for(int i = 0; i < cols; ++i) c.col(qr.colsPermutation().coeff(i)) = b.col(i);
|
||||
VERIFY_IS_APPROX(m1, c);
|
||||
|
||||
MatrixType m2 = MatrixType::Random(cols,cols2);
|
||||
MatrixType m3 = m1*m2;
|
||||
m2 = MatrixType::Random(cols,cols2);
|
||||
qr.solve(m3, &m2);
|
||||
VERIFY_IS_APPROX(m3, m1*m2);
|
||||
}
|
||||
|
||||
template<typename MatrixType> void qr_invertible()
|
||||
{
|
||||
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
|
||||
int size = ei_random<int>(10,50);
|
||||
|
||||
MatrixType m1(size, size), m2(size, size), m3(size, size);
|
||||
m1 = MatrixType::Random(size,size);
|
||||
|
||||
if (ei_is_same_type<RealScalar,float>::ret)
|
||||
{
|
||||
// let's build a matrix more stable to inverse
|
||||
MatrixType a = MatrixType::Random(size,size*2);
|
||||
m1 += a * a.adjoint();
|
||||
}
|
||||
|
||||
FullPivotingHouseholderQR<MatrixType> qr(m1);
|
||||
m3 = MatrixType::Random(size,size);
|
||||
qr.solve(m3, &m2);
|
||||
VERIFY_IS_APPROX(m3, m1*m2);
|
||||
}
|
||||
|
||||
template<typename MatrixType> void qr_verify_assert()
|
||||
{
|
||||
MatrixType tmp;
|
||||
|
||||
FullPivotingHouseholderQR<MatrixType> qr;
|
||||
VERIFY_RAISES_ASSERT(qr.matrixR())
|
||||
VERIFY_RAISES_ASSERT(qr.solve(tmp,&tmp))
|
||||
VERIFY_RAISES_ASSERT(qr.matrixQ())
|
||||
}
|
||||
|
||||
void test_qr_fullpivoting()
|
||||
{
|
||||
for(int i = 0; i < 1; i++) {
|
||||
// FIXME : very weird bug here
|
||||
// CALL_SUBTEST( qr(Matrix2f()) );
|
||||
CALL_SUBTEST( qr<MatrixXf>() );
|
||||
CALL_SUBTEST( qr<MatrixXd>() );
|
||||
CALL_SUBTEST( qr<MatrixXcd>() );
|
||||
}
|
||||
|
||||
for(int i = 0; i < g_repeat; i++) {
|
||||
CALL_SUBTEST( qr_invertible<MatrixXf>() );
|
||||
CALL_SUBTEST( qr_invertible<MatrixXd>() );
|
||||
CALL_SUBTEST( qr_invertible<MatrixXcf>() );
|
||||
CALL_SUBTEST( qr_invertible<MatrixXcd>() );
|
||||
}
|
||||
|
||||
CALL_SUBTEST(qr_verify_assert<Matrix3f>());
|
||||
CALL_SUBTEST(qr_verify_assert<Matrix3d>());
|
||||
CALL_SUBTEST(qr_verify_assert<MatrixXf>());
|
||||
CALL_SUBTEST(qr_verify_assert<MatrixXd>());
|
||||
CALL_SUBTEST(qr_verify_assert<MatrixXcf>());
|
||||
CALL_SUBTEST(qr_verify_assert<MatrixXcd>());
|
||||
}
|
Loading…
Reference in New Issue
Block a user