CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TrackWithVertexSelector.cc
Go to the documentation of this file.
3 //
4 // constructors and destructor
5 //
6 
7 namespace {
8  constexpr float fakeBeamSpotTimeWidth = 0.175f;
9 }
10 
12  : numberOfValidHits_(iConfig.getParameter<uint32_t>("numberOfValidHits")),
13  numberOfValidPixelHits_(iConfig.getParameter<uint32_t>("numberOfValidPixelHits")),
14  numberOfLostHits_(iConfig.getParameter<uint32_t>("numberOfLostHits")),
15  normalizedChi2_(iConfig.getParameter<double>("normalizedChi2")),
16  ptMin_(iConfig.getParameter<double>("ptMin")),
17  ptMax_(iConfig.getParameter<double>("ptMax")),
18  etaMin_(iConfig.getParameter<double>("etaMin")),
19  etaMax_(iConfig.getParameter<double>("etaMax")),
20  dzMax_(iConfig.getParameter<double>("dzMax")),
21  d0Max_(iConfig.getParameter<double>("d0Max")),
22  ptErrorCut_(iConfig.getParameter<double>("ptErrorCut")),
23  quality_(iConfig.getParameter<std::string>("quality")),
24  nVertices_(iConfig.getParameter<bool>("useVtx") ? iConfig.getParameter<uint32_t>("nVertices") : 0),
25  vertexToken_(iC.consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexTag"))),
26  timesToken_(iC.consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("timesTag"))),
27  timeResosToken_(iC.consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("timeResosTag"))),
28  vtxFallback_(iConfig.getParameter<bool>("vtxFallback")),
29  zetaVtx_(iConfig.getParameter<double>("zetaVtx")),
30  rhoVtx_(iConfig.getParameter<double>("rhoVtx")),
31  nSigmaDtVertex_(iConfig.getParameter<double>("nSigmaDtVertex")) {}
32 
34 
37  event.getByToken(vertexToken_, hVtx);
38  vcoll_ = hVtx.product();
39 
41  event.getByToken(timesToken_, hTimes);
42  timescoll_ = hTimes.isValid() ? hTimes.product() : nullptr;
43 
45  event.getByToken(timeResosToken_, hTimeResos);
46  timeresoscoll_ = hTimeResos.isValid() ? hTimeResos.product() : nullptr;
47 }
48 
50  using std::abs;
52  (static_cast<unsigned int>(t.hitPattern().numberOfValidPixelHits()) >= numberOfValidPixelHits_) &&
54  (t.ptError() / t.pt() * std::max(1., t.normalizedChi2()) <= ptErrorCut_) &&
55  (t.quality(t.qualityByName(quality_))) && (t.pt() >= ptMin_) && (t.pt() <= ptMax_) && (abs(t.eta()) <= etaMax_) &&
56  (abs(t.eta()) >= etaMin_) && (abs(t.dz()) <= dzMax_) && (abs(t.d0()) <= d0Max_)) {
57  return true;
58  }
59  return false;
60 }
61 
62 bool TrackWithVertexSelector::testTrack(const reco::TrackRef &tref) const { return testTrack(*tref); }
63 
65  bool ok = false;
66  if (!vtxs.empty()) {
67  unsigned int tested = 1;
68  for (reco::VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end(); it != ed; ++it) {
69  if ((std::abs(t.dxy(it->position())) < rhoVtx_) && (std::abs(t.dz(it->position())) < zetaVtx_)) {
70  ok = true;
71  break;
72  }
73  if (tested++ >= nVertices_)
74  break;
75  }
76  } else if (vtxFallback_) {
77  return ((std::abs(t.vertex().z()) < 15.9) && (t.vertex().Rho() < 0.2));
78  }
79  return ok;
80 }
81 
83  const auto &t = *tref;
84  const bool timeAvailable = timescoll_ != nullptr && timeresoscoll_ != nullptr;
85  bool ok = false;
86  if (!vtxs.empty()) {
87  unsigned int tested = 1;
88  for (reco::VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end(); it != ed; ++it) {
89  const bool useTime = timeAvailable && it->t() != 0.;
90  float time = useTime ? (*timescoll_)[tref] : -1.f;
91  float timeReso = useTime ? (*timeresoscoll_)[tref] : -1.f;
92  timeReso = (timeReso > 1e-6 ? timeReso : fakeBeamSpotTimeWidth);
93 
94  if (edm::isNotFinite(time)) {
95  time = 0.0;
96  timeReso = 2.0 * fakeBeamSpotTimeWidth;
97  }
98 
99  const double vtxSigmaT2 = it->tError() * it->tError();
100  const double vtxTrackErr = std::sqrt(vtxSigmaT2 + timeReso * timeReso);
101 
102  if ((std::abs(t.dxy(it->position())) < rhoVtx_) && (std::abs(t.dz(it->position())) < zetaVtx_) &&
103  (!useTime || (std::abs(time - it->t()) / vtxTrackErr < nSigmaDtVertex_))) {
104  ok = true;
105  break;
106  }
107  if (tested++ >= nVertices_)
108  break;
109  }
110  } else if (vtxFallback_) {
111  return ((std::abs(t.vertex().z()) < 15.9) && (t.vertex().Rho() < 0.2));
112  }
113  return ok;
114 }
115 
117  if (!testTrack(t))
118  return false;
119  if (nVertices_ == 0)
120  return true;
121  return testVertices(t, *vcoll_);
122 }
123 
125  if (!testTrack(tref))
126  return false;
127  if (nVertices_ == 0)
128  return true;
129  return testVertices(tref, *vcoll_);
130 }
bool testTrack(const reco::Track &t) const
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:611
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:593
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
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:801
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
etaMax_(conf.getParameter< double >("etaMax"))
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
std::vector< Vertex > VertexCollection
Definition: Vertex.h:12
void init(const edm::Event &event, const edm::EventSetup &)
bool operator()(const reco::Track &t) const
ptMin_(conf.getParameter< double >("ptMin"))
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:676
edm::EDGetTokenT< edm::ValueMap< float > > timesToken_
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
T sqrt(T t)
Definition: SSEVec.h:19
edm::EDGetTokenT< edm::ValueMap< float > > timeResosToken_
double pt() const
track transverse momentum
Definition: TrackBase.h:637
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:754
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:798
bool isValid() const
Definition: HandleBase.h:70
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:622
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:126
T const * product() const
Definition: Handle.h:70
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:504
edm::ValueMap< float > const * timeresoscoll_
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:552
int numberOfValidPixelHits() const
Definition: HitPattern.h:831
TrackWithVertexSelector(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
edm::ValueMap< float > const * timescoll_
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:608
etaMin_(conf.getParameter< double >("etaMin"))