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