CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // $Id: EgammaHLTTrackIsolation.h,v 1.5 2009/01/20 11:31:29 covarell Exp $
19 //
20 
21 
30 
32 
36 
41 
42 
43 
45 
46 
47 
49 {
50 
51  public:
52 
53 
54  EgammaHLTTrackIsolation(double egTrkIso_PtMin,
55  double egTrkIso_ConeSize,
56  double egTrkIso_ZSpan,
57  double egTrkIso_RSpan,
58  double egTrkIso_VetoConeSize,
59  double egTrkIso_stripBarrel=0,
60  double egTrkIso_stripEndcap=0) :
61  ptMin(egTrkIso_PtMin),
62  conesize(egTrkIso_ConeSize),
63  zspan(egTrkIso_ZSpan),
64  rspan(egTrkIso_RSpan),
65  vetoConesize(egTrkIso_VetoConeSize),
66  stripBarrel(egTrkIso_stripBarrel),
67  stripEndcap(egTrkIso_stripEndcap)
68  {
69 
70  /*
71  std::cout << "EgammaHLTTrackIsolation instance:"
72  << " ptMin=" << ptMin << " "
73  << " conesize="<< conesize << " "
74  << " zspan=" << zspan << " "
75  << " rspan=" << rspan << " "
76  << " vetoConesize="<< vetoConesize
77  << std::endl;
78  */
79  }
80 
81 
83  std::pair<int,float> electronIsolation(const reco::Track * const tr, const reco::TrackCollection* isoTracks);
84  std::pair<int,float> electronIsolation(const reco::Track * const tr, const reco::ElectronCollection* allEle, const reco::TrackCollection* isoTracks);
85  std::pair<int,float> electronIsolation(const reco::Track * const tr, const reco::TrackCollection* isoTracks, GlobalPoint vertex);
86 
90  std::pair<int,float> photonIsolation(const reco::RecoCandidate * const recocand, const reco::TrackCollection* isoTracks, bool useVertex);
91  std::pair<int,float> photonIsolation(const reco::RecoCandidate * const recocand, const reco::TrackCollection* isoTracks, GlobalPoint vertex);
92  std::pair<int,float> photonIsolation(const reco::RecoCandidate * const recocand, const reco::ElectronCollection* allEle, const reco::TrackCollection* isoTracks);
93 
95  int electronTrackCount(const reco::Track * const tr, const reco::TrackCollection* isoTracks)
96  {return electronIsolation(tr,isoTracks).first;}
97  int electronTrackCount(const reco::Track * const tr, const reco::TrackCollection* isoTracks, GlobalPoint vertex)
98  {return electronIsolation(tr,isoTracks,vertex).first;}
99 
103  int photonTrackCount(const reco::RecoCandidate * const recocand, const reco::TrackCollection* isoTracks, bool useVertex)
104  {return photonIsolation(recocand,isoTracks,useVertex).first;}
105  int photonTrackCount(const reco::RecoCandidate * const recocand, const reco::TrackCollection* isoTracks, GlobalPoint vertex)
106  {return photonIsolation(recocand,isoTracks,vertex).first;}
107  int photonTrackCount(const reco::RecoCandidate * const recocand, const reco::ElectronCollection* allEle, const reco::TrackCollection* isoTracks)
108  {return photonIsolation(recocand,allEle,isoTracks).first;}
109 
111  float electronPtSum(const reco::Track * const tr, const reco::TrackCollection* isoTracks)
112  {return electronIsolation(tr,isoTracks).second;}
113  float electronPtSum(const reco::Track * const tr, const reco::TrackCollection* isoTracks, GlobalPoint vertex)
114  {return electronIsolation(tr,isoTracks,vertex).second;}
115  float electronPtSum(const reco::Track * const tr, const reco::ElectronCollection* allEle ,const reco::TrackCollection* isoTracks)
116  {return electronIsolation(tr,allEle,isoTracks).second;}
117 
121  float photonPtSum(const reco::RecoCandidate * const recocand, const reco::TrackCollection* isoTracks, bool useVertex)
122  {return photonIsolation(recocand,isoTracks, useVertex).second;}
123  float photonPtSum(const reco::RecoCandidate * const recocand, const reco::TrackCollection* isoTracks, GlobalPoint vertex)
124  {return photonIsolation(recocand,isoTracks, vertex).second;}
125  float photonPtSum(const reco::RecoCandidate * const recocand, const reco::ElectronCollection* allEle, const reco::TrackCollection* isoTracks)
126  {return photonIsolation(recocand,allEle,isoTracks).second;}
127 
128 
130  double getPtMin() { return ptMin;}
132  double getConeSize() { return conesize; }
134  double getZspan() {return zspan; }
136  double getRspan() { return rspan; }
138  double getvetoConesize() { return vetoConesize; }
139 
140  private:
141  // Call track reconstruction
142  std::pair<int,float> findIsoTracks(GlobalVector mom, GlobalPoint vtx, const reco::TrackCollection* isoTracks, bool isElectron, bool useVertex=true);
143  std::pair<int,float> findIsoTracksWithoutEle(GlobalVector mom, GlobalPoint vtx, const reco::ElectronCollection* allEle, const reco::TrackCollection* isoTracks);
144 
145  // Parameters of isolation cone geometry.
146  double ptMin;
147  double conesize;
148  double zspan;
149  double rspan;
150  double vetoConesize;
151 
152  //added for inner eta strip veto (I'll keep the violation of CMS naming conventions to be consistant)
153  double stripBarrel;
154  double stripEndcap;
155 
156 
157 
158 };
159 
160 
161 #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:10
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)