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 Attributes
RecoTrackSelectorBase Class Reference

#include <RecoTrackSelectorBase.h>

Inheritance diagram for RecoTrackSelectorBase:
RecoTrackRefSelector RecoTrackSelector RecoTrackViewRefSelector

Public Member Functions

void init (const edm::Event &event, const edm::EventSetup &es)
 
bool operator() (const reco::Track &t) const
 
 RecoTrackSelectorBase ()
 
 RecoTrackSelectorBase (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
 

Private Attributes

std::vector
< reco::TrackBase::TrackAlgorithm
algorithm_
 
std::vector
< reco::TrackBase::TrackAlgorithm
algorithmMask_
 
edm::EDGetTokenT< reco::BeamSpotbsSrcToken_
 
double lip_
 
double maxChi2_
 
double maxRapidity_
 
int min3DLayer_
 
int minHit_
 
int minLayer_
 
int minPixelHit_
 
double minRapidity_
 
std::vector
< reco::TrackBase::TrackAlgorithm
originalAlgorithm_
 
double ptMin_
 
std::vector
< reco::TrackBase::TrackQuality
quality_
 
double tip_
 
bool usePV_
 
reco::Track::Point vertex_
 
edm::EDGetTokenT
< reco::VertexCollection
vertexToken_
 

Detailed Description

Definition at line 14 of file RecoTrackSelectorBase.h.

Constructor & Destructor Documentation

RecoTrackSelectorBase::RecoTrackSelectorBase ( )
inline

Definition at line 16 of file RecoTrackSelectorBase.h.

16 {}
RecoTrackSelectorBase::RecoTrackSelectorBase ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)
inline

Definition at line 17 of file RecoTrackSelectorBase.h.

References reco::TrackBase::algoByName(), HLT_25ns14e33_v1_cff::algorithm, algorithm_, algorithmMask_, edm::ConsumesCollector::consumes(), edm::ParameterSet::getParameter(), originalAlgorithm_, HLT_25ns14e33_v1_cff::quality, quality_, reco::TrackBase::qualityByName(), AlCaHLTBitMon_QueryRunRegistry::string, usePV_, and vertexToken_.

17  :
18  ptMin_(cfg.getParameter<double>("ptMin")),
19  minRapidity_(cfg.getParameter<double>("minRapidity")),
20  maxRapidity_(cfg.getParameter<double>("maxRapidity")),
21  tip_(cfg.getParameter<double>("tip")),
22  lip_(cfg.getParameter<double>("lip")),
23  maxChi2_(cfg.getParameter<double>("maxChi2")),
24  minHit_(cfg.getParameter<int>("minHit")),
25  minPixelHit_(cfg.getParameter<int>("minPixelHit")),
26  minLayer_(cfg.getParameter<int>("minLayer")),
27  min3DLayer_(cfg.getParameter<int>("min3DLayer")),
28  usePV_(cfg.getParameter<bool>("usePV")),
30  if (usePV_)
32  for(const std::string& quality: cfg.getParameter<std::vector<std::string> >("quality"))
34  for(const std::string& algorithm: cfg.getParameter<std::vector<std::string> >("algorithm"))
36  for(const std::string& algorithm: cfg.getParameter<std::vector<std::string> >("originalAlgorithm"))
38  for(const std::string& algorithm: cfg.getParameter<std::vector<std::string> >("algorithmMaskContains"))
40  }
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
std::vector< reco::TrackBase::TrackAlgorithm > algorithmMask_
std::vector< reco::TrackBase::TrackQuality > quality_
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:123
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
std::vector< reco::TrackBase::TrackAlgorithm > algorithm_
edm::EDGetTokenT< reco::BeamSpot > bsSrcToken_
static TrackAlgorithm algoByName(const std::string &name)
Definition: TrackBase.cc:135
std::vector< reco::TrackBase::TrackAlgorithm > originalAlgorithm_

Member Function Documentation

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

Definition at line 42 of file RecoTrackSelectorBase.h.

References SiPixelRawToDigiRegional_cfi::beamSpot, bsSrcToken_, usePV_, vertex_, and vertexToken_.

Referenced by RecoTrackViewRefSelector::select(), RecoTrackRefSelector::select(), and RecoTrackSelector::select().

42  {
44  event.getByToken(bsSrcToken_,beamSpot);
45  vertex_ = beamSpot->position();
46  if (!usePV_) return;
48  event.getByToken(vertexToken_, hVtx);
49  if (hVtx->empty()) return;
50  vertex_ = (*hVtx)[0].position();
51  }
reco::Track::Point vertex_
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
edm::EDGetTokenT< reco::BeamSpot > bsSrcToken_
bool RecoTrackSelectorBase::operator() ( const reco::Track t) const
inline

Definition at line 53 of file RecoTrackSelectorBase.h.

References ecalcalib_dqm_sourceclient-live_cfg::algo, reco::TrackBase::algo(), reco::TrackBase::algoMask(), algorithm_, algorithmMask_, reco::TrackBase::dsz(), reco::TrackBase::dxy(), reco::TrackBase::eta(), spr::find(), reco::TrackBase::hitPattern(), i, lip_, maxChi2_, maxRapidity_, min3DLayer_, minHit_, minLayer_, minPixelHit_, minRapidity_, reco::TrackBase::normalizedChi2(), reco::HitPattern::numberOfValidHits(), reco::HitPattern::numberOfValidPixelHits(), reco::HitPattern::numberOfValidStripLayersWithMonoAndStereo(), reco::TrackBase::originalAlgo(), originalAlgorithm_, reco::HitPattern::pixelLayersWithMeasurement(), reco::TrackBase::pt(), ptMin_, reco::TrackBase::quality(), quality_, tip_, reco::HitPattern::trackerLayersWithMeasurement(), and vertex_.

53  {
54  bool quality_ok = true;
55  if (quality_.size()!=0) {
56  quality_ok = false;
57  for (unsigned int i = 0; i<quality_.size();++i) {
58  if (t.quality(quality_[i])){
59  quality_ok = true;
60  break;
61  }
62  }
63  }
64 
65  bool algo_ok = true;
66  if (algorithm_.size()!=0) {
67  if (std::find(algorithm_.begin(),algorithm_.end(),t.algo())==algorithm_.end()) algo_ok = false;
68  }
69  if (!originalAlgorithm_.empty() && algo_ok) {
70  if (std::find(originalAlgorithm_.begin(), originalAlgorithm_.end(), t.originalAlgo()) == originalAlgorithm_.end()) algo_ok = false;
71  }
72  if(!algorithmMask_.empty() && algo_ok) {
73  if(std::find_if(algorithmMask_.begin(), algorithmMask_.end(), [&](reco::TrackBase::TrackAlgorithm algo) -> bool { // for some reason I have to either explicitly give the return type, or use static_cast<bool>()
74  return t.algoMask()[algo];
75  }) == algorithmMask_.end()) algo_ok = false;
76  }
77  return
78  (
79  (algo_ok & quality_ok) &&
85  fabs(t.pt()) >= ptMin_ &&
86  t.eta() >= minRapidity_ && t.eta() <= maxRapidity_ &&
87  fabs(t.dxy(vertex_)) <= tip_ &&
88  fabs(t.dsz(vertex_)) <= lip_ &&
90  );
91  }
int i
Definition: DBlmapReader.cc:9
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:548
int numberOfValidHits() const
Definition: HitPattern.h:801
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:497
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
double dsz() const
dsz parameter (THIS IS NOT the SZ impact parameter to (0,0,0) if refPoint is far from (0...
Definition: TrackBase.h:590
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:516
TrackAlgorithm
track algorithm
Definition: TrackBase.h:99
reco::Track::Point vertex_
TrackAlgorithm algo() const
Definition: TrackBase.h:484
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:638
int numberOfValidStripLayersWithMonoAndStereo(uint16_t stripdet, uint16_t layer) const
Definition: HitPattern.cc:339
double pt() const
track transverse momentum
Definition: TrackBase.h:608
std::vector< reco::TrackBase::TrackAlgorithm > algorithmMask_
std::vector< reco::TrackBase::TrackQuality > quality_
TrackAlgorithm originalAlgo() const
Definition: TrackBase.h:488
AlgoMask algoMask() const
Definition: TrackBase.h:361
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:437
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:497
std::vector< reco::TrackBase::TrackAlgorithm > algorithm_
int numberOfValidPixelHits() const
Definition: HitPattern.h:816
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:578
std::vector< reco::TrackBase::TrackAlgorithm > originalAlgorithm_

Member Data Documentation

std::vector<reco::TrackBase::TrackAlgorithm> RecoTrackSelectorBase::algorithm_
private

Definition at line 111 of file RecoTrackSelectorBase.h.

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

std::vector<reco::TrackBase::TrackAlgorithm> RecoTrackSelectorBase::algorithmMask_
private

Definition at line 113 of file RecoTrackSelectorBase.h.

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

edm::EDGetTokenT<reco::BeamSpot> RecoTrackSelectorBase::bsSrcToken_
private

Definition at line 107 of file RecoTrackSelectorBase.h.

Referenced by init().

double RecoTrackSelectorBase::lip_
private

Definition at line 99 of file RecoTrackSelectorBase.h.

Referenced by operator()().

double RecoTrackSelectorBase::maxChi2_
private

Definition at line 100 of file RecoTrackSelectorBase.h.

Referenced by operator()().

double RecoTrackSelectorBase::maxRapidity_
private

Definition at line 97 of file RecoTrackSelectorBase.h.

Referenced by operator()().

int RecoTrackSelectorBase::min3DLayer_
private

Definition at line 104 of file RecoTrackSelectorBase.h.

Referenced by operator()().

int RecoTrackSelectorBase::minHit_
private

Definition at line 101 of file RecoTrackSelectorBase.h.

Referenced by operator()().

int RecoTrackSelectorBase::minLayer_
private

Definition at line 103 of file RecoTrackSelectorBase.h.

Referenced by operator()().

int RecoTrackSelectorBase::minPixelHit_
private

Definition at line 102 of file RecoTrackSelectorBase.h.

Referenced by operator()().

double RecoTrackSelectorBase::minRapidity_
private

Definition at line 96 of file RecoTrackSelectorBase.h.

Referenced by operator()().

std::vector<reco::TrackBase::TrackAlgorithm> RecoTrackSelectorBase::originalAlgorithm_
private

Definition at line 112 of file RecoTrackSelectorBase.h.

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

double RecoTrackSelectorBase::ptMin_
private

Definition at line 95 of file RecoTrackSelectorBase.h.

Referenced by operator()().

std::vector<reco::TrackBase::TrackQuality> RecoTrackSelectorBase::quality_
private

Definition at line 110 of file RecoTrackSelectorBase.h.

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

double RecoTrackSelectorBase::tip_
private

Definition at line 98 of file RecoTrackSelectorBase.h.

Referenced by operator()().

bool RecoTrackSelectorBase::usePV_
private

Definition at line 105 of file RecoTrackSelectorBase.h.

Referenced by init(), and RecoTrackSelectorBase().

reco::Track::Point RecoTrackSelectorBase::vertex_
private

Definition at line 115 of file RecoTrackSelectorBase.h.

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

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

Definition at line 108 of file RecoTrackSelectorBase.h.

Referenced by init(), and RecoTrackSelectorBase().