// This file is part of Eigen, a lightweight C++ template library // for linear algebra. Eigen itself is part of the KDE project. // // Copyright (C) 2008 Benoit Jacob // Copyright (C) 2008 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_PART_H #define EIGEN_PART_H /** \nonstableyet * \class Part * * \brief Expression of a triangular matrix extracted from a given matrix * * \param MatrixType the type of the object in which we are taking the triangular part * \param Mode the kind of triangular matrix expression to construct. Can be UpperTriangular, StrictlyUpperTriangular, * UnitUpperTriangular, LowerTriangular, StrictlyLowerTriangular, UnitLowerTriangular. This is in fact a bit field; it must have either * UpperTriangularBit or LowerTriangularBit, and additionnaly it may have either ZeroDiagBit or * UnitDiagBit. * * This class represents an expression of the upper or lower triangular part of * a square matrix, possibly with a further assumption on the diagonal. It is the return type * of MatrixBase::part() and most of the time this is the only way it is used. * * \sa MatrixBase::part() */ template struct ei_traits > : ei_traits { typedef typename ei_nested::type MatrixTypeNested; typedef typename ei_unref::type _MatrixTypeNested; enum { Flags = (_MatrixTypeNested::Flags & (HereditaryBits) & (~(PacketAccessBit | DirectAccessBit | LinearAccessBit))) | Mode, CoeffReadCost = _MatrixTypeNested::CoeffReadCost }; }; template class Part : public MatrixBase > { public: EIGEN_GENERIC_PUBLIC_INTERFACE(Part) inline Part(const MatrixType& matrix) : m_matrix(matrix) { ei_assert(ei_are_flags_consistent::ret); } /** \sa MatrixBase::operator+=() */ template Part& operator+=(const Other& other); /** \sa MatrixBase::operator-=() */ template Part& operator-=(const Other& other); /** \sa MatrixBase::operator*=() */ Part& operator*=(const typename ei_traits::Scalar& other); /** \sa MatrixBase::operator/=() */ Part& operator/=(const typename ei_traits::Scalar& other); /** \sa operator=(), MatrixBase::lazyAssign() */ template void lazyAssign(const Other& other); /** \sa MatrixBase::operator=() */ template Part& operator=(const Other& other); inline int rows() const { return m_matrix.rows(); } inline int cols() const { return m_matrix.cols(); } inline int stride() const { return m_matrix.stride(); } inline Scalar coeff(int row, int col) const { // SelfAdjointBit doesn't play any role here: just because a matrix is selfadjoint doesn't say anything about // each individual coefficient, except for the not-very-useful-here fact that diagonal coefficients are real. if( ((Flags & LowerTriangularBit) && (col>row)) || ((Flags & UpperTriangularBit) && (row>col)) ) return (Scalar)0; if(Flags & UnitDiagBit) return col==row ? (Scalar)1 : m_matrix.coeff(row, col); else if(Flags & ZeroDiagBit) return col==row ? (Scalar)0 : m_matrix.coeff(row, col); else return m_matrix.coeff(row, col); } inline Scalar& coeffRef(int row, int col) { EIGEN_STATIC_ASSERT(!(Flags & UnitDiagBit), WRITING_TO_TRIANGULAR_PART_WITH_UNIT_DIAGONAL_IS_NOT_SUPPORTED) EIGEN_STATIC_ASSERT(!(Flags & SelfAdjointBit), COEFFICIENT_WRITE_ACCESS_TO_SELFADJOINT_NOT_SUPPORTED) ei_assert( (Mode==UpperTriangular && col>=row) || (Mode==LowerTriangular && col<=row) || (Mode==StrictlyUpperTriangular && col>row) || (Mode==StrictlyLowerTriangular && col row(int i) { return Base::row(i); } const Block row(int i) const { return Base::row(i); } /** discard any writes to a column */ const Block col(int i) { return Base::col(i); } const Block col(int i) const { return Base::col(i); } template void swap(const MatrixBase& other) { Part,Mode>(const_cast(m_matrix)).lazyAssign(other.derived()); } protected: const typename MatrixType::Nested m_matrix; private: Part& operator=(const Part&); }; /** \nonstableyet * \returns an expression of a triangular matrix 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. * * \addexample PartExample \label How to extract a triangular part of an arbitrary matrix * * Example: \include MatrixBase_extract.cpp * Output: \verbinclude MatrixBase_extract.out * * \sa class Part, part(), marked() */ template template const Part MatrixBase::part() const { return derived(); } template template inline Part& Part::operator=(const Other& other) { if(Other::Flags & EvalBeforeAssigningBit) { typename MatrixBase::PlainMatrixType other_evaluated(other.rows(), other.cols()); other_evaluated.template part().lazyAssign(other); lazyAssign(other_evaluated); } else lazyAssign(other.derived()); return *this; } template struct ei_part_assignment_impl { enum { col = (UnrollCount-1) / Derived1::RowsAtCompileTime, row = (UnrollCount-1) % Derived1::RowsAtCompileTime }; inline static void run(Derived1 &dst, const Derived2 &src) { ei_part_assignment_impl::run(dst, src); if(Mode == SelfAdjoint) { if(row == col) dst.coeffRef(row, col) = ei_real(src.coeff(row, col)); else if(row < col) dst.coeffRef(col, row) = ei_conj(dst.coeffRef(row, col) = src.coeff(row, col)); } else { ei_assert(Mode == UpperTriangular || Mode == LowerTriangular || Mode == StrictlyUpperTriangular || Mode == StrictlyLowerTriangular); if((Mode == UpperTriangular && row <= col) || (Mode == LowerTriangular && row >= col) || (Mode == StrictlyUpperTriangular && row < col) || (Mode == StrictlyLowerTriangular && row > col)) dst.copyCoeff(row, col, src); } } }; template struct ei_part_assignment_impl { inline static void run(Derived1 &dst, const Derived2 &src) { if(!(Mode & ZeroDiagBit)) dst.copyCoeff(0, 0, src); } }; // prevent buggy user code from causing an infinite recursion template struct ei_part_assignment_impl { inline static void run(Derived1 &, const Derived2 &) {} }; template struct ei_part_assignment_impl { inline static void run(Derived1 &dst, const Derived2 &src) { for(int j = 0; j < dst.cols(); ++j) for(int i = 0; i <= j; ++i) dst.copyCoeff(i, j, src); } }; template struct ei_part_assignment_impl { 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); } }; template struct ei_part_assignment_impl { inline static void run(Derived1 &dst, const Derived2 &src) { for(int j = 0; j < dst.cols(); ++j) for(int i = 0; i < j; ++i) dst.copyCoeff(i, j, src); } }; template struct ei_part_assignment_impl { 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); } }; template struct ei_part_assignment_impl { inline static void run(Derived1 &dst, const Derived2 &src) { for(int j = 0; j < dst.cols(); ++j) { for(int i = 0; i < j; ++i) dst.coeffRef(j, i) = ei_conj(dst.coeffRef(i, j) = src.coeff(i, j)); dst.coeffRef(j, j) = ei_real(src.coeff(j, j)); } } }; template template void Part::lazyAssign(const Other& other) { const bool unroll = MatrixType::SizeAtCompileTime * Other::CoeffReadCost / 2 <= EIGEN_UNROLLING_LIMIT; ei_assert(m_matrix.rows() == other.rows() && m_matrix.cols() == other.cols()); ei_part_assignment_impl ::run(m_matrix.const_cast_derived(), other.derived()); } /** \nonstableyet * \returns a lvalue pseudo-expression allowing to perform special operations on \c *this. * * The \a Mode parameter can have the following values: \c UpperTriangular, \c StrictlyUpperTriangular, \c LowerTriangular, * \c StrictlyLowerTriangular, \c SelfAdjoint. * * \addexample PartExample \label How to write to a triangular part of a matrix * * Example: \include MatrixBase_part.cpp * Output: \verbinclude MatrixBase_part.out * * \sa class Part, MatrixBase::extract(), MatrixBase::marked() */ template template inline Part MatrixBase::part() { return Part(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 { if(cols() != rows()) return false; RealScalar maxAbsOnUpperTriangularPart = static_cast(-1); for(int j = 0; j < cols(); ++j) for(int i = 0; i <= j; ++i) { RealScalar absValue = ei_abs(coeff(i,j)); if(absValue > maxAbsOnUpperTriangularPart) maxAbsOnUpperTriangularPart = absValue; } for(int j = 0; j < cols()-1; ++j) for(int i = j+1; i < rows(); ++i) if(!ei_isMuchSmallerThan(coeff(i, j), maxAbsOnUpperTriangularPart, prec)) 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 { if(cols() != rows()) return false; 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; } for(int j = 1; j < cols(); ++j) for(int i = 0; i < j; ++i) if(!ei_isMuchSmallerThan(coeff(i, j), maxAbsOnLowerTriangularPart, prec)) return false; return true; } template template inline Part& Part::operator+=(const Other& other) { return *this = m_matrix + other; } template template inline Part& Part::operator-=(const Other& other) { return *this = m_matrix - other; } template inline Part& Part::operator*= (const typename ei_traits::Scalar& other) { return *this = m_matrix * other; } template inline Part& Part::operator/= (const typename ei_traits::Scalar& other) { return *this = m_matrix / other; } #endif // EIGEN_PART_H