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