CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions
ConversionTools Class Reference

#include <ConversionTools.h>

Public Member Functions

 ConversionTools ()
 

Static Public Member Functions

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)
 
static bool hasMatchedConversion (const reco::TrackRef &trk, const reco::ConversionCollection &convCol, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)
 
static bool hasMatchedConversion (const reco::SuperCluster &sc, const reco::ConversionCollection &convCol, const math::XYZPoint &beamspot, float dRMax=0.1, float dEtaMax=999., float dPhiMax=999., float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)
 
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)
 
static bool isGoodConversion (const reco::Conversion &conv, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)
 
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)
 
static reco::Conversion const * matchedConversion (const reco::GsfElectronCore &eleCore, const reco::ConversionCollection &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
 
static reco::Conversion const * matchedConversion (const reco::TrackRef &trk, const reco::ConversionCollection &convCol, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)
 
static reco::Conversion const * matchedConversion (const reco::SuperCluster &sc, const reco::ConversionCollection &convCol, const math::XYZPoint &beamspot, float dRMax=0.1, float dEtaMax=999., float dPhiMax=999., float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)
 
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)
 
static bool matchesConversion (const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
 
static bool matchesConversion (const reco::GsfElectronCore &eleCore, const reco::Conversion &conv, bool allowCkfMatch=true)
 
static bool matchesConversion (const reco::SuperCluster &sc, const reco::Conversion &conv, float dRMax=0.1, float dEtaMax=999., float dPhiMax=999.)
 
static bool matchesConversion (const edm::RefToBase< reco::Track > &trk, const reco::Conversion &conv)
 
static bool matchesConversion (const reco::TrackRef &trk, const reco::Conversion &conv)
 
static bool matchesConversion (const reco::GsfTrackRef &trk, const reco::Conversion &conv)
 

Detailed Description

Definition at line 31 of file ConversionTools.h.

Constructor & Destructor Documentation

◆ ConversionTools()

ConversionTools::ConversionTools ( )
inline

Definition at line 33 of file ConversionTools.h.

33 {}

Member Function Documentation

◆ getVtxFitProb()

float ConversionTools::getVtxFitProb ( const reco::Conversion conv)
static

Definition at line 438 of file ConversionTools.cc.

References conv, and L1BJetProducer_cff::vtx.

Referenced by GsfElectronAlgo::createElectron(), and EG9X105XObjectUpdateModifier::modifyObject().

438  {
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 }
EPOS::IO_EPOS conv

◆ hasMatchedConversion() [1/3]

bool ConversionTools::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 
)
static

Definition at line 181 of file ConversionTools.cc.

References vertexSelectForHeavyFlavorDQM_cfi::probMin.

Referenced by LeptonSkimming::filter(), DQMExample_Step1::MediumEle(), GsfEleConversionVetoCut::operator()(), ElectronIdentifier::passID(), SoftPFElectronTagInfoProducer::produce(), pat::PATElectronProducer::produce(), EgammaCutBasedEleId::TestWP(), and GsfEleConversionVetoCut::value().

187  {
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 }
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
static bool isGoodConversion(const reco::Conversion &conv, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)

◆ hasMatchedConversion() [2/3]

bool ConversionTools::hasMatchedConversion ( const reco::TrackRef trk,
const reco::ConversionCollection convCol,
const math::XYZPoint beamspot,
float  lxyMin = 2.0,
float  probMin = 1e-6,
unsigned int  nHitsBeforeVtxMax = 1 
)
static

Definition at line 205 of file ConversionTools.cc.

References edm::Ref< C, T, F >::isNull(), and vertexSelectForHeavyFlavorDQM_cfi::probMin.

210  {
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 }
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
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)

◆ hasMatchedConversion() [3/3]

bool ConversionTools::hasMatchedConversion ( const reco::SuperCluster sc,
const reco::ConversionCollection convCol,
const math::XYZPoint beamspot,
float  dRMax = 0.1,
float  dEtaMax = 999.,
float  dPhiMax = 999.,
float  lxyMin = 2.0,
float  probMin = 1e-6,
unsigned int  nHitsBeforeVtxMax = 1 
)
static

Definition at line 230 of file ConversionTools.cc.

References vertexSelectForHeavyFlavorDQM_cfi::probMin.

238  {
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 }
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
static bool isGoodConversion(const reco::Conversion &conv, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)

◆ hasMatchedPromptElectron()

bool ConversionTools::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 
)
static

Definition at line 389 of file ConversionTools.cc.

References vertexSelectForHeavyFlavorDQM_cfi::probMin.

Referenced by pat::PATPhotonProducer::produce().

396  {
397  return !(matchedPromptElectron(sc, eleCol, convCol, beamspot, allowCkfMatch, lxyMin, probMin, nHitsBeforeVtxMax) ==
398  nullptr);
399 }
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)

◆ isGoodConversion()

bool ConversionTools::isGoodConversion ( const reco::Conversion conv,
const math::XYZPoint beamspot,
float  lxyMin = 2.0,
float  probMin = 1e-6,
unsigned int  nHitsBeforeVtxMax = 1 
)
static

Definition at line 15 of file ConversionTools.cc.

References conv, vertexSelectForHeavyFlavorDQM_cfi::probMin, and L1BJetProducer_cff::vtx.

19  {
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 }
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
EPOS::IO_EPOS conv

◆ matchedConversion() [1/4]

reco::Conversion const * ConversionTools::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 
)
static

Definition at line 255 of file ConversionTools.cc.

References edm::match(), and vertexSelectForHeavyFlavorDQM_cfi::probMin.

Referenced by ElectronConversionRejectionValidator::analyze(), GsfElectronAlgo::createElectron(), and EG9X105XObjectUpdateModifier::modifyObject().

261  {
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 }
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
static bool isGoodConversion(const reco::Conversion &conv, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10

◆ matchedConversion() [2/4]

reco::Conversion const * ConversionTools::matchedConversion ( const reco::GsfElectronCore eleCore,
const reco::ConversionCollection convCol,
const math::XYZPoint beamspot,
bool  allowCkfMatch = true,
float  lxyMin = 2.0,
float  probMin = 1e-6,
unsigned int  nHitsBeforeVtxMax = 0 
)
static

Definition at line 288 of file ConversionTools.cc.

References edm::match(), and vertexSelectForHeavyFlavorDQM_cfi::probMin.

294  {
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 }
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
static bool isGoodConversion(const reco::Conversion &conv, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10

◆ matchedConversion() [3/4]

reco::Conversion const * ConversionTools::matchedConversion ( const reco::TrackRef trk,
const reco::ConversionCollection convCol,
const math::XYZPoint beamspot,
float  lxyMin = 2.0,
float  probMin = 1e-6,
unsigned int  nHitsBeforeVtxMax = 1 
)
static

Definition at line 321 of file ConversionTools.cc.

References edm::Ref< C, T, F >::isNull(), edm::match(), and vertexSelectForHeavyFlavorDQM_cfi::probMin.

326  {
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 }
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
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)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10

◆ matchedConversion() [4/4]

reco::Conversion const * ConversionTools::matchedConversion ( const reco::SuperCluster sc,
const reco::ConversionCollection convCol,
const math::XYZPoint beamspot,
float  dRMax = 0.1,
float  dEtaMax = 999.,
float  dPhiMax = 999.,
float  lxyMin = 2.0,
float  probMin = 1e-6,
unsigned int  nHitsBeforeVtxMax = 1 
)
static

Definition at line 355 of file ConversionTools.cc.

References photonValidator_cfi::dEtaMax, HLT_2023v12_cff::dPhiMax, metBenchmark_cfi::dRMax, edm::match(), and vertexSelectForHeavyFlavorDQM_cfi::probMin.

363  {
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 }
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
static bool isGoodConversion(const reco::Conversion &conv, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10

◆ matchedPromptElectron()

reco::GsfElectron const * ConversionTools::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 
)
static

Definition at line 402 of file ConversionTools.cc.

References edm::Ref< C, T, F >::isNull(), edm::match(), reco::HitPattern::MISSING_INNER_HITS, and vertexSelectForHeavyFlavorDQM_cfi::probMin.

409  {
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 }
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)
bool isNull() const
Checks for null.
Definition: Ref.h:235
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10

◆ matchesConversion() [1/6]

bool ConversionTools::matchesConversion ( const reco::GsfElectron ele,
const reco::Conversion conv,
bool  allowCkfMatch = true,
bool  allowAmbiguousGsfMatch = false 
)
static

Definition at line 53 of file ConversionTools.cc.

References reco::GsfElectron::ambiguousGsfTracks(), conv, and trackerHitRTTI::vector.

Referenced by PF_PU_AssoMapAlgos::ComesFromConversion(), ReducedEGProducer::linkConversionsByTrackRef(), and pat::PATConversionProducer::produce().

56  {
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 }
auto const & ambiguousGsfTracks() const
Definition: GsfElectron.h:767
EPOS::IO_EPOS conv

◆ matchesConversion() [2/6]

bool ConversionTools::matchesConversion ( const reco::GsfElectronCore eleCore,
const reco::Conversion conv,
bool  allowCkfMatch = true 
)
static

Definition at line 83 of file ConversionTools.cc.

References conv, reco::GsfElectronCore::ctfTrack(), reco::GsfElectronCore::gsfTrack(), edm::ProductID::id(), edm::Ref< C, T, F >::id(), edm::Ref< C, T, F >::isNonnull(), and edm::Ref< C, T, F >::key().

85  {
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 }
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
TrackRef ctfTrack() const
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
key_type key() const
Accessor for product key.
Definition: Ref.h:250
ProductIndex id() const
Definition: ProductID.h:35
EPOS::IO_EPOS conv
const GsfTrackRef & gsfTrack() const

◆ matchesConversion() [3/6]

bool ConversionTools::matchesConversion ( const reco::SuperCluster sc,
const reco::Conversion conv,
float  dRMax = 0.1,
float  dEtaMax = 999.,
float  dPhiMax = 999. 
)
static

Definition at line 103 of file ConversionTools.cc.

References conv, reco::deltaPhi(), reco::deltaR(), HLT_2023v12_cff::dEta, photonValidator_cfi::dEtaMax, HLT_2023v12_cff::dPhiMax, HGC3DClusterGenMatchSelector_cfi::dR, metBenchmark_cfi::dRMax, and reco::CaloCluster::position().

104  {
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 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:153
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
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

◆ matchesConversion() [4/6]

bool ConversionTools::matchesConversion ( const edm::RefToBase< reco::Track > &  trk,
const reco::Conversion conv 
)
static

Definition at line 130 of file ConversionTools.cc.

References conv, edm::ProductID::id(), edm::RefToBase< T >::id(), edm::RefToBase< T >::isNull(), edm::RefToBase< T >::key(), and trackerHitRTTI::vector.

130  {
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 }
ProductID id() const
Definition: RefToBase.h:216
bool isNull() const
Checks for null.
Definition: RefToBase.h:297
ProductIndex id() const
Definition: ProductID.h:35
size_t key() const
Definition: RefToBase.h:221
EPOS::IO_EPOS conv

◆ matchesConversion() [5/6]

bool ConversionTools::matchesConversion ( const reco::TrackRef trk,
const reco::Conversion conv 
)
static

Definition at line 147 of file ConversionTools.cc.

References conv, edm::ProductID::id(), edm::Ref< C, T, F >::id(), edm::Ref< C, T, F >::isNull(), edm::Ref< C, T, F >::key(), and trackerHitRTTI::vector.

147  {
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 }
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
key_type key() const
Accessor for product key.
Definition: Ref.h:250
bool isNull() const
Checks for null.
Definition: Ref.h:235
ProductIndex id() const
Definition: ProductID.h:35
EPOS::IO_EPOS conv

◆ matchesConversion() [6/6]

bool ConversionTools::matchesConversion ( const reco::GsfTrackRef trk,
const reco::Conversion conv 
)
static

Definition at line 164 of file ConversionTools.cc.

References conv, edm::ProductID::id(), edm::Ref< C, T, F >::id(), edm::Ref< C, T, F >::isNull(), edm::Ref< C, T, F >::key(), and trackerHitRTTI::vector.

164  {
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 }
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
key_type key() const
Accessor for product key.
Definition: Ref.h:250
bool isNull() const
Checks for null.
Definition: Ref.h:235
ProductIndex id() const
Definition: ProductID.h:35
EPOS::IO_EPOS conv