CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Attributes
reweight::LumiReWeighting Class Reference

#include <LumiReweightingStandAlone.h>

Public Member Functions

double ITweight (int npv)
 
double ITweight3BX (float ave_int)
 
 LumiReWeighting ()
 
 LumiReWeighting (std::string generatedFile, std::string dataFile, std::string GenHistName, std::string DataHistName)
 
 LumiReWeighting (const std::vector< float > &MC_distr, const std::vector< float > &Lumi_distr)
 
double weight (float n_int)
 
double weight3D (int pv1, int pv2, int pv3)
 
void weight3D_init (float ScaleFactor, std::string WeightOutputFile="")
 
void weight3D_set (std::string WeightFileName)
 
double weightOOT (int npv_in_time, int npv_m50nsBX)
 
void weightOOT_init ()
 

Protected Attributes

TH1F * Data_distr_
 
TFile * dataFile_
 
std::string dataFileName_
 
std::string DataHistName_
 
bool FirstWarning_
 
TFile * generatedFile_
 
std::string generatedFileName_
 
std::string GenHistName_
 
TH1F * MC_distr_
 
double Weight3D_ [50][50][50]
 
double WeightOOTPU_ [25][25]
 
TH1F * weights_
 

Detailed Description

Definition at line 228 of file LumiReweightingStandAlone.h.

Constructor & Destructor Documentation

reweight::LumiReWeighting::LumiReWeighting ( )
inline

Definition at line 231 of file LumiReweightingStandAlone.h.

231 { } ;
reweight::LumiReWeighting::LumiReWeighting ( std::string  generatedFile,
std::string  dataFile,
std::string  GenHistName,
std::string  DataHistName 
)
inline

Definition at line 233 of file LumiReweightingStandAlone.h.

References gather_cfg::cout.

236  :
237  generatedFileName_( generatedFile),
238  dataFileName_ ( dataFile ),
239  GenHistName_ ( GenHistName ),
240  DataHistName_ ( DataHistName )
241  {
242  generatedFile_ = new TFile( generatedFileName_.c_str() ) ; //MC distribution
243  dataFile_ = new TFile( dataFileName_.c_str() ); //Data distribution
244 
245  Data_distr_ = new TH1F( *(static_cast<TH1F*>(dataFile_->Get( DataHistName_.c_str() )->Clone() )) );
246  MC_distr_ = new TH1F( *(static_cast<TH1F*>(generatedFile_->Get( GenHistName_.c_str() )->Clone() )) );
247 
248  // normalize both histograms first
249 
250  Data_distr_->Scale( 1.0/ Data_distr_->Integral() );
251  MC_distr_->Scale( 1.0/ MC_distr_->Integral() );
252 
253  weights_ = new TH1F( *(Data_distr_)) ;
254 
255  // MC * data/MC = data, so the weights are data/MC:
256 
257  weights_->SetName("lumiWeights");
258 
259  TH1F* den = new TH1F(*(MC_distr_));
260 
261  weights_->Divide( den ); // so now the average weight should be 1.0
262 
263  std::cout << " Lumi/Pileup Reweighting: Computed Weights per In-Time Nint " << std::endl;
264 
265  int NBins = weights_->GetNbinsX();
266 
267  for(int ibin = 1; ibin<NBins+1; ++ibin){
268  std::cout << " " << ibin-1 << " " << weights_->GetBinContent(ibin) << std::endl;
269  }
270 
271  weightOOT_init();
272 
273  FirstWarning_ = true;
274 
275  }
reweight::LumiReWeighting::LumiReWeighting ( const std::vector< float > &  MC_distr,
const std::vector< float > &  Lumi_distr 
)
inline

Definition at line 278 of file LumiReweightingStandAlone.h.

References MessageLogger_cfi::cerr, and gather_cfg::cout.

278  {
279  // no histograms for input: use vectors
280 
281  // now, make histograms out of them:
282 
283  // first, check they are the same size...
284 
285  if( MC_distr.size() != Lumi_distr.size() ){
286 
287  std::cerr <<"ERROR: LumiReWeighting: input vectors have different sizes. Quitting... \n";
288  return;
289 
290  }
291 
292  Int_t NBins = MC_distr.size();
293 
294  MC_distr_ = new TH1F("MC_distr","MC dist",NBins,-0.5, float(NBins)-0.5);
295  Data_distr_ = new TH1F("Data_distr","Data dist",NBins,-0.5, float(NBins)-0.5);
296 
297  weights_ = new TH1F("luminumer","luminumer",NBins,-0.5, float(NBins)-0.5);
298  TH1F* den = new TH1F("lumidenom","lumidenom",NBins,-0.5, float(NBins)-0.5);
299 
300  for(int ibin = 1; ibin<NBins+1; ++ibin ) {
301  weights_->SetBinContent(ibin, Lumi_distr[ibin-1]);
302  Data_distr_->SetBinContent(ibin, Lumi_distr[ibin-1]);
303  den->SetBinContent(ibin,MC_distr[ibin-1]);
304  MC_distr_->SetBinContent(ibin,MC_distr[ibin-1]);
305  }
306 
307  // check integrals, make sure things are normalized
308 
309  float deltaH = weights_->Integral();
310  if(fabs(1.0 - deltaH) > 0.02 ) { //*OOPS*...
311  weights_->Scale( 1.0/ deltaH );
312  Data_distr_->Scale( 1.0/ deltaH );
313  }
314  float deltaMC = den->Integral();
315  if(fabs(1.0 - deltaMC) > 0.02 ) {
316  den->Scale(1.0/ deltaMC );
317  MC_distr_->Scale(1.0/ deltaMC );
318  }
319 
320  weights_->Divide( den ); // so now the average weight should be 1.0
321 
322  std::cout << " Lumi/Pileup Reweighting: Computed Weights per In-Time Nint " << std::endl;
323 
324  for(int ibin = 1; ibin<NBins+1; ++ibin){
325  std::cout << " " << ibin-1 << " " << weights_->GetBinContent(ibin) << std::endl;
326  }
327 
328  weightOOT_init();
329 
330  FirstWarning_ = true;
331 
332  }

Member Function Documentation

double reweight::LumiReWeighting::ITweight ( int  npv)
inline

Definition at line 1298 of file LumiReweightingStandAlone.h.

References stringResolutionProvider_cfi::bin.

1298  {
1299  int bin = weights_->GetXaxis()->FindBin( npv );
1300  return weights_->GetBinContent( bin );
1301  }
bin
set the eta bin as selection string.
double reweight::LumiReWeighting::ITweight3BX ( float  ave_int)
inline

Definition at line 1303 of file LumiReweightingStandAlone.h.

References stringResolutionProvider_cfi::bin.

1303  {
1304  int bin = weights_->GetXaxis()->FindBin( ave_int );
1305  return weights_->GetBinContent( bin );
1306  }
bin
set the eta bin as selection string.
double reweight::LumiReWeighting::weight ( float  n_int)
inline

Definition at line 1308 of file LumiReweightingStandAlone.h.

References stringResolutionProvider_cfi::bin.

1308  {
1309  int bin = weights_->GetXaxis()->FindBin( n_int );
1310  return weights_->GetBinContent( bin );
1311  }
bin
set the eta bin as selection string.
double reweight::LumiReWeighting::weight3D ( int  pv1,
int  pv2,
int  pv3 
)
inline

Definition at line 1314 of file LumiReweightingStandAlone.h.

References min().

1314  {
1315 
1316  using std::min;
1317 
1318  int npm1 = min(pv1,34);
1319  int np0 = min(pv2,34);
1320  int npp1 = min(pv3,34);
1321 
1322  return Weight3D_[npm1][np0][npp1];
1323 
1324  }
T min(T a, T b)
Definition: MathUtil.h:58
void reweight::LumiReWeighting::weight3D_init ( float  ScaleFactor,
std::string  WeightOutputFile = "" 
)
inline

Definition at line 334 of file LumiReweightingStandAlone.h.

References runEdmFileComparison::base, gather_cfg::cout, JetChargeProducer_cfi::exp, factorial(), objects.autophobj::float, mps_fire::i, createfilelist::int, gen::k, SiStripPI::mean, min(), lumiCalc2::outfile, and hybridSuperClusters_cfi::xi.

334  {
335 
336  //create histogram to write output weights, save pain of generating them again...
337 
338  TH3D* WHist = new TH3D("WHist","3D weights",50,0.,50.,50,0.,50.,50,0.,50. );
339  TH3D* DHist = new TH3D("DHist","3D weights",50,0.,50.,50,0.,50.,50,0.,50. );
340  TH3D* MHist = new TH3D("MHist","3D weights",50,0.,50.,50,0.,50.,50,0.,50. );
341 
342 
343  using std::min;
344 
345  if( MC_distr_->GetEntries() == 0 ) {
346  std::cout << " MC and Data distributions are not initialized! You must call the LumiReWeighting constructor. " << std::endl;
347  }
348 
349  // arrays for storing number of interactions
350 
351  double MC_ints[50][50][50];
352  double Data_ints[50][50][50];
353 
354  for (int i=0; i<50; i++) {
355  for(int j=0; j<50; j++) {
356  for(int k=0; k<50; k++) {
357  MC_ints[i][j][k] = 0.;
358  Data_ints[i][j][k] = 0.;
359  }
360  }
361  }
362 
363  double factorial[50];
364  double PowerSer[50];
365  double base = 1.;
366 
367  factorial[0] = 1.;
368  PowerSer[0]=1.;
369 
370  for (int i = 1; i<50; ++i) {
371  base = base*float(i);
372  factorial[i] = base;
373  }
374 
375 
376  double x;
377  double xweight;
378  double probi, probj, probk;
379  double Expval, mean;
380  int xi;
381 
382  // Get entries for Data, MC, fill arrays:
383  int NMCbin = MC_distr_->GetNbinsX();
384 
385  for (int jbin=1;jbin<NMCbin+1;jbin++) {
386  x = MC_distr_->GetBinCenter(jbin);
387  xweight = MC_distr_->GetBinContent(jbin); //use as weight for matrix
388 
389  //for Summer 11, we have this int feature:
390  xi = int(x);
391 
392  // Generate Poisson distribution for each value of the mean
393  mean = double(xi);
394 
395  if(mean<0.) {
396  std::cout << "LumiReweighting:BadInputValue" << " Your histogram generates MC luminosity values less than zero!"
397  << " Please Check. Terminating." << std::endl;
398  }
399 
400 
401  if(mean==0.){
402  Expval = 1.;
403  }
404  else {
405  Expval = exp(-1.*mean);
406  }
407 
408  base = 1.;
409 
410  for (int i = 1; i<50; ++i) {
411  base = base*mean;
412  PowerSer[i] = base; // PowerSer is mean^i
413  }
414 
415  // compute poisson probability for each Nvtx in weight matrix
416  for (int i=0; i<50; i++) {
417  probi = PowerSer[i]/factorial[i]*Expval;
418  for(int j=0; j<50; j++) {
419  probj = PowerSer[j]/factorial[j]*Expval;
420  for(int k=0; k<50; k++) {
421  probk = PowerSer[k]/factorial[k]*Expval;
422  // joint probability is product of event weights multiplied by weight of input distribution bin
423  MC_ints[i][j][k] = MC_ints[i][j][k]+probi*probj*probk*xweight;
424  }
425  }
426  }
427 
428  }
429 
430  int NDatabin = Data_distr_->GetNbinsX();
431 
432  for (int jbin=1;jbin<NDatabin+1;jbin++) {
433  mean = (Data_distr_->GetBinCenter(jbin))*ScaleFactor;
434  xweight = Data_distr_->GetBinContent(jbin);
435 
436  // Generate poisson distribution for each value of the mean
437  if(mean<0.) {
438  std::cout << "LumiReweighting:BadInputValue" << " Your histogram generates MC luminosity values less than zero!"
439  << " Please Check. Terminating." << std::endl;
440  }
441 
442  if(mean==0.){
443  Expval = 1.;
444  }
445  else {
446  Expval = exp(-1.*mean);
447  }
448 
449  base = 1.;
450 
451  for (int i = 1; i<50; ++i) {
452  base = base*mean;
453  PowerSer[i] = base;
454  }
455 
456  // compute poisson probability for each Nvtx in weight matrix
457 
458  for (int i=0; i<50; i++) {
459  probi = PowerSer[i]/factorial[i]*Expval;
460  for(int j=0; j<50; j++) {
461  probj = PowerSer[j]/factorial[j]*Expval;
462  for(int k=0; k<50; k++) {
463  probk = PowerSer[k]/factorial[k]*Expval;
464  // joint probability is product of event weights multiplied by weight of input distribution bin
465  Data_ints[i][j][k] = Data_ints[i][j][k]+probi*probj*probk*xweight;
466  }
467  }
468  }
469 
470  }
471 
472 
473  for (int i=0; i<50; i++) {
474  //if(i<5) std::cout << "i = " << i << std::endl;
475  for(int j=0; j<50; j++) {
476  for(int k=0; k<50; k++) {
477  if( (MC_ints[i][j][k])>0.) {
478  Weight3D_[i][j][k] = Data_ints[i][j][k]/MC_ints[i][j][k];
479  }
480  else {
481  Weight3D_[i][j][k] = 0.;
482  }
483  WHist->SetBinContent( i+1,j+1,k+1,Weight3D_[i][j][k] );
484  DHist->SetBinContent( i+1,j+1,k+1,Data_ints[i][j][k] );
485  MHist->SetBinContent( i+1,j+1,k+1,MC_ints[i][j][k] );
486  // if(i<5 && j<5 && k<5) std::cout << Weight3D_[i][j][k] << " " ;
487  }
488  // if(i<5 && j<5) std::cout << std::endl;
489  }
490  }
491 
492  if(! WeightOutputFile.empty() ) {
493  std::cout << " 3D Weight Matrix initialized! " << std::endl;
494  std::cout << " Writing weights to file " << WeightOutputFile << " for re-use... " << std::endl;
495 
496 
497  TFile * outfile = new TFile(WeightOutputFile.c_str(),"RECREATE");
498  WHist->Write();
499  MHist->Write();
500  DHist->Write();
501  outfile->Write();
502  outfile->Close();
503  outfile->Delete();
504  }
505 
506  return;
507  }
T min(T a, T b)
Definition: MathUtil.h:58
base
Make Sure CMSSW is Setup ##.
int k[5][pyjets_maxn]
int factorial(int n)
factorial function
void reweight::LumiReWeighting::weight3D_set ( std::string  WeightFileName)
inline

Definition at line 510 of file LumiReweightingStandAlone.h.

References gather_cfg::cout, mps_fire::i, and gen::k.

510  {
511 
512  TFile *infile = new TFile(WeightFileName.c_str());
513  TH1F *WHist = (TH1F*)infile->Get("WHist");
514 
515  // Check if the histogram exists
516  if (!WHist) {
517  std::cout << " Could not find the histogram WHist in the file "
518  << "in the file " << WeightFileName << "." << std::endl;
519  return;
520  }
521 
522  for (int i=0; i<50; i++) {
523  for(int j=0; j<50; j++) {
524  for(int k=0; k<50; k++) {
525  Weight3D_[i][j][k] = WHist->GetBinContent(i,j,k);
526  }
527  }
528  }
529 
530  std::cout << " 3D Weight Matrix initialized! " << std::endl;
531 
532  return;
533 
534 
535  }
int k[5][pyjets_maxn]
double reweight::LumiReWeighting::weightOOT ( int  npv_in_time,
int  npv_m50nsBX 
)
inline

Definition at line 1328 of file LumiReweightingStandAlone.h.

References stringResolutionProvider_cfi::bin, MessageLogger_cfi::cerr, and gather_cfg::cout.

1328  {
1329 
1330  static const double Correct_Weights2011[25] = { // residual correction to match lumi spectrum
1331  5.30031,
1332  2.07903,
1333  1.40729,
1334  1.27687,
1335  1.0702,
1336  0.902094,
1337  0.902345,
1338  0.931449,
1339  0.78202,
1340  0.824686,
1341  0.837735,
1342  0.910261,
1343  1.01394,
1344  1.1599,
1345  1.12778,
1346  1.58423,
1347  1.78868,
1348  1.58296,
1349  2.3291,
1350  3.86641,
1351  0,
1352  0,
1353  0,
1354  0,
1355  0
1356  };
1357 
1358 
1359  if(FirstWarning_) {
1360 
1361  std::cout << " **** Warning: Out-of-time pileup reweighting appropriate only for PU_S3 **** " << std::endl;
1362  std::cout << " **** will be applied **** " << std::endl;
1363 
1364  FirstWarning_ = false;
1365 
1366  }
1367 
1368 
1369  // Note: for the "uncorrelated" out-of-time pileup, reweighting is only done on the 50ns
1370  // "late" bunch (BX=+1), since that is basically the only one that matters in terms of
1371  // energy deposition.
1372 
1373  if(npv_in_time < 0) {
1374  std::cerr << " no in-time beam crossing found\n! " ;
1375  std::cerr << " Returning event weight=0\n! ";
1376  return 0.;
1377  }
1378  if(npv_m50nsBX < 0) {
1379  std::cerr << " no out-of-time beam crossing found\n! " ;
1380  std::cerr << " Returning event weight=0\n! ";
1381  return 0.;
1382  }
1383 
1384  int bin = weights_->GetXaxis()->FindBin( npv_in_time );
1385 
1386  double inTimeWeight = weights_->GetBinContent( bin );
1387 
1388  double TotalWeight = 1.0;
1389 
1390 
1391  TotalWeight = inTimeWeight * WeightOOTPU_[bin-1][npv_m50nsBX] * Correct_Weights2011[bin-1];
1392 
1393 
1394  return TotalWeight;
1395 
1396  }
bin
set the eta bin as selection string.
void reweight::LumiReWeighting::weightOOT_init ( )
inline

Definition at line 539 of file LumiReweightingStandAlone.h.

References iint.

539  {
540 
541  // The following are poisson distributions with different means, where the maximum
542  // of the function has been normalized to weight 1.0
543  // These are used to reweight the out-of-time pileup to match the in-time distribution.
544  // The total event weight is the product of the in-time weight, the out-of-time weight,
545  // and a residual correction to fix the distortions caused by the fact that the out-of-time
546  // distribution is not flat.
547 
548  static const double weight_24[25] = {
549  0,
550  0,
551  0,
552  0,
553  2.46277e-06,
554  2.95532e-05,
555  0.000104668,
556  0.000401431,
557  0.00130034,
558  0.00342202,
559  0.00818132,
560  0.0175534,
561  0.035784,
562  0.0650836,
563  0.112232,
564  0.178699,
565  0.268934,
566  0.380868,
567  0.507505,
568  0.640922,
569  0.768551,
570  0.877829,
571  0.958624,
572  0.99939,
573  1
574  };
575 
576  static const double weight_23[25] = {
577  0,
578  1.20628e-06,
579  1.20628e-06,
580  2.41255e-06,
581  1.20628e-05,
582  6.39326e-05,
583  0.000252112,
584  0.000862487,
585  0.00244995,
586  0.00616527,
587  0.0140821,
588  0.0293342,
589  0.0564501,
590  0.100602,
591  0.164479,
592  0.252659,
593  0.36268,
594  0.491427,
595  0.627979,
596  0.75918,
597  0.873185,
598  0.957934,
599  0.999381,
600  1,
601  0.957738
602  };
603 
604  static const double weight_22[25] = {
605  0,
606  0,
607  0,
608  5.88636e-06,
609  3.0609e-05,
610  0.000143627,
611  0.000561558,
612  0.00173059,
613  0.00460078,
614  0.0110616,
615  0.0238974,
616  0.0475406,
617  0.0875077,
618  0.148682,
619  0.235752,
620  0.343591,
621  0.473146,
622  0.611897,
623  0.748345,
624  0.865978,
625  0.953199,
626  0.997848,
627  1,
628  0.954245,
629  0.873688
630  };
631 
632  static const double weight_21[25] = {
633  0,
634  0,
635  1.15381e-06,
636  8.07665e-06,
637  7.1536e-05,
638  0.000280375,
639  0.00107189,
640  0.00327104,
641  0.00809396,
642  0.0190978,
643  0.0401894,
644  0.0761028,
645  0.13472,
646  0.216315,
647  0.324649,
648  0.455125,
649  0.598241,
650  0.739215,
651  0.861866,
652  0.953911,
653  0.998918,
654  1,
655  0.956683,
656  0.872272,
657  0.76399
658  };
659 
660 
661  static const double weight_20[25] = {
662  0,
663  0,
664  1.12532e-06,
665  2.58822e-05,
666  0.000145166,
667  0.000633552,
668  0.00215048,
669  0.00592816,
670  0.0145605,
671  0.0328367,
672  0.0652649,
673  0.11893,
674  0.19803,
675  0.305525,
676  0.436588,
677  0.581566,
678  0.727048,
679  0.8534,
680  0.949419,
681  0.999785,
682  1,
683  0.953008,
684  0.865689,
685  0.753288,
686  0.62765
687  };
688  static const double weight_19[25] = {
689  0,
690  0,
691  1.20714e-05,
692  5.92596e-05,
693  0.000364337,
694  0.00124994,
695  0.00403953,
696  0.0108149,
697  0.025824,
698  0.0544969,
699  0.103567,
700  0.17936,
701  0.283532,
702  0.416091,
703  0.562078,
704  0.714714,
705  0.846523,
706  0.947875,
707  1,
708  0.999448,
709  0.951404,
710  0.859717,
711  0.742319,
712  0.613601,
713  0.48552
714  };
715 
716  static const double weight_18[25] = {
717  0,
718  3.20101e-06,
719  2.88091e-05,
720  0.000164319,
721  0.000719161,
722  0.00250106,
723  0.00773685,
724  0.0197513,
725  0.0443693,
726  0.0885998,
727  0.159891,
728  0.262607,
729  0.392327,
730  0.543125,
731  0.69924,
732  0.837474,
733  0.943486,
734  0.998029,
735  1,
736  0.945937,
737  0.851807,
738  0.729309,
739  0.596332,
740  0.467818,
741  0.350434
742  };
743 
744 
745  static const double weight_17[25] = {
746  1.03634e-06,
747  7.25437e-06,
748  4.97443e-05,
749  0.000340956,
750  0.00148715,
751  0.00501485,
752  0.0143067,
753  0.034679,
754  0.0742009,
755  0.140287,
756  0.238288,
757  0.369416,
758  0.521637,
759  0.682368,
760  0.828634,
761  0.939655,
762  1,
763  0.996829,
764  0.94062,
765  0.841575,
766  0.716664,
767  0.582053,
768  0.449595,
769  0.331336,
770  0.234332
771  };
772 
773 
774  static const double weight_16[25] = {
775  4.03159e-06,
776  2.41895e-05,
777  0.000141106,
778  0.00081942,
779  0.00314565,
780  0.00990662,
781  0.026293,
782  0.0603881,
783  0.120973,
784  0.214532,
785  0.343708,
786  0.501141,
787  0.665978,
788  0.820107,
789  0.938149,
790  1,
791  0.99941,
792  0.940768,
793  0.837813,
794  0.703086,
795  0.564023,
796  0.42928,
797  0.312515,
798  0.216251,
799  0.14561
800  };
801 
802 
803  static const double weight_15[25] = {
804  9.76084e-07,
805  5.07564e-05,
806  0.000303562,
807  0.00174036,
808  0.00617959,
809  0.0188579,
810  0.047465,
811  0.101656,
812  0.189492,
813  0.315673,
814  0.474383,
815  0.646828,
816  0.809462,
817  0.934107,
818  0.998874,
819  1,
820  0.936163,
821  0.827473,
822  0.689675,
823  0.544384,
824  0.40907,
825  0.290648,
826  0.198861,
827  0.12951,
828  0.0808051
829  };
830 
831 
832  static const double weight_14[25] = {
833  1.13288e-05,
834  0.000124617,
835  0.000753365,
836  0.00345056,
837  0.0123909,
838  0.0352712,
839  0.0825463,
840  0.16413,
841  0.287213,
842  0.44615,
843  0.625826,
844  0.796365,
845  0.930624,
846  0.999958,
847  1,
848  0.934414,
849  0.816456,
850  0.672939,
851  0.523033,
852  0.386068,
853  0.269824,
854  0.180342,
855  0.114669,
856  0.0698288,
857  0.0406496
858  };
859 
860 
861  static const double weight_13[25] = {
862  2.54296e-05,
863  0.000261561,
864  0.00167018,
865  0.00748083,
866  0.0241308,
867  0.0636801,
868  0.138222,
869  0.255814,
870  0.414275,
871  0.600244,
872  0.779958,
873  0.92256,
874  0.999155,
875  1,
876  0.927126,
877  0.804504,
878  0.651803,
879  0.497534,
880  0.35976,
881  0.245834,
882  0.160904,
883  0.0991589,
884  0.0585434,
885  0.0332437,
886  0.0180159
887  };
888 
889  static const double weight_12[25] = {
890  5.85742e-05,
891  0.000627706,
892  0.00386677,
893  0.0154068,
894  0.0465892,
895  0.111683,
896  0.222487,
897  0.381677,
898  0.5719,
899  0.765001,
900  0.915916,
901  1,
902  0.999717,
903  0.921443,
904  0.791958,
905  0.632344,
906  0.475195,
907  0.334982,
908  0.223666,
909  0.141781,
910  0.0851538,
911  0.048433,
912  0.0263287,
913  0.0133969,
914  0.00696683
915  };
916 
917 
918  static const double weight_11[25] = {
919  0.00015238,
920  0.00156064,
921  0.00846044,
922  0.0310939,
923  0.0856225,
924  0.187589,
925  0.343579,
926  0.541892,
927  0.74224,
928  0.909269,
929  0.998711,
930  1,
931  0.916889,
932  0.77485,
933  0.608819,
934  0.447016,
935  0.307375,
936  0.198444,
937  0.121208,
938  0.070222,
939  0.0386492,
940  0.0201108,
941  0.0100922,
942  0.00484937,
943  0.00222458
944  };
945 
946  static const double weight_10[25] = {
947  0.000393044,
948  0.00367001,
949  0.0179474,
950  0.060389,
951  0.151477,
952  0.302077,
953  0.503113,
954  0.720373,
955  0.899568,
956  1,
957  0.997739,
958  0.909409,
959  0.75728,
960  0.582031,
961  0.415322,
962  0.277663,
963  0.174147,
964  0.102154,
965  0.0566719,
966  0.0298642,
967  0.0147751,
968  0.00710995,
969  0.00319628,
970  0.00140601,
971  0.000568796
972  };
973 
974 
975  static const double weight_9[25] = {
976  0.00093396,
977  0.00854448,
978  0.0380306,
979  0.113181,
980  0.256614,
981  0.460894,
982  0.690242,
983  0.888781,
984  1,
985  0.998756,
986  0.899872,
987  0.735642,
988  0.552532,
989  0.382726,
990  0.246114,
991  0.147497,
992  0.0825541,
993  0.0441199,
994  0.0218157,
995  0.0103578,
996  0.00462959,
997  0.0019142,
998  0.000771598,
999  0.000295893,
1000  0.000111529
1001  };
1002 
1003 
1004  static const double weight_8[25] = {
1005  0.00240233,
1006  0.0192688,
1007  0.0768653,
1008  0.205008,
1009  0.410958,
1010  0.65758,
1011  0.875657,
1012  0.999886,
1013  1,
1014  0.889476,
1015  0.711446,
1016  0.517781,
1017  0.345774,
1018  0.212028,
1019  0.121208,
1020  0.0644629,
1021  0.0324928,
1022  0.0152492,
1023  0.00673527,
1024  0.0028547,
1025  0.00117213,
1026  0.000440177,
1027  0.000168471,
1028  5.80689e-05,
1029  1.93563e-05
1030  };
1031 
1032  static const double weight_7[25] = {
1033  0.00617233,
1034  0.0428714,
1035  0.150018,
1036  0.350317,
1037  0.612535,
1038  0.856525,
1039  0.999923,
1040  1,
1041  0.87544,
1042  0.679383,
1043  0.478345,
1044  0.303378,
1045  0.176923,
1046  0.0950103,
1047  0.0476253,
1048  0.0222211,
1049  0.00972738,
1050  0.00392962,
1051  0.0015258,
1052  0.000559168,
1053  0.000183928,
1054  6.77983e-05,
1055  1.67818e-05,
1056  7.38398e-06,
1057  6.71271e-07
1058  };
1059 
1060  static const double weight_6[25] = {
1061  0.0154465,
1062  0.0923472,
1063  0.277322,
1064  0.55552,
1065  0.833099,
1066  0.999035,
1067  1,
1068  0.855183,
1069  0.641976,
1070  0.428277,
1071  0.256804,
1072  0.139798,
1073  0.0700072,
1074  0.0321586,
1075  0.0137971,
1076  0.00544756,
1077  0.00202316,
1078  0.000766228,
1079  0.000259348,
1080  8.45836e-05,
1081  1.80362e-05,
1082  8.70713e-06,
1083  3.73163e-06,
1084  6.21938e-07,
1085  0
1086  };
1087 
1088 
1089  static const double weight_5[25] = {
1090  0.0382845,
1091  0.191122,
1092  0.478782,
1093  0.797314,
1094  1,
1095  0.997148,
1096  0.831144,
1097  0.59461,
1098  0.371293,
1099  0.205903,
1100  0.103102,
1101  0.0471424,
1102  0.0194997,
1103  0.00749415,
1104  0.00273709,
1105  0.000879189,
1106  0.000286049,
1107  0.000102364,
1108  1.70606e-05,
1109  3.98081e-06,
1110  2.27475e-06,
1111  0,
1112  0,
1113  0,
1114  0
1115  };
1116 
1117 
1118  static const double weight_4[25] = {
1119  0.0941305,
1120  0.373824,
1121  0.750094,
1122  1,
1123  0.997698,
1124  0.800956,
1125  0.532306,
1126  0.304597,
1127  0.152207,
1128  0.0676275,
1129  0.0270646,
1130  0.00975365,
1131  0.00326077,
1132  0.00101071,
1133  0.000301781,
1134  7.41664e-05,
1135  1.58563e-05,
1136  3.58045e-06,
1137  1.02299e-06,
1138  0,
1139  5.11493e-07,
1140  0,
1141  0,
1142  0,
1143  0
1144  };
1145 
1146 
1147  static const double weight_3[25] = {
1148  0.222714,
1149  0.667015,
1150  1,
1151  0.999208,
1152  0.750609,
1153  0.449854,
1154  0.224968,
1155  0.0965185,
1156  0.0361225,
1157  0.012084,
1158  0.00359618,
1159  0.000977166,
1160  0.000239269,
1161  6.29422e-05,
1162  1.16064e-05,
1163  1.78559e-06,
1164  0,
1165  4.46398e-07,
1166  0,
1167  0,
1168  0,
1169  0,
1170  0,
1171  0,
1172  0
1173  };
1174 
1175  static const double weight_2[25] = {
1176  0.499541,
1177  0.999607,
1178  1,
1179  0.666607,
1180  0.333301,
1181  0.13279,
1182  0.0441871,
1183  0.0127455,
1184  0.00318434,
1185  0.00071752,
1186  0.000132204,
1187  2.69578e-05,
1188  5.16999e-06,
1189  2.21571e-06,
1190  0,
1191  0,
1192  0,
1193  0,
1194  0,
1195  0,
1196  0,
1197  0,
1198  0,
1199  0,
1200  0
1201  };
1202 
1203  static const double weight_1[25] = {
1204  0.999165,
1205  1,
1206  0.499996,
1207  0.166868,
1208  0.0414266,
1209  0.00831053,
1210  0.00137472,
1211  0.000198911,
1212  2.66302e-05,
1213  2.44563e-06,
1214  2.71737e-07,
1215  2.71737e-07,
1216  0,
1217  0,
1218  0,
1219  0,
1220  0,
1221  0,
1222  0,
1223  0,
1224  0,
1225  0,
1226  0,
1227  0,
1228  0
1229  };
1230 
1231  static const double weight_0[25] = {
1232  1,
1233  0,
1234  0,
1235  0,
1236  0,
1237  0,
1238  0,
1239  0,
1240  0,
1241  0,
1242  0,
1243  0,
1244  0,
1245  0,
1246  0,
1247  0,
1248  0,
1249  0,
1250  0,
1251  0,
1252  0,
1253  0,
1254  0,
1255  0,
1256  0
1257  };
1258 
1259  //WeightOOTPU_ = {0};
1260 
1261  const double* WeightPtr = 0;
1262 
1263  for(int iint = 0; iint<25; ++iint){
1264  if(iint ==0) WeightPtr = weight_0;
1265  if(iint ==1) WeightPtr = weight_1;
1266  if(iint ==2) WeightPtr = weight_2;
1267  if(iint ==3) WeightPtr = weight_3;
1268  if(iint ==4) WeightPtr = weight_4;
1269  if(iint ==5) WeightPtr = weight_5;
1270  if(iint ==6) WeightPtr = weight_6;
1271  if(iint ==7) WeightPtr = weight_7;
1272  if(iint ==8) WeightPtr = weight_8;
1273  if(iint ==9) WeightPtr = weight_9;
1274  if(iint ==10) WeightPtr = weight_10;
1275  if(iint ==11) WeightPtr = weight_11;
1276  if(iint ==12) WeightPtr = weight_12;
1277  if(iint ==13) WeightPtr = weight_13;
1278  if(iint ==14) WeightPtr = weight_14;
1279  if(iint ==15) WeightPtr = weight_15;
1280  if(iint ==16) WeightPtr = weight_16;
1281  if(iint ==17) WeightPtr = weight_17;
1282  if(iint ==18) WeightPtr = weight_18;
1283  if(iint ==19) WeightPtr = weight_19;
1284  if(iint ==20) WeightPtr = weight_20;
1285  if(iint ==21) WeightPtr = weight_21;
1286  if(iint ==22) WeightPtr = weight_22;
1287  if(iint ==23) WeightPtr = weight_23;
1288  if(iint ==24) WeightPtr = weight_24;
1289 
1290  for(int ibin = 0; ibin<25; ++ibin){
1291  WeightOOTPU_[iint][ibin] = *(WeightPtr+ibin);
1292  }
1293  }
1294 
1295  }
int iint

Member Data Documentation

TH1F* reweight::LumiReWeighting::Data_distr_
protected

Definition at line 1410 of file LumiReweightingStandAlone.h.

TFile* reweight::LumiReWeighting::dataFile_
protected

Definition at line 1405 of file LumiReweightingStandAlone.h.

std::string reweight::LumiReWeighting::dataFileName_
protected

Definition at line 1401 of file LumiReweightingStandAlone.h.

std::string reweight::LumiReWeighting::DataHistName_
protected

Definition at line 1403 of file LumiReweightingStandAlone.h.

bool reweight::LumiReWeighting::FirstWarning_
protected

Definition at line 1415 of file LumiReweightingStandAlone.h.

TFile* reweight::LumiReWeighting::generatedFile_
protected

Definition at line 1404 of file LumiReweightingStandAlone.h.

std::string reweight::LumiReWeighting::generatedFileName_
protected

Definition at line 1400 of file LumiReweightingStandAlone.h.

std::string reweight::LumiReWeighting::GenHistName_
protected

Definition at line 1402 of file LumiReweightingStandAlone.h.

TH1F* reweight::LumiReWeighting::MC_distr_
protected

Definition at line 1409 of file LumiReweightingStandAlone.h.

double reweight::LumiReWeighting::Weight3D_[50][50][50]
protected

Definition at line 1413 of file LumiReweightingStandAlone.h.

double reweight::LumiReWeighting::WeightOOTPU_[25][25]
protected

Definition at line 1412 of file LumiReweightingStandAlone.h.

TH1F* reweight::LumiReWeighting::weights_
protected

Definition at line 1406 of file LumiReweightingStandAlone.h.