CMS 3D CMS Logo

RecoTrackSelectorBase.h
Go to the documentation of this file.
1 #ifndef CommonTools_RecoAlgos_RecoTrackSelectorBase_h
2 #define CommonTools_RecoAlgos_RecoTrackSelectorBase_h
3 
8 
14 
15 
17 public:
20  ptMin_(cfg.getParameter<double>("ptMin")),
21  minRapidity_(cfg.getParameter<double>("minRapidity")),
22  maxRapidity_(cfg.getParameter<double>("maxRapidity")),
23  tip_(cfg.getParameter<double>("tip")),
24  lip_(cfg.getParameter<double>("lip")),
25  maxChi2_(cfg.getParameter<double>("maxChi2")),
26  minHit_(cfg.getParameter<int>("minHit")),
27  minPixelHit_(cfg.getParameter<int>("minPixelHit")),
28  minLayer_(cfg.getParameter<int>("minLayer")),
29  min3DLayer_(cfg.getParameter<int>("min3DLayer")),
30  usePV_(cfg.getParameter<bool>("usePV")),
31  bsSrcToken_(iC.consumes<reco::BeamSpot>(cfg.getParameter<edm::InputTag>("beamSpot"))) {
32  if (usePV_)
34  for(const std::string& quality: cfg.getParameter<std::vector<std::string> >("quality"))
36  for(const std::string& algorithm: cfg.getParameter<std::vector<std::string> >("algorithm"))
38  for(const std::string& algorithm: cfg.getParameter<std::vector<std::string> >("originalAlgorithm"))
40  for(const std::string& algorithm: cfg.getParameter<std::vector<std::string> >("algorithmMaskContains"))
42  }
43 
44  void init(const edm::Event& event, const edm::EventSetup& es) {
46  event.getByToken(bsSrcToken_,beamSpot);
47  vertex_ = beamSpot->position();
48  if (!usePV_) return;
50  event.getByToken(vertexToken_, hVtx);
51  if (hVtx->empty()) return;
52  vertex_ = (*hVtx)[0].position();
53  }
54 
55  bool operator()( const reco::TrackRef& tref ) const {
56  return (*this)(*tref);
57  }
58 
59  bool operator()( const reco::Track & t) const {
60  bool quality_ok = true;
61  if (quality_.size()!=0) {
62  quality_ok = false;
63  for (unsigned int i = 0; i<quality_.size();++i) {
64  if (t.quality(quality_[i])){
65  quality_ok = true;
66  break;
67  }
68  }
69  }
70 
71  bool algo_ok = true;
72  if (algorithm_.size()!=0) {
73  if (std::find(algorithm_.begin(),algorithm_.end(),t.algo())==algorithm_.end()) algo_ok = false;
74  }
75  if (!originalAlgorithm_.empty() && algo_ok) {
76  if (std::find(originalAlgorithm_.begin(), originalAlgorithm_.end(), t.originalAlgo()) == originalAlgorithm_.end()) algo_ok = false;
77  }
78  if(!algorithmMask_.empty() && algo_ok) {
79  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>()
80  return t.algoMask()[algo];
81  }) == algorithmMask_.end()) algo_ok = false;
82  }
83  return
84  (
85  (algo_ok & quality_ok) &&
91  fabs(t.pt()) >= ptMin_ &&
92  t.eta() >= minRapidity_ && t.eta() <= maxRapidity_ &&
93  fabs(t.dxy(vertex_)) <= tip_ &&
94  fabs(t.dsz(vertex_)) <= lip_ &&
96  );
97  }
98 
99 
100 private:
101  double ptMin_;
102  double minRapidity_;
103  double maxRapidity_;
104  double tip_;
105  double lip_;
106  double maxChi2_;
107  int minHit_;
111  bool usePV_;
112 
115 
116  std::vector<reco::TrackBase::TrackQuality> quality_;
117  std::vector<reco::TrackBase::TrackAlgorithm> algorithm_;
118  std::vector<reco::TrackBase::TrackAlgorithm> originalAlgorithm_;
119  std::vector<reco::TrackBase::TrackAlgorithm> algorithmMask_;
120 
122 };
123 
124 #endif
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
void init(const edm::Event &event, const edm::EventSetup &es)
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
int numberOfValidHits() const
Definition: HitPattern.h:823
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:493
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
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:598
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:512
TrackAlgorithm
track algorithm
Definition: TrackBase.h:99
reco::Track::Point vertex_
TrackAlgorithm algo() const
Definition: TrackBase.h:492
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:646
int numberOfValidStripLayersWithMonoAndStereo(uint16_t stripdet, uint16_t layer) const
Definition: HitPattern.cc:335
double pt() const
track transverse momentum
Definition: TrackBase.h:616
std::vector< reco::TrackBase::TrackAlgorithm > algorithmMask_
math::XYZPoint Point
point in the space
Definition: TrackBase.h:83
std::vector< reco::TrackBase::TrackQuality > quality_
RecoTrackSelectorBase(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
TrackAlgorithm originalAlgo() const
Definition: TrackBase.h:496
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:125
AlgoMask algoMask() const
Definition: TrackBase.h:363
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:445
bool operator()(const reco::Track &t) const
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:505
std::vector< reco::TrackBase::TrackAlgorithm > algorithm_
fixed size matrix
HLT enums.
edm::EDGetTokenT< reco::BeamSpot > bsSrcToken_
static TrackAlgorithm algoByName(const std::string &name)
Definition: TrackBase.cc:137
int numberOfValidPixelHits() const
Definition: HitPattern.h:838
bool operator()(const reco::TrackRef &tref) const
const Point & position() const
position
Definition: BeamSpot.h:62
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
std::vector< reco::TrackBase::TrackAlgorithm > originalAlgorithm_
Definition: event.py:1