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 
7 
13 
14 
16 public:
19  ptMin_(cfg.getParameter<double>("ptMin")),
20  minRapidity_(cfg.getParameter<double>("minRapidity")),
21  maxRapidity_(cfg.getParameter<double>("maxRapidity")),
22  tip_(cfg.getParameter<double>("tip")),
23  lip_(cfg.getParameter<double>("lip")),
24  maxChi2_(cfg.getParameter<double>("maxChi2")),
25  minHit_(cfg.getParameter<int>("minHit")),
26  minPixelHit_(cfg.getParameter<int>("minPixelHit")),
27  minLayer_(cfg.getParameter<int>("minLayer")),
28  min3DLayer_(cfg.getParameter<int>("min3DLayer")),
29  usePV_(cfg.getParameter<bool>("usePV")),
30  bsSrcToken_(iC.consumes<reco::BeamSpot>(cfg.getParameter<edm::InputTag>("beamSpot"))) {
31  if (usePV_)
33  for(const std::string& quality: cfg.getParameter<std::vector<std::string> >("quality"))
35  for(const std::string& algorithm: cfg.getParameter<std::vector<std::string> >("algorithm"))
37  for(const std::string& algorithm: cfg.getParameter<std::vector<std::string> >("originalAlgorithm"))
39  for(const std::string& algorithm: cfg.getParameter<std::vector<std::string> >("algorithmMaskContains"))
41  }
42 
43  void init(const edm::Event& event, const edm::EventSetup& es) {
45  event.getByToken(bsSrcToken_,beamSpot);
46  vertex_ = beamSpot->position();
47  if (!usePV_) return;
49  event.getByToken(vertexToken_, hVtx);
50  if (hVtx->empty()) return;
51  vertex_ = (*hVtx)[0].position();
52  }
53 
54  bool operator()( const reco::TrackRef& tref ) const {
55  return (*this)(*tref);
56  }
57 
58  bool operator()( const reco::Track & t) const {
59  bool quality_ok = true;
60  if (quality_.size()!=0) {
61  quality_ok = false;
62  for (unsigned int i = 0; i<quality_.size();++i) {
63  if (t.quality(quality_[i])){
64  quality_ok = true;
65  break;
66  }
67  }
68  }
69 
70  bool algo_ok = true;
71  if (algorithm_.size()!=0) {
72  if (std::find(algorithm_.begin(),algorithm_.end(),t.algo())==algorithm_.end()) algo_ok = false;
73  }
74  if (!originalAlgorithm_.empty() && algo_ok) {
75  if (std::find(originalAlgorithm_.begin(), originalAlgorithm_.end(), t.originalAlgo()) == originalAlgorithm_.end()) algo_ok = false;
76  }
77  if(!algorithmMask_.empty() && algo_ok) {
78  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>()
79  return t.algoMask()[algo];
80  }) == algorithmMask_.end()) algo_ok = false;
81  }
82  return
83  (
84  (algo_ok & quality_ok) &&
90  fabs(t.pt()) >= ptMin_ &&
91  t.eta() >= minRapidity_ && t.eta() <= maxRapidity_ &&
92  fabs(t.dxy(vertex_)) <= tip_ &&
93  fabs(t.dsz(vertex_)) <= lip_ &&
95  );
96  }
97 
98 
99 private:
100  double ptMin_;
101  double minRapidity_;
102  double maxRapidity_;
103  double tip_;
104  double lip_;
105  double maxChi2_;
106  int minHit_;
110  bool usePV_;
111 
114 
115  std::vector<reco::TrackBase::TrackQuality> quality_;
116  std::vector<reco::TrackBase::TrackAlgorithm> algorithm_;
117  std::vector<reco::TrackBase::TrackAlgorithm> originalAlgorithm_;
118  std::vector<reco::TrackBase::TrackAlgorithm> algorithmMask_;
119 
121 };
122 
123 #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:499
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:518
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:341
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