CMS 3D CMS Logo

EgammaHLTTrackIsolation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: EgammaHLTAlgos
4 // Class : EgammaHLTTrackIsolation
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author: Monica Vazquez Acosta
10 // Created: Tue Jun 13 12:17:19 CEST 2006
11 //
12 
13 // system include files
14 
15 // user include files
17 
18 std::pair<int, float> EgammaHLTTrackIsolation::electronIsolation(const reco::Track* const tr,
19  const reco::TrackCollection* isoTracks) {
20  GlobalPoint vtx(0, 0, tr->vertex().z());
21  const reco::Track::Vector& p = tr->momentum();
22  GlobalVector mom(p.x(), p.y(), p.z());
23  return findIsoTracks(mom, vtx, isoTracks, true);
24 }
25 
26 std::pair<int, float> EgammaHLTTrackIsolation::electronIsolation(const reco::Track* const tr,
27  const reco::TrackCollection* isoTracks,
28  GlobalPoint zvtx) {
29  // Just to insure consistency with no-vertex-code
30  GlobalPoint vtx(0, 0, zvtx.z());
31  const reco::Track::Vector& p = tr->momentum();
32  GlobalVector mom(p.x(), p.y(), p.z());
33  return findIsoTracks(mom, vtx, isoTracks, true);
34 }
35 
36 std::pair<int, float> EgammaHLTTrackIsolation::electronIsolation(const reco::Track* const tr,
37  const reco::ElectronCollection* allEle,
38  const reco::TrackCollection* isoTracks) {
39  GlobalPoint vtx(0, 0, tr->vertex().z());
40  const reco::Track::Vector& p = tr->momentum();
41  GlobalVector mom(p.x(), p.y(), p.z());
42  return findIsoTracksWithoutEle(mom, vtx, allEle, isoTracks);
43 }
44 
45 std::pair<int, float> EgammaHLTTrackIsolation::photonIsolation(const reco::RecoCandidate* const recocandidate,
46  const reco::TrackCollection* isoTracks,
47  bool useVertex) {
48  if (useVertex) {
49  GlobalPoint vtx(0, 0, recocandidate->vertex().z());
50  return photonIsolation(recocandidate, isoTracks, vtx);
51  } else {
52  reco::RecoCandidate::Point pos = recocandidate->superCluster()->position();
53  GlobalVector mom(pos.x(), pos.y(), pos.z());
54  return findIsoTracks(mom, GlobalPoint(), isoTracks, false, false);
55  }
56 }
57 
58 std::pair<int, float> EgammaHLTTrackIsolation::photonIsolation(const reco::RecoCandidate* const recocandidate,
59  const reco::TrackCollection* isoTracks,
60  GlobalPoint zvtx) {
61  // to insure consistency with no-free-vertex-code
62  GlobalPoint vtx(0, 0, zvtx.z());
63 
64  reco::RecoCandidate::Point pos = recocandidate->superCluster()->position();
65  GlobalVector mom(pos.x() - vtx.x(), pos.y() - vtx.y(), pos.z() - vtx.z());
66 
67  return findIsoTracks(mom, vtx, isoTracks, false);
68 }
69 
70 std::pair<int, float> EgammaHLTTrackIsolation::photonIsolation(const reco::RecoCandidate* const recocandidate,
71  const reco::ElectronCollection* allEle,
72  const reco::TrackCollection* isoTracks) {
73  reco::RecoCandidate::Point pos = recocandidate->superCluster()->position();
74  GlobalVector mom(pos.x(), pos.y(), pos.z());
75  return findIsoTracksWithoutEle(mom, GlobalPoint(), allEle, isoTracks);
76 }
77 
79  GlobalVector mom, GlobalPoint vtx, const reco::TrackCollection* isoTracks, bool isElectron, bool useVertex) {
80  // Check that reconstructed tracks fit within cone boundaries,
81  // (Note: tracks will not always stay within boundaries)
82  int ntrack = 0;
83  float ptSum = 0.;
84 
85  for (reco::TrackCollection::const_iterator trItr = isoTracks->begin(); trItr != isoTracks->end(); ++trItr) {
86  GlobalPoint ivtx(trItr->vertex().x(), trItr->vertex().y(), trItr->vertex().z());
87  reco::Track::Vector ip = trItr->momentum();
88  GlobalVector imom(ip.x(), ip.y(), ip.z());
89 
90  float pt = imom.perp();
91  float dperp = 0.;
92  float dz = 0.;
93  float deta = 0.;
94  float dphi = 0.;
95  if (useVertex) {
96  dperp = ivtx.perp() - vtx.perp();
97  dz = ivtx.z() - vtx.z();
98  deta = imom.eta() - mom.eta();
99  dphi = imom.phi() - mom.phi();
100  } else {
101  //in case of unkown photon vertex, modify direction of photon to point from
102  //current track vertex to sc instead of from (0.,0.,0.) to sc. In this
103  //way, all tracks are considered based on direction alone.
104  GlobalVector mom_temp = mom - GlobalVector(ivtx.x(), ivtx.y(), ivtx.z());
105  deta = imom.eta() - mom_temp.eta();
106  dphi = imom.phi() - mom_temp.phi();
107  }
108  // Correct dmom_phi's from [-2pi,2pi] to [-pi,pi]
109  if (dphi > M_PI)
110  dphi = dphi - 2 * M_PI;
111  else if (dphi < -M_PI)
112  dphi = dphi + 2 * M_PI;
113 
114  float R = sqrt(dphi * dphi + deta * deta);
115 
116  // Apply boundary cut
117  // bool selected=false;
118 
119  // if (pt > ptMin && R < conesize &&
120  // fabs(dperp) < rspan && fabs(dz) < zspan) selected=true;
121 
122  // if (selected) {
123  // ntrack++;
124  // if (!isElectron || R > vetoConesize) ptSum+=pt; //to exclude electron track
125  // }
126  // float theVetoVar = R;
127  // if (isElectron) theVetoVar = R;
128 
129  //hmm how do I figure out if this is barrel or endcap?
130  //abs(mom.eta())<1.5 is the obvious choice but that will be electron not detector eta for electrons
131  //well lets leave it as that for now, its what reco does (well with eta=1.479)
132  double innerStrip = fabs(mom.eta()) < 1.479 ? stripBarrel : stripEndcap;
133 
134  if (pt > ptMin && R < conesize && R > vetoConesize && fabs(dperp) < rspan && fabs(dz) < zspan &&
135  fabs(deta) >= innerStrip) {
136  ntrack++;
137  ptSum += pt;
138  }
139  }
140 
141  // if (isElectron) ntrack-=1; //to exclude electron track
142 
143  return (std::pair<int, float>(ntrack, ptSum));
144 }
145 
148  const reco::ElectronCollection* allEle,
149  const reco::TrackCollection* isoTracks) {
150  // Check that reconstructed tracks fit within cone boundaries,
151  // (Note: tracks will not always stay within boundaries)
152  int iele = 0;
153  int ntrack = 0;
154  float ptSum = 0.;
155  std::vector<float> etaele;
156  std::vector<float> phiele;
157 
158  // std::cout << "allEle.size() = " << allEle->size() << std::endl;
159 
160  // Store ALL electrons eta and phi
161  for (reco::ElectronCollection::const_iterator iElectron = allEle->begin(); iElectron != allEle->end(); iElectron++) {
162  iele++;
163  reco::TrackRef anothereletrackref = iElectron->track();
164  etaele.push_back(anothereletrackref->momentum().eta());
165  phiele.push_back(anothereletrackref->momentum().phi());
166  // std::cout << "Electron " << iele << ": phi = " << anothereletrackref->momentum().phi() << ", eta = " << anothereletrackref->momentum().eta() << std::endl;
167  }
168 
169  for (reco::TrackCollection::const_iterator trItr = isoTracks->begin(); trItr != isoTracks->end(); ++trItr) {
170  GlobalPoint ivtx(trItr->vertex().x(), trItr->vertex().y(), trItr->vertex().z());
171  reco::Track::Vector ip = trItr->momentum();
172  GlobalVector imom(ip.x(), ip.y(), ip.z());
173 
174  float pt = imom.perp();
175  float dperp = ivtx.perp() - vtx.perp();
176  float dz = ivtx.z() - vtx.z();
177  float deta = imom.eta() - mom.eta();
178  float dphi = imom.phi() - mom.phi();
179 
180  // Correct dmom_phi's from [-2pi,2pi] to [-pi,pi]
181  if (dphi > M_PI)
182  dphi = dphi - 2 * M_PI;
183  else if (dphi < -M_PI)
184  dphi = dphi + 2 * M_PI;
185 
186  float R = sqrt(dphi * dphi + deta * deta);
187 
188  // Apply boundary cut
189  bool selected = false;
190  bool passedconeveto = true;
191 
192  //hmm how do I figure out if this is barrel or endcap?
193  //abs(mom.eta())<1.5 is the obvious choice but that will be electron not detector eta for electrons
194  //well lets leave it as that for now, its what reco does (well with eta=1.479)
195  double innerStrip = fabs(mom.eta()) < 1.479 ? stripBarrel : stripEndcap;
196 
197  if (pt > ptMin && R < conesize && fabs(dperp) < rspan && fabs(dz) < zspan && fabs(deta) >= innerStrip)
198  selected = true;
199 
200  // Check that NO electron is counted in the isolation
201  for (unsigned int eleItr = 0; eleItr < etaele.size(); ++eleItr) {
202  deta = etaele[eleItr] - imom.eta();
203  dphi = phiele[eleItr] - imom.phi();
204 
205  // Correct dmom_phi's from [-2pi,2pi] to [-pi,pi]
206  if (dphi > M_PI)
207  dphi = dphi - 2 * M_PI;
208  else if (dphi < -M_PI)
209  dphi = dphi + 2 * M_PI;
210 
211  R = sqrt(dphi * dphi + deta * deta);
212  if (R < vetoConesize)
213  passedconeveto = false;
214  }
215 
216  if (selected && passedconeveto) {
217  ntrack++;
218  ptSum += pt; //to exclude electron tracks
219  }
220  }
221 
222  // ntrack-=1; //to exclude electron track
223 
224  return (std::pair<int, float>(ntrack, ptSum));
225 }
Vector3DBase
Definition: Vector3DBase.h:8
ntrack
#define ntrack
reco::isElectron
bool isElectron(const Candidate &part)
Definition: pdgIdUtils.h:7
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
pos
Definition: PixelAliasList.h:18
reco::RecoCandidate::superCluster
virtual reco::SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: RecoCandidate.cc:25
EgammaHLTTrackIsolation::stripBarrel
double stripBarrel
Definition: EgammaHLTTrackIsolation.h:179
GlobalVector
Global3DVector GlobalVector
Definition: GlobalVector.h:10
edm::Ref< TrackCollection >
EgammaHLTTrackIsolation::electronIsolation
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.
Definition: EgammaHLTTrackIsolation.cc:18
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
EgammaHLTTrackIsolation::conesize
double conesize
Definition: EgammaHLTTrackIsolation.h:173
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
reco::Track
Definition: Track.h:27
EgammaHLTTrackIsolation.h
GlobalPoint
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Point3DBase< float, GlobalTag >
EgammaHLTTrackIsolation::zspan
double zspan
Definition: EgammaHLTTrackIsolation.h:174
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
PV3DBase::eta
T eta() const
Definition: PV3DBase.h:73
HLT_FULL_cff.useVertex
useVertex
Definition: HLT_FULL_cff.py:26537
EgammaHLTTrackIsolation::findIsoTracks
std::pair< int, float > findIsoTracks(GlobalVector mom, GlobalPoint vtx, const reco::TrackCollection *isoTracks, bool isElectron, bool useVertex=true)
Definition: EgammaHLTTrackIsolation.cc:78
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
EgammaHLTTrackIsolation::photonIsolation
std::pair< int, float > photonIsolation(const reco::RecoCandidate *const recocand, const reco::TrackCollection *isoTracks, bool useVertex)
Definition: EgammaHLTTrackIsolation.cc:45
EgammaHLTTrackIsolation::findIsoTracksWithoutEle
std::pair< int, float > findIsoTracksWithoutEle(GlobalVector mom, GlobalPoint vtx, const reco::ElectronCollection *allEle, const reco::TrackCollection *isoTracks)
Definition: EgammaHLTTrackIsolation.cc:146
reco::RecoCandidate
Definition: RecoCandidate.h:20
reco::ElectronCollection
std::vector< Electron > ElectronCollection
collectin of Electron objects
Definition: ElectronFwd.h:9
reco::TrackBase::vertex
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead
Definition: TrackBase.h:676
reco::LeafCandidate::vertex
const Point & vertex() const override
vertex position (overwritten by PF...)
Definition: LeafCandidate.h:165
EgammaHLTTrackIsolation::vetoConesize
double vetoConesize
Definition: EgammaHLTTrackIsolation.h:176
EgammaHLTTrackIsolation::stripEndcap
double stripEndcap
Definition: EgammaHLTTrackIsolation.h:180
extraflags_cff.vtx
vtx
Definition: extraflags_cff.py:19
PVValHelper::dz
Definition: PVValidationHelpers.h:51
EgammaHLTTrackIsolation::ptMin
double ptMin
Definition: EgammaHLTTrackIsolation.h:172
reco::Candidate::Point
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
reco::TrackBase::momentum
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:664
EgammaHLTTrackIsolation::rspan
double rspan
Definition: EgammaHLTTrackIsolation.h:175
dttmaxenums::R
Definition: DTTMax.h:29
reco::TrackCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
PV3DBase::phi
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
reco::TrackBase::Vector
math::XYZVector Vector
spatial vector
Definition: TrackBase.h:77