Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
NumericalDerivator.cxx
Go to the documentation of this file.
1// @(#)root/mathcore:$Id$
2// Authors: L. Moneta, J.T. Offermann, E.G.P. Bos 2013-2018
3//
4/**********************************************************************
5 * *
6 * Copyright (c) 2013 , LCG ROOT MathLib Team *
7 * *
8 **********************************************************************/
9/**
10 * \class NumericalDerivator
11 *
12 * Original version created on: Aug 14, 2013
13 * Authors: L. Moneta, J. T. Offermann
14 * Modified version created on: Sep 27, 2017
15 * Author: E. G. P. Bos
16 *
17 * NumericalDerivator was essentially a slightly modified copy of code
18 * written by M. Winkler, F. James, L. Moneta, and A. Zsenei for Minuit2,
19 * Copyright (c) 2005 LCG ROOT Math team, CERN/PH-SFT. Original version:
20 * https://github.com/lmoneta/root/blob/lvmini/math/mathcore/src/NumericalDerivator.cxx
21 *
22 * This class attempts to more closely follow the Minuit2 implementation.
23 * Modified things (w.r.t. original) are indicated by MODIFIED.
24 */
25
27#include <cmath>
28#include <algorithm>
29#include <iostream>
30#include <TMath.h>
31#include <cassert>
33
34#include <Math/Minimizer.h> // needed here because in Fitter is only a forward declaration
35
36namespace ROOT {
37namespace Minuit2 {
38
43
46 : fStepTolerance(step_tolerance),
47 fGradTolerance(grad_tolerance),
48 fUp(error_level),
49 fNCycles(ncycles),
50 fAlwaysExactlyMimicMinuit2(always_exactly_mimic_minuit2)
51{
52}
53
55
56/// This function sets internal state based on input parameters. This state
57/// setup is used in the actual (partial) derivative calculations.
58void NumericalDerivator::SetupDifferentiate(unsigned int nDim, const FCNBase *function, const double *cx,
59 std::span<const ROOT::Fit::ParameterSettings> parameters)
60{
61 assert(function != nullptr && "function is a nullptr");
62
63 fVx.resize(nDim);
64 fVxExternal.resize(nDim);
65 fVxFValCache.resize(nDim);
66 std::copy(cx, cx + nDim, fVx.data());
67
68 // convert to Minuit external parameters
69 for (unsigned i = 0; i < nDim; i++) {
70 fVxExternal[i] = Int2ext(parameters[i], fVx[i]);
71 }
72
73 if (fVx != fVxFValCache) {
75 fVal = (*function)(fVxExternal); // value of function at given points
76 }
77
78 fDfmin = 8. * fPrecision.Eps2() * (std::abs(fVal) + fUp);
79 fVrysml = 8. * fPrecision.Eps() * fPrecision.Eps();
80}
81
82DerivatorElement NumericalDerivator::PartialDerivative(unsigned int nDim, const FCNBase *function, const double *x,
83 std::span<const ROOT::Fit::ParameterSettings> parameters,
85{
86 SetupDifferentiate(nDim, function, x, parameters);
87 return FastPartialDerivative(function, parameters, i_component, previous);
88}
89
90// leaves the parameter setup to the caller
92 std::span<const ROOT::Fit::ParameterSettings> parameters,
93 unsigned int i_component, const DerivatorElement &previous)
94{
96
97 double xtf = fVx[i_component];
98 double epspri = fPrecision.Eps2() + std::abs(deriv.derivative * fPrecision.Eps2());
99 double step_old = 0.;
100
101 for (unsigned int j = 0; j < fNCycles; ++j) {
102 double optstp = std::sqrt(fDfmin / (std::abs(deriv.second_derivative) + epspri));
103 double step = std::max(optstp, std::abs(0.1 * deriv.step_size));
104
105 if (parameters[i_component].IsBound()) {
106 if (step > 0.5)
107 step = 0.5;
108 }
109
110 double stpmax = 10. * std::abs(deriv.step_size);
111 if (step > stpmax)
112 step = stpmax;
113
114 double stpmin = std::max(fVrysml, 8. * std::abs(fPrecision.Eps2() * fVx[i_component]));
115 if (step < stpmin)
116 step = stpmin;
117 if (std::abs((step - step_old) / step) < fStepTolerance) {
118 break;
119 }
120 deriv.step_size = step;
121 step_old = step;
122 fVx[i_component] = xtf + step;
124 double fs1 = (*function)(fVxExternal);
125
126 fVx[i_component] = xtf - step;
128 double fs2 = (*function)(fVxExternal);
129
132
133 double fGrd_old = deriv.derivative;
134 deriv.derivative = 0.5 * (fs1 - fs2) / step;
135
136 deriv.second_derivative = (fs1 + fs2 - 2. * fVal) / step / step;
137
138 if (std::abs(fGrd_old - deriv.derivative) / (std::abs(deriv.derivative) + fDfmin / step) < fGradTolerance) {
139 break;
140 }
141 }
142 return deriv;
143}
144
145DerivatorElement NumericalDerivator::operator()(unsigned int nDim, const FCNBase *function, const double *x,
146 std::span<const ROOT::Fit::ParameterSettings> parameters,
147 unsigned int i_component, const DerivatorElement &previous)
148{
149 return PartialDerivative(nDim, function, x, parameters, i_component, previous);
150}
151
152std::vector<DerivatorElement>
153NumericalDerivator::Differentiate(unsigned int nDim, const FCNBase *function, const double *cx,
154 std::span<const ROOT::Fit::ParameterSettings> parameters,
155 std::span<const DerivatorElement> previous_gradient)
156{
157 SetupDifferentiate(nDim, function, cx, parameters);
158
159 std::vector<DerivatorElement> gradient;
160 gradient.reserve(nDim);
161
162 for (unsigned int ix = 0; ix < nDim; ++ix) {
163 gradient.emplace_back(FastPartialDerivative(function, parameters, ix, previous_gradient[ix]));
164 }
165
166 return gradient;
167}
168
170{
171 // return external value from internal value for parameter i
172 if (parameter.IsBound()) {
173 if (parameter.IsDoubleBound()) {
174 return fDoubleLimTrafo.Int2ext(val, parameter.UpperLimit(), parameter.LowerLimit());
175 } else if (parameter.HasUpperLimit() && !parameter.HasLowerLimit()) {
176 return fUpperLimTrafo.Int2ext(val, parameter.UpperLimit());
177 } else {
178 return fLowerLimTrafo.Int2ext(val, parameter.LowerLimit());
179 }
180 }
181
182 return val;
183}
184
186{
187 // return the internal value for parameter i with external value val
188
189 if (parameter.IsBound()) {
190 if (parameter.IsDoubleBound()) {
191 return fDoubleLimTrafo.Ext2int(val, parameter.UpperLimit(), parameter.LowerLimit(), fPrecision);
192 } else if (parameter.HasUpperLimit() && !parameter.HasLowerLimit()) {
193 return fUpperLimTrafo.Ext2int(val, parameter.UpperLimit(), fPrecision);
194 } else {
195 return fLowerLimTrafo.Ext2int(val, parameter.LowerLimit(), fPrecision);
196 }
197 }
198
199 return val;
200}
201
203{
204 // return the derivative of the int->ext transformation: dPext(i) / dPint(i)
205 // for the parameter i with value val
206
207 double dd = 1.;
208 if (parameter.IsBound()) {
209 if (parameter.IsDoubleBound()) {
210 dd = fDoubleLimTrafo.DInt2Ext(val, parameter.UpperLimit(), parameter.LowerLimit());
211 } else if (parameter.HasUpperLimit() && !parameter.HasLowerLimit()) {
212 dd = fUpperLimTrafo.DInt2Ext(val, parameter.UpperLimit());
213 } else {
214 dd = fLowerLimTrafo.DInt2Ext(val, parameter.LowerLimit());
215 }
216 }
217
218 return dd;
219}
220
221// MODIFIED:
222/// This function was not implemented as in Minuit2. Now it copies the behavior
223/// of InitialGradientCalculator. See https://github.com/roofit-dev/root/issues/10
224void NumericalDerivator::SetInitialGradient(std::span<const ROOT::Fit::ParameterSettings> parameters,
225 std::vector<DerivatorElement> &gradient)
226{
227 // set an initial gradient using some given steps
228 // (used in the first iteration)
229
230 double eps2 = fPrecision.Eps2();
231
232 unsigned ix = 0;
233 for (auto parameter = parameters.begin(); parameter != parameters.end(); ++parameter, ++ix) {
234 // What Minuit2 calls "Error" is stepsize on the ROOT side.
235 double werr = parameter->StepSize();
236
237 // Actually, extParVal in Minuit2 is the external parameter value, so that is
238 // what we called var before and var is unnecessary here.
239 double extParVal = parameter->Value();
240
241 // However, we do need var below, so let's calculate it using Ext2int:
242 double var = Ext2int(*parameter, extParVal);
243
245 // this transformation can lose a few bits, but Minuit2 does it too
246 extParVal = Int2ext(*parameter, var);
247 }
248
249 double extParVal2 = extParVal + werr;
250
251 if (parameter->HasUpperLimit() && extParVal2 > parameter->UpperLimit()) {
252 extParVal2 = parameter->UpperLimit();
253 }
254
255 double var2 = Ext2int(*parameter, extParVal2);
256 double vplu = var2 - var;
257
259
260 if (parameter->HasLowerLimit() && extParVal2 < parameter->LowerLimit()) {
261 extParVal2 = parameter->LowerLimit();
262 }
263
265 double vmin = var2 - var;
266 double gsmin = 8. * eps2 * (std::abs(var) + eps2);
267 // protect against very small step sizes which can cause dirin to zero and then nan values in grd
268 double dirin = std::max(0.5 * (std::abs(vplu) + std::abs(vmin)), gsmin);
269 double g2 = 2.0 * fUp / (dirin * dirin);
270 double gstep = std::max(gsmin, 0.1 * dirin);
271 double grd = g2 * dirin;
272
273 if (parameter->IsBound()) {
274 gstep = std::min(gstep, 0.5);
275 }
276
277 gradient[ix].derivative = grd;
278 gradient[ix].second_derivative = g2;
279 gradient[ix].step_size = gstep;
280 }
281}
282
283std::ostream &operator<<(std::ostream &out, const DerivatorElement &value)
284{
285 return out << "(derivative: " << value.derivative << ", second_derivative: " << value.second_derivative
286 << ", step_size: " << value.step_size << ")";
287}
288
289} // namespace Minuit2
290} // namespace ROOT
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Class, describing value, limits and step size of the parameters Provides functionality also to set/re...
Interface (abstract class) defining the function to be minimized, which has to be implemented by the ...
Definition FCNBase.h:37
double Eps() const
eps returns the smallest possible number so that 1.+eps > 1.
double Eps2() const
eps2 returns 2*sqrt(eps)
void SetInitialGradient(std::span< const ROOT::Fit::ParameterSettings > parameters, std::vector< DerivatorElement > &gradient)
This function was not implemented as in Minuit2.
double Int2ext(const ROOT::Fit::ParameterSettings &parameter, double val) const
ROOT::Minuit2::SqrtUpParameterTransformation fUpperLimTrafo
ROOT::Minuit2::MnMachinePrecision fPrecision
NumericalDerivator(bool always_exactly_mimic_minuit2=true)
void SetupDifferentiate(unsigned int nDim, const FCNBase *function, const double *cx, std::span< const ROOT::Fit::ParameterSettings > parameters)
This function sets internal state based on input parameters.
double DInt2Ext(const ROOT::Fit::ParameterSettings &parameter, double val) const
double Ext2int(const ROOT::Fit::ParameterSettings &parameter, double val) const
DerivatorElement operator()(unsigned int nDim, const FCNBase *function, const double *x, std::span< const ROOT::Fit::ParameterSettings > parameters, unsigned int i_component, const DerivatorElement &previous)
std::vector< DerivatorElement > Differentiate(unsigned int nDim, const FCNBase *function, const double *x, std::span< const ROOT::Fit::ParameterSettings > parameters, std::span< const DerivatorElement > previous_gradient)
ROOT::Minuit2::SinParameterTransformation fDoubleLimTrafo
DerivatorElement FastPartialDerivative(const FCNBase *function, std::span< const ROOT::Fit::ParameterSettings > parameters, unsigned int i_component, const DerivatorElement &previous)
ROOT::Minuit2::SqrtLowParameterTransformation fLowerLimTrafo
DerivatorElement PartialDerivative(unsigned int nDim, const FCNBase *function, const double *x, std::span< const ROOT::Fit::ParameterSettings > parameters, unsigned int i_component, DerivatorElement previous)
long double DInt2Ext(long double Value, long double Upper, long double Lower) const
long double Int2ext(long double Value, long double Upper, long double Lower) const
long double Ext2int(long double Value, long double Upper, long double Lower, const MnMachinePrecision &) const
long double Ext2int(long double Value, long double Lower, const MnMachinePrecision &) const
long double DInt2Ext(long double Value, long double Lower) const
long double Int2ext(long double Value, long double Lower) const
long double Int2ext(long double Value, long double Upper) const
long double DInt2Ext(long double Value, long double Upper) const
long double Ext2int(long double Value, long double Upper, const MnMachinePrecision &) const
Double_t x[n]
Definition legend1.C:17
std::ostream & operator<<(std::ostream &, const LAVector &)
Definition MnMatrix.cxx:691