CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LumiReWeighting.cc
Go to the documentation of this file.
1 #ifndef PhysicsTools_Utilities_interface_LumiReWeighting_cc
2 #define PhysicsTools_Utilities_interface_LumiReWeighting_cc
3 
4 
20 #include "TRandom1.h"
21 #include "TRandom2.h"
22 #include "TRandom3.h"
23 #include "TStopwatch.h"
24 #include "TH1.h"
25 #include "TFile.h"
26 #include <string>
27 #include <algorithm>
28 #include <boost/shared_ptr.hpp>
35 
36 using namespace edm;
37 
39  std::string dataFile,
40  std::string GenHistName = "pileup",
41  std::string DataHistName = "pileup" ) :
42  generatedFileName_( generatedFile),
43  dataFileName_ ( dataFile ),
44  GenHistName_ ( GenHistName ),
45  DataHistName_ ( DataHistName )
46  {
47  generatedFile_ = boost::shared_ptr<TFile>( new TFile(generatedFileName_.c_str()) ); //MC distribution
48  dataFile_ = boost::shared_ptr<TFile>( new TFile(dataFileName_.c_str()) ); //Data distribution
49 
50  Data_distr_ = boost::shared_ptr<TH1>( (static_cast<TH1*>(dataFile_->Get( DataHistName_.c_str() )->Clone() )) );
51  MC_distr_ = boost::shared_ptr<TH1>( (static_cast<TH1*>(generatedFile_->Get( GenHistName_.c_str() )->Clone() )) );
52 
53  // MC * data/MC = data, so the weights are data/MC:
54 
55  // normalize both histograms first
56 
57  Data_distr_->Scale( 1.0/ Data_distr_->Integral() );
58  MC_distr_->Scale( 1.0/ MC_distr_->Integral() );
59 
60  weights_ = boost::shared_ptr<TH1>( static_cast<TH1*>(Data_distr_->Clone()) );
61 
62  weights_->SetName("lumiWeights");
63 
64  TH1* den = dynamic_cast<TH1*>(MC_distr_->Clone());
65 
66  //den->Scale(1.0/ den->Integral());
67 
68  weights_->Divide( den ); // so now the average weight should be 1.0
69 
70  std::cout << " Lumi/Pileup Reweighting: Computed Weights per In-Time Nint " << std::endl;
71 
72  int NBins = weights_->GetNbinsX();
73 
74  for(int ibin = 1; ibin<NBins+1; ++ibin){
75  std::cout << " " << ibin-1 << " " << weights_->GetBinContent(ibin) << std::endl;
76  }
77 
78 
79  FirstWarning_ = true;
80  OldLumiSection_ = -1;
81 }
82 
83 LumiReWeighting::LumiReWeighting( std::vector< float > MC_distr, std::vector< float > Lumi_distr) {
84  // no histograms for input: use vectors
85 
86  // now, make histograms out of them:
87 
88  // first, check they are the same size...
89 
90  if( MC_distr.size() != Lumi_distr.size() ){
91 
92  std::cerr <<"ERROR: LumiReWeighting: input vectors have different sizes. Quitting... \n";
93  return;
94 
95  }
96 
97  Int_t NBins = MC_distr.size();
98 
99  MC_distr_ = boost::shared_ptr<TH1> ( new TH1F("MC_distr","MC dist",NBins,-0.5, float(NBins)-0.5) );
100  Data_distr_ = boost::shared_ptr<TH1> ( new TH1F("Data_distr","Data dist",NBins,-0.5, float(NBins)-0.5) );
101 
102  weights_ = boost::shared_ptr<TH1> ( new TH1F("luminumer","luminumer",NBins,-0.5, float(NBins)-0.5) );
103  TH1* den = new TH1F("lumidenom","lumidenom",NBins,-0.5, float(NBins)-0.5) ;
104 
105  for(int ibin = 1; ibin<NBins+1; ++ibin ) {
106  weights_->SetBinContent(ibin, Lumi_distr[ibin-1]);
107  Data_distr_->SetBinContent(ibin, Lumi_distr[ibin-1]);
108  den->SetBinContent(ibin,MC_distr[ibin-1]);
109  MC_distr_->SetBinContent(ibin,MC_distr[ibin-1]);
110  }
111 
112  // check integrals, make sure things are normalized
113 
114  float deltaH = weights_->Integral();
115  if(fabs(1.0 - deltaH) > 0.02 ) { //*OOPS*...
116  weights_->Scale( 1.0/ weights_->Integral() );
117  Data_distr_->Scale( 1.0/ Data_distr_->Integral() );
118  }
119  float deltaMC = den->Integral();
120  if(fabs(1.0 - deltaMC) > 0.02 ) {
121  den->Scale(1.0/ den->Integral());
122  MC_distr_->Scale(1.0/ MC_distr_->Integral());
123  }
124 
125  weights_->Divide( den ); // so now the average weight should be 1.0
126 
127  std::cout << " Lumi/Pileup Reweighting: Computed Weights per In-Time Nint " << std::endl;
128 
129  for(int ibin = 1; ibin<NBins+1; ++ibin){
130  std::cout << " " << ibin-1 << " " << weights_->GetBinContent(ibin) << std::endl;
131  }
132 
133  FirstWarning_ = true;
134  OldLumiSection_ = -1;
135 }
136 
137 double LumiReWeighting::weight( int npv ) {
138  int bin = weights_->GetXaxis()->FindBin( npv );
139  return weights_->GetBinContent( bin );
140 }
141 
142 double LumiReWeighting::weight( float npv ) {
143  int bin = weights_->GetXaxis()->FindBin( npv );
144  return weights_->GetBinContent( bin );
145 }
146 
147 
148 
149 // This version of weight does all of the work for you, assuming you want to re-weight
150 // using the true number of interactions in the in-time beam crossing.
151 
152 
154 
155  // find provenance of event objects, just to check at the job beginning if there might be an issue
156 
157  if(FirstWarning_) {
158 
160  edm::ProcessHistory::const_iterator PHist_iter = PHist.begin();
161 
162  for(; PHist_iter<PHist.end() ;++PHist_iter) {
163  edm::ProcessConfiguration PConf = *(PHist_iter);
164  edm::ReleaseVersion Release = PConf.releaseVersion() ;
165  const std::string Process = PConf.processName();
166 
167  }
168  // SetFirstFalse();
169  FirstWarning_ = false;
170  }
171 
172  // get pileup summary information
173 
175  e.getByLabel(edm::InputTag("addPileupInfo"), PupInfo);
176 
177  std::vector<PileupSummaryInfo>::const_iterator PVI;
178 
179  int npv = -1;
180  for(PVI = PupInfo->begin(); PVI != PupInfo->end(); ++PVI) {
181 
182  int BX = PVI->getBunchCrossing();
183 
184  if(BX == 0) {
185  npv = PVI->getPU_NumInteractions();
186  continue;
187  }
188 
189  }
190 
191  if(npv < 0) std::cerr << " no in-time beam crossing found\n! " ;
192 
193  return weight(npv);
194 
195 }
196 
197 // Use this routine to re-weight out-of-time pileup to match the in-time distribution
198 // As of May 2011, CMS is only sensitive to a bunch that is 50ns "late", which corresponds to
199 // BunchCrossing +1. So, we use that here for re-weighting.
200 
202 
203 
204  //int Run = e.run();
205  int LumiSection = e.luminosityBlock();
206 
207  // do some caching here, attempt to catch file boundaries
208 
209  if(LumiSection != OldLumiSection_) {
210 
212  edm::ProcessHistory::const_iterator PHist_iter = PHist.begin();
213 
214  for(; PHist_iter<PHist.end() ;++PHist_iter) {
215  edm::ProcessConfiguration PConf = *(PHist_iter);
216  edm::ReleaseVersion Release = PConf.releaseVersion() ;
217  const std::string Process = PConf.processName();
218 
219  }
220  OldLumiSection_ = LumiSection;
221  }
222 
223  // find the pileup summary information
224 
226  e.getByLabel(edm::InputTag("addPileupInfo"), PupInfo);
227 
228  std::vector<PileupSummaryInfo>::const_iterator PVI;
229 
230  int npv = -1;
231  int npv50ns = -1;
232 
233  for(PVI = PupInfo->begin(); PVI != PupInfo->end(); ++PVI) {
234 
235  int BX = PVI->getBunchCrossing();
236 
237  if(BX == 0) {
238  npv = PVI->getPU_NumInteractions();
239  }
240 
241  if(BX == 1) {
242  npv50ns = PVI->getPU_NumInteractions();
243  }
244 
245  }
246 
247  // Note: for the "uncorrelated" out-of-time pileup, reweighting is only done on the 50ns
248  // "late" bunch (BX=+1), since that is basically the only one that matters in terms of
249  // energy deposition.
250 
251  if(npv < 0) {
252  std::cerr << " no in-time beam crossing found\n! " ;
253  std::cerr << " Returning event weight=0\n! ";
254  return 0.;
255  }
256  if(npv50ns < 0) {
257  std::cerr << " no out-of-time beam crossing found\n! " ;
258  std::cerr << " Returning event weight=0\n! ";
259  return 0.;
260  }
261 
262  int bin = weights_->GetXaxis()->FindBin( npv );
263 
264  double inTimeWeight = weights_->GetBinContent( bin );
265 
266  double TotalWeight = 1.0;
267 
268  TotalWeight = inTimeWeight;
269 
270  return TotalWeight;
271 
272 }
273 
274 #endif
collection_type::const_iterator const_iterator
const_iterator begin() const
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:59
boost::shared_ptr< TH1 > MC_distr_
std::string const & processName() const
double weight(int npv)
double weightOOT(const edm::EventBase &e)
boost::shared_ptr< TFile > dataFile_
boost::shared_ptr< TFile > generatedFile_
virtual ProcessHistory const & processHistory() const =0
ReleaseVersion const & releaseVersion() const
boost::shared_ptr< TH1 > Data_distr_
std::string ReleaseVersion
Definition: ReleaseVersion.h:7
const_iterator end() const
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:86
boost::shared_ptr< TH1 > weights_
tuple cout
Definition: gather_cfg.py:121
std::string generatedFileName_