#include <ConversionTools.h>
Public Member Functions | |
ConversionTools () | |
Static Public Member Functions | |
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 | hasMatchedConversion (const reco::TrackRef &trk, const edm::Handle< 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 edm::Handle< 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 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) |
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::ConversionRef | matchedConversion (const reco::SuperCluster &sc, const edm::Handle< 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::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) |
static reco::ConversionRef | matchedConversion (const reco::TrackRef &trk, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1) |
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) |
static bool | matchesConversion (const reco::GsfTrackRef &trk, const reco::Conversion &conv) |
static bool | matchesConversion (const reco::TrackRef &trk, const reco::Conversion &conv) |
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::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true) |
Definition at line 37 of file ConversionTools.h.
ConversionTools::ConversionTools | ( | ) | [inline] |
Definition at line 40 of file ConversionTools.h.
{}
bool ConversionTools::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] |
Definition at line 145 of file ConversionTools.cc.
Referenced by pat::PATElectronProducer::produce(), SoftPFElectronTagInfoProducer::produce(), and EgammaCutBasedEleId::TestWP().
{ //check if a given electron candidate matches to at least one conversion candidate in the //collection which also passes the selection cuts, optionally match with the closestckf track in //in addition to just the gsf track (enabled in default arguments) for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) { if (!matchesConversion(ele, *it, allowCkfMatch)) continue; if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue; return true; } return false; }
bool ConversionTools::hasMatchedConversion | ( | const reco::TrackRef & | trk, |
const edm::Handle< reco::ConversionCollection > & | convCol, | ||
const math::XYZPoint & | beamspot, | ||
float | lxyMin = 2.0 , |
||
float | probMin = 1e-6 , |
||
unsigned int | nHitsBeforeVtxMax = 1 |
||
) | [static] |
Definition at line 165 of file ConversionTools.cc.
References edm::Ref< C, T, F >::isNull().
{ //check if a given track matches to at least one conversion candidate in the //collection which also passes the selection cuts if (trk.isNull()) return false; for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) { if (!matchesConversion(trk, *it)) continue; if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue; return true; } return false; }
bool ConversionTools::hasMatchedConversion | ( | const reco::SuperCluster & | sc, |
const edm::Handle< 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 186 of file ConversionTools.cc.
{ //check if a given SuperCluster matches to at least one conversion candidate in the //collection which also passes the selection cuts for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) { if (!matchesConversion(sc, *it)) continue; if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue; return true; } return false; }
bool ConversionTools::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 |
||
) | [static] |
Definition at line 292 of file ConversionTools.cc.
References edm::Ref< C, T, F >::isNull().
Referenced by ggPFPhotons::PFElectronVeto().
{ //check if a given SuperCluster matches to at least one GsfElectron having zero expected inner hits //and not matching any conversion in the collection passing the quality cuts if (sc.isNull()) return false; for (GsfElectronCollection::const_iterator it = eleCol->begin(); it!=eleCol->end(); ++it) { //match electron to supercluster if (it->superCluster()!=sc) continue; //check expected inner hits if (it->gsfTrack()->trackerExpectedHitsInner().numberOfHits()>0) continue; //check if electron is matching to a conversion if (hasMatchedConversion(*it,convCol,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue; return true; } return false; }
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 19 of file ConversionTools.cc.
References reco::Vertex::chi2(), reco::Conversion::conversionVertex(), reco::Vertex::isValid(), reco::Vertex::ndof(), reco::Conversion::nHitsBeforeVtx(), reco::Conversion::refittedPairMomentum(), reco::Vertex::x(), and reco::Vertex::y().
{ //Check if a given conversion candidate passes the conversion selection cuts const reco::Vertex &vtx = conv.conversionVertex(); //vertex validity if (!vtx.isValid()) return false; //fit probability if (TMath::Prob( vtx.chi2(), vtx.ndof() )<probMin) return false; //compute transverse decay length math::XYZVector mom(conv.refittedPairMomentum()); double dbsx = vtx.x() - beamspot.x(); double dbsy = vtx.y() - beamspot.y(); double lxy = (mom.x()*dbsx + mom.y()*dbsy)/mom.rho(); //transverse decay length if ( lxy<lxyMin ) return false; //loop through daughters to check nhitsbeforevtx for (std::vector<uint8_t>::const_iterator it = conv.nHitsBeforeVtx().begin(); it!=conv.nHitsBeforeVtx().end(); ++it) { if ( (*it)>nHitsBeforeVtxMax ) return false; } return true; }
reco::ConversionRef ConversionTools::matchedConversion | ( | const reco::SuperCluster & | sc, |
const edm::Handle< 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 264 of file ConversionTools.cc.
References edm::match(), and rho.
{ //check if a given SuperCluster matches to at least one conversion candidate in the //collection which also passes the selection cuts //If multiple conversions are found, returned reference corresponds to minimum //conversion radius ConversionRef match; double minRho = 999.; for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) { float rho = it->conversionVertex().position().rho(); if (rho>minRho) continue; if (!matchesConversion(sc, *it, dRMax,dEtaMax,dPhiMax)) continue; if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue; minRho = rho; match = ConversionRef(convCol,it-convCol->begin()); } return match; }
reco::ConversionRef ConversionTools::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 |
||
) | [static] |
Definition at line 207 of file ConversionTools.cc.
References edm::match(), and rho.
Referenced by ElectronConversionRejectionValidator::analyze().
{ //check if a given electron candidate matches to at least one conversion candidate in the //collection which also passes the selection cuts, optionally match with the closestckf track in //in addition to just the gsf track (enabled in default arguments) //If multiple conversions are found, returned reference corresponds to minimum //conversion radius ConversionRef match; double minRho = 999.; for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) { float rho = it->conversionVertex().position().rho(); if (rho>minRho) continue; if (!matchesConversion(ele, *it, allowCkfMatch)) continue; if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue; minRho = rho; match = ConversionRef(convCol,it-convCol->begin()); } return match; }
reco::ConversionRef ConversionTools::matchedConversion | ( | const reco::TrackRef & | trk, |
const edm::Handle< reco::ConversionCollection > & | convCol, | ||
const math::XYZPoint & | beamspot, | ||
float | lxyMin = 2.0 , |
||
float | probMin = 1e-6 , |
||
unsigned int | nHitsBeforeVtxMax = 1 |
||
) | [static] |
Definition at line 235 of file ConversionTools.cc.
References edm::Ref< C, T, F >::isNull(), edm::match(), and rho.
{ //check if a given track matches to at least one conversion candidate in the //collection which also passes the selection cuts //If multiple conversions are found, returned reference corresponds to minimum //conversion radius ConversionRef match; if (trk.isNull()) return match; double minRho = 999.; for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) { float rho = it->conversionVertex().position().rho(); if (rho>minRho) continue; if (!matchesConversion(trk, *it)) continue; if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue; minRho = rho; match = ConversionRef(convCol,it-convCol->begin()); } return match; }
reco::GsfElectronRef ConversionTools::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 |
||
) | [static] |
Definition at line 322 of file ConversionTools.cc.
References edm::Ref< C, T, F >::isNull(), and edm::match().
{ //check if a given SuperCluster matches to at least one GsfElectron having zero expected inner hits //and not matching any conversion in the collection passing the quality cuts GsfElectronRef match; if (sc.isNull()) return match; for (GsfElectronCollection::const_iterator it = eleCol->begin(); it!=eleCol->end(); ++it) { //match electron to supercluster if (it->superCluster()!=sc) continue; //check expected inner hits if (it->gsfTrack()->trackerExpectedHitsInner().numberOfHits()>0) continue; //check if electron is matching to a conversion if (hasMatchedConversion(*it,convCol,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue; match = GsfElectronRef(eleCol,it-eleCol->begin()); } return match; }
bool ConversionTools::matchesConversion | ( | const reco::GsfTrackRef & | trk, |
const reco::Conversion & | conv | ||
) | [static] |
Definition at line 128 of file ConversionTools.cc.
References convBrem_cff::convTracks, edm::Ref< C, T, F >::id(), edm::ProductID::id(), edm::Ref< C, T, F >::isNull(), edm::Ref< C, T, F >::key(), and reco::Conversion::tracks().
{ //check if given track matches given conversion (matching by ref) if (trk.isNull()) return false; const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks(); for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) { if (trk.id()==it->id() && trk.key()==it->key()) return true; } return false; }
bool ConversionTools::matchesConversion | ( | const reco::TrackRef & | trk, |
const reco::Conversion & | conv | ||
) | [static] |
Definition at line 112 of file ConversionTools.cc.
References convBrem_cff::convTracks, edm::Ref< C, T, F >::id(), edm::ProductID::id(), edm::Ref< C, T, F >::isNull(), edm::Ref< C, T, F >::key(), and reco::Conversion::tracks().
{ //check if given track matches given conversion (matching by ref) if (trk.isNull()) return false; const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks(); for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) { if (trk.id()==it->id() && trk.key()==it->key()) return true; } return false; }
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 68 of file ConversionTools.cc.
References reco::Conversion::conversionVertex(), SiPixelRawToDigiRegional_cfi::deltaPhi, deltaR(), dPhi(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, reco::Vertex::position(), reco::CaloCluster::position(), and reco::Conversion::refittedPairMomentum().
{ //check if a given SuperCluster matches a given conversion (no quality cuts applied) //matching is geometric between conversion momentum and vector joining conversion vertex //to supercluster position math::XYZVector mom(conv.refittedPairMomentum()); math::XYZPoint scpos(sc.position()); math::XYZPoint cvtx(conv.conversionVertex().position()); math::XYZVector cscvector = scpos - cvtx; float dR = reco::deltaR(mom,cscvector); float dEta = mom.eta() - cscvector.eta(); float dPhi = reco::deltaPhi(mom.phi(),cscvector.phi()); if (dR>dRMax) return false; if (dEta>dEtaMax) return false; if (dPhi>dPhiMax) return false; return true; }
bool ConversionTools::matchesConversion | ( | const edm::RefToBase< reco::Track > & | trk, |
const reco::Conversion & | conv | ||
) | [static] |
Definition at line 96 of file ConversionTools.cc.
References convBrem_cff::convTracks, edm::RefToBase< T >::id(), edm::ProductID::id(), edm::RefToBase< T >::isNull(), edm::RefToBase< T >::key(), and reco::Conversion::tracks().
{ //check if given track matches given conversion (matching by ref) if (trk.isNull()) return false; const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks(); for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) { if (trk.id()==it->id() && trk.key()==it->key()) return true; } return false; }
bool ConversionTools::matchesConversion | ( | const reco::GsfElectron & | ele, |
const reco::Conversion & | conv, | ||
bool | allowCkfMatch = true |
||
) | [static] |
Definition at line 51 of file ConversionTools.cc.
References reco::GsfElectron::closestCtfTrackRef(), convBrem_cff::convTracks, reco::GsfElectron::gsfTrack(), edm::Ref< C, T, F >::id(), edm::ProductID::id(), edm::Ref< C, T, F >::isNonnull(), edm::Ref< C, T, F >::key(), and reco::Conversion::tracks().
Referenced by PF_PU_AssoMapAlgos::ComesFromConversion(), and pat::PATConversionProducer::produce().
{ //check if a given GsfElectron matches a given conversion (no quality cuts applied) //matching is always attempted through the gsf track ref, and optionally attempted through the //closest ctf track ref const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks(); for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) { if ( ele.gsfTrack().isNonnull() && ele.gsfTrack().id()==it->id() && ele.gsfTrack().key()==it->key()) return true; else if ( allowCkfMatch && ele.closestCtfTrackRef().isNonnull() && ele.closestCtfTrackRef().id()==it->id() && ele.closestCtfTrackRef().key()==it->key() ) return true; } return false; }