SparseTriangularView.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009-2015 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_SPARSE_TRIANGULARVIEW_H
12 #define EIGEN_SPARSE_TRIANGULARVIEW_H
13 
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
27 template<typename MatrixType, unsigned int Mode> class TriangularViewImpl<MatrixType,Mode,Sparse>
28  : public SparseMatrixBase<TriangularView<MatrixType,Mode> >
29 {
30  enum { SkipFirst = ((Mode&Lower) && !(MatrixType::Flags&RowMajorBit))
31  || ((Mode&Upper) && (MatrixType::Flags&RowMajorBit)),
32  SkipLast = !SkipFirst,
33  SkipDiag = (Mode&ZeroDiag) ? 1 : 0,
34  HasUnitDiag = (Mode&UnitDiag) ? 1 : 0
35  };
36 
38 
39  protected:
40  // dummy solve function to make TriangularView happy.
41  void solve() const;
42 
44  public:
45 
47 
48  typedef typename MatrixType::Nested MatrixTypeNested;
49  typedef std::remove_reference_t<MatrixTypeNested> MatrixTypeNestedNonRef;
51 
52  template<typename RhsType, typename DstType>
54  EIGEN_STRONG_INLINE void _solve_impl(const RhsType &rhs, DstType &dst) const {
55  if(!(internal::is_same<RhsType,DstType>::value && internal::extract_data(dst) == internal::extract_data(rhs)))
56  dst = rhs;
57  this->solveInPlace(dst);
58  }
59 
61  template<typename OtherDerived> void solveInPlace(MatrixBase<OtherDerived>& other) const;
62 
64  template<typename OtherDerived> void solveInPlace(SparseMatrixBase<OtherDerived>& other) const;
65 
66 };
67 
68 namespace internal {
69 
70 template<typename ArgType, unsigned int Mode>
71 struct unary_evaluator<TriangularView<ArgType,Mode>, IteratorBased>
72  : evaluator_base<TriangularView<ArgType,Mode> >
73 {
74  typedef TriangularView<ArgType,Mode> XprType;
75 
76 protected:
77 
78  typedef typename XprType::Scalar Scalar;
79  typedef typename XprType::StorageIndex StorageIndex;
80  typedef typename evaluator<ArgType>::InnerIterator EvalIterator;
81 
82  enum { SkipFirst = ((Mode&Lower) && !(ArgType::Flags&RowMajorBit))
83  || ((Mode&Upper) && (ArgType::Flags&RowMajorBit)),
84  SkipLast = !SkipFirst,
85  SkipDiag = (Mode&ZeroDiag) ? 1 : 0,
86  HasUnitDiag = (Mode&UnitDiag) ? 1 : 0
87  };
88 
89 public:
90 
91  enum {
92  CoeffReadCost = evaluator<ArgType>::CoeffReadCost,
93  Flags = XprType::Flags
94  };
95 
96  explicit unary_evaluator(const XprType &xpr) : m_argImpl(xpr.nestedExpression()), m_arg(xpr.nestedExpression()) {}
97 
98  inline Index nonZerosEstimate() const {
99  return m_argImpl.nonZerosEstimate();
100  }
101 
102  class InnerIterator : public EvalIterator
103  {
104  typedef EvalIterator Base;
105  public:
106 
107  EIGEN_STRONG_INLINE InnerIterator(const unary_evaluator& xprEval, Index outer)
108  : Base(xprEval.m_argImpl,outer), m_returnOne(false), m_containsDiag(Base::outer()<xprEval.m_arg.innerSize())
109  {
110  if(SkipFirst)
111  {
112  while((*this) && ((HasUnitDiag||SkipDiag) ? this->index()<=outer : this->index()<outer))
114  if(HasUnitDiag)
115  m_returnOne = m_containsDiag;
116  }
117  else if(HasUnitDiag && ((!Base::operator bool()) || Base::index()>=Base::outer()))
118  {
119  if((!SkipFirst) && Base::operator bool())
121  m_returnOne = m_containsDiag;
122  }
123  }
124 
125  EIGEN_STRONG_INLINE InnerIterator& operator++()
126  {
127  if(HasUnitDiag && m_returnOne)
128  m_returnOne = false;
129  else
130  {
132  if(HasUnitDiag && (!SkipFirst) && ((!Base::operator bool()) || Base::index()>=Base::outer()))
133  {
134  if((!SkipFirst) && Base::operator bool())
136  m_returnOne = m_containsDiag;
137  }
138  }
139  return *this;
140  }
141 
142  EIGEN_STRONG_INLINE operator bool() const
143  {
144  if(HasUnitDiag && m_returnOne)
145  return true;
146  if(SkipFirst) return Base::operator bool();
147  else
148  {
149  if (SkipDiag) return (Base::operator bool() && this->index() < this->outer());
150  else return (Base::operator bool() && this->index() <= this->outer());
151  }
152  }
153 
154 // inline Index row() const { return (ArgType::Flags&RowMajorBit ? Base::outer() : this->index()); }
155 // inline Index col() const { return (ArgType::Flags&RowMajorBit ? this->index() : Base::outer()); }
156  inline StorageIndex index() const
157  {
158  if(HasUnitDiag && m_returnOne) return internal::convert_index<StorageIndex>(Base::outer());
159  else return Base::index();
160  }
161  inline Scalar value() const
162  {
163  if(HasUnitDiag && m_returnOne) return Scalar(1);
164  else return Base::value();
165  }
166 
167  protected:
168  bool m_returnOne;
169  bool m_containsDiag;
170  private:
171  Scalar& valueRef();
172  };
173 
174 protected:
175  evaluator<ArgType> m_argImpl;
176  const ArgType& m_arg;
177 };
178 
179 } // end namespace internal
180 
181 template<typename Derived>
182 template<int Mode>
183 inline const TriangularView<const Derived, Mode>
185 {
186  return TriangularView<const Derived, Mode>(derived());
187 }
188 
189 } // end namespace Eigen
190 
191 #endif // EIGEN_SPARSE_TRIANGULARVIEW_H
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:883
#define EIGEN_SPARSE_PUBLIC_INTERFACE(Derived)
Definition: SparseUtil.h:45
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182
Base class of any sparse matrices or sparse expressions.
const TriangularView< const Derived, Mode > triangularView() const
internal::remove_all_t< MatrixTypeNested > MatrixTypeNestedCleaned
void solveInPlace(MatrixBase< OtherDerived > &other) const
void solveInPlace(SparseMatrixBase< OtherDerived > &other) const
std::remove_reference_t< MatrixTypeNested > MatrixTypeNestedNonRef
Expression of a triangular part in a matrix.
@ UnitDiag
Definition: Constants.h:215
@ ZeroDiag
Definition: Constants.h:217
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
const unsigned int RowMajorBit
Definition: Constants.h:68
bfloat16 operator++(bfloat16 &a)
Definition: BFloat16.h:298
typename remove_all< T >::type remove_all_t
Definition: Meta.h:119
EIGEN_ALWAYS_INLINE const T::Scalar * extract_data(const T &m)
Definition: BlasUtil.h:593
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Definition: BFloat16.h:222