CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PFEnergyResolution.cc
Go to the documentation of this file.
2 #include <TMath.h>
3 #include <cmath>
4 /*
5 PFEnergyResolution::PFEnergyResolution(const edm::ParameterSet& parameters)
6 {
7 //--- nothing to be done yet
8 }
9 */
11  //--- nothing to be done yet
12 }
13 
15  //--- nothing to be done yet
16 }
17 
18 double PFEnergyResolution::getEnergyResolutionEm(double CorrectedEnergy, double eta) const {
19  //The parameters S,N,C has been determined with the Full Sim on CMSSW_2_1_0_pre4.
20  //The resolution must be a function of the corrected energy available in PFEnergyCalibration
21  //Jonathan Biteau July 2008
22 
23  double C;
24  double S;
25  double N;
26  if (TMath::Abs(eta) < 1.48) {
27  C = 0.35 / 100;
28  S = 5.51 / 100;
29  N = 98. / 1000.;
30  } else {
31  C = 0;
32  S = 12.8 / 100;
33  N = 440. / 1000.;
34  }
35  double result = TMath::Sqrt(C * C * CorrectedEnergy * CorrectedEnergy + S * S * CorrectedEnergy + N * N);
36  return result;
37 }
38 
39 double PFEnergyResolution::getEnergyResolutionHad(double energyHCAL, double eta, double phi) const {
40  //--- estimate **relative** resolution of energy measurement (sigmaE/E)
41  // for hadrons in depositing energy in HCAL
42  // (eta and phi dependence not implemented yet)
43 
44  return 1.49356 / sqrt(energyHCAL) + 6.62527e-03 * sqrt(energyHCAL) - 6.33966e-02;
45 }
46 /*
47 double PFEnergyResolution::getEnergyResolutionHad(double energyECAL, double energyHCAL, double eta, double phi) const
48 {
49 //--- estimate **relative** resolution of energy measurement (sigmaE/E)
50 // for hadrons depositing energy in ECAL and HCAL
51 // (currently, the resolution function for hadrons
52 // is assumed to be the same in ECAL and HCAL)
53 
54  return getEnergyResolutionHad(energyECAL + energyHCAL, theta, phi);
55 }
56 */
double getEnergyResolutionHad(double energyHCAL, double eta, double phi) const
tuple result
Definition: mps_fire.py:311
double getEnergyResolutionEm(double CorrectedEnergy, double eta) const
T sqrt(T t)
Definition: SSEVec.h:19
T Abs(T a)
Definition: MathUtil.h:49
#define N
Definition: blowfish.cc:9
double energyHCAL(std::vector< DetId > &vdets, edm::Handle< T > &hits, double hbThr=-100, double heThr=-100, double hfThr=-100, double hoThr=-100, double tMin=-500, double tMax=500, int useRaw=0, bool debug=false)
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:97