CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonTrackCut.cc
Go to the documentation of this file.
4 
6 {
7 public:
9 
10  result_type operator()(const reco::MuonPtr&) const override final;
11  CandidateType candidateType() const override final { return MUON; }
12  double value(const reco::CandidatePtr&) const override final;
13 
14 private:
15  // inner track selection cuts
20 
22 
23 };
25  MuonTrackCut, "MuonTrackCut");
26 
27 // Define constructors and initialization routines
28 MuonTrackCut::MuonTrackCut(const edm::ParameterSet& c):
30  minTrackerLayersWithMeasurement_(-1),
31  minPixelLayersWithMeasurement_(-1),
32  minNumberOfValidPixelHits_(-1),
33  minValidFraction_(-1),
34  trackQuality_(reco::Track::undefQuality),
35  minNumberOfValidMuonHits_(-1),
36  doInnerTrack_(false), doGlobalTrack_(false)
37 {
38 
39  if ( c.existsAs<edm::ParameterSet>("innerTrack") )
40  {
41  doInnerTrack_ = true;
42 
43  const edm::ParameterSet cc = c.getParameter<edm::ParameterSet>("innerTrack");
44  if ( cc.exists("minTrackerLayersWithMeasurement") ) minTrackerLayersWithMeasurement_ = cc.getParameter<int>("minTrackerLayersWithMeasurement");
45  if ( cc.exists("minPixelLayersWithMeasurement") ) minPixelLayersWithMeasurement_ = cc.getParameter<int>("minPixelLayersWithMeasurement");
46  if ( cc.exists("minNumberOfValidPixelHits") ) minNumberOfValidPixelHits_ = cc.getParameter<int>("minNumberOfValidPixelHits");
47  if ( cc.exists("minValidFraction") ) minValidFraction_ = cc.getParameter<double>("minValidFraction");
48  const std::string trackQualityStr = cc.exists("trackQuality") ? cc.getParameter<std::string>("trackQuality") : "";
49  trackQuality_ = reco::Track::qualityByName(trackQualityStr);
50  }
51  if ( c.existsAs<edm::ParameterSet>("globalTrack") )
52  {
53  doGlobalTrack_ = true;
54 
55  const edm::ParameterSet cc = c.getParameter<edm::ParameterSet>("globalTrack");
56  //if ( cc.exists("maxNormalizedChi2") ) maxNormalizedChi2_ = cc.getParameter<double>("maxNormalizedChi2");
57  if ( cc.exists("minNumberOfValidMuonHits") ) minNumberOfValidMuonHits_ = cc.getParameter<int>("minNumberOfValidMuonHits");
58  }
59 
60 }
61 
62 // Functors for evaluation
63 CutApplicatorBase::result_type MuonTrackCut::operator()(const reco::MuonPtr& muon) const
64 {
65  if ( doInnerTrack_ )
66  {
67  const reco::TrackRef t = muon->innerTrack();
68  if ( t.isNull() ) return false;
69  const auto& h = t->hitPattern();
70  if ( trackQuality_ != reco::Track::undefQuality and !t->quality(trackQuality_) ) return false;
71  if ( h.trackerLayersWithMeasurement() < minTrackerLayersWithMeasurement_ ) return false;
72  if ( h.pixelLayersWithMeasurement() < minPixelLayersWithMeasurement_ ) return false;
73  if ( h.numberOfValidPixelHits() < minNumberOfValidPixelHits_ ) return false;
74  if ( t->validFraction() <= minValidFraction_ ) return false;
75  }
76  if ( doGlobalTrack_ )
77  {
78  const reco::TrackRef t = muon->globalTrack();
79  if ( t.isNull() ) return false;
80  const auto& h = t->hitPattern();
81  if ( h.numberOfValidMuonHits() < minNumberOfValidMuonHits_ ) return false;
82  // if ( t->normalizedChi2() > maxNormalizedChi2_ ) return false; Not used for
83  }
84 
85  return true;
86 }
87 
88 double MuonTrackCut::value(const reco::CandidatePtr& cand) const
89 {
90  const reco::MuonPtr muon(cand);
91  if ( doInnerTrack_ )
92  {
93  const reco::TrackRef t = muon->innerTrack();
94  if ( t.isNull() ) return 0;
95  const auto& h = t->hitPattern();
96  if ( trackQuality_ != reco::Track::undefQuality and !t->quality(trackQuality_) ) return t->quality(trackQuality_);
97  if ( h.trackerLayersWithMeasurement() < minTrackerLayersWithMeasurement_ ) return h.trackerLayersWithMeasurement();
98  if ( h.pixelLayersWithMeasurement() < minPixelLayersWithMeasurement_ ) return h.pixelLayersWithMeasurement();
99  if ( h.numberOfValidPixelHits() < minNumberOfValidPixelHits_ ) return h.numberOfValidPixelHits();
100  if ( t->validFraction() <= minValidFraction_ ) return t->validFraction();
101 
102  return t->validFraction();
103  }
104  if ( doGlobalTrack_ )
105  {
106  const reco::TrackRef t = muon->globalTrack();
107  if ( t.isNull() ) return 0;
108  const auto& h = t->hitPattern();
109  if ( h.numberOfValidMuonHits() < minNumberOfValidMuonHits_ ) return h.numberOfValidMuonHits();
110  // if ( t->normalizedChi2() > maxNormalizedChi2_ ) return false; Not used for
111 
112  return h.numberOfValidMuonHits();
113  }
114 
115  return 0;
116 }
T getParameter(std::string const &) const
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
bool exists(std::string const &parameterName) const
checks if a parameter exists
int minNumberOfValidPixelHits_
Definition: MuonTrackCut.cc:16
result_type operator()(const reco::MuonPtr &) const overridefinal
Definition: MuonTrackCut.cc:63
MuonTrackCut(const edm::ParameterSet &c)
Definition: MuonTrackCut.cc:28
bool isNull() const
Checks for null.
Definition: Ref.h:249
double value(const reco::CandidatePtr &) const overridefinal
Definition: MuonTrackCut.cc:88
int minNumberOfValidMuonHits_
Definition: MuonTrackCut.cc:19
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:125
CandidateType candidateType() const overridefinal
Definition: MuonTrackCut.cc:11
int minTrackerLayersWithMeasurement_
Definition: MuonTrackCut.cc:16
This class analyses the reconstruction quality for a given track.
Definition: TrackQuality.h:28
string const
Definition: compareJSON.py:14
reco::Track::TrackQuality trackQuality_
Definition: MuonTrackCut.cc:18
double minValidFraction_
Definition: MuonTrackCut.cc:17
#define DEFINE_EDM_PLUGIN(factory, type, name)
volatile std::atomic< bool > shutdown_flag false
bool doInnerTrack_
Definition: MuonTrackCut.cc:21
int minPixelLayersWithMeasurement_
Definition: MuonTrackCut.cc:16
bool doGlobalTrack_
Definition: MuonTrackCut.cc:21