Port LU module to evaluators (except image() and kernel())

This commit is contained in:
Gael Guennebaud 2014-02-20 15:26:15 +01:00
parent b2e1453e1e
commit 93125e372d
3 changed files with 178 additions and 61 deletions

View File

@ -12,6 +12,19 @@
namespace Eigen {
namespace internal {
template<typename _MatrixType> struct traits<FullPivLU<_MatrixType> >
: traits<_MatrixType>
{
typedef traits<_MatrixType> BaseTraits;
enum {
Flags = BaseTraits::Flags & RowMajorBit,
CoeffReadCost = Dynamic
};
};
} // end namespace internal
/** \ingroup LU_Module
*
* \class FullPivLU
@ -61,6 +74,7 @@ template<typename _MatrixType> class FullPivLU
typedef typename internal::plain_col_type<MatrixType, Index>::type IntColVectorType;
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationQType;
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationPType;
typedef typename MatrixType::PlainObject PlainObject;
/**
* \brief Default Constructor.
@ -209,6 +223,15 @@ template<typename _MatrixType> class FullPivLU
*
* \sa TriangularView::solve(), kernel(), inverse()
*/
#ifdef EIGEN_TEST_EVALUATORS
template<typename Rhs>
inline const Solve<FullPivLU, Rhs>
solve(const MatrixBase<Rhs>& b) const
{
eigen_assert(m_isInitialized && "LU is not initialized.");
return Solve<FullPivLU, Rhs>(*this, b.derived());
}
#else
template<typename Rhs>
inline const internal::solve_retval<FullPivLU, Rhs>
solve(const MatrixBase<Rhs>& b) const
@ -216,6 +239,7 @@ template<typename _MatrixType> class FullPivLU
eigen_assert(m_isInitialized && "LU is not initialized.");
return internal::solve_retval<FullPivLU, Rhs>(*this, b.derived());
}
#endif
/** \returns the determinant of the matrix of which
* *this is the LU decomposition. It has only linear complexity
@ -359,6 +383,14 @@ template<typename _MatrixType> class FullPivLU
*
* \sa MatrixBase::inverse()
*/
#ifdef EIGEN_TEST_EVALUATORS
inline const Inverse<FullPivLU> inverse() const
{
eigen_assert(m_isInitialized && "LU is not initialized.");
eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the inverse of a non-square matrix!");
return Inverse<FullPivLU>(*this);
}
#else
inline const internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType> inverse() const
{
eigen_assert(m_isInitialized && "LU is not initialized.");
@ -366,12 +398,19 @@ template<typename _MatrixType> class FullPivLU
return internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType>
(*this, MatrixType::Identity(m_lu.rows(), m_lu.cols()));
}
#endif
MatrixType reconstructedMatrix() const;
inline Index rows() const { return m_lu.rows(); }
inline Index cols() const { return m_lu.cols(); }
#ifndef EIGEN_PARSED_BY_DOXYGEN
template<typename RhsType, typename DstType>
EIGEN_DEVICE_FUNC
void _solve_impl(const RhsType &rhs, DstType &dst) const;
#endif
protected:
MatrixType m_lu;
PermutationPType m_p;
@ -662,6 +701,61 @@ struct image_retval<FullPivLU<_MatrixType> >
/***** Implementation of solve() *****************************************************/
} // end namespace internal
#ifndef EIGEN_PARSED_BY_DOXYGEN
template<typename _MatrixType>
template<typename RhsType, typename DstType>
void FullPivLU<_MatrixType>::_solve_impl(const RhsType &rhs, DstType &dst) const {
/* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1}.
* So we proceed as follows:
* Step 1: compute c = P * rhs.
* Step 2: replace c by the solution x to Lx = c. Exists because L is invertible.
* Step 3: replace c by the solution x to Ux = c. May or may not exist.
* Step 4: result = Q * c;
*/
const Index rows = this->rows(),
cols = this->cols(),
nonzero_pivots = this->nonzeroPivots();
eigen_assert(rhs.rows() == rows);
const Index smalldim = (std::min)(rows, cols);
if(nonzero_pivots == 0)
{
dst.setZero();
return;
}
typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
// Step 1
c = permutationP() * rhs;
// Step 2
m_lu.topLeftCorner(smalldim,smalldim)
.template triangularView<UnitLower>()
.solveInPlace(c.topRows(smalldim));
if(rows>cols)
c.bottomRows(rows-cols) -= m_lu.bottomRows(rows-cols) * c.topRows(cols);
// Step 3
m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
.template triangularView<Upper>()
.solveInPlace(c.topRows(nonzero_pivots));
// Step 4
for(Index i = 0; i < nonzero_pivots; ++i)
dst.row(permutationQ().indices().coeff(i)) = c.row(i);
for(Index i = nonzero_pivots; i < m_lu.cols(); ++i)
dst.row(permutationQ().indices().coeff(i)).setZero();
}
#endif
namespace internal {
#ifdef EIGEN_TEST_EVALUATORS
template<typename _MatrixType, typename Rhs>
struct solve_retval<FullPivLU<_MatrixType>, Rhs>
: solve_retval_base<FullPivLU<_MatrixType>, Rhs>
@ -670,53 +764,21 @@ struct solve_retval<FullPivLU<_MatrixType>, Rhs>
template<typename Dest> void evalTo(Dest& dst) const
{
/* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1}.
* So we proceed as follows:
* Step 1: compute c = P * rhs.
* Step 2: replace c by the solution x to Lx = c. Exists because L is invertible.
* Step 3: replace c by the solution x to Ux = c. May or may not exist.
* Step 4: result = Q * c;
*/
dec()._solve_impl(rhs(), dst);
}
};
#endif
const Index rows = dec().rows(), cols = dec().cols(),
nonzero_pivots = dec().nonzeroPivots();
eigen_assert(rhs().rows() == rows);
const Index smalldim = (std::min)(rows, cols);
/***** Implementation of inverse() *****************************************************/
if(nonzero_pivots == 0)
{
dst.setZero();
return;
}
typename Rhs::PlainObject c(rhs().rows(), rhs().cols());
// Step 1
c = dec().permutationP() * rhs();
// Step 2
dec().matrixLU()
.topLeftCorner(smalldim,smalldim)
.template triangularView<UnitLower>()
.solveInPlace(c.topRows(smalldim));
if(rows>cols)
{
c.bottomRows(rows-cols)
-= dec().matrixLU().bottomRows(rows-cols)
* c.topRows(cols);
}
// Step 3
dec().matrixLU()
.topLeftCorner(nonzero_pivots, nonzero_pivots)
.template triangularView<Upper>()
.solveInPlace(c.topRows(nonzero_pivots));
// Step 4
for(Index i = 0; i < nonzero_pivots; ++i)
dst.row(dec().permutationQ().indices().coeff(i)) = c.row(i);
for(Index i = nonzero_pivots; i < dec().matrixLU().cols(); ++i)
dst.row(dec().permutationQ().indices().coeff(i)).setZero();
template<typename DstXprType, typename MatrixType, typename Scalar>
struct Assignment<DstXprType, Inverse<FullPivLU<MatrixType> >, internal::assign_op<Scalar>, Dense2Dense, Scalar>
{
typedef FullPivLU<MatrixType> LuType;
typedef Inverse<LuType> SrcXprType;
static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<Scalar> &)
{
dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
}
};

View File

@ -13,6 +13,19 @@
namespace Eigen {
namespace internal {
template<typename _MatrixType> struct traits<PartialPivLU<_MatrixType> >
: traits<_MatrixType>
{
typedef traits<_MatrixType> BaseTraits;
enum {
Flags = BaseTraits::Flags & RowMajorBit,
CoeffReadCost = Dynamic
};
};
} // end namespace internal
/** \ingroup LU_Module
*
* \class PartialPivLU
@ -62,6 +75,7 @@ template<typename _MatrixType> class PartialPivLU
typedef typename MatrixType::Index Index;
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;
typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType;
typedef typename MatrixType::PlainObject PlainObject;
/**
@ -128,6 +142,15 @@ template<typename _MatrixType> class PartialPivLU
*
* \sa TriangularView::solve(), inverse(), computeInverse()
*/
#ifdef EIGEN_TEST_EVALUATORS
template<typename Rhs>
inline const Solve<PartialPivLU, Rhs>
solve(const MatrixBase<Rhs>& b) const
{
eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
return Solve<PartialPivLU, Rhs>(*this, b.derived());
}
#else
template<typename Rhs>
inline const internal::solve_retval<PartialPivLU, Rhs>
solve(const MatrixBase<Rhs>& b) const
@ -135,6 +158,7 @@ template<typename _MatrixType> class PartialPivLU
eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
return internal::solve_retval<PartialPivLU, Rhs>(*this, b.derived());
}
#endif
/** \returns the inverse of the matrix of which *this is the LU decomposition.
*
@ -143,12 +167,20 @@ template<typename _MatrixType> class PartialPivLU
*
* \sa MatrixBase::inverse(), LU::inverse()
*/
#ifdef EIGEN_TEST_EVALUATORS
inline const Inverse<PartialPivLU> inverse() const
{
eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
return Inverse<PartialPivLU>(*this);
}
#else
inline const internal::solve_retval<PartialPivLU,typename MatrixType::IdentityReturnType> inverse() const
{
eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
return internal::solve_retval<PartialPivLU,typename MatrixType::IdentityReturnType>
(*this, MatrixType::Identity(m_lu.rows(), m_lu.cols()));
}
#endif
/** \returns the determinant of the matrix of which
* *this is the LU decomposition. It has only linear complexity
@ -170,6 +202,30 @@ template<typename _MatrixType> class PartialPivLU
inline Index rows() const { return m_lu.rows(); }
inline Index cols() const { return m_lu.cols(); }
#ifndef EIGEN_PARSED_BY_DOXYGEN
template<typename RhsType, typename DstType>
EIGEN_DEVICE_FUNC
void _solve_impl(const RhsType &rhs, DstType &dst) const {
/* The decomposition PA = LU can be rewritten as A = P^{-1} L U.
* So we proceed as follows:
* Step 1: compute c = Pb.
* Step 2: replace c by the solution x to Lx = c.
* Step 3: replace c by the solution x to Ux = c.
*/
eigen_assert(rhs.rows() == m_lu.rows());
// Step 1
dst = permutationP() * rhs;
// Step 2
m_lu.template triangularView<UnitLower>().solveInPlace(dst);
// Step 3
m_lu.template triangularView<Upper>().solveInPlace(dst);
}
#endif
protected:
MatrixType m_lu;
PermutationType m_p;
@ -434,6 +490,7 @@ MatrixType PartialPivLU<MatrixType>::reconstructedMatrix() const
namespace internal {
#ifndef EIGEN_TEST_EVALUATORS
template<typename _MatrixType, typename Rhs>
struct solve_retval<PartialPivLU<_MatrixType>, Rhs>
: solve_retval_base<PartialPivLU<_MatrixType>, Rhs>
@ -442,23 +499,21 @@ struct solve_retval<PartialPivLU<_MatrixType>, Rhs>
template<typename Dest> void evalTo(Dest& dst) const
{
/* The decomposition PA = LU can be rewritten as A = P^{-1} L U.
* So we proceed as follows:
* Step 1: compute c = Pb.
* Step 2: replace c by the solution x to Lx = c.
* Step 3: replace c by the solution x to Ux = c.
*/
dec()._solve_impl(rhs(), dst);
}
};
#endif
eigen_assert(rhs().rows() == dec().matrixLU().rows());
/***** Implementation of inverse() *****************************************************/
// Step 1
dst = dec().permutationP() * rhs();
// Step 2
dec().matrixLU().template triangularView<UnitLower>().solveInPlace(dst);
// Step 3
dec().matrixLU().template triangularView<Upper>().solveInPlace(dst);
template<typename DstXprType, typename MatrixType, typename Scalar>
struct Assignment<DstXprType, Inverse<PartialPivLU<MatrixType> >, internal::assign_op<Scalar>, Dense2Dense, Scalar>
{
typedef PartialPivLU<MatrixType> LuType;
typedef Inverse<LuType> SrcXprType;
static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<Scalar> &)
{
dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
}
};