Mantid
Loading...
Searching...
No Matches
LevenbergMarquardtMinimizer.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//----------------------------------------------------------------------
12
15
17#include "MantidKernel/Logger.h"
18
19#include <gsl/gsl_blas.h>
20#include <gsl/gsl_version.h>
21
23namespace {
24// Get a reference to the logger
25Kernel::Logger g_log("LevenbergMarquardtMinimizer");
26
27bool cannotReachSpecifiedToleranceInF(int errorCode) { return errorCode == GSL_ETOLF; }
28bool cannotReachSpecifiedToleranceInX(int errorCode) { return errorCode == GSL_ETOLX; }
29} // namespace
30
31// clang-format off
32DECLARE_FUNCMINIMIZER(LevenbergMarquardtMinimizer, Levenberg-Marquardt)
33// clang-format on
34
36 : m_data(nullptr), gslContainer(), m_gslSolver(nullptr), m_function(), m_absError(1e-4), m_relError(1e-4) {
37 declareProperty("AbsError", m_absError,
38 "Absolute error allowed for "
39 "parameters - a stopping parameter "
40 "in success.");
41 declareProperty("RelError", m_relError,
42 "Relative error allowed for "
43 "parameters - a stopping parameter "
44 "in success.");
45}
46
47void LevenbergMarquardtMinimizer::initialize(API::ICostFunction_sptr costFunction, size_t /*maxIterations*/) {
48 // set-up GSL container to be used with GSL simplex algorithm
49 auto leastSquares = std::dynamic_pointer_cast<CostFunctions::CostFuncLeastSquares>(costFunction);
50 if (leastSquares) {
51 m_data = std::make_unique<GSL_FitData>(leastSquares);
52 } else {
53 throw std::runtime_error("LevenbergMarquardt can only be used with Least "
54 "squares cost function.");
55 }
56
57 // specify the type of GSL solver to use
58 const gsl_multifit_fdfsolver_type *T = gsl_multifit_fdfsolver_lmsder;
59
60 // setup GSL container
62 gslContainer.df = &gsl_df;
63 gslContainer.fdf = &gsl_fdf;
64
65 gslContainer.n = m_data->n;
66 gslContainer.p = m_data->p;
67 gslContainer.params = m_data.get();
68
69 // setup GSL solver
70 m_gslSolver = gsl_multifit_fdfsolver_alloc(T, m_data->n, m_data->p);
71 if (!m_gslSolver) {
72 throw std::runtime_error("Levenberg-Marquardt minimizer failed to initialize. \n" + std::to_string(m_data->n) +
73 " data points, " + std::to_string(m_data->p) + " fitting parameters. ");
74 }
75 gsl_multifit_fdfsolver_set(m_gslSolver, &gslContainer, m_data->initFuncParams);
76
77 m_function = leastSquares->getFittingFunction();
78}
79
81 if (m_gslSolver) {
82 gsl_multifit_fdfsolver_free(m_gslSolver);
83 }
84}
85
86bool LevenbergMarquardtMinimizer::iterate(size_t /*iteration*/) {
87 m_absError = getProperty("AbsError");
88 m_relError = getProperty("RelError");
89
90 int retVal = gsl_multifit_fdfsolver_iterate(m_gslSolver);
91
92 // From experience it is found that gsl_multifit_fdfsolver_iterate
93 // occasionally get
94 // stock - even after having achieved a sensible fit. This seem in particular
95 // to be a
96 // problem on Linux. However, to force GSL not to return ga ga have to do
97 // stuff in the
98 // if statement below
99 // GSL 1.14 changed return value from GSL_CONTINUE->GSL_ENOPROG for
100 // non-converging fits at 10 iterations
101 if (retVal == GSL_CONTINUE || retVal == GSL_ENOPROG) {
102 size_t ia = 0;
103 for (size_t i = 0; i < m_function->nParams(); i++) {
104 if (m_function->isActive(i)) {
105 m_function->setActiveParameter(i, gsl_vector_get(m_gslSolver->x, ia));
106 ++ia;
107 }
108 }
109 m_function->applyTies();
110 retVal = GSL_CONTINUE;
111 }
112
113 if (retVal && retVal != GSL_CONTINUE) {
114 m_errorString = gsl_strerror(retVal);
115 if (cannotReachSpecifiedToleranceInF(retVal)) {
116 m_errorString = "Changes in function value are too small";
117 } else if (cannotReachSpecifiedToleranceInX(retVal)) {
118 m_errorString = "Changes in parameter value are too small";
119 }
120 return false;
121 }
122
123 retVal = hasConverged();
124 return retVal != GSL_SUCCESS;
125}
126
128 return gsl_multifit_test_delta(m_gslSolver->dx, m_gslSolver->x, m_absError, m_relError);
129}
130
132 double chi = gsl_blas_dnrm2(m_gslSolver->f);
133 return chi * chi;
134}
135
136/* Calculates covariance matrix
137 *
138 * @param epsrel :: Is used to remove linear-dependent columns
139 * @param covar :: Returned covariance matrix, here as
140 */
141void LevenbergMarquardtMinimizer::calCovarianceMatrix(double epsrel, gsl_matrix *covar) {
142#if GSL_MAJOR_VERSION < 2
143 gsl_multifit_covar(m_gslSolver->J, epsrel, covar);
144#else
145 gsl_matrix *J = gsl_matrix_alloc(gslContainer.n, gslContainer.p);
146 gsl_multifit_fdfsolver_jac(m_gslSolver, J);
147 gsl_multifit_covar(J, epsrel, covar);
148 gsl_matrix_free(J);
149#endif
150}
151
152} // namespace Mantid::CurveFitting::FuncMinimisers
#define DECLARE_FUNCMINIMIZER(classname, username)
Macro for declaring a new type of minimizers to be used with the FuncMinimizerFactory.
const std::vector< Type > & m_data
std::string m_errorString
Error string.
Implementing Levenberg-Marquardt by wrapping the IFuncMinimizer interface around the GSL implementati...
gsl_multifit_fdfsolver * m_gslSolver
pointer to the GSL solver doing the work
double costFunctionVal() override
Return current value of the cost function.
API::IFunction_sptr m_function
Stored to access IFunction interface in iterate()
void initialize(API::ICostFunction_sptr costFunction, size_t maxIterations=0) override
Initialize minimizer, i.e. pass a function to minimize.
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
int gsl_fdf(const gsl_vector *x, void *params, gsl_vector *f, gsl_matrix *J)
Fit derivatives and function GSL wrapper.
int gsl_f(const gsl_vector *x, void *params, gsl_vector *f)
Fit GSL function wrapper.
int gsl_df(const gsl_vector *x, void *params, gsl_matrix *J)
Fit GSL derivative function wrapper.
std::string to_string(const wide_integer< Bits, Signed > &n)