compbio
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 namespace Eigen {
14 
35 template <typename _Scalar>
37 {
38  typedef _Scalar Scalar;
40  public:
41  typedef typename Vector::StorageIndex StorageIndex;
42  enum {
43  ColsAtCompileTime = Dynamic,
44  MaxColsAtCompileTime = Dynamic
45  };
46 
47  DiagonalPreconditioner() : m_isInitialized(false) {}
48 
49  template<typename MatType>
50  explicit DiagonalPreconditioner(const MatType& mat) : m_invdiag(mat.cols())
51  {
52  compute(mat);
53  }
54 
55  Index rows() const { return m_invdiag.size(); }
56  Index cols() const { return m_invdiag.size(); }
57 
58  template<typename MatType>
59  DiagonalPreconditioner& analyzePattern(const MatType& )
60  {
61  return *this;
62  }
63 
64  template<typename MatType>
65  DiagonalPreconditioner& factorize(const MatType& mat)
66  {
67  m_invdiag.resize(mat.cols());
68  for(int j=0; j<mat.outerSize(); ++j)
69  {
70  typename MatType::InnerIterator it(mat,j);
71  while(it && it.index()!=j) ++it;
72  if(it && it.index()==j && it.value()!=Scalar(0))
73  m_invdiag(j) = Scalar(1)/it.value();
74  else
75  m_invdiag(j) = Scalar(1);
76  }
77  m_isInitialized = true;
78  return *this;
79  }
80 
81  template<typename MatType>
82  DiagonalPreconditioner& compute(const MatType& mat)
83  {
84  return factorize(mat);
85  }
86 
88  template<typename Rhs, typename Dest>
89  void _solve_impl(const Rhs& b, Dest& x) const
90  {
91  x = m_invdiag.array() * b.array() ;
92  }
93 
94  template<typename Rhs> inline const Solve<DiagonalPreconditioner, Rhs>
95  solve(const MatrixBase<Rhs>& b) const
96  {
97  eigen_assert(m_isInitialized && "DiagonalPreconditioner is not initialized.");
98  eigen_assert(m_invdiag.size()==b.rows()
99  && "DiagonalPreconditioner::solve(): invalid number of rows of the right hand side matrix b");
100  return Solve<DiagonalPreconditioner, Rhs>(*this, b.derived());
101  }
102 
103  ComputationInfo info() { return Success; }
104 
105  protected:
106  Vector m_invdiag;
107  bool m_isInitialized;
108 };
109 
127 template <typename _Scalar>
129 {
130  typedef _Scalar Scalar;
131  typedef typename NumTraits<Scalar>::Real RealScalar;
133  using Base::m_invdiag;
134  public:
135 
137 
138  template<typename MatType>
139  explicit LeastSquareDiagonalPreconditioner(const MatType& mat) : Base()
140  {
141  compute(mat);
142  }
143 
144  template<typename MatType>
145  LeastSquareDiagonalPreconditioner& analyzePattern(const MatType& )
146  {
147  return *this;
148  }
149 
150  template<typename MatType>
151  LeastSquareDiagonalPreconditioner& factorize(const MatType& mat)
152  {
153  // Compute the inverse squared-norm of each column of mat
154  m_invdiag.resize(mat.cols());
155  for(Index j=0; j<mat.outerSize(); ++j)
156  {
157  RealScalar sum = mat.innerVector(j).squaredNorm();
158  if(sum>0)
159  m_invdiag(j) = RealScalar(1)/sum;
160  else
161  m_invdiag(j) = RealScalar(1);
162  }
163  Base::m_isInitialized = true;
164  return *this;
165  }
166 
167  template<typename MatType>
168  LeastSquareDiagonalPreconditioner& compute(const MatType& mat)
169  {
170  return factorize(mat);
171  }
172 
173  ComputationInfo info() { return Success; }
174 
175  protected:
176 };
177 
186 {
187  public:
188 
190 
191  template<typename MatrixType>
192  explicit IdentityPreconditioner(const MatrixType& ) {}
193 
194  template<typename MatrixType>
195  IdentityPreconditioner& analyzePattern(const MatrixType& ) { return *this; }
196 
197  template<typename MatrixType>
198  IdentityPreconditioner& factorize(const MatrixType& ) { return *this; }
199 
200  template<typename MatrixType>
201  IdentityPreconditioner& compute(const MatrixType& ) { return *this; }
202 
203  template<typename Rhs>
204  inline const Rhs& solve(const Rhs& b) const { return b; }
205 
206  ComputationInfo info() { return Success; }
207 };
208 
209 } // end namespace Eigen
210 
211 #endif // EIGEN_BASIC_PRECONDITIONERS_H
A preconditioner based on the digonal entries.
Definition: BasicPreconditioners.h:36
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:85
Holds information about the various numeric (i.e.
Definition: NumTraits.h:150
Jacobi preconditioner for LeastSquaresConjugateGradient.
Definition: BasicPreconditioners.h:128
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Resizes *this to a rows x cols matrix.
Definition: PlainObjectBase.h:273
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
Computation was successful.
Definition: Constants.h:432
A naive preconditioner which approximates any matrix as the identity matrix.
Definition: BasicPreconditioners.h:185
const int Dynamic
This value means that a positive quantity (e.g., a size) is not known at compile-time, and that instead the value is stored in some runtime variable.
Definition: Constants.h:21
Pseudo expression representing a solving operation.
Definition: Solve.h:62
ComputationInfo
Enum for reporting the status of a computation.
Definition: Constants.h:430
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48