CMS 3D CMS Logo

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;
51  if ((t.numberOfValidHits() >= numberOfValidHits_) &&
52  (static_cast<unsigned int>(t.hitPattern().numberOfValidPixelHits()) >= numberOfValidPixelHits_) &&
53  (t.numberOfLostHits() <= numberOfLostHits_) && (t.normalizedChi2() <= normalizedChi2_) &&
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 }
electrons_cff.bool
bool
Definition: electrons_cff.py:372
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
edm::Handle::product
T const * product() const
Definition: Handle.h:70
edm::isNotFinite
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
edm
HLT enums.
Definition: AlignableModifier.h:19
reco::VertexCollection
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
TrackWithVertexSelector::testVertices
bool testVertices(const reco::Track &t, const reco::VertexCollection &vtxs) const
Definition: TrackWithVertexSelector.cc:64
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
TrackWithVertexSelector::init
void init(const edm::Event &event, const edm::EventSetup &)
Definition: TrackWithVertexSelector.h:24
TrackWithVertexSelector::vtxFallback_
bool vtxFallback_
Definition: TrackWithVertexSelector.h:57
convertSQLiteXML.ok
bool ok
Definition: convertSQLiteXML.py:98
edm::Handle< reco::VertexCollection >
TrackWithVertexSelector::ptMax_
double ptMax_
Definition: TrackWithVertexSelector.h:49
TrackWithVertexSelector::timeresoscoll_
const edm::ValueMap< float > * timeresoscoll_
Definition: TrackWithVertexSelector.h:62
edm::Ref< TrackCollection >
TrackWithVertexSelector::vcoll_
const reco::VertexCollection * vcoll_
Definition: TrackWithVertexSelector.h:60
TrackWithVertexSelector::numberOfValidHits_
uint32_t numberOfValidHits_
Definition: TrackWithVertexSelector.h:45
TrackWithVertexSelector::vertexToken_
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
Definition: TrackWithVertexSelector.h:55
TrackWithVertexSelector::timesToken_
edm::EDGetTokenT< edm::ValueMap< float > > timesToken_
Definition: TrackWithVertexSelector.h:56
TrackWithVertexSelector::~TrackWithVertexSelector
~TrackWithVertexSelector()
Definition: TrackWithVertexSelector.cc:33
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
TrackWithVertexSelector::rhoVtx_
double rhoVtx_
Definition: TrackWithVertexSelector.h:58
TrackWithVertexSelector::timeResosToken_
edm::EDGetTokenT< edm::ValueMap< float > > timeResosToken_
Definition: TrackWithVertexSelector.h:56
TrackWithVertexSelector::zetaVtx_
double zetaVtx_
Definition: TrackWithVertexSelector.h:58
reco::Track
Definition: Track.h:27
TrackWithVertexSelector::ptMin_
double ptMin_
Definition: TrackWithVertexSelector.h:49
OrderedSet.t
t
Definition: OrderedSet.py:90
TrackWithVertexSelector::ptErrorCut_
double ptErrorCut_
Definition: TrackWithVertexSelector.h:51
TrackWithVertexSelector::d0Max_
double d0Max_
Definition: TrackWithVertexSelector.h:50
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
edm::ParameterSet
Definition: ParameterSet.h:36
TrackWithVertexSelector::numberOfValidPixelHits_
uint32_t numberOfValidPixelHits_
Definition: TrackWithVertexSelector.h:46
TrackWithVertexSelector::testTrack
bool testTrack(const reco::Track &t) const
Definition: TrackWithVertexSelector.cc:49
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
TrackWithVertexSelector::nVertices_
uint32_t nVertices_
Definition: TrackWithVertexSelector.h:54
TrackWithVertexSelector::quality_
std::string quality_
Definition: TrackWithVertexSelector.h:52
TrackWithVertexSelector::operator()
bool operator()(const reco::Track &t) const
Definition: TrackWithVertexSelector.cc:116
TrackWithVertexSelector::nSigmaDtVertex_
double nSigmaDtVertex_
Definition: TrackWithVertexSelector.h:58
TrackWithVertexSelector::dzMax_
double dzMax_
Definition: TrackWithVertexSelector.h:50
std
Definition: JetResolutionObject.h:76
isFinite.h
TrackWithVertexSelector::etaMin_
double etaMin_
Definition: TrackWithVertexSelector.h:49
TrackWithVertexSelector::timescoll_
const edm::ValueMap< float > * timescoll_
Definition: TrackWithVertexSelector.h:61
TrackWithVertexSelector.h
TrackWithVertexSelector::numberOfLostHits_
uint32_t numberOfLostHits_
Definition: TrackWithVertexSelector.h:47
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
ntuplemaker.time
time
Definition: ntuplemaker.py:310
event
Definition: event.py:1
edm::Event
Definition: Event.h:73
TrackWithVertexSelector::etaMax_
double etaMax_
Definition: TrackWithVertexSelector.h:49
TrackWithVertexSelector::TrackWithVertexSelector
TrackWithVertexSelector(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
Definition: TrackWithVertexSelector.h:19
edm::ConsumesCollector
Definition: ConsumesCollector.h:39
TrackWithVertexSelector::normalizedChi2_
double normalizedChi2_
Definition: TrackWithVertexSelector.h:48
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37