48 ClassImp(RooMCIntegrator);
58 void RooMCIntegrator::registerIntegrator(RooNumIntFactory& fact)
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) ;
66 RooCategory genType(
"genType",
"Generator Type") ;
67 genType.defineType(
"QuasiRandom",RooMCIntegrator::QuasiRandom) ;
68 genType.defineType(
"PseudoRandom",RooMCIntegrator::PseudoRandom) ;
69 genType.setIndex(RooMCIntegrator::QuasiRandom) ;
71 RooCategory verbose(
"verbose",
"Verbose flag") ;
72 verbose.defineType(
"true",1) ;
73 verbose.defineType(
"false",0) ;
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) ;
82 RooMCIntegrator* proto =
new RooMCIntegrator() ;
85 fact.storeProtoIntegrator(proto,RooArgSet(samplingMode,genType,verbose,alpha,nRefineIter,nRefinePerDim,nIntPerDim)) ;
88 RooNumIntConfig::defaultConfig().methodND().setLabel(proto->IsA()->GetName()) ;
98 RooMCIntegrator::RooMCIntegrator()
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)
118 if(!(_valid= _grid.isValid()))
return;
119 if(_verbose) _grid.Print();
128 RooMCIntegrator::RooMCIntegrator(
const RooAbsFunc&
function,
const RooNumIntConfig& config) :
129 RooAbsIntegrator(function), _grid(function)
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) ;
141 if(!(_valid= _grid.isValid()))
return;
142 if(_verbose) _grid.Print();
151 RooAbsIntegrator* RooMCIntegrator::clone(
const RooAbsFunc&
function,
const RooNumIntConfig& config)
const
153 return new RooMCIntegrator(
function,config) ;
161 RooMCIntegrator::~RooMCIntegrator()
172 Bool_t RooMCIntegrator::checkLimits()
const
174 return _grid.initialize(*integrand());
185 Double_t RooMCIntegrator::integral(
const Double_t* )
188 vegas(AllStages,_nRefinePerDim*_grid.getDimension(),_nRefineIter);
189 Double_t ret = vegas(ReuseGrid,_nIntegratePerDim*_grid.getDimension(),1);
202 Double_t RooMCIntegrator::vegas(Stage stage, UInt_t calls, UInt_t iterations, Double_t *absError)
207 if(stage == AllStages) _grid.initialize(*_function);
210 if(stage <= ReuseGrid) {
219 if(stage <= RefineGrid) {
220 UInt_t bins = RooGrid::maxBins;
222 UInt_t dim(_grid.getDimension());
225 if(_mode != ImportanceOnly) {
228 boxes = (UInt_t)floor(TMath::Power(calls/2.0,1.0/dim));
232 if (2*boxes >= RooGrid::maxBins) {
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;
243 oocxcoutD((TObject*)0,Integration) <<
"RooMCIntegrator: using importance sampling with " << bins <<
" bins and "
244 << boxes <<
" boxes" << endl;
249 Double_t tot_boxes = TMath::Power((Double_t)boxes,(Double_t)dim);
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);
257 _jac = _grid.getVolume()*TMath::Power((Double_t)bins,(Double_t)dim)/calls;
260 _grid.setNBoxes(boxes);
261 if(bins != _grid.getNBins()) _grid.resize(bins);
265 UInt_t *box= _grid.createIndexVector();
266 UInt_t *bin= _grid.createIndexVector();
267 Double_t *x= _grid.createPoint();
270 Double_t cum_int(0),cum_sig(0);
273 for (UInt_t it = 0; it < iterations; it++) {
274 Double_t intgrl(0),intgrl_sq(0),sig(0),jacbin(_jac);
276 _it_num = _it_start + it;
286 for(UInt_t k = 0; k < _calls_per_box; k++) {
289 _grid.generatePoint(box, x, bin, bin_vol, _genType == QuasiRandom ? kTRUE : kFALSE);
291 Double_t fval= jacbin*bin_vol*integrand(x);
293 Double_t d = fval - m;
295 q+= d * d * (k / (k + 1.0));
297 if (_mode != Stratified) _grid.accumulate(bin, fval*fval);
299 intgrl += m * _calls_per_box;
300 Double_t f_sq_sum = q * _calls_per_box ;
304 if (_mode == Stratified) _grid.accumulate(bin, f_sq_sum);
307 if(_timer.RealTime() > 1) {
308 oocoutW((TObject*)0,Integration) <<
"RooMCIntegrator: still working..." << endl;
312 _timer.Start(kFALSE);
315 }
while(_grid.nextBox(box));
319 sig = sig / (_calls_per_box - 1.0) ;
323 else if (_sum_wgts > 0) {
324 wgt = _sum_wgts / _samples;
329 intgrl_sq = intgrl * intgrl;
336 _wtd_int_sum += intgrl * wgt;
337 _chi_sum += intgrl_sq * wgt;
339 cum_int = _wtd_int_sum / _sum_wgts;
340 cum_sig = sqrt (1 / _sum_wgts);
343 _chisq = (_chi_sum - _wtd_int_sum * cum_int)/(_samples - 1.0);
347 cum_int += (intgrl - cum_int) / (it + 1.0);
350 oocxcoutD((TObject*)0,Integration) <<
"=== Iteration " << _it_num <<
" : I = " << intgrl <<
" +/- " << sqrt(sig) << endl
351 <<
" Cumulative : I = " << cum_int <<
" +/- " << cum_sig <<
"( chi2 = " << _chisq
354 if (oodologD((TObject*)0,Integration)) {
355 if(it + 1 == iterations) _grid.Print(
"V");
357 _grid.refine(_alpha);
365 if(absError) *absError = cum_sig;