CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFEnergyCalibration.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterTools_PFEnergyCalibration_h
2 #define RecoParticleFlow_PFClusterTools_PFEnergyCalibration_h
3 
5 
6 class TF1;
7 
8 // -*- C++ -*-
9 //
10 // Package: PFClusterTools
11 // Class: PFEnergyCalibration
12 //
28 //
29 // Original Author: Christian Veelken
30 // Created: Tue Aug 8 16:26:18 CDT 2006
31 // $Id: PFEnergyCalibration.h,v 1.10 2009/05/19 18:59:33 pjanot Exp $
32 //
33 //
34 
35 #include <iostream>
36 
37 //#include "FWCore/ParameterSet/interface/ParameterSet.h"
38 
40 {
41  public:
42  PFEnergyCalibration(); // default constructor;
43  // needed by PFRootEvent
44 
45  PFEnergyCalibration( double e_slope,
46  double e_offset,
47  double eh_eslope,
48  double eh_hslope,
49  double eh_offset,
50  double h_slope,
51  double h_offset,
52  double h_damping,
53  unsigned newCalib = 0);
54 
56 
57  // ecal calibration
58  double energyEm(double uncalibratedEnergyECAL,
59  double eta=0, double phi=0) const;
60 
61  double energyEm(const reco::PFCluster& clusterEcal,
62  std::vector<double> &EclustersPS1,
63  std::vector<double> &EclustersPS2,
64  bool crackCorrection = true);
65 
66  double energyEm(const reco::PFCluster& clusterEcal,
67  std::vector<double> &EclustersPS1,
68  std::vector<double> &EclustersPS2,
69  double &ps1,double&ps2,
70  bool crackCorrection=true);
71 
72  // HCAL only calibration
73  double energyHad(double uncalibratedEnergyHCAL,
74  double eta=0, double phi=0) const;
75 
76 
77  // ECAL+HCAL (abc) calibration
78  double energyEmHad(double uncalibratedEnergyECAL,
79  double uncalibratedEnergyHCAL,
80  double eta=0, double phi=0) const;
81 
82  // ECAL+HCAL (abc-alpha-beta) calibration, with E and eta dependent coefficients
83  void energyEmHad(double t, double& e, double&h, double eta, double phi) const;
84 
85  // set calibration parameters for energy deposits of electrons and photons in ECAL; this member function is needed by PFRootEvent
87  double paramECAL_offset);
88 
89  // Initialize E- and eta-dependent coefficient functional form
91 
92  double paramECAL_slope() const {return paramECAL_slope_;}
93 
94  double paramECAL_offset() const {return paramECAL_offset_;}
95 
96  double paramECALplusHCAL_slopeECAL() const {
98  }
99 
102  }
103 
105 
106  double paramHCAL_slope() const {return paramHCAL_slope_;}
107  double paramHCAL_offset() const {return paramHCAL_offset_;}
108  double paramHCAL_damping() const {return paramHCAL_damping_;}
109 
110 
111  friend std::ostream& operator<<(std::ostream& out,
112  const PFEnergyCalibration& calib);
113 
114  protected:
117 
121 
125 
126  private:
127 
128  double minimum(double a,double b);
129  double dCrackPhi(double phi, double eta);
130  double CorrPhi(double phi, double eta);
131  double CorrEta(double eta);
132  double CorrBarrel(double E, double eta);
133  double Alpha(double eta);
134  double Beta(double E, double eta);
135  double Gamma(double etaEcal);
136  double EcorrBarrel(double E, double eta, double phi, bool crackCorrection=true);
137  double EcorrZoneBeforePS(double E, double eta);
138  double EcorrPS(double eEcal,double ePS1,double ePS2,double etaEcal);
139  double EcorrPS(double eEcal,double ePS1,double ePS2,double etaEcal,double&, double&);
140  double EcorrPS_ePSNil(double eEcal,double eta);
141  double EcorrZoneAfterPS(double E, double eta);
142  double Ecorr(double eEcal,double ePS1,double ePS2,double eta,double phi,bool crackCorrection=true);
143  double Ecorr(double eEcal,double ePS1,double ePS2,double eta,double phi,double&,double&,bool crackCorrection=true);
144 
145  // Barrel calibration (eta 0.00 -> 1.48)
146  TF1* faBarrel;
147  TF1* fbBarrel;
148  TF1* fcBarrel;
149  TF1* faEtaBarrel;
150  TF1* fbEtaBarrel;
151 
152  // Endcap calibration (eta 1.48 -> 3.xx)
153  TF1* faEndcap;
154  TF1* fbEndcap;
155  TF1* fcEndcap;
156  TF1* faEtaEndcap;
157  TF1* fbEtaEndcap;
158 
159  // Threshold correction (offset)
160  double threshE, threshH;
161 
162 };
163 
164 #endif
165 
166 
double Beta(double E, double eta)
double EcorrZoneAfterPS(double E, double eta)
double Alpha(double eta)
double paramHCAL_slope() const
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
double EcorrBarrel(double E, double eta, double phi, bool crackCorrection=true)
double EcorrZoneBeforePS(double E, double eta)
double EcorrPS_ePSNil(double eEcal, double eta)
double CorrEta(double eta)
double paramECAL_slope() const
double paramECAL_offset() const
T eta() const
double Gamma(double etaEcal)
MVATrainerComputer * calib
Definition: MVATrainer.cc:64
double CorrPhi(double phi, double eta)
double CorrBarrel(double E, double eta)
friend std::ostream & operator<<(std::ostream &out, const PFEnergyCalibration &calib)
double paramECALplusHCAL_offset() const
double paramECALplusHCAL_slopeHCAL() const
double paramHCAL_damping() const
double energyEmHad(double uncalibratedEnergyECAL, double uncalibratedEnergyHCAL, double eta=0, double phi=0) const
double dCrackPhi(double phi, double eta)
double energyEm(double uncalibratedEnergyECAL, double eta=0, double phi=0) const
tuple out
Definition: dbtoconf.py:99
double b
Definition: hdecay.h:120
void setCalibrationParametersEm(double paramECAL_slope, double paramECAL_offset)
double EcorrPS(double eEcal, double ePS1, double ePS2, double etaEcal)
double energyHad(double uncalibratedEnergyHCAL, double eta=0, double phi=0) const
double Ecorr(double eEcal, double ePS1, double ePS2, double eta, double phi, bool crackCorrection=true)
double paramECALplusHCAL_slopeECAL() const
double a
Definition: hdecay.h:121
double paramHCAL_offset() const
double minimum(double a, double b)
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
Definition: DDAxes.h:10