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