CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PFEnergyCalibration.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterTools_PFEnergyCalibration_h
2 #define RecoParticleFlow_PFClusterTools_PFEnergyCalibration_h
3 
6 
9 
10 #include <TF1.h>
11 
12 // -*- C++ -*-
13 //
14 // Package: PFClusterTools
15 // Class: PFEnergyCalibration
16 //
32 //
33 // Original Author: Christian Veelken
34 // Created: Tue Aug 8 16:26:18 CDT 2006
35 //
36 //
37 
38 #include <iostream>
39 
40 //#include "FWCore/ParameterSet/interface/ParameterSet.h"
41 
43 public:
45 
46  // ecal calibration for photons
47  double energyEm(const reco::PFCluster& clusterEcal, double ePS1, double ePS2, bool crackCorrection = true) const;
48 
50  double clusterEnergy = 0.;
51  double ps1Energy = 0.;
52  double ps2Energy = 0.;
53  };
54 
56  reco::PFCluster const& eeCluster,
57  std::vector<reco::PFCluster const*> const& psClusterPointers,
58  ESChannelStatus const& channelStatus,
59  bool applyCrackCorrections) const;
60 
61  // ECAL+HCAL (abc) calibration, with E and eta dependent coefficients, for hadrons
62  void energyEmHad(double t, double& e, double& h, double eta, double phi) const;
63 
64  // Set the run-dependent calibration functions from the global tag
66 
68  esEEInterCalib_ = esEEInterCalib;
69  }
70 
71  friend std::ostream& operator<<(std::ostream& out, const PFEnergyCalibration& calib);
72 
73 private:
74  // ecal calibration for photons
75  double energyEm(const reco::PFCluster& clusterEcal,
76  double ePS1,
77  double ePS2,
78  double& ps1,
79  double& ps2,
80  bool crackCorrection = true) const;
81 
82  // Calibration functions from global tag
85 
86  // Barrel calibration (eta 0.00 -> 1.48)
87  std::unique_ptr<TF1> faBarrel;
88  std::unique_ptr<TF1> fbBarrel;
89  std::unique_ptr<TF1> fcBarrel;
90  std::unique_ptr<TF1> faEtaBarrelEH;
91  std::unique_ptr<TF1> fbEtaBarrelEH;
92  std::unique_ptr<TF1> faEtaBarrelH;
93  std::unique_ptr<TF1> fbEtaBarrelH;
94 
95  // Endcap calibration (eta 1.48 -> 3.xx)
96  std::unique_ptr<TF1> faEndcap;
97  std::unique_ptr<TF1> fbEndcap;
98  std::unique_ptr<TF1> fcEndcap;
99  std::unique_ptr<TF1> faEtaEndcapEH;
100  std::unique_ptr<TF1> fbEtaEndcapEH;
101  std::unique_ptr<TF1> faEtaEndcapH;
102  std::unique_ptr<TF1> fbEtaEndcapH;
103 
104  //added by Bhumika on 2 august 2018
105  std::unique_ptr<TF1> fcEtaBarrelEH;
106  std::unique_ptr<TF1> fcEtaEndcapEH;
107  std::unique_ptr<TF1> fdEtaEndcapEH;
108  std::unique_ptr<TF1> fcEtaBarrelH;
109  std::unique_ptr<TF1> fcEtaEndcapH;
110  std::unique_ptr<TF1> fdEtaEndcapH;
111 
112  double minimum(double a, double b) const;
113  double dCrackPhi(double phi, double eta) const;
114  double CorrPhi(double phi, double eta) const;
115  double CorrEta(double eta) const;
116  double CorrBarrel(double E, double eta) const;
117  double Alpha(double eta) const;
118  double Beta(double E, double eta) const;
119  double Gamma(double etaEcal) const;
120  double EcorrBarrel(double E, double eta, double phi, bool crackCorrection = true) const;
121  double EcorrZoneBeforePS(double E, double eta) const;
122  double EcorrPS(double eEcal, double ePS1, double ePS2, double etaEcal) const;
123  double EcorrPS(double eEcal, double ePS1, double ePS2, double etaEcal, double&, double&) const;
124  double EcorrPS_ePSNil(double eEcal, double eta) const;
125  double EcorrZoneAfterPS(double E, double eta) const;
126  double Ecorr(double eEcal, double ePS1, double ePS2, double eta, double phi, bool crackCorrection = true) const;
127  double Ecorr(double eEcal,
128  double ePS1,
129  double ePS2,
130  double eta,
131  double phi,
132  double&,
133  double&,
134  bool crackCorrection = true) const;
135 
136  // The calibration functions
137  double aBarrel(double x) const;
138  double bBarrel(double x) const;
139  double cBarrel(double x) const;
140  double aEtaBarrelEH(double x) const;
141  double bEtaBarrelEH(double x) const;
142  double aEtaBarrelH(double x) const;
143  double bEtaBarrelH(double x) const;
144  double aEndcap(double x) const;
145  double bEndcap(double x) const;
146  double cEndcap(double x) const;
147  double aEtaEndcapEH(double x) const;
148  double bEtaEndcapEH(double x) const;
149  double aEtaEndcapH(double x) const;
150  double bEtaEndcapH(double x) const;
151  //added by Bhumika on 3 august 2018
152  double cEtaBarrelEH(double x) const;
153  double cEtaEndcapEH(double x) const;
154  double dEtaEndcapEH(double x) const;
155  double cEtaBarrelH(double x) const;
156  double cEtaEndcapH(double x) const;
157  double dEtaEndcapH(double x) const;
158 
159  // Threshold correction (offset)
160  const double threshE = 3.5;
161  const double threshH = 2.5;
162 };
163 
164 #endif
std::unique_ptr< TF1 > faEtaBarrelH
std::unique_ptr< TF1 > fcEtaEndcapEH
const PerformancePayloadFromTFormula * pfCalibrations
double aEtaEndcapEH(double x) const
std::unique_ptr< TF1 > faEtaBarrelEH
CalibratedEndcapPFClusterEnergies calibrateEndcapClusterEnergies(reco::PFCluster const &eeCluster, std::vector< reco::PFCluster const * > const &psClusterPointers, ESChannelStatus const &channelStatus, bool applyCrackCorrections) const
double bBarrel(double x) const
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
double cEtaEndcapEH(double x) const
double energyEm(const reco::PFCluster &clusterEcal, double ePS1, double ePS2, bool crackCorrection=true) const
double aEtaBarrelEH(double x) const
double aEndcap(double x) const
double bEndcap(double x) const
const ESEEIntercalibConstants * esEEInterCalib_
double bEtaEndcapH(double x) const
double Ecorr(double eEcal, double ePS1, double ePS2, double eta, double phi, bool crackCorrection=true) const
void initAlphaGamma_ESplanes_fromDB(const ESEEIntercalibConstants *esEEInterCalib)
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
double cEtaBarrelEH(double x) const
std::unique_ptr< TF1 > fcEtaBarrelEH
double bEtaBarrelEH(double x) const
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
double Beta(double E, double eta) const
std::unique_ptr< TF1 > faBarrel
std::unique_ptr< TF1 > faEndcap
double dEtaEndcapEH(double x) const
std::unique_ptr< TF1 > fcEtaBarrelH
double bEtaEndcapEH(double x) const
double EcorrZoneAfterPS(double E, double eta) const
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 dCrackPhi(double phi, double eta) const
double a
Definition: hdecay.h:119
std::unique_ptr< TF1 > fcEtaEndcapH
std::unique_ptr< TF1 > fbEndcap
std::unique_ptr< TF1 > fcBarrel
void setCalibrationFunctions(const PerformancePayloadFromTFormula *thePFCal)
tuple applyCrackCorrections
std::unique_ptr< TF1 > fdEtaEndcapH
double CorrEta(double eta) const
double aEtaEndcapH(double x) const
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
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
std::unique_ptr< TF1 > fbEtaEndcapH
double EcorrPS(double eEcal, double ePS1, double ePS2, double etaEcal) const
double minimum(double a, double b) const