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: 2010/02/11 00:10:50 $
8  * $Revision: 1.3 $
9  *
10  */
15 
17  public:
19  typedef std::vector<const reco::Track *> container;
20  typedef container::const_iterator const_iterator;
21 
25  ptMin_(cfg.getParameter<double>("ptMin")),
26  minRapidity_(cfg.getParameter<double>("minRapidity")),
27  maxRapidity_(cfg.getParameter<double>("maxRapidity")),
28  tip_(cfg.getParameter<double>("tip")),
29  lip_(cfg.getParameter<double>("lip")),
30  minHit_(cfg.getParameter<int>("minHit")),
31  min3DHit_(cfg.getParameter<int>("min3DHit")),
32  maxChi2_(cfg.getParameter<double>("maxChi2")),
33  bsSrc_(cfg.getParameter<edm::InputTag>("beamSpot"))
34  {
35  std::vector<std::string> quality = cfg.getParameter<std::vector<std::string> >("quality");
36  for (unsigned int j=0;j<quality.size();j++) quality_.push_back(reco::TrackBase::qualityByName(quality[j]));
37  std::vector<std::string> algorithm = cfg.getParameter<std::vector<std::string> >("algorithm");
38  for (unsigned int j=0;j<algorithm.size();j++) algorithm_.push_back(reco::TrackBase::algoByName(algorithm[j]));
39  }
40 
41  RecoTrackSelector ( double ptMin, double minRapidity, double maxRapidity,
42  double tip, double lip, int minHit, int min3DHit, double maxChi2,
43  std::vector<std::string> quality , std::vector<std::string> algorithm ) :
44  ptMin_( ptMin ), minRapidity_( minRapidity ), maxRapidity_( maxRapidity ),
45  tip_( tip ), lip_( lip ), minHit_( minHit ), min3DHit_( min3DHit), maxChi2_( maxChi2 )
46  {
47  for (unsigned int j=0;j<quality.size();j++) quality_.push_back(reco::TrackBase::qualityByName(quality[j]));
48  for (unsigned int j=0;j<algorithm.size();j++) algorithm_.push_back(reco::TrackBase::algoByName(algorithm[j]));
49  }
50 
51  const_iterator begin() const { return selected_.begin(); }
52  const_iterator end() const { return selected_.end(); }
53 
55  selected_.clear();
57  event.getByLabel(bsSrc_,beamSpot);
58  for( reco::TrackCollection::const_iterator trk = c->begin();
59  trk != c->end(); ++ trk )
60  if ( operator()(*trk,beamSpot.product()) ) {
61  selected_.push_back( & * trk );
62  }
63  }
64 
66  bool operator()( const reco::Track & t, const reco::BeamSpot* bs) {
67  bool quality_ok = true;
68  if (quality_.size()!=0) {
69  quality_ok = false;
70  for (unsigned int i = 0; i<quality_.size();++i) {
71  if (t.quality(quality_[i])){
72  quality_ok = true;
73  break;
74  }
75  }
76  }
77  bool algo_ok = true;
78  if (algorithm_.size()!=0) {
79  if (std::find(algorithm_.begin(),algorithm_.end(),t.algo())==algorithm_.end()) algo_ok = false;
80  }
81  return
85  fabs(t.pt()) >= ptMin_ &&
86  t.eta() >= minRapidity_ && t.eta() <= maxRapidity_ &&
87  fabs(t.dxy(bs->position())) <= tip_ &&
88  fabs(t.dsz(bs->position())) <= lip_ &&
89  t.normalizedChi2()<=maxChi2_ &&
90  quality_ok &&
91  algo_ok);
92  }
93 
94  size_t size() const { return selected_.size(); }
95 
96  protected:
97  double ptMin_;
98  double minRapidity_;
99  double maxRapidity_;
100  double tip_;
101  double lip_;
102  int minHit_;
104  double maxChi2_;
105  std::vector<reco::TrackBase::TrackQuality> quality_;
106  std::vector<reco::TrackBase::TrackAlgorithm> algorithm_;
109 };
110 
111 #endif
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
< trclass="colgroup">< tdclass="colgroup"colspan=5 > Ecal cluster collections</td ></tr >< tr >< td >< ahref="classreco_1_1BasicCluster.html"> reco::BasicCluster</a ></td >< td >< ahref="DataFormats_EgammaReco.html"> reco::BasicClusterCollection</a ></td >< td >< ahref="#"> hybridSuperClusters</a ></td >< tdclass="description"> Basic clusters reconstructed with hybrid algorithm(barrel only)</td >< td >S.Rahatlou</td ></tr >< tr >< td >< a href
edm::InputTag bsSrc_
container::const_iterator const_iterator
RecoTrackSelector(double ptMin, double minRapidity, double maxRapidity, double tip, double lip, int minHit, int min3DHit, double maxChi2, std::vector< std::string > quality, std::vector< std::string > algorithm)
bool operator()(const reco::Track &t, const reco::BeamSpot *bs)
Operator() performs the selection: e.g. if (recoTrackSelector(track)) {...}.
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:111
std::vector< reco::TrackBase::TrackQuality > quality_
const_iterator end() const
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
int pixelLayersWithMeasurement() const
Definition: HitPattern.h:710
RecoTrackSelector()
Constructors.
int numberOfValidStripLayersWithMonoAndStereo() const
Definition: HitPattern.cc:224
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:125
TrackAlgorithm algo() const
Definition: TrackBase.h:332
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:141
const_iterator begin() const
int trackerLayersWithMeasurement() const
Definition: HitPattern.h:705
double pt() const
track transverse momentum
Definition: TrackBase.h:131
int j
Definition: DBlmapReader.cc:9
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:223
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
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:46
std::vector< reco::TrackBase::TrackAlgorithm > algorithm_
RecoTrackSelector(const edm::ParameterSet &cfg)
T const * product() const
Definition: Handle.h:74
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:377
void select(const edm::Handle< collection > &c, const edm::Event &event, const edm::EventSetup &)
size_t size() const
std::vector< const reco::Track * > container
static TrackAlgorithm algoByName(const std::string &name)
Definition: TrackBase.cc:55
const Point & position() const
position
Definition: BeamSpot.h:63
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:121
reco::TrackCollection collection