TriangularMatrix.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
00005 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
00006 //
00007 // Eigen is free software; you can redistribute it and/or
00008 // modify it under the terms of the GNU Lesser General Public
00009 // License as published by the Free Software Foundation; either
00010 // version 3 of the License, or (at your option) any later version.
00011 //
00012 // Alternatively, you can redistribute it and/or
00013 // modify it under the terms of the GNU General Public License as
00014 // published by the Free Software Foundation; either version 2 of
00015 // the License, or (at your option) any later version.
00016 //
00017 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00018 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00019 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00020 // GNU General Public License for more details.
00021 //
00022 // You should have received a copy of the GNU Lesser General Public
00023 // License and a copy of the GNU General Public License along with
00024 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00025 
00026 #ifndef EIGEN_TRIANGULARMATRIX_H
00027 #define EIGEN_TRIANGULARMATRIX_H
00028 
00029 namespace internal {
00030   
00031 template<int Side, typename TriangularType, typename Rhs> struct triangular_solve_retval;
00032   
00033 }
00034 
00042 template<typename Derived> class TriangularBase : public EigenBase<Derived>
00043 {
00044   public:
00045 
00046     enum {
00047       Mode = internal::traits<Derived>::Mode,
00048       CoeffReadCost = internal::traits<Derived>::CoeffReadCost,
00049       RowsAtCompileTime = internal::traits<Derived>::RowsAtCompileTime,
00050       ColsAtCompileTime = internal::traits<Derived>::ColsAtCompileTime,
00051       MaxRowsAtCompileTime = internal::traits<Derived>::MaxRowsAtCompileTime,
00052       MaxColsAtCompileTime = internal::traits<Derived>::MaxColsAtCompileTime
00053     };
00054     typedef typename internal::traits<Derived>::Scalar Scalar;
00055     typedef typename internal::traits<Derived>::StorageKind StorageKind;
00056     typedef typename internal::traits<Derived>::Index Index;
00057     typedef typename internal::traits<Derived>::DenseMatrixType DenseMatrixType;
00058     typedef DenseMatrixType DenseType;
00059 
00060     inline TriangularBase() { eigen_assert(!((Mode&UnitDiag) && (Mode&ZeroDiag))); }
00061 
00062     inline Index rows() const { return derived().rows(); }
00063     inline Index cols() const { return derived().cols(); }
00064     inline Index outerStride() const { return derived().outerStride(); }
00065     inline Index innerStride() const { return derived().innerStride(); }
00066 
00067     inline Scalar coeff(Index row, Index col) const  { return derived().coeff(row,col); }
00068     inline Scalar& coeffRef(Index row, Index col) { return derived().coeffRef(row,col); }
00069 
00072     template<typename Other>
00073     EIGEN_STRONG_INLINE void copyCoeff(Index row, Index col, Other& other)
00074     {
00075       derived().coeffRef(row, col) = other.coeff(row, col);
00076     }
00077 
00078     inline Scalar operator()(Index row, Index col) const
00079     {
00080       check_coordinates(row, col);
00081       return coeff(row,col);
00082     }
00083     inline Scalar& operator()(Index row, Index col)
00084     {
00085       check_coordinates(row, col);
00086       return coeffRef(row,col);
00087     }
00088 
00089     #ifndef EIGEN_PARSED_BY_DOXYGEN
00090     inline const Derived& derived() const { return *static_cast<const Derived*>(this); }
00091     inline Derived& derived() { return *static_cast<Derived*>(this); }
00092     #endif // not EIGEN_PARSED_BY_DOXYGEN
00093 
00094     template<typename DenseDerived>
00095     void evalTo(MatrixBase<DenseDerived> &other) const;
00096     template<typename DenseDerived>
00097     void evalToLazy(MatrixBase<DenseDerived> &other) const;
00098 
00099     DenseMatrixType toDenseMatrix() const
00100     {
00101       DenseMatrixType res(rows(), cols());
00102       evalToLazy(res);
00103       return res;
00104     }
00105 
00106   protected:
00107 
00108     void check_coordinates(Index row, Index col) const
00109     {
00110       EIGEN_ONLY_USED_FOR_DEBUG(row);
00111       EIGEN_ONLY_USED_FOR_DEBUG(col);
00112       eigen_assert(col>=0 && col<cols() && row>=0 && row<rows());
00113       const int mode = int(Mode) & ~SelfAdjoint;
00114       EIGEN_ONLY_USED_FOR_DEBUG(mode);
00115       eigen_assert((mode==Upper && col>=row)
00116                 || (mode==Lower && col<=row)
00117                 || ((mode==StrictlyUpper || mode==UnitUpper) && col>row)
00118                 || ((mode==StrictlyLower || mode==UnitLower) && col<row));
00119     }
00120 
00121     #ifdef EIGEN_INTERNAL_DEBUGGING
00122     void check_coordinates_internal(Index row, Index col) const
00123     {
00124       check_coordinates(row, col);
00125     }
00126     #else
00127     void check_coordinates_internal(Index , Index ) const {}
00128     #endif
00129 
00130 };
00131 
00149 namespace internal {
00150 template<typename MatrixType, unsigned int _Mode>
00151 struct traits<TriangularView<MatrixType, _Mode> > : traits<MatrixType>
00152 {
00153   typedef typename nested<MatrixType>::type MatrixTypeNested;
00154   typedef typename remove_reference<MatrixTypeNested>::type MatrixTypeNestedNonRef;
00155   typedef typename remove_all<MatrixTypeNested>::type MatrixTypeNestedCleaned;
00156   typedef MatrixType ExpressionType;
00157   typedef typename MatrixType::PlainObject DenseMatrixType;
00158   enum {
00159     Mode = _Mode,
00160     Flags = (MatrixTypeNestedCleaned::Flags & (HereditaryBits) & (~(PacketAccessBit | DirectAccessBit | LinearAccessBit))) | Mode,
00161     CoeffReadCost = MatrixTypeNestedCleaned::CoeffReadCost
00162   };
00163 };
00164 }
00165 
00166 template<int Mode, bool LhsIsTriangular,
00167          typename Lhs, bool LhsIsVector,
00168          typename Rhs, bool RhsIsVector>
00169 struct TriangularProduct;
00170 
00171 template<typename _MatrixType, unsigned int _Mode> class TriangularView
00172   : public TriangularBase<TriangularView<_MatrixType, _Mode> >
00173 {
00174   public:
00175 
00176     typedef TriangularBase<TriangularView> Base;
00177     typedef typename internal::traits<TriangularView>::Scalar Scalar;
00178 
00179     typedef _MatrixType MatrixType;
00180     typedef typename internal::traits<TriangularView>::DenseMatrixType DenseMatrixType;
00181     typedef DenseMatrixType PlainObject;
00182 
00183   protected:
00184     typedef typename internal::traits<TriangularView>::MatrixTypeNested MatrixTypeNested;
00185     typedef typename internal::traits<TriangularView>::MatrixTypeNestedNonRef MatrixTypeNestedNonRef;
00186     typedef typename internal::traits<TriangularView>::MatrixTypeNestedCleaned MatrixTypeNestedCleaned;
00187 
00188     typedef typename internal::remove_all<typename MatrixType::ConjugateReturnType>::type MatrixConjugateReturnType;
00189     
00190   public:
00191     using Base::evalToLazy;
00192   
00193 
00194     typedef typename internal::traits<TriangularView>::StorageKind StorageKind;
00195     typedef typename internal::traits<TriangularView>::Index Index;
00196 
00197     enum {
00198       Mode = _Mode,
00199       TransposeMode = (Mode & Upper ? Lower : 0)
00200                     | (Mode & Lower ? Upper : 0)
00201                     | (Mode & (UnitDiag))
00202                     | (Mode & (ZeroDiag))
00203     };
00204 
00205     inline TriangularView(const MatrixType& matrix) : m_matrix(matrix)
00206     {}
00207 
00208     inline Index rows() const { return m_matrix.rows(); }
00209     inline Index cols() const { return m_matrix.cols(); }
00210     inline Index outerStride() const { return m_matrix.outerStride(); }
00211     inline Index innerStride() const { return m_matrix.innerStride(); }
00212 
00214     template<typename Other> TriangularView&  operator+=(const DenseBase<Other>& other) { return *this = m_matrix + other.derived(); }
00216     template<typename Other> TriangularView&  operator-=(const DenseBase<Other>& other) { return *this = m_matrix - other.derived(); }
00218     TriangularView&  operator*=(const typename internal::traits<MatrixType>::Scalar& other) { return *this = m_matrix * other; }
00220     TriangularView&  operator/=(const typename internal::traits<MatrixType>::Scalar& other) { return *this = m_matrix / other; }
00221 
00223     void fill(const Scalar& value) { setConstant(value); }
00225     TriangularView& setConstant(const Scalar& value)
00226     { return *this = MatrixType::Constant(rows(), cols(), value); }
00228     TriangularView& setZero() { return setConstant(Scalar(0)); }
00230     TriangularView& setOnes() { return setConstant(Scalar(1)); }
00231 
00235     inline Scalar coeff(Index row, Index col) const
00236     {
00237       Base::check_coordinates_internal(row, col);
00238       return m_matrix.coeff(row, col);
00239     }
00240 
00244     inline Scalar& coeffRef(Index row, Index col)
00245     {
00246       Base::check_coordinates_internal(row, col);
00247       return m_matrix.const_cast_derived().coeffRef(row, col);
00248     }
00249 
00250     const MatrixTypeNestedCleaned& nestedExpression() const { return m_matrix; }
00251     MatrixTypeNestedCleaned& nestedExpression() { return *const_cast<MatrixTypeNestedCleaned*>(&m_matrix); }
00252 
00254     template<typename OtherDerived>
00255     TriangularView& operator=(const TriangularBase<OtherDerived>& other);
00256 
00257     template<typename OtherDerived>
00258     TriangularView& operator=(const MatrixBase<OtherDerived>& other);
00259 
00260     TriangularView& operator=(const TriangularView& other)
00261     { return *this = other.nestedExpression(); }
00262 
00263     template<typename OtherDerived>
00264     void lazyAssign(const TriangularBase<OtherDerived>& other);
00265 
00266     template<typename OtherDerived>
00267     void lazyAssign(const MatrixBase<OtherDerived>& other);
00268 
00270     inline TriangularView<MatrixConjugateReturnType,Mode> conjugate()
00271     { return m_matrix.conjugate(); }
00273     inline const TriangularView<MatrixConjugateReturnType,Mode> conjugate() const
00274     { return m_matrix.conjugate(); }
00275 
00277     inline TriangularView<typename MatrixType::AdjointReturnType,TransposeMode> adjoint()
00278     { return m_matrix.adjoint(); }
00280     inline const TriangularView<typename MatrixType::AdjointReturnType,TransposeMode> adjoint() const
00281     { return m_matrix.adjoint(); }
00282 
00284     inline TriangularView<Transpose<MatrixType>,TransposeMode> transpose()
00285     {
00286       EIGEN_STATIC_ASSERT_LVALUE(MatrixType)
00287       return m_matrix.const_cast_derived().transpose();
00288     }
00290     inline const TriangularView<Transpose<MatrixType>,TransposeMode> transpose() const
00291     { return m_matrix.transpose(); }
00292 
00294     template<typename OtherDerived>
00295     TriangularProduct<Mode,true,MatrixType,false,OtherDerived,OtherDerived::IsVectorAtCompileTime>
00296     operator*(const MatrixBase<OtherDerived>& rhs) const
00297     {
00298       return TriangularProduct
00299               <Mode,true,MatrixType,false,OtherDerived,OtherDerived::IsVectorAtCompileTime>
00300               (m_matrix, rhs.derived());
00301     }
00302 
00304     template<typename OtherDerived> friend
00305     TriangularProduct<Mode,false,OtherDerived,OtherDerived::IsVectorAtCompileTime,MatrixType,false>
00306     operator*(const MatrixBase<OtherDerived>& lhs, const TriangularView& rhs)
00307     {
00308       return TriangularProduct
00309               <Mode,false,OtherDerived,OtherDerived::IsVectorAtCompileTime,MatrixType,false>
00310               (lhs.derived(),rhs.m_matrix);
00311     }
00312 
00313     #ifdef EIGEN2_SUPPORT
00314     template<typename OtherDerived>
00315     struct eigen2_product_return_type
00316     {
00317       typedef typename TriangularView<MatrixType,Mode>::DenseMatrixType DenseMatrixType;
00318       typedef typename OtherDerived::PlainObject::DenseType OtherPlainObject;
00319       typedef typename ProductReturnType<DenseMatrixType, OtherPlainObject>::Type ProdRetType;
00320       typedef typename ProdRetType::PlainObject type;
00321     };
00322     template<typename OtherDerived>
00323     const typename eigen2_product_return_type<OtherDerived>::type
00324     operator*(const EigenBase<OtherDerived>& rhs) const
00325     {
00326       typename OtherDerived::PlainObject::DenseType rhsPlainObject;
00327       rhs.evalTo(rhsPlainObject);
00328       return this->toDenseMatrix() * rhsPlainObject;
00329     }
00330     template<typename OtherMatrixType>
00331     bool isApprox(const TriangularView<OtherMatrixType, Mode>& other, typename NumTraits<Scalar>::Real precision = NumTraits<Scalar>::dummy_precision()) const
00332     {
00333       return this->toDenseMatrix().isApprox(other.toDenseMatrix(), precision);
00334     }
00335     template<typename OtherDerived>
00336     bool isApprox(const MatrixBase<OtherDerived>& other, typename NumTraits<Scalar>::Real precision = NumTraits<Scalar>::dummy_precision()) const
00337     {
00338       return this->toDenseMatrix().isApprox(other, precision);
00339     }
00340     #endif // EIGEN2_SUPPORT
00341 
00342     template<int Side, typename Other>
00343     inline const internal::triangular_solve_retval<Side,TriangularView, Other>
00344     solve(const MatrixBase<Other>& other) const;
00345 
00346     template<int Side, typename OtherDerived>
00347     void solveInPlace(const MatrixBase<OtherDerived>& other) const;
00348 
00349     template<typename Other>
00350     inline const internal::triangular_solve_retval<OnTheLeft,TriangularView, Other> 
00351     solve(const MatrixBase<Other>& other) const
00352     { return solve<OnTheLeft>(other); }
00353 
00354     template<typename OtherDerived>
00355     void solveInPlace(const MatrixBase<OtherDerived>& other) const
00356     { return solveInPlace<OnTheLeft>(other); }
00357 
00358     const SelfAdjointView<MatrixTypeNestedNonRef,Mode> selfadjointView() const
00359     {
00360       EIGEN_STATIC_ASSERT((Mode&UnitDiag)==0,PROGRAMMING_ERROR);
00361       return SelfAdjointView<MatrixTypeNestedNonRef,Mode>(m_matrix);
00362     }
00363     SelfAdjointView<MatrixTypeNestedNonRef,Mode> selfadjointView()
00364     {
00365       EIGEN_STATIC_ASSERT((Mode&UnitDiag)==0,PROGRAMMING_ERROR);
00366       return SelfAdjointView<MatrixTypeNestedNonRef,Mode>(m_matrix);
00367     }
00368 
00369     template<typename OtherDerived>
00370     void swap(TriangularBase<OtherDerived> const & other)
00371     {
00372       TriangularView<SwapWrapper<MatrixType>,Mode>(const_cast<MatrixType&>(m_matrix)).lazyAssign(other.derived());
00373     }
00374 
00375     template<typename OtherDerived>
00376     void swap(MatrixBase<OtherDerived> const & other)
00377     {
00378       TriangularView<SwapWrapper<MatrixType>,Mode>(const_cast<MatrixType&>(m_matrix)).lazyAssign(other.derived());
00379     }
00380 
00381     Scalar determinant() const
00382     {
00383       if (Mode & UnitDiag)
00384         return 1;
00385       else if (Mode & ZeroDiag)
00386         return 0;
00387       else
00388         return m_matrix.diagonal().prod();
00389     }
00390     
00391     // TODO simplify the following:
00392     template<typename ProductDerived, typename Lhs, typename Rhs>
00393     EIGEN_STRONG_INLINE TriangularView& operator=(const ProductBase<ProductDerived, Lhs,Rhs>& other)
00394     {
00395       setZero();
00396       return assignProduct(other,1);
00397     }
00398     
00399     template<typename ProductDerived, typename Lhs, typename Rhs>
00400     EIGEN_STRONG_INLINE TriangularView& operator+=(const ProductBase<ProductDerived, Lhs,Rhs>& other)
00401     {
00402       return assignProduct(other,1);
00403     }
00404     
00405     template<typename ProductDerived, typename Lhs, typename Rhs>
00406     EIGEN_STRONG_INLINE TriangularView& operator-=(const ProductBase<ProductDerived, Lhs,Rhs>& other)
00407     {
00408       return assignProduct(other,-1);
00409     }
00410     
00411     
00412     template<typename ProductDerived>
00413     EIGEN_STRONG_INLINE TriangularView& operator=(const ScaledProduct<ProductDerived>& other)
00414     {
00415       setZero();
00416       return assignProduct(other,other.alpha());
00417     }
00418     
00419     template<typename ProductDerived>
00420     EIGEN_STRONG_INLINE TriangularView& operator+=(const ScaledProduct<ProductDerived>& other)
00421     {
00422       return assignProduct(other,other.alpha());
00423     }
00424     
00425     template<typename ProductDerived>
00426     EIGEN_STRONG_INLINE TriangularView& operator-=(const ScaledProduct<ProductDerived>& other)
00427     {
00428       return assignProduct(other,-other.alpha());
00429     }
00430     
00431   protected:
00432     
00433     template<typename ProductDerived, typename Lhs, typename Rhs>
00434     EIGEN_STRONG_INLINE TriangularView& assignProduct(const ProductBase<ProductDerived, Lhs,Rhs>& prod, const Scalar& alpha);
00435 
00436     const MatrixTypeNested m_matrix;
00437 };
00438 
00439 /***************************************************************************
00440 * Implementation of triangular evaluation/assignment
00441 ***************************************************************************/
00442 
00443 namespace internal {
00444 
00445 template<typename Derived1, typename Derived2, unsigned int Mode, int UnrollCount, bool ClearOpposite>
00446 struct triangular_assignment_selector
00447 {
00448   enum {
00449     col = (UnrollCount-1) / Derived1::RowsAtCompileTime,
00450     row = (UnrollCount-1) % Derived1::RowsAtCompileTime
00451   };
00452   
00453   typedef typename Derived1::Scalar Scalar;
00454 
00455   inline static void run(Derived1 &dst, const Derived2 &src)
00456   {
00457     triangular_assignment_selector<Derived1, Derived2, Mode, UnrollCount-1, ClearOpposite>::run(dst, src);
00458 
00459     eigen_assert( Mode == Upper || Mode == Lower
00460             || Mode == StrictlyUpper || Mode == StrictlyLower
00461             || Mode == UnitUpper || Mode == UnitLower);
00462     if((Mode == Upper && row <= col)
00463     || (Mode == Lower && row >= col)
00464     || (Mode == StrictlyUpper && row < col)
00465     || (Mode == StrictlyLower && row > col)
00466     || (Mode == UnitUpper && row < col)
00467     || (Mode == UnitLower && row > col))
00468       dst.copyCoeff(row, col, src);
00469     else if(ClearOpposite)
00470     {
00471       if (Mode&UnitDiag && row==col)
00472         dst.coeffRef(row, col) = Scalar(1);
00473       else
00474         dst.coeffRef(row, col) = Scalar(0);
00475     }
00476   }
00477 };
00478 
00479 // prevent buggy user code from causing an infinite recursion
00480 template<typename Derived1, typename Derived2, unsigned int Mode, bool ClearOpposite>
00481 struct triangular_assignment_selector<Derived1, Derived2, Mode, 0, ClearOpposite>
00482 {
00483   inline static void run(Derived1 &, const Derived2 &) {}
00484 };
00485 
00486 template<typename Derived1, typename Derived2, bool ClearOpposite>
00487 struct triangular_assignment_selector<Derived1, Derived2, Upper, Dynamic, ClearOpposite>
00488 {
00489   typedef typename Derived1::Index Index;
00490   typedef typename Derived1::Scalar Scalar;
00491   inline static void run(Derived1 &dst, const Derived2 &src)
00492   {
00493     for(Index j = 0; j < dst.cols(); ++j)
00494     {
00495       Index maxi = (std::min)(j, dst.rows()-1);
00496       for(Index i = 0; i <= maxi; ++i)
00497         dst.copyCoeff(i, j, src);
00498       if (ClearOpposite)
00499         for(Index i = maxi+1; i < dst.rows(); ++i)
00500           dst.coeffRef(i, j) = Scalar(0);
00501     }
00502   }
00503 };
00504 
00505 template<typename Derived1, typename Derived2, bool ClearOpposite>
00506 struct triangular_assignment_selector<Derived1, Derived2, Lower, Dynamic, ClearOpposite>
00507 {
00508   typedef typename Derived1::Index Index;
00509   inline static void run(Derived1 &dst, const Derived2 &src)
00510   {
00511     for(Index j = 0; j < dst.cols(); ++j)
00512     {
00513       for(Index i = j; i < dst.rows(); ++i)
00514         dst.copyCoeff(i, j, src);
00515       Index maxi = (std::min)(j, dst.rows());
00516       if (ClearOpposite)
00517         for(Index i = 0; i < maxi; ++i)
00518           dst.coeffRef(i, j) = static_cast<typename Derived1::Scalar>(0);
00519     }
00520   }
00521 };
00522 
00523 template<typename Derived1, typename Derived2, bool ClearOpposite>
00524 struct triangular_assignment_selector<Derived1, Derived2, StrictlyUpper, Dynamic, ClearOpposite>
00525 {
00526   typedef typename Derived1::Index Index;
00527   inline static void run(Derived1 &dst, const Derived2 &src)
00528   {
00529     for(Index j = 0; j < dst.cols(); ++j)
00530     {
00531       Index maxi = (std::min)(j, dst.rows());
00532       for(Index i = 0; i < maxi; ++i)
00533         dst.copyCoeff(i, j, src);
00534       if (ClearOpposite)
00535         for(Index i = maxi; i < dst.rows(); ++i)
00536           dst.coeffRef(i, j) = 0;
00537     }
00538   }
00539 };
00540 
00541 template<typename Derived1, typename Derived2, bool ClearOpposite>
00542 struct triangular_assignment_selector<Derived1, Derived2, StrictlyLower, Dynamic, ClearOpposite>
00543 {
00544   typedef typename Derived1::Index Index;
00545   inline static void run(Derived1 &dst, const Derived2 &src)
00546   {
00547     for(Index j = 0; j < dst.cols(); ++j)
00548     {
00549       for(Index i = j+1; i < dst.rows(); ++i)
00550         dst.copyCoeff(i, j, src);
00551       Index maxi = (std::min)(j, dst.rows()-1);
00552       if (ClearOpposite)
00553         for(Index i = 0; i <= maxi; ++i)
00554           dst.coeffRef(i, j) = static_cast<typename Derived1::Scalar>(0);
00555     }
00556   }
00557 };
00558 
00559 template<typename Derived1, typename Derived2, bool ClearOpposite>
00560 struct triangular_assignment_selector<Derived1, Derived2, UnitUpper, Dynamic, ClearOpposite>
00561 {
00562   typedef typename Derived1::Index Index;
00563   inline static void run(Derived1 &dst, const Derived2 &src)
00564   {
00565     for(Index j = 0; j < dst.cols(); ++j)
00566     {
00567       Index maxi = (std::min)(j, dst.rows());
00568       for(Index i = 0; i < maxi; ++i)
00569         dst.copyCoeff(i, j, src);
00570       if (ClearOpposite)
00571       {
00572         for(Index i = maxi+1; i < dst.rows(); ++i)
00573           dst.coeffRef(i, j) = 0;
00574       }
00575     }
00576     dst.diagonal().setOnes();
00577   }
00578 };
00579 template<typename Derived1, typename Derived2, bool ClearOpposite>
00580 struct triangular_assignment_selector<Derived1, Derived2, UnitLower, Dynamic, ClearOpposite>
00581 {
00582   typedef typename Derived1::Index Index;
00583   inline static void run(Derived1 &dst, const Derived2 &src)
00584   {
00585     for(Index j = 0; j < dst.cols(); ++j)
00586     {
00587       Index maxi = (std::min)(j, dst.rows());
00588       for(Index i = maxi+1; i < dst.rows(); ++i)
00589         dst.copyCoeff(i, j, src);
00590       if (ClearOpposite)
00591       {
00592         for(Index i = 0; i < maxi; ++i)
00593           dst.coeffRef(i, j) = 0;
00594       }
00595     }
00596     dst.diagonal().setOnes();
00597   }
00598 };
00599 
00600 } // end namespace internal
00601 
00602 // FIXME should we keep that possibility
00603 template<typename MatrixType, unsigned int Mode>
00604 template<typename OtherDerived>
00605 inline TriangularView<MatrixType, Mode>&
00606 TriangularView<MatrixType, Mode>::operator=(const MatrixBase<OtherDerived>& other)
00607 {
00608   if(OtherDerived::Flags & EvalBeforeAssigningBit)
00609   {
00610     typename internal::plain_matrix_type<OtherDerived>::type other_evaluated(other.rows(), other.cols());
00611     other_evaluated.template triangularView<Mode>().lazyAssign(other.derived());
00612     lazyAssign(other_evaluated);
00613   }
00614   else
00615     lazyAssign(other.derived());
00616   return *this;
00617 }
00618 
00619 // FIXME should we keep that possibility
00620 template<typename MatrixType, unsigned int Mode>
00621 template<typename OtherDerived>
00622 void TriangularView<MatrixType, Mode>::lazyAssign(const MatrixBase<OtherDerived>& other)
00623 {
00624   enum {
00625     unroll = MatrixType::SizeAtCompileTime != Dynamic
00626           && internal::traits<OtherDerived>::CoeffReadCost != Dynamic
00627           && MatrixType::SizeAtCompileTime*internal::traits<OtherDerived>::CoeffReadCost/2 <= EIGEN_UNROLLING_LIMIT
00628   };
00629   eigen_assert(m_matrix.rows() == other.rows() && m_matrix.cols() == other.cols());
00630 
00631   internal::triangular_assignment_selector
00632     <MatrixType, OtherDerived, int(Mode),
00633     unroll ? int(MatrixType::SizeAtCompileTime) : Dynamic,
00634     false // do not change the opposite triangular part
00635     >::run(m_matrix.const_cast_derived(), other.derived());
00636 }
00637 
00638 
00639 
00640 template<typename MatrixType, unsigned int Mode>
00641 template<typename OtherDerived>
00642 inline TriangularView<MatrixType, Mode>&
00643 TriangularView<MatrixType, Mode>::operator=(const TriangularBase<OtherDerived>& other)
00644 {
00645   eigen_assert(Mode == int(OtherDerived::Mode));
00646   if(internal::traits<OtherDerived>::Flags & EvalBeforeAssigningBit)
00647   {
00648     typename OtherDerived::DenseMatrixType other_evaluated(other.rows(), other.cols());
00649     other_evaluated.template triangularView<Mode>().lazyAssign(other.derived().nestedExpression());
00650     lazyAssign(other_evaluated);
00651   }
00652   else
00653     lazyAssign(other.derived().nestedExpression());
00654   return *this;
00655 }
00656 
00657 template<typename MatrixType, unsigned int Mode>
00658 template<typename OtherDerived>
00659 void TriangularView<MatrixType, Mode>::lazyAssign(const TriangularBase<OtherDerived>& other)
00660 {
00661   enum {
00662     unroll = MatrixType::SizeAtCompileTime != Dynamic
00663                    && internal::traits<OtherDerived>::CoeffReadCost != Dynamic
00664                    && MatrixType::SizeAtCompileTime * internal::traits<OtherDerived>::CoeffReadCost / 2
00665                         <= EIGEN_UNROLLING_LIMIT
00666   };
00667   eigen_assert(m_matrix.rows() == other.rows() && m_matrix.cols() == other.cols());
00668 
00669   internal::triangular_assignment_selector
00670     <MatrixType, OtherDerived, int(Mode),
00671     unroll ? int(MatrixType::SizeAtCompileTime) : Dynamic,
00672     false // preserve the opposite triangular part
00673     >::run(m_matrix.const_cast_derived(), other.derived().nestedExpression());
00674 }
00675 
00676 /***************************************************************************
00677 * Implementation of TriangularBase methods
00678 ***************************************************************************/
00679 
00682 template<typename Derived>
00683 template<typename DenseDerived>
00684 void TriangularBase<Derived>::evalTo(MatrixBase<DenseDerived> &other) const
00685 {
00686   if(internal::traits<Derived>::Flags & EvalBeforeAssigningBit)
00687   {
00688     typename internal::plain_matrix_type<Derived>::type other_evaluated(rows(), cols());
00689     evalToLazy(other_evaluated);
00690     other.derived().swap(other_evaluated);
00691   }
00692   else
00693     evalToLazy(other.derived());
00694 }
00695 
00698 template<typename Derived>
00699 template<typename DenseDerived>
00700 void TriangularBase<Derived>::evalToLazy(MatrixBase<DenseDerived> &other) const
00701 {
00702   enum {
00703     unroll = DenseDerived::SizeAtCompileTime != Dynamic
00704                    && internal::traits<Derived>::CoeffReadCost != Dynamic
00705                    && DenseDerived::SizeAtCompileTime * internal::traits<Derived>::CoeffReadCost / 2
00706                         <= EIGEN_UNROLLING_LIMIT
00707   };
00708   other.derived().resize(this->rows(), this->cols());
00709 
00710   internal::triangular_assignment_selector
00711     <DenseDerived, typename internal::traits<Derived>::MatrixTypeNestedCleaned, Derived::Mode,
00712     unroll ? int(DenseDerived::SizeAtCompileTime) : Dynamic,
00713     true // clear the opposite triangular part
00714     >::run(other.derived(), derived().nestedExpression());
00715 }
00716 
00717 /***************************************************************************
00718 * Implementation of TriangularView methods
00719 ***************************************************************************/
00720 
00721 /***************************************************************************
00722 * Implementation of MatrixBase methods
00723 ***************************************************************************/
00724 
00725 #ifdef EIGEN2_SUPPORT
00726 
00727 // implementation of part<>(), including the SelfAdjoint case.
00728 
00729 namespace internal {
00730 template<typename MatrixType, unsigned int Mode>
00731 struct eigen2_part_return_type
00732 {
00733   typedef TriangularView<MatrixType, Mode> type;
00734 };
00735 
00736 template<typename MatrixType>
00737 struct eigen2_part_return_type<MatrixType, SelfAdjoint>
00738 {
00739   typedef SelfAdjointView<MatrixType, Upper> type;
00740 };
00741 }
00742 
00744 template<typename Derived>
00745 template<unsigned int Mode>
00746 const typename internal::eigen2_part_return_type<Derived, Mode>::type MatrixBase<Derived>::part() const
00747 {
00748   return derived();
00749 }
00750 
00752 template<typename Derived>
00753 template<unsigned int Mode>
00754 typename internal::eigen2_part_return_type<Derived, Mode>::type MatrixBase<Derived>::part()
00755 {
00756   return derived();
00757 }
00758 #endif
00759 
00771 template<typename Derived>
00772 template<unsigned int Mode>
00773 typename MatrixBase<Derived>::template TriangularViewReturnType<Mode>::Type
00774 MatrixBase<Derived>::triangularView()
00775 {
00776   return derived();
00777 }
00778 
00780 template<typename Derived>
00781 template<unsigned int Mode>
00782 typename MatrixBase<Derived>::template ConstTriangularViewReturnType<Mode>::Type
00783 MatrixBase<Derived>::triangularView() const
00784 {
00785   return derived();
00786 }
00787 
00793 template<typename Derived>
00794 bool MatrixBase<Derived>::isUpperTriangular(RealScalar prec) const
00795 {
00796   RealScalar maxAbsOnUpperPart = static_cast<RealScalar>(-1);
00797   for(Index j = 0; j < cols(); ++j)
00798   {
00799     Index maxi = (std::min)(j, rows()-1);
00800     for(Index i = 0; i <= maxi; ++i)
00801     {
00802       RealScalar absValue = internal::abs(coeff(i,j));
00803       if(absValue > maxAbsOnUpperPart) maxAbsOnUpperPart = absValue;
00804     }
00805   }
00806   RealScalar threshold = maxAbsOnUpperPart * prec;
00807   for(Index j = 0; j < cols(); ++j)
00808     for(Index i = j+1; i < rows(); ++i)
00809       if(internal::abs(coeff(i, j)) > threshold) return false;
00810   return true;
00811 }
00812 
00818 template<typename Derived>
00819 bool MatrixBase<Derived>::isLowerTriangular(RealScalar prec) const
00820 {
00821   RealScalar maxAbsOnLowerPart = static_cast<RealScalar>(-1);
00822   for(Index j = 0; j < cols(); ++j)
00823     for(Index i = j; i < rows(); ++i)
00824     {
00825       RealScalar absValue = internal::abs(coeff(i,j));
00826       if(absValue > maxAbsOnLowerPart) maxAbsOnLowerPart = absValue;
00827     }
00828   RealScalar threshold = maxAbsOnLowerPart * prec;
00829   for(Index j = 1; j < cols(); ++j)
00830   {
00831     Index maxi = (std::min)(j, rows()-1);
00832     for(Index i = 0; i < maxi; ++i)
00833       if(internal::abs(coeff(i, j)) > threshold) return false;
00834   }
00835   return true;
00836 }
00837 
00838 #endif // EIGEN_TRIANGULARMATRIX_H


libicr
Author(s): Robert Krug
autogenerated on Mon Jan 6 2014 11:33:49