CMS 3D CMS Logo

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

#include <L1JPTOffsetCorrector.h>

Inheritance diagram for L1JPTOffsetCorrector:
JetCorrector

Public Member Functions

double correction (const LorentzVector &fJet) const override
 get 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...
 
 L1JPTOffsetCorrector (const JetCorrectorParameters &fConfig, const edm::ParameterSet &fParameters)
 
bool refRequired () const override
 if correction needs the jet reference More...
 
 ~L1JPTOffsetCorrector () 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

FactorizedJetCorrectorCalculatormCorrector
 
bool mIsOffsetSet
 
std::string mOffsetService
 

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 16 of file L1JPTOffsetCorrector.h.

Constructor & Destructor Documentation

◆ L1JPTOffsetCorrector()

L1JPTOffsetCorrector::L1JPTOffsetCorrector ( const JetCorrectorParameters fConfig,
const edm::ParameterSet fParameters 
)

Definition at line 23 of file L1JPTOffsetCorrector.cc.

References JetCorrectorParameters::definitions(), edm::ParameterSet::getParameter(), JetCorrectorParameters::Definitions::level(), and AlCaHLTBitMon_QueryRunRegistry::string.

23  {
24  mOffsetService = fConfig.getParameter<std::string>("offsetService");
25  mIsOffsetSet = false;
26  if (!mOffsetService.empty())
27  mIsOffsetSet = true;
28  if (fParam.definitions().level() != "L1JPTOffset")
29  throw cms::Exception("L1OffsetCorrector")
30  << " correction level: " << fParam.definitions().level() << " is not L1JPTOffset";
31  vector<JetCorrectorParameters> vParam;
32  vParam.push_back(fParam);
34 }
FactorizedJetCorrectorCalculator * mCorrector

◆ ~L1JPTOffsetCorrector()

L1JPTOffsetCorrector::~L1JPTOffsetCorrector ( )
override

Definition at line 38 of file L1JPTOffsetCorrector.cc.

38 { delete mCorrector; }
FactorizedJetCorrectorCalculator * mCorrector

Member Function Documentation

◆ correction() [1/3]

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

get correction using Jet information only

Implements JetCorrector.

Definition at line 42 of file L1JPTOffsetCorrector.cc.

References Exception.

42  {
43  throw cms::Exception("EventRequired") << "Wrong interface correction(LorentzVector), event required!";
44  return 1.0;
45 }

◆ correction() [2/3]

double L1JPTOffsetCorrector::correction ( const reco::Jet fJet) const
overridevirtual

apply correction using Jet information only

Implements JetCorrector.

Definition at line 49 of file L1JPTOffsetCorrector.cc.

References Exception.

49  {
50  throw cms::Exception("EventRequired") << "Wrong interface correction(reco::Jet), event required!";
51  return 1.0;
52 }

◆ correction() [3/3]

double L1JPTOffsetCorrector::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 56 of file L1JPTOffsetCorrector.cc.

References JetCorrector::correction(), reco::LeafCandidate::energy(), hcaldqm::fEvent, reco::JPTJet::getCaloJetRef(), JetCorrector::getJetCorrector(), hltrates_dqm_sourceclient-live_cfg::offset, reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::LeafCandidate::pz(), mps_fire::result, and contentValuesCheck::values.

58  {
59  double result = 1.;
60  const reco::JPTJet& jptjet = dynamic_cast<const reco::JPTJet&>(fJet);
61  const edm::RefToBase<reco::Jet>& rawcalojet = jptjet.getCaloJetRef();
62  //------ access the offset correction service ----------------
63  double offset = 1.0;
64  if (mIsOffsetSet) {
65  const JetCorrector* OffsetCorrector = JetCorrector::getJetCorrector(mOffsetService, fSetup);
66  offset = OffsetCorrector->correction(*rawcalojet, fEvent, fSetup);
67  }
68  //------ calculate the correction for the JPT jet ------------
69  TLorentzVector JPTrawP4(rawcalojet->px(), rawcalojet->py(), rawcalojet->pz(), rawcalojet->energy());
71  values.setJPTrawP4(JPTrawP4);
72  values.setJPTrawOff(offset);
73  values.setJetE(fJet.energy());
75  return result;
76 }
double pz() const final
z coordinate of momentum vector
virtual double correction(const LorentzVector &fJet) const =0
get correction using Jet information only
double px() const final
x coordinate of momentum vector
Jets made from CaloJets corrected for ZSP and tracks.
Definition: JPTJet.h:28
FactorizedJetCorrectorCalculator * mCorrector
const edm::RefToBase< reco::Jet > & getCaloJetRef() const
Definition: JPTJet.h:118
double py() const final
y coordinate of momentum vector
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: JetCorrector.cc:48
double energy() const final
energy

◆ eventRequired()

bool L1JPTOffsetCorrector::eventRequired ( ) const
inlineoverridevirtual

if correction needs event information

Implements JetCorrector.

Definition at line 33 of file L1JPTOffsetCorrector.h.

33 { return true; }

◆ refRequired()

bool L1JPTOffsetCorrector::refRequired ( ) const
inlineoverridevirtual

if correction needs the jet reference

Implements JetCorrector.

Definition at line 34 of file L1JPTOffsetCorrector.h.

34 { return false; }

Member Data Documentation

◆ mCorrector

FactorizedJetCorrectorCalculator* L1JPTOffsetCorrector::mCorrector
private

Definition at line 40 of file L1JPTOffsetCorrector.h.

◆ mIsOffsetSet

bool L1JPTOffsetCorrector::mIsOffsetSet
private

Definition at line 39 of file L1JPTOffsetCorrector.h.

◆ mOffsetService

std::string L1JPTOffsetCorrector::mOffsetService
private

Definition at line 38 of file L1JPTOffsetCorrector.h.