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_ = 0
 
edm::EDGetTokenT< edm::ValueMap< float > > timeResosToken_
 
edm::ValueMap< float > const * timescoll_ = 0
 
edm::EDGetTokenT< edm::ValueMap< float > > timesToken_
 
reco::VertexCollection const * vcoll_ = 0
 
edm::EDGetTokenT< reco::VertexCollectionvertexToken_
 
bool vtxFallback_
 
double zetaVtx_
 

Detailed Description

Definition at line 18 of file TrackWithVertexSelector.h.

Member Typedef Documentation

Definition at line 61 of file TrackWithVertexSelector.h.

Constructor & Destructor Documentation

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

Definition at line 20 of file TrackWithVertexSelector.h.

References ~TrackWithVertexSelector().

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

Definition at line 11 of file TrackWithVertexSelector.cc.

11  :
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")) {
32 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
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 ( )

Definition at line 34 of file TrackWithVertexSelector.cc.

Referenced by TrackWithVertexSelector().

34 { }

Member Function Documentation

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

Definition at line 25 of file TrackWithVertexSelector.h.

References init(), operator()(), and lumiQTWidget::t.

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

25 {init(event);}
void init(const edm::Event &event, const edm::EventSetup &)
void TrackWithVertexSelector::init ( const edm::Event event)

Definition at line 37 of file TrackWithVertexSelector.cc.

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

37  {
39  event.getByToken(vertexToken_, hVtx);
40  vcoll_ = hVtx.product();
41 
43  event.getByToken(timesToken_, hTimes);
44  timescoll_ = hTimes.isValid() ? hTimes.product() : nullptr;
45 
47  event.getByToken(timeResosToken_, hTimeResos);
48  timeresoscoll_ = hTimeResos.isValid() ? hTimeResos.product() : nullptr;
49 }
reco::VertexCollection const * vcoll_
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
edm::EDGetTokenT< edm::ValueMap< float > > timesToken_
edm::EDGetTokenT< edm::ValueMap< float > > timeResosToken_
bool isValid() const
Definition: HandleBase.h:74
T const * product() const
Definition: Handle.h:81
edm::ValueMap< float > const * timeresoscoll_
edm::ValueMap< float > const * timescoll_
bool TrackWithVertexSelector::operator() ( const reco::Track t) const

Definition at line 126 of file TrackWithVertexSelector.cc.

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

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

126  {
127  if (!testTrack(t)) return false;
128  if (nVertices_ == 0) return true;
129  return testVertices(t, *vcoll_);
130 }
bool testTrack(const reco::Track &t) const
bool testVertices(const reco::Track &t, const reco::VertexCollection &vtxs) const
reco::VertexCollection const * vcoll_
bool TrackWithVertexSelector::operator() ( const reco::Track t,
const edm::Event iEvent 
)
inline

Definition at line 29 of file TrackWithVertexSelector.h.

References init(), operator()(), and lumiQTWidget::t.

29  {
30  init(iEvent); return (*this)(t);
31  }
void init(const edm::Event &event, const edm::EventSetup &)
bool TrackWithVertexSelector::operator() ( const reco::TrackRef t) const

Definition at line 132 of file TrackWithVertexSelector.cc.

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

132  {
133  if (!testTrack(tref)) return false;
134  if (nVertices_ == 0) return true;
135  return testVertices(tref, *vcoll_);
136 }
bool testTrack(const reco::Track &t) const
bool testVertices(const reco::Track &t, const reco::VertexCollection &vtxs) const
reco::VertexCollection const * vcoll_
bool TrackWithVertexSelector::operator() ( const reco::TrackRef t,
const edm::Event iEvent 
)
inline

Definition at line 34 of file TrackWithVertexSelector.h.

References init(), lumiQTWidget::t, testTrack(), and testVertices().

34  {
35  init(iEvent); return (*this)(t);
36  }
void init(const edm::Event &event, const edm::EventSetup &)
bool TrackWithVertexSelector::testTrack ( const reco::Track t) const

Definition at line 51 of file TrackWithVertexSelector.cc.

References funct::abs(), reco::TrackBase::d0(), d0Max_, reco::TrackBase::dz(), dzMax_, reco::TrackBase::eta(), etaMax_, etaMin_, reco::TrackBase::hitPattern(), SiStripPI::max, reco::TrackBase::normalizedChi2(), normalizedChi2_, reco::TrackBase::numberOfLostHits(), numberOfLostHits_, reco::TrackBase::numberOfValidHits(), numberOfValidHits_, reco::HitPattern::numberOfValidPixelHits(), numberOfValidPixelHits_, reco::TrackBase::pt(), reco::TrackBase::ptError(), ptErrorCut_, ptMax_, ptMin_, reco::TrackBase::quality(), quality_, and reco::TrackBase::qualityByName().

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

51  {
52  using std::abs;
54  (static_cast<unsigned int>(t.hitPattern().numberOfValidPixelHits()) >= numberOfValidPixelHits_) &&
57  (t.ptError()/t.pt()*std::max(1.,t.normalizedChi2()) <= ptErrorCut_) &&
58  (t.quality(t.qualityByName(quality_))) &&
59  (t.pt() >= ptMin_) &&
60  (t.pt() <= ptMax_) &&
61  (abs(t.eta()) <= etaMax_) &&
62  (abs(t.eta()) >= etaMin_) &&
63  (abs(t.dz()) <= dzMax_) &&
64  (abs(t.d0()) <= d0Max_) ) {
65  return true;
66  }
67  return false;
68 }
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:597
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:561
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:826
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:651
double pt() const
track transverse momentum
Definition: TrackBase.h:621
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:763
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:820
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:609
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:125
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:446
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:510
int numberOfValidPixelHits() const
Definition: HitPattern.h:838
bool TrackWithVertexSelector::testTrack ( const reco::TrackRef t) const

Definition at line 70 of file TrackWithVertexSelector.cc.

References testTrack().

70  {
71  return testTrack(*tref);
72 }
bool testTrack(const reco::Track &t) const
bool TrackWithVertexSelector::testVertices ( const reco::Track t,
const reco::VertexCollection vtxs 
) const

Definition at line 74 of file TrackWithVertexSelector.cc.

References funct::abs(), reco::TrackBase::dxy(), reco::TrackBase::dz(), nVertices_, convertSQLiteXML::ok, rhoVtx_, reco::TrackBase::vertex(), vtxFallback_, and zetaVtx_.

Referenced by operator()().

74  {
75  bool ok = false;
76  if (!vtxs.empty()) {
77  unsigned int tested = 1;
78  for (reco::VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end();
79  it != ed; ++it) {
80  if ( (std::abs(t.dxy(it->position())) < rhoVtx_) &&
81  (std::abs(t.dz(it->position())) < zetaVtx_) ) {
82  ok = true; break;
83  }
84  if (tested++ >= nVertices_) break;
85  }
86  } else if (vtxFallback_) {
87  return ( (std::abs(t.vertex().z()) < 15.9) && (t.vertex().Rho() < 0.2) );
88  }
89  return ok;
90 }
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:687
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
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:609
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:591
bool TrackWithVertexSelector::testVertices ( const reco::TrackRef t,
const reco::VertexCollection vtxs 
) const

Definition at line 92 of file TrackWithVertexSelector.cc.

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

92  {
93  const auto& t = *tref;
94  const bool timeAvailable = timescoll_ != nullptr && timeresoscoll_ != nullptr;
95  bool ok = false;
96  if (!vtxs.empty()) {
97  unsigned int tested = 1;
98  for (reco::VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end();
99  it != ed; ++it) {
100  const bool useTime = timeAvailable && it->t() != 0.;
101  float time = useTime ? (*timescoll_)[tref] : -1.f;
102  float timeReso = useTime ? (*timeresoscoll_)[tref] : -1.f;
103  timeReso = ( timeReso > 1e-6 ? timeReso : fakeBeamSpotTimeWidth );
104 
105  if( edm::isNotFinite(time) ) {
106  time = 0.0;
107  timeReso = 2.0*fakeBeamSpotTimeWidth;
108  }
109 
110  const double vtxSigmaT2 = it->tError() * it->tError();
111  const double vtxTrackErr = std::sqrt( vtxSigmaT2 + timeReso*timeReso );
112 
113  if ( (std::abs(t.dxy(it->position())) < rhoVtx_) &&
114  (std::abs(t.dz(it->position())) < zetaVtx_) &&
115  ( !useTime || (std::abs(time - it->t())/vtxTrackErr < nSigmaDtVertex_) ) ) {
116  ok = true; break;
117  }
118  if (tested++ >= nVertices_) break;
119  }
120  } else if (vtxFallback_) {
121  return ( (std::abs(t.vertex().z()) < 15.9) && (t.vertex().Rho() < 0.2) );
122  }
123  return ok;
124 }
bool isNotFinite(T x)
Definition: isFinite.h:10
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::ValueMap< float > const * timeresoscoll_
edm::ValueMap< float > const * timescoll_

Member Data Documentation

double TrackWithVertexSelector::d0Max_
private

Definition at line 48 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::dzMax_
private

Definition at line 48 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::etaMax_
private

Definition at line 47 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::etaMin_
private

Definition at line 47 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::normalizedChi2_
private

Definition at line 46 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::nSigmaDtVertex_
private

Definition at line 56 of file TrackWithVertexSelector.h.

Referenced by testVertices().

uint32_t TrackWithVertexSelector::numberOfLostHits_
private

Definition at line 45 of file TrackWithVertexSelector.h.

Referenced by testTrack().

uint32_t TrackWithVertexSelector::numberOfValidHits_
private

Definition at line 43 of file TrackWithVertexSelector.h.

Referenced by testTrack().

uint32_t TrackWithVertexSelector::numberOfValidPixelHits_
private

Definition at line 44 of file TrackWithVertexSelector.h.

Referenced by testTrack().

uint32_t TrackWithVertexSelector::nVertices_
private

Definition at line 52 of file TrackWithVertexSelector.h.

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

double TrackWithVertexSelector::ptErrorCut_
private

Definition at line 49 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::ptMax_
private

Definition at line 47 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::ptMin_
private

Definition at line 47 of file TrackWithVertexSelector.h.

Referenced by testTrack().

std::string TrackWithVertexSelector::quality_
private

Definition at line 50 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::rhoVtx_
private

Definition at line 56 of file TrackWithVertexSelector.h.

Referenced by testVertices().

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

Definition at line 60 of file TrackWithVertexSelector.h.

Referenced by init(), and testVertices().

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

Definition at line 54 of file TrackWithVertexSelector.h.

Referenced by init().

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

Definition at line 59 of file TrackWithVertexSelector.h.

Referenced by init(), and testVertices().

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

Definition at line 54 of file TrackWithVertexSelector.h.

Referenced by init().

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

Definition at line 58 of file TrackWithVertexSelector.h.

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

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

Definition at line 53 of file TrackWithVertexSelector.h.

Referenced by init().

bool TrackWithVertexSelector::vtxFallback_
private

Definition at line 55 of file TrackWithVertexSelector.h.

Referenced by testVertices().

double TrackWithVertexSelector::zetaVtx_
private

Definition at line 56 of file TrackWithVertexSelector.h.

Referenced by testVertices().