CMS 3D CMS Logo

EgammaHLTTrackIsolation.h
Go to the documentation of this file.
1 #ifndef EgammaHLTAlgos_EgammaHLTTrackIsolation_h
2 #define EgammaHLTAlgos_EgammaHLTTrackIsolation_h
3 // -*- C++ -*-
4 //
5 // Package: EgammaHLTAlgos
6 // Class : EgammaHLTTrackIsolation
7 //
15 //
16 // Original Author: Monica Vazquez Acosta - CERN
17 // Created: Tue Jun 13 12:19:32 CEST 2006
18 //
19 
28 
30 
34 
39 
41 
43 public:
44  EgammaHLTTrackIsolation(double egTrkIso_PtMin,
45  double egTrkIso_ConeSize,
46  double egTrkIso_ZSpan,
47  double egTrkIso_RSpan,
48  double egTrkIso_VetoConeSize,
49  double egTrkIso_stripBarrel = 0,
50  double egTrkIso_stripEndcap = 0)
51  : ptMin(egTrkIso_PtMin),
52  conesize(egTrkIso_ConeSize),
53  zspan(egTrkIso_ZSpan),
54  rspan(egTrkIso_RSpan),
55  vetoConesize(egTrkIso_VetoConeSize),
56  stripBarrel(egTrkIso_stripBarrel),
57  stripEndcap(egTrkIso_stripEndcap) {
58  /*
59  std::cout << "EgammaHLTTrackIsolation instance:"
60  << " ptMin=" << ptMin << " "
61  << " conesize="<< conesize << " "
62  << " zspan=" << zspan << " "
63  << " rspan=" << rspan << " "
64  << " vetoConesize="<< vetoConesize
65  << std::endl;
66  */
67  }
68 
70  std::pair<int, float> electronIsolation(const reco::Track* const tr, const reco::TrackCollection* isoTracks);
71  std::pair<int, float> electronIsolation(const reco::Track* const tr,
72  const reco::ElectronCollection* allEle,
73  const reco::TrackCollection* isoTracks);
74  std::pair<int, float> electronIsolation(const reco::Track* const tr,
75  const reco::TrackCollection* isoTracks,
77 
81  std::pair<int, float> photonIsolation(const reco::RecoCandidate* const recocand,
82  const reco::TrackCollection* isoTracks,
83  bool useVertex);
84  std::pair<int, float> photonIsolation(const reco::RecoCandidate* const recocand,
85  const reco::TrackCollection* isoTracks,
87  std::pair<int, float> photonIsolation(const reco::RecoCandidate* const recocand,
88  const reco::ElectronCollection* allEle,
89  const reco::TrackCollection* isoTracks);
90 
92  int electronTrackCount(const reco::Track* const tr, const reco::TrackCollection* isoTracks) {
93  return electronIsolation(tr, isoTracks).first;
94  }
95  int electronTrackCount(const reco::Track* const tr, const reco::TrackCollection* isoTracks, GlobalPoint vertex) {
96  return electronIsolation(tr, isoTracks, vertex).first;
97  }
98 
102  int photonTrackCount(const reco::RecoCandidate* const recocand,
103  const reco::TrackCollection* isoTracks,
104  bool useVertex) {
105  return photonIsolation(recocand, isoTracks, useVertex).first;
106  }
107  int photonTrackCount(const reco::RecoCandidate* const recocand,
108  const reco::TrackCollection* isoTracks,
110  return photonIsolation(recocand, isoTracks, vertex).first;
111  }
112  int photonTrackCount(const reco::RecoCandidate* const recocand,
113  const reco::ElectronCollection* allEle,
114  const reco::TrackCollection* isoTracks) {
115  return photonIsolation(recocand, allEle, isoTracks).first;
116  }
117 
119  float electronPtSum(const reco::Track* const tr, const reco::TrackCollection* isoTracks) {
120  return electronIsolation(tr, isoTracks).second;
121  }
122  float electronPtSum(const reco::Track* const tr, const reco::TrackCollection* isoTracks, GlobalPoint vertex) {
123  return electronIsolation(tr, isoTracks, vertex).second;
124  }
125  float electronPtSum(const reco::Track* const tr,
126  const reco::ElectronCollection* allEle,
127  const reco::TrackCollection* isoTracks) {
128  return electronIsolation(tr, allEle, isoTracks).second;
129  }
130 
134  float photonPtSum(const reco::RecoCandidate* const recocand, const reco::TrackCollection* isoTracks, bool useVertex) {
135  return photonIsolation(recocand, isoTracks, useVertex).second;
136  }
137  float photonPtSum(const reco::RecoCandidate* const recocand,
138  const reco::TrackCollection* isoTracks,
140  return photonIsolation(recocand, isoTracks, vertex).second;
141  }
142  float photonPtSum(const reco::RecoCandidate* const recocand,
143  const reco::ElectronCollection* allEle,
144  const reco::TrackCollection* isoTracks) {
145  return photonIsolation(recocand, allEle, isoTracks).second;
146  }
147 
149  double getPtMin() { return ptMin; }
151  double getConeSize() { return conesize; }
153  double getZspan() { return zspan; }
155  double getRspan() { return rspan; }
157  double getvetoConesize() { return vetoConesize; }
158 
159 private:
160  // Call track reconstruction
161  std::pair<int, float> findIsoTracks(GlobalVector mom,
163  const reco::TrackCollection* isoTracks,
164  bool isElectron,
165  bool useVertex = true);
166  std::pair<int, float> findIsoTracksWithoutEle(GlobalVector mom,
168  const reco::ElectronCollection* allEle,
169  const reco::TrackCollection* isoTracks);
170 
171  // Parameters of isolation cone geometry.
172  double ptMin;
173  double conesize;
174  double zspan;
175  double rspan;
176  double vetoConesize;
177 
178  //added for inner eta strip veto (I'll keep the violation of CMS naming conventions to be consistant)
179  double stripBarrel;
180  double stripEndcap;
181 };
182 
183 #endif
float photonPtSum(const reco::RecoCandidate *const recocand, const reco::TrackCollection *isoTracks, GlobalPoint vertex)
double getConeSize()
Get isolation cone size.
double getZspan()
Get maximum ivertex z-coordinate spread.
int photonTrackCount(const reco::RecoCandidate *const recocand, const reco::TrackCollection *isoTracks, GlobalPoint vertex)
double getRspan()
Get maximum transverse distance of ivertex from beam line.
std::pair< int, float > electronIsolation(const reco::Track *const tr, const reco::TrackCollection *isoTracks)
Get number of tracks and Pt sum of tracks inside an isolation cone for electrons. ...
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
int electronTrackCount(const reco::Track *const tr, const reco::TrackCollection *isoTracks)
Get number of tracks inside an isolation cone for electrons.
int photonTrackCount(const reco::RecoCandidate *const recocand, const reco::TrackCollection *isoTracks, bool useVertex)
float electronPtSum(const reco::Track *const tr, const reco::ElectronCollection *allEle, const reco::TrackCollection *isoTracks)
int photonTrackCount(const reco::RecoCandidate *const recocand, const reco::ElectronCollection *allEle, const reco::TrackCollection *isoTracks)
bool isElectron(const Candidate &part)
Definition: pdgIdUtils.h:7
float electronPtSum(const reco::Track *const tr, const reco::TrackCollection *isoTracks)
Get Pt sum of tracks inside an isolation cone for electrons.
float photonPtSum(const reco::RecoCandidate *const recocand, const reco::ElectronCollection *allEle, const reco::TrackCollection *isoTracks)
std::vector< Electron > ElectronCollection
collectin of Electron objects
Definition: ElectronFwd.h:9
std::pair< int, float > photonIsolation(const reco::RecoCandidate *const recocand, const reco::TrackCollection *isoTracks, bool useVertex)
EgammaHLTTrackIsolation(double egTrkIso_PtMin, double egTrkIso_ConeSize, double egTrkIso_ZSpan, double egTrkIso_RSpan, double egTrkIso_VetoConeSize, double egTrkIso_stripBarrel=0, double egTrkIso_stripEndcap=0)
std::pair< int, float > findIsoTracks(GlobalVector mom, GlobalPoint vtx, const reco::TrackCollection *isoTracks, bool isElectron, bool useVertex=true)
double getPtMin()
Get pt cut for itracks.
float electronPtSum(const reco::Track *const tr, const reco::TrackCollection *isoTracks, GlobalPoint vertex)
std::pair< int, float > findIsoTracksWithoutEle(GlobalVector mom, GlobalPoint vtx, const reco::ElectronCollection *allEle, const reco::TrackCollection *isoTracks)
double getvetoConesize()
Get veto cone size.
float photonPtSum(const reco::RecoCandidate *const recocand, const reco::TrackCollection *isoTracks, bool useVertex)
int electronTrackCount(const reco::Track *const tr, const reco::TrackCollection *isoTracks, GlobalPoint vertex)