Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooGradMinimizerFcn.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * AL, Alfio Lazzaro, INFN Milan, alfio.lazzaro@mi.infn.it *
7 * PB, Patrick Bos, NL eScience Center, p.bos@esciencecenter.nl *
8 * VC, Vince Croft, DIANA / NYU, vincent.croft@cern.ch *
9 * *
10 * *
11 * Redistribution and use in source and binary forms, *
12 * with or without modification, are permitted according to the terms *
13 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
14 *****************************************************************************/
15
16//////////////////////////////////////////////////////////////////////////////
17//
18// RooGradMinimizerFcn is am interface class to the ROOT::Math function
19// for minimization. See RooGradMinimizer.cxx for more information.
20//
21
22#include "RooGradMinimizerFcn.h"
23
24#include "RooAbsArg.h"
25#include "RooAbsPdf.h"
26#include "RooArgSet.h"
27#include "RooRealVar.h"
28#include "RooMsgService.h"
29#include "RooMinimizer.h"
30
31#include "Riostream.h"
32#include "TIterator.h"
33#include "TClass.h"
34#include "Fit/Fitter.h"
35#include "Math/Minimizer.h"
36
37#include <algorithm> // std::equal
38#include <iostream>
39
41 : RooAbsMinimizerFcn(RooArgList( * std::unique_ptr<RooArgSet>(funct->getParameters(RooArgSet{})) ), context, verbose),
42 _grad(getNDim()), _grad_params(getNDim()), _funct(funct),
43 has_been_calculated(getNDim())
44{
45 // TODO: added "parameters" after rewrite in april 2020, check if correct
46 auto parameters = _context->fitter()->Config().ParamsSettings();
47 synchronizeParameterSettings(parameters, kTRUE, verbose);
51}
52
54 : RooAbsMinimizerFcn(other), _grad(other._grad), _grad_params(other._grad_params), _gradf(other._gradf), _funct(other._funct),
55 has_been_calculated(other.has_been_calculated), none_have_been_calculated(other.none_have_been_calculated)
56{
57}
58
60{
61 return new RooGradMinimizerFcn(*this);
62}
63
65 std::vector<ROOT::Fit::ParameterSettings> &parameter_settings) const
66{
67 _gradf.SetInitialGradient(nullptr, parameter_settings, _grad);
68}
69
70////////////////////////////////////////////////////////////////////////////////
71
72double RooGradMinimizerFcn::DoEval(const double *x) const
73{
74 Bool_t parameters_changed = kFALSE;
75
76 // Set the parameter values for this iteration
77 for (unsigned index = 0; index < NDim(); index++) {
78 // also check whether the function was already evaluated for this set of parameters
79 parameters_changed |= SetPdfParamVal(index, x[index]);
80 }
81
82 // Calculate the function for these parameters
84 double fvalue = _funct->getVal();
86
87 if (!parameters_changed) {
88 return fvalue;
89 }
90
91 if (!std::isfinite(fvalue) || RooAbsReal::numEvalErrors() > 0 || fvalue > 1e30) {
92
93 if (_printEvalErrors >= 0) {
94
95 if (_doEvalErrorWall) {
96 oocoutW(static_cast<RooAbsArg *>(nullptr), Eval)
97 << "RooGradMinimizerFcn: Minimized function has error status." << std::endl
98 << "Returning maximum FCN so far (" << _maxFCN
99 << ") to force MIGRAD to back out of this region. Error log follows" << std::endl;
100 } else {
101 oocoutW(static_cast<RooAbsArg *>(nullptr), Eval)
102 << "RooGradMinimizerFcn: Minimized function has error status but is ignored" << std::endl;
103 }
104
105 TIterator *iter = _floatParamList->createIterator();
106 RooRealVar *var;
108 ooccoutW(static_cast<RooAbsArg *>(nullptr), Eval) << "Parameter values: ";
109 while ((var = (RooRealVar *)iter->Next())) {
110 if (first) {
111 first = kFALSE;
112 } else
113 ooccoutW(static_cast<RooAbsArg *>(nullptr), Eval) << ", ";
114 ooccoutW(static_cast<RooAbsArg *>(nullptr), Eval) << var->GetName() << "=" << var->getVal();
115 }
116 delete iter;
117 ooccoutW(static_cast<RooAbsArg *>(nullptr), Eval) << std::endl;
118
119 RooAbsReal::printEvalErrors(ooccoutW(static_cast<RooAbsArg *>(nullptr), Eval), _printEvalErrors);
120 ooccoutW(static_cast<RooAbsArg *>(nullptr), Eval) << std::endl;
121 }
122
123 if (_doEvalErrorWall) {
124 fvalue = _maxFCN + 1;
125 }
126
128 _numBadNLL++;
129 } else if (fvalue > _maxFCN) {
130 _maxFCN = fvalue;
131 }
132
133 // Optional logging
134 if (_verbose) {
135 std::cout << "\nprevFCN" << (_funct->isOffsetting() ? "-offset" : "") << " = " << std::setprecision(10) << fvalue
136 << std::setprecision(4) << " ";
137 std::cout.flush();
138 }
139
140 _evalCounter++;
141 //#ifndef NDEBUG
142 // std::cout << "RooGradMinimizerFcn " << this << " evaluations (in DoEval): " << _evalCounter <<
143 // std::endl;
144 //#endif
145 return fvalue;
146}
147
149{
150 for (auto it = has_been_calculated.begin(); it != has_been_calculated.end(); ++it) {
151 *it = false;
152 }
154}
155
156bool RooGradMinimizerFcn::syncParameter(double x, std::size_t ix) const
157{
158 bool parameter_has_changed = (_grad_params[ix] != x);
159
160 if (parameter_has_changed) {
161 _grad_params[ix] = x;
162 // Set the parameter values for this iteration
163 // TODO: this is already done in DoEval as well; find efficient way to do only once
164 SetPdfParamVal(ix, x);
165
168 }
169 }
170
171 return parameter_has_changed;
172}
173
174bool RooGradMinimizerFcn::syncParameters(const double *x) const
175{
176 bool has_been_synced = false;
177
178 for (std::size_t ix = 0; ix < NDim(); ++ix) {
179 bool parameter_has_changed = (_grad_params[ix] != x[ix]);
180
181 if (parameter_has_changed) {
182 _grad_params[ix] = x[ix];
183 // Set the parameter values for this iteration
184 // TODO: this is already done in DoEval as well; find efficient way to do only once
185 SetPdfParamVal(ix, x[ix]);
186 }
187
188 has_been_synced |= parameter_has_changed;
189 }
190
191 if (has_been_synced) {
193 }
194
195 return has_been_synced;
196}
197
198void RooGradMinimizerFcn::runDerivator(unsigned int i_component) const
199{
200 // check whether the derivative was already calculated for this set of parameters
201 if (!has_been_calculated[i_component]) {
202 // Calculate the derivative etc for these parameters
203 _grad[i_component] =
205 _context->fitter()->Config().ParamsSettings(), i_component, _grad[i_component]);
206 has_been_calculated[i_component] = true;
208 }
209}
210
211double RooGradMinimizerFcn::DoDerivative(const double *x, unsigned int i_component) const
212{
214 runDerivator(i_component);
215 return _grad[i_component].derivative;
216}
217
218double RooGradMinimizerFcn::DoDerivativeWithPrevResult(const double *x, unsigned int i_component,
219 double *previous_grad, double *previous_g2,
220 double *previous_gstep) const
221{
223 _grad[i_component] = {previous_grad[i_component], previous_g2[i_component], previous_gstep[i_component]};
224 runDerivator(i_component);
225 previous_grad[i_component] = _grad[i_component].derivative;
226 previous_g2[i_component] = _grad[i_component].second_derivative;
227 previous_gstep[i_component] = _grad[i_component].step_size;
228 return _grad[i_component].derivative;
229}
230
231////////////////////////////////////////////////////////////////////////////////
232
234{
235 assert(istrat >= 0);
236 ROOT::Minuit2::MnStrategy strategy(static_cast<unsigned int>(istrat));
237
240 setNcycles(strategy.GradientNCycles());
241}
242
243Bool_t
244RooGradMinimizerFcn::Synchronize(std::vector<ROOT::Fit::ParameterSettings> &parameters, Bool_t optConst, Bool_t verbose)
245{
246 Bool_t returnee = synchronizeParameterSettings(parameters, optConst, verbose);
250 return returnee;
251}
#define oocoutW(o, a)
#define ooccoutW(o, a)
const Bool_t kFALSE
Definition RtypesCore.h:101
bool Bool_t
Definition RtypesCore.h:63
const Bool_t kTRUE
Definition RtypesCore.h:100
const std::vector< ROOT::Fit::ParameterSettings > & ParamsSettings() const
get the vector of parameter settings (const method)
Definition FitConfig.h:86
ROOT::Math::MinimizerOptions & MinimizerOptions()
access to the minimizer control parameter (non const method)
Definition FitConfig.h:167
const FitConfig & Config() const
access to the fit configuration (const method)
Definition Fitter.h:412
Interface (abstract class) for multi-dimensional functions providing a gradient calculation.
Definition IFunction.h:343
double ErrorDef() const
error definition
API class for defining three levels of strategies: low (0), medium (1), high (>=2); acts on: Migrad (...
Definition MnStrategy.h:27
double GradientStepTolerance() const
Definition MnStrategy.h:41
double GradientTolerance() const
Definition MnStrategy.h:42
unsigned int GradientNCycles() const
Definition MnStrategy.h:40
DerivatorElement PartialDerivative(const ROOT::Math::IBaseFunctionMultiDim *function, const double *x, const std::vector< ROOT::Fit::ParameterSettings > &parameters, unsigned int i_component, DerivatorElement previous)
void SetInitialGradient(const ROOT::Math::IBaseFunctionMultiDim *function, const std::vector< ROOT::Fit::ParameterSettings > &parameters, std::vector< DerivatorElement > &gradient)
This function was not implemented as in Minuit2.
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:69
Bool_t SetPdfParamVal(int index, double value) const
Set value of parameter i.
std::unique_ptr< RooArgList > _floatParamList
Bool_t synchronizeParameterSettings(std::vector< ROOT::Fit::ParameterSettings > &parameters, Bool_t optConst, Bool_t verbose)
Informs Minuit through its parameter_settings vector of RooFit parameter properties.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
virtual Bool_t isOffsetting() const
Definition RooAbsReal.h:372
static void setHideOffset(Bool_t flag)
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
static Int_t numEvalErrors()
Return the number of logged evaluation errors since the last clearing.
static void printEvalErrors(std::ostream &os=std::cout, Int_t maxPerNode=10000000)
Print all outstanding logged evaluation error on the given ostream.
static void clearEvalErrorLog()
Clear the stack of evaluation error messages.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:35
void synchronizeGradientParameterSettings(std::vector< ROOT::Fit::ParameterSettings > &parameter_settings) const
void setErrorLevel(double error_level) const
double DoDerivativeWithPrevResult(const double *x, unsigned int i_component, double *previous_grad, double *previous_g2, double *previous_gstep) const override
Bool_t Synchronize(std::vector< ROOT::Fit::ParameterSettings > &parameter_settings, Bool_t optConst, Bool_t verbose=kFALSE) override
Like synchronizeParameterSettings, Synchronize informs Minuit through its parameter_settings vector o...
unsigned int NDim() const override
Retrieve the dimension of the function.
void resetHasBeenCalculatedFlags() const
void setNcycles(unsigned int ncycles) const
double DoEval(const double *x) const override
void setGradTolerance(double grad_tolerance) const
std::vector< bool > has_been_calculated
RooGradMinimizerFcn(RooAbsReal *funct, RooMinimizer *context, bool verbose=false)
void runDerivator(unsigned int i_component) const
std::vector< double > _grad_params
void setStepTolerance(double step_tolerance) const
ROOT::Math::IMultiGradFunction * Clone() const override
Clone a function.
std::vector< ROOT::Minuit2::DerivatorElement > _grad
bool syncParameter(double x, std::size_t ix) const
ROOT::Minuit2::NumericalDerivator _gradf
bool syncParameters(const double *x) const
double DoDerivative(const double *x, unsigned int icoord) const override
RooMinimizer is a wrapper class around ROOT::Fit:Fitter that provides a seamless interface between th...
ROOT::Math::IMultiGenFunction * getMultiGenFcn() const
ROOT::Fit::Fitter * fitter()
Return underlying ROOT fitter object.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
Iterator abstract base class.
Definition TIterator.h:30
virtual TObject * Next()=0
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Double_t x[n]
Definition legend1.C:17
Definition first.py:1