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