CMS 3D CMS Logo

PFEnergyCalibration.cc
Go to the documentation of this file.
1 
4 
6 
7 #include <TMath.h>
8 #include <cmath>
9 #include <vector>
10 #include <TF1.h>
11 #include <map>
12 #include <algorithm>
13 #include <numeric>
14 
15 using namespace std;
16 
17 PFEnergyCalibration::PFEnergyCalibration() : pfCalibrations(nullptr), esEEInterCalib_(nullptr) {
19 }
20 
22 
24  threshE = 3.5;
25  threshH = 2.5;
26 
27  //calibChrisClean.C calibration parameters bhumika Nov, 2018
28  faBarrel = std::make_unique<TF1>(
29  "faBarrel", "[0]+((([1]+([2]/sqrt(x)))*exp(-(x^[6]/[3])))-([4]*exp(-(x^[7]/[5]))))", 0., 1000.);
30  faBarrel->SetParameter(0, -30.7141);
31  faBarrel->SetParameter(1, 31.7583);
32  faBarrel->SetParameter(2, 4.40594);
33  faBarrel->SetParameter(3, 1.70914);
34  faBarrel->SetParameter(4, 0.0613696);
35  faBarrel->SetParameter(5, 0.000104857);
36  faBarrel->SetParameter(6, -1.38927);
37  faBarrel->SetParameter(7, -0.743082);
38  fbBarrel = std::make_unique<TF1>(
39  "fbBarrel", "[0]+((([1]+([2]/sqrt(x)))*exp(-(x^[6]/[3])))-([4]*exp(-(x^[7]/[5]))))", 0., 1000.);
40  fbBarrel->SetParameter(0, 2.25366);
41  fbBarrel->SetParameter(1, 0.537715);
42  fbBarrel->SetParameter(2, -4.81375);
43  fbBarrel->SetParameter(3, 12.109);
44  fbBarrel->SetParameter(4, 1.80577);
45  fbBarrel->SetParameter(5, 0.187919);
46  fbBarrel->SetParameter(6, -6.26234);
47  fbBarrel->SetParameter(7, -0.607392);
48  fcBarrel = std::make_unique<TF1>(
49  "fcBarrel", "[0]+((([1]+([2]/sqrt(x)))*exp(-(x^[6]/[3])))-([4]*exp(-(x^[7]/[5]))))", 0., 1000.);
50  fcBarrel->SetParameter(0, 1.5125962);
51  fcBarrel->SetParameter(1, 0.855057);
52  fcBarrel->SetParameter(2, -6.04199);
53  fcBarrel->SetParameter(3, 2.08229);
54  fcBarrel->SetParameter(4, 0.592266);
55  fcBarrel->SetParameter(5, 0.0291232);
56  fcBarrel->SetParameter(6, 0.364802);
57  fcBarrel->SetParameter(7, -1.50142);
58  faEtaBarrelEH = std::make_unique<TF1>("faEtaBarrelEH", "[0]+[1]*exp(-x/[2])", 0., 1000.);
59  faEtaBarrelEH->SetParameter(0, 0.0185555);
60  faEtaBarrelEH->SetParameter(1, -0.0470674);
61  faEtaBarrelEH->SetParameter(2, 396.959);
62  fbEtaBarrelEH = std::make_unique<TF1>("fbEtaBarrelEH", "[0]+[1]*exp(-x/[2])", 0., 1000.);
63  fbEtaBarrelEH->SetParameter(0, 0.0396458);
64  fbEtaBarrelEH->SetParameter(1, 0.114128);
65  fbEtaBarrelEH->SetParameter(2, 251.405);
66  faEtaBarrelH = std::make_unique<TF1>("faEtaBarrelH", "[0]+[1]*x", 0., 1000.);
67  faEtaBarrelH->SetParameter(0, 0.00434994);
68  faEtaBarrelH->SetParameter(1, -5.16564e-06);
69  fbEtaBarrelH = std::make_unique<TF1>("fbEtaBarrelH", "[0]+[1]*exp(-x/[2])", 0., 1000.);
70  fbEtaBarrelH->SetParameter(0, -0.0232604);
71  fbEtaBarrelH->SetParameter(1, 0.0937525);
72  fbEtaBarrelH->SetParameter(2, 34.9935);
73 
74  faEndcap = std::make_unique<TF1>(
75  "faEndcap", "[0]+((([1]+([2]/sqrt(x)))*exp(-(x^[6]/[3])))-([4]*exp(-(x^[7]/[5]))))", 0., 1000.);
76  faEndcap->SetParameter(0, 1.17227);
77  faEndcap->SetParameter(1, 13.1489);
78  faEndcap->SetParameter(2, -29.1672);
79  faEndcap->SetParameter(3, 0.604223);
80  faEndcap->SetParameter(4, 0.0426363);
81  faEndcap->SetParameter(5, 3.30898e-15);
82  faEndcap->SetParameter(6, 0.165293);
83  faEndcap->SetParameter(7, -7.56786);
84  fbEndcap = std::make_unique<TF1>(
85  "fbEndcap", "[0]+((([1]+([2]/sqrt(x)))*exp(-(x^[6]/[3])))-([4]*exp(-(x^[7]/[5]))))", 0., 1000.);
86  fbEndcap->SetParameter(0, -0.974251);
87  fbEndcap->SetParameter(1, 1.61733);
88  fbEndcap->SetParameter(2, 0.0629183);
89  fbEndcap->SetParameter(3, 7.78495);
90  fbEndcap->SetParameter(4, -0.774289);
91  fbEndcap->SetParameter(5, 7.81399e-05);
92  fbEndcap->SetParameter(6, 0.139116);
93  fbEndcap->SetParameter(7, -4.25551);
94  fcEndcap = std::make_unique<TF1>(
95  "fcEndcap", "[0]+((([1]+([2]/sqrt(x)))*exp(-(x^[6]/[3])))-([4]*exp(-(x^[7]/[5]))))", 0., 1000.);
96  fcEndcap->SetParameter(0, 1.01863);
97  fcEndcap->SetParameter(1, 1.29787);
98  fcEndcap->SetParameter(2, -3.97293);
99  fcEndcap->SetParameter(3, 21.7805);
100  fcEndcap->SetParameter(4, 0.810195);
101  fcEndcap->SetParameter(5, 0.234134);
102  fcEndcap->SetParameter(6, 1.42226);
103  fcEndcap->SetParameter(7, -0.0997326);
104  faEtaEndcapEH = std::make_unique<TF1>("faEtaEndcapEH", "[0]+[1]*exp(-x/[2])", 0., 1000.);
105  faEtaEndcapEH->SetParameter(0, 0.0112692);
106  faEtaEndcapEH->SetParameter(1, -2.68063);
107  faEtaEndcapEH->SetParameter(2, 2.90973);
108  fbEtaEndcapEH = std::make_unique<TF1>("fbEtaEndcapEH", "[0]+[1]*exp(-x/[2])", 0., 1000.);
109  fbEtaEndcapEH->SetParameter(0, -0.0192991);
110  fbEtaEndcapEH->SetParameter(1, -0.265);
111  fbEtaEndcapEH->SetParameter(2, 80.5502);
112  faEtaEndcapH = std::make_unique<TF1>("faEtaEndcapH", "[0]+[1]*exp(-x/[2])+[3]*[3]*exp(-x*x/([4]*[4]))", 0., 1000.);
113  faEtaEndcapH->SetParameter(0, -0.0106029);
114  faEtaEndcapH->SetParameter(1, -0.692207);
115  faEtaEndcapH->SetParameter(2, 0.0542991);
116  faEtaEndcapH->SetParameter(3, -0.171435);
117  faEtaEndcapH->SetParameter(4, -61.2277);
118  fbEtaEndcapH = std::make_unique<TF1>("fbEtaEndcapH", "[0]+[1]*exp(-x/[2])+[3]*[3]*exp(-x*x/([4]*[4]))", 0., 1000.);
119  fbEtaEndcapH->SetParameter(0, 0.0214894);
120  fbEtaEndcapH->SetParameter(1, -0.266704);
121  fbEtaEndcapH->SetParameter(2, 5.2112);
122  fbEtaEndcapH->SetParameter(3, 0.303578);
123  fbEtaEndcapH->SetParameter(4, -104.367);
124 
125  //added by Bhumika on 2 august 2018
126 
127  fcEtaBarrelH = std::make_unique<TF1>("fcEtaBarrelH", "[3]*((x-[0])^[1])+[2]", 0., 1000.);
128  fcEtaBarrelH->SetParameter(0, 0);
129  fcEtaBarrelH->SetParameter(1, 2);
130  fcEtaBarrelH->SetParameter(2, 0);
131  fcEtaBarrelH->SetParameter(3, 1);
132 
133  fcEtaEndcapH = std::make_unique<TF1>("fcEtaEndcapH", "[3]*((x-[0])^[1])+[2]", 0., 1000.);
134  fcEtaEndcapH->SetParameter(0, 0);
135  fcEtaEndcapH->SetParameter(1, 0);
136  fcEtaEndcapH->SetParameter(2, 0.05);
137  fcEtaEndcapH->SetParameter(3, 0);
138 
139  fdEtaEndcapH = std::make_unique<TF1>("fdEtaEndcapH", "[3]*((x-[0])^[1])+[2]", 0., 1000.);
140  fdEtaEndcapH->SetParameter(0, 1.5);
141  fdEtaEndcapH->SetParameter(1, 4);
142  fdEtaEndcapH->SetParameter(2, -1.1);
143  fdEtaEndcapH->SetParameter(3, 1.0);
144 
145  fcEtaBarrelEH = std::make_unique<TF1>("fcEtaBarrelEH", "[3]*((x-[0])^[1])+[2]", 0., 1000.);
146  fcEtaBarrelEH->SetParameter(0, 0);
147  fcEtaBarrelEH->SetParameter(1, 2);
148  fcEtaBarrelEH->SetParameter(2, 0);
149  fcEtaBarrelEH->SetParameter(3, 1);
150 
151  fcEtaEndcapEH = std::make_unique<TF1>("fcEtaEndcapEH", "[3]*((x-[0])^[1])+[2]", 0., 1000.);
152  fcEtaEndcapEH->SetParameter(0, 0);
153  fcEtaEndcapEH->SetParameter(1, 0);
154  fcEtaEndcapEH->SetParameter(2, 0);
155  fcEtaEndcapEH->SetParameter(3, 0);
156 
157  fdEtaEndcapEH = std::make_unique<TF1>("fdEtaEndcapEH", "[3]*((x-[0])^[1])+[2]", 0., 1000.);
158  fdEtaEndcapEH->SetParameter(0, 1.5);
159  fdEtaEndcapEH->SetParameter(1, 2.0);
160  fdEtaEndcapEH->SetParameter(2, 0.6);
161  fdEtaEndcapEH->SetParameter(3, 1.0);
162 }
163 
164 void PFEnergyCalibration::energyEmHad(double t, double& e, double& h, double eta, double phi) const {
165  // Use calorimetric energy as true energy for neutral particles
166  double tt = t;
167  double ee = e;
168  double hh = h;
169  double a = 1.;
170  double b = 1.;
171  double etaCorrE = 1.;
172  double etaCorrH = 1.;
173  auto absEta = std::abs(eta);
174  t = min(999.9, max(tt, e + h));
175  if (t < 1.)
176  return;
177 
178  // Barrel calibration
179  if (absEta < 1.48) {
180  // The energy correction
181  a = e > 0. ? aBarrel(t) : 1.;
182  b = e > 0. ? bBarrel(t) : cBarrel(t);
183  double thresh = e > 0. ? threshE : threshH;
184 
185  // Protection against negative calibration
186  if (a < -0.25 || b < -0.25) {
187  a = 1.;
188  b = 1.;
189  thresh = 0.;
190  }
191 
192  // The new estimate of the true energy
193  t = min(999.9, max(tt, thresh + a * e + b * h));
194 
195  // The angular correction
196  if (e > 0. && thresh > 0.) {
197  etaCorrE = 1.0 + aEtaBarrelEH(t) + 1.3 * bEtaBarrelEH(t) * cEtaBarrelEH(absEta);
198  etaCorrH = 1.0;
199  } else {
200  etaCorrE = 1.0 + aEtaBarrelH(t) + 1.3 * bEtaBarrelH(t) * cEtaBarrelH(absEta);
201  etaCorrH = 1.0 + aEtaBarrelH(t) + bEtaBarrelH(t) * cEtaBarrelH(absEta);
202  }
203  if (e > 0. && thresh > 0.)
204  e = h > 0. ? threshE - threshH + etaCorrE * a * e : threshE + etaCorrE * a * e;
205  if (h > 0. && thresh > 0.) {
206  h = threshH + etaCorrH * b * h;
207  }
208 
209  // Endcap calibration
210  } else {
211  // The energy correction
212  a = e > 0. ? aEndcap(t) : 1.;
213  b = e > 0. ? bEndcap(t) : cEndcap(t);
214  double thresh = e > 0. ? threshE : threshH;
215 
216  // Protection against negative calibration
217  if (a < -0.25 || b < -0.25) {
218  a = 1.;
219  b = 1.;
220  thresh = 0.;
221  }
222 
223  // The new estimate of the true energy
224  t = min(999.9, max(tt, thresh + a * e + b * h));
225 
226  // The angular correction
227  double dEta = std::abs(absEta - 1.5);
228  double etaPow = dEta * dEta * dEta * dEta;
229 
230  if (e > 0. && thresh > 0.) {
231  if (absEta < 2.5) {
232  etaCorrE = 1. + aEtaEndcapEH(t) + bEtaEndcapEH(t) * cEtaEndcapEH(absEta);
233  } else {
234  etaCorrE = 1. + aEtaEndcapEH(t) + 1.3 * bEtaEndcapEH(t) * dEtaEndcapEH(absEta);
235  }
236 
237  etaCorrH = 1. + aEtaEndcapEH(t) + bEtaEndcapEH(t) * (0.04 + etaPow);
238  } else {
239  etaCorrE = 1.;
240  if (absEta < 2.5) {
241  etaCorrH = 1. + aEtaEndcapH(t) + bEtaEndcapH(t) * cEtaEndcapH(absEta);
242  } else {
243  etaCorrH = 1. + aEtaEndcapH(t) + bEtaEndcapH(t) * dEtaEndcapH(absEta);
244  }
245  }
246 
247  //t = min(999.9,max(tt, thresh + etaCorrE*a*e + etaCorrH*b*h));
248 
249  if (e > 0. && thresh > 0.)
250  e = h > 0. ? threshE - threshH + etaCorrE * a * e : threshE + etaCorrE * a * e;
251  if (h > 0. && thresh > 0.) {
252  h = threshH + etaCorrH * b * h;
253  }
254  }
255 
256  // Protection
257  if (e < 0. || h < 0.) {
258  // Some protection against crazy calibration
259  if (e < 0.)
260  e = ee;
261  if (h < 0.)
262  h = hh;
263  }
264 
265  // And that's it !
266 }
267 
268 // The calibration functions
269 double PFEnergyCalibration::aBarrel(double x) const {
270  if (pfCalibrations) {
274 
275  } else {
276  return faBarrel->Eval(x);
277  }
278 }
279 
280 double PFEnergyCalibration::bBarrel(double x) const {
281  if (pfCalibrations) {
285 
286  } else {
287  return fbBarrel->Eval(x);
288  }
289 }
290 
291 double PFEnergyCalibration::cBarrel(double x) const {
292  if (pfCalibrations) {
296 
297  } else {
298  return fcBarrel->Eval(x);
299  }
300 }
301 
302 double PFEnergyCalibration::aEtaBarrelEH(double x) const {
303  if (pfCalibrations) {
307 
308  } else {
309  return faEtaBarrelEH->Eval(x);
310  }
311 }
312 
313 double PFEnergyCalibration::bEtaBarrelEH(double x) const {
314  if (pfCalibrations) {
318 
319  } else {
320  return fbEtaBarrelEH->Eval(x);
321  }
322 }
323 
324 double PFEnergyCalibration::aEtaBarrelH(double x) const {
325  if (pfCalibrations) {
329 
330  } else {
331  return faEtaBarrelH->Eval(x);
332  }
333 }
334 
335 double PFEnergyCalibration::bEtaBarrelH(double x) const {
336  if (pfCalibrations) {
340 
341  } else {
342  return fbEtaBarrelH->Eval(x);
343  }
344 }
345 
346 double PFEnergyCalibration::aEndcap(double x) const {
347  if (pfCalibrations) {
351 
352  } else {
353  return faEndcap->Eval(x);
354  }
355 }
356 
357 double PFEnergyCalibration::bEndcap(double x) const {
358  if (pfCalibrations) {
362 
363  } else {
364  return fbEndcap->Eval(x);
365  }
366 }
367 
368 double PFEnergyCalibration::cEndcap(double x) const {
369  if (pfCalibrations) {
373 
374  } else {
375  return fcEndcap->Eval(x);
376  }
377 }
378 
379 double PFEnergyCalibration::aEtaEndcapEH(double x) const {
380  if (pfCalibrations) {
384 
385  } else {
386  return faEtaEndcapEH->Eval(x);
387  }
388 }
389 
390 double PFEnergyCalibration::bEtaEndcapEH(double x) const {
391  if (pfCalibrations) {
395 
396  } else {
397  return fbEtaEndcapEH->Eval(x);
398  }
399 }
400 
401 double PFEnergyCalibration::aEtaEndcapH(double x) const {
402  if (pfCalibrations) {
406 
407  } else {
408  return faEtaEndcapH->Eval(x);
409  }
410 }
411 
412 double PFEnergyCalibration::bEtaEndcapH(double x) const {
413  if (pfCalibrations) {
417  } else {
418  return fbEtaEndcapH->Eval(x);
419  }
420 }
421 
422 //added by Bhumika Kansal on 3 august 2018
423 
424 double PFEnergyCalibration::cEtaBarrelH(double x) const {
425  if (pfCalibrations) {
429 
430  } else {
431  return fcEtaBarrelH->Eval(x);
432  }
433 }
434 double PFEnergyCalibration::cEtaEndcapH(double x) const {
435  if (pfCalibrations) {
439 
440  } else {
441  return fcEtaEndcapH->Eval(x);
442  }
443 }
444 
445 double PFEnergyCalibration::dEtaEndcapH(double x) const {
446  if (pfCalibrations) {
450 
451  } else {
452  return fdEtaEndcapH->Eval(x);
453  }
454 }
455 
456 double PFEnergyCalibration::cEtaBarrelEH(double x) const {
457  if (pfCalibrations) {
461 
462  } else {
463  return fcEtaBarrelEH->Eval(x);
464  }
465 }
466 
467 double PFEnergyCalibration::cEtaEndcapEH(double x) const {
468  if (pfCalibrations) {
472 
473  } else {
474  return fcEtaEndcapEH->Eval(x);
475  }
476 }
477 
478 double PFEnergyCalibration::dEtaEndcapEH(double x) const {
479  if (pfCalibrations) {
483 
484  } else {
485  return fdEtaEndcapEH->Eval(x);
486  }
487 }
488 //
490  std::vector<double>& EclustersPS1,
491  std::vector<double>& EclustersPS2,
492  bool crackCorrection) const {
493  double ePS1(std::accumulate(EclustersPS1.begin(), EclustersPS1.end(), 0.0));
494  double ePS2(std::accumulate(EclustersPS2.begin(), EclustersPS2.end(), 0.0));
495  return energyEm(clusterEcal, ePS1, ePS2, crackCorrection);
496 }
497 
499  double ePS1,
500  double ePS2,
501  bool crackCorrection) const {
502  double eEcal = clusterEcal.energy();
503  //temporaty ugly fix
504  reco::PFCluster myPFCluster = clusterEcal;
505  myPFCluster.calculatePositionREP();
506  double eta = myPFCluster.positionREP().eta();
507  double phi = myPFCluster.positionREP().phi();
508 
509  double calibrated = Ecorr(eEcal, ePS1, ePS2, eta, phi, crackCorrection);
510  // if(eEcal!=0 && calibrated==0) std::cout<<"Eecal = "<<eEcal<<" eta = "<<eta<<" phi = "<<phi<<std::endl;
511  return calibrated;
512 }
513 
515  std::vector<double>& EclustersPS1,
516  std::vector<double>& EclustersPS2,
517  double& ps1,
518  double& ps2,
519  bool crackCorrection) const {
520  double ePS1(std::accumulate(EclustersPS1.begin(), EclustersPS1.end(), 0.0));
521  double ePS2(std::accumulate(EclustersPS2.begin(), EclustersPS2.end(), 0.0));
522  return energyEm(clusterEcal, ePS1, ePS2, ps1, ps2, crackCorrection);
523 }
525  double ePS1,
526  double ePS2,
527  double& ps1,
528  double& ps2,
529  bool crackCorrection) const {
530  double eEcal = clusterEcal.energy();
531  //temporaty ugly fix
532  reco::PFCluster myPFCluster = clusterEcal;
533  myPFCluster.calculatePositionREP();
534  double eta = myPFCluster.positionREP().eta();
535  double phi = myPFCluster.positionREP().phi();
536 
537  double calibrated = Ecorr(eEcal, ePS1, ePS2, eta, phi, ps1, ps2, crackCorrection);
538  // if(eEcal!=0 && calibrated==0) std::cout<<"Eecal = "<<eEcal<<" eta = "<<eta<<" phi = "<<phi<<std::endl;
539  return calibrated;
540 }
541 
542 std::ostream& operator<<(std::ostream& out, const PFEnergyCalibration& calib) {
543  if (!out)
544  return out;
545 
546  out << "PFEnergyCalibration -- " << endl;
547 
548  if (calib.pfCalibrations) {
549  static const std::map<std::string, PerformanceResult::ResultType> functType = {
550  {"PFfa_BARREL", PerformanceResult::PFfa_BARREL},
551  {"PFfa_ENDCAP", PerformanceResult::PFfa_ENDCAP},
552  {"PFfb_BARREL", PerformanceResult::PFfb_BARREL},
553  {"PFfb_ENDCAP", PerformanceResult::PFfb_ENDCAP},
554  {"PFfc_BARREL", PerformanceResult::PFfc_BARREL},
555  {"PFfc_ENDCAP", PerformanceResult::PFfc_ENDCAP},
556  {"PFfaEta_BARRELH", PerformanceResult::PFfaEta_BARRELH},
557  {"PFfaEta_ENDCAPH", PerformanceResult::PFfaEta_ENDCAPH},
558  {"PFfbEta_BARRELH", PerformanceResult::PFfbEta_BARRELH},
559  {"PFfbEta_ENDCAPH", PerformanceResult::PFfbEta_ENDCAPH},
560  {"PFfaEta_BARRELEH", PerformanceResult::PFfaEta_BARRELEH},
561  {"PFfaEta_ENDCAPEH", PerformanceResult::PFfaEta_ENDCAPEH},
562  {"PFfbEta_BARRELEH", PerformanceResult::PFfbEta_BARRELEH},
563  {"PFfbEta_ENDCAPEH", PerformanceResult::PFfbEta_ENDCAPEH},
564  {"PFfcEta_BARRELH", PerformanceResult::PFfcEta_BARRELH},
565  {"PFfcEta_ENDCAPH", PerformanceResult::PFfcEta_ENDCAPH},
566  {"PFfdEta_ENDCAPH", PerformanceResult::PFfdEta_ENDCAPH},
567  {"PFfcEta_BARRELEH", PerformanceResult::PFfcEta_BARRELEH},
568  {"PFfcEta_ENDCAPEH", PerformanceResult::PFfcEta_ENDCAPEH},
569  {"PFfdEta_ENDCAPEH", PerformanceResult::PFfdEta_ENDCAPEH}
570 
571  };
572 
573  for (std::map<std::string, PerformanceResult::ResultType>::const_iterator func = functType.begin();
574  func != functType.end();
575  ++func) {
576  cout << "Function: " << func->first << endl;
577  PerformanceResult::ResultType fType = func->second;
578  calib.pfCalibrations->printFormula(fType);
579  }
580 
581  } else {
582  std::cout << "Default calibration functions : " << std::endl;
583 
584  calib.faBarrel->Print();
585  calib.fbBarrel->Print();
586  calib.fcBarrel->Print();
587  calib.faEtaBarrelEH->Print();
588  calib.fbEtaBarrelEH->Print();
589  calib.faEtaBarrelH->Print();
590  calib.fbEtaBarrelH->Print();
591  calib.faEndcap->Print();
592  calib.fbEndcap->Print();
593  calib.fcEndcap->Print();
594  calib.faEtaEndcapEH->Print();
595  calib.fbEtaEndcapEH->Print();
596  calib.faEtaEndcapH->Print();
597  calib.fbEtaEndcapH->Print();
598  //
599  }
600 
601  return out;
602 }
603 
616 
622 
623 //useful to compute the signed distance to the closest crack in the barrel
624 double PFEnergyCalibration::minimum(double a, double b) const {
625  if (TMath::Abs(b) < TMath::Abs(a))
626  a = b;
627  return a;
628 }
629 
630 namespace {
631  constexpr double pi = M_PI; // 3.14159265358979323846;
632 
633  std::vector<double> fillcPhi() {
634  std::vector<double> retValue;
635  retValue.resize(18, 0);
636  retValue[0] = 2.97025;
637  for (unsigned i = 1; i <= 17; ++i)
638  retValue[i] = retValue[0] - 2 * i * pi / 18;
639 
640  return retValue;
641  }
642 
643  //Location of the 18 phi-cracks
644  const std::vector<double> cPhi = fillcPhi();
645 } // namespace
646 
647 //compute the unsigned distance to the closest phi-crack in the barrel
648 double PFEnergyCalibration::dCrackPhi(double phi, double eta) const {
649  //Shift of this location if eta<0
650  constexpr double delta_cPhi = 0.00638;
651 
652  double m; //the result
653 
654  //the location is shifted
655  if (eta < 0)
656  phi += delta_cPhi;
657 
658  if (phi >= -pi && phi <= pi) {
659  //the problem of the extrema
660  if (phi < cPhi[17] || phi >= cPhi[0]) {
661  if (phi < 0)
662  phi += 2 * pi;
663  m = minimum(phi - cPhi[0], phi - cPhi[17] - 2 * pi);
664  }
665 
666  //between these extrema...
667  else {
668  bool OK = false;
669  unsigned i = 16;
670  while (!OK) {
671  if (phi < cPhi[i]) {
672  m = minimum(phi - cPhi[i + 1], phi - cPhi[i]);
673  OK = true;
674  } else
675  i -= 1;
676  }
677  }
678  } else {
679  m = 0.; //if there is a problem, we assum that we are in a crack
680  std::cout << "Problem in dminphi" << std::endl;
681  }
682  if (eta < 0)
683  m = -m; //because of the disymetry
684  return m;
685 }
686 
687 // corrects the effect of phi-cracks
688 double PFEnergyCalibration::CorrPhi(double phi, double eta) const {
689  // we use 3 gaussians to correct the phi-cracks effect
690  constexpr double p1 = 5.59379e-01;
691  constexpr double p2 = -1.26607e-03;
692  constexpr double p3 = 9.61133e-04;
693 
694  constexpr double p4 = 1.81691e-01;
695  constexpr double p5 = -4.97535e-03;
696  constexpr double p6 = 1.31006e-03;
697 
698  constexpr double p7 = 1.38498e-01;
699  constexpr double p8 = 1.18599e-04;
700  constexpr double p9 = 2.01858e-03;
701 
702  double dminphi = dCrackPhi(phi, eta);
703 
704  double result =
705  (1 + p1 * TMath::Gaus(dminphi, p2, p3) + p4 * TMath::Gaus(dminphi, p5, p6) + p7 * TMath::Gaus(dminphi, p8, p9));
706 
707  return result;
708 }
709 
710 // corrects the effect of |eta|-cracks
711 double PFEnergyCalibration::CorrEta(double eta) const {
712  // we use a gaussian with a screwness for each of the 5 |eta|-cracks
713  constexpr double a[] = {6.13349e-01, 5.08146e-01, 4.44480e-01, 3.3487e-01, 7.65627e-01}; // amplitude
714  constexpr double m[] = {-1.79514e-02, 4.44747e-01, 7.92824e-01, 1.14090e+00, 1.47464e+00}; // mean
715  constexpr double s[] = {7.92382e-03, 3.06028e-03, 3.36139e-03, 3.94521e-03, 8.63950e-04}; // sigma
716  constexpr double sa[] = {1.27228e+01, 3.81517e-02, 1.63507e-01, -6.56480e-02, 1.87160e-01}; // screwness amplitude
717  constexpr double ss[] = {5.48753e-02, -1.00223e-02, 2.22866e-03, 4.26288e-04, 2.67937e-03}; // screwness sigma
718  double result = 1;
719 
720  for (unsigned i = 0; i <= 4; i++)
721  result += a[i] * TMath::Gaus(eta, m[i], s[i]) *
722  (1 + sa[i] * TMath::Sign(1., eta - m[i]) * TMath::Exp(-TMath::Abs(eta - m[i]) / ss[i]));
723 
724  return result;
725 }
726 
727 //corrects the global behaviour in the barrel
728 double PFEnergyCalibration::CorrBarrel(double E, double eta) const {
729  //Energy dependency
730  /*
731  //YM Parameters 52XX:
732  constexpr double p0=1.00000e+00;
733  constexpr double p1=3.27753e+01;
734  constexpr double p2=2.28552e-02;
735  constexpr double p3=3.06139e+00;
736  constexpr double p4=2.25135e-01;
737  constexpr double p5=1.47824e+00;
738  constexpr double p6=1.09e-02;
739  constexpr double p7=4.19343e+01;
740  */
741  constexpr double p0 = 0.9944;
742  constexpr double p1 = 9.827;
743  constexpr double p2 = 1.503;
744  constexpr double p3 = 1.196;
745  constexpr double p4 = 0.3349;
746  constexpr double p5 = 0.89;
747  constexpr double p6 = 0.004361;
748  constexpr double p7 = 51.51;
749  //Eta dependency
750  constexpr double p8 = 2.705593e-03;
751 
752  double result =
753  (p0 + 1 / (p1 + p2 * TMath::Power(E, p3)) + p4 * TMath::Exp(-E / p5) + p6 * TMath::Exp(-E * E / (p7 * p7))) *
754  (1 + p8 * eta * eta);
755 
756  return result;
757 }
758 
767 
768 //Alpha, Beta, Gamma give the weight of each sub-detector (PS layer1, PS layer2 and Ecal) in the areas of the endcaps where there is a PS
769 // Etot = Beta*eEcal + Gamma*(ePS1 + Alpha*ePS2)
770 
771 double PFEnergyCalibration::Alpha(double eta) const {
772  //Energy dependency
773  constexpr double p0 = 5.97621e-01;
774 
775  //Eta dependency
776  constexpr double p1 = -1.86407e-01;
777  constexpr double p2 = 3.85197e-01;
778 
779  //so that <feta()> = 1
780  constexpr double norm = (p1 + p2 * (2.6 + 1.656) / 2);
781 
782  double result = p0 * (p1 + p2 * eta) / norm;
783 
784  return result;
785 }
786 
787 double PFEnergyCalibration::Beta(double E, double eta) const {
788  //Energy dependency
789  constexpr double p0 = 0.032;
790  constexpr double p1 = 9.70394e-02;
791  constexpr double p2 = 2.23072e+01;
792  constexpr double p3 = 100;
793 
794  //Eta dependency
795  constexpr double p4 = 1.02496e+00;
796  constexpr double p5 = -4.40176e-03;
797 
798  //so that <feta()> = 1
799  constexpr double norm = (p4 + p5 * (2.6 + 1.656) / 2);
800 
801  double result = (1.0012 + p0 * TMath::Exp(-E / p3) + p1 * TMath::Exp(-E / p2)) * (p4 + p5 * eta) / norm;
802  return result;
803 }
804 
805 double PFEnergyCalibration::Gamma(double etaEcal) const {
806  //Energy dependency
807  constexpr double p0 = 2.49752e-02;
808 
809  //Eta dependency
810  constexpr double p1 = 6.48816e-02;
811  constexpr double p2 = -1.59517e-02;
812 
813  //so that <feta()> = 1
814  constexpr double norm = (p1 + p2 * (2.6 + 1.656) / 2);
815 
816  double result = p0 * (p1 + p2 * etaEcal) / norm;
817 
818  return result;
819 }
820 
826 
827 // returns the corrected energy in the barrel (0,1.48)
828 // Global Behaviour, phi and eta cracks are taken into account
829 double PFEnergyCalibration::EcorrBarrel(double E, double eta, double phi, bool crackCorrection) const {
830  // double result = E*CorrBarrel(E,eta)*CorrEta(eta)*CorrPhi(phi,eta);
831  double correction = crackCorrection ? std::max(CorrEta(eta), CorrPhi(phi, eta)) : 1.;
832  double result = E * CorrBarrel(E, eta) * correction;
833 
834  return result;
835 }
836 
837 // returns the corrected energy in the area between the barrel and the PS (1.48,1.65)
838 double PFEnergyCalibration::EcorrZoneBeforePS(double E, double eta) const {
839  //Energy dependency
840  constexpr double p0 = 1;
841  constexpr double p1 = 0.18;
842  constexpr double p2 = 8.;
843 
844  //Eta dependency
845  constexpr double p3 = 0.3;
846  constexpr double p4 = 1.11;
847  constexpr double p5 = 0.025;
848  constexpr double p6 = 1.49;
849  constexpr double p7 = 0.6;
850 
851  //so that <feta()> = 1
852  constexpr double norm = 1.21;
853 
854  double result = E * (p0 + p1 * TMath::Exp(-E / p2)) * (p3 + p4 * TMath::Gaus(eta, p6, p5) + p7 * eta) / norm;
855 
856  return result;
857 }
858 
859 // returns the corrected energy in the PS (1.65,2.6)
860 // only when (ePS1>0)||(ePS2>0)
861 double PFEnergyCalibration::EcorrPS(double eEcal, double ePS1, double ePS2, double etaEcal) const {
862  // gives the good weights to each subdetector
863  double E = Beta(1.0155 * eEcal + 0.025 * (ePS1 + 0.5976 * ePS2) / 9e-5, etaEcal) * eEcal +
864  Gamma(etaEcal) * (ePS1 + Alpha(etaEcal) * ePS2) / 9e-5;
865 
866  //Correction of the residual energy dependency
867  constexpr double p0 = 1.00;
868  constexpr double p1 = 2.18;
869  constexpr double p2 = 1.94;
870  constexpr double p3 = 4.13;
871  constexpr double p4 = 1.127;
872 
873  double result = E * (p0 + p1 * TMath::Exp(-E / p2) - p3 * TMath::Exp(-E / p4));
874 
875  return result;
876 }
877 
878 // returns the corrected energy in the PS (1.65,2.6)
879 // only when (ePS1>0)||(ePS2>0)
881  double eEcal, double ePS1, double ePS2, double etaEcal, double& outputPS1, double& outputPS2) const {
882  // gives the good weights to each subdetector
883  double gammaprime = Gamma(etaEcal) / 9e-5;
884 
885  if (outputPS1 == 0 && outputPS2 == 0 && esEEInterCalib_ != nullptr) {
886  // both ES planes working
887  // scaling factor accounting for data-mc
888  outputPS1 = gammaprime * ePS1 * esEEInterCalib_->getGammaLow0();
889  outputPS2 = gammaprime * Alpha(etaEcal) * ePS2 * esEEInterCalib_->getGammaLow3();
890  } else if (outputPS1 == 0 && outputPS2 == -1 && esEEInterCalib_ != nullptr) {
891  // ESP1 only working
892  double corrTotES = gammaprime * ePS1 * esEEInterCalib_->getGammaLow0() * esEEInterCalib_->getGammaLow1();
893  outputPS1 = gammaprime * ePS1 * esEEInterCalib_->getGammaLow0();
894  outputPS2 = corrTotES - outputPS1;
895  } else if (outputPS1 == -1 && outputPS2 == 0 && esEEInterCalib_ != nullptr) {
896  // ESP2 only working
897  double corrTotES =
898  gammaprime * Alpha(etaEcal) * ePS2 * esEEInterCalib_->getGammaLow3() * esEEInterCalib_->getGammaLow2();
899  outputPS2 = gammaprime * Alpha(etaEcal) * ePS2 * esEEInterCalib_->getGammaLow3();
900  outputPS1 = corrTotES - outputPS2;
901  } else {
902  // none working
903  outputPS1 = gammaprime * ePS1;
904  outputPS2 = gammaprime * Alpha(etaEcal) * ePS2;
905  }
906 
907  double E = Beta(1.0155 * eEcal + 0.025 * (ePS1 + 0.5976 * ePS2) / 9e-5, etaEcal) * eEcal + outputPS1 + outputPS2;
908 
909  //Correction of the residual energy dependency
910  constexpr double p0 = 1.00;
911  constexpr double p1 = 2.18;
912  constexpr double p2 = 1.94;
913  constexpr double p3 = 4.13;
914  constexpr double p4 = 1.127;
915 
916  double corrfac = (p0 + p1 * TMath::Exp(-E / p2) - p3 * TMath::Exp(-E / p4));
917  outputPS1 *= corrfac;
918  outputPS2 *= corrfac;
919  double result = E * corrfac;
920 
921  return result;
922 }
923 
924 // returns the corrected energy in the PS (1.65,2.6)
925 // only when (ePS1=0)&&(ePS2=0)
926 double PFEnergyCalibration::EcorrPS_ePSNil(double eEcal, double eta) const {
927  //Energy dependency
928  constexpr double p0 = 1.02;
929  constexpr double p1 = 0.165;
930  constexpr double p2 = 6.5;
931  constexpr double p3 = 2.1;
932 
933  //Eta dependency
934  constexpr double p4 = 1.02496e+00;
935  constexpr double p5 = -4.40176e-03;
936 
937  //so that <feta()> = 1
938  constexpr double norm = (p4 + p5 * (2.6 + 1.656) / 2);
939 
940  double result = eEcal * (p0 + p1 * TMath::Exp(-TMath::Abs(eEcal - p3) / p2)) * (p4 + p5 * eta) / norm;
941 
942  return result;
943 }
944 
945 // returns the corrected energy in the area between the end of the PS and the end of the endcap (2.6,2.98)
946 double PFEnergyCalibration::EcorrZoneAfterPS(double E, double eta) const {
947  //Energy dependency
948  constexpr double p0 = 1;
949  constexpr double p1 = 0.058;
950  constexpr double p2 = 12.5;
951  constexpr double p3 = -1.05444e+00;
952  constexpr double p4 = -5.39557e+00;
953  constexpr double p5 = 8.38444e+00;
954  constexpr double p6 = 6.10998e-01;
955 
956  //Eta dependency
957  constexpr double p7 = 1.06161e+00;
958  constexpr double p8 = 0.41;
959  constexpr double p9 = 2.918;
960  constexpr double p10 = 0.0181;
961  constexpr double p11 = 2.05;
962  constexpr double p12 = 2.99;
963  constexpr double p13 = 0.0287;
964 
965  //so that <feta()> = 1
966  constexpr double norm = 1.045;
967 
968  double result = E * (p0 + p1 * TMath::Exp(-(E - p3) / p2) + 1 / (p4 + p5 * TMath::Power(E, p6))) *
969  (p7 + p8 * TMath::Gaus(eta, p9, p10) + p11 * TMath::Gaus(eta, p12, p13)) / norm;
970  return result;
971 }
972 
973 // returns the corrected energy everywhere
974 // this work should be improved between 1.479 and 1.52 (junction barrel-endcap)
976  double eEcal, double ePS1, double ePS2, double eta, double phi, bool crackCorrection) const {
977  constexpr double endBarrel = 1.48;
978  constexpr double beginingPS = 1.65;
979  constexpr double endPS = 2.6;
980  constexpr double endEndCap = 2.98;
981 
982  double result = 0;
983 
984  eta = TMath::Abs(eta);
985 
986  if (eEcal > 0) {
987  if (eta <= endBarrel)
988  result = EcorrBarrel(eEcal, eta, phi, crackCorrection);
989  else if (eta <= beginingPS)
990  result = EcorrZoneBeforePS(eEcal, eta);
991  else if ((eta < endPS) && ePS1 == 0 && ePS2 == 0)
992  result = EcorrPS_ePSNil(eEcal, eta);
993  else if (eta < endPS)
994  result = EcorrPS(eEcal, ePS1, ePS2, eta);
995  else if (eta < endEndCap)
996  result = EcorrZoneAfterPS(eEcal, eta);
997  else
998  result = eEcal;
999  } else
1000  result = eEcal; // useful if eEcal=0 or eta>2.98
1001  //protection
1002  if (result < eEcal)
1003  result = eEcal;
1004  return result;
1005 }
1006 
1007 // returns the corrected energy everywhere
1008 // this work should be improved between 1.479 and 1.52 (junction barrel-endcap)
1009 double PFEnergyCalibration::Ecorr(double eEcal,
1010  double ePS1,
1011  double ePS2,
1012  double eta,
1013  double phi,
1014  double& ps1,
1015  double& ps2,
1016  bool crackCorrection) const {
1017  constexpr double endBarrel = 1.48;
1018  constexpr double beginingPS = 1.65;
1019  constexpr double endPS = 2.6;
1020  constexpr double endEndCap = 2.98;
1021 
1022  double result = 0;
1023 
1024  eta = TMath::Abs(eta);
1025 
1026  if (eEcal > 0) {
1027  if (eta <= endBarrel)
1028  result = EcorrBarrel(eEcal, eta, phi, crackCorrection);
1029  else if (eta <= beginingPS)
1030  result = EcorrZoneBeforePS(eEcal, eta);
1031  else if ((eta < endPS) && ePS1 == 0 && ePS2 == 0)
1032  result = EcorrPS_ePSNil(eEcal, eta);
1033  else if (eta < endPS)
1034  result = EcorrPS(eEcal, ePS1, ePS2, eta, ps1, ps2);
1035  else if (eta < endEndCap)
1036  result = EcorrZoneAfterPS(eEcal, eta);
1037  else
1038  result = eEcal;
1039  } else
1040  result = eEcal; // useful if eEcal=0 or eta>2.98
1041  // protection
1042  if (result < eEcal)
1043  result = eEcal;
1044  return result;
1045 }
std::unique_ptr< TF1 > faEtaBarrelH
std::unique_ptr< TF1 > fcEtaEndcapEH
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
const PerformancePayloadFromTFormula * pfCalibrations
double aEtaEndcapEH(double x) const
float getResult(PerformanceResult::ResultType, const BinningPointByMap &) const override
std::unique_ptr< TF1 > faEtaBarrelEH
double bBarrel(double x) const
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:46
double cEtaEndcapEH(double x) const
T Sign(T A, T B)
Definition: MathUtil.h:54
double aEtaBarrelEH(double x) const
double aEndcap(double x) const
double bEndcap(double x) const
const ESEEIntercalibConstants * esEEInterCalib_
double bEtaEndcapH(double x) const
#define nullptr
double Ecorr(double eEcal, double ePS1, double ePS2, double eta, double phi, bool crackCorrection=true) const
std::unique_ptr< TF1 > fbEtaBarrelEH
double aEtaBarrelH(double x) const
double EcorrZoneBeforePS(double E, double eta) const
double cEtaBarrelH(double x) const
double Alpha(double eta) const
void energyEmHad(double t, double &e, double &h, double eta, double phi) const
const Double_t pi
double cEtaBarrelEH(double x) const
std::unique_ptr< TF1 > fcEtaBarrelEH
void calculatePositionREP()
computes posrep_ once and for all
Definition: PFCluster.h:99
double bEtaBarrelEH(double x) const
double p4[4]
Definition: TauolaWrapper.h:92
std::unique_ptr< TF1 > fbBarrel
std::unique_ptr< TF1 > faEtaEndcapEH
friend std::ostream & operator<<(std::ostream &out, const PFEnergyCalibration &calib)
double cEndcap(double x) const
const REPPoint & positionREP() const
cluster position: rho, eta, phi
Definition: PFCluster.h:96
T Abs(T a)
Definition: MathUtil.h:49
double Beta(double E, double eta) const
std::unique_ptr< TF1 > faBarrel
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::unique_ptr< TF1 > faEndcap
T min(T a, T b)
Definition: MathUtil.h:58
double dEtaEndcapEH(double x) const
std::unique_ptr< TF1 > fcEtaBarrelH
double p2[4]
Definition: TauolaWrapper.h:90
double energy() const
cluster energy
Definition: PFCluster.h:78
#define M_PI
double bEtaEndcapEH(double x) const
double EcorrZoneAfterPS(double E, double eta) const
std::pair< int, edm::FunctionWithDict > OK
Definition: findMethod.cc:126
bool insert(BinningVariables::BinningVariablesType, float)
double bEtaBarrelH(double x) const
double dEtaEndcapH(double x) const
std::unique_ptr< TF1 > fbEtaBarrelH
std::unique_ptr< TF1 > fdEtaEndcapEH
double b
Definition: hdecay.h:118
std::unique_ptr< TF1 > faEtaEndcapH
double CorrBarrel(double E, double eta) const
double Gamma(double etaEcal) const
double energyEm(const reco::PFCluster &clusterEcal, std::vector< double > &EclustersPS1, std::vector< double > &EclustersPS2, bool crackCorrection=true) const
double dCrackPhi(double phi, double eta) const
double p1[4]
Definition: TauolaWrapper.h:89
double a
Definition: hdecay.h:119
std::unique_ptr< TF1 > fcEtaEndcapH
std::unique_ptr< TF1 > fbEndcap
std::unique_ptr< TF1 > fcBarrel
void printFormula(PerformanceResult::ResultType res) const
std::unique_ptr< TF1 > fdEtaEndcapH
double CorrEta(double eta) const
double aEtaEndcapH(double x) const
double aBarrel(double x) const
double CorrPhi(double phi, double eta) const
double EcorrPS_ePSNil(double eEcal, double eta) const
std::unique_ptr< TF1 > fbEtaEndcapEH
double cBarrel(double x) const
double cEtaEndcapH(double x) const
double EcorrBarrel(double E, double eta, double phi, bool crackCorrection=true) const
std::unique_ptr< TF1 > fcEndcap
#define constexpr
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
std::unique_ptr< TF1 > fbEtaEndcapH
double EcorrPS(double eEcal, double ePS1, double ePS2, double etaEcal) const
double p3[4]
Definition: TauolaWrapper.h:91
double minimum(double a, double b) const