CMS 3D CMS Logo

MuonMvaEstimator.cc
Go to the documentation of this file.
2 
13 
14 using namespace pat;
15 
16 MuonMvaEstimator::MuonMvaEstimator(const edm::FileInPath& weightsfile, float dRmax):
17  dRmax_(dRmax)
18 {
19  gbrForest_ = createGBRForest( weightsfile );
20 }
21 
23 
24 namespace {
25 
26  enum inputIndexes {
27  kPt,
28  kEta,
29  kJetNDauCharged,
30  kMiniRelIsoCharged,
31  kMiniRelIsoNeutral,
32  kJetPtRel,
33  kJetBTagCSV,
34  kJetPtRatio,
35  kLog_abs_dxyBS,
36  kSip,
37  kLog_abs_dzPV,
38  kSegmentCompatibility,
39  kLast
40  };
41 
43  bool subtractMuon=true)
44  {
46  if (subtractMuon) jp4-=muP4;
47  float dot = muP4.Vect().Dot( jp4.Vect() );
48  float ptrel = muP4.P2() - dot*dot/jp4.P2();
49  ptrel = ptrel>0 ? sqrt(ptrel) : 0.0;
50  return ptrel;
51  }
52 }
53 
55  const reco::Vertex& vertex,
57  float& jetPtRatio,
58  float& jetPtRel,
59  float& miniIsoValue,
60  const reco::JetCorrector* correctorL1,
61  const reco::JetCorrector* correctorL1L2L3Res
62  ) const
63 {
64  float var[kLast]{};
65 
66  var[kPt] = muon.pt();
67  var[kEta] = muon.eta();
68  var[kSegmentCompatibility] = muon.segmentCompatibility();
69  var[kMiniRelIsoCharged] = muon.miniPFIsolation().chargedHadronIso() / muon.pt();
70  var[kMiniRelIsoNeutral] = miniIsoValue - var[kMiniRelIsoCharged];
71 
72  double dB2D = fabs(muon.dB(pat::Muon::PV2D));
73  double dB3D = muon.dB(pat::Muon::PV3D);
74  double edB3D = muon.edB(pat::Muon::PV3D);
75  double dz = fabs(muon.muonBestTrack()->dz(vertex.position()));
76  var[kSip] = edB3D>0?fabs(dB3D/edB3D):0.0;
77  var[kLog_abs_dxyBS] = dB2D>0?log(dB2D):0;
78  var[kLog_abs_dzPV] = dz>0?log(dz):0;
79 
80  //Initialise loop variables
81  double minDr = 9999;
82  double jecL1L2L3Res = 1.;
83  double jecL1 = 1.;
84 
85  // Compute corrected isolation variables
86  double chIso = muon.pfIsolationR04().sumChargedHadronPt;
87  double nIso = muon.pfIsolationR04().sumNeutralHadronEt;
88  double phoIso = muon.pfIsolationR04().sumPhotonEt;
89  double puIso = muon.pfIsolationR04().sumPUPt;
90  double dbCorrectedIsolation = chIso + std::max( nIso + phoIso - .5*puIso, 0. ) ;
91  double dbCorrectedRelIso = dbCorrectedIsolation/muon.pt();
92 
93  var[kJetPtRatio] = 1./(1+dbCorrectedRelIso);
94  var[kJetPtRel] = 0;
95  var[kJetBTagCSV] = -999;
96  var[kJetNDauCharged] = -1;
97 
98 
99  for (const auto& tagI: bTags){
100  // for each muon with the lepton
101  double dr = deltaR(*(tagI.first), muon);
102  if(dr > minDr) continue;
103  minDr = dr;
104 
105  const reco::Candidate::LorentzVector& muP4(muon.p4());
106  reco::Candidate::LorentzVector jetP4(tagI.first->p4());
107 
108  if (correctorL1 && correctorL1L2L3Res){
109  jecL1L2L3Res = correctorL1L2L3Res->correction(*(tagI.first));
110  jecL1 = correctorL1->correction(*(tagI.first));
111  }
112 
113  // Get b-jet info
114  var[kJetBTagCSV] = tagI.second;
115  var[kJetNDauCharged] = 0;
116  for (auto jet: tagI.first->getJetConstituentsQuick()){
117  const reco::PFCandidate *pfcand = dynamic_cast<const reco::PFCandidate*>(jet);
118  if (pfcand==nullptr) throw cms::Exception("ConfigurationError") << "Cannot get jet constituents";
119  if (pfcand->charge()==0) continue;
120  auto bestTrackPtr = pfcand->bestTrack();
121  if (!bestTrackPtr) continue;
122  if (!bestTrackPtr->quality(reco::Track::highPurity)) continue;
123  if (bestTrackPtr->pt()<1.) continue;
124  if (bestTrackPtr->hitPattern().numberOfValidHits()<8) continue;
125  if (bestTrackPtr->hitPattern().numberOfValidPixelHits()<2) continue;
126  if (bestTrackPtr->normalizedChi2()>=5) continue;
127 
128  if (std::fabs(bestTrackPtr->dxy(vertex.position())) > 0.2) continue;
129  if (std::fabs(bestTrackPtr->dz(vertex.position())) > 17) continue;
130  var[kJetNDauCharged]++;
131  }
132 
133  if(minDr < dRmax_){
134  if ((jetP4-muP4).Rho()<0.0001){
135  var[kJetPtRel] = 0;
136  var[kJetPtRatio] = 1;
137  } else {
138  jetP4 -= muP4/jecL1;
139  jetP4 *= jecL1L2L3Res;
140  jetP4 += muP4;
141 
142  var[kJetPtRatio] = muP4.pt()/jetP4.pt();
143  var[kJetPtRel] = ptRel(muP4,jetP4);
144  }
145  }
146  }
147 
148  if (var[kJetPtRatio] > 1.5) var[kJetPtRatio] = 1.5;
149  if (var[kJetBTagCSV] < 0) var[kJetBTagCSV] = 0;
150  jetPtRatio = var[kJetPtRatio];
151  jetPtRel = var[kJetPtRel];
152  return gbrForest_->GetClassifier(var);
153 };
double eta() const final
momentum pseudorapidity
example_track example_track const char *const kPt
Definition: TSelector.cc:34
float computeMva(const pat::Muon &imuon, const reco::Vertex &vertex, const reco::JetTagCollection &bTags, float &jetPtRatio, float &jetPtRel, float &miniIsoValue, const reco::JetCorrector *correctorL1=0, const reco::JetCorrector *correctorL1L2L3Res=0) const
double dB(IPTYPE type) const
double pt() const final
transverse momentum
int charge() const final
electric charge
Definition: LeafCandidate.h:91
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:49
MuonMvaEstimator(const edm::FileInPath &weightsfile, float dRmax)
float sumPhotonEt
sum pt of PF photons
const Point & position() const
position
Definition: Vertex.h:109
float sumNeutralHadronEt
sum pt of neutral hadrons
Definition: HeavyIon.h:7
const PFIsolation & miniPFIsolation() const
Definition: Lepton.h:196
inputIndexes
double segmentCompatibility(reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration) const
Returns the segment compatibility, using muon::segmentCompatibility (DataFormats/MuonReco/interface/M...
float sumPUPt
sum pt of charged Particles not from PV (for Pu corrections)
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:660
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
double edB(IPTYPE type) const
const MuonPFIsolation & pfIsolationR04() const
Definition: Muon.h:168
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
T dot(const Basic3DVector &v) const
Scalar product, or "dot" product, with a vector of same type.
std::unique_ptr< const GBRForest > gbrForest_
reco::TrackRef muonBestTrack() const override
Track selected to be the best measurement of the muon parameters (including PFlow global information)...
float chargedHadronIso() const
Definition: PFIsolation.h:33
Analysis-level muon class.
Definition: Muon.h:51
const reco::Track * bestTrack() const override
Definition: PFCandidate.h:162
std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightsFile)
float sumChargedHadronPt
sum-pt of charged Hadron