CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackWithVertexSelector.cc
Go to the documentation of this file.
2 //
3 // constructors and destructor
4 //
5 
7  numberOfValidHits_(iConfig.getParameter<uint32_t>("numberOfValidHits")),
8  numberOfValidPixelHits_(iConfig.getParameter<uint32_t>("numberOfValidPixelHits")),
9  numberOfLostHits_(iConfig.getParameter<uint32_t>("numberOfLostHits")),
10  normalizedChi2_(iConfig.getParameter<double>("normalizedChi2")),
11  ptMin_(iConfig.getParameter<double>("ptMin")),
12  ptMax_(iConfig.getParameter<double>("ptMax")),
13  etaMin_(iConfig.getParameter<double>("etaMin")),
14  etaMax_(iConfig.getParameter<double>("etaMax")),
15  dzMax_(iConfig.getParameter<double>("dzMax")),
16  d0Max_(iConfig.getParameter<double>("d0Max")),
17  ptErrorCut_(iConfig.getParameter<double>("ptErrorCut")),
18  quality_(iConfig.getParameter<std::string>("quality")),
19  nVertices_(iConfig.getParameter<bool>("useVtx") ? iConfig.getParameter<uint32_t>("nVertices") : 0),
20  vertexToken_(iC.consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexTag"))),
21  vtxFallback_(iConfig.getParameter<bool>("vtxFallback")),
22  zetaVtx_(iConfig.getParameter<double>("zetaVtx")),
23  rhoVtx_(iConfig.getParameter<double>("rhoVtx")) {
24 }
25 
27 
29  using std::abs;
31  (static_cast<unsigned int>(t.hitPattern().numberOfValidPixelHits()) >= numberOfValidPixelHits_) &&
34  (t.ptError()/t.pt()*std::max(1.,t.normalizedChi2()) <= ptErrorCut_) &&
35  (t.quality(t.qualityByName(quality_))) &&
36  (t.pt() >= ptMin_) &&
37  (t.pt() <= ptMax_) &&
38  (abs(t.eta()) <= etaMax_) &&
39  (abs(t.eta()) >= etaMin_) &&
40  (abs(t.dz()) <= dzMax_) &&
41  (abs(t.d0()) <= d0Max_) ) {
42  return true;
43  }
44  return false;
45 }
46 
48  bool ok = false;
49  if (vtxs.size() > 0) {
50  unsigned int tested = 1;
51  for (reco::VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end();
52  it != ed; ++it) {
53  if ((std::abs(t.dxy(it->position())) < rhoVtx_) &&
54  (std::abs(t.dz(it->position())) < zetaVtx_)) {
55  ok = true; break;
56  }
57  if (tested++ >= nVertices_) break;
58  }
59  } else if (vtxFallback_) {
60  return ( (std::abs(t.vertex().z()) < 15.9) && (t.vertex().Rho() < 0.2) );
61  }
62  return ok;
63 }
64 
66  if (!testTrack(t)) return false;
67  if (nVertices_ == 0) return true;
69  evt.getByToken(vertexToken_, hVtx);
70  return testVertices(t, *hVtx);
71 }
72 
74  if (!testTrack(t)) return false;
75  if (nVertices_ == 0) return true;
76  return testVertices(t, vtxs);
77 }
bool testTrack(const reco::Track &t) const
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:573
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:537
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
bool testVertices(const reco::Track &t, const reco::VertexCollection &vtxs) const
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:802
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:663
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:627
const T & max(const T &a, const T &b)
double pt() const
track transverse momentum
Definition: TrackBase.h:597
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:739
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:796
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:585
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:102
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:364
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:479
bool operator()(const reco::Track &t, const edm::Event &iEvent) const
int numberOfValidPixelHits() const
Definition: HitPattern.h:752
TrackWithVertexSelector(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:567