// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2008 Benoit Jacob // Copyright (C) 2008-2009 Gael Guennebaud // // 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 . #ifndef EIGEN_TRIANGULARMATRIX_H #define EIGEN_TRIANGULARMATRIX_H /** \internal * * \class TriangularBase * * \brief Base class for triangular part in a matrix */ template class TriangularBase : public AnyMatrixBase { public: enum { Mode = ei_traits::Mode, CoeffReadCost = ei_traits::CoeffReadCost, RowsAtCompileTime = ei_traits::RowsAtCompileTime, ColsAtCompileTime = ei_traits::ColsAtCompileTime, MaxRowsAtCompileTime = ei_traits::MaxRowsAtCompileTime, MaxColsAtCompileTime = ei_traits::MaxColsAtCompileTime }; typedef typename ei_traits::Scalar Scalar; inline TriangularBase() { ei_assert(ei_are_flags_consistent::ret); } inline int rows() const { return derived().rows(); } inline int cols() const { return derived().cols(); } inline int stride() const { return derived().stride(); } inline Scalar coeff(int row, int col) const { return derived().coeff(row,col); } inline Scalar& coeffRef(int row, int col) { return derived().coeffRef(row,col); } /** \see MatrixBase::copyCoeff(row,col) */ template EIGEN_STRONG_INLINE void copyCoeff(int row, int col, Other& other) { derived().coeffRef(row, col) = other.coeff(row, col); } inline Scalar operator()(int row, int col) const { check_coordinates(row, col); return coeff(row,col); } inline Scalar& operator()(int row, int col) { check_coordinates(row, col); return coeffRef(row,col); } #ifndef EIGEN_PARSED_BY_DOXYGEN inline const Derived& derived() const { return *static_cast(this); } inline Derived& derived() { return *static_cast(this); } #endif // not EIGEN_PARSED_BY_DOXYGEN template void evalTo(MatrixBase &other) const; template void evalToLazy(MatrixBase &other) const; protected: void check_coordinates(int row, int col) { ei_assert(col>=0 && col=0 && row=row) || (Mode==LowerTriangular && col<=row) || ((Mode==StrictlyUpperTriangular || Mode==UnitUpperTriangular) && col>row) || ((Mode==StrictlyLowerTriangular || Mode==UnitLowerTriangular) && col struct ei_traits > : ei_traits { typedef typename ei_nested::type MatrixTypeNested; typedef typename ei_unref::type _MatrixTypeNested; typedef MatrixType ExpressionType; enum { Mode = _Mode, Flags = (_MatrixTypeNested::Flags & (HereditaryBits) & (~(PacketAccessBit | DirectAccessBit | LinearAccessBit))) | Mode, CoeffReadCost = _MatrixTypeNested::CoeffReadCost }; }; template struct TriangularProduct; template class TriangularView : public TriangularBase > { public: typedef TriangularBase Base; typedef typename ei_traits::Scalar Scalar; typedef _MatrixType MatrixType; typedef typename MatrixType::PlainMatrixType DenseMatrixType; typedef typename MatrixType::Nested MatrixTypeNested; typedef typename ei_cleantype::type _MatrixTypeNested; enum { Mode = _Mode, TransposeMode = (Mode & UpperTriangularBit ? LowerTriangularBit : 0) | (Mode & LowerTriangularBit ? UpperTriangularBit : 0) | (Mode & (ZeroDiagBit | UnitDiagBit)) }; inline TriangularView(const MatrixType& matrix) : m_matrix(matrix) { ei_assert(ei_are_flags_consistent::ret); } inline int rows() const { return m_matrix.rows(); } inline int cols() const { return m_matrix.cols(); } inline int stride() const { return m_matrix.stride(); } /** \sa MatrixBase::operator+=() */ template TriangularView& operator+=(const Other& other) { return *this = m_matrix + other; } /** \sa MatrixBase::operator-=() */ template TriangularView& operator-=(const Other& other) { return *this = m_matrix - other; } /** \sa MatrixBase::operator*=() */ TriangularView& operator*=(const typename ei_traits::Scalar& other) { return *this = m_matrix * other; } /** \sa MatrixBase::operator/=() */ TriangularView& operator/=(const typename ei_traits::Scalar& other) { return *this = m_matrix / other; } /** \sa MatrixBase::fill() */ void fill(const Scalar& value) { setConstant(value); } /** \sa MatrixBase::setConstant() */ TriangularView& setConstant(const Scalar& value) { return *this = MatrixType::Constant(rows(), cols(), value); } /** \sa MatrixBase::setZero() */ TriangularView& setZero() { return setConstant(Scalar(0)); } /** \sa MatrixBase::setOnes() */ TriangularView& setOnes() { return setConstant(Scalar(1)); } /** \sa MatrixBase::coeff() * \warning the coordinates must fit into the referenced triangular part */ inline Scalar coeff(int row, int col) const { Base::check_coordinates_internal(row, col); return m_matrix.coeff(row, col); } /** \sa MatrixBase::coeffRef() * \warning the coordinates must fit into the referenced triangular part */ inline Scalar& coeffRef(int row, int col) { Base::check_coordinates_internal(row, col); return m_matrix.const_cast_derived().coeffRef(row, col); } /** \internal */ const MatrixType& _expression() const { return m_matrix; } /** Assigns a triangular matrix to a triangular part of a dense matrix */ template TriangularView& operator=(const TriangularBase& other); template TriangularView& operator=(const MatrixBase& other); TriangularView& operator=(const TriangularView& other) { return *this = other._expression(); } template void lazyAssign(const TriangularBase& other); template void lazyAssign(const MatrixBase& other); /** \sa MatrixBase::adjoint() */ inline TriangularView,TransposeMode> adjoint() { return m_matrix.adjoint().nestByValue(); } /** \sa MatrixBase::adjoint() const */ inline const TriangularView,TransposeMode> adjoint() const { return m_matrix.adjoint().nestByValue(); } /** \sa MatrixBase::transpose() */ inline TriangularView >,TransposeMode> transpose() { return m_matrix.transpose().nestByValue(); } /** \sa MatrixBase::transpose() const */ inline const TriangularView >,TransposeMode> transpose() const { return m_matrix.transpose().nestByValue(); } DenseMatrixType toDenseMatrix() const { DenseMatrixType res(rows(), cols()); evalToLazy(res); return res; } /** Efficient triangular matrix times vector/matrix product */ template TriangularProduct operator*(const MatrixBase& rhs) const { return TriangularProduct (m_matrix, rhs.derived()); } /** Efficient vector/matrix times triangular matrix product */ template friend TriangularProduct operator*(const MatrixBase& lhs, const TriangularView& rhs) { return TriangularProduct (lhs.derived(),rhs.m_matrix); } template typename ei_plain_matrix_type_column_major::type solve(const MatrixBase& other) const; template void solveInPlace(const MatrixBase& other) const; template typename ei_plain_matrix_type_column_major::type solve(const MatrixBase& other) const { return solve(other); } template void solveInPlace(const MatrixBase& other) const { return solveInPlace(other); } const SelfAdjointView<_MatrixTypeNested,Mode> selfadjointView() const { EIGEN_STATIC_ASSERT((Mode&UnitDiagBit)==0,PROGRAMMING_ERROR); return SelfAdjointView<_MatrixTypeNested,Mode>(m_matrix); } SelfAdjointView<_MatrixTypeNested,Mode> selfadjointView() { EIGEN_STATIC_ASSERT((Mode&UnitDiagBit)==0,PROGRAMMING_ERROR); return SelfAdjointView<_MatrixTypeNested,Mode>(m_matrix); } template void swap(TriangularBase EIGEN_REF_TO_TEMPORARY other) { TriangularView,Mode>(const_cast(m_matrix)).lazyAssign(other.derived()); } template void swap(MatrixBase EIGEN_REF_TO_TEMPORARY other) { TriangularView,Mode>(const_cast(m_matrix)).lazyAssign(other.derived()); } protected: const MatrixTypeNested m_matrix; }; /*************************************************************************** * Implementation of triangular evaluation/assignment ***************************************************************************/ template struct ei_triangular_assignment_selector { enum { col = (UnrollCount-1) / Derived1::RowsAtCompileTime, row = (UnrollCount-1) % Derived1::RowsAtCompileTime }; inline static void run(Derived1 &dst, const Derived2 &src) { ei_triangular_assignment_selector::run(dst, src); ei_assert( Mode == UpperTriangular || Mode == LowerTriangular || Mode == StrictlyUpperTriangular || Mode == StrictlyLowerTriangular || Mode == UnitUpperTriangular || Mode == UnitLowerTriangular); if((Mode == UpperTriangular && row <= col) || (Mode == LowerTriangular && row >= col) || (Mode == StrictlyUpperTriangular && row < col) || (Mode == StrictlyLowerTriangular && row > col) || (Mode == UnitUpperTriangular && row < col) || (Mode == UnitLowerTriangular && row > col)) dst.copyCoeff(row, col, src); else if(ClearOpposite) { if (Mode&UnitDiagBit && row==col) dst.coeffRef(row, col) = 1; else dst.coeffRef(row, col) = 0; } } }; // prevent buggy user code from causing an infinite recursion template struct ei_triangular_assignment_selector { inline static void run(Derived1 &, const Derived2 &) {} }; template struct ei_triangular_assignment_selector { inline static void run(Derived1 &dst, const Derived2 &src) { for(int j = 0; j < dst.cols(); ++j) { int maxi = std::min(j, dst.rows()-1); for(int i = 0; i <= maxi; ++i) dst.copyCoeff(i, j, src); if (ClearOpposite) for(int i = maxi+1; i < dst.rows(); ++i) dst.coeffRef(i, j) = 0; } } }; template struct ei_triangular_assignment_selector { inline static void run(Derived1 &dst, const Derived2 &src) { for(int j = 0; j < dst.cols(); ++j) { for(int i = j; i < dst.rows(); ++i) dst.copyCoeff(i, j, src); int maxi = std::min(j, dst.rows()); if (ClearOpposite) for(int i = 0; i < maxi; ++i) dst.coeffRef(i, j) = 0; } } }; template struct ei_triangular_assignment_selector { inline static void run(Derived1 &dst, const Derived2 &src) { for(int j = 0; j < dst.cols(); ++j) { int maxi = std::min(j, dst.rows()); for(int i = 0; i < maxi; ++i) dst.copyCoeff(i, j, src); if (ClearOpposite) for(int i = maxi; i < dst.rows(); ++i) dst.coeffRef(i, j) = 0; } } }; template struct ei_triangular_assignment_selector { inline static void run(Derived1 &dst, const Derived2 &src) { for(int j = 0; j < dst.cols(); ++j) { for(int i = j+1; i < dst.rows(); ++i) dst.copyCoeff(i, j, src); int maxi = std::min(j, dst.rows()-1); if (ClearOpposite) for(int i = 0; i <= maxi; ++i) dst.coeffRef(i, j) = 0; } } }; template struct ei_triangular_assignment_selector { inline static void run(Derived1 &dst, const Derived2 &src) { for(int j = 0; j < dst.cols(); ++j) { int maxi = std::min(j, dst.rows()); for(int i = 0; i < maxi; ++i) dst.copyCoeff(i, j, src); if (ClearOpposite) { for(int i = maxi+1; i < dst.rows(); ++i) dst.coeffRef(i, j) = 0; } } dst.diagonal().setOnes(); } }; template struct ei_triangular_assignment_selector { inline static void run(Derived1 &dst, const Derived2 &src) { for(int j = 0; j < dst.cols(); ++j) { int maxi = std::min(j, dst.rows()); for(int i = maxi+1; i < dst.rows(); ++i) dst.copyCoeff(i, j, src); if (ClearOpposite) { for(int i = 0; i < maxi; ++i) dst.coeffRef(i, j) = 0; } } dst.diagonal().setOnes(); } }; // FIXME should we keep that possibility template template inline TriangularView& TriangularView::operator=(const MatrixBase& other) { if(OtherDerived::Flags & EvalBeforeAssigningBit) { typename OtherDerived::PlainMatrixType other_evaluated(other.rows(), other.cols()); other_evaluated.template triangularView().lazyAssign(other.derived()); lazyAssign(other_evaluated); } else lazyAssign(other.derived()); return *this; } // FIXME should we keep that possibility template template void TriangularView::lazyAssign(const MatrixBase& other) { const bool unroll = MatrixType::SizeAtCompileTime * ei_traits::CoeffReadCost / 2 <= EIGEN_UNROLLING_LIMIT; ei_assert(m_matrix.rows() == other.rows() && m_matrix.cols() == other.cols()); ei_triangular_assignment_selector ::run(m_matrix.const_cast_derived(), other.derived()); } template template inline TriangularView& TriangularView::operator=(const TriangularBase& other) { ei_assert(Mode == OtherDerived::Mode); if(ei_traits::Flags & EvalBeforeAssigningBit) { typename OtherDerived::DenseMatrixType other_evaluated(other.rows(), other.cols()); other_evaluated.template triangularView().lazyAssign(other.derived()); lazyAssign(other_evaluated); } else lazyAssign(other.derived()); return *this; } template template void TriangularView::lazyAssign(const TriangularBase& other) { const bool unroll = MatrixType::SizeAtCompileTime * ei_traits::CoeffReadCost / 2 <= EIGEN_UNROLLING_LIMIT; ei_assert(m_matrix.rows() == other.rows() && m_matrix.cols() == other.cols()); ei_triangular_assignment_selector ::run(m_matrix.const_cast_derived(), other.derived()._expression()); } /*************************************************************************** * Implementation of TriangularBase methods ***************************************************************************/ /** Assigns a triangular or selfadjoint matrix to a dense matrix. * If the matrix is triangular, the opposite part is set to zero. */ template template void TriangularBase::evalTo(MatrixBase &other) const { if(ei_traits::Flags & EvalBeforeAssigningBit) { typename Derived::PlainMatrixType other_evaluated(rows(), cols()); evalToLazy(other_evaluated); other.derived().swap(other_evaluated); } else evalToLazy(other.derived()); } /** Assigns a triangular or selfadjoint matrix to a dense matrix. * If the matrix is triangular, the opposite part is set to zero. */ template template void TriangularBase::evalToLazy(MatrixBase &other) const { const bool unroll = DenseDerived::SizeAtCompileTime * Derived::CoeffReadCost / 2 <= EIGEN_UNROLLING_LIMIT; ei_assert(this->rows() == other.rows() && this->cols() == other.cols()); ei_triangular_assignment_selector ::ExpressionType, Derived::Mode, unroll ? int(DenseDerived::SizeAtCompileTime) : Dynamic, true // clear the opposite triangular part >::run(other.derived(), derived()._expression()); } /*************************************************************************** * Implementation of TriangularView methods ***************************************************************************/ /*************************************************************************** * Implementation of MatrixBase methods ***************************************************************************/ /** \deprecated use MatrixBase::triangularView() */ template template EIGEN_DEPRECATED const TriangularView MatrixBase::part() const { return derived(); } /** \deprecated use MatrixBase::triangularView() */ template template EIGEN_DEPRECATED TriangularView MatrixBase::part() { return derived(); } /** \nonstableyet * \returns an expression of a triangular view extracted from the current matrix * * The parameter \a Mode can have the following values: \c UpperTriangular, \c StrictlyUpperTriangular, \c UnitUpperTriangular, * \c LowerTriangular, \c StrictlyLowerTriangular, \c UnitLowerTriangular. * * Example: \include MatrixBase_extract.cpp * Output: \verbinclude MatrixBase_extract.out * * \sa class TriangularView */ template template TriangularView MatrixBase::triangularView() { return derived(); } /** This is the const version of MatrixBase::triangularView() */ template template const TriangularView MatrixBase::triangularView() const { return derived(); } /** \returns true if *this is approximately equal to an upper triangular matrix, * within the precision given by \a prec. * * \sa isLowerTriangular(), extract(), part(), marked() */ template bool MatrixBase::isUpperTriangular(RealScalar prec) const { RealScalar maxAbsOnUpperTriangularPart = static_cast(-1); for(int j = 0; j < cols(); ++j) { int maxi = std::min(j, rows()-1); for(int i = 0; i <= maxi; ++i) { RealScalar absValue = ei_abs(coeff(i,j)); if(absValue > maxAbsOnUpperTriangularPart) maxAbsOnUpperTriangularPart = absValue; } } RealScalar threshold = maxAbsOnUpperTriangularPart * prec; for(int j = 0; j < cols(); ++j) for(int i = j+1; i < rows(); ++i) if(ei_abs(coeff(i, j)) > threshold) return false; return true; } /** \returns true if *this is approximately equal to a lower triangular matrix, * within the precision given by \a prec. * * \sa isUpperTriangular(), extract(), part(), marked() */ template bool MatrixBase::isLowerTriangular(RealScalar prec) const { RealScalar maxAbsOnLowerTriangularPart = static_cast(-1); for(int j = 0; j < cols(); ++j) for(int i = j; i < rows(); ++i) { RealScalar absValue = ei_abs(coeff(i,j)); if(absValue > maxAbsOnLowerTriangularPart) maxAbsOnLowerTriangularPart = absValue; } RealScalar threshold = maxAbsOnLowerTriangularPart * prec; for(int j = 1; j < cols(); ++j) { int maxi = std::min(j, rows()-1); for(int i = 0; i < maxi; ++i) if(ei_abs(coeff(i, j)) > threshold) return false; } return true; } #endif // EIGEN_TRIANGULARMATRIX_H