CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MuonSegmentCompatibilityCut.cc
Go to the documentation of this file.
4 
6 public:
8 
9  result_type operator()(const reco::MuonPtr&) const final;
10  CandidateType candidateType() const final { return MUON; }
11  double value(const reco::CandidatePtr&) const final;
12 
13 private:
16 };
18 
19 // Define constructors and initialization routines
20 MuonSegmentCompatibilityCut::MuonSegmentCompatibilityCut(const edm::ParameterSet& c)
21  : CutApplicatorBase(c),
22  minCompatGlb_(c.getParameter<double>("minCompatGlb")),
23  minCompatNonGlb_(c.getParameter<double>("minCompatNonGlb")) {
24  const edm::ParameterSet cc = c.getParameter<edm::ParameterSet>("goodGLB");
25  maxGlbNormChi2_ = cc.getParameter<double>("maxGlbNormChi2");
26  maxChi2LocalPos_ = cc.getParameter<double>("maxChi2LocalPos");
27  maxTrkKink_ = cc.getParameter<double>("maxTrkKink");
28 }
29 
30 // Functors for evaluation
32  const bool isGoodGlb =
33  (muon->isGlobalMuon() and muon->globalTrack()->normalizedChi2() < maxGlbNormChi2_ and
34  muon->combinedQuality().chi2LocalPosition < maxChi2LocalPos_ and muon->combinedQuality().trkKink < maxTrkKink_);
35 
36  const double compat = muon::segmentCompatibility(*muon);
37 
38  return compat > (isGoodGlb ? minCompatGlb_ : minCompatNonGlb_);
39 }
40 
42  const reco::MuonPtr muon(cand);
43  return muon::segmentCompatibility(*muon);
44 }
const edm::EventSetup & c
MuonSegmentCompatibilityCut(const edm::ParameterSet &c)
float segmentCompatibility(const reco::Muon &muon, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
double value(const reco::CandidatePtr &) const final
CandidateType candidateType() const final
result_type operator()(const reco::MuonPtr &) const final
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
#define DEFINE_EDM_PLUGIN(factory, type, name)