CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ConversionTools.cc
Go to the documentation of this file.
1 
2 #include <TMath.h>
12 
13 using namespace edm;
14 using namespace reco;
15 
16 
17 //--------------------------------------------------------------------------------------------------
18 bool ConversionTools::isGoodConversion(const Conversion &conv, const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
19 {
20 
21  //Check if a given conversion candidate passes the conversion selection cuts
22 
23  const reco::Vertex &vtx = conv.conversionVertex();
24 
25  //vertex validity
26  if (!vtx.isValid()) return false;
27 
28  //fit probability
29  if (TMath::Prob( vtx.chi2(), vtx.ndof() )<probMin) return false;
30 
31  //compute transverse decay length
33  double dbsx = vtx.x() - beamspot.x();
34  double dbsy = vtx.y() - beamspot.y();
35  double lxy = (mom.x()*dbsx + mom.y()*dbsy)/mom.rho();
36 
37  //transverse decay length
38  if ( lxy<lxyMin )
39  return false;
40 
41  //loop through daughters to check nhitsbeforevtx
42  for (std::vector<uint8_t>::const_iterator it = conv.nHitsBeforeVtx().begin(); it!=conv.nHitsBeforeVtx().end(); ++it) {
43  if ( (*it)>nHitsBeforeVtxMax ) return false;
44  }
45 
46  return true;
47 }
48 
49 //--------------------------------------------------------------------------------------------------
50 bool ConversionTools::matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch)
51 {
52 
53  //check if a given GsfElectron matches a given conversion (no quality cuts applied)
54  //matching is always attempted through the gsf track ref, and optionally attempted through the
55  //closest ctf track ref
56 
57  const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
58  for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) {
59  if ( ele.gsfTrack().isNonnull() && ele.gsfTrack().id()==it->id() && ele.gsfTrack().key()==it->key()) return true;
60  else if ( allowCkfMatch && ele.closestCtfTrackRef().isNonnull() && ele.closestCtfTrackRef().id()==it->id() && ele.closestCtfTrackRef().key()==it->key() ) return true;
61  }
62 
63  return false;
64 }
65 
66 //--------------------------------------------------------------------------------------------------
67 bool ConversionTools::matchesConversion(const reco::SuperCluster &sc, const reco::Conversion &conv, float dRMax, float dEtaMax, float dPhiMax) {
68 
69  //check if a given SuperCluster matches a given conversion (no quality cuts applied)
70  //matching is geometric between conversion momentum and vector joining conversion vertex
71  //to supercluster position
72 
73 
75 
76  math::XYZPoint scpos(sc.position());
78 
79 
80  math::XYZVector cscvector = scpos - cvtx;
81  float dR = reco::deltaR(mom,cscvector);
82  float dEta = mom.eta() - cscvector.eta();
83  float dPhi = reco::deltaPhi(mom.phi(),cscvector.phi());
84 
85  if (dR>dRMax) return false;
86  if (dEta>dEtaMax) return false;
87  if (dPhi>dPhiMax) return false;
88 
89  return true;
90 
91 }
92 
93 
94 //--------------------------------------------------------------------------------------------------
96 {
97 
98  //check if given track matches given conversion (matching by ref)
99 
100  if (trk.isNull()) return false;
101 
102  const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
103  for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) {
104  if (trk.id()==it->id() && trk.key()==it->key()) return true;
105  }
106 
107  return false;
108 }
109 
110 //--------------------------------------------------------------------------------------------------
112 {
113 
114  //check if given track matches given conversion (matching by ref)
115 
116  if (trk.isNull()) return false;
117 
118  const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
119  for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) {
120  if (trk.id()==it->id() && trk.key()==it->key()) return true;
121  }
122 
123  return false;
124 }
125 
126 //--------------------------------------------------------------------------------------------------
128 {
129 
130  //check if given track matches given conversion (matching by ref)
131 
132  if (trk.isNull()) return false;
133 
134  const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
135  for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) {
136  if (trk.id()==it->id() && trk.key()==it->key()) return true;
137  }
138 
139  return false;
140 }
141 
142 
143 //--------------------------------------------------------------------------------------------------
146  const math::XYZPoint &beamspot, bool allowCkfMatch, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
147 {
148  //check if a given electron candidate matches to at least one conversion candidate in the
149  //collection which also passes the selection cuts, optionally match with the closestckf track in
150  //in addition to just the gsf track (enabled in default arguments)
151 
152  for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
153  if (!matchesConversion(ele, *it, allowCkfMatch)) continue;
154  if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
155 
156  return true;
157  }
158 
159  return false;
160 
161 }
162 
163 //--------------------------------------------------------------------------------------------------
166  const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
167 {
168  //check if a given track matches to at least one conversion candidate in the
169  //collection which also passes the selection cuts
170 
171  if (trk.isNull()) return false;
172 
173  for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
174  if (!matchesConversion(trk, *it)) continue;
175  if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
176 
177  return true;
178  }
179 
180  return false;
181 
182 }
183 
184 //--------------------------------------------------------------------------------------------------
187  const math::XYZPoint &beamspot, float dRMax, float dEtaMax, float dPhiMax, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
188 {
189 
190  //check if a given SuperCluster matches to at least one conversion candidate in the
191  //collection which also passes the selection cuts
192 
193  for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
194  if (!matchesConversion(sc, *it)) continue;
195  if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
196 
197  return true;
198  }
199 
200  return false;
201 
202 }
203 
204 
205 //--------------------------------------------------------------------------------------------------
208  const math::XYZPoint &beamspot, bool allowCkfMatch, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
209 {
210  //check if a given electron candidate matches to at least one conversion candidate in the
211  //collection which also passes the selection cuts, optionally match with the closestckf track in
212  //in addition to just the gsf track (enabled in default arguments)
213  //If multiple conversions are found, returned reference corresponds to minimum
214  //conversion radius
215 
217 
218  double minRho = 999.;
219  for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
220  float rho = it->conversionVertex().position().rho();
221  if (rho>minRho) continue;
222  if (!matchesConversion(ele, *it, allowCkfMatch)) continue;
223  if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
224 
225  minRho = rho;
226  match = ConversionRef(convCol,it-convCol->begin());
227  }
228 
229  return match;
230 
231 }
232 
233 //--------------------------------------------------------------------------------------------------
236  const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
237 {
238  //check if a given track matches to at least one conversion candidate in the
239  //collection which also passes the selection cuts
240  //If multiple conversions are found, returned reference corresponds to minimum
241  //conversion radius
242 
244 
245  if (trk.isNull()) return match;
246 
247  double minRho = 999.;
248  for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
249  float rho = it->conversionVertex().position().rho();
250  if (rho>minRho) continue;
251  if (!matchesConversion(trk, *it)) continue;
252  if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
253 
254  minRho = rho;
255  match = ConversionRef(convCol,it-convCol->begin());
256  }
257 
258  return match;
259 
260 }
261 
262 //--------------------------------------------------------------------------------------------------
265  const math::XYZPoint &beamspot, float dRMax, float dEtaMax, float dPhiMax, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
266 {
267 
268  //check if a given SuperCluster matches to at least one conversion candidate in the
269  //collection which also passes the selection cuts
270  //If multiple conversions are found, returned reference corresponds to minimum
271  //conversion radius
272 
274 
275  double minRho = 999.;
276  for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
277  float rho = it->conversionVertex().position().rho();
278  if (rho>minRho) continue;
279  if (!matchesConversion(sc, *it, dRMax,dEtaMax,dPhiMax)) continue;
280  if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
281 
282  minRho = rho;
283  match = ConversionRef(convCol,it-convCol->begin());
284  }
285 
286  return match;
287 
288 }
289 
290 //--------------------------------------------------------------------------------------------------
292  const edm::Handle<reco::ConversionCollection> &convCol, const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
293 {
294 
295  //check if a given SuperCluster matches to at least one GsfElectron having zero expected inner hits
296  //and not matching any conversion in the collection passing the quality cuts
297 
298  if (sc.isNull()) return false;
299 
300  for (GsfElectronCollection::const_iterator it = eleCol->begin(); it!=eleCol->end(); ++it) {
301  //match electron to supercluster
302  if (it->superCluster()!=sc) continue;
303 
304  //check expected inner hits
305  if (it->gsfTrack()->trackerExpectedHitsInner().numberOfHits()>0) continue;
306 
307  //check if electron is matching to a conversion
308  if (hasMatchedConversion(*it,convCol,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
309 
310 
311  return true;
312  }
313 
314  return false;
315 
316 
317 }
318 
319 
320 //--------------------------------------------------------------------------------------------------
322  const edm::Handle<reco::ConversionCollection> &convCol, const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
323 {
324 
325  //check if a given SuperCluster matches to at least one GsfElectron having zero expected inner hits
326  //and not matching any conversion in the collection passing the quality cuts
327 
329 
330  if (sc.isNull()) return match;
331 
332  for (GsfElectronCollection::const_iterator it = eleCol->begin(); it!=eleCol->end(); ++it) {
333  //match electron to supercluster
334  if (it->superCluster()!=sc) continue;
335 
336  //check expected inner hits
337  if (it->gsfTrack()->trackerExpectedHitsInner().numberOfHits()>0) continue;
338 
339  //check if electron is matching to a conversion
340  if (hasMatchedConversion(*it,convCol,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
341 
342 
343  match = GsfElectronRef(eleCol,it-eleCol->begin());
344  }
345 
346  return match;
347 
348 
349 }
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:97
tuple convTracks
Definition: convBrem_cff.py:35
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:124
static HepMC::IO_HEPEVT conv
double y() const
y coordinate
Definition: Vertex.h:96
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:60
Definition: DDAxes.h:10
static bool hasMatchedPromptElectron(const reco::SuperClusterRef &sc, const edm::Handle< reco::GsfElectronCollection > &eleCol, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
bool isNull() const
Checks for null.
Definition: RefToBase.h:270
ProductID id() const
Definition: RefToBase.h:220
const Point & position() const
position
Definition: Vertex.h:92
math::XYZVectorF refittedPairMomentum() const
Conversion tracks momentum from the tracks refitted with vertex constraint.
Definition: Conversion.cc:248
TrackRef closestCtfTrackRef() const
Definition: GsfElectron.h:185
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
bool isNull() const
Checks for null.
Definition: Ref.h:247
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
double chi2() const
chi-squares
Definition: Vertex.h:81
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
size_t key() const
Definition: RefToBase.h:228
static bool hasMatchedConversion(const reco::GsfElectron &ele, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
static bool isGoodConversion(const reco::Conversion &conv, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)
edm::Ref< GsfElectronCollection > GsfElectronRef
reference to an object in a collection of GsfElectron objects
double ndof() const
Definition: Vertex.h:88
static reco::GsfElectronRef matchedPromptElectron(const reco::SuperClusterRef &sc, const edm::Handle< reco::GsfElectronCollection > &eleCol, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:12
double x() const
x coordinate
Definition: Vertex.h:94
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
key_type key() const
Accessor for product key.
Definition: Ref.h:266
const std::vector< uint8_t > & nHitsBeforeVtx() const
Vector of the number of hits before the vertex along each track trajector.
Definition: Conversion.h:161
edm::Ref< ConversionCollection > ConversionRef
reference to an object in a collection of Conversion objects
Definition: ConversionFwd.h:15
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
ProductIndex id() const
Definition: ProductID.h:38
static reco::ConversionRef matchedConversion(const reco::GsfElectron &ele, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
std::vector< edm::RefToBase< reco::Track > > const & tracks() const
vector of track to base references
Definition: Conversion.cc:176
virtual GsfTrackRef gsfTrack() const
reference to a GsfTrack
Definition: GsfElectron.h:169
std::string match(BranchDescription const &a, BranchDescription const &b, std::string const &fileName)