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) ;
63 genType.defineType(
"QuasiRandom",RooMCIntegrator::QuasiRandom) ;
64 genType.defineType(
"PseudoRandom",RooMCIntegrator::PseudoRandom) ;
65 genType.setIndex(RooMCIntegrator::QuasiRandom) ;
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) ;
79 return std::make_unique<RooMCIntegrator>(function,config);
83 std::string
name =
"RooMCIntegrator";
102RooMCIntegrator::RooMCIntegrator(
const RooAbsFunc& function, SamplingMode mode,
103 GeneratorType genType,
bool verbose) :
105 _alpha(1.5), _mode(mode), _genType(genType),
106 _nRefineIter(5),_nRefinePerDim(1000),_nIntegratePerDim(5000)
109 if(!(_valid= _grid.isValid()))
return;
110 if(_verbose) _grid.print(std::cout);
125 _mode = (SamplingMode) configSet.
getCatIndex(
"samplingMode",Importance) ;
126 _genType = (GeneratorType) configSet.
getCatIndex(
"genType",QuasiRandom) ;
132 if(!(_valid= _grid.isValid()))
return;
133 if(_verbose) _grid.print(std::cout);
141bool RooMCIntegrator::checkLimits()
const
143 return _grid.initialize(*integrand());
154double RooMCIntegrator::integral(
const double* )
157 vegas(AllStages,_nRefinePerDim*_grid.getDimension(),_nRefineIter);
158 double ret = vegas(ReuseGrid,_nIntegratePerDim*_grid.getDimension(),1);
171double RooMCIntegrator::vegas(Stage stage,
UInt_t calls,
UInt_t iterations,
double *absError)
176 if(stage == AllStages) _grid.initialize(*_function);
179 if(stage <= ReuseGrid) {
188 if(stage <= RefineGrid) {
189 UInt_t bins = RooGrid::maxBins;
194 if(_mode != ImportanceOnly) {
197 boxes = (
UInt_t)
floor(std::pow(calls/2.0,1.0/dim));
201 if (2*boxes >= RooGrid::maxBins) {
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;
212 oocxcoutD((
TObject*)
nullptr,Integration) <<
"RooMCIntegrator: using importance sampling with " << bins <<
" bins and "
213 << boxes <<
" boxes" << std::endl;
218 double tot_boxes = std::pow((
double)boxes,(
double)dim);
221 _calls_per_box = (
UInt_t)(calls/tot_boxes);
222 if(_calls_per_box < 2) _calls_per_box= 2;
223 calls= (
UInt_t)(_calls_per_box*tot_boxes);
226 _jac = _grid.getVolume()*std::pow((
double)bins,(
double)dim)/calls;
229 _grid.setNBoxes(boxes);
230 if(bins != _grid.getNBins()) _grid.resize(bins);
234 std::vector<UInt_t>
box(_grid.getDimension());
235 std::vector<UInt_t> bin(_grid.getDimension());
236 std::vector<double>
x(_grid.getDimension());
244 for (
UInt_t it = 0; it < iterations; it++) {
250 _it_num = _it_start + it;
256 _grid.firstBox(
box.data());
261 for(
UInt_t k = 0; k < _calls_per_box; k++) {
264 _grid.generatePoint(
box.data(),
x.data(), bin.data(), bin_vol, _genType == QuasiRandom ?
true :
false);
266 double fval= jacbin*bin_vol*integrand(
x.data());
270 q+=
d *
d * (k / (k + 1.0));
272 if (_mode != Stratified) _grid.accumulate(bin.data(), fval*fval);
274 intgrl +=
m * _calls_per_box;
275 double f_sq_sum =
q * _calls_per_box ;
279 if (_mode == Stratified) _grid.accumulate(bin.data(), f_sq_sum);
282 if(_timer.RealTime() > 30) {
283 std::size_t index = 0;
284 std::size_t sizeOfDim = 1;
286 for (
unsigned int i=0; i < _grid.getDimension(); ++i) {
287 index +=
box[i] * sizeOfDim;
288 sizeOfDim *= _grid.getNBoxes();
290 oocoutP(
nullptr, Integration) <<
"RooMCIntegrator: still working ... iteration "
291 << it <<
'/' << iterations <<
" box " << index <<
"/"<< std::pow(_grid.getNBoxes(), _grid.getDimension()) << std::endl;
298 }
while(_grid.nextBox(
box.data()));
302 sig = sig / (_calls_per_box - 1.0) ;
306 else if (_sum_wgts > 0) {
307 wgt = _sum_wgts / _samples;
312 intgrl_sq = intgrl * intgrl;
319 _wtd_int_sum += intgrl * wgt;
320 _chi_sum += intgrl_sq * wgt;
322 cum_int = _wtd_int_sum / _sum_wgts;
323 cum_sig =
sqrt (1 / _sum_wgts);
326 _chisq = (_chi_sum - _wtd_int_sum * cum_int)/(_samples - 1.0);
330 cum_int += (intgrl - cum_int) / (it + 1.0);
333 oocxcoutD((
TObject*)
nullptr,Integration) <<
"=== Iteration " << _it_num <<
" : I = " << intgrl <<
" +/- " <<
sqrt(sig) << std::endl
334 <<
" Cumulative : I = " << cum_int <<
" +/- " << cum_sig <<
"( chi2 = " << _chisq
338 if(it + 1 == iterations) _grid.print(std::cout,
true);
340 _grid.refine(_alpha);
343 if(absError) *absError = cum_sig;
int Int_t
Signed integer 4 bytes (int).
unsigned int UInt_t
Unsigned integer 4 bytes (unsigned int).
double getRealValue(const char *name, double defVal=0.0, bool verbose=false) const
Get value of a RooAbsReal stored in set with given name.
Int_t getCatIndex(const char *name, Int_t defVal=0, bool verbose=false) const
Get index value of a RooAbsCategory stored in set with given name.
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
Abstract interface for integrators of real-valued functions that implement the RooAbsFunc interface.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Object to represent discrete states.
bool setLabel(const char *label, bool printError=true) override
Set value by specifying the name of the desired state.
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.
bool registerPlugin(std::string const &name, Creator const &creator, const RooArgSet &defConfig, bool canIntegrate1D, bool canIntegrate2D, bool canIntegrateND, bool canIntegrateOpenEnded, const char *depName="")
Method accepting registration of a prototype numeric integrator along with a RooArgSet of its default...
Variable that can be changed from the outside.
Mother of all ROOT objects.
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
RVec< PromoteType< T > > floor(const RVec< T > &v)
VecExpr< UnaryOp< Sqrt< T >, VecExpr< A, T, D >, T >, T, D > sqrt(const VecExpr< A, T, D > &rhs)
void function(const Char_t *name_, T fun, const Char_t *docstring=0)