diff --git a/Eigen/src/Core/AssignEvaluator.h b/Eigen/src/Core/AssignEvaluator.h index 084d650d8..73f20c1a1 100644 --- a/Eigen/src/Core/AssignEvaluator.h +++ b/Eigen/src/Core/AssignEvaluator.h @@ -721,7 +721,7 @@ void call_assignment_no_alias(Dst& dst, const Src& src, const Func& func) (int(Dst::ColsAtCompileTime) == 1 && int(Src::RowsAtCompileTime) == 1)) && int(Dst::SizeAtCompileTime) != 1 }; - + dst.resize(NeedToTranspose ? src.cols() : src.rows(), NeedToTranspose ? src.rows() : src.cols()); diff --git a/Eigen/src/LU/FullPivLU.h b/Eigen/src/LU/FullPivLU.h index 44699b763..86de91ccb 100644 --- a/Eigen/src/LU/FullPivLU.h +++ b/Eigen/src/LU/FullPivLU.h @@ -12,6 +12,19 @@ namespace Eigen { +namespace internal { +template struct traits > + : 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 class FullPivLU typedef typename internal::plain_col_type::type IntColVectorType; typedef PermutationMatrix PermutationQType; typedef PermutationMatrix PermutationPType; + typedef typename MatrixType::PlainObject PlainObject; /** * \brief Default Constructor. @@ -209,6 +223,15 @@ template class FullPivLU * * \sa TriangularView::solve(), kernel(), inverse() */ +#ifdef EIGEN_TEST_EVALUATORS + template + inline const Solve + solve(const MatrixBase& b) const + { + eigen_assert(m_isInitialized && "LU is not initialized."); + return Solve(*this, b.derived()); + } +#else template inline const internal::solve_retval solve(const MatrixBase& b) const @@ -216,6 +239,7 @@ template class FullPivLU eigen_assert(m_isInitialized && "LU is not initialized."); return internal::solve_retval(*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 class FullPivLU * * \sa MatrixBase::inverse() */ +#ifdef EIGEN_TEST_EVALUATORS + inline const Inverse 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(*this); + } +#else inline const internal::solve_retval inverse() const { eigen_assert(m_isInitialized && "LU is not initialized."); @@ -366,11 +398,18 @@ template class FullPivLU return internal::solve_retval (*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 + EIGEN_DEVICE_FUNC + void _solve_impl(const RhsType &rhs, DstType &dst) const; + #endif protected: MatrixType m_lu; @@ -662,6 +701,61 @@ struct image_retval > /***** Implementation of solve() *****************************************************/ +} // end namespace internal + +#ifndef EIGEN_PARSED_BY_DOXYGEN +template +template +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() + .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() + .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 struct solve_retval, Rhs> : solve_retval_base, Rhs> @@ -670,53 +764,21 @@ struct solve_retval, Rhs> template 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() - .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() - .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 +struct Assignment >, internal::assign_op, Dense2Dense, Scalar> +{ + typedef FullPivLU LuType; + typedef Inverse SrcXprType; + static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op &) + { + dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols())); } }; diff --git a/Eigen/src/LU/PartialPivLU.h b/Eigen/src/LU/PartialPivLU.h index 1d389ecac..ac53f7ab5 100644 --- a/Eigen/src/LU/PartialPivLU.h +++ b/Eigen/src/LU/PartialPivLU.h @@ -13,6 +13,19 @@ namespace Eigen { +namespace internal { +template struct traits > + : 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 class PartialPivLU typedef typename MatrixType::Index Index; typedef PermutationMatrix PermutationType; typedef Transpositions TranspositionType; + typedef typename MatrixType::PlainObject PlainObject; /** @@ -128,6 +142,15 @@ template class PartialPivLU * * \sa TriangularView::solve(), inverse(), computeInverse() */ +#ifdef EIGEN_TEST_EVALUATORS + template + inline const Solve + solve(const MatrixBase& b) const + { + eigen_assert(m_isInitialized && "PartialPivLU is not initialized."); + return Solve(*this, b.derived()); + } +#else template inline const internal::solve_retval solve(const MatrixBase& b) const @@ -135,6 +158,7 @@ template class PartialPivLU eigen_assert(m_isInitialized && "PartialPivLU is not initialized."); return internal::solve_retval(*this, b.derived()); } +#endif /** \returns the inverse of the matrix of which *this is the LU decomposition. * @@ -143,12 +167,20 @@ template class PartialPivLU * * \sa MatrixBase::inverse(), LU::inverse() */ +#ifdef EIGEN_TEST_EVALUATORS + inline const Inverse inverse() const + { + eigen_assert(m_isInitialized && "PartialPivLU is not initialized."); + return Inverse(*this); + } +#else inline const internal::solve_retval inverse() const { eigen_assert(m_isInitialized && "PartialPivLU is not initialized."); return internal::solve_retval (*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 @@ -169,6 +201,30 @@ template class PartialPivLU inline Index rows() const { return m_lu.rows(); } inline Index cols() const { return m_lu.cols(); } + + #ifndef EIGEN_PARSED_BY_DOXYGEN + template + 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().solveInPlace(dst); + + // Step 3 + m_lu.template triangularView().solveInPlace(dst); + } + #endif protected: MatrixType m_lu; @@ -434,6 +490,7 @@ MatrixType PartialPivLU::reconstructedMatrix() const namespace internal { +#ifndef EIGEN_TEST_EVALUATORS template struct solve_retval, Rhs> : solve_retval_base, Rhs> @@ -442,23 +499,21 @@ struct solve_retval, Rhs> template 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().solveInPlace(dst); - - // Step 3 - dec().matrixLU().template triangularView().solveInPlace(dst); +template +struct Assignment >, internal::assign_op, Dense2Dense, Scalar> +{ + typedef PartialPivLU LuType; + typedef Inverse SrcXprType; + static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op &) + { + dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols())); } };