CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 testTrack (const reco::Track &t) const
 
bool testVertices (const reco::Track &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_
 
uint32_t numberOfLostHits_
 
uint32_t numberOfValidHits_
 
uint32_t numberOfValidPixelHits_
 
uint32_t nVertices_
 
double ptErrorCut_
 
double ptMax_
 
double ptMin_
 
std::string quality_
 
double rhoVtx_
 
reco::VertexCollection const * vcoll_ = 0
 
edm::EDGetTokenT
< reco::VertexCollection
vertexToken_
 
bool vtxFallback_
 
double zetaVtx_
 

Detailed Description

Definition at line 17 of file TrackWithVertexSelector.h.

Member Typedef Documentation

Definition at line 50 of file TrackWithVertexSelector.h.

Constructor & Destructor Documentation

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

Definition at line 19 of file TrackWithVertexSelector.h.

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

Definition at line 6 of file TrackWithVertexSelector.cc.

6  :
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),
21  vtxFallback_(iConfig.getParameter<bool>("vtxFallback")),
22  zetaVtx_(iConfig.getParameter<double>("zetaVtx")),
23  rhoVtx_(iConfig.getParameter<double>("rhoVtx")) {
24 }
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_
TrackWithVertexSelector::~TrackWithVertexSelector ( )

Definition at line 26 of file TrackWithVertexSelector.cc.

26 { }

Member Function Documentation

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 &)
void TrackWithVertexSelector::init ( const edm::Event event)

Definition at line 29 of file TrackWithVertexSelector.cc.

References edm::Handle< T >::product(), vcoll_, and vertexToken_.

29  {
31  event.getByToken(vertexToken_, hVtx);
32  vcoll_ = hVtx.product();
33 }
reco::VertexCollection const * vcoll_
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
T const * product() const
Definition: Handle.h:81
bool TrackWithVertexSelector::operator() ( const reco::Track t) const

Definition at line 73 of file TrackWithVertexSelector.cc.

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

73  {
74  if (!testTrack(t)) return false;
75  if (nVertices_ == 0) return true;
76  return testVertices(t, *vcoll_);
77 }
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(), and lumiQTWidget::t.

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

Definition at line 36 of file TrackWithVertexSelector.cc.

References funct::abs(), reco::TrackBase::d0(), d0Max_, reco::TrackBase::dz(), dzMax_, reco::TrackBase::eta(), etaMax_, etaMin_, reco::TrackBase::hitPattern(), bookConverter::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()().

36  {
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 }
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
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:821
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
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
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
bool TrackWithVertexSelector::testVertices ( const reco::Track t,
const reco::VertexCollection vtxs 
) const

Definition at line 55 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()().

55  {
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 }
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:682
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:604
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

Member Data Documentation

double TrackWithVertexSelector::d0Max_
private

Definition at line 40 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::dzMax_
private

Definition at line 40 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::etaMax_
private

Definition at line 39 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::etaMin_
private

Definition at line 39 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::normalizedChi2_
private

Definition at line 38 of file TrackWithVertexSelector.h.

Referenced by testTrack().

uint32_t TrackWithVertexSelector::numberOfLostHits_
private

Definition at line 37 of file TrackWithVertexSelector.h.

Referenced by testTrack().

uint32_t TrackWithVertexSelector::numberOfValidHits_
private

Definition at line 35 of file TrackWithVertexSelector.h.

Referenced by testTrack().

uint32_t TrackWithVertexSelector::numberOfValidPixelHits_
private

Definition at line 36 of file TrackWithVertexSelector.h.

Referenced by testTrack().

uint32_t TrackWithVertexSelector::nVertices_
private

Definition at line 44 of file TrackWithVertexSelector.h.

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

double TrackWithVertexSelector::ptErrorCut_
private

Definition at line 41 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::ptMax_
private

Definition at line 39 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::ptMin_
private

Definition at line 39 of file TrackWithVertexSelector.h.

Referenced by testTrack().

std::string TrackWithVertexSelector::quality_
private

Definition at line 42 of file TrackWithVertexSelector.h.

Referenced by testTrack().

double TrackWithVertexSelector::rhoVtx_
private

Definition at line 47 of file TrackWithVertexSelector.h.

Referenced by testVertices().

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

Definition at line 49 of file TrackWithVertexSelector.h.

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

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

Definition at line 45 of file TrackWithVertexSelector.h.

Referenced by init().

bool TrackWithVertexSelector::vtxFallback_
private

Definition at line 46 of file TrackWithVertexSelector.h.

Referenced by testVertices().

double TrackWithVertexSelector::zetaVtx_
private

Definition at line 47 of file TrackWithVertexSelector.h.

Referenced by testVertices().