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 52 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 56 of file PATJetCorrExtractor.h.

References patCaloMETCorrections_cff::jetCorrEtaMax.

59  {
60  JetCorrExtractorT<pat::Jet> jetCorrExtractor;
61  return jetCorrExtractor(rawJet, jetCorr, jetCorrEtaMax, rawJetP4_specified);
62  }
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 64 of file PATJetCorrExtractor.h.

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

67  {
69 
70  try {
71  corrJetP4 = jet.correctedP4(jetCorrLabel);
72  if(rawJetP4_specified != nullptr) {
73  //MM: compensate for potential removal of constituents (as muons)
74  //similar effect in JetMETCorrection/Type1MET/interface/JetCorrExtractor.h
75  reco::Candidate::LorentzVector rawJetP4 = jet.correctedP4("Uncorrected");
76  double corrFactor = corrJetP4.pt()/rawJetP4.pt();
77  corrJetP4 = (*rawJetP4_specified);
78  corrJetP4 *= corrFactor;
79  if(corrFactor<0) {
80  edm::LogWarning("PATJetCorrExtractor") << "Negative jet energy scale correction noticed" << ".\n";
81  }
82  }
83  } catch( cms::Exception e ) {
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:156
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