Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
SparseTriangularView.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2009-2015 Gael Guennebaud <[email protected]>
5// Copyright (C) 2012 Désiré Nuentsa-Wakam <[email protected]>
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// IWYU pragma: private
15#include "./InternalHeaderCheck.h"
16
17namespace Eigen {
18
29template <typename MatrixType, unsigned int Mode>
30class TriangularViewImpl<MatrixType, Mode, Sparse> : public SparseMatrixBase<TriangularView<MatrixType, Mode> > {
31 enum {
32 SkipFirst =
33 ((Mode & Lower) && !(MatrixType::Flags & RowMajorBit)) || ((Mode & Upper) && (MatrixType::Flags & RowMajorBit)),
34 SkipLast = !SkipFirst,
35 SkipDiag = (Mode & ZeroDiag) ? 1 : 0,
36 HasUnitDiag = (Mode & UnitDiag) ? 1 : 0
37 };
38
40
41 protected:
42 // dummy solve function to make TriangularView happy.
43 void solve() const;
44
46
47 public:
48 EIGEN_SPARSE_PUBLIC_INTERFACE(TriangularViewType)
49
50 typedef typename MatrixType::Nested MatrixTypeNested;
51 typedef std::remove_reference_t<MatrixTypeNested> MatrixTypeNestedNonRef;
52 typedef internal::remove_all_t<MatrixTypeNested> MatrixTypeNestedCleaned;
53
54 template <typename RhsType, typename DstType>
55 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void _solve_impl(const RhsType& rhs, DstType& dst) const {
56 if (!(internal::is_same<RhsType, DstType>::value && internal::extract_data(dst) == internal::extract_data(rhs)))
57 dst = rhs;
58 this->solveInPlace(dst);
59 }
60
62 template <typename OtherDerived>
64
66 template <typename OtherDerived>
68};
69
70namespace internal {
71
72template <typename ArgType, unsigned int Mode>
73struct unary_evaluator<TriangularView<ArgType, Mode>, IteratorBased> : evaluator_base<TriangularView<ArgType, Mode> > {
74 typedef TriangularView<ArgType, Mode> XprType;
75
76 protected:
77 typedef typename XprType::Scalar Scalar;
78 typedef typename XprType::StorageIndex StorageIndex;
79 typedef typename evaluator<ArgType>::InnerIterator EvalIterator;
80
81 enum {
82 SkipFirst =
83 ((Mode & Lower) && !(ArgType::Flags & RowMajorBit)) || ((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 enum { CoeffReadCost = evaluator<ArgType>::CoeffReadCost, Flags = XprType::Flags };
91
92 explicit unary_evaluator(const XprType& xpr) : m_argImpl(xpr.nestedExpression()), m_arg(xpr.nestedExpression()) {}
93
94 inline Index nonZerosEstimate() const { return m_argImpl.nonZerosEstimate(); }
95
96 class InnerIterator : public EvalIterator {
97 typedef EvalIterator Base;
98
99 public:
100 EIGEN_STRONG_INLINE InnerIterator(const unary_evaluator& xprEval, Index outer)
101 : Base(xprEval.m_argImpl, outer),
102 m_returnOne(false),
103 m_containsDiag(Base::outer() < xprEval.m_arg.innerSize()) {
104 if (SkipFirst) {
105 while ((*this) && ((HasUnitDiag || SkipDiag) ? this->index() <= outer : this->index() < outer))
106 Base::operator++();
107 if (HasUnitDiag) m_returnOne = m_containsDiag;
108 } else if (HasUnitDiag && ((!Base::operator bool()) || Base::index() >= Base::outer())) {
109 if ((!SkipFirst) && Base::operator bool()) Base::operator++();
110 m_returnOne = m_containsDiag;
111 }
112 }
113
114 EIGEN_STRONG_INLINE InnerIterator& operator++() {
115 if (HasUnitDiag && m_returnOne)
116 m_returnOne = false;
117 else {
118 Base::operator++();
119 if (HasUnitDiag && (!SkipFirst) && ((!Base::operator bool()) || Base::index() >= Base::outer())) {
120 if ((!SkipFirst) && Base::operator bool()) Base::operator++();
121 m_returnOne = m_containsDiag;
122 }
123 }
124 return *this;
125 }
126
127 EIGEN_STRONG_INLINE operator bool() const {
128 if (HasUnitDiag && m_returnOne) return true;
129 if (SkipFirst)
130 return Base::operator bool();
131 else {
132 if (SkipDiag)
133 return (Base::operator bool() && this->index() < this->outer());
134 else
135 return (Base::operator bool() && this->index() <= this->outer());
136 }
137 }
138
139 inline Index row() const { return (ArgType::Flags & RowMajorBit ? Base::outer() : this->index()); }
140 inline Index col() const { return (ArgType::Flags & RowMajorBit ? this->index() : Base::outer()); }
141 inline StorageIndex index() const {
142 if (HasUnitDiag && m_returnOne)
143 return internal::convert_index<StorageIndex>(Base::outer());
144 else
145 return Base::index();
146 }
147 inline Scalar value() const {
148 if (HasUnitDiag && m_returnOne)
149 return Scalar(1);
150 else
151 return Base::value();
152 }
153
154 protected:
155 bool m_returnOne;
156 bool m_containsDiag;
157
158 private:
159 Scalar& valueRef();
160 };
161
162 protected:
163 evaluator<ArgType> m_argImpl;
164 const ArgType& m_arg;
165};
166
167} // end namespace internal
168
169template <typename Derived>
170template <int Mode>
171inline const TriangularView<const Derived, Mode> SparseMatrixBase<Derived>::triangularView() const {
172 return TriangularView<const Derived, Mode>(derived());
173}
174
175} // end namespace Eigen
176
177#endif // EIGEN_SPARSE_TRIANGULARVIEW_H
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
Base class of any sparse matrices or sparse expressions.
Definition SparseMatrixBase.h:30
void solveInPlace(MatrixBase< OtherDerived > &other) const
void solveInPlace(SparseMatrixBase< OtherDerived > &other) const
Expression of a triangular part in a matrix.
Definition TriangularMatrix.h:167
@ 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:70
Namespace containing all symbols from the Eigen library.
Definition Core:137
Definition Constants.h:519