CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
L1FastjetCorrector Class Reference

#include <L1FastjetCorrector.h>

Inheritance diagram for L1FastjetCorrector:
JetCorrector

Public Member Functions

double correction (const LorentzVector &fJet) const override
 apply correction using Jet information only More...
 
double correction (const reco::Jet &fJet) const override
 apply correction using Jet information only More...
 
double correction (const reco::Jet &fJet, const edm::Event &fEvent, const edm::EventSetup &fSetup) const override
 apply correction using all event information More...
 
bool eventRequired () const override
 if correction needs event information More...
 
 L1FastjetCorrector (const JetCorrectorParameters &fParam, const edm::ParameterSet &fConfig)
 
bool refRequired () const override
 if correction needs the jet reference More...
 
 ~L1FastjetCorrector () override
 
- Public Member Functions inherited from JetCorrector
virtual double correction (const reco::Jet &fJet, const edm::RefToBase< reco::Jet > &fJetRef, const edm::Event &fEvent, const edm::EventSetup &fSetup) const
 apply correction using all event information More...
 
virtual double correction (const reco::Jet &fJet, const edm::RefToBase< reco::Jet > &fJetRef, const edm::Event &fEvent, const edm::EventSetup &fSetup, LorentzVector &corrected) const
 Apply vectorial correction using all event information. More...
 
 JetCorrector ()
 
virtual bool vectorialCorrection () const
 if vectorial correction is provided More...
 
virtual ~JetCorrector ()
 

Private Attributes

FactorizedJetCorrectorCalculator const * mCorrector
 
edm::InputTag srcRho_
 

Additional Inherited Members

- Public Types inherited from JetCorrector
typedef reco::Particle::LorentzVector LorentzVector
 
- Static Public Member Functions inherited from JetCorrector
static const JetCorrectorgetJetCorrector (const std::string &fName, const edm::EventSetup &fSetup)
 retrieve corrector from the event setup. troughs exception if something is missing More...
 

Detailed Description

Definition at line 17 of file L1FastjetCorrector.h.

Constructor & Destructor Documentation

L1FastjetCorrector::L1FastjetCorrector ( const JetCorrectorParameters fParam,
const edm::ParameterSet fConfig 
)

Definition at line 21 of file L1FastjetCorrector.cc.

References JetCorrectorParameters::definitions(), JetCorrectorParameters::Definitions::level(), and mCorrector.

22  : srcRho_(fConfig.getParameter<edm::InputTag>("srcRho")) {
23  if (fParam.definitions().level() != "L1FastJet")
24  throw cms::Exception("L1FastjetCorrector")
25  << " correction level: " << fParam.definitions().level() << " is not L1FastJet";
26  vector<JetCorrectorParameters> vParam;
27  vParam.push_back(fParam);
29 }
T getParameter(std::string const &) const
const Definitions & definitions() const
FactorizedJetCorrectorCalculator const * mCorrector
L1FastjetCorrector::~L1FastjetCorrector ( )
override

Definition at line 32 of file L1FastjetCorrector.cc.

References mCorrector.

32 { delete mCorrector; }
FactorizedJetCorrectorCalculator const * mCorrector

Member Function Documentation

double L1FastjetCorrector::correction ( const LorentzVector fJet) const
overridevirtual

apply correction using Jet information only

Implements JetCorrector.

Definition at line 39 of file L1FastjetCorrector.cc.

References Exception.

39  {
40  throw cms::Exception("EventRequired") << "Wrong interface correction(LorentzVector), event required!";
41  return 1.0;
42 }
double L1FastjetCorrector::correction ( const reco::Jet fJet) const
overridevirtual

apply correction using Jet information only

Implements JetCorrector.

Definition at line 45 of file L1FastjetCorrector.cc.

References Exception.

45  {
46  throw cms::Exception("EventRequired") << "Wrong interface correction(reco::Jet), event required!";
47  return 1.0;
48 }
double L1FastjetCorrector::correction ( const reco::Jet fJet,
const edm::Event fEvent,
const edm::EventSetup fSetup 
) const
overridevirtual

apply correction using all event information

Reimplemented from JetCorrector.

Definition at line 51 of file L1FastjetCorrector.cc.

References reco::LeafCandidate::energy(), reco::LeafCandidate::eta(), edm::Event::getByLabel(), FactorizedJetCorrectorCalculator::getCorrection(), reco::Jet::jetArea(), mCorrector, reco::LeafCandidate::pt(), mps_fire::result, rho, FactorizedJetCorrectorCalculator::VariableValues::setJetA(), FactorizedJetCorrectorCalculator::VariableValues::setJetE(), FactorizedJetCorrectorCalculator::VariableValues::setJetEta(), FactorizedJetCorrectorCalculator::VariableValues::setJetPt(), FactorizedJetCorrectorCalculator::VariableValues::setRho(), srcRho_, and contentValuesCheck::values.

53  {
55  fEvent.getByLabel(srcRho_, rho);
56  double result(1.0);
58  values.setJetEta(fJet.eta());
59  values.setJetPt(fJet.pt());
60  values.setJetE(fJet.energy());
61  values.setJetA(fJet.jetArea());
62  values.setRho(*rho);
63  result = mCorrector->getCorrection(values);
64  return result;
65 }
double eta() const final
momentum pseudorapidity
double pt() const final
transverse momentum
double energy() const final
energy
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:488
virtual float jetArea() const
get jet area
Definition: Jet.h:103
FactorizedJetCorrectorCalculator const * mCorrector
bool L1FastjetCorrector::eventRequired ( ) const
inlineoverridevirtual

if correction needs event information

Implements JetCorrector.

Definition at line 33 of file L1FastjetCorrector.h.

33 { return true; }
bool L1FastjetCorrector::refRequired ( ) const
inlineoverridevirtual

if correction needs the jet reference

Implements JetCorrector.

Definition at line 36 of file L1FastjetCorrector.h.

36 { return false; }

Member Data Documentation

FactorizedJetCorrectorCalculator const* L1FastjetCorrector::mCorrector
private

Definition at line 41 of file L1FastjetCorrector.h.

Referenced by correction(), L1FastjetCorrector(), and ~L1FastjetCorrector().

edm::InputTag L1FastjetCorrector::srcRho_
private

Definition at line 40 of file L1FastjetCorrector.h.

Referenced by correction().