CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/PhysicsTools/Utilities/src/LumiReWeighting.cc

Go to the documentation of this file.
00001 #ifndef PhysicsTools_Utilities_interface_LumiReWeighting_cc
00002 #define PhysicsTools_Utilities_interface_LumiReWeighting_cc
00003 
00004 
00020 #include "TRandom1.h"
00021 #include "TRandom2.h"
00022 #include "TRandom3.h"
00023 #include "TStopwatch.h"
00024 #include "TH1.h"
00025 #include "TFile.h"
00026 #include <string>
00027 #include <algorithm>
00028 #include <boost/shared_ptr.hpp>
00029 #include "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h" 
00030 #include "PhysicsTools/Utilities/interface/LumiReWeighting.h"
00031 #include "DataFormats/Common/interface/Handle.h"
00032 #include "DataFormats/Provenance/interface/ProcessHistory.h"
00033 #include "DataFormats/Provenance/interface/Provenance.h"
00034 #include "FWCore/Common/interface/EventBase.h"
00035 
00036 using namespace edm;
00037 
00038 LumiReWeighting::LumiReWeighting( std::string generatedFile,
00039                    std::string dataFile,
00040                    std::string GenHistName = "pileup",
00041                    std::string DataHistName = "pileup" ) :
00042       generatedFileName_( generatedFile), 
00043       dataFileName_     ( dataFile ), 
00044       GenHistName_        ( GenHistName ), 
00045       DataHistName_        ( DataHistName )
00046       {
00047         generatedFile_ = boost::shared_ptr<TFile>( new TFile(generatedFileName_.c_str()) ); //MC distribution
00048         dataFile_      = boost::shared_ptr<TFile>( new TFile(dataFileName_.c_str()) );      //Data distribution
00049 
00050         Data_distr_ = boost::shared_ptr<TH1>(  (static_cast<TH1*>(dataFile_->Get( DataHistName_.c_str() )->Clone() )) );
00051         MC_distr_ = boost::shared_ptr<TH1>(  (static_cast<TH1*>(generatedFile_->Get( GenHistName_.c_str() )->Clone() )) );
00052 
00053         // MC * data/MC = data, so the weights are data/MC:
00054 
00055         // normalize both histograms first
00056 
00057         Data_distr_->Scale( 1.0/ Data_distr_->Integral() );
00058         MC_distr_->Scale( 1.0/ MC_distr_->Integral() );
00059 
00060         weights_ = boost::shared_ptr<TH1>( static_cast<TH1*>(Data_distr_->Clone()) );
00061 
00062         weights_->SetName("lumiWeights");
00063 
00064         TH1* den = dynamic_cast<TH1*>(MC_distr_->Clone());
00065 
00066         //den->Scale(1.0/ den->Integral());
00067 
00068         weights_->Divide( den );  // so now the average weight should be 1.0
00069 
00070         std::cout << " Lumi/Pileup Reweighting: Computed Weights per In-Time Nint " << std::endl;
00071 
00072         int NBins = weights_->GetNbinsX();
00073 
00074         for(int ibin = 1; ibin<NBins+1; ++ibin){
00075           std::cout << "   " << ibin-1 << " " << weights_->GetBinContent(ibin) << std::endl;
00076         }
00077 
00078 
00079         weightOOT_init();
00080 
00081         FirstWarning_ = true;
00082         OldLumiSection_ = -1;
00083 }
00084 
00085 LumiReWeighting::LumiReWeighting( std::vector< float > MC_distr, std::vector< float > Lumi_distr) {
00086   // no histograms for input: use vectors
00087   
00088   // now, make histograms out of them:
00089 
00090   // first, check they are the same size...
00091 
00092   if( MC_distr.size() != Lumi_distr.size() ){   
00093 
00094     std::cerr <<"ERROR: LumiReWeighting: input vectors have different sizes. Quitting... \n";
00095     return;
00096 
00097   }
00098 
00099   Int_t NBins = MC_distr.size();
00100 
00101   MC_distr_ = boost::shared_ptr<TH1> ( new TH1F("MC_distr","MC dist",NBins,-0.5, float(NBins)-0.5) );
00102   Data_distr_ = boost::shared_ptr<TH1> ( new TH1F("Data_distr","Data dist",NBins,-0.5, float(NBins)-0.5) );
00103 
00104   weights_ = boost::shared_ptr<TH1> ( new TH1F("luminumer","luminumer",NBins,-0.5, float(NBins)-0.5) );
00105   TH1* den = new TH1F("lumidenom","lumidenom",NBins,-0.5, float(NBins)-0.5) ;
00106 
00107   for(int ibin = 1; ibin<NBins+1; ++ibin ) {
00108     weights_->SetBinContent(ibin, Lumi_distr[ibin-1]);
00109     Data_distr_->SetBinContent(ibin, Lumi_distr[ibin-1]);
00110     den->SetBinContent(ibin,MC_distr[ibin-1]);
00111     MC_distr_->SetBinContent(ibin,MC_distr[ibin-1]);
00112   }
00113 
00114   // check integrals, make sure things are normalized
00115 
00116   float deltaH = weights_->Integral();
00117   if(fabs(1.0 - deltaH) > 0.02 ) { //*OOPS*...
00118     weights_->Scale( 1.0/ weights_->Integral() );
00119     Data_distr_->Scale( 1.0/ Data_distr_->Integral() );
00120   }
00121   float deltaMC = den->Integral();
00122   if(fabs(1.0 - deltaMC) > 0.02 ) {
00123     den->Scale(1.0/ den->Integral());
00124     MC_distr_->Scale(1.0/ MC_distr_->Integral());
00125   }
00126 
00127   weights_->Divide( den );  // so now the average weight should be 1.0    
00128 
00129   std::cout << " Lumi/Pileup Reweighting: Computed Weights per In-Time Nint " << std::endl;
00130 
00131   for(int ibin = 1; ibin<NBins+1; ++ibin){
00132     std::cout << "   " << ibin-1 << " " << weights_->GetBinContent(ibin) << std::endl;
00133   }
00134 
00135   weightOOT_init();
00136 
00137   FirstWarning_ = true;
00138   OldLumiSection_ = -1;
00139 }
00140 
00141 double LumiReWeighting::weight( int npv ) {
00142   int bin = weights_->GetXaxis()->FindBin( npv );
00143   return weights_->GetBinContent( bin );
00144 }
00145 
00146 double LumiReWeighting::weight( float npv ) {
00147   int bin = weights_->GetXaxis()->FindBin( npv );
00148   return weights_->GetBinContent( bin );
00149 }
00150 
00151 
00152 
00153 // This version of weight does all of the work for you, assuming you want to re-weight
00154 // using the true number of interactions in the in-time beam crossing.
00155 
00156 
00157 double LumiReWeighting::weight( const edm::EventBase &e ) {
00158 
00159   // find provenance of event objects, just to check at the job beginning if there might be an issue  
00160 
00161   if(FirstWarning_) {
00162 
00163     edm::ReleaseVersion TargetRelease = edm::ReleaseVersion("\"CMSSW_4_2_2_patch2\"");
00164     edm::ProcessHistory PHist = e.processHistory();
00165     edm::ProcessHistory::const_iterator PHist_iter = PHist.begin();
00166 
00167     for(; PHist_iter<PHist.end() ;++PHist_iter) {
00168       edm::ProcessConfiguration PConf = *(PHist_iter);
00169       edm::ReleaseVersion Release =  PConf.releaseVersion() ;
00170       const std::string Process =  PConf.processName();
00171 
00172       if((Process=="HLT") && (Release==TargetRelease)) {
00173         std::cout << " **** Warning! You are using in-time-only pileup reweighting for Release " << Release << " **** " << std::endl;
00174         std::cout << " **** There is a weightOOT() function available if needed ****  " << std::endl;
00175       }
00176     }
00177     //    SetFirstFalse();
00178     FirstWarning_ = false;
00179   }
00180 
00181   // get pileup summary information
00182 
00183   Handle<std::vector< PileupSummaryInfo > >  PupInfo;
00184   e.getByLabel(edm::InputTag("addPileupInfo"), PupInfo);
00185 
00186   std::vector<PileupSummaryInfo>::const_iterator PVI;
00187 
00188   int npv = -1;
00189   for(PVI = PupInfo->begin(); PVI != PupInfo->end(); ++PVI) {
00190 
00191     int BX = PVI->getBunchCrossing();
00192 
00193     if(BX == 0) { 
00194       npv = PVI->getPU_NumInteractions();
00195       continue;
00196     }
00197 
00198   }
00199 
00200   if(npv < 0) std::cerr << " no in-time beam crossing found\n! " ;
00201 
00202   return weight(npv);
00203 
00204 }
00205 
00206 void LumiReWeighting::weightOOT_init() {
00207 
00208   // The following are poisson distributions with different means, where the maximum
00209   // of the function has been normalized to weight 1.0
00210   // These are used to reweight the out-of-time pileup to match the in-time distribution.
00211   // The total event weight is the product of the in-time weight, the out-of-time weight,
00212   // and a residual correction to fix the distortions caused by the fact that the out-of-time
00213   // distribution is not flat.
00214 
00215   static double weight_24[25] = {
00216     0,
00217     0,
00218     0,
00219     0,
00220     2.46277e-06,
00221     2.95532e-05,
00222     0.000104668,
00223     0.000401431,
00224     0.00130034,
00225     0.00342202,
00226     0.00818132,
00227     0.0175534,
00228     0.035784,
00229     0.0650836,
00230     0.112232,
00231     0.178699,
00232     0.268934,
00233     0.380868,
00234     0.507505,
00235     0.640922,
00236     0.768551,
00237     0.877829,
00238     0.958624,
00239     0.99939,
00240     1
00241   };
00242 
00243   static double weight_23[25] = {
00244     0,
00245     1.20628e-06,
00246     1.20628e-06,
00247     2.41255e-06,
00248     1.20628e-05,
00249     6.39326e-05,
00250     0.000252112,
00251     0.000862487,
00252     0.00244995,
00253     0.00616527,
00254     0.0140821,
00255     0.0293342,
00256     0.0564501,
00257     0.100602,
00258     0.164479,
00259     0.252659,
00260     0.36268,
00261     0.491427,
00262     0.627979,
00263     0.75918,
00264     0.873185,
00265     0.957934,
00266     0.999381,
00267     1,
00268     0.957738
00269   };
00270 
00271   static double weight_22[25] = {
00272     0,
00273     0,
00274     0,
00275     5.88636e-06,
00276     3.0609e-05,
00277     0.000143627,
00278     0.000561558,
00279     0.00173059,
00280     0.00460078,
00281     0.0110616,
00282     0.0238974,
00283     0.0475406,
00284     0.0875077,
00285     0.148682,
00286     0.235752,
00287     0.343591,
00288     0.473146,
00289     0.611897,
00290     0.748345,
00291     0.865978,
00292     0.953199,
00293     0.997848,
00294     1,
00295     0.954245,
00296     0.873688
00297   };
00298 
00299   static double weight_21[25] = {
00300     0,
00301     0,
00302     1.15381e-06,
00303     8.07665e-06,
00304     7.1536e-05,
00305     0.000280375,
00306     0.00107189,
00307     0.00327104,
00308     0.00809396,
00309     0.0190978,
00310     0.0401894,
00311     0.0761028,
00312     0.13472,
00313     0.216315,
00314     0.324649,
00315     0.455125,
00316     0.598241,
00317     0.739215,
00318     0.861866,
00319     0.953911,
00320     0.998918,
00321     1,
00322     0.956683,
00323     0.872272,
00324     0.76399
00325   };
00326  
00327  
00328   static double weight_20[25] = {
00329     0,
00330     0,
00331     1.12532e-06,
00332     2.58822e-05,
00333     0.000145166,
00334     0.000633552,
00335     0.00215048,
00336     0.00592816,
00337     0.0145605,
00338     0.0328367,
00339     0.0652649,
00340     0.11893,
00341     0.19803,
00342     0.305525,
00343     0.436588,
00344     0.581566,
00345     0.727048,
00346     0.8534,
00347     0.949419,
00348     0.999785,
00349     1,
00350     0.953008,
00351     0.865689,
00352     0.753288,
00353     0.62765
00354   }; 
00355   static double weight_19[25] = {
00356     0,
00357     0,
00358     1.20714e-05,
00359     5.92596e-05,
00360     0.000364337,
00361     0.00124994,
00362     0.00403953,
00363     0.0108149,
00364     0.025824,
00365     0.0544969,
00366     0.103567,
00367     0.17936,
00368     0.283532,
00369     0.416091,
00370     0.562078,
00371     0.714714,
00372     0.846523,
00373     0.947875,
00374     1,
00375     0.999448,
00376     0.951404,
00377     0.859717,
00378     0.742319,
00379     0.613601,
00380     0.48552
00381   };
00382 
00383   static double weight_18[25] = {
00384     0,
00385     3.20101e-06,
00386     2.88091e-05,
00387     0.000164319,
00388     0.000719161,
00389     0.00250106,
00390     0.00773685,
00391     0.0197513,
00392     0.0443693,
00393     0.0885998,
00394     0.159891,
00395     0.262607,
00396     0.392327,
00397     0.543125,
00398     0.69924,
00399     0.837474,
00400     0.943486,
00401     0.998029,
00402     1,
00403     0.945937,
00404     0.851807,
00405     0.729309,
00406     0.596332,
00407     0.467818,
00408     0.350434
00409   };
00410 
00411  
00412   static double weight_17[25] = {
00413     1.03634e-06,
00414     7.25437e-06,
00415     4.97443e-05,
00416     0.000340956,
00417     0.00148715,
00418     0.00501485,
00419     0.0143067,
00420     0.034679,
00421     0.0742009,
00422     0.140287,
00423     0.238288,
00424     0.369416,
00425     0.521637,
00426     0.682368,
00427     0.828634,
00428     0.939655,
00429     1,
00430     0.996829,
00431     0.94062,
00432     0.841575,
00433     0.716664,
00434     0.582053,
00435     0.449595,
00436     0.331336,
00437     0.234332
00438   };
00439 
00440  
00441   static double weight_16[25] = {
00442     4.03159e-06,
00443     2.41895e-05,
00444     0.000141106,
00445     0.00081942,
00446     0.00314565,
00447     0.00990662,
00448     0.026293,
00449     0.0603881,
00450     0.120973,
00451     0.214532,
00452     0.343708,
00453     0.501141,
00454     0.665978,
00455     0.820107,
00456     0.938149,
00457     1,
00458     0.99941,
00459     0.940768,
00460     0.837813,
00461     0.703086,
00462     0.564023,
00463     0.42928,
00464     0.312515,
00465     0.216251,
00466     0.14561
00467   };
00468  
00469  
00470   static double weight_15[25] = {
00471     9.76084e-07,
00472     5.07564e-05,
00473     0.000303562,
00474     0.00174036,
00475     0.00617959,
00476     0.0188579,
00477     0.047465,
00478     0.101656,
00479     0.189492,
00480     0.315673,
00481     0.474383,
00482     0.646828,
00483     0.809462,
00484     0.934107,
00485     0.998874,
00486     1,
00487     0.936163,
00488     0.827473,
00489     0.689675,
00490     0.544384,
00491     0.40907,
00492     0.290648,
00493     0.198861,
00494     0.12951,
00495     0.0808051
00496   };
00497  
00498  
00499   static double weight_14[25] = {
00500     1.13288e-05,
00501     0.000124617,
00502     0.000753365,
00503     0.00345056,
00504     0.0123909,
00505     0.0352712,
00506     0.0825463,
00507     0.16413,
00508     0.287213,
00509     0.44615,
00510     0.625826,
00511     0.796365,
00512     0.930624,
00513     0.999958,
00514     1,
00515     0.934414,
00516     0.816456,
00517     0.672939,
00518     0.523033,
00519     0.386068,
00520     0.269824,
00521     0.180342,
00522     0.114669,
00523     0.0698288,
00524     0.0406496
00525   };
00526 
00527  
00528   static double weight_13[25] = {
00529     2.54296e-05,
00530     0.000261561,
00531     0.00167018,
00532     0.00748083,
00533     0.0241308,
00534     0.0636801,
00535     0.138222,
00536     0.255814,
00537     0.414275,
00538     0.600244,
00539     0.779958,
00540     0.92256,
00541     0.999155,
00542     1,
00543     0.927126,
00544     0.804504,
00545     0.651803,
00546     0.497534,
00547     0.35976,
00548     0.245834,
00549     0.160904,
00550     0.0991589,
00551     0.0585434,
00552     0.0332437,
00553     0.0180159
00554   };
00555 
00556   static double weight_12[25] = {
00557     5.85742e-05,
00558     0.000627706,
00559     0.00386677,
00560     0.0154068,
00561     0.0465892,
00562     0.111683,
00563     0.222487,
00564     0.381677,
00565     0.5719,
00566     0.765001,
00567     0.915916,
00568     1,
00569     0.999717,
00570     0.921443,
00571     0.791958,
00572     0.632344,
00573     0.475195,
00574     0.334982,
00575     0.223666,
00576     0.141781,
00577     0.0851538,
00578     0.048433,
00579     0.0263287,
00580     0.0133969,
00581     0.00696683
00582   };
00583 
00584  
00585   static double weight_11[25] = {
00586     0.00015238,
00587     0.00156064,
00588     0.00846044,
00589     0.0310939,
00590     0.0856225,
00591     0.187589,
00592     0.343579,
00593     0.541892,
00594     0.74224,
00595     0.909269,
00596     0.998711,
00597     1,
00598     0.916889,
00599     0.77485,
00600     0.608819,
00601     0.447016,
00602     0.307375,
00603     0.198444,
00604     0.121208,
00605     0.070222,
00606     0.0386492,
00607     0.0201108,
00608     0.0100922,
00609     0.00484937,
00610     0.00222458
00611   };
00612 
00613   static double weight_10[25] = {
00614     0.000393044,
00615     0.00367001,
00616     0.0179474,
00617     0.060389,
00618     0.151477,
00619     0.302077,
00620     0.503113,
00621     0.720373,
00622     0.899568,
00623     1,
00624     0.997739,
00625     0.909409,
00626     0.75728,
00627     0.582031,
00628     0.415322,
00629     0.277663,
00630     0.174147,
00631     0.102154,
00632     0.0566719,
00633     0.0298642,
00634     0.0147751,
00635     0.00710995,
00636     0.00319628,
00637     0.00140601,
00638     0.000568796
00639   };
00640 
00641  
00642   static double weight_9[25] = {
00643     0.00093396,
00644     0.00854448,
00645     0.0380306,
00646     0.113181,
00647     0.256614,
00648     0.460894,
00649     0.690242,
00650     0.888781,
00651     1,
00652     0.998756,
00653     0.899872,
00654     0.735642,
00655     0.552532,
00656     0.382726,
00657     0.246114,
00658     0.147497,
00659     0.0825541,
00660     0.0441199,
00661     0.0218157,
00662     0.0103578,
00663     0.00462959,
00664     0.0019142,
00665     0.000771598,
00666     0.000295893,
00667     0.000111529
00668   };
00669 
00670  
00671   static double weight_8[25] = {
00672     0.00240233,
00673     0.0192688,
00674     0.0768653,
00675     0.205008,
00676     0.410958,
00677     0.65758,
00678     0.875657,
00679     0.999886,
00680     1,
00681     0.889476,
00682     0.711446,
00683     0.517781,
00684     0.345774,
00685     0.212028,
00686     0.121208,
00687     0.0644629,
00688     0.0324928,
00689     0.0152492,
00690     0.00673527,
00691     0.0028547,
00692     0.00117213,
00693     0.000440177,
00694     0.000168471,
00695     5.80689e-05,
00696     1.93563e-05
00697   };
00698 
00699   static double weight_7[25] = {
00700     0.00617233,
00701     0.0428714,
00702     0.150018,
00703     0.350317,
00704     0.612535,
00705     0.856525,
00706     0.999923,
00707     1,
00708     0.87544,
00709     0.679383,
00710     0.478345,
00711     0.303378,
00712     0.176923,
00713     0.0950103,
00714     0.0476253,
00715     0.0222211,
00716     0.00972738,
00717     0.00392962,
00718     0.0015258,
00719     0.000559168,
00720     0.000183928,
00721     6.77983e-05,
00722     1.67818e-05,
00723     7.38398e-06,
00724     6.71271e-07
00725   };
00726  
00727   static double weight_6[25] = {
00728     0.0154465,
00729     0.0923472,
00730     0.277322,
00731     0.55552,
00732     0.833099,
00733     0.999035,
00734     1,
00735     0.855183,
00736     0.641976,
00737     0.428277,
00738     0.256804,
00739     0.139798,
00740     0.0700072,
00741     0.0321586,
00742     0.0137971,
00743     0.00544756,
00744     0.00202316,
00745     0.000766228,
00746     0.000259348,
00747     8.45836e-05,
00748     1.80362e-05,
00749     8.70713e-06,
00750     3.73163e-06,
00751     6.21938e-07,
00752     0
00753   };
00754  
00755  
00756   static double weight_5[25] = {
00757     0.0382845,
00758     0.191122,
00759     0.478782,
00760     0.797314,
00761     1,
00762     0.997148,
00763     0.831144,
00764     0.59461,
00765     0.371293,
00766     0.205903,
00767     0.103102,
00768     0.0471424,
00769     0.0194997,
00770     0.00749415,
00771     0.00273709,
00772     0.000879189,
00773     0.000286049,
00774     0.000102364,
00775     1.70606e-05,
00776     3.98081e-06,
00777     2.27475e-06,
00778     0,
00779     0,
00780     0,
00781     0
00782   };
00783  
00784  
00785   static double weight_4[25] = {
00786     0.0941305,
00787     0.373824,
00788     0.750094,
00789     1,
00790     0.997698,
00791     0.800956,
00792     0.532306,
00793     0.304597,
00794     0.152207,
00795     0.0676275,
00796     0.0270646,
00797     0.00975365,
00798     0.00326077,
00799     0.00101071,
00800     0.000301781,
00801     7.41664e-05,
00802     1.58563e-05,
00803     3.58045e-06,
00804     1.02299e-06,
00805     0,
00806     5.11493e-07,
00807     0,
00808     0,
00809     0,
00810     0
00811   };
00812  
00813  
00814   static double weight_3[25] = {
00815     0.222714,
00816     0.667015,
00817     1,
00818     0.999208,
00819     0.750609,
00820     0.449854,
00821     0.224968,
00822     0.0965185,
00823     0.0361225,
00824     0.012084,
00825     0.00359618,
00826     0.000977166,
00827     0.000239269,
00828     6.29422e-05,
00829     1.16064e-05,
00830     1.78559e-06,
00831     0,
00832     4.46398e-07,
00833     0,
00834     0,
00835     0,
00836     0,
00837     0,
00838     0,
00839     0
00840   };
00841  
00842   static double weight_2[25] = {
00843     0.499541,
00844     0.999607,
00845     1,
00846     0.666607,
00847     0.333301,
00848     0.13279,
00849     0.0441871,
00850     0.0127455,
00851     0.00318434,
00852     0.00071752,
00853     0.000132204,
00854     2.69578e-05,
00855     5.16999e-06,
00856     2.21571e-06,
00857     0,
00858     0,
00859     0,
00860     0,
00861     0,
00862     0,
00863     0,
00864     0,
00865     0,
00866     0,
00867     0
00868   };
00869  
00870   static double weight_1[25] = {
00871     0.999165,
00872     1,
00873     0.499996,
00874     0.166868,
00875     0.0414266,
00876     0.00831053,
00877     0.00137472,
00878     0.000198911,
00879     2.66302e-05,
00880     2.44563e-06,
00881     2.71737e-07,
00882     2.71737e-07,
00883     0,
00884     0,
00885     0,
00886     0,
00887     0,
00888     0,
00889     0,
00890     0,
00891     0,
00892     0,
00893     0,
00894     0,
00895     0
00896   };
00897  
00898   static double weight_0[25] = {
00899     1,
00900     0,
00901     0,
00902     0,
00903     0,
00904     0,
00905     0,
00906     0,
00907     0,
00908     0,
00909     0,
00910     0,
00911     0,
00912     0,
00913     0,
00914     0,
00915     0,
00916     0,
00917     0,
00918     0,
00919     0,
00920     0,
00921     0,
00922     0,
00923     0
00924   };
00925 
00926   //WeightOOTPU_ = {0};
00927 
00928   double* WeightPtr = 0;
00929 
00930   for(int iint = 0; iint<25; ++iint){
00931     if(iint ==0) WeightPtr = weight_0;
00932     if(iint ==1) WeightPtr = weight_1;
00933     if(iint ==2) WeightPtr = weight_2;
00934     if(iint ==3) WeightPtr = weight_3;
00935     if(iint ==4) WeightPtr = weight_4;
00936     if(iint ==5) WeightPtr = weight_5;
00937     if(iint ==6) WeightPtr = weight_6;
00938     if(iint ==7) WeightPtr = weight_7;
00939     if(iint ==8) WeightPtr = weight_8;
00940     if(iint ==9) WeightPtr = weight_9;
00941     if(iint ==10) WeightPtr = weight_10;
00942     if(iint ==11) WeightPtr = weight_11;
00943     if(iint ==12) WeightPtr = weight_12;
00944     if(iint ==13) WeightPtr = weight_13;
00945     if(iint ==14) WeightPtr = weight_14;
00946     if(iint ==15) WeightPtr = weight_15;
00947     if(iint ==16) WeightPtr = weight_16;
00948     if(iint ==17) WeightPtr = weight_17;
00949     if(iint ==18) WeightPtr = weight_18;
00950     if(iint ==19) WeightPtr = weight_19;
00951     if(iint ==20) WeightPtr = weight_20;
00952     if(iint ==21) WeightPtr = weight_21;
00953     if(iint ==22) WeightPtr = weight_22;
00954     if(iint ==23) WeightPtr = weight_23;
00955     if(iint ==24) WeightPtr = weight_24;
00956 
00957     for(int ibin = 0; ibin<25; ++ibin){
00958       WeightOOTPU_[iint][ibin] = *(WeightPtr+ibin);
00959     }
00960   }
00961 
00962 }
00963 
00964 // Use this routine to re-weight out-of-time pileup to match the in-time distribution
00965 // As of May 2011, CMS is only sensitive to a bunch that is 50ns "late", which corresponds to
00966 // BunchCrossing +1.  So, we use that here for re-weighting.
00967 
00968 double LumiReWeighting::weightOOT( const edm::EventBase &e ) {
00969 
00970 
00971   static double Correct_Weights2011[25] = { // residual correction to match lumi spectrum
00972     5.30031,
00973     2.07903,
00974     1.40729,
00975     1.27687,
00976     1.0702,
00977     0.902094,
00978     0.902345,
00979     0.931449,
00980     0.78202,
00981     0.824686,
00982     0.837735,
00983     0.910261,
00984     1.01394,
00985     1.1599,
00986     1.12778,
00987     1.58423,
00988     1.78868,
00989     1.58296,
00990     2.3291,
00991     3.86641,
00992     0,
00993     0,
00994     0,
00995     0,
00996     0
00997   };                        
00998 
00999   //int Run = e.run();
01000   int LumiSection = e.luminosityBlock();
01001   
01002   // do some caching here, attempt to catch file boundaries
01003 
01004   if(LumiSection != OldLumiSection_) {
01005 
01006     Reweight_4_2_2p2_ = false;
01007 
01008     static edm::ReleaseVersion TargetRelease = edm::ReleaseVersion("\"CMSSW_4_2_2_patch2\"");
01009     edm::ProcessHistory PHist = e.processHistory();
01010     edm::ProcessHistory::const_iterator PHist_iter = PHist.begin();
01011 
01012     // check to see if we need to do the special out-of-time poisson weighting for CMSSW_4_2_2_patch2 MC:
01013 
01014     for(; PHist_iter<PHist.end() ;++PHist_iter) {
01015       edm::ProcessConfiguration PConf = *(PHist_iter);
01016       edm::ReleaseVersion Release =  PConf.releaseVersion() ;
01017       const std::string Process =  PConf.processName();
01018 
01019       if((Process=="HLT") && (Release==TargetRelease)) {
01020 
01021         Reweight_4_2_2p2_ = true;
01022 
01023         if(FirstWarning_) {
01024 
01025           std::cout << " **** Warning: Out-of-time pileup reweighting appropriate for Release " << Release << " **** " << std::endl;
01026           std::cout << " **** will be applied  ****  " << std::endl;
01027 
01028           FirstWarning_ = false;
01029         }
01030       }
01031     }
01032     OldLumiSection_ = LumiSection;
01033   }
01034 
01035   // find the pileup summary information
01036 
01037   Handle<std::vector< PileupSummaryInfo > >  PupInfo;
01038   e.getByLabel(edm::InputTag("addPileupInfo"), PupInfo);
01039 
01040   std::vector<PileupSummaryInfo>::const_iterator PVI;
01041 
01042   int npv = -1;
01043   int npv50ns = -1;
01044 
01045   for(PVI = PupInfo->begin(); PVI != PupInfo->end(); ++PVI) {
01046 
01047     int BX = PVI->getBunchCrossing();
01048 
01049     if(BX == 0) { 
01050       npv = PVI->getPU_NumInteractions();
01051     }
01052 
01053     if(BX == 1) { 
01054       npv50ns = PVI->getPU_NumInteractions();
01055     }
01056 
01057   }
01058 
01059   // Note: for the "uncorrelated" out-of-time pileup, reweighting is only done on the 50ns
01060   // "late" bunch (BX=+1), since that is basically the only one that matters in terms of 
01061   // energy deposition.  
01062 
01063   if(npv < 0) {
01064     std::cerr << " no in-time beam crossing found\n! " ;
01065     std::cerr << " Returning event weight=0\n! ";
01066     return 0.;
01067   }
01068   if(npv50ns < 0) {
01069     std::cerr << " no out-of-time beam crossing found\n! " ;
01070     std::cerr << " Returning event weight=0\n! ";
01071     return 0.;
01072   }
01073 
01074   int bin = weights_->GetXaxis()->FindBin( npv );
01075 
01076   double inTimeWeight = weights_->GetBinContent( bin );
01077 
01078   double TotalWeight = 1.0;
01079 
01080   if(Reweight_4_2_2p2_) {
01081     TotalWeight = inTimeWeight * WeightOOTPU_[bin-1][npv50ns] * Correct_Weights2011[bin-1];
01082   }
01083   else {
01084     TotalWeight = inTimeWeight;
01085   }
01086 
01087   return TotalWeight;
01088  
01089 }
01090 
01091 #endif