#include <PhysicsTools/RecoAlgos/interface/TrackWithVertexSelector.h>
Public Member Functions | |
bool | operator() (const reco::Track &t, const reco::VertexCollection &vtxs) const |
bool | operator() (const reco::Track &t, const edm::Event &iEvent) const |
bool | testTrack (const reco::Track &t) const |
bool | testVertices (const reco::Track &t, const reco::VertexCollection &vtxs) const |
TrackWithVertexSelector (const edm::ParameterSet &) | |
~TrackWithVertexSelector () | |
Private Types | |
typedef math::XYZPoint | Point |
Private Member Functions | |
bool | testPoint (const Point &point, const Point &vtx) const |
Private Attributes | |
double | d0Max_ |
double | dzMax_ |
double | etaMax_ |
double | etaMin_ |
double | normalizedChi2_ |
uint32_t | numberOfLostHits_ |
uint32_t | numberOfValidHits_ |
uint32_t | numberOfValidPixelHits_ |
uint32_t | nVertices_ |
double | ptMax_ |
double | ptMin_ |
double | rhoVtx_ |
edm::InputTag | vertexTag_ |
bool | vtxFallback_ |
double | zetaVtx_ |
Definition at line 16 of file TrackWithVertexSelector.h.
typedef math::XYZPoint TrackWithVertexSelector::Point [private] |
Definition at line 37 of file TrackWithVertexSelector.h.
TrackWithVertexSelector::TrackWithVertexSelector | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 6 of file TrackWithVertexSelector.cc.
00006 : 00007 numberOfValidHits_(iConfig.getParameter<uint32_t>("numberOfValidHits")), 00008 numberOfValidPixelHits_(iConfig.getParameter<uint32_t>("numberOfValidPixelHits")), 00009 numberOfLostHits_(iConfig.getParameter<uint32_t>("numberOfLostHits")), 00010 normalizedChi2_(iConfig.getParameter<double>("normalizedChi2")), 00011 ptMin_(iConfig.getParameter<double>("ptMin")), 00012 ptMax_(iConfig.getParameter<double>("ptMax")), 00013 etaMin_(iConfig.getParameter<double>("etaMin")), 00014 etaMax_(iConfig.getParameter<double>("etaMax")), 00015 dzMax_(iConfig.getParameter<double>("dzMax")), 00016 d0Max_(iConfig.getParameter<double>("d0Max")), 00017 nVertices_(iConfig.getParameter<bool>("useVtx") ? iConfig.getParameter<uint32_t>("nVertices") : 0), 00018 vertexTag_(iConfig.getParameter<edm::InputTag>("vertexTag")), 00019 vtxFallback_(iConfig.getParameter<bool>("vtxFallback")), 00020 zetaVtx_(iConfig.getParameter<double>("zetaVtx")), 00021 rhoVtx_(iConfig.getParameter<double>("rhoVtx")) { 00022 }
TrackWithVertexSelector::~TrackWithVertexSelector | ( | ) |
bool TrackWithVertexSelector::operator() | ( | const reco::Track & | t, | |
const reco::VertexCollection & | vtxs | |||
) | const |
Definition at line 68 of file TrackWithVertexSelector.cc.
References funct::abs(), nVertices_, testTrack(), and testVertices().
00068 { 00069 using std::abs; 00070 if (!testTrack(t)) return false; 00071 if (nVertices_ == 0) return true; 00072 return testVertices(t, vtxs); 00073 }
bool TrackWithVertexSelector::operator() | ( | const reco::Track & | t, | |
const edm::Event & | iEvent | |||
) | const |
Definition at line 59 of file TrackWithVertexSelector.cc.
References funct::abs(), edm::Event::getByLabel(), nVertices_, testTrack(), testVertices(), and vertexTag_.
00059 { 00060 using std::abs; 00061 if (!testTrack(t)) return false; 00062 if (nVertices_ == 0) return true; 00063 edm::Handle<reco::VertexCollection> hVtx; 00064 evt.getByLabel(vertexTag_, hVtx); 00065 return testVertices(t, *hVtx); 00066 }
Definition at line 75 of file TrackWithVertexSelector.cc.
References funct::abs(), d, rhoVtx_, and zetaVtx_.
Referenced by testVertices().
00075 { 00076 using std::abs; 00077 math::XYZVector d = point - vtx; 00078 return ((abs(d.z()) < zetaVtx_) && (abs(d.Rho()) < rhoVtx_)); 00079 }
bool TrackWithVertexSelector::testTrack | ( | const reco::Track & | t | ) | const |
Definition at line 26 of file TrackWithVertexSelector.cc.
References funct::abs(), reco::TrackBase::d0(), d0Max_, reco::TrackBase::dz(), dzMax_, reco::TrackBase::eta(), etaMax_, etaMin_, reco::TrackBase::hitPattern(), reco::TrackBase::normalizedChi2(), normalizedChi2_, reco::TrackBase::numberOfLostHits(), numberOfLostHits_, reco::TrackBase::numberOfValidHits(), numberOfValidHits_, numberOfValidPixelHits_, reco::TrackBase::pt(), ptMax_, and ptMin_.
Referenced by operator()().
00026 { 00027 using std::abs; 00028 if ((t.numberOfValidHits() >= numberOfValidHits_) && 00029 (static_cast<unsigned int>(t.hitPattern().numberOfValidPixelHits()) >= numberOfValidPixelHits_) && 00030 (t.numberOfLostHits() <= numberOfLostHits_) && 00031 (t.normalizedChi2() <= normalizedChi2_) && 00032 (t.pt() >= ptMin_) && 00033 (t.pt() <= ptMax_) && 00034 (abs(t.eta()) <= etaMax_) && 00035 (abs(t.eta()) >= etaMin_) && 00036 (abs(t.dz()) <= dzMax_) && 00037 (abs(t.d0()) <= d0Max_) ) { 00038 return true; 00039 } 00040 return false; 00041 }
bool TrackWithVertexSelector::testVertices | ( | const reco::Track & | t, | |
const reco::VertexCollection & | vtxs | |||
) | const |
Definition at line 43 of file TrackWithVertexSelector.cc.
References funct::abs(), it, nVertices_, testPoint(), reco::TrackBase::vertex(), and vtxFallback_.
Referenced by operator()().
00043 { 00044 bool ok = false; 00045 const Point &pca = t.vertex(); 00046 if (vtxs.size() > 0) { 00047 unsigned int tested = 1; 00048 for (reco::VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end(); 00049 it != ed; ++it) { 00050 if (testPoint(pca, it->position())) { ok = true; break; } 00051 if (tested++ >= nVertices_) break; 00052 } 00053 } else if (vtxFallback_) { 00054 return ( (std::abs(pca.z()) < 15.9) && (pca.Rho() < 0.2) ); 00055 } 00056 return ok; 00057 }
double TrackWithVertexSelector::d0Max_ [private] |
double TrackWithVertexSelector::dzMax_ [private] |
double TrackWithVertexSelector::etaMax_ [private] |
double TrackWithVertexSelector::etaMin_ [private] |
double TrackWithVertexSelector::normalizedChi2_ [private] |
uint32_t TrackWithVertexSelector::numberOfLostHits_ [private] |
uint32_t TrackWithVertexSelector::numberOfValidHits_ [private] |
uint32_t TrackWithVertexSelector::numberOfValidPixelHits_ [private] |
uint32_t TrackWithVertexSelector::nVertices_ [private] |
Definition at line 32 of file TrackWithVertexSelector.h.
Referenced by operator()(), and testVertices().
double TrackWithVertexSelector::ptMax_ [private] |
double TrackWithVertexSelector::ptMin_ [private] |
double TrackWithVertexSelector::rhoVtx_ [private] |
bool TrackWithVertexSelector::vtxFallback_ [private] |
double TrackWithVertexSelector::zetaVtx_ [private] |