CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LumiReWeighting.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_Utilities_interface_LumiReWeighting_h
2 #define PhysicsTools_Utilities_interface_LumiReWeighting_h
3 
4 
18 #include "TH1.h"
19 #include "TFile.h"
20 #include <cmath>
21 #include <string>
22 #include <boost/shared_ptr.hpp>
23 #include <vector>
24 
25 namespace reweight{
26 
27  // add a class to shift the mean of a poisson-like luminosity distribution by an arbitrary amount.
28  // Only valid for small (<1.5) shifts about the 2011 lumi distribution for now, because shifts are non-linear
29  // Each PoissonMeanShifter does one shift, so defining multiples can give you an arbitrary collection
30 
32 
33  public:
34 
36 
37  PoissonMeanShifter( float Shift ){
38 
39  // these are the polynomial or exponential coefficients for each bin of a 25-bin sequence that
40  // convert the Distribution of the 2011 luminosity to something with a lower or higher peak luminosity.
41  // The distributions aren't quite poisson because they model luminosity decreasing during a fill. This implies that
42  // they do get wider as the mean increases, so the weights are not linear with increasing mean.
43 
44  constexpr double p0_minus[20] = { 1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1. };
45  constexpr double p1_minus[20] = {
46  -0.677786,
47  -0.619614,
48  -0.49465,
49  -0.357963,
50  -0.238359,
51  -0.110002,
52  0.0348629,
53  0.191263,
54  0.347648,
55  0.516615,
56  0.679646,
57  0.836673,
58  0.97764,
59  1.135,
60  1.29922,
61  1.42467,
62  1.55901,
63  1.61762,
64  1.67275,
65  1.96008
66  };
67  constexpr double p2_minus[20] = {
68  0.526164,
69  0.251816,
70  0.11049,
71  0.026917,
72  -0.0464692,
73  -0.087022,
74  -0.0931581,
75  -0.0714295,
76  -0.0331772,
77  0.0347473,
78  0.108658,
79  0.193048,
80  0.272314,
81  0.376357,
82  0.4964,
83  0.58854,
84  0.684959,
85  0.731063,
86  0.760044,
87  1.02386
88  };
89 
90  constexpr double p1_expoM[5] = {
91  1.63363e-03,
92  6.79290e-04,
93  3.69900e-04,
94  2.24349e-04,
95  9.87156e-06
96  };
97 
98  constexpr double p2_expoM[5] = {
99  2.64692,
100  3.26585,
101  3.53229,
102  4.18035,
103  5.64027
104  };
105 
106 
107  constexpr double p0_plus[20] = { 1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1.,1. };
108  constexpr double p1_plus[20] = {
109  -0.739059,
110  -0.594445,
111  -0.477276,
112  -0.359707,
113  -0.233573,
114  -0.103458,
115  0.0373401,
116  0.176571,
117  0.337617,
118  0.499074,
119  0.675126,
120  0.840522,
121  1.00917,
122  1.15847,
123  1.23816,
124  1.44271,
125  1.52982,
126  1.46385,
127  1.5802,
128  0.988689
129  };
130  constexpr double p2_plus[20] = {
131  0.208068,
132  0.130033,
133  0.0850356,
134  0.0448344,
135  0.000749832,
136  -0.0331347,
137  -0.0653281,
138  -0.0746009,
139  -0.0800667,
140  -0.0527636,
141  -0.00402649,
142  0.103338,
143  0.261261,
144  0.491084,
145  0.857966,
146  1.19495,
147  1.75071,
148  2.65559,
149  3.35433,
150  5.48835
151  };
152 
153  constexpr double p1_expoP[5] = {
154  1.42463e-01,
155  4.18966e-02,
156  1.12697e-01,
157  1.66197e-01,
158  1.50768e-01
159  };
160 
161  constexpr double p2_expoP[5] = {
162  1.98758,
163  2.27217,
164  2.26799,
165  2.38455,
166  2.52428
167  };
168 
169  // initialize weights based on desired Shift
170 
171 
172 
173  for (int ibin=0;ibin<20;ibin++) {
174 
175  if( Shift < .0) {
176  Pweight_[ibin] = p0_minus[ibin] + p1_minus[ibin]*Shift + p2_minus[ibin]*Shift*Shift;
177  }
178  else {
179  Pweight_[ibin] = p0_plus[ibin] + p1_plus[ibin]*Shift + p2_plus[ibin]*Shift*Shift;
180  }
181  }
182 
183  // last few bins fit better to an exponential...
184 
185  for (int ibin=20;ibin<25;ibin++) {
186  if( Shift < 0.) {
187  Pweight_[ibin] = p1_expoM[ibin-20]*exp(p2_expoM[ibin-20]*Shift);
188  }
189  else {
190  Pweight_[ibin] = p1_expoP[ibin-20]*exp(p2_expoP[ibin-20]*Shift);
191  }
192  }
193 
194  };
195 
196  double ShiftWeight( int ibin ) {
197  if(ibin<25 && ibin>=0) { return Pweight_[ibin]; }
198  else { return 0;}
199 
200  };
201 
202  double ShiftWeight( float pvnum ) {
203  return ShiftWeight(int(pvnum));
204  };
205 
206  private:
207 
208  double Pweight_[25];
209 
210  };
211 
212 
213 
214 } // end of reweight namespace
215 
216 
217 namespace edm {
218  class EventBase;
220  public:
221  LumiReWeighting( std::string generatedFile,
222  std::string dataFile,
223  std::string GenHistName,
224  std::string DataHistName);
225 
226  LumiReWeighting( const std::vector< float >& MC_distr, const std::vector< float >& Lumi_distr);
227 
228  LumiReWeighting ( ) { } ;
229 
230  double weight( int npv ) ;
231 
232  double weight( float npv ) ;
233 
234  double weight( const edm::EventBase &e ) ;
235 
236  double weightOOT( const edm::EventBase &e ) ;
237 
238  protected:
239 
244  boost::shared_ptr<TFile> generatedFile_;
245  boost::shared_ptr<TFile> dataFile_;
246  boost::shared_ptr<TH1> weights_;
247 
248  //keep copies of normalized distributions:
249 
250  boost::shared_ptr<TH1> MC_distr_;
251  boost::shared_ptr<TH1> Data_distr_;
252 
255 
256 
257  };
258 }
259 
260 
261 
262 #endif
boost::shared_ptr< TH1 > MC_distr_
double weight(int npv)
double weightOOT(const edm::EventBase &e)
double ShiftWeight(float pvnum)
boost::shared_ptr< TFile > dataFile_
boost::shared_ptr< TFile > generatedFile_
boost::shared_ptr< TH1 > Data_distr_
boost::shared_ptr< TH1 > weights_
std::string generatedFileName_
#define constexpr