CMS 3D CMS Logo

EgammaHcalIsolation.h
Go to the documentation of this file.
1 #ifndef EgammaIsolationAlgos_EgammaHcalIsolation_h
2 #define EgammaIsolationAlgos_EgammaHcalIsolation_h
3 //*****************************************************************************
4 // File: EgammaHcalIsolation.h
5 // ----------------------------------------------------------------------------
6 // OrigAuth: Matthias Mozer
7 // Institute: IIHE-VUB
8 //=============================================================================
9 //*****************************************************************************
10 
11 //C++ includes
12 #include <vector>
13 #include <functional>
14 
15 //CMSSW includes
25 
26 //Sum helper functions
27 double scaleToE(const double &eta);
28 double scaleToEt(const double &eta);
29 
31 public:
32  enum HcalDepth { AllDepths = 0, Depth1 = 1, Depth2 = 2 };
33 
34  //constructors
36  double intRadius,
37  double eLowB,
38  double eLowE,
39  double etLowB,
40  double etLowE,
41  edm::ESHandle<CaloGeometry> theCaloGeom,
42  const HBHERecHitCollection &mhbhe);
43 
44  //destructor
46 
47  //AllDepths
48  double getHcalESum(const reco::Candidate *c) const { return getHcalESum(c->get<reco::SuperClusterRef>().get()); }
49  double getHcalEtSum(const reco::Candidate *c) const { return getHcalEtSum(c->get<reco::SuperClusterRef>().get()); }
50  double getHcalESum(const reco::SuperCluster *sc) const { return getHcalESum(sc->position()); }
51  double getHcalEtSum(const reco::SuperCluster *sc) const { return getHcalEtSum(sc->position()); }
52  double getHcalESum(const math::XYZPoint &p) const { return getHcalESum(GlobalPoint(p.x(), p.y(), p.z())); }
53  double getHcalEtSum(const math::XYZPoint &p) const { return getHcalEtSum(GlobalPoint(p.x(), p.y(), p.z())); }
54  double getHcalESum(const GlobalPoint &pclu) const { return getHcalSum(pclu, AllDepths, &scaleToE); }
55  double getHcalEtSum(const GlobalPoint &pclu) const { return getHcalSum(pclu, AllDepths, &scaleToEt); }
56 
57  //Depth1
58  double getHcalESumDepth1(const reco::Candidate *c) const {
60  }
61  double getHcalEtSumDepth1(const reco::Candidate *c) const {
63  }
64  double getHcalESumDepth1(const reco::SuperCluster *sc) const { return getHcalESumDepth1(sc->position()); }
65  double getHcalEtSumDepth1(const reco::SuperCluster *sc) const { return getHcalEtSumDepth1(sc->position()); }
66  double getHcalESumDepth1(const math::XYZPoint &p) const {
67  return getHcalESumDepth1(GlobalPoint(p.x(), p.y(), p.z()));
68  }
69  double getHcalEtSumDepth1(const math::XYZPoint &p) const {
70  return getHcalEtSumDepth1(GlobalPoint(p.x(), p.y(), p.z()));
71  }
72  double getHcalESumDepth1(const GlobalPoint &pclu) const { return getHcalSum(pclu, Depth1, &scaleToE); }
73  double getHcalEtSumDepth1(const GlobalPoint &pclu) const { return getHcalSum(pclu, Depth1, &scaleToEt); }
74 
75  //Depth2
76  double getHcalESumDepth2(const reco::Candidate *c) const {
78  }
79  double getHcalEtSumDepth2(const reco::Candidate *c) const {
81  }
82  double getHcalESumDepth2(const reco::SuperCluster *sc) const { return getHcalESumDepth2(sc->position()); }
83  double getHcalEtSumDepth2(const reco::SuperCluster *sc) const { return getHcalEtSumDepth2(sc->position()); }
84  double getHcalESumDepth2(const math::XYZPoint &p) const {
85  return getHcalESumDepth2(GlobalPoint(p.x(), p.y(), p.z()));
86  }
87  double getHcalEtSumDepth2(const math::XYZPoint &p) const {
88  return getHcalEtSumDepth2(GlobalPoint(p.x(), p.y(), p.z()));
89  }
90  double getHcalESumDepth2(const GlobalPoint &pclu) const { return getHcalSum(pclu, Depth2, &scaleToE); }
91  double getHcalEtSumDepth2(const GlobalPoint &pclu) const { return getHcalSum(pclu, Depth2, &scaleToEt); }
92 
93 private:
94  bool isDepth2(const DetId &) const;
95  double getHcalSum(const GlobalPoint &, const HcalDepth &, double (*)(const double &)) const;
96 
97  double extRadius_;
98  double intRadius_;
99  double eLowB_;
100  double eLowE_;
101  double etLowB_;
102  double etLowE_;
103 
106 
108 };
109 
110 #endif
bool isDepth2(const DetId &) const
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:153
double getHcalESumDepth2(const math::XYZPoint &p) const
edm::ESHandle< CaloGeometry > theCaloGeom_
double getHcalEtSumDepth1(const reco::SuperCluster *sc) const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
double getHcalEtSumDepth2(const reco::SuperCluster *sc) const
double getHcalEtSumDepth2(const reco::Candidate *c) const
CaloDualConeSelector< HBHERecHit > * doubleConeSel_
double getHcalEtSum(const reco::Candidate *c) const
double getHcalESum(const math::XYZPoint &p) const
double getHcalESumDepth1(const GlobalPoint &pclu) const
double getHcalESumDepth2(const reco::SuperCluster *sc) const
double getHcalESumDepth1(const math::XYZPoint &p) const
const HBHERecHitCollection & mhbhe_
double getHcalEtSum(const reco::SuperCluster *sc) const
double getHcalEtSum(const math::XYZPoint &p) const
double getHcalESum(const reco::SuperCluster *sc) const
double getHcalEtSumDepth2(const math::XYZPoint &p) const
double getHcalESum(const GlobalPoint &pclu) const
double getHcalSum(const GlobalPoint &, const HcalDepth &, double(*)(const double &)) const
double getHcalESumDepth2(const reco::Candidate *c) const
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
double getHcalEtSumDepth1(const math::XYZPoint &p) const
double getHcalESumDepth2(const GlobalPoint &pclu) const
EgammaHcalIsolation(double extRadius, double intRadius, double eLowB, double eLowE, double etLowB, double etLowE, edm::ESHandle< CaloGeometry > theCaloGeom, const HBHERecHitCollection &mhbhe)
Definition: DetId.h:17
double getHcalEtSum(const GlobalPoint &pclu) const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
double getHcalESumDepth1(const reco::SuperCluster *sc) const
double getHcalEtSumDepth1(const GlobalPoint &pclu) const
double getHcalESum(const reco::Candidate *c) const
double scaleToEt(const double &eta)
double scaleToE(const double &eta)
double getHcalEtSumDepth1(const reco::Candidate *c) const
double getHcalEtSumDepth2(const GlobalPoint &pclu) const
T get() const
get a component
Definition: Candidate.h:222
double getHcalESumDepth1(const reco::Candidate *c) const