Eigen-unsupported  3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
IterationController.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
5 
6 /* NOTE The class IterationController has been adapted from the iteration
7  * class of the GMM++ and ITL libraries.
8  */
9 
10 //=======================================================================
11 // Copyright (C) 1997-2001
12 // Authors: Andrew Lumsdaine <lums@osl.iu.edu>
13 // Lie-Quan Lee <llee@osl.iu.edu>
14 //
15 // This file is part of the Iterative Template Library
16 //
17 // You should have received a copy of the License Agreement for the
18 // Iterative Template Library along with the software; see the
19 // file LICENSE.
20 //
21 // Permission to modify the code and to distribute modified code is
22 // granted, provided the text of this NOTICE is retained, a notice that
23 // the code was modified is included with the above COPYRIGHT NOTICE and
24 // with the COPYRIGHT NOTICE in the LICENSE file, and that the LICENSE
25 // file is distributed with the modified code.
26 //
27 // LICENSOR MAKES NO REPRESENTATIONS OR WARRANTIES, EXPRESS OR IMPLIED.
28 // By way of example, but not limitation, Licensor MAKES NO
29 // REPRESENTATIONS OR WARRANTIES OF MERCHANTABILITY OR FITNESS FOR ANY
30 // PARTICULAR PURPOSE OR THAT THE USE OF THE LICENSED SOFTWARE COMPONENTS
31 // OR DOCUMENTATION WILL NOT INFRINGE ANY PATENTS, COPYRIGHTS, TRADEMARKS
32 // OR OTHER RIGHTS.
33 //=======================================================================
34 
35 //========================================================================
36 //
37 // Copyright (C) 2002-2007 Yves Renard
38 //
39 // This file is a part of GETFEM++
40 //
41 // Getfem++ is free software; you can redistribute it and/or modify
42 // it under the terms of the GNU Lesser General Public License as
43 // published by the Free Software Foundation; version 2.1 of the License.
44 //
45 // This program is distributed in the hope that it will be useful,
46 // but WITHOUT ANY WARRANTY; without even the implied warranty of
47 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
48 // GNU Lesser General Public License for more details.
49 // You should have received a copy of the GNU Lesser General Public
50 // License along with this program; if not, write to the Free Software
51 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301,
52 // USA.
53 //
54 //========================================================================
55 
56 #include "../../../../Eigen/src/Core/util/NonMPL2.h"
57 
58 #ifndef EIGEN_ITERATION_CONTROLLER_H
59 #define EIGEN_ITERATION_CONTROLLER_H
60 
61 #include "./InternalHeaderCheck.h"
62 
63 namespace Eigen {
64 
74 {
75  protected :
76  double m_rhsn;
77  size_t m_maxiter;
78  int m_noise;
79  double m_resmax;
80  double m_resminreach, m_resadd;
81  size_t m_nit;
82  double m_res;
83  bool m_written;
84  void (*m_callback)(const IterationController&);
85  public :
86 
87  void init()
88  {
89  m_nit = 0; m_res = 0.0; m_written = false;
90  m_resminreach = 1E50; m_resadd = 0.0;
91  m_callback = 0;
92  }
93 
94  IterationController(double r = 1.0E-8, int noi = 0, size_t mit = size_t(-1))
95  : m_rhsn(1.0), m_maxiter(mit), m_noise(noi), m_resmax(r) { init(); }
96 
97  void operator ++(int) { m_nit++; m_written = false; m_resadd += m_res; }
98  void operator ++() { (*this)++; }
99 
100  bool first() { return m_nit == 0; }
101 
102  /* get/set the "noisyness" (verbosity) of the solvers */
103  int noiseLevel() const { return m_noise; }
104  void setNoiseLevel(int n) { m_noise = n; }
105  void reduceNoiseLevel() { if (m_noise > 0) m_noise--; }
106 
107  double maxResidual() const { return m_resmax; }
108  void setMaxResidual(double r) { m_resmax = r; }
109 
110  double residual() const { return m_res; }
111 
112  /* change the user-definable callback, called after each iteration */
113  void setCallback(void (*t)(const IterationController&))
114  {
115  m_callback = t;
116  }
117 
118  size_t iteration() const { return m_nit; }
119  void setIteration(size_t i) { m_nit = i; }
120 
121  size_t maxIterarions() const { return m_maxiter; }
122  void setMaxIterations(size_t i) { m_maxiter = i; }
123 
124  double rhsNorm() const { return m_rhsn; }
125  void setRhsNorm(double r) { m_rhsn = r; }
126 
127  bool converged() const { return m_res <= m_rhsn * m_resmax; }
128  bool converged(double nr)
129  {
130  using std::abs;
131  m_res = abs(nr);
132  m_resminreach = (std::min)(m_resminreach, m_res);
133  return converged();
134  }
135  template<typename VectorType> bool converged(const VectorType &v)
136  { return converged(v.squaredNorm()); }
137 
138  bool finished(double nr)
139  {
140  if (m_callback) m_callback(*this);
141  if (m_noise > 0 && !m_written)
142  {
143  converged(nr);
144  m_written = true;
145  }
146  return (m_nit >= m_maxiter || converged(nr));
147  }
148  template <typename VectorType>
149  bool finished(const MatrixBase<VectorType> &v)
150  { return finished(double(v.squaredNorm())); }
151 
152 };
153 
154 } // end namespace Eigen
155 
156 #endif // EIGEN_ITERATION_CONTROLLER_H
Controls the iterations of the iterative solvers.
Definition: IterationController.h:74
size_t m_maxiter
Max. number of iterations.
Definition: IterationController.h:77
double m_rhsn
Right hand side norm.
Definition: IterationController.h:76
double m_resmax
maximum residual
Definition: IterationController.h:79
size_t m_nit
iteration number
Definition: IterationController.h:81
double m_res
last computed residual
Definition: IterationController.h:82
int m_noise
if noise > 0 iterations are printed
Definition: IterationController.h:78
Namespace containing all symbols from the Eigen library.
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)