CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Attributes
TrackWithVertexSelector Class Reference

#include <TrackWithVertexSelector.h>

Public Member Functions

void init (const edm::Event &event, const edm::EventSetup &)
 
void init (const edm::Event &event)
 
bool operator() (const reco::Track &t) const
 
bool operator() (const reco::Track &t, const edm::Event &iEvent)
 
bool operator() (const reco::TrackRef &t) const
 
bool operator() (const reco::TrackRef &t, const edm::Event &iEvent)
 
bool testTrack (const reco::Track &t) const
 
bool testTrack (const reco::TrackRef &t) const
 
bool testVertices (const reco::Track &t, const reco::VertexCollection &vtxs) const
 
bool testVertices (const reco::TrackRef &t, const reco::VertexCollection &vtxs) const
 
 TrackWithVertexSelector (const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
 
 TrackWithVertexSelector (const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
 
 ~TrackWithVertexSelector ()
 

Private Types

typedef math::XYZPoint Point
 

Private Attributes

double d0Max_
 
double dzMax_
 
double etaMax_
 
double etaMin_
 
double normalizedChi2_
 
double nSigmaDtVertex_
 
uint32_t numberOfLostHits_
 
uint32_t numberOfValidHits_
 
uint32_t numberOfValidPixelHits_
 
uint32_t nVertices_
 
double ptErrorCut_
 
double ptMax_
 
double ptMin_
 
std::string quality_
 
double rhoVtx_
 
edm::ValueMap< float > const * timeresoscoll_ = nullptr
 
edm::EDGetTokenT< edm::ValueMap< float > > timeResosToken_
 
edm::ValueMap< float > const * timescoll_ = nullptr
 
edm::EDGetTokenT< edm::ValueMap< float > > timesToken_
 
reco::VertexCollection const * vcoll_ = nullptr
 
edm::EDGetTokenT< reco::VertexCollectionvertexToken_
 
bool vtxFallback_
 
double zetaVtx_
 

Detailed Description

Definition at line 17 of file TrackWithVertexSelector.h.

Member Typedef Documentation

◆ Point

Definition at line 63 of file TrackWithVertexSelector.h.

Constructor & Destructor Documentation

◆ TrackWithVertexSelector() [1/2]

TrackWithVertexSelector::TrackWithVertexSelector ( const edm::ParameterSet iConfig,
edm::ConsumesCollector &&  iC 
)
inlineexplicit

Definition at line 19 of file TrackWithVertexSelector.h.

20  : TrackWithVertexSelector(iConfig, iC) {}
TrackWithVertexSelector(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)

◆ TrackWithVertexSelector() [2/2]

TrackWithVertexSelector::TrackWithVertexSelector ( const edm::ParameterSet iConfig,
edm::ConsumesCollector iC 
)
explicit

Definition at line 11 of file TrackWithVertexSelector.cc.

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),
28  vtxFallback_(iConfig.getParameter<bool>("vtxFallback")),
29  zetaVtx_(iConfig.getParameter<double>("zetaVtx")),
30  rhoVtx_(iConfig.getParameter<double>("rhoVtx")),
31  nSigmaDtVertex_(iConfig.getParameter<double>("nSigmaDtVertex")) {}
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
edm::EDGetTokenT< edm::ValueMap< float > > timesToken_
edm::EDGetTokenT< edm::ValueMap< float > > timeResosToken_

◆ ~TrackWithVertexSelector()

TrackWithVertexSelector::~TrackWithVertexSelector ( )

Definition at line 33 of file TrackWithVertexSelector.cc.

33 {}

Member Function Documentation

◆ init() [1/2]

void TrackWithVertexSelector::init ( const edm::Event event,
const edm::EventSetup  
)
inline

Definition at line 24 of file TrackWithVertexSelector.h.

References init().

Referenced by init(), and operator()().

24 { init(event); }
void init(const edm::Event &event, const edm::EventSetup &)
Definition: event.py:1

◆ init() [2/2]

void TrackWithVertexSelector::init ( const edm::Event event)

Definition at line 35 of file TrackWithVertexSelector.cc.

References edm::HandleBase::isValid(), nVertices_, edm::Handle< T >::product(), timeresoscoll_, timeResosToken_, timescoll_, timesToken_, vcoll_, and vertexToken_.

35  {
37  event.getByToken(vertexToken_, hVtx);
38  vcoll_ = nVertices_ > 0 ? hVtx.product() : nullptr;
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 }
reco::VertexCollection const * vcoll_
T const * product() const
Definition: Handle.h:70
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
edm::EDGetTokenT< edm::ValueMap< float > > timesToken_
edm::EDGetTokenT< edm::ValueMap< float > > timeResosToken_
edm::ValueMap< float > const * timeresoscoll_
bool isValid() const
Definition: HandleBase.h:70
edm::ValueMap< float > const * timescoll_

◆ operator()() [1/4]

bool TrackWithVertexSelector::operator() ( const reco::Track t) const

Definition at line 116 of file TrackWithVertexSelector.cc.

References nVertices_, submitPVValidationJobs::t, testTrack(), testVertices(), and vcoll_.

116  {
117  if (!testTrack(t))
118  return false;
119  if (nVertices_ == 0)
120  return true;
121  return testVertices(t, *vcoll_);
122 }
reco::VertexCollection const * vcoll_
bool testVertices(const reco::Track &t, const reco::VertexCollection &vtxs) const
bool testTrack(const reco::Track &t) const

◆ operator()() [2/4]

bool TrackWithVertexSelector::operator() ( const reco::Track t,
const edm::Event iEvent 
)
inline

Definition at line 28 of file TrackWithVertexSelector.h.

References iEvent, init(), and submitPVValidationJobs::t.

28  {
29  init(iEvent);
30  return (*this)(t);
31  }
void init(const edm::Event &event, const edm::EventSetup &)
int iEvent
Definition: GenABIO.cc:224

◆ operator()() [3/4]

bool TrackWithVertexSelector::operator() ( const reco::TrackRef t) const

Definition at line 124 of file TrackWithVertexSelector.cc.

References nVertices_, testTrack(), testVertices(), and vcoll_.

124  {
125  if (!testTrack(tref))
126  return false;
127  if (nVertices_ == 0)
128  return true;
129  return testVertices(tref, *vcoll_);
130 }
reco::VertexCollection const * vcoll_
bool testVertices(const reco::Track &t, const reco::VertexCollection &vtxs) const
bool testTrack(const reco::Track &t) const

◆ operator()() [4/4]

bool TrackWithVertexSelector::operator() ( const reco::TrackRef t,
const edm::Event iEvent 
)
inline

Definition at line 34 of file TrackWithVertexSelector.h.

References iEvent, init(), and submitPVValidationJobs::t.

34  {
35  init(iEvent);
36  return (*this)(t);
37  }
void init(const edm::Event &event, const edm::EventSetup &)
int iEvent
Definition: GenABIO.cc:224

◆ testTrack() [1/2]

bool TrackWithVertexSelector::testTrack ( const reco::Track t) const

Definition at line 49 of file TrackWithVertexSelector.cc.

References funct::abs(), d0Max_, dzMax_, etaMax_, etaMin_, SiStripPI::max, normalizedChi2_, numberOfLostHits_, numberOfValidHits_, numberOfValidPixelHits_, ptErrorCut_, ptMax_, ptMin_, quality_, and submitPVValidationJobs::t.

Referenced by operator()(), and testTrack().

49  {
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 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ testTrack() [2/2]

bool TrackWithVertexSelector::testTrack ( const reco::TrackRef t) const

Definition at line 62 of file TrackWithVertexSelector.cc.

References testTrack().

62 { return testTrack(*tref); }
bool testTrack(const reco::Track &t) const

◆ testVertices() [1/2]

bool TrackWithVertexSelector::testVertices ( const reco::Track t,
const reco::VertexCollection vtxs 
) const

Definition at line 64 of file TrackWithVertexSelector.cc.

References funct::abs(), nVertices_, convertSQLiteXML::ok, rhoVtx_, submitPVValidationJobs::t, vtxFallback_, and zetaVtx_.

Referenced by operator()().

64  {
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 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ testVertices() [2/2]

bool TrackWithVertexSelector::testVertices ( const reco::TrackRef t,
const reco::VertexCollection vtxs 
) const

Definition at line 82 of file TrackWithVertexSelector.cc.

References funct::abs(), MillePedeFileConverter_cfg::e, edm::isNotFinite(), nSigmaDtVertex_, nVertices_, convertSQLiteXML::ok, rhoVtx_, mathSSE::sqrt(), submitPVValidationJobs::t, hcalRecHitTable_cff::time, timeresoscoll_, timescoll_, vtxFallback_, and zetaVtx_.

82  {
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 }
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::ValueMap< float > const * timeresoscoll_
edm::ValueMap< float > const * timescoll_

Member Data Documentation

◆ d0Max_

double TrackWithVertexSelector::d0Max_
private

Definition at line 50 of file TrackWithVertexSelector.h.

Referenced by testTrack().

◆ dzMax_

double TrackWithVertexSelector::dzMax_
private

Definition at line 50 of file TrackWithVertexSelector.h.

Referenced by testTrack().

◆ etaMax_

double TrackWithVertexSelector::etaMax_
private

Definition at line 49 of file TrackWithVertexSelector.h.

Referenced by testTrack().

◆ etaMin_

double TrackWithVertexSelector::etaMin_
private

Definition at line 49 of file TrackWithVertexSelector.h.

Referenced by testTrack().

◆ normalizedChi2_

double TrackWithVertexSelector::normalizedChi2_
private

Definition at line 48 of file TrackWithVertexSelector.h.

Referenced by testTrack().

◆ nSigmaDtVertex_

double TrackWithVertexSelector::nSigmaDtVertex_
private

Definition at line 58 of file TrackWithVertexSelector.h.

Referenced by testVertices().

◆ numberOfLostHits_

uint32_t TrackWithVertexSelector::numberOfLostHits_
private

Definition at line 47 of file TrackWithVertexSelector.h.

Referenced by testTrack().

◆ numberOfValidHits_

uint32_t TrackWithVertexSelector::numberOfValidHits_
private

Definition at line 45 of file TrackWithVertexSelector.h.

Referenced by testTrack().

◆ numberOfValidPixelHits_

uint32_t TrackWithVertexSelector::numberOfValidPixelHits_
private

Definition at line 46 of file TrackWithVertexSelector.h.

Referenced by testTrack().

◆ nVertices_

uint32_t TrackWithVertexSelector::nVertices_
private

Definition at line 54 of file TrackWithVertexSelector.h.

Referenced by init(), operator()(), and testVertices().

◆ ptErrorCut_

double TrackWithVertexSelector::ptErrorCut_
private

Definition at line 51 of file TrackWithVertexSelector.h.

Referenced by testTrack().

◆ ptMax_

double TrackWithVertexSelector::ptMax_
private

Definition at line 49 of file TrackWithVertexSelector.h.

Referenced by testTrack().

◆ ptMin_

double TrackWithVertexSelector::ptMin_
private

Definition at line 49 of file TrackWithVertexSelector.h.

Referenced by testTrack().

◆ quality_

std::string TrackWithVertexSelector::quality_
private

Definition at line 52 of file TrackWithVertexSelector.h.

Referenced by testTrack().

◆ rhoVtx_

double TrackWithVertexSelector::rhoVtx_
private

Definition at line 58 of file TrackWithVertexSelector.h.

Referenced by testVertices().

◆ timeresoscoll_

edm::ValueMap<float> const* TrackWithVertexSelector::timeresoscoll_ = nullptr
private

Definition at line 62 of file TrackWithVertexSelector.h.

Referenced by init(), and testVertices().

◆ timeResosToken_

edm::EDGetTokenT<edm::ValueMap<float> > TrackWithVertexSelector::timeResosToken_
private

Definition at line 56 of file TrackWithVertexSelector.h.

Referenced by init().

◆ timescoll_

edm::ValueMap<float> const* TrackWithVertexSelector::timescoll_ = nullptr
private

Definition at line 61 of file TrackWithVertexSelector.h.

Referenced by init(), and testVertices().

◆ timesToken_

edm::EDGetTokenT<edm::ValueMap<float> > TrackWithVertexSelector::timesToken_
private

Definition at line 56 of file TrackWithVertexSelector.h.

Referenced by init().

◆ vcoll_

reco::VertexCollection const* TrackWithVertexSelector::vcoll_ = nullptr
private

Definition at line 60 of file TrackWithVertexSelector.h.

Referenced by init(), and operator()().

◆ vertexToken_

edm::EDGetTokenT<reco::VertexCollection> TrackWithVertexSelector::vertexToken_
private

Definition at line 55 of file TrackWithVertexSelector.h.

Referenced by init().

◆ vtxFallback_

bool TrackWithVertexSelector::vtxFallback_
private

Definition at line 57 of file TrackWithVertexSelector.h.

Referenced by testVertices().

◆ zetaVtx_

double TrackWithVertexSelector::zetaVtx_
private

Definition at line 58 of file TrackWithVertexSelector.h.

Referenced by testVertices().