TriangularSolver.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) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_SPARSETRIANGULARSOLVER_H
11 #define EIGEN_SPARSETRIANGULARSOLVER_H
12 
13 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
17 namespace internal {
18 
19 template<typename Lhs, typename Rhs, int Mode,
20  int UpLo = (Mode & Lower)
21  ? Lower
22  : (Mode & Upper)
23  ? Upper
24  : -1,
25  int StorageOrder = int(traits<Lhs>::Flags) & RowMajorBit>
26 struct sparse_solve_triangular_selector;
27 
28 // forward substitution, row-major
29 template<typename Lhs, typename Rhs, int Mode>
30 struct sparse_solve_triangular_selector<Lhs,Rhs,Mode,Lower,RowMajor>
31 {
32  typedef typename Rhs::Scalar Scalar;
33  typedef evaluator<Lhs> LhsEval;
34  typedef typename evaluator<Lhs>::InnerIterator LhsIterator;
35  static void run(const Lhs& lhs, Rhs& other)
36  {
37  LhsEval lhsEval(lhs);
38  for(Index col=0 ; col<other.cols() ; ++col)
39  {
40  for(Index i=0; i<lhs.rows(); ++i)
41  {
42  Scalar tmp = other.coeff(i,col);
43  Scalar lastVal(0);
44  Index lastIndex = 0;
45  for(LhsIterator it(lhsEval, i); it; ++it)
46  {
47  lastVal = it.value();
48  lastIndex = it.index();
49  if(lastIndex==i)
50  break;
51  tmp -= lastVal * other.coeff(lastIndex,col);
52  }
53  if (Mode & UnitDiag)
54  other.coeffRef(i,col) = tmp;
55  else
56  {
57  eigen_assert(lastIndex==i);
58  other.coeffRef(i,col) = tmp/lastVal;
59  }
60  }
61  }
62  }
63 };
64 
65 // backward substitution, row-major
66 template<typename Lhs, typename Rhs, int Mode>
67 struct sparse_solve_triangular_selector<Lhs,Rhs,Mode,Upper,RowMajor>
68 {
69  typedef typename Rhs::Scalar Scalar;
70  typedef evaluator<Lhs> LhsEval;
71  typedef typename evaluator<Lhs>::InnerIterator LhsIterator;
72  static void run(const Lhs& lhs, Rhs& other)
73  {
74  LhsEval lhsEval(lhs);
75  for(Index col=0 ; col<other.cols() ; ++col)
76  {
77  for(Index i=lhs.rows()-1 ; i>=0 ; --i)
78  {
79  Scalar tmp = other.coeff(i,col);
80  Scalar l_ii(0);
81  LhsIterator it(lhsEval, i);
82  while(it && it.index()<i)
83  ++it;
84  if(!(Mode & UnitDiag))
85  {
86  eigen_assert(it && it.index()==i);
87  l_ii = it.value();
88  ++it;
89  }
90  else if (it && it.index() == i)
91  ++it;
92  for(; it; ++it)
93  {
94  tmp -= it.value() * other.coeff(it.index(),col);
95  }
96 
97  if (Mode & UnitDiag) other.coeffRef(i,col) = tmp;
98  else other.coeffRef(i,col) = tmp/l_ii;
99  }
100  }
101  }
102 };
103 
104 // forward substitution, col-major
105 template<typename Lhs, typename Rhs, int Mode>
106 struct sparse_solve_triangular_selector<Lhs,Rhs,Mode,Lower,ColMajor>
107 {
108  typedef typename Rhs::Scalar Scalar;
109  typedef evaluator<Lhs> LhsEval;
110  typedef typename evaluator<Lhs>::InnerIterator LhsIterator;
111  static void run(const Lhs& lhs, Rhs& other)
112  {
113  LhsEval lhsEval(lhs);
114  for(Index col=0 ; col<other.cols() ; ++col)
115  {
116  for(Index i=0; i<lhs.cols(); ++i)
117  {
118  Scalar& tmp = other.coeffRef(i,col);
119  if (!numext::is_exactly_zero(tmp)) // optimization when other is actually sparse
120  {
121  LhsIterator it(lhsEval, i);
122  while(it && it.index()<i)
123  ++it;
124  if(!(Mode & UnitDiag))
125  {
126  eigen_assert(it && it.index()==i);
127  tmp /= it.value();
128  }
129  if (it && it.index()==i)
130  ++it;
131  for(; it; ++it)
132  other.coeffRef(it.index(), col) -= tmp * it.value();
133  }
134  }
135  }
136  }
137 };
138 
139 // backward substitution, col-major
140 template<typename Lhs, typename Rhs, int Mode>
141 struct sparse_solve_triangular_selector<Lhs,Rhs,Mode,Upper,ColMajor>
142 {
143  typedef typename Rhs::Scalar Scalar;
144  typedef evaluator<Lhs> LhsEval;
145  typedef typename evaluator<Lhs>::InnerIterator LhsIterator;
146  static void run(const Lhs& lhs, Rhs& other)
147  {
148  LhsEval lhsEval(lhs);
149  for(Index col=0 ; col<other.cols() ; ++col)
150  {
151  for(Index i=lhs.cols()-1; i>=0; --i)
152  {
153  Scalar& tmp = other.coeffRef(i,col);
154  if (!numext::is_exactly_zero(tmp)) // optimization when other is actually sparse
155  {
156  if(!(Mode & UnitDiag))
157  {
158  // TODO replace this by a binary search. make sure the binary search is safe for partially sorted elements
159  LhsIterator it(lhsEval, i);
160  while(it && it.index()!=i)
161  ++it;
162  eigen_assert(it && it.index()==i);
163  other.coeffRef(i,col) /= it.value();
164  }
165  LhsIterator it(lhsEval, i);
166  for(; it && it.index()<i; ++it)
167  other.coeffRef(it.index(), col) -= tmp * it.value();
168  }
169  }
170  }
171  }
172 };
173 
174 } // end namespace internal
175 
176 #ifndef EIGEN_PARSED_BY_DOXYGEN
177 
178 template<typename ExpressionType,unsigned int Mode>
179 template<typename OtherDerived>
180 void TriangularViewImpl<ExpressionType,Mode,Sparse>::solveInPlace(MatrixBase<OtherDerived>& other) const
181 {
182  eigen_assert(derived().cols() == derived().rows() && derived().cols() == other.rows());
183  eigen_assert((!(Mode & ZeroDiag)) && bool(Mode & (Upper|Lower)));
184 
185  enum { copy = internal::traits<OtherDerived>::Flags & RowMajorBit };
186 
187  typedef std::conditional_t<copy,
188  typename internal::plain_matrix_type_column_major<OtherDerived>::type, OtherDerived&> OtherCopy;
189  OtherCopy otherCopy(other.derived());
190 
191  internal::sparse_solve_triangular_selector<ExpressionType, std::remove_reference_t<OtherCopy>, Mode>::run(derived().nestedExpression(), otherCopy);
192 
193  if (copy)
194  other = otherCopy;
195 }
196 #endif
197 
198 // pure sparse path
199 
200 namespace internal {
201 
202 template<typename Lhs, typename Rhs, int Mode,
203  int UpLo = (Mode & Lower)
204  ? Lower
205  : (Mode & Upper)
206  ? Upper
207  : -1,
208  int StorageOrder = int(Lhs::Flags) & (RowMajorBit)>
209 struct sparse_solve_triangular_sparse_selector;
210 
211 // forward substitution, col-major
212 template<typename Lhs, typename Rhs, int Mode, int UpLo>
213 struct sparse_solve_triangular_sparse_selector<Lhs,Rhs,Mode,UpLo,ColMajor>
214 {
215  typedef typename Rhs::Scalar Scalar;
216  typedef typename promote_index_type<typename traits<Lhs>::StorageIndex,
217  typename traits<Rhs>::StorageIndex>::type StorageIndex;
218  static void run(const Lhs& lhs, Rhs& other)
219  {
220  const bool IsLower = (UpLo==Lower);
221  AmbiVector<Scalar,StorageIndex> tempVector(other.rows()*2);
222  tempVector.setBounds(0,other.rows());
223 
224  Rhs res(other.rows(), other.cols());
225  res.reserve(other.nonZeros());
226 
227  for(Index col=0 ; col<other.cols() ; ++col)
228  {
229  // FIXME estimate number of non zeros
230  tempVector.init(.99/*float(other.col(col).nonZeros())/float(other.rows())*/);
231  tempVector.setZero();
232  tempVector.restart();
233  for (typename Rhs::InnerIterator rhsIt(other, col); rhsIt; ++rhsIt)
234  {
235  tempVector.coeffRef(rhsIt.index()) = rhsIt.value();
236  }
237 
238  for(Index i=IsLower?0:lhs.cols()-1;
239  IsLower?i<lhs.cols():i>=0;
240  i+=IsLower?1:-1)
241  {
242  tempVector.restart();
243  Scalar& ci = tempVector.coeffRef(i);
244  if (!numext::is_exactly_zero(ci))
245  {
246  // find
247  typename Lhs::InnerIterator it(lhs, i);
248  if(!(Mode & UnitDiag))
249  {
250  if (IsLower)
251  {
252  eigen_assert(it.index()==i);
253  ci /= it.value();
254  }
255  else
256  ci /= lhs.coeff(i,i);
257  }
258  tempVector.restart();
259  if (IsLower)
260  {
261  if (it.index()==i)
262  ++it;
263  for(; it; ++it)
264  tempVector.coeffRef(it.index()) -= ci * it.value();
265  }
266  else
267  {
268  for(; it && it.index()<i; ++it)
269  tempVector.coeffRef(it.index()) -= ci * it.value();
270  }
271  }
272  }
273 
274 
275 // Index count = 0;
276  // FIXME compute a reference value to filter zeros
277  for (typename AmbiVector<Scalar,StorageIndex>::Iterator it(tempVector/*,1e-12*/); it; ++it)
278  {
279 // ++ count;
280 // std::cerr << "fill " << it.index() << ", " << col << "\n";
281 // std::cout << it.value() << " ";
282  // FIXME use insertBack
283  res.insert(it.index(), col) = it.value();
284  }
285 // std::cout << "tempVector.nonZeros() == " << int(count) << " / " << (other.rows()) << "\n";
286  }
287  res.finalize();
288  other = res.markAsRValue();
289  }
290 };
291 
292 } // end namespace internal
293 
294 #ifndef EIGEN_PARSED_BY_DOXYGEN
295 template<typename ExpressionType,unsigned int Mode>
296 template<typename OtherDerived>
297 void TriangularViewImpl<ExpressionType,Mode,Sparse>::solveInPlace(SparseMatrixBase<OtherDerived>& other) const
298 {
299  eigen_assert(derived().cols() == derived().rows() && derived().cols() == other.rows());
300  eigen_assert( (!(Mode & ZeroDiag)) && bool(Mode & (Upper|Lower)));
301 
302 // enum { copy = internal::traits<OtherDerived>::Flags & RowMajorBit };
303 
304 // typedef std::conditional_t<copy,
305 // typename internal::plain_matrix_type_column_major<OtherDerived>::type, OtherDerived&> OtherCopy;
306 // OtherCopy otherCopy(other.derived());
307 
308  internal::sparse_solve_triangular_sparse_selector<ExpressionType, OtherDerived, Mode>::run(derived().nestedExpression(), other.derived());
309 
310 // if (copy)
311 // other = otherCopy;
312 }
313 #endif
314 
315 } // end namespace Eigen
316 
317 #endif // EIGEN_SPARSETRIANGULARSOLVER_H
ColXpr col(Index i)
This is the const version of col().
#define eigen_assert(x)
Definition: Macros.h:902
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
@ UnitDiag
Definition: Constants.h:215
@ ZeroDiag
Definition: Constants.h:217
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ ColMajor
Definition: Constants.h:321
@ RowMajor
Definition: Constants.h:323
const unsigned int RowMajorBit
Definition: Constants.h:68
bool is_exactly_zero(const X &x)
Definition: Meta.h:475
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82