00001 #include "PhysicsTools/RecoAlgos/interface/TrackWithVertexSelector.h"
00002
00003
00004
00005
00006 TrackWithVertexSelector::TrackWithVertexSelector(const edm::ParameterSet& iConfig) :
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 }
00023
00024 TrackWithVertexSelector::~TrackWithVertexSelector() { }
00025
00026 bool TrackWithVertexSelector::testTrack(const reco::Track &t) const {
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 }
00042
00043 bool TrackWithVertexSelector::testVertices(const reco::Track &t, const reco::VertexCollection &vtxs) const {
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 }
00058
00059 bool TrackWithVertexSelector::operator()(const reco::Track &t, const edm::Event &evt) const {
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 }
00067
00068 bool TrackWithVertexSelector::operator()(const reco::Track &t, const reco::VertexCollection &vtxs) const {
00069 using std::abs;
00070 if (!testTrack(t)) return false;
00071 if (nVertices_ == 0) return true;
00072 return testVertices(t, vtxs);
00073 }
00074
00075 bool TrackWithVertexSelector::testPoint(const Point &point, const Point &vtx) const {
00076 using std::abs;
00077 math::XYZVector d = point - vtx;
00078 return ((abs(d.z()) < zetaVtx_) && (abs(d.Rho()) < rhoVtx_));
00079 }