CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LXXXCorrectorImpl.cc
Go to the documentation of this file.
1 // Implementation of class LXXXCorrectorImpl.
2 // Generic LX jet corrector class.
3 
6 
15 
16 using namespace std;
17 
18 
21 {
22 }
23 
24 std::unique_ptr<reco::JetCorrectorImpl>
26 {
27  unsigned int level =0;
28  auto calculator = getCalculator(fSetup,
29  [&level](std::string const& levelName)
30  {
31  if (levelName == "L2Relative")
32  level = 2;
33  else if (levelName == "L3Absolute")
34  level = 3;
35  else if (levelName == "L4EMF")
36  level = 4;
37  else if (levelName == "L5Flavor")
38  level = 5;
39  else if (levelName == "L7Parton")
40  level = 7;
41  else
42  throw cms::Exception("LXXXCorrectorImpl")<<" unknown correction level "<<levelName;
43  });
44  return std::unique_ptr<reco::JetCorrectorImpl>(new LXXXCorrectorImpl(calculator,level));
45 }
46 
47 void
49 {
51  addToDescription(desc);
52 
53  iDescriptions.addDefault(desc);
54 }
55 
56 //------------------------------------------------------------------------
57 //--- LXXXCorrectorImpl constructor ------------------------------------------
58 //------------------------------------------------------------------------
59 LXXXCorrectorImpl::LXXXCorrectorImpl(std::shared_ptr<FactorizedJetCorrectorCalculator const> calculator, unsigned int level):
60  mLevel(level),
61  mCorrector(calculator)
62 {
63 }
64 //------------------------------------------------------------------------
65 //--- Returns correction for a given 4-vector ----------------------------
66 //------------------------------------------------------------------------
67 double LXXXCorrectorImpl::correction(const LorentzVector& fJet) const
68 {
69  // L4 correction requires more information than a simple 4-vector
70  if (mLevel == 4) {
71  throw cms::Exception("Invalid jet type") << "L4EMFCorrection is applicable to CaloJets only";
72  return 1;
73  }
74 
76  values.setJetEta(fJet.eta());
77  values.setJetE(fJet.energy());
78  values.setJetPt(fJet.pt());
79  values.setJetPhi(fJet.phi());
80 
81  return mCorrector->getCorrection(values);
82 }
83 //------------------------------------------------------------------------
84 //--- Returns correction for a given jet ---------------------------------
85 //------------------------------------------------------------------------
86 double LXXXCorrectorImpl::correction(const reco::Jet& fJet) const
87 {
88  double result = 1.;
89  // L4 correction applies to Calojets only
90  if (mLevel == 4) {
91  const reco::CaloJet& caloJet = dynamic_cast <const reco::CaloJet&> (fJet);
93  values.setJetEta(fJet.eta());
94  values.setJetPt(fJet.pt());
95  values.setJetEMF(caloJet.emEnergyFraction());
96  result = mCorrector->getCorrection(values);
97  }
98  else
99  result = correction(fJet.p4());
100  return result;
101 }
Jets made from CaloTowers.
Definition: CaloJet.h:29
LXXXCorrectorImpl(std::shared_ptr< FactorizedJetCorrectorCalculator const > calculator, unsigned int level)
Base class for all types of Jets.
Definition: Jet.h:20
reco::Particle::LorentzVector LorentzVector
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
LXXXCorrectorImplMaker(edm::ParameterSet const &, edm::ConsumesCollector)
static void fillDescriptions(edm::ConfigurationDescriptions &iDescriptions)
void addDefault(ParameterSetDescription const &psetDescription)
tuple result
Definition: query.py:137
const char * levelName(LogLevel)
Definition: fwLog.cc:34
std::unique_ptr< reco::JetCorrectorImpl > make(edm::Event const &, edm::EventSetup const &)
virtual double correction(const LorentzVector &fJet) const override
get correction using Jet information only
static void addToDescription(edm::ParameterSetDescription &iDescription)
tuple level
Definition: testEve_cfg.py:34
std::shared_ptr< FactorizedJetCorrectorCalculator const > mCorrector
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
std::shared_ptr< FactorizedJetCorrectorCalculator const > getCalculator(edm::EventSetup const &, std::function< void(std::string const &)> levelCheck)
float emEnergyFraction() const
Definition: CaloJet.h:102