CMS 3D CMS Logo

SoftMuonMvaEstimator.cc
Go to the documentation of this file.
2 
9 
10 #include "TMath.h"
11 
12 using namespace pat;
13 
15  gbrForest_ = createGBRForest(weightsfile);
16 }
17 
19 
20 namespace {
21  enum inputIndexes {
22  kSegmentCompatibility,
23  kChi2LocalMomentum,
24  kChi2LocalPosition,
25  kGlbTrackProbability,
26  kIValidFraction,
27  kLayersWithMeasurement,
28  kTrkKink,
29  kLog2PlusGlbKink,
30  kTimeAtIpInOutErr,
31  kOuterChi2,
32  kInnerChi2,
33  kTrkRelChi2,
34  kVMuonHitComb,
35  kQProd,
36  kPID,
37  kPt,
38  kEta,
39  kMomID,
40  kLast
41  };
42 }
43 
45  float var[kLast]{};
46 
47  reco::TrackRef gTrack = muon.globalTrack();
48  reco::TrackRef iTrack = muon.innerTrack();
49  reco::TrackRef oTrack = muon.outerTrack();
50 
51  if (!(muon.innerTrack().isNonnull() and muon.outerTrack().isNonnull() and muon.globalTrack().isNonnull())) {
52  return -1;
53  }
54 
55  //VARIABLE EXTRACTION
56  var[kPt] = muon.pt();
57  var[kEta] = muon.eta();
58  var[kMomID] = -1;
59  var[kPID] = -1;
60 
61  var[kChi2LocalMomentum] = muon.combinedQuality().chi2LocalMomentum;
62  var[kChi2LocalPosition] = muon.combinedQuality().chi2LocalPosition;
63  var[kGlbTrackProbability] = muon.combinedQuality().glbTrackProbability;
64  var[kTrkRelChi2] = muon.combinedQuality().trkRelChi2;
65 
66  var[kTrkKink] = muon.combinedQuality().trkKink;
67  var[kLog2PlusGlbKink] = TMath::Log(2 + muon.combinedQuality().glbKink);
68  var[kSegmentCompatibility] = muon.segmentCompatibility();
69 
70  var[kTimeAtIpInOutErr] = muon.time().timeAtIpInOutErr;
71 
72  //TRACK RELATED VARIABLES
73 
74  var[kIValidFraction] = iTrack->validFraction();
75  var[kInnerChi2] = iTrack->normalizedChi2();
76  var[kLayersWithMeasurement] = iTrack->hitPattern().trackerLayersWithMeasurement();
77 
78  var[kOuterChi2] = oTrack->normalizedChi2();
79 
80  var[kQProd] = iTrack->charge() * oTrack->charge();
81 
82  //vComb Calculation
83 
84  const reco::HitPattern& gMpattern = gTrack->hitPattern();
85 
86  std::vector<int> fvDThits{0, 0, 0, 0};
87  std::vector<int> fvRPChits{0, 0, 0, 0};
88  std::vector<int> fvCSChits{0, 0, 0, 0};
89 
90  var[kVMuonHitComb] = 0;
91 
92  for (int i = 0; i < gMpattern.numberOfAllHits(reco::HitPattern::TRACK_HITS); i++) {
93  uint32_t hit = gMpattern.getHitPattern(reco::HitPattern::TRACK_HITS, i);
94  if (!gMpattern.validHitFilter(hit))
95  continue;
96 
97  int muStation0 = gMpattern.getMuonStation(hit) - 1;
98  if (muStation0 >= 0 && muStation0 < 4) {
99  if (gMpattern.muonDTHitFilter(hit))
100  fvDThits[muStation0]++;
101  if (gMpattern.muonRPCHitFilter(hit))
102  fvRPChits[muStation0]++;
103  if (gMpattern.muonCSCHitFilter(hit))
104  fvCSChits[muStation0]++;
105  }
106  }
107 
108  for (unsigned int station = 0; station < 4; ++station) {
109  var[kVMuonHitComb] += (fvDThits[station]) / 2.;
110  var[kVMuonHitComb] += fvRPChits[station];
111 
112  if (fvCSChits[station] > 6) {
113  var[kVMuonHitComb] += 6;
114  } else {
115  var[kVMuonHitComb] += fvCSChits[station];
116  }
117  }
118 
119  if (var[kChi2LocalMomentum] < 5000 and var[kChi2LocalPosition] < 2000 and var[kGlbTrackProbability] < 5000 and
120  var[kTrkKink] < 900 and var[kLog2PlusGlbKink] < 50 and var[kTimeAtIpInOutErr] < 4 and var[kOuterChi2] < 1000 and
121  var[kInnerChi2] < 10 and var[kTrkRelChi2] < 3) {
122  return gbrForest_->GetAdaBoostClassifier(var);
123  } else {
124  return -1;
125  }
126 }
kPt
example_track example_trackconst char *const kPt
Definition: TSelector.cc:33
reco::HitPattern::getHitPattern
uint16_t getHitPattern(HitCategory category, int position) const
Definition: HitPattern.h:530
mps_fire.i
i
Definition: mps_fire.py:428
Muon.h
muon
Definition: MuonCocktails.h:17
reco::HitPattern::muonRPCHitFilter
static bool muonRPCHitFilter(uint16_t pattern)
Definition: HitPattern.h:640
relativeConstraints.station
station
Definition: relativeConstraints.py:67
Muon.h
GBRForestTools.h
pat::SoftMuonMvaEstimator::gbrForest_
std::unique_ptr< const GBRForest > gbrForest_
Definition: SoftMuonMvaEstimator.h:27
reco::HitPattern::validHitFilter
static bool validHitFilter(uint16_t pattern)
Definition: HitPattern.h:765
createGBRForest
std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightsFile)
Definition: GBRForestTools.cc:244
reco::HitPattern::muonCSCHitFilter
static bool muonCSCHitFilter(uint16_t pattern)
Definition: HitPattern.h:632
pat::Muon
Analysis-level muon class.
Definition: Muon.h:51
edm::Ref< TrackCollection >
FileInPath.h
trigObjTnPSource_cfi.var
var
Definition: trigObjTnPSource_cfi.py:21
edm::FileInPath
Definition: FileInPath.h:64
reco::HitPattern
Definition: HitPattern.h:147
MuonSelectors.h
SoftMuonMvaEstimator.h
pat::SoftMuonMvaEstimator::~SoftMuonMvaEstimator
~SoftMuonMvaEstimator()
Definition: SoftMuonMvaEstimator.cc:18
reco::HitPattern::getMuonStation
static uint16_t getMuonStation(uint16_t pattern)
Muon station (1-4). Only valid for muon patterns, of course. only for patterns from muon,...
Definition: HitPattern.h:732
pat
Definition: HeavyIon.h:7
reco::HitPattern::TRACK_HITS
Definition: HitPattern.h:155
pat::SoftMuonMvaEstimator::computeMva
float computeMva(const pat::Muon &imuon) const
Definition: SoftMuonMvaEstimator.cc:44
reco::HitPattern::numberOfAllHits
int numberOfAllHits(HitCategory category) const
Definition: HitPattern.h:773
pat::SoftMuonMvaEstimator::SoftMuonMvaEstimator
SoftMuonMvaEstimator(const edm::FileInPath &weightsfile)
Definition: SoftMuonMvaEstimator.cc:14
inputIndexes
inputIndexes
Definition: MuonMvaEstimator.cc:24
Candidate.h
JetComb::kEta
Definition: TtSemiLepJetComb.h:17
reco::HitPattern::muonDTHitFilter
static bool muonDTHitFilter(uint16_t pattern)
Definition: HitPattern.h:624
hit
Definition: SiStripHitEffFromCalibTree.cc:88