CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/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         FirstWarning_ = true;
00080         OldLumiSection_ = -1;
00081 }
00082 
00083 LumiReWeighting::LumiReWeighting( std::vector< float > MC_distr, std::vector< float > Lumi_distr) {
00084   // no histograms for input: use vectors
00085   
00086   // now, make histograms out of them:
00087 
00088   // first, check they are the same size...
00089 
00090   if( MC_distr.size() != Lumi_distr.size() ){   
00091 
00092     std::cerr <<"ERROR: LumiReWeighting: input vectors have different sizes. Quitting... \n";
00093     return;
00094 
00095   }
00096 
00097   Int_t NBins = MC_distr.size();
00098 
00099   MC_distr_ = boost::shared_ptr<TH1> ( new TH1F("MC_distr","MC dist",NBins,-0.5, float(NBins)-0.5) );
00100   Data_distr_ = boost::shared_ptr<TH1> ( new TH1F("Data_distr","Data dist",NBins,-0.5, float(NBins)-0.5) );
00101 
00102   weights_ = boost::shared_ptr<TH1> ( new TH1F("luminumer","luminumer",NBins,-0.5, float(NBins)-0.5) );
00103   TH1* den = new TH1F("lumidenom","lumidenom",NBins,-0.5, float(NBins)-0.5) ;
00104 
00105   for(int ibin = 1; ibin<NBins+1; ++ibin ) {
00106     weights_->SetBinContent(ibin, Lumi_distr[ibin-1]);
00107     Data_distr_->SetBinContent(ibin, Lumi_distr[ibin-1]);
00108     den->SetBinContent(ibin,MC_distr[ibin-1]);
00109     MC_distr_->SetBinContent(ibin,MC_distr[ibin-1]);
00110   }
00111 
00112   // check integrals, make sure things are normalized
00113 
00114   float deltaH = weights_->Integral();
00115   if(fabs(1.0 - deltaH) > 0.02 ) { //*OOPS*...
00116     weights_->Scale( 1.0/ weights_->Integral() );
00117     Data_distr_->Scale( 1.0/ Data_distr_->Integral() );
00118   }
00119   float deltaMC = den->Integral();
00120   if(fabs(1.0 - deltaMC) > 0.02 ) {
00121     den->Scale(1.0/ den->Integral());
00122     MC_distr_->Scale(1.0/ MC_distr_->Integral());
00123   }
00124 
00125   weights_->Divide( den );  // so now the average weight should be 1.0    
00126 
00127   std::cout << " Lumi/Pileup Reweighting: Computed Weights per In-Time Nint " << std::endl;
00128 
00129   for(int ibin = 1; ibin<NBins+1; ++ibin){
00130     std::cout << "   " << ibin-1 << " " << weights_->GetBinContent(ibin) << std::endl;
00131   }
00132 
00133   FirstWarning_ = true;
00134   OldLumiSection_ = -1;
00135 }
00136 
00137 double LumiReWeighting::weight( int npv ) {
00138   int bin = weights_->GetXaxis()->FindBin( npv );
00139   return weights_->GetBinContent( bin );
00140 }
00141 
00142 double LumiReWeighting::weight( float npv ) {
00143   int bin = weights_->GetXaxis()->FindBin( npv );
00144   return weights_->GetBinContent( bin );
00145 }
00146 
00147 
00148 
00149 // This version of weight does all of the work for you, assuming you want to re-weight
00150 // using the true number of interactions in the in-time beam crossing.
00151 
00152 
00153 double LumiReWeighting::weight( const edm::EventBase &e ) {
00154 
00155   // find provenance of event objects, just to check at the job beginning if there might be an issue  
00156 
00157   if(FirstWarning_) {
00158 
00159     edm::ProcessHistory PHist = e.processHistory();
00160     edm::ProcessHistory::const_iterator PHist_iter = PHist.begin();
00161 
00162     for(; PHist_iter<PHist.end() ;++PHist_iter) {
00163       edm::ProcessConfiguration PConf = *(PHist_iter);
00164       edm::ReleaseVersion Release =  PConf.releaseVersion() ;
00165       const std::string Process =  PConf.processName();
00166 
00167     }
00168     //    SetFirstFalse();
00169     FirstWarning_ = false;
00170   }
00171 
00172   // get pileup summary information
00173 
00174   Handle<std::vector< PileupSummaryInfo > >  PupInfo;
00175   e.getByLabel(edm::InputTag("addPileupInfo"), PupInfo);
00176 
00177   std::vector<PileupSummaryInfo>::const_iterator PVI;
00178 
00179   int npv = -1;
00180   for(PVI = PupInfo->begin(); PVI != PupInfo->end(); ++PVI) {
00181 
00182     int BX = PVI->getBunchCrossing();
00183 
00184     if(BX == 0) { 
00185       npv = PVI->getPU_NumInteractions();
00186       continue;
00187     }
00188 
00189   }
00190 
00191   if(npv < 0) std::cerr << " no in-time beam crossing found\n! " ;
00192 
00193   return weight(npv);
00194 
00195 }
00196 
00197 // Use this routine to re-weight out-of-time pileup to match the in-time distribution
00198 // As of May 2011, CMS is only sensitive to a bunch that is 50ns "late", which corresponds to
00199 // BunchCrossing +1.  So, we use that here for re-weighting.
00200 
00201 double LumiReWeighting::weightOOT( const edm::EventBase &e ) {
00202 
00203 
00204   //int Run = e.run();
00205   int LumiSection = e.luminosityBlock();
00206   
00207   // do some caching here, attempt to catch file boundaries
00208 
00209   if(LumiSection != OldLumiSection_) {
00210 
00211     edm::ProcessHistory PHist = e.processHistory();
00212     edm::ProcessHistory::const_iterator PHist_iter = PHist.begin();
00213 
00214     for(; PHist_iter<PHist.end() ;++PHist_iter) {
00215       edm::ProcessConfiguration PConf = *(PHist_iter);
00216       edm::ReleaseVersion Release =  PConf.releaseVersion() ;
00217       const std::string Process =  PConf.processName();
00218 
00219     }
00220     OldLumiSection_ = LumiSection;
00221   }
00222 
00223   // find the pileup summary information
00224 
00225   Handle<std::vector< PileupSummaryInfo > >  PupInfo;
00226   e.getByLabel(edm::InputTag("addPileupInfo"), PupInfo);
00227 
00228   std::vector<PileupSummaryInfo>::const_iterator PVI;
00229 
00230   int npv = -1;
00231   int npv50ns = -1;
00232 
00233   for(PVI = PupInfo->begin(); PVI != PupInfo->end(); ++PVI) {
00234 
00235     int BX = PVI->getBunchCrossing();
00236 
00237     if(BX == 0) { 
00238       npv = PVI->getPU_NumInteractions();
00239     }
00240 
00241     if(BX == 1) { 
00242       npv50ns = PVI->getPU_NumInteractions();
00243     }
00244 
00245   }
00246 
00247   // Note: for the "uncorrelated" out-of-time pileup, reweighting is only done on the 50ns
00248   // "late" bunch (BX=+1), since that is basically the only one that matters in terms of 
00249   // energy deposition.  
00250 
00251   if(npv < 0) {
00252     std::cerr << " no in-time beam crossing found\n! " ;
00253     std::cerr << " Returning event weight=0\n! ";
00254     return 0.;
00255   }
00256   if(npv50ns < 0) {
00257     std::cerr << " no out-of-time beam crossing found\n! " ;
00258     std::cerr << " Returning event weight=0\n! ";
00259     return 0.;
00260   }
00261 
00262   int bin = weights_->GetXaxis()->FindBin( npv );
00263 
00264   double inTimeWeight = weights_->GetBinContent( bin );
00265 
00266   double TotalWeight = 1.0;
00267 
00268   TotalWeight = inTimeWeight;
00269 
00270   return TotalWeight;
00271  
00272 }
00273 
00274 #endif