CMS 3D CMS Logo

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