CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
pat::MuonMvaEstimator Class Reference

#include <MuonMvaEstimator.h>

Public Member Functions

float computeMva (const pat::Muon &imuon, const reco::Vertex &vertex, const reco::JetTagCollection &bTags, float &jetPtRatio, float &jetPtRel, const reco::JetCorrector *correctorL1=0, const reco::JetCorrector *correctorL1L2L3Res=0) const
 
 MuonMvaEstimator (const edm::FileInPath &weightsfile, float dRmax)
 
 ~MuonMvaEstimator ()
 

Private Attributes

float dRmax_
 
std::unique_ptr< const GBRForestgbrForest_
 

Detailed Description

Definition at line 25 of file MuonMvaEstimator.h.

Constructor & Destructor Documentation

MuonMvaEstimator::MuonMvaEstimator ( const edm::FileInPath weightsfile,
float  dRmax 
)

Definition at line 16 of file MuonMvaEstimator.cc.

References createGBRForest(), and gbrForest_.

16  :
17  dRmax_(dRmax)
18 {
19  gbrForest_ = createGBRForest( weightsfile );
20 }
std::unique_ptr< const GBRForest > gbrForest_
std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightsFile)
MuonMvaEstimator::~MuonMvaEstimator ( )

Definition at line 22 of file MuonMvaEstimator.cc.

22 { }

Member Function Documentation

float MuonMvaEstimator::computeMva ( const pat::Muon imuon,
const reco::Vertex vertex,
const reco::JetTagCollection bTags,
float &  jetPtRatio,
float &  jetPtRel,
const reco::JetCorrector correctorL1 = 0,
const reco::JetCorrector correctorL1L2L3Res = 0 
) const

Definition at line 54 of file MuonMvaEstimator.cc.

References reco::PFCandidate::bestTrack(), pat::Muon::BS2D, reco::LeafCandidate::charge(), pat::PFIsolation::chargedHadronIso(), reco::JetCorrector::correction(), pat::Muon::dB(), deltaR(), runTauDisplay::dr, dRmax_, PVValHelper::dz, pat::Muon::edB(), reco::LeafCandidate::eta(), Exception, gbrForest_, reco::TrackBase::highPurity, metsig::jet, JetComb::kEta, kPt, cmsBatch::log, SiStripPI::max, pat::Lepton< LeptonType >::miniPFIsolation(), pat::Muon::muonBestTrack(), pat::PFIsolation::neutralHadronIso(), reco::LeafCandidate::p4(), pfDeepBoostedJetPreprocessParams_cfi::pfcand, reco::Muon::pfIsolationR04(), reco::Vertex::position(), reco::LeafCandidate::pt(), electrons_cff::ptRel, pat::Muon::PV3D, pat::Muon::segmentCompatibility(), reco::MuonPFIsolation::sumChargedHadronPt, reco::MuonPFIsolation::sumNeutralHadronEt, reco::MuonPFIsolation::sumPhotonEt, reco::MuonPFIsolation::sumPUPt, and JetChargeProducer_cfi::var.

61 {
62  float var[kLast]{};
63 
64  var[kPt] = muon.pt();
65  var[kEta] = muon.eta();
66  var[kSegmentCompatibility] = muon.segmentCompatibility();
67  var[kMiniRelIsoCharged] = muon.miniPFIsolation().chargedHadronIso();
68  var[kMiniRelIsoNeutral] = muon.miniPFIsolation().neutralHadronIso();
69 
70  double dB2D = fabs(muon.dB(pat::Muon::BS2D));
71  double dB3D = muon.dB(pat::Muon::PV3D);
72  double edB3D = muon.edB(pat::Muon::PV3D);
73  double dz = fabs(muon.muonBestTrack()->dz(vertex.position()));
74  var[kSip] = edB3D>0?fabs(dB3D/edB3D):0.0;
75  var[kLog_abs_dxyBS] = dB2D>0?log(dB2D):0;
76  var[kLog_abs_dzPV] = dz>0?log(dz):0;
77 
78  //Initialise loop variables
79  double minDr = 9999;
80  double jecL1L2L3Res = 1.;
81  double jecL1 = 1.;
82 
83  // Compute corrected isolation variables
84  double chIso = muon.pfIsolationR04().sumChargedHadronPt;
85  double nIso = muon.pfIsolationR04().sumNeutralHadronEt;
86  double phoIso = muon.pfIsolationR04().sumPhotonEt;
87  double puIso = muon.pfIsolationR04().sumPUPt;
88  double dbCorrectedIsolation = chIso + std::max( nIso + phoIso - .5*puIso, 0. ) ;
89  double dbCorrectedRelIso = dbCorrectedIsolation/muon.pt();
90 
91  var[kJetPtRatio] = 1./(1+dbCorrectedRelIso);
92  var[kJetPtRel] = 0;
93  var[kJetBTagCSV] = -999;
94  var[kJetNDauCharged] = -1;
95 
96 
97  for (const auto& tagI: bTags){
98  // for each muon with the lepton
99  double dr = deltaR(*(tagI.first), muon);
100  if(dr > minDr) continue;
101  minDr = dr;
102 
103  const reco::Candidate::LorentzVector& muP4(muon.p4());
104  reco::Candidate::LorentzVector jetP4(tagI.first->p4());
105 
106  if (correctorL1 && correctorL1L2L3Res){
107  jecL1L2L3Res = correctorL1L2L3Res->correction(*(tagI.first));
108  jecL1 = correctorL1->correction(*(tagI.first));
109  }
110 
111  // Get b-jet info
112  var[kJetBTagCSV] = tagI.second;
113  var[kJetNDauCharged] = 0;
114  for (auto jet: tagI.first->getJetConstituentsQuick()){
115  const reco::PFCandidate *pfcand = dynamic_cast<const reco::PFCandidate*>(jet);
116  if (pfcand==nullptr) throw cms::Exception("ConfigurationError") << "Cannot get jet constituents";
117  if (pfcand->charge()==0) continue;
118  auto bestTrackPtr = pfcand->bestTrack();
119  if (!bestTrackPtr) continue;
120  if (!bestTrackPtr->quality(reco::Track::highPurity)) continue;
121  if (bestTrackPtr->pt()<1.) continue;
122  if (bestTrackPtr->hitPattern().numberOfValidHits()<8) continue;
123  if (bestTrackPtr->hitPattern().numberOfValidPixelHits()<2) continue;
124  if (bestTrackPtr->normalizedChi2()>=5) continue;
125 
126  if (std::fabs(bestTrackPtr->dxy(vertex.position())) > 0.2) continue;
127  if (std::fabs(bestTrackPtr->dz(vertex.position())) > 17) continue;
128  var[kJetNDauCharged]++;
129  }
130 
131  if(minDr < dRmax_){
132  if ((jetP4-muP4).Rho()<0.0001){
133  var[kJetPtRel] = 0;
134  var[kJetPtRatio] = 1;
135  } else {
136  jetP4 -= muP4/jecL1;
137  jetP4 *= jecL1L2L3Res;
138  jetP4 += muP4;
139 
140  var[kJetPtRatio] = muP4.pt()/jetP4.pt();
141  var[kJetPtRel] = ptRel(muP4,jetP4);
142  }
143  }
144  }
145 
146  if (var[kJetPtRatio] > 1.5) var[kJetPtRatio] = 1.5;
147  if (var[kJetBTagCSV] < 0) var[kJetBTagCSV] = 0;
148  jetPtRatio = var[kJetPtRatio];
149  jetPtRel = var[kJetPtRel];
150  return gbrForest_->GetClassifier(var);
151 };
example_track example_track const char *const kPt
Definition: TSelector.cc:34
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
const Point & position() const
position
Definition: Vertex.h:109
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
std::unique_ptr< const GBRForest > gbrForest_
const reco::Track * bestTrack() const override
Definition: PFCandidate.h:162

Member Data Documentation

float pat::MuonMvaEstimator::dRmax_
private

Definition at line 43 of file MuonMvaEstimator.h.

Referenced by computeMva().

std::unique_ptr<const GBRForest> pat::MuonMvaEstimator::gbrForest_
private

Definition at line 42 of file MuonMvaEstimator.h.

Referenced by computeMva(), and MuonMvaEstimator().