CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTrackSelector.h
Go to the documentation of this file.
1 #ifndef RecoSelectors_RecoTrackSelector_h
2 #define RecoSelectors_RecoTrackSelector_h
3 /* \class RecoTrackSelector
4  *
5  * \author Giuseppe Cerati, INFN
6  *
7  * $Date: 2013/06/05 14:49:09 $
8  * $Revision: 1.3.12.2 $
9  *
10  */
16 
18  public:
20  typedef std::vector<const reco::Track *> container;
21  typedef container::const_iterator const_iterator;
22 
26  RecoTrackSelector( cfg, iC ) {}
28  ptMin_(cfg.getParameter<double>("ptMin")),
29  minRapidity_(cfg.getParameter<double>("minRapidity")),
30  maxRapidity_(cfg.getParameter<double>("maxRapidity")),
31  tip_(cfg.getParameter<double>("tip")),
32  lip_(cfg.getParameter<double>("lip")),
33  minHit_(cfg.getParameter<int>("minHit")),
34  min3DHit_(cfg.getParameter<int>("min3DHit")),
35  maxChi2_(cfg.getParameter<double>("maxChi2")),
36  bsSrcToken_(iC.consumes<reco::BeamSpot>(cfg.getParameter<edm::InputTag>("beamSpot"))),
37  bs(0)
38  {
39  std::vector<std::string> quality = cfg.getParameter<std::vector<std::string> >("quality");
40  for (unsigned int j=0;j<quality.size();j++) quality_.push_back(reco::TrackBase::qualityByName(quality[j]));
41  std::vector<std::string> algorithm = cfg.getParameter<std::vector<std::string> >("algorithm");
42  for (unsigned int j=0;j<algorithm.size();j++) algorithm_.push_back(reco::TrackBase::algoByName(algorithm[j]));
43  }
44 
45  RecoTrackSelector ( double ptMin, double minRapidity, double maxRapidity,
46  double tip, double lip, int minHit, int min3DHit, double maxChi2,
47  const std::vector<std::string>& quality , const std::vector<std::string>& algorithm ) :
48  ptMin_( ptMin ), minRapidity_( minRapidity ), maxRapidity_( maxRapidity ),
49  tip_( tip ), lip_( lip ), minHit_( minHit ), min3DHit_( min3DHit), maxChi2_( maxChi2 ), bs(0)
50  {
51  for (unsigned int j=0;j<quality.size();j++) quality_.push_back(reco::TrackBase::qualityByName(quality[j]));
52  for (unsigned int j=0;j<algorithm.size();j++) algorithm_.push_back(reco::TrackBase::algoByName(algorithm[j]));
53  }
54 
55  const_iterator begin() const { return selected_.begin(); }
56  const_iterator end() const { return selected_.end(); }
57 
59  selected_.clear();
61  event.getByToken(bsSrcToken_,beamSpot);
62  bs = beamSpot.product();
63  for( reco::TrackCollection::const_iterator trk = c->begin();
64  trk != c->end(); ++ trk )
65  if ( operator()(*trk) ) {
66  selected_.push_back( & * trk );
67  }
68  }
69 
72 
73  if ((bs==0)|| (previousEvent != event.id())) {
75  event.getByToken(bsSrcToken_,beamSpot);
76  bs = beamSpot.product();
77  previousEvent = event.id();
78  }
79  return operator()(t);
80  }
81 
82  bool operator()( const reco::Track & t, const reco::BeamSpot* bs_) {
83  bs = bs_;
84  return operator()(t);
85  }
86 
87  bool operator()( const reco::Track & t) {
88  bool quality_ok = true;
89  if (quality_.size()!=0) {
90  quality_ok = false;
91  for (unsigned int i = 0; i<quality_.size();++i) {
92  if (t.quality(quality_[i])){
93  quality_ok = true;
94  break;
95  }
96  }
97  }
98  bool algo_ok = true;
99  if (algorithm_.size()!=0) {
100  if (std::find(algorithm_.begin(),algorithm_.end(),t.algo())==algorithm_.end()) algo_ok = false;
101  }
102  return
106  fabs(t.pt()) >= ptMin_ &&
107  t.eta() >= minRapidity_ && t.eta() <= maxRapidity_ &&
108  fabs(t.dxy(bs->position())) <= tip_ &&
109  fabs(t.dsz(bs->position())) <= lip_ &&
110  t.normalizedChi2()<=maxChi2_ &&
111  quality_ok &&
112  algo_ok);
113  }
114 
115  size_t size() const { return selected_.size(); }
116 
117  protected:
118  double ptMin_;
119  double minRapidity_;
120  double maxRapidity_;
121  double tip_;
122  double lip_;
123  int minHit_;
125  double maxChi2_;
126  std::vector<reco::TrackBase::TrackQuality> quality_;
127  std::vector<reco::TrackBase::TrackAlgorithm> algorithm_;
132 };
133 
134 #endif
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
tuple cfg
Definition: looper.py:259
container::const_iterator const_iterator
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:609
std::vector< reco::TrackBase::TrackQuality > quality_
const_iterator end() const
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
RecoTrackSelector()
Constructors.
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:458
edm::EDGetTokenT< reco::BeamSpot > bsSrcToken_
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:651
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:477
edm::EventID previousEvent
TrackAlgorithm algo() const
Definition: TrackBase.h:423
RecoTrackSelector(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:699
int numberOfValidStripLayersWithMonoAndStereo(uint16_t stripdet, uint16_t layer) const
Definition: HitPattern.cc:300
const_iterator begin() const
bool operator()(const reco::Track &t, edm::Event &event)
Operator() performs the selection: e.g. if (recoTrackSelector(track)) {...}.
double pt() const
track transverse momentum
Definition: TrackBase.h:669
RecoTrackSelector(double ptMin, double minRapidity, double maxRapidity, double tip, double lip, int minHit, int min3DHit, double maxChi2, const std::vector< std::string > &quality, const std::vector< std::string > &algorithm)
int j
Definition: DBlmapReader.cc:9
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
RecoTrackSelector(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:108
T const * product() const
Definition: Handle.h:81
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:384
std::vector< reco::TrackBase::TrackAlgorithm > algorithm_
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:551
void select(const edm::Handle< collection > &c, const edm::Event &event, const edm::EventSetup &)
edm::EventID id() const
Definition: EventBase.h:60
bool operator()(const reco::Track &t, const reco::BeamSpot *bs_)
size_t size() const
std::vector< const reco::Track * > container
const reco::BeamSpot * bs
bool operator()(const reco::Track &t)
static TrackAlgorithm algoByName(const std::string &name)
Definition: TrackBase.cc:120
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:639
reco::TrackCollection collection