#include <JetCorrector.h>
Public Types | |
typedef reco::Particle::LorentzVector | LorentzVector |
Public Member Functions | |
virtual double | correction (const LorentzVector &fJet) const =0 |
get correction using Jet information only | |
virtual double | correction (const reco::Jet &fJet) const =0 |
apply correction using Jet information only | |
virtual double | correction (const reco::Jet &fJet, const edm::Event &fEvent, const edm::EventSetup &fSetup) const |
apply correction using all event information | |
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 | |
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. | |
virtual bool | eventRequired () const =0 |
if correction needs event information | |
JetCorrector () | |
virtual bool | refRequired () const =0 |
if correction needs the jet reference | |
virtual bool | vectorialCorrection () const |
if vectorial correction is provided | |
virtual | ~JetCorrector () |
Static Public Member Functions | |
static const JetCorrector * | getJetCorrector (const std::string &fName, const edm::EventSetup &fSetup) |
retrieve corrector from the event setup. troughs exception if something is missing |
Definition at line 20 of file JetCorrector.h.
Definition at line 24 of file JetCorrector.h.
JetCorrector::JetCorrector | ( | ) | [inline] |
Definition at line 26 of file JetCorrector.h.
{};
virtual JetCorrector::~JetCorrector | ( | ) | [inline, virtual] |
Definition at line 27 of file JetCorrector.h.
{};
virtual double JetCorrector::correction | ( | const LorentzVector & | fJet | ) | const [pure virtual] |
get correction using Jet information only
Implemented in L1FastjetCorrector, L1JPTOffsetCorrector, L1OffsetCorrector, L6SLBCorrector, LXXXCorrector, JetPartonCorrector, ChainedJetCorrector, TauJetCorrector, and TCTauCorrector.
Referenced by TauJetCorrectorExample::analyze(), calcTopMass::analyze(), JetCorrectorOnTheFly< Jet >::analyze(), JetCorrectorDemo::analyze(), JPTJetTester::analyze(), JetMETHLTOfflineSource::analyze(), JetCorExample< Jet >::analyze(), FlavorJetCorrectionExample::analyze(), PartonJetCorrectionExample::analyze(), CMSDAS11DijetAnalyzer::analyze(), PFJetTester::analyze(), CaloJetTester::analyze(), L1JPTOffsetCorrector::correction(), correction(), TopSingleLepton::MonitorEnsemble::fill(), TopHLTSingleLepton::MonitorEnsemble::fill(), CorrectJet::operator()(), EcalTDigitizer< ESDigitizerTraits >::run(), TauMETAlgo::run(), SelectionStep< Object >::select(), and SelectionStepHLT< Object >::select().
virtual double JetCorrector::correction | ( | const reco::Jet & | fJet | ) | const [pure virtual] |
apply correction using Jet information only
Implemented in L1FastjetCorrector, L1JPTOffsetCorrector, L1OffsetCorrector, L6SLBCorrector, LXXXCorrector, ChainedJetCorrector, and TauJetCorrector.
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 26 of file JetCorrector.cc.
References correction(), eventRequired(), and refRequired().
{ if (eventRequired () && refRequired()) { edm::LogError ("Missing Jet Correction Method") << "Undefined Jet Correction method requiring event data and jet reference is called" << std::endl; return 0; } return correction (fJet); }
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 38 of file JetCorrector.cc.
References correction(), and vectorialCorrection().
{ if ( vectorialCorrection() ) { edm::LogError ("Missing Jet Correction Method") << "Undefined Jet (vectorial) correction method requiring event data is called" << std::endl; return 0; } return correction (fJet); }
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 15 of file JetCorrector.cc.
References correction(), eventRequired(), and refRequired().
{ if (eventRequired () && !refRequired()) { edm::LogError ("Missing Jet Correction Method") << "Undefined Jet Correction method requiring event data is called" << std::endl; return 0; } return correction (fJet); }
virtual bool JetCorrector::eventRequired | ( | ) | const [pure virtual] |
if correction needs event information
Implemented in L1FastjetCorrector, L1JPTOffsetCorrector, L1OffsetCorrector, L6SLBCorrector, LXXXCorrector, JetPartonCorrector, ChainedJetCorrector, TauJetCorrector, and TCTauCorrector.
Referenced by correction().
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 51 of file JetCorrector.cc.
References edm::EventSetup::get(), edm::eventsetup::EventSetupRecord::get(), patZpeak::handle, and record.
Referenced by CMSDAS11DijetTestAnalyzer::analyze(), TauJetCorrectorExample::analyze(), calcTopMass::analyze(), JetCorrectorDemo::analyze(), JetCorrectorOnTheFly< Jet >::analyze(), JetMETHLTOfflineSource::analyze(), JPTJetTester::analyze(), JetCorExample< Jet >::analyze(), FlavorJetCorrectionExample::analyze(), PartonJetCorrectionExample::analyze(), CMSDAS11DijetAnalyzer::analyze(), PFJetTester::analyze(), CaloJetTester::analyze(), L1JPTOffsetCorrector::correction(), TopSingleLepton::MonitorEnsemble::fill(), TopHLTSingleLepton::MonitorEnsemble::fill(), TopDiLeptonOffline::MonitorEnsemble::fill(), Type1PFMET::produce(), cms::TauMET::produce(), cms::Type1MET::produce(), SelectionStep< Object >::select(), SelectionStepHLT< Object >::select(), and CorrectJet::setEventSetup().
{ const JetCorrectionsRecord& record = fSetup.get <JetCorrectionsRecord> (); edm::ESHandle <JetCorrector> handle; record.get (fName, handle); return &*handle; }
virtual bool JetCorrector::refRequired | ( | ) | const [pure virtual] |
if correction needs the jet reference
Implemented in L1FastjetCorrector, L1JPTOffsetCorrector, L1OffsetCorrector, L6SLBCorrector, LXXXCorrector, and ChainedJetCorrector.
Referenced by correction().
bool JetCorrector::vectorialCorrection | ( | ) | const [inline, virtual] |
if vectorial correction is provided
Definition at line 67 of file JetCorrector.h.
Referenced by CMSDAS11DijetTestAnalyzer::analyze(), and correction().
{ return false; }