CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
L1FastjetCorrector.cc
Go to the documentation of this file.
1 //
3 // L1FastjetCorrector
4 // ------------------
5 //
6 // 08/09/2009 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
8 
10 
13 
14 
15 using namespace std;
16 
17 
19 // construction / destruction
21 
22 //______________________________________________________________________________
24  : srcRho_(fConfig.getParameter<edm::InputTag>("srcRho"))
25 {
26  if (fParam.definitions().level() != "L1FastJet")
27  throw cms::Exception("L1FastjetCorrector")<<" correction level: "<<fParam.definitions().level()<<" is not L1FastJet";
28  vector<JetCorrectorParameters> vParam;
29  vParam.push_back(fParam);
30  mCorrector = new FactorizedJetCorrector(vParam);
31 }
32 
33 //______________________________________________________________________________
35 {
36  delete mCorrector;
37 }
38 
39 
41 // implementation of member functions
43 
44 //______________________________________________________________________________
46 {
47  throw cms::Exception("EventRequired")
48  <<"Wrong interface correction(LorentzVector), event required!";
49  return 1.0;
50 }
51 
52 
53 //______________________________________________________________________________
54 double L1FastjetCorrector::correction (const reco::Jet& fJet) const
55 {
56  throw cms::Exception("EventRequired")
57  <<"Wrong interface correction(reco::Jet), event required!";
58  return 1.0;
59 }
60 
61 
62 //______________________________________________________________________________
64  const edm::Event& fEvent,
65  const edm::EventSetup& fSetup) const
66 {
68  fEvent.getByLabel(srcRho_,rho);
69  double result(1.0);
70  mCorrector->setJetEta(fJet.eta());
71  mCorrector->setJetPt(fJet.pt());
72  mCorrector->setJetA(fJet.jetArea());
73  mCorrector->setRho(*rho);
74  result = mCorrector->getCorrection();
75  return result;
76 }
77 
78 
79 
80 
virtual double correction(const LorentzVector &fJet) const
apply correction using Jet information only
FactorizedJetCorrector * mCorrector
Base class for all types of Jets.
Definition: Jet.h:21
const Definitions & definitions() const
Definition: DDAxes.h:10
virtual double eta() const
momentum pseudorapidity
tuple result
Definition: query.py:137
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
L1FastjetCorrector(const JetCorrectorParameters &fParam, const edm::ParameterSet &fConfig)
virtual double pt() const
transverse momentum
virtual float jetArea() const
get jet area
Definition: Jet.h:106
reco::Particle::LorentzVector LorentzVector
Definition: JetCorrector.h:24