CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Member Functions
JetCorrector Class Referenceabstract

#include <JetCorrector.h>

Inheritance diagram for JetCorrector:
ChainedJetCorrector L1FastjetCorrector L1JPTOffsetCorrector L1OffsetCorrector L6SLBCorrector LXXXCorrector

Public Types

typedef reco::Particle::LorentzVector LorentzVector
 

Public Member Functions

virtual double correction (const LorentzVector &fJet) const =0
 get correction using Jet information only More...
 
virtual double correction (const reco::Jet &fJet) const =0
 apply correction using Jet information only More...
 
virtual double correction (const reco::Jet &fJet, 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) 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...
 
virtual bool eventRequired () const =0
 if correction needs event information More...
 
 JetCorrector ()
 
virtual bool refRequired () const =0
 if correction needs the jet reference More...
 
virtual bool vectorialCorrection () const
 if vectorial correction is provided More...
 
virtual ~JetCorrector ()
 

Static Public Member Functions

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 19 of file JetCorrector.h.

Member Typedef Documentation

◆ LorentzVector

Definition at line 21 of file JetCorrector.h.

Constructor & Destructor Documentation

◆ JetCorrector()

JetCorrector::JetCorrector ( )
inline

Definition at line 23 of file JetCorrector.h.

23 {};

◆ ~JetCorrector()

virtual JetCorrector::~JetCorrector ( )
inlinevirtual

Definition at line 24 of file JetCorrector.h.

24 {};

Member Function Documentation

◆ correction() [1/5]

virtual double JetCorrector::correction ( const LorentzVector fJet) const
pure virtual

◆ correction() [2/5]

virtual double JetCorrector::correction ( const reco::Jet fJet) const
pure virtual

apply correction using Jet information only

Implemented in L6SLBCorrector, L1FastjetCorrector, L1JPTOffsetCorrector, L1OffsetCorrector, LXXXCorrector, and ChainedJetCorrector.

◆ correction() [3/5]

double JetCorrector::correction ( const reco::Jet fJet,
const edm::Event fEvent,
const edm::EventSetup fSetup 
) const
virtual

apply correction using all event information

Reimplemented in L1FastjetCorrector, L1JPTOffsetCorrector, L1OffsetCorrector, and ChainedJetCorrector.

Definition at line 14 of file JetCorrector.cc.

References correction(), eventRequired(), and refRequired().

14  {
15  if (eventRequired() && !refRequired()) {
16  edm::LogError("Missing Jet Correction Method")
17  << "Undefined Jet Correction method requiring event data is called" << std::endl;
18  return 0;
19  }
20  return correction(fJet);
21 }
virtual double correction(const LorentzVector &fJet) const =0
get correction using Jet information only
Log< level::Error, false > LogError
virtual bool eventRequired() const =0
if correction needs event information
virtual bool refRequired() const =0
if correction needs the jet reference

◆ correction() [4/5]

double JetCorrector::correction ( const reco::Jet fJet,
const edm::RefToBase< reco::Jet > &  fJetRef,
const edm::Event fEvent,
const edm::EventSetup fSetup 
) const
virtual

apply correction using all event information

Reimplemented in L6SLBCorrector, and ChainedJetCorrector.

Definition at line 23 of file JetCorrector.cc.

References correction(), eventRequired(), and refRequired().

26  {
27  if (eventRequired() && refRequired()) {
28  edm::LogError("Missing Jet Correction Method")
29  << "Undefined Jet Correction method requiring event data and jet reference is called" << std::endl;
30  return 0;
31  }
32  return correction(fJet);
33 }
virtual double correction(const LorentzVector &fJet) const =0
get correction using Jet information only
Log< level::Error, false > LogError
virtual bool eventRequired() const =0
if correction needs event information
virtual bool refRequired() const =0
if correction needs the jet reference

◆ correction() [5/5]

double JetCorrector::correction ( const reco::Jet fJet,
const edm::RefToBase< reco::Jet > &  fJetRef,
const edm::Event fEvent,
const edm::EventSetup fSetup,
LorentzVector corrected 
) const
virtual

Apply vectorial correction using all event information.

Definition at line 35 of file JetCorrector.cc.

References correction(), and vectorialCorrection().

39  {
40  if (vectorialCorrection()) {
41  edm::LogError("Missing Jet Correction Method")
42  << "Undefined Jet (vectorial) correction method requiring event data is called" << std::endl;
43  return 0;
44  }
45  return correction(fJet);
46 }
virtual bool vectorialCorrection() const
if vectorial correction is provided
Definition: JetCorrector.h:62
virtual double correction(const LorentzVector &fJet) const =0
get correction using Jet information only
Log< level::Error, false > LogError

◆ eventRequired()

virtual bool JetCorrector::eventRequired ( ) const
pure virtual

if correction needs event information

Implemented in L6SLBCorrector, L1FastjetCorrector, L1JPTOffsetCorrector, L1OffsetCorrector, LXXXCorrector, and ChainedJetCorrector.

Referenced by correction().

◆ getJetCorrector()

const JetCorrector * JetCorrector::getJetCorrector ( const std::string &  fName,
const edm::EventSetup fSetup 
)
static

retrieve corrector from the event setup. troughs exception if something is missing

Definition at line 48 of file JetCorrector.cc.

References MainPageGenerator::fName, edm::EventSetup::get(), patZpeak::handle, and AlCaHarvesting_cff::record.

Referenced by L1JPTOffsetCorrector::correction().

◆ refRequired()

virtual bool JetCorrector::refRequired ( ) const
pure virtual

if correction needs the jet reference

Implemented in L6SLBCorrector, L1FastjetCorrector, L1OffsetCorrector, L1JPTOffsetCorrector, LXXXCorrector, and ChainedJetCorrector.

Referenced by correction().

◆ vectorialCorrection()

bool JetCorrector::vectorialCorrection ( ) const
inlinevirtual

if vectorial correction is provided

Definition at line 62 of file JetCorrector.h.

Referenced by correction().

62 { return false; }