00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027 #include "RooFit.h"
00028 #include "Riostream.h"
00029
00030 #include "TMath.h"
00031 #include "TClass.h"
00032 #include "RooMCIntegrator.h"
00033 #include "RooArgSet.h"
00034 #include "RooNumber.h"
00035 #include "RooAbsArg.h"
00036 #include "RooNumIntFactory.h"
00037 #include "RooRealVar.h"
00038 #include "RooCategory.h"
00039 #include "RooMsgService.h"
00040
00041 #include <math.h>
00042 #include <assert.h>
00043
00044
00045
00046 ClassImp(RooMCIntegrator)
00047 ;
00048
00049
00050
00051
00052
00053 void RooMCIntegrator::registerIntegrator(RooNumIntFactory& fact)
00054 {
00055
00056
00057
00058 RooCategory samplingMode("samplingMode","Sampling Mode") ;
00059 samplingMode.defineType("Importance",RooMCIntegrator::Importance) ;
00060 samplingMode.defineType("ImportanceOnly",RooMCIntegrator::ImportanceOnly) ;
00061 samplingMode.defineType("Stratified",RooMCIntegrator::Stratified) ;
00062 samplingMode.setIndex(RooMCIntegrator::Importance) ;
00063
00064 RooCategory genType("genType","Generator Type") ;
00065 genType.defineType("QuasiRandom",RooMCIntegrator::QuasiRandom) ;
00066 genType.defineType("PseudoRandom",RooMCIntegrator::PseudoRandom) ;
00067 genType.setIndex(RooMCIntegrator::QuasiRandom) ;
00068
00069 RooCategory verbose("verbose","Verbose flag") ;
00070 verbose.defineType("true",1) ;
00071 verbose.defineType("false",0) ;
00072 verbose.setIndex(0) ;
00073
00074 RooRealVar alpha("alpha","Grid structure constant",1.5) ;
00075 RooRealVar nRefineIter("nRefineIter","Number of refining iterations",5) ;
00076 RooRealVar nRefinePerDim("nRefinePerDim","Number of refining samples (per dimension)",1000) ;
00077 RooRealVar nIntPerDim("nIntPerDim","Number of integration samples (per dimension)",5000) ;
00078
00079
00080 RooMCIntegrator* proto = new RooMCIntegrator() ;
00081
00082
00083 fact.storeProtoIntegrator(proto,RooArgSet(samplingMode,genType,verbose,alpha,nRefineIter,nRefinePerDim,nIntPerDim)) ;
00084
00085
00086 RooNumIntConfig::defaultConfig().methodND().setLabel(proto->IsA()->GetName()) ;
00087 }
00088
00089
00090
00091
00092 RooMCIntegrator::RooMCIntegrator()
00093 {
00094
00095
00096
00097 }
00098
00099
00100
00101
00102 RooMCIntegrator::RooMCIntegrator(const RooAbsFunc& function, SamplingMode mode,
00103 GeneratorType genType, Bool_t verbose) :
00104 RooAbsIntegrator(function), _grid(function), _verbose(verbose),
00105 _alpha(1.5), _mode(mode), _genType(genType),
00106 _nRefineIter(5),_nRefinePerDim(1000),_nIntegratePerDim(5000)
00107 {
00108
00109
00110
00111
00112
00113
00114
00115 if(!(_valid= _grid.isValid())) return;
00116 if(_verbose) _grid.Print();
00117 }
00118
00119
00120
00121
00122 RooMCIntegrator::RooMCIntegrator(const RooAbsFunc& function, const RooNumIntConfig& config) :
00123 RooAbsIntegrator(function), _grid(function)
00124 {
00125
00126
00127
00128 const RooArgSet& configSet = config.getConfigSection(IsA()->GetName()) ;
00129 _verbose = (Bool_t) configSet.getCatIndex("verbose",0) ;
00130 _alpha = configSet.getRealValue("alpha",1.5) ;
00131 _mode = (SamplingMode) configSet.getCatIndex("samplingMode",Importance) ;
00132 _genType = (GeneratorType) configSet.getCatIndex("genType",QuasiRandom) ;
00133 _nRefineIter = (Int_t) configSet.getRealValue("nRefineIter",5) ;
00134 _nRefinePerDim = (Int_t) configSet.getRealValue("nRefinePerDim",1000) ;
00135 _nIntegratePerDim = (Int_t) configSet.getRealValue("nIntPerDim",5000) ;
00136
00137
00138 if(!(_valid= _grid.isValid())) return;
00139 if(_verbose) _grid.Print();
00140 }
00141
00142
00143
00144
00145 RooAbsIntegrator* RooMCIntegrator::clone(const RooAbsFunc& function, const RooNumIntConfig& config) const
00146 {
00147
00148
00149
00150 return new RooMCIntegrator(function,config) ;
00151 }
00152
00153
00154
00155
00156 RooMCIntegrator::~RooMCIntegrator()
00157 {
00158
00159 }
00160
00161
00162
00163
00164 Bool_t RooMCIntegrator::checkLimits() const
00165 {
00166
00167
00168
00169
00170 return _grid.initialize(*integrand());
00171 }
00172
00173
00174
00175
00176 Double_t RooMCIntegrator::integral(const Double_t* )
00177 {
00178
00179
00180
00181
00182
00183 _timer.Start(kTRUE);
00184 vegas(AllStages,_nRefinePerDim*_grid.getDimension(),_nRefineIter);
00185 Double_t ret = vegas(ReuseGrid,_nIntegratePerDim*_grid.getDimension(),1);
00186 return ret ;
00187 }
00188
00189
00190
00191
00192 Double_t RooMCIntegrator::vegas(Stage stage, UInt_t calls, UInt_t iterations, Double_t *absError)
00193 {
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203 if(stage == AllStages) _grid.initialize(*_function);
00204
00205
00206 if(stage <= ReuseGrid) {
00207 _wtd_int_sum = 0;
00208 _sum_wgts = 0;
00209 _chi_sum = 0;
00210 _it_num = 1;
00211 _samples = 0;
00212 }
00213
00214
00215 if(stage <= RefineGrid) {
00216 UInt_t bins = RooGrid::maxBins;
00217 UInt_t boxes = 1;
00218 UInt_t dim(_grid.getDimension());
00219
00220
00221 if(_mode != ImportanceOnly) {
00222
00223
00224 boxes = (UInt_t)floor(TMath::Power(calls/2.0,1.0/dim));
00225
00226
00227 _mode = Importance;
00228 if (2*boxes >= RooGrid::maxBins) {
00229 _mode = Stratified;
00230
00231 Int_t box_per_bin= (boxes > RooGrid::maxBins) ? boxes/RooGrid::maxBins : 1;
00232 bins= boxes/box_per_bin;
00233 if(bins > RooGrid::maxBins) bins= RooGrid::maxBins;
00234 boxes = box_per_bin * bins;
00235 oocxcoutD((TObject*)0,Integration) << "RooMCIntegrator: using stratified sampling with " << bins << " bins and "
00236 << box_per_bin << " boxes/bin" << endl;
00237 }
00238 else {
00239 oocxcoutD((TObject*)0,Integration) << "RooMCIntegrator: using importance sampling with " << bins << " bins and "
00240 << boxes << " boxes" << endl;
00241 }
00242 }
00243
00244
00245 Double_t tot_boxes = TMath::Power((Double_t)boxes,(Double_t)dim);
00246
00247
00248 _calls_per_box = (UInt_t)(calls/tot_boxes);
00249 if(_calls_per_box < 2) _calls_per_box= 2;
00250 calls= (UInt_t)(_calls_per_box*tot_boxes);
00251
00252
00253 _jac = _grid.getVolume()*TMath::Power((Double_t)bins,(Double_t)dim)/calls;
00254
00255
00256 _grid.setNBoxes(boxes);
00257 if(bins != _grid.getNBins()) _grid.resize(bins);
00258 }
00259
00260
00261 UInt_t *box= _grid.createIndexVector();
00262 UInt_t *bin= _grid.createIndexVector();
00263 Double_t *x= _grid.createPoint();
00264
00265
00266 Double_t cum_int(0),cum_sig(0);
00267 _it_start = _it_num;
00268 _chisq = 0.0;
00269 for (UInt_t it = 0; it < iterations; it++) {
00270 Double_t intgrl(0),intgrl_sq(0),sig(0),jacbin(_jac);
00271
00272 _it_num = _it_start + it;
00273
00274
00275 _grid.resetValues();
00276
00277
00278 _grid.firstBox(box);
00279 do {
00280 Double_t m(0),q(0);
00281
00282 for(UInt_t k = 0; k < _calls_per_box; k++) {
00283
00284 Double_t bin_vol(0);
00285 _grid.generatePoint(box, x, bin, bin_vol, _genType == QuasiRandom ? kTRUE : kFALSE);
00286
00287 Double_t fval= jacbin*bin_vol*integrand(x);
00288
00289 Double_t d = fval - m;
00290 m+= d / (k + 1.0);
00291 q+= d * d * (k / (k + 1.0));
00292
00293 if (_mode != Stratified) _grid.accumulate(bin, fval*fval);
00294 }
00295 intgrl += m * _calls_per_box;
00296 Double_t f_sq_sum = q * _calls_per_box ;
00297 sig += f_sq_sum ;
00298
00299
00300 if (_mode == Stratified) _grid.accumulate(bin, f_sq_sum);
00301
00302
00303 if(_timer.RealTime() > 1) {
00304 oocoutW((TObject*)0,Integration) << "RooMCIntegrator: still working..." << endl;
00305 _timer.Start(kTRUE);
00306 }
00307 else {
00308 _timer.Start(kFALSE);
00309 }
00310
00311 } while(_grid.nextBox(box));
00312
00313
00314 Double_t wgt;
00315 sig = sig / (_calls_per_box - 1.0) ;
00316 if (sig > 0) {
00317 wgt = 1.0 / sig;
00318 }
00319 else if (_sum_wgts > 0) {
00320 wgt = _sum_wgts / _samples;
00321 }
00322 else {
00323 wgt = 0.0;
00324 }
00325 intgrl_sq = intgrl * intgrl;
00326 _result = intgrl;
00327 _sigma = sqrt(sig);
00328
00329 if (wgt > 0.0) {
00330 _samples++ ;
00331 _sum_wgts += wgt;
00332 _wtd_int_sum += intgrl * wgt;
00333 _chi_sum += intgrl_sq * wgt;
00334
00335 cum_int = _wtd_int_sum / _sum_wgts;
00336 cum_sig = sqrt (1 / _sum_wgts);
00337
00338 if (_samples > 1) {
00339 _chisq = (_chi_sum - _wtd_int_sum * cum_int)/(_samples - 1.0);
00340 }
00341 }
00342 else {
00343 cum_int += (intgrl - cum_int) / (it + 1.0);
00344 cum_sig = 0.0;
00345 }
00346 oocxcoutD((TObject*)0,Integration) << "=== Iteration " << _it_num << " : I = " << intgrl << " +/- " << sqrt(sig) << endl
00347 << " Cumulative : I = " << cum_int << " +/- " << cum_sig << "( chi2 = " << _chisq
00348 << ")" << endl;
00349
00350 if (oodologD((TObject*)0,Integration)) {
00351 if(it + 1 == iterations) _grid.Print("V");
00352 }
00353 _grid.refine(_alpha);
00354 }
00355
00356
00357 delete[] bin;
00358 delete[] box;
00359 delete[] x;
00360
00361 if(absError) *absError = cum_sig;
00362 return cum_int;
00363 }