Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #ifndef EIGEN_SPARSE_TRIANGULARVIEW_H
00011 #define EIGEN_SPARSE_TRIANGULARVIEW_H
00012
00013 namespace Eigen {
00014
00015 namespace internal {
00016
00017 template<typename MatrixType, int Mode>
00018 struct traits<SparseTriangularView<MatrixType,Mode> >
00019 : public traits<MatrixType>
00020 {};
00021
00022 }
00023
00024 template<typename MatrixType, int Mode> class SparseTriangularView
00025 : public SparseMatrixBase<SparseTriangularView<MatrixType,Mode> >
00026 {
00027 enum { SkipFirst = ((Mode&Lower) && !(MatrixType::Flags&RowMajorBit))
00028 || ((Mode&Upper) && (MatrixType::Flags&RowMajorBit)),
00029 SkipLast = !SkipFirst,
00030 HasUnitDiag = (Mode&UnitDiag) ? 1 : 0
00031 };
00032
00033 public:
00034
00035 EIGEN_SPARSE_PUBLIC_INTERFACE(SparseTriangularView)
00036
00037 class InnerIterator;
00038 class ReverseInnerIterator;
00039
00040 inline Index rows() const { return m_matrix.rows(); }
00041 inline Index cols() const { return m_matrix.cols(); }
00042
00043 typedef typename MatrixType::Nested MatrixTypeNested;
00044 typedef typename internal::remove_reference<MatrixTypeNested>::type MatrixTypeNestedNonRef;
00045 typedef typename internal::remove_all<MatrixTypeNested>::type MatrixTypeNestedCleaned;
00046
00047 inline SparseTriangularView(const MatrixType& matrix) : m_matrix(matrix) {}
00048
00050 inline const MatrixTypeNestedCleaned& nestedExpression() const { return m_matrix; }
00051
00052 template<typename OtherDerived>
00053 typename internal::plain_matrix_type_column_major<OtherDerived>::type
00054 solve(const MatrixBase<OtherDerived>& other) const;
00055
00056 template<typename OtherDerived> void solveInPlace(MatrixBase<OtherDerived>& other) const;
00057 template<typename OtherDerived> void solveInPlace(SparseMatrixBase<OtherDerived>& other) const;
00058
00059 protected:
00060 MatrixTypeNested m_matrix;
00061 };
00062
00063 template<typename MatrixType, int Mode>
00064 class SparseTriangularView<MatrixType,Mode>::InnerIterator : public MatrixTypeNestedCleaned::InnerIterator
00065 {
00066 typedef typename MatrixTypeNestedCleaned::InnerIterator Base;
00067 public:
00068
00069 EIGEN_STRONG_INLINE InnerIterator(const SparseTriangularView& view, Index outer)
00070 : Base(view.nestedExpression(), outer), m_returnOne(false)
00071 {
00072 if(SkipFirst)
00073 {
00074 while((*this) && (HasUnitDiag ? this->index()<=outer : this->index()<outer))
00075 Base::operator++();
00076 if(HasUnitDiag)
00077 m_returnOne = true;
00078 }
00079 else if(HasUnitDiag && ((!Base::operator bool()) || Base::index()>=Base::outer()))
00080 {
00081 if((!SkipFirst) && Base::operator bool())
00082 Base::operator++();
00083 m_returnOne = true;
00084 }
00085 }
00086
00087 EIGEN_STRONG_INLINE InnerIterator& operator++()
00088 {
00089 if(HasUnitDiag && m_returnOne)
00090 m_returnOne = false;
00091 else
00092 {
00093 Base::operator++();
00094 if(HasUnitDiag && (!SkipFirst) && ((!Base::operator bool()) || Base::index()>=Base::outer()))
00095 {
00096 if((!SkipFirst) && Base::operator bool())
00097 Base::operator++();
00098 m_returnOne = true;
00099 }
00100 }
00101 return *this;
00102 }
00103
00104 inline Index row() const { return Base::row(); }
00105 inline Index col() const { return Base::col(); }
00106 inline Index index() const
00107 {
00108 if(HasUnitDiag && m_returnOne) return Base::outer();
00109 else return Base::index();
00110 }
00111 inline Scalar value() const
00112 {
00113 if(HasUnitDiag && m_returnOne) return Scalar(1);
00114 else return Base::value();
00115 }
00116
00117 EIGEN_STRONG_INLINE operator bool() const
00118 {
00119 if(HasUnitDiag && m_returnOne)
00120 return true;
00121 return (SkipFirst ? Base::operator bool() : (Base::operator bool() && this->index() <= this->outer()));
00122 }
00123 protected:
00124 bool m_returnOne;
00125 };
00126
00127 template<typename MatrixType, int Mode>
00128 class SparseTriangularView<MatrixType,Mode>::ReverseInnerIterator : public MatrixTypeNestedCleaned::ReverseInnerIterator
00129 {
00130 typedef typename MatrixTypeNestedCleaned::ReverseInnerIterator Base;
00131 public:
00132
00133 EIGEN_STRONG_INLINE ReverseInnerIterator(const SparseTriangularView& view, Index outer)
00134 : Base(view.nestedExpression(), outer)
00135 {
00136 eigen_assert((!HasUnitDiag) && "ReverseInnerIterator does not support yet triangular views with a unit diagonal");
00137 if(SkipLast)
00138 while((*this) && this->index()>outer)
00139 --(*this);
00140 }
00141
00142 EIGEN_STRONG_INLINE InnerIterator& operator--()
00143 { Base::operator--(); return *this; }
00144
00145 inline Index row() const { return Base::row(); }
00146 inline Index col() const { return Base::col(); }
00147
00148 EIGEN_STRONG_INLINE operator bool() const
00149 {
00150 return SkipLast ? Base::operator bool() : (Base::operator bool() && this->index() >= this->outer());
00151 }
00152 };
00153
00154 template<typename Derived>
00155 template<int Mode>
00156 inline const SparseTriangularView<Derived, Mode>
00157 SparseMatrixBase<Derived>::triangularView() const
00158 {
00159 return derived();
00160 }
00161
00162 }
00163
00164 #endif // EIGEN_SPARSE_TRIANGULARVIEW_H