CMS 3D CMS Logo

List of all members | Public Member Functions
PATJetCorrExtractor Class Reference

#include <PATJetCorrExtractor.h>

Public Member Functions

reco::Candidate::LorentzVector operator() (const pat::Jet rawJet, const reco::JetCorrector *jetCorr, double jetCorrEtaMax=9.9, const reco::Candidate::LorentzVector *const rawJetP4_specified=0) const
 
reco::Candidate::LorentzVector operator() (const pat::Jet &jet, const std::string &jetCorrLabel, double jetCorrEtaMax=9.9, const reco::Candidate::LorentzVector *const rawJetP4_specified=0) const
 

Detailed Description

Retrieve jet energy correction factor for pat::Jets (of either PF-type or Calo-type)

NOTE: this specialization of the "generic" template defined in JetMETCorrections/Type1MET/interface/JetCorrExtractorT.h is to be used for pat::Jets only

Author
Christian Veelken, LLR

Definition at line 51 of file PATJetCorrExtractor.h.

Member Function Documentation

reco::Candidate::LorentzVector PATJetCorrExtractor::operator() ( const pat::Jet  rawJet,
const reco::JetCorrector jetCorr,
double  jetCorrEtaMax = 9.9,
const reco::Candidate::LorentzVector *const  rawJetP4_specified = 0 
) const
inline

Definition at line 53 of file PATJetCorrExtractor.h.

References HLT_2018_cff::jetCorrEtaMax.

57  {
58  JetCorrExtractorT<pat::Jet> jetCorrExtractor;
59  return jetCorrExtractor(rawJet, jetCorr, jetCorrEtaMax, rawJetP4_specified);
60  }
reco::Candidate::LorentzVector PATJetCorrExtractor::operator() ( const pat::Jet jet,
const std::string &  jetCorrLabel,
double  jetCorrEtaMax = 9.9,
const reco::Candidate::LorentzVector *const  rawJetP4_specified = 0 
) const
inline

Definition at line 62 of file PATJetCorrExtractor.h.

References pat::Jet::availableJECLevels(), pat::Jet::correctedP4(), Exception, and format_vstring().

66  {
68 
69  try {
70  corrJetP4 = jet.correctedP4(jetCorrLabel);
71  if (rawJetP4_specified != nullptr) {
72  //MM: compensate for potential removal of constituents (as muons)
73  //similar effect in JetMETCorrection/Type1MET/interface/JetCorrExtractor.h
74  reco::Candidate::LorentzVector rawJetP4 = jet.correctedP4("Uncorrected");
75  double corrFactor = corrJetP4.pt() / rawJetP4.pt();
76  corrJetP4 = (*rawJetP4_specified);
77  corrJetP4 *= corrFactor;
78  if (corrFactor < 0) {
79  edm::LogWarning("PATJetCorrExtractor") << "Negative jet energy scale correction noticed"
80  << ".\n";
81  }
82  }
83  } catch (cms::Exception const&) {
84  throw cms::Exception("InvalidRequest")
85  << "The JEC level " << jetCorrLabel << " does not exist !!\n"
86  << "Available levels = " << format_vstring(jet.availableJECLevels()) << ".\n";
87  }
88 
89  return corrJetP4;
90  }
const LorentzVector correctedP4(const std::string &level, const std::string &flavor="none", const std::string &set="") const
Definition: Jet.h:164
const std::vector< std::string > availableJECLevels(const int &set=0) const
std::string format_vstring(const std::vector< std::string > &vs)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37