CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTrackSelectorBase.h
Go to the documentation of this file.
1 #ifndef CommonTools_RecoAlgos_RecoTrackSelectorBase_h
2 #define CommonTools_RecoAlgos_RecoTrackSelectorBase_h
3 
6 
12 
13 
15 public:
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")),
29  bsSrcToken_(iC.consumes<reco::BeamSpot>(cfg.getParameter<edm::InputTag>("beamSpot"))) {
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  }
41 
42  void init(const edm::Event& event, const edm::EventSetup& es) {
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  }
52 
53  bool operator()( const reco::TrackRef& tref ) const {
54  return (*this)(*tref);
55  }
56 
57  bool operator()( const reco::Track & t) const {
58  bool quality_ok = true;
59  if (quality_.size()!=0) {
60  quality_ok = false;
61  for (unsigned int i = 0; i<quality_.size();++i) {
62  if (t.quality(quality_[i])){
63  quality_ok = true;
64  break;
65  }
66  }
67  }
68 
69  bool algo_ok = true;
70  if (algorithm_.size()!=0) {
71  if (std::find(algorithm_.begin(),algorithm_.end(),t.algo())==algorithm_.end()) algo_ok = false;
72  }
73  if (!originalAlgorithm_.empty() && algo_ok) {
74  if (std::find(originalAlgorithm_.begin(), originalAlgorithm_.end(), t.originalAlgo()) == originalAlgorithm_.end()) algo_ok = false;
75  }
76  if(!algorithmMask_.empty() && algo_ok) {
77  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>()
78  return t.algoMask()[algo];
79  }) == algorithmMask_.end()) algo_ok = false;
80  }
81  return
82  (
83  (algo_ok & quality_ok) &&
89  fabs(t.pt()) >= ptMin_ &&
90  t.eta() >= minRapidity_ && t.eta() <= maxRapidity_ &&
91  fabs(t.dxy(vertex_)) <= tip_ &&
92  fabs(t.dsz(vertex_)) <= lip_ &&
94  );
95  }
96 
97 
98 private:
99  double ptMin_;
100  double minRapidity_;
101  double maxRapidity_;
102  double tip_;
103  double lip_;
104  double maxChi2_;
105  int minHit_;
109  bool usePV_;
110 
113 
114  std::vector<reco::TrackBase::TrackQuality> quality_;
115  std::vector<reco::TrackBase::TrackAlgorithm> algorithm_;
116  std::vector<reco::TrackBase::TrackAlgorithm> originalAlgorithm_;
117  std::vector<reco::TrackBase::TrackAlgorithm> algorithmMask_;
118 
120 };
121 
122 #endif
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
void init(const edm::Event &event, const edm::EventSetup &es)
tuple cfg
Definition: looper.py:293
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: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: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
tuple quality
[pTError/pT]*max(1,normChi2) &lt;= ptErrorCut
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
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
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_
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
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_