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 
28 
31  event.getByToken(vertexToken_, hVtx);
32  vcoll_ = hVtx.product();
33 }
34 
35 
37  using std::abs;
39  (static_cast<unsigned int>(t.hitPattern().numberOfValidPixelHits()) >= numberOfValidPixelHits_) &&
42  (t.ptError()/t.pt()*std::max(1.,t.normalizedChi2()) <= ptErrorCut_) &&
43  (t.quality(t.qualityByName(quality_))) &&
44  (t.pt() >= ptMin_) &&
45  (t.pt() <= ptMax_) &&
46  (abs(t.eta()) <= etaMax_) &&
47  (abs(t.eta()) >= etaMin_) &&
48  (abs(t.dz()) <= dzMax_) &&
49  (abs(t.d0()) <= d0Max_) ) {
50  return true;
51  }
52  return false;
53 }
54 
56  bool ok = false;
57  if (vtxs.size() > 0) {
58  unsigned int tested = 1;
59  for (reco::VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end();
60  it != ed; ++it) {
61  if ((std::abs(t.dxy(it->position())) < rhoVtx_) &&
62  (std::abs(t.dz(it->position())) < zetaVtx_)) {
63  ok = true; break;
64  }
65  if (tested++ >= nVertices_) break;
66  }
67  } else if (vtxFallback_) {
68  return ( (std::abs(t.vertex().z()) < 15.9) && (t.vertex().Rho() < 0.2) );
69  }
70  return ok;
71 }
72 
74  if (!testTrack(t)) return false;
75  if (nVertices_ == 0) return true;
76  return testVertices(t, *vcoll_);
77 }
78 
bool testTrack(const reco::Track &t) const
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:592
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:556
bool testVertices(const reco::Track &t, const reco::VertexCollection &vtxs) const
reco::VertexCollection const * vcoll_
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:821
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
void init(const edm::Event &event, const edm::EventSetup &)
bool operator()(const reco::Track &t) const
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:682
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:646
double pt() const
track transverse momentum
Definition: TrackBase.h:616
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:758
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:815
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
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:604
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:125
T const * product() const
Definition: Handle.h:81
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:445
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:505
int numberOfValidPixelHits() const
Definition: HitPattern.h:838
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:586