SkylineInplaceLU.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 Guillaume Saupin <guillaume.saupin@cea.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_SKYLINEINPLACELU_H
11 #define EIGEN_SKYLINEINPLACELU_H
12 
13 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
26 template<typename MatrixType>
28 protected:
29  typedef typename MatrixType::Scalar Scalar;
30  typedef typename MatrixType::Index Index;
31 
33 
34 public:
35 
38  SkylineInplaceLU(MatrixType& matrix, int flags = 0)
39  : /*m_matrix(matrix.rows(), matrix.cols()),*/ m_flags(flags), m_status(0), m_lu(matrix) {
40  m_precision = RealScalar(0.1) * Eigen::dummy_precision<RealScalar > ();
41  m_lu.IsRowMajor ? computeRowMajor() : compute();
42  }
43 
55  m_precision = v;
56  }
57 
62  return m_precision;
63  }
64 
73  void setFlags(int f) {
74  m_flags = f;
75  }
76 
78  int flags() const {
79  return m_flags;
80  }
81 
82  void setOrderingMethod(int m) {
83  m_flags = m;
84  }
85 
86  int orderingMethod() const {
87  return m_flags;
88  }
89 
91  void compute();
92  void computeRowMajor();
93 
95  //inline const MatrixType& matrixL() const { return m_matrixL; }
96 
98  //inline const MatrixType& matrixU() const { return m_matrixU; }
99 
100  template<typename BDerived, typename XDerived>
102  const int transposed = 0) const;
103 
105  inline bool succeeded(void) const {
106  return m_succeeded;
107  }
108 
109 protected:
111  int m_flags;
112  mutable int m_status;
115 };
116 
120 template<typename MatrixType>
121 //template<typename Scalar_>
123  const size_t rows = m_lu.rows();
124  const size_t cols = m_lu.cols();
125 
126  eigen_assert(rows == cols && "We do not (yet) support rectangular LU.");
127  eigen_assert(!m_lu.IsRowMajor && "LU decomposition does not work with rowMajor Storage");
128 
129  for (Index row = 0; row < rows; row++) {
130  const double pivot = m_lu.coeffDiag(row);
131 
132  //Lower matrix Columns update
133  const Index& col = row;
134  for (typename MatrixType::InnerLowerIterator lIt(m_lu, col); lIt; ++lIt) {
135  lIt.valueRef() /= pivot;
136  }
137 
138  //Upper matrix update -> contiguous memory access
139  typename MatrixType::InnerLowerIterator lIt(m_lu, col);
140  for (Index rrow = row + 1; rrow < m_lu.rows(); rrow++) {
141  typename MatrixType::InnerUpperIterator uItPivot(m_lu, row);
142  typename MatrixType::InnerUpperIterator uIt(m_lu, rrow);
143  const double coef = lIt.value();
144 
145  uItPivot += (rrow - row - 1);
146 
147  //update upper part -> contiguous memory access
148  for (++uItPivot; uIt && uItPivot;) {
149  uIt.valueRef() -= uItPivot.value() * coef;
150 
151  ++uIt;
152  ++uItPivot;
153  }
154  ++lIt;
155  }
156 
157  //Upper matrix update -> non contiguous memory access
158  typename MatrixType::InnerLowerIterator lIt3(m_lu, col);
159  for (Index rrow = row + 1; rrow < m_lu.rows(); rrow++) {
160  typename MatrixType::InnerUpperIterator uItPivot(m_lu, row);
161  const double coef = lIt3.value();
162 
163  //update lower part -> non contiguous memory access
164  for (Index i = 0; i < rrow - row - 1; i++) {
165  m_lu.coeffRefLower(rrow, row + i + 1) -= uItPivot.value() * coef;
166  ++uItPivot;
167  }
168  ++lIt3;
169  }
170  //update diag -> contiguous
171  typename MatrixType::InnerLowerIterator lIt2(m_lu, col);
172  for (Index rrow = row + 1; rrow < m_lu.rows(); rrow++) {
173 
174  typename MatrixType::InnerUpperIterator uItPivot(m_lu, row);
175  typename MatrixType::InnerUpperIterator uIt(m_lu, rrow);
176  const double coef = lIt2.value();
177 
178  uItPivot += (rrow - row - 1);
179  m_lu.coeffRefDiag(rrow) -= uItPivot.value() * coef;
180  ++lIt2;
181  }
182  }
183 }
184 
185 template<typename MatrixType>
187  const size_t rows = m_lu.rows();
188  const size_t cols = m_lu.cols();
189 
190  eigen_assert(rows == cols && "We do not (yet) support rectangular LU.");
191  eigen_assert(m_lu.IsRowMajor && "You're trying to apply rowMajor decomposition on a ColMajor matrix !");
192 
193  for (Index row = 0; row < rows; row++) {
194  typename MatrixType::InnerLowerIterator llIt(m_lu, row);
195 
196 
197  for (Index col = llIt.col(); col < row; col++) {
198  if (m_lu.coeffExistLower(row, col)) {
199  const double diag = m_lu.coeffDiag(col);
200 
201  typename MatrixType::InnerLowerIterator lIt(m_lu, row);
202  typename MatrixType::InnerUpperIterator uIt(m_lu, col);
203 
204 
205  const Index offset = lIt.col() - uIt.row();
206 
207 
208  Index stop = offset > 0 ? col - lIt.col() : col - uIt.row();
209 
210  //#define VECTORIZE
211 #ifdef VECTORIZE
212  Map<VectorXd > rowVal(lIt.valuePtr() + (offset > 0 ? 0 : -offset), stop);
213  Map<VectorXd > colVal(uIt.valuePtr() + (offset > 0 ? offset : 0), stop);
214 
215 
216  Scalar newCoeff = m_lu.coeffLower(row, col) - rowVal.dot(colVal);
217 #else
218  if (offset > 0) //Skip zero value of lIt
219  uIt += offset;
220  else //Skip zero values of uIt
221  lIt += -offset;
222  Scalar newCoeff = m_lu.coeffLower(row, col);
223 
224  for (Index k = 0; k < stop; ++k) {
225  const Scalar tmp = newCoeff;
226  newCoeff = tmp - lIt.value() * uIt.value();
227  ++lIt;
228  ++uIt;
229  }
230 #endif
231 
232  m_lu.coeffRefLower(row, col) = newCoeff / diag;
233  }
234  }
235 
236  //Upper matrix update
237  const Index col = row;
238  typename MatrixType::InnerUpperIterator uuIt(m_lu, col);
239  for (Index rrow = uuIt.row(); rrow < col; rrow++) {
240 
241  typename MatrixType::InnerLowerIterator lIt(m_lu, rrow);
242  typename MatrixType::InnerUpperIterator uIt(m_lu, col);
243  const Index offset = lIt.col() - uIt.row();
244 
245  Index stop = offset > 0 ? rrow - lIt.col() : rrow - uIt.row();
246 
247 #ifdef VECTORIZE
248  Map<VectorXd > rowVal(lIt.valuePtr() + (offset > 0 ? 0 : -offset), stop);
249  Map<VectorXd > colVal(uIt.valuePtr() + (offset > 0 ? offset : 0), stop);
250 
251  Scalar newCoeff = m_lu.coeffUpper(rrow, col) - rowVal.dot(colVal);
252 #else
253  if (offset > 0) //Skip zero value of lIt
254  uIt += offset;
255  else //Skip zero values of uIt
256  lIt += -offset;
257  Scalar newCoeff = m_lu.coeffUpper(rrow, col);
258  for (Index k = 0; k < stop; ++k) {
259  const Scalar tmp = newCoeff;
260  newCoeff = tmp - lIt.value() * uIt.value();
261 
262  ++lIt;
263  ++uIt;
264  }
265 #endif
266  m_lu.coeffRefUpper(rrow, col) = newCoeff;
267  }
268 
269 
270  //Diag matrix update
271  typename MatrixType::InnerLowerIterator lIt(m_lu, row);
272  typename MatrixType::InnerUpperIterator uIt(m_lu, row);
273 
274  const Index offset = lIt.col() - uIt.row();
275 
276 
277  Index stop = offset > 0 ? lIt.size() : uIt.size();
278 #ifdef VECTORIZE
279  Map<VectorXd > rowVal(lIt.valuePtr() + (offset > 0 ? 0 : -offset), stop);
280  Map<VectorXd > colVal(uIt.valuePtr() + (offset > 0 ? offset : 0), stop);
281  Scalar newCoeff = m_lu.coeffDiag(row) - rowVal.dot(colVal);
282 #else
283  if (offset > 0) //Skip zero value of lIt
284  uIt += offset;
285  else //Skip zero values of uIt
286  lIt += -offset;
287  Scalar newCoeff = m_lu.coeffDiag(row);
288  for (Index k = 0; k < stop; ++k) {
289  const Scalar tmp = newCoeff;
290  newCoeff = tmp - lIt.value() * uIt.value();
291  ++lIt;
292  ++uIt;
293  }
294 #endif
295  m_lu.coeffRefDiag(row) = newCoeff;
296  }
297 }
298 
307 template<typename MatrixType>
308 template<typename BDerived, typename XDerived>
310  const size_t rows = m_lu.rows();
311  const size_t cols = m_lu.cols();
312 
313 
314  for (Index row = 0; row < rows; row++) {
315  x->coeffRef(row) = b.coeff(row);
316  Scalar newVal = x->coeff(row);
317  typename MatrixType::InnerLowerIterator lIt(m_lu, row);
318 
319  Index col = lIt.col();
320  while (lIt.col() < row) {
321 
322  newVal -= x->coeff(col++) * lIt.value();
323  ++lIt;
324  }
325 
326  x->coeffRef(row) = newVal;
327  }
328 
329 
330  for (Index col = rows - 1; col > 0; col--) {
331  x->coeffRef(col) = x->coeff(col) / m_lu.coeffDiag(col);
332 
333  const Scalar x_col = x->coeff(col);
334 
335  typename MatrixType::InnerUpperIterator uIt(m_lu, col);
336  uIt += uIt.size()-1;
337 
338 
339  while (uIt) {
340  x->coeffRef(uIt.row()) -= x_col * uIt.value();
341  //TODO : introduce --operator
342  uIt += -1;
343  }
344 
345 
346  }
347  x->coeffRef(0) = x->coeff(0) / m_lu.coeffDiag(0);
348 
349  return true;
350 }
351 
352 } // end namespace Eigen
353 
354 #endif // EIGEN_SKYLINEINPLACELU_H
Matrix3f m
Array< int, Dynamic, 1 > v
int i
RowXpr row(Index i) const
ColXpr col(Index i) const
#define eigen_assert(x)
Matrix< float, 1, Dynamic > MatrixType
Inplace LU decomposition of a skyline matrix and associated features.
RealScalar precision() const
void setPrecision(RealScalar v)
MatrixType::Scalar Scalar
bool solve(const MatrixBase< BDerived > &b, MatrixBase< XDerived > *x, const int transposed=0) const
MatrixType::Index Index
bool succeeded(void) const
SkylineInplaceLU(MatrixType &matrix, int flags=0)
NumTraits< typename MatrixType::Scalar >::Real RealScalar
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend
Derived::Index cols
Derived::Index rows