CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LXXXCorrector.cc
Go to the documentation of this file.
1 // Implementation of class LXXXCorrector.
2 // Generic LX jet corrector class.
3 
6 
13 
14 using namespace std;
15 
16 
17 //------------------------------------------------------------------------
18 //--- LXXXCorrector constructor ------------------------------------------
19 //------------------------------------------------------------------------
21 {
22  string level = fParam.definitions().level();
23  if (level == "L2Relative")
24  mLevel = 2;
25  else if (level == "L3Absolute")
26  mLevel = 3;
27  else if (level == "L4EMF")
28  mLevel = 4;
29  else if (level == "L5Flavor")
30  mLevel = 5;
31  else if (level == "L7Parton")
32  mLevel = 7;
33  else
34  throw cms::Exception("LXXXCorrector")<<" unknown correction level "<<level;
35  vector<JetCorrectorParameters> vParam;
36  vParam.push_back(fParam);
37  mCorrector = new FactorizedJetCorrectorCalculator(vParam);
38 }
39 //------------------------------------------------------------------------
40 //--- LXXXCorrector destructor -------------------------------------------
41 //------------------------------------------------------------------------
43 {
44  delete mCorrector;
45 }
46 //------------------------------------------------------------------------
47 //--- Returns correction for a given 4-vector ----------------------------
48 //------------------------------------------------------------------------
49 double LXXXCorrector::correction(const LorentzVector& fJet) const
50 {
51  // L4 correction requires more information than a simple 4-vector
52  if (mLevel == 4) {
53  throw cms::Exception("Invalid jet type") << "L4EMFCorrection is applicable to CaloJets only";
54  return 1;
55  }
56 
58  values.setJetEta(fJet.eta());
59  values.setJetE(fJet.energy());
60  values.setJetPt(fJet.pt());
61  values.setJetPhi(fJet.phi());
62 
63  return mCorrector->getCorrection(values);
64 }
65 //------------------------------------------------------------------------
66 //--- Returns correction for a given jet ---------------------------------
67 //------------------------------------------------------------------------
68 double LXXXCorrector::correction(const reco::Jet& fJet) const
69 {
70  double result = 1.;
71  // L4 correction applies to Calojets only
72  if (mLevel == 4) {
73  const reco::CaloJet& caloJet = dynamic_cast <const reco::CaloJet&> (fJet);
75  values.setJetEta(fJet.eta());
76  values.setJetPt(fJet.pt());
77  values.setJetEMF(caloJet.emEnergyFraction());
78  result = mCorrector->getCorrection(values);
79  }
80  else
81  result = correction(fJet.p4());
82  return result;
83 }
LXXXCorrector(const JetCorrectorParameters &fConfig, const edm::ParameterSet &fParameters)
Jets made from CaloTowers.
Definition: CaloJet.h:29
Base class for all types of Jets.
Definition: Jet.h:20
const Definitions & definitions() const
virtual double correction(const LorentzVector &fJet) const
get correction using Jet information only
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
tuple result
Definition: query.py:137
virtual ~LXXXCorrector()
tuple level
Definition: testEve_cfg.py:34
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
float emEnergyFraction() const
Definition: CaloJet.h:102
reco::Particle::LorentzVector LorentzVector
Definition: JetCorrector.h:23