BasicPreconditioners.h 6.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226
  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. #ifndef EIGEN_BASIC_PRECONDITIONERS_H
  10. #define EIGEN_BASIC_PRECONDITIONERS_H
  11. namespace Eigen {
  12. /** \ingroup IterativeLinearSolvers_Module
  13. * \brief A preconditioner based on the digonal entries
  14. *
  15. * This class allows to approximately solve for A.x = b problems assuming A is a diagonal matrix.
  16. * In other words, this preconditioner neglects all off diagonal entries and, in Eigen's language, solves for:
  17. \code
  18. A.diagonal().asDiagonal() . x = b
  19. \endcode
  20. *
  21. * \tparam _Scalar the type of the scalar.
  22. *
  23. * \implsparsesolverconcept
  24. *
  25. * This preconditioner is suitable for both selfadjoint and general problems.
  26. * The diagonal entries are pre-inverted and stored into a dense vector.
  27. *
  28. * \note A variant that has yet to be implemented would attempt to preserve the norm of each column.
  29. *
  30. * \sa class LeastSquareDiagonalPreconditioner, class ConjugateGradient
  31. */
  32. template <typename _Scalar>
  33. class DiagonalPreconditioner
  34. {
  35. typedef _Scalar Scalar;
  36. typedef Matrix<Scalar,Dynamic,1> Vector;
  37. public:
  38. typedef typename Vector::StorageIndex StorageIndex;
  39. enum {
  40. ColsAtCompileTime = Dynamic,
  41. MaxColsAtCompileTime = Dynamic
  42. };
  43. DiagonalPreconditioner() : m_isInitialized(false) {}
  44. template<typename MatType>
  45. explicit DiagonalPreconditioner(const MatType& mat) : m_invdiag(mat.cols())
  46. {
  47. compute(mat);
  48. }
  49. Index rows() const { return m_invdiag.size(); }
  50. Index cols() const { return m_invdiag.size(); }
  51. template<typename MatType>
  52. DiagonalPreconditioner& analyzePattern(const MatType& )
  53. {
  54. return *this;
  55. }
  56. template<typename MatType>
  57. DiagonalPreconditioner& factorize(const MatType& mat)
  58. {
  59. m_invdiag.resize(mat.cols());
  60. for(int j=0; j<mat.outerSize(); ++j)
  61. {
  62. typename MatType::InnerIterator it(mat,j);
  63. while(it && it.index()!=j) ++it;
  64. if(it && it.index()==j && it.value()!=Scalar(0))
  65. m_invdiag(j) = Scalar(1)/it.value();
  66. else
  67. m_invdiag(j) = Scalar(1);
  68. }
  69. m_isInitialized = true;
  70. return *this;
  71. }
  72. template<typename MatType>
  73. DiagonalPreconditioner& compute(const MatType& mat)
  74. {
  75. return factorize(mat);
  76. }
  77. /** \internal */
  78. template<typename Rhs, typename Dest>
  79. void _solve_impl(const Rhs& b, Dest& x) const
  80. {
  81. x = m_invdiag.array() * b.array() ;
  82. }
  83. template<typename Rhs> inline const Solve<DiagonalPreconditioner, Rhs>
  84. solve(const MatrixBase<Rhs>& b) const
  85. {
  86. eigen_assert(m_isInitialized && "DiagonalPreconditioner is not initialized.");
  87. eigen_assert(m_invdiag.size()==b.rows()
  88. && "DiagonalPreconditioner::solve(): invalid number of rows of the right hand side matrix b");
  89. return Solve<DiagonalPreconditioner, Rhs>(*this, b.derived());
  90. }
  91. ComputationInfo info() { return Success; }
  92. protected:
  93. Vector m_invdiag;
  94. bool m_isInitialized;
  95. };
  96. /** \ingroup IterativeLinearSolvers_Module
  97. * \brief Jacobi preconditioner for LeastSquaresConjugateGradient
  98. *
  99. * This class allows to approximately solve for A' A x = A' b problems assuming A' A is a diagonal matrix.
  100. * In other words, this preconditioner neglects all off diagonal entries and, in Eigen's language, solves for:
  101. \code
  102. (A.adjoint() * A).diagonal().asDiagonal() * x = b
  103. \endcode
  104. *
  105. * \tparam _Scalar the type of the scalar.
  106. *
  107. * \implsparsesolverconcept
  108. *
  109. * The diagonal entries are pre-inverted and stored into a dense vector.
  110. *
  111. * \sa class LeastSquaresConjugateGradient, class DiagonalPreconditioner
  112. */
  113. template <typename _Scalar>
  114. class LeastSquareDiagonalPreconditioner : public DiagonalPreconditioner<_Scalar>
  115. {
  116. typedef _Scalar Scalar;
  117. typedef typename NumTraits<Scalar>::Real RealScalar;
  118. typedef DiagonalPreconditioner<_Scalar> Base;
  119. using Base::m_invdiag;
  120. public:
  121. LeastSquareDiagonalPreconditioner() : Base() {}
  122. template<typename MatType>
  123. explicit LeastSquareDiagonalPreconditioner(const MatType& mat) : Base()
  124. {
  125. compute(mat);
  126. }
  127. template<typename MatType>
  128. LeastSquareDiagonalPreconditioner& analyzePattern(const MatType& )
  129. {
  130. return *this;
  131. }
  132. template<typename MatType>
  133. LeastSquareDiagonalPreconditioner& factorize(const MatType& mat)
  134. {
  135. // Compute the inverse squared-norm of each column of mat
  136. m_invdiag.resize(mat.cols());
  137. if(MatType::IsRowMajor)
  138. {
  139. m_invdiag.setZero();
  140. for(Index j=0; j<mat.outerSize(); ++j)
  141. {
  142. for(typename MatType::InnerIterator it(mat,j); it; ++it)
  143. m_invdiag(it.index()) += numext::abs2(it.value());
  144. }
  145. for(Index j=0; j<mat.cols(); ++j)
  146. if(numext::real(m_invdiag(j))>RealScalar(0))
  147. m_invdiag(j) = RealScalar(1)/numext::real(m_invdiag(j));
  148. }
  149. else
  150. {
  151. for(Index j=0; j<mat.outerSize(); ++j)
  152. {
  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. template<typename MatType>
  164. LeastSquareDiagonalPreconditioner& compute(const MatType& mat)
  165. {
  166. return factorize(mat);
  167. }
  168. ComputationInfo info() { return Success; }
  169. protected:
  170. };
  171. /** \ingroup IterativeLinearSolvers_Module
  172. * \brief A naive preconditioner which approximates any matrix as the identity matrix
  173. *
  174. * \implsparsesolverconcept
  175. *
  176. * \sa class DiagonalPreconditioner
  177. */
  178. class IdentityPreconditioner
  179. {
  180. public:
  181. IdentityPreconditioner() {}
  182. template<typename MatrixType>
  183. explicit IdentityPreconditioner(const MatrixType& ) {}
  184. template<typename MatrixType>
  185. IdentityPreconditioner& analyzePattern(const MatrixType& ) { return *this; }
  186. template<typename MatrixType>
  187. IdentityPreconditioner& factorize(const MatrixType& ) { return *this; }
  188. template<typename MatrixType>
  189. IdentityPreconditioner& compute(const MatrixType& ) { return *this; }
  190. template<typename Rhs>
  191. inline const Rhs& solve(const Rhs& b) const { return b; }
  192. ComputationInfo info() { return Success; }
  193. };
  194. } // end namespace Eigen
  195. #endif // EIGEN_BASIC_PRECONDITIONERS_H