CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PATJetCorrExtractor.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_PatUtils_PATJetCorrExtractor_h
2 #define PhysicsTools_PatUtils_PATJetCorrExtractor_h
3 
21 
23 
28 
29 #include <string>
30 #include <vector>
31 
32 namespace
33 {
34  std::string format_vstring(const std::vector<std::string>& v)
35  {
36  std::string retVal;
37 
38  retVal.append("{ ");
39 
40  unsigned numEntries = v.size();
41  for ( unsigned iEntry = 0; iEntry < numEntries; ++iEntry ) {
42  retVal.append(v[iEntry]);
43  if ( iEntry < (numEntries - 1) ) retVal.append(", ");
44  }
45 
46  retVal.append(" }");
47 
48  return retVal;
49  }
50 }
51 
53 {
54  public:
55 
57  double jetCorrEtaMax = 9.9,
58  const reco::Candidate::LorentzVector * const rawJetP4_specified = nullptr) const
59  {
60  JetCorrExtractorT<pat::Jet> jetCorrExtractor;
61  return jetCorrExtractor(rawJet, jetCorr, jetCorrEtaMax, rawJetP4_specified);
62  }
63 
65  double jetCorrEtaMax = 9.9,
66  const reco::Candidate::LorentzVector* const rawJetP4_specified = nullptr) const
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  }
91 };
92 
93 #endif
94 
95 
const LorentzVector correctedP4(const std::string &level, const std::string &flavor="none", const std::string &set="") const
Definition: Jet.h:155
Long64_t numEntries(TFile *hdl, std::string const &trname)
Definition: CollUtil.cc:50
const std::vector< std::string > availableJECLevels(const int &set=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
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
std::string format_vstring(const std::vector< std::string > &vs)
Analysis-level calorimeter jet class.
Definition: Jet.h:77
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37