$darkmode
Eigen  5.0.1-dev
BasicPreconditioners.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2011-2014 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_BASIC_PRECONDITIONERS_H
11 #define EIGEN_BASIC_PRECONDITIONERS_H
12 
13 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
38 template <typename Scalar_>
40  typedef Scalar_ Scalar;
42 
43  public:
44  typedef typename Vector::StorageIndex StorageIndex;
45  enum { ColsAtCompileTime = Dynamic, MaxColsAtCompileTime = Dynamic };
46 
47  DiagonalPreconditioner() : m_isInitialized(false) {}
48 
49  template <typename MatType>
50  explicit DiagonalPreconditioner(const MatType& mat) : m_invdiag(mat.cols()) {
51  compute(mat);
52  }
53 
54  constexpr Index rows() const noexcept { return m_invdiag.size(); }
55  constexpr Index cols() const noexcept { return m_invdiag.size(); }
56 
57  template <typename MatType>
58  DiagonalPreconditioner& analyzePattern(const MatType&) {
59  return *this;
60  }
61 
62  template <typename MatType>
63  DiagonalPreconditioner& factorize(const MatType& mat) {
64  m_invdiag.resize(mat.cols());
65  for (int j = 0; j < mat.outerSize(); ++j) {
66  typename MatType::InnerIterator it(mat, j);
67  while (it && it.index() != j) ++it;
68  if (it && it.index() == j && it.value() != Scalar(0))
69  m_invdiag(j) = Scalar(1) / it.value();
70  else
71  m_invdiag(j) = Scalar(1);
72  }
73  m_isInitialized = true;
74  return *this;
75  }
76 
77  template <typename MatType>
78  DiagonalPreconditioner& compute(const MatType& mat) {
79  return factorize(mat);
80  }
81 
83  template <typename Rhs, typename Dest>
84  void _solve_impl(const Rhs& b, Dest& x) const {
85  x = m_invdiag.array() * b.array();
86  }
87 
88  template <typename Rhs>
89  inline const Solve<DiagonalPreconditioner, Rhs> solve(const MatrixBase<Rhs>& b) const {
90  eigen_assert(m_isInitialized && "DiagonalPreconditioner is not initialized.");
91  eigen_assert(m_invdiag.size() == b.rows() &&
92  "DiagonalPreconditioner::solve(): invalid number of rows of the right hand side matrix b");
94  }
95 
96  ComputationInfo info() { return Success; }
97 
98  protected:
99  Vector m_invdiag;
100  bool m_isInitialized;
101 };
102 
120 template <typename Scalar_>
122  typedef Scalar_ Scalar;
123  typedef typename NumTraits<Scalar>::Real RealScalar;
125  using Base::m_invdiag;
126 
127  public:
129 
130  template <typename MatType>
131  explicit LeastSquareDiagonalPreconditioner(const MatType& mat) : Base() {
132  compute(mat);
133  }
134 
135  template <typename MatType>
136  LeastSquareDiagonalPreconditioner& analyzePattern(const MatType&) {
137  return *this;
138  }
139 
140  template <typename MatType>
141  LeastSquareDiagonalPreconditioner& factorize(const MatType& mat) {
142  // Compute the inverse squared-norm of each column of mat
143  m_invdiag.resize(mat.cols());
144  if (MatType::IsRowMajor) {
145  m_invdiag.setZero();
146  for (Index j = 0; j < mat.outerSize(); ++j) {
147  for (typename MatType::InnerIterator it(mat, j); it; ++it) m_invdiag(it.index()) += numext::abs2(it.value());
148  }
149  for (Index j = 0; j < mat.cols(); ++j)
150  if (numext::real(m_invdiag(j)) > RealScalar(0)) m_invdiag(j) = RealScalar(1) / numext::real(m_invdiag(j));
151  } else {
152  for (Index j = 0; j < mat.outerSize(); ++j) {
153  RealScalar sum = mat.col(j).squaredNorm();
154  if (sum > RealScalar(0))
155  m_invdiag(j) = RealScalar(1) / sum;
156  else
157  m_invdiag(j) = RealScalar(1);
158  }
159  }
160  Base::m_isInitialized = true;
161  return *this;
162  }
163 
164  template <typename MatType>
165  LeastSquareDiagonalPreconditioner& compute(const MatType& mat) {
166  return factorize(mat);
167  }
168 
169  ComputationInfo info() { return Success; }
170 
171  protected:
172 };
173 
182  public:
184 
185  template <typename MatrixType>
186  explicit IdentityPreconditioner(const MatrixType&) {}
187 
188  template <typename MatrixType>
189  IdentityPreconditioner& analyzePattern(const MatrixType&) {
190  return *this;
191  }
192 
193  template <typename MatrixType>
194  IdentityPreconditioner& factorize(const MatrixType&) {
195  return *this;
196  }
197 
198  template <typename MatrixType>
199  IdentityPreconditioner& compute(const MatrixType&) {
200  return *this;
201  }
202 
203  template <typename Rhs>
204  inline const Rhs& solve(const Rhs& b) const {
205  return b;
206  }
207 
208  ComputationInfo info() { return Success; }
209 };
210 
211 } // end namespace Eigen
212 
213 #endif // EIGEN_BASIC_PRECONDITIONERS_H
A preconditioner based on the digonal entries.
Definition: BasicPreconditioners.h:39
constexpr Derived & derived()
Definition: EigenBase.h:49
constexpr void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:268
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:232
Jacobi preconditioner for LeastSquaresConjugateGradient.
Definition: BasicPreconditioners.h:121
Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:567
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Definition: Constants.h:440
A naive preconditioner which approximates any matrix as the identity matrix.
Definition: BasicPreconditioners.h:181
const int Dynamic
Definition: Constants.h:25
Pseudo expression representing a solving operation.
Definition: Solve.h:62
constexpr Index rows() const noexcept
Definition: EigenBase.h:59
ComputationInfo
Definition: Constants.h:438
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52