Mantid
Loading...
Searching...
No Matches
LevenbergMarquardtMDMinimizer.cpp
Go to the documentation of this file.
1// Mantid Repository : https://github.com/mantidproject/mantid
2//
3// Copyright © 2018 ISIS Rutherford Appleton Laboratory UKRI,
4// NScD Oak Ridge National Laboratory, European Spallation Source,
5// Institut Laue - Langevin & CSNS, Institute of High Energy Physics, CAS
6// SPDX - License - Identifier: GPL - 3.0 +
7//----------------------------------------------------------------------
8// Includes
9//----------------------------------------------------------------------
13
16#include "MantidAPI/IFunction.h"
17
18#include "MantidKernel/Logger.h"
19
20#include <cmath>
21
23namespace {
25Kernel::Logger g_log("LevenbergMarquardMD");
26} // namespace
27
28// clang-format off
29DECLARE_FUNCMINIMIZER(LevenbergMarquardtMDMinimizer, Levenberg-MarquardtMD)
30// clang-format on
31
32
34 : IFuncMinimizer(), m_tau(1e-6), m_mu(1e-6), m_nu(2.0), m_rho(1.0), m_F(0.0) {
35 declareProperty("MuMax", 1e6, "Maximum value of mu - a stopping parameter in failure.");
36 declareProperty("AbsError", 0.0001,
37 "Absolute error allowed for parameters - "
38 "a stopping parameter in success.");
39 declareProperty("Verbose", false, "Make output more verbose.");
40}
41
44 m_costFunction = std::dynamic_pointer_cast<CostFunctions::CostFuncFitting>(function);
45 if (!m_costFunction) {
46 throw std::invalid_argument("Levenberg-Marquardt minimizer works only with "
47 "functions which define the Hessian. Different function was given.");
48 }
49 m_mu = 0;
50 m_nu = 2.0;
51 m_rho = 1.0;
52}
53
55bool LevenbergMarquardtMDMinimizer::iterate(size_t /*iteration*/) {
56 const bool verbose = getProperty("Verbose");
57 const double muMax = getProperty("MuMax");
58 const double absError = getProperty("AbsError");
59
60 if (!m_costFunction) {
61 throw std::runtime_error("Cost function isn't set up.");
62 }
63 size_t n = m_costFunction->nParams();
64
65 if (n == 0) {
66 m_errorString = "No parameters to fit.";
67 return false;
68 }
69
70 if (m_mu > muMax) {
71 m_errorString = "Failed to converge, maximum mu reached.";
72 return false;
73 }
74
75 // calculate the first and second derivatives of the cost function.
76 if (m_mu == 0.0 || m_rho > 0) {
77 // calculate everything first time or
78 // if last iteration was good
79 m_F = m_costFunction->valDerivHessian();
80 }
81 // else if m_rho < 0 last iteration was bad: reuse m_der and m_hessian
82
83 // Calculate damping to hessian
84 if (m_mu == 0) // first iteration or accidental zero
85 {
86 m_mu = m_tau;
87 m_nu = 2.0;
88 }
89
90 if (verbose) {
91 g_log.warning() << "===========================================================\n";
92 g_log.warning() << "mu=" << m_mu << "\n\n";
93 }
94
95 if (m_D.empty()) {
96 m_D.resize(n);
97 }
98
99 // copy the hessian
100 EigenMatrix H(m_costFunction->getHessian());
101 EigenVector dd(m_costFunction->getDeriv());
102
103 // scaling factors
104 std::vector<double> sf(n);
105
106 for (size_t i = 0; i < n; ++i) {
107 double d = fabs(dd.get(i));
108 if (m_D[i] > d)
109 d = m_D[i];
110 m_D[i] = d;
111 double tmp = H.get(i, i) + m_mu * d;
112 H.set(i, i, tmp);
113 sf[i] = sqrt(tmp);
114 if (tmp == 0.0) {
115 m_errorString = "Function doesn't depend on parameter " + m_costFunction->parameterName(i);
116 return false;
117 }
118 }
119
120 // apply scaling
121 for (size_t i = 0; i < n; ++i) {
122 double d = dd.get(i);
123 dd.set(i, d / sf[i]);
124 for (size_t j = i; j < n; ++j) {
125 const double f = sf[i] * sf[j];
126 double tmp = H.get(i, j);
127 H.set(i, j, tmp / f);
128 if (i != j) {
129 tmp = H.get(j, i);
130 H.set(j, i, tmp / f);
131 }
132 }
133 }
134
135 if (verbose && m_rho > 0) {
136 g_log.warning() << "Hessian:\n" << H;
137 g_log.warning() << "Right-hand side:\n";
138 for (size_t j = 0; j < n; ++j) {
139 g_log.warning() << dd.get(j) << ' ';
140 }
141 g_log.warning() << '\n';
142 g_log.warning() << "Determinant=" << H.det() << '\n';
143 }
144
145 // Parameter corrections
146 EigenVector dx(n);
147 // To find dx solve the system of linear equations H * dx == -m_der
148 dd *= -1.0;
149 try {
150 H.solve(dd, dx);
151 } catch (std::runtime_error &error) {
152 m_errorString = error.what();
153 return false;
154 }
155
156 if (verbose) {
157 g_log.warning() << "\nScaling factors:\n";
158 for (size_t j = 0; j < n; ++j) {
159 g_log.warning() << sf[j] << ' ';
160 }
161 g_log.warning() << '\n';
162 g_log.warning() << "Corrections:\n";
163 for (size_t j = 0; j < n; ++j) {
164 g_log.warning() << dx.get(j) << ' ';
165 }
166 g_log.warning() << "\n\n";
167 }
168
169 // restore scaling
170 for (size_t i = 0; i < n; ++i) {
171 double d = dx.get(i);
172 dx.set(i, d / sf[i]);
173 d = dd.get(i);
174 dd.set(i, d * sf[i]);
175 }
176
177 // save previous state
178 m_costFunction->push();
179 // Update the parameters of the cost function.
180 EigenVector parameters(n);
181 m_costFunction->getParameters(parameters);
182 parameters += dx;
183 m_costFunction->setParameters(parameters);
184 if (verbose) {
185 for (size_t i = 0; i < n; ++i) {
186 g_log.warning() << "Parameter(" << i << ")=" << parameters[i] << '\n';
187 }
188 }
189 m_costFunction->getFittingFunction()->applyTies();
190
191 // --- prepare for the next iteration --- //
192
193 double dL;
194 // der -> - der - 0.5 * hessian * dx
195
196 EigenVector Trdx = m_costFunction->getHessian().tr() * dx;
197 Trdx *= -0.5;
198 dd += Trdx;
199 // calculate the linear part of the change in cost function
200 dL = dd.dot(dx);
201
202 double F1 = m_costFunction->val();
203 if (verbose) {
204 g_log.warning() << '\n';
205 g_log.warning() << "Old cost function " << m_F << '\n';
206 g_log.warning() << "New cost function " << F1 << '\n';
207 g_log.warning() << "Linear part " << dL << '\n';
208 }
209
210 // Try the stop condition
211 if (m_rho >= 0) {
212 EigenVector p(n);
213 m_costFunction->getParameters(p);
214 double dx_norm = dx.norm();
215 if (dx_norm < absError) {
216 if (verbose) {
217 g_log.warning() << "Successful fit, parameters changed by less than " << absError << '\n';
218 }
219 return false;
220 }
221 if (m_rho == 0) {
222 if (m_F != F1) {
223 this->m_errorString = "Failed to converge, rho == 0";
224 }
225 if (verbose) {
226 g_log.warning() << "Successful fit, cost function didn't change.\n";
227 }
228 return false;
229 }
230 }
231
232 if (fabs(dL) == 0.0) {
233 if (m_F == F1)
234 m_rho = 1.0;
235 else
236 m_rho = 0;
237 } else {
238 m_rho = (m_F - F1) / dL;
239 if (m_rho == 0) {
240 return false;
241 }
242 }
243 if (verbose) {
244 g_log.warning() << "rho=" << m_rho << '\n';
245 }
246
247 if (m_rho > 0) { // good progress, decrease m_mu but no more than by 1/3
248 // rho = 1 - (2*rho - 1)^3
249 m_rho = 2.0 * m_rho - 1.0;
250 m_rho = 1.0 - m_rho * m_rho * m_rho;
251 const double I3 = 1.0 / 3.0;
252 if (m_rho > I3)
253 m_rho = I3;
254 if (m_rho < 0.0001)
255 m_rho = 0.1;
256 m_mu *= m_rho;
257 m_nu = 2.0;
258 m_F = F1;
259 if (verbose) {
260 g_log.warning() << "Good iteration, accept new parameters.\n";
261 g_log.warning() << "rho=" << m_rho << '\n';
262 }
263 // drop saved state, accept new parameters
264 m_costFunction->drop();
265 } else { // bad iteration. increase m_mu and revert changes to parameters
266 m_mu *= m_nu;
267 m_nu *= 2.0;
268 // undo parameter update
269 m_costFunction->pop();
270 m_F = m_costFunction->val();
271 if (verbose) {
272 g_log.warning() << "Bad iteration, increase mu and revert changes to parameters.\n";
273 }
274 }
275
276 return true;
277}
278
281 if (!m_costFunction) {
282 throw std::runtime_error("Cost function isn't set up.");
283 }
284 return m_costFunction->val();
285}
286
287} // namespace Mantid::CurveFitting::FuncMinimisers
gsl_vector * tmp
#define DECLARE_FUNCMINIMIZER(classname, username)
Macro for declaring a new type of minimizers to be used with the FuncMinimizerFactory.
double error
#define fabs(x)
Definition Matrix.cpp:22
An interface for function minimizers.
std::string m_errorString
Error string.
A wrapper around Eigen::Matrix.
Definition EigenMatrix.h:33
A wrapper around Eigen::Vector.
Definition EigenVector.h:27
void set(const size_t i, const double value)
Set an element.
double dot(const EigenVector &v) const
Calculate the dot product.
double get(const size_t i) const
Get an element.
double norm() const
Get vector norm (length)
double m_mu
The damping mu parameter in the Levenberg-Marquardt method.
std::shared_ptr< CostFunctions::CostFuncFitting > m_costFunction
Pointer to the cost function.
double m_rho
The rho parameter in the Levenberg-Marquardt method.
void initialize(API::ICostFunction_sptr function, size_t maxIterations=0) override
Initialize minimizer, i.e. pass a function to minimize.
double m_nu
The nu parameter in the Levenberg-Marquardt method.
double costFunctionVal() override
Return current value of the cost function.
double m_tau
The tau parameter in the Levenberg-Marquardt method.
void warning(const std::string &msg)
Logs at warning level.
Definition Logger.cpp:117
TypedValue getProperty(const std::string &name) const override
Get the value of a property.
Kernel::Logger g_log("ExperimentInfo")
static logger object
std::shared_ptr< ICostFunction > ICostFunction_sptr
define a shared pointer to a cost function