CMS 3D CMS Logo

ConversionTools.cc
Go to the documentation of this file.
8 
9 #include <TMath.h>
10 
11 using namespace edm;
12 using namespace reco;
13 
14 //--------------------------------------------------------------------------------------------------
16  const math::XYZPoint &beamspot,
17  float lxyMin,
18  float probMin,
19  unsigned int nHitsBeforeVtxMax) {
20  //Check if a given conversion candidate passes the conversion selection cuts
21 
22  const reco::Vertex &vtx = conv.conversionVertex();
23 
24  //vertex validity
25  if (!vtx.isValid())
26  return false;
27 
28  //fit probability
29  if (TMath::Prob(vtx.chi2(), vtx.ndof()) < probMin)
30  return false;
31 
32  //compute transverse decay length
33  math::XYZVector mom(conv.refittedPairMomentum());
34  double dbsx = vtx.x() - beamspot.x();
35  double dbsy = vtx.y() - beamspot.y();
36  double lxy = (mom.x() * dbsx + mom.y() * dbsy) / mom.rho();
37 
38  //transverse decay length
39  if (lxy < lxyMin)
40  return false;
41 
42  //loop through daughters to check nhitsbeforevtx
43  for (std::vector<uint8_t>::const_iterator it = conv.nHitsBeforeVtx().begin(); it != conv.nHitsBeforeVtx().end();
44  ++it) {
45  if ((*it) > nHitsBeforeVtxMax)
46  return false;
47  }
48 
49  return true;
50 }
51 
52 //--------------------------------------------------------------------------------------------------
54  const reco::Conversion &conv,
55  bool allowCkfMatch,
56  bool allowAmbiguousGsfMatch) {
57  //check if a given GsfElectron matches a given conversion (no quality cuts applied)
58  //matching is always attempted through the gsf track ref, and optionally attempted through the
59  //closest ctf track ref
60 
61  const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
62  for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it = convTracks.begin(); it != convTracks.end();
63  ++it) {
64  if (ele.reco::GsfElectron::gsfTrack().isNonnull() && ele.reco::GsfElectron::gsfTrack().id() == it->id() &&
65  ele.reco::GsfElectron::gsfTrack().key() == it->key())
66  return true;
67  else if (allowCkfMatch && ele.reco::GsfElectron::closestCtfTrackRef().isNonnull() &&
68  ele.reco::GsfElectron::closestCtfTrackRef().id() == it->id() &&
69  ele.reco::GsfElectron::closestCtfTrackRef().key() == it->key())
70  return true;
71  if (allowAmbiguousGsfMatch) {
72  for (auto const &tk : ele.ambiguousGsfTracks()) {
73  if (tk.isNonnull() && tk.id() == it->id() && tk.key() == it->key())
74  return true;
75  }
76  }
77  }
78 
79  return false;
80 }
81 
82 //--------------------------------------------------------------------------------------------------
84  const reco::Conversion &conv,
85  bool allowCkfMatch) {
86  //check if a given GsfElectronCore matches a given conversion (no quality cuts applied)
87  //matching is always attempted through the gsf track ref, and optionally attempted through the
88  //closest ctf track ref
89 
90  for (const auto &trkRef : conv.tracks()) {
91  if (eleCore.gsfTrack().isNonnull() && eleCore.gsfTrack().id() == trkRef.id() &&
92  eleCore.gsfTrack().key() == trkRef.key())
93  return true;
94  else if (allowCkfMatch && eleCore.ctfTrack().isNonnull() && eleCore.ctfTrack().id() == trkRef.id() &&
95  eleCore.ctfTrack().key() == trkRef.key())
96  return true;
97  }
98 
99  return false;
100 }
101 
102 //--------------------------------------------------------------------------------------------------
104  const reco::SuperCluster &sc, const reco::Conversion &conv, float dRMax, float dEtaMax, float dPhiMax) {
105  //check if a given SuperCluster matches a given conversion (no quality cuts applied)
106  //matching is geometric between conversion momentum and vector joining conversion vertex
107  //to supercluster position
108 
109  math::XYZVector mom(conv.refittedPairMomentum());
110 
111  const math::XYZPoint &scpos(sc.position());
112  math::XYZPoint cvtx(conv.conversionVertex().position());
113 
114  math::XYZVector cscvector = scpos - cvtx;
115  float dR = reco::deltaR(mom, cscvector);
116  float dEta = mom.eta() - cscvector.eta();
117  float dPhi = reco::deltaPhi(mom.phi(), cscvector.phi());
118 
119  if (dR > dRMax)
120  return false;
121  if (dEta > dEtaMax)
122  return false;
123  if (dPhi > dPhiMax)
124  return false;
125 
126  return true;
127 }
128 
129 //--------------------------------------------------------------------------------------------------
131  //check if given track matches given conversion (matching by ref)
132 
133  if (trk.isNull())
134  return false;
135 
136  const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
137  for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it = convTracks.begin(); it != convTracks.end();
138  ++it) {
139  if (trk.id() == it->id() && trk.key() == it->key())
140  return true;
141  }
142 
143  return false;
144 }
145 
146 //--------------------------------------------------------------------------------------------------
148  //check if given track matches given conversion (matching by ref)
149 
150  if (trk.isNull())
151  return false;
152 
153  const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
154  for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it = convTracks.begin(); it != convTracks.end();
155  ++it) {
156  if (trk.id() == it->id() && trk.key() == it->key())
157  return true;
158  }
159 
160  return false;
161 }
162 
163 //--------------------------------------------------------------------------------------------------
165  //check if given track matches given conversion (matching by ref)
166 
167  if (trk.isNull())
168  return false;
169 
170  const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
171  for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it = convTracks.begin(); it != convTracks.end();
172  ++it) {
173  if (trk.id() == it->id() && trk.key() == it->key())
174  return true;
175  }
176 
177  return false;
178 }
179 
180 //--------------------------------------------------------------------------------------------------
182  const reco::ConversionCollection &convCol,
183  const math::XYZPoint &beamspot,
184  bool allowCkfMatch,
185  float lxyMin,
186  float probMin,
187  unsigned int nHitsBeforeVtxMax) {
188  //check if a given electron candidate matches to at least one conversion candidate in the
189  //collection which also passes the selection cuts, optionally match with the closestckf track in
190  //in addition to just the gsf track (enabled in default arguments)
191 
192  for (auto const &it : convCol) {
193  if (!matchesConversion(ele, it, allowCkfMatch))
194  continue;
195  if (!isGoodConversion(it, beamspot, lxyMin, probMin, nHitsBeforeVtxMax))
196  continue;
197 
198  return true;
199  }
200 
201  return false;
202 }
203 
204 //--------------------------------------------------------------------------------------------------
206  const reco::ConversionCollection &convCol,
207  const math::XYZPoint &beamspot,
208  float lxyMin,
209  float probMin,
210  unsigned int nHitsBeforeVtxMax) {
211  //check if a given track matches to at least one conversion candidate in the
212  //collection which also passes the selection cuts
213 
214  if (trk.isNull())
215  return false;
216 
217  for (auto const &it : convCol) {
218  if (!matchesConversion(trk, it))
219  continue;
220  if (!isGoodConversion(it, beamspot, lxyMin, probMin, nHitsBeforeVtxMax))
221  continue;
222 
223  return true;
224  }
225 
226  return false;
227 }
228 
229 //--------------------------------------------------------------------------------------------------
231  const reco::ConversionCollection &convCol,
232  const math::XYZPoint &beamspot,
233  float dRMax,
234  float dEtaMax,
235  float dPhiMax,
236  float lxyMin,
237  float probMin,
238  unsigned int nHitsBeforeVtxMax) {
239  //check if a given SuperCluster matches to at least one conversion candidate in the
240  //collection which also passes the selection cuts
241 
242  for (auto const &it : convCol) {
243  if (!matchesConversion(sc, it))
244  continue;
245  if (!isGoodConversion(it, beamspot, lxyMin, probMin, nHitsBeforeVtxMax))
246  continue;
247 
248  return true;
249  }
250 
251  return false;
252 }
253 
254 //--------------------------------------------------------------------------------------------------
256  const reco::ConversionCollection &convCol,
257  const math::XYZPoint &beamspot,
258  bool allowCkfMatch,
259  float lxyMin,
260  float probMin,
261  unsigned int nHitsBeforeVtxMax) {
262  //check if a given electron candidate matches to at least one conversion candidate in the
263  //collection which also passes the selection cuts, optionally match with the closestckf track in
264  //in addition to just the gsf track (enabled in default arguments)
265  //If multiple conversions are found, returned reference corresponds to minimum
266  //conversion radius
267 
268  reco::Conversion const *match = nullptr;
269 
270  double minRho = 999.;
271  for (auto const &it : convCol) {
272  float rho = it.conversionVertex().position().rho();
273  if (rho > minRho)
274  continue;
275  if (!matchesConversion(ele, it, allowCkfMatch))
276  continue;
277  if (!isGoodConversion(it, beamspot, lxyMin, probMin, nHitsBeforeVtxMax))
278  continue;
279 
280  minRho = rho;
281  match = &it;
282  }
283 
284  return match;
285 }
286 
287 //--------------------------------------------------------------------------------------------------
289  const reco::ConversionCollection &convCol,
290  const math::XYZPoint &beamspot,
291  bool allowCkfMatch,
292  float lxyMin,
293  float probMin,
294  unsigned int nHitsBeforeVtxMax) {
295  //check if a given electron candidate matches to at least one conversion candidate in the
296  //collection which also passes the selection cuts, optionally match with the closestckf track in
297  //in addition to just the gsf track (enabled in default arguments)
298  //If multiple conversions are found, returned reference corresponds to minimum
299  //conversion radius
300 
301  reco::Conversion const *match = nullptr;
302 
303  double minRho = 999.;
304  for (auto const &it : convCol) {
305  float rho = it.conversionVertex().position().rho();
306  if (rho > minRho)
307  continue;
308  if (!matchesConversion(eleCore, it, allowCkfMatch))
309  continue;
310  if (!isGoodConversion(it, beamspot, lxyMin, probMin, nHitsBeforeVtxMax))
311  continue;
312 
313  minRho = rho;
314  match = &it;
315  }
316 
317  return match;
318 }
319 
320 //--------------------------------------------------------------------------------------------------
322  const reco::ConversionCollection &convCol,
323  const math::XYZPoint &beamspot,
324  float lxyMin,
325  float probMin,
326  unsigned int nHitsBeforeVtxMax) {
327  //check if a given track matches to at least one conversion candidate in the
328  //collection which also passes the selection cuts
329  //If multiple conversions are found, returned reference corresponds to minimum
330  //conversion radius
331 
332  reco::Conversion const *match = nullptr;
333 
334  if (trk.isNull())
335  return match;
336 
337  double minRho = 999.;
338  for (auto const &it : convCol) {
339  float rho = it.conversionVertex().position().rho();
340  if (rho > minRho)
341  continue;
342  if (!matchesConversion(trk, it))
343  continue;
344  if (!isGoodConversion(it, beamspot, lxyMin, probMin, nHitsBeforeVtxMax))
345  continue;
346 
347  minRho = rho;
348  match = &it;
349  }
350 
351  return match;
352 }
353 
354 //--------------------------------------------------------------------------------------------------
356  const reco::ConversionCollection &convCol,
357  const math::XYZPoint &beamspot,
358  float dRMax,
359  float dEtaMax,
360  float dPhiMax,
361  float lxyMin,
362  float probMin,
363  unsigned int nHitsBeforeVtxMax) {
364  //check if a given SuperCluster matches to at least one conversion candidate in the
365  //collection which also passes the selection cuts
366  //If multiple conversions are found, returned reference corresponds to minimum
367  //conversion radius
368 
369  reco::Conversion const *match = nullptr;
370 
371  double minRho = 999.;
372  for (auto const &it : convCol) {
373  float rho = it.conversionVertex().position().rho();
374  if (rho > minRho)
375  continue;
376  if (!matchesConversion(sc, it, dRMax, dEtaMax, dPhiMax))
377  continue;
378  if (!isGoodConversion(it, beamspot, lxyMin, probMin, nHitsBeforeVtxMax))
379  continue;
380 
381  minRho = rho;
382  match = &it;
383  }
384 
385  return match;
386 }
387 
388 //--------------------------------------------------------------------------------------------------
390  const reco::GsfElectronCollection &eleCol,
391  const reco::ConversionCollection &convCol,
392  const math::XYZPoint &beamspot,
393  bool allowCkfMatch,
394  float lxyMin,
395  float probMin,
396  unsigned int nHitsBeforeVtxMax) {
397  return !(matchedPromptElectron(sc, eleCol, convCol, beamspot, allowCkfMatch, lxyMin, probMin, nHitsBeforeVtxMax) ==
398  nullptr);
399 }
400 
401 //--------------------------------------------------------------------------------------------------
403  const reco::GsfElectronCollection &eleCol,
404  const reco::ConversionCollection &convCol,
405  const math::XYZPoint &beamspot,
406  bool allowCkfMatch,
407  float lxyMin,
408  float probMin,
409  unsigned int nHitsBeforeVtxMax) {
410  //check if a given SuperCluster matches to at least one GsfElectron having zero expected inner hits
411  //and not matching any conversion in the collection passing the quality cuts
412 
413  reco::GsfElectron const *match = nullptr;
414 
415  if (sc.isNull())
416  return match;
417 
418  for (auto const &it : eleCol) {
419  //match electron to supercluster
420  if (it.superCluster() != sc)
421  continue;
422 
423  //check expected inner hits
424  if (it.gsfTrack()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) > 0)
425  continue;
426 
427  //check if electron is matching to a conversion
428  if (hasMatchedConversion(it, convCol, beamspot, allowCkfMatch, lxyMin, probMin, nHitsBeforeVtxMax))
429  continue;
430 
431  match = &it;
432  }
433 
434  return match;
435 }
436 
437 //--------------------------------------------------------------------------------------------------
439  if (conv != nullptr) {
440  const reco::Vertex &vtx = conv->conversionVertex();
441  if (vtx.isValid()) {
442  return TMath::Prob(vtx.chi2(), vtx.ndof());
443  }
444  }
445  return -1;
446 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:154
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
TrackRef ctfTrack() const
static bool hasMatchedPromptElectron(const reco::SuperClusterRef &sc, const reco::GsfElectronCollection &eleCol, const reco::ConversionCollection &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
static reco::Conversion const * matchedConversion(const reco::GsfElectron &ele, const reco::ConversionCollection &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
key_type key() const
Accessor for product key.
Definition: Ref.h:250
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
ProductID id() const
Definition: RefToBase.h:216
static float getVtxFitProb(const reco::Conversion *conv)
static bool hasMatchedConversion(const reco::GsfElectron &ele, const reco::ConversionCollection &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
auto const & ambiguousGsfTracks() const
Definition: GsfElectron.h:767
bool isNull() const
Checks for null.
Definition: RefToBase.h:297
bool isNull() const
Checks for null.
Definition: Ref.h:235
static bool isGoodConversion(const reco::Conversion &conv, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
ProductIndex id() const
Definition: ProductID.h:35
size_t key() const
Definition: RefToBase.h:221
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
EPOS::IO_EPOS conv
const GsfTrackRef & gsfTrack() const
fixed size matrix
HLT enums.
static reco::GsfElectron const * matchedPromptElectron(const reco::SuperClusterRef &sc, const reco::GsfElectronCollection &eleCol, const reco::ConversionCollection &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
std::string match(BranchDescription const &a, BranchDescription const &b, std::string const &fileName)