Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooMCIntegrator.cxx
Go to the documentation of this file.
1/// \cond ROOFIT_INTERNAL
2
3/*****************************************************************************
4 * Project: RooFit *
5 * Package: RooFitCore *
6 * @(#)root/roofitcore:$Id$
7 * Authors: *
8 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
9 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
10 * *
11 * Copyright (c) 2000-2005, Regents of the University of California *
12 * and Stanford University. All rights reserved. *
13 * *
14 * Redistribution and use in source and binary forms, *
15 * with or without modification, are permitted according to the terms *
16 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
17 *****************************************************************************/
18
19/**
20\file RooMCIntegrator.cxx
21\class RooMCIntegrator
22\ingroup Roofitcore
23
24Implements an adaptive multi-dimensional Monte Carlo
25numerical integration, following the VEGAS algorithm originally described
26in G. P. Lepage, J. Comp. Phys. 27, 192(1978). This implementation is
27based on a C version from the 0.9 beta release of the GNU scientific library.
28**/
29
30#include "Riostream.h"
31
32#include "TMath.h"
33#include "RooMCIntegrator.h"
34#include "RooArgSet.h"
35#include "RooNumber.h"
36#include "RooNumIntFactory.h"
37#include "RooRealVar.h"
38#include "RooCategory.h"
39#include "RooMsgService.h"
40
41#include <cmath>
42
43
44using std::endl;
45
46
47// Register this class with RooNumIntFactory
48
49
50////////////////////////////////////////////////////////////////////////////////
51/// This function registers class RooMCIntegrator, its configuration options
52/// and its capabilities with RooNumIntFactory
53
54void RooMCIntegrator::registerIntegrator(RooNumIntFactory& fact)
55{
56 RooCategory samplingMode("samplingMode","Sampling Mode") ;
57 samplingMode.defineType("Importance",RooMCIntegrator::Importance) ;
58 samplingMode.defineType("ImportanceOnly",RooMCIntegrator::ImportanceOnly) ;
59 samplingMode.defineType("Stratified",RooMCIntegrator::Stratified) ;
60 samplingMode.setIndex(RooMCIntegrator::Importance) ;
61
62 RooCategory genType("genType","Generator Type") ;
63 genType.defineType("QuasiRandom",RooMCIntegrator::QuasiRandom) ;
64 genType.defineType("PseudoRandom",RooMCIntegrator::PseudoRandom) ;
65 genType.setIndex(RooMCIntegrator::QuasiRandom) ;
66
67 RooCategory verbose("verbose","Verbose flag") ;
68 verbose.defineType("true",1) ;
69 verbose.defineType("false",0) ;
70 verbose.setIndex(0) ;
71
72 RooRealVar alpha("alpha","Grid structure constant",1.5) ;
73 RooRealVar nRefineIter("nRefineIter","Number of refining iterations",5) ;
74 RooRealVar nRefinePerDim("nRefinePerDim","Number of refining samples (per dimension)",1000) ;
75 RooRealVar nIntPerDim("nIntPerDim","Number of integration samples (per dimension)",5000) ;
76
77 // Create prototype integrator
78 auto creator = [](const RooAbsFunc& function, const RooNumIntConfig& config) {
79 return std::make_unique<RooMCIntegrator>(function,config);
80 };
81
82 // Register prototype and default config with factory
83 std::string name = "RooMCIntegrator";
85 /*canIntegrate1D=*/true,
86 /*canIntegrate2D=*/true,
87 /*canIntegrateND=*/true,
88 /*canIntegrateOpenEnded=*/false);
89
90 // Make this method the default for all N>2-dim integrals
91 RooNumIntConfig::defaultConfig().methodND().setLabel(name) ;
92}
93
94
95////////////////////////////////////////////////////////////////////////////////
96/// Construct an integrator over 'function' with given sampling mode
97/// and generator type. The sampling mode can be 'Importance'
98/// (default), 'ImportanceOnly' and 'Stratified'. The generator type
99/// can be 'QuasiRandom' (default) and 'PseudoRandom'. Consult the original
100/// VEGAS documentation on details of the mode and type parameters.
101
102RooMCIntegrator::RooMCIntegrator(const RooAbsFunc& function, SamplingMode mode,
103 GeneratorType genType, bool verbose) :
104 RooAbsIntegrator(function), _grid(function), _verbose(verbose),
105 _alpha(1.5), _mode(mode), _genType(genType),
107{
108 // coverity[UNINIT_CTOR]
109 if(!(_valid= _grid.isValid())) return;
110 if(_verbose) _grid.print(std::cout);
111}
112
113
114
115////////////////////////////////////////////////////////////////////////////////
116/// Construct an integrator over 'function' where the configuration details
117/// are taken from 'config'
118
119RooMCIntegrator::RooMCIntegrator(const RooAbsFunc& function, const RooNumIntConfig& config) :
121{
122 const RooArgSet& configSet = config.getConfigSection("RooMCIntegrator") ;
123 _verbose = (bool) configSet.getCatIndex("verbose",0) ;
124 _alpha = configSet.getRealValue("alpha",1.5) ;
125 _mode = (SamplingMode) configSet.getCatIndex("samplingMode",Importance) ;
126 _genType = (GeneratorType) configSet.getCatIndex("genType",QuasiRandom) ;
127 _nRefineIter = (Int_t) configSet.getRealValue("nRefineIter",5) ;
128 _nRefinePerDim = (Int_t) configSet.getRealValue("nRefinePerDim",1000) ;
129 _nIntegratePerDim = (Int_t) configSet.getRealValue("nIntPerDim",5000) ;
130
131 // check that our grid initialized without errors
132 if(!(_valid= _grid.isValid())) return;
133 if(_verbose) _grid.print(std::cout);
134}
135
136////////////////////////////////////////////////////////////////////////////////
137/// Check if we can integrate over the current domain. If return value
138/// is true we cannot handle the current limits (e.g. where the domain
139/// of one or more observables is open ended.
140
141bool RooMCIntegrator::checkLimits() const
142{
143 return _grid.initialize(*integrand());
144}
145
146
147
148////////////////////////////////////////////////////////////////////////////////
149/// Evaluate the integral using a fixed number of calls to evaluate the integrand
150/// equal to about 10k per dimension. Use the first 5k calls to refine the grid
151/// over 5 iterations of 1k calls each, and the remaining 5k calls for a single
152/// high statistics integration.
153
154double RooMCIntegrator::integral(const double* /*yvec*/)
155{
156 _timer.Start(true);
157 vegas(AllStages,_nRefinePerDim*_grid.getDimension(),_nRefineIter);
158 double ret = vegas(ReuseGrid,_nIntegratePerDim*_grid.getDimension(),1);
159 return ret ;
160}
161
162
163
164////////////////////////////////////////////////////////////////////////////////
165/// Perform one step of Monte Carlo integration using the specified number of iterations
166/// with (approximately) the specified number of integrand evaluation calls per iteration.
167/// Use the VEGAS algorithm, starting from the specified stage. Returns the best estimate
168/// of the integral. Also sets *absError to the estimated absolute error of the integral
169/// estimate if absError is non-zero.
170
171double RooMCIntegrator::vegas(Stage stage, UInt_t calls, UInt_t iterations, double *absError)
172{
173 //cout << "VEGAS stage = " << stage << " calls = " << calls << " iterations = " << iterations << std::endl ;
174
175 // reset the grid to its initial state if we are starting from scratch
176 if(stage == AllStages) _grid.initialize(*_function);
177
178 // reset the results of previous calculations on this grid, but reuse the grid itself.
179 if(stage <= ReuseGrid) {
180 _wtd_int_sum = 0;
181 _sum_wgts = 0;
182 _chi_sum = 0;
183 _it_num = 1;
184 _samples = 0;
185 }
186
187 // refine the results of previous calculations on the current grid.
188 if(stage <= RefineGrid) {
189 UInt_t bins = RooGrid::maxBins;
190 UInt_t boxes = 1;
191 UInt_t dim(_grid.getDimension());
192
193 // select the sampling mode for the next step
194 if(_mode != ImportanceOnly) {
195 // calculate the largest number of equal subdivisions ("boxes") along each
196 // axis that results in an average of no more than 2 integrand calls per cell
197 boxes = (UInt_t)floor(std::pow(calls/2.0,1.0/dim));
198 // use stratified sampling if we are allowed enough calls (or equivalently,
199 // if the dimension is low enough)
200 _mode = Importance;
201 if (2*boxes >= RooGrid::maxBins) {
202 _mode = Stratified;
203 // adjust the number of bins and boxes to give an integral number >= 1 of boxes per bin
204 Int_t box_per_bin= (boxes > RooGrid::maxBins) ? boxes/RooGrid::maxBins : 1;
205 bins= boxes/box_per_bin;
206 if(bins > RooGrid::maxBins) bins= RooGrid::maxBins;
207 boxes = box_per_bin * bins;
208 oocxcoutD((TObject*)nullptr,Integration) << "RooMCIntegrator: using stratified sampling with " << bins << " bins and "
209 << box_per_bin << " boxes/bin" << std::endl;
210 }
211 else {
212 oocxcoutD((TObject*)nullptr,Integration) << "RooMCIntegrator: using importance sampling with " << bins << " bins and "
213 << boxes << " boxes" << std::endl;
214 }
215 }
216
217 // calculate the total number of n-dim boxes for this step
218 double tot_boxes = std::pow((double)boxes,(double)dim);
219
220 // increase the total number of calls to get at least 2 calls per box, if necessary
224
225 // calculate the Jacobean factor: volume/(avg # of calls/bin)
226 _jac = _grid.getVolume()*std::pow((double)bins,(double)dim)/calls;
227
228 // setup our grid to use the calculated number of boxes and bins
229 _grid.setNBoxes(boxes);
230 if(bins != _grid.getNBins()) _grid.resize(bins);
231 }
232
233 // allocate memory for some book-keeping arrays
234 std::vector<UInt_t> box(_grid.getDimension());
235 std::vector<UInt_t> bin(_grid.getDimension());
236 std::vector<double> x(_grid.getDimension());
237
238
239 // loop over iterations for this step
240 double cum_int(0);
241 double cum_sig(0);
243 _chisq = 0.0;
244 for (UInt_t it = 0; it < iterations; it++) {
245 double intgrl(0);
246 double intgrl_sq(0);
247 double sig(0);
248 double jacbin(_jac);
249
250 _it_num = _it_start + it;
251
252 // reset the values associated with each grid cell
253 _grid.resetValues();
254
255 // loop over grid boxes
256 _grid.firstBox(box.data());
257 do {
258 double m(0);
259 double q(0);
260 // loop over integrand evaluations within this grid box
261 for(UInt_t k = 0; k < _calls_per_box; k++) {
262 // generate a random point in this box
263 double bin_vol(0);
264 _grid.generatePoint(box.data(), x.data(), bin.data(), bin_vol, _genType == QuasiRandom ? true : false);
265 // evaluate the integrand at the generated point
266 double fval= jacbin*bin_vol*integrand(x.data());
267 // update mean and variance calculations
268 double d = fval - m;
269 m+= d / (k + 1.0);
270 q+= d * d * (k / (k + 1.0));
271 // accumulate the results of this evaluation (importance sampling only)
272 if (_mode != Stratified) _grid.accumulate(bin.data(), fval*fval);
273 }
275 double f_sq_sum = q * _calls_per_box ;
276 sig += f_sq_sum ;
277
278 // accumulate the results for this grid box (stratified sampling only)
279 if (_mode == Stratified) _grid.accumulate(bin.data(), f_sq_sum);
280
281 // print occasional progress messages
282 if(_timer.RealTime() > 30) {
283 std::size_t index = 0;
284 std::size_t sizeOfDim = 1;
285
286 for (unsigned int i=0; i < _grid.getDimension(); ++i) {
287 index += box[i] * sizeOfDim;
288 sizeOfDim *= _grid.getNBoxes();
289 }
290 oocoutP(nullptr, Integration) << "RooMCIntegrator: still working ... iteration "
291 << it << '/' << iterations << " box " << index << "/"<< std::pow(_grid.getNBoxes(), _grid.getDimension()) << std::endl;
292 _timer.Start(true);
293 }
294 else {
295 _timer.Start(false);
296 }
297
298 } while(_grid.nextBox(box.data()));
299
300 // compute final results for this iteration
301 double wgt;
302 sig = sig / (_calls_per_box - 1.0) ;
303 if (sig > 0) {
304 wgt = 1.0 / sig;
305 }
306 else if (_sum_wgts > 0) {
308 }
309 else {
310 wgt = 0.0;
311 }
313 _result = intgrl;
314 _sigma = sqrt(sig);
315
316 if (wgt > 0.0) {
317 _samples++ ;
318 _sum_wgts += wgt;
321
323 cum_sig = sqrt (1 / _sum_wgts);
324
325 if (_samples > 1) {
327 }
328 }
329 else {
330 cum_int += (intgrl - cum_int) / (it + 1.0);
331 cum_sig = 0.0;
332 }
333 oocxcoutD((TObject*)nullptr,Integration) << "=== Iteration " << _it_num << " : I = " << intgrl << " +/- " << sqrt(sig) << std::endl
334 << " Cumulative : I = " << cum_int << " +/- " << cum_sig << "( chi2 = " << _chisq
335 << ")" << std::endl;
336 // print the grid after the final iteration
337 if (oodologD((TObject*)nullptr,Integration)) {
338 if(it + 1 == iterations) _grid.print(std::cout, true);
339 }
340 _grid.refine(_alpha);
341 }
342
344 return cum_int;
345}
346
347/// \endcond
#define d(i)
Definition RSha256.hxx:102
#define oocxcoutD(o, a)
#define oodologD(o, a)
#define oocoutP(o, a)
int Int_t
Definition RtypesCore.h:45
unsigned int UInt_t
Definition RtypesCore.h:46
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 char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char mode
char name[80]
Definition TGX11.cxx:110
float * q
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
Definition RooAbsFunc.h:27
Abstract interface for integrators of real-valued functions that implement the RooAbsFunc interface.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Object to represent discrete states.
Definition RooCategory.h:28
Holds the configuration parameters of the various numeric integrators used by RooRealIntegral.
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
static RooNumIntConfig & defaultConfig()
Return reference to instance of default numeric integrator configuration object.
Factory to instantiate numeric integrators from a given function binding and a given configuration.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Mother of all ROOT objects.
Definition TObject.h:41
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition fillpatterns.C:1
Double_t x[n]
Definition legend1.C:17
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Definition RExports.h:167
TMarker m
Definition textangle.C:8