CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetCorrExtractorT.h
Go to the documentation of this file.
1 #ifndef JetMETCorrections_Type1MET_JetCorrExtractorT_h
2 #define JetMETCorrections_Type1MET_JetCorrExtractorT_h
3 
23 
28 
29 namespace jetcorrextractor
30 {
31  // never heard of copysign?
32  inline double sign(double x)
33  {
34  if ( x > 0. ) return +1.;
35  else if ( x < 0. ) return -1.;
36  else return 0.;
37  }
38 }
39 
40 template <typename T>
42 {
43  public:
44 
46  double jetCorrEtaMax = 9.9,
47  const reco::Candidate::LorentzVector * const rawJetP4_specified = nullptr) const
48  {
49 
50  // allow to specify four-vector to be used as "raw" (uncorrected) jet momentum,
51  // call 'rawJet.p4()' in case four-vector not specified explicitely
52  reco::Candidate::LorentzVector rawJetP4 = ( rawJetP4_specified ) ?
53  (*rawJetP4_specified) : rawJet.p4();
54 
55  double jetCorrFactor = 1.;
56  if ( fabs(rawJetP4.eta()) < jetCorrEtaMax ) {
57  jetCorrFactor = getCorrection(rawJet, jetCorr);
58  } else {
59  reco::Candidate::PolarLorentzVector modJetPolarP4(rawJetP4);
60  modJetPolarP4.SetEta(jetcorrextractor::sign(rawJetP4.eta())*jetCorrEtaMax);
61 
62  reco::Candidate::LorentzVector modJetP4(modJetPolarP4);
63 
64  T modJet(rawJet);
65  modJet.setP4(modJetP4);
66 
67  jetCorrFactor = getCorrection(modJet, jetCorr);
68  if(jetCorrFactor<0) {
69  edm::LogWarning("JetCorrExtractor") << "Negative jet energy scale correction noticed" << ".\n";
70  }
71  }
72 
73  reco::Candidate::LorentzVector corrJetP4 = rawJetP4;
74  corrJetP4 *= jetCorrFactor;
75 
76  return corrJetP4;
77  }
78 
79  reco::Candidate::LorentzVector operator()(const T& rawJet, const std::string& jetCorrLabel,
80  double jetCorrEtaMax = 9.9,
81  const reco::Candidate::LorentzVector * const rawJetP4_specified = nullptr) const
82  {
83  edm::LogWarning("JetCorrExtractor") << "JetCorrExtractorT<T>::operator(const T&, const std::string&, ...) is deprecated.\n"
84  << "Please use JetCorrExtractorT<T>::operator(const T&, const reco::JetCorrector*, ...) instead.\n"
85  << "Jet remains uncorrected!";
86  reco::Candidate::LorentzVector rawJetP4 = ( rawJetP4_specified ) ?
87  (*rawJetP4_specified) : rawJet.p4();
88  return rawJetP4;
89  }
90 
91  private:
92 
93  static double getCorrection(const T& jet, const reco::JetCorrector* jetCorr)
94  {
95  return jetCorr->correction(jet);
96  }
97 
98 
99 };
100 
101 #endif
reco::Candidate::LorentzVector operator()(const T &rawJet, const reco::JetCorrector *jetCorr, double jetCorrEtaMax=9.9, const reco::Candidate::LorentzVector *const rawJetP4_specified=0) const
double sign(double x)
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:47
static double getCorrection(const T &jet, const reco::JetCorrector *jetCorr)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
reco::Candidate::LorentzVector operator()(const T &rawJet, const std::string &jetCorrLabel, double jetCorrEtaMax=9.9, const reco::Candidate::LorentzVector *const rawJetP4_specified=0) const
long double T
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:39