Eigen-unsupported  5.0.1-dev+7c7d8473
 
Loading...
Searching...
No Matches
IncompleteLU.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2011 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_INCOMPLETE_LU_H
11#define EIGEN_INCOMPLETE_LU_H
12
13// IWYU pragma: private
14#include "./InternalHeaderCheck.h"
15
16namespace Eigen {
17
18template <typename Scalar_>
19class IncompleteLU : public SparseSolverBase<IncompleteLU<Scalar_> > {
20 protected:
22 using Base::m_isInitialized;
23
24 typedef Scalar_ Scalar;
25 typedef Matrix<Scalar, Dynamic, 1> Vector;
26 typedef typename Vector::Index Index;
27 typedef SparseMatrix<Scalar, RowMajor> FactorType;
28
29 public:
30 typedef Matrix<Scalar, Dynamic, Dynamic> MatrixType;
31
32 IncompleteLU() {}
33
34 template <typename MatrixType>
35 IncompleteLU(const MatrixType& mat) {
36 compute(mat);
37 }
38
39 Index rows() const { return m_lu.rows(); }
40 Index cols() const { return m_lu.cols(); }
41
42 template <typename MatrixType>
43 IncompleteLU& compute(const MatrixType& mat) {
44 m_lu = mat;
45 int size = mat.cols();
46 Vector diag(size);
47 for (int i = 0; i < size; ++i) {
48 typename FactorType::InnerIterator k_it(m_lu, i);
49 for (; k_it && k_it.index() < i; ++k_it) {
50 int k = k_it.index();
51 k_it.valueRef() /= diag(k);
52
53 typename FactorType::InnerIterator j_it(k_it);
54 typename FactorType::InnerIterator kj_it(m_lu, k);
55 while (kj_it && kj_it.index() <= k) ++kj_it;
56 for (++j_it; j_it;) {
57 if (kj_it.index() == j_it.index()) {
58 j_it.valueRef() -= k_it.value() * kj_it.value();
59 ++j_it;
60 ++kj_it;
61 } else if (kj_it.index() < j_it.index())
62 ++kj_it;
63 else
64 ++j_it;
65 }
66 }
67 if (k_it && k_it.index() == i)
68 diag(i) = k_it.value();
69 else
70 diag(i) = 1;
71 }
72 m_isInitialized = true;
73 return *this;
74 }
75
76 template <typename Rhs, typename Dest>
77 void _solve_impl(const Rhs& b, Dest& x) const {
78 x = m_lu.template triangularView<UnitLower>().solve(b);
79 x = m_lu.template triangularView<Upper>().solve(x);
80 }
81
82 protected:
83 FactorType m_lu;
84};
85
86} // end namespace Eigen
87
88#endif // EIGEN_INCOMPLETE_LU_H
Namespace containing all symbols from the Eigen library.