CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Private Attributes
L1JPTOffsetCorrectorImpl Class Reference

#include <L1JPTOffsetCorrectorImpl.h>

Inheritance diagram for L1JPTOffsetCorrectorImpl:
reco::JetCorrectorImpl

Public Types

typedef
L1JPTOffsetCorrectorImplMaker 
Maker
 
- Public Types inherited from reco::JetCorrectorImpl
typedef
reco::Particle::LorentzVector 
LorentzVector
 

Public Member Functions

virtual double correction (const LorentzVector &fJet) const override
 get correction using Jet information only More...
 
virtual double correction (const reco::Jet &fJet) const override
 apply correction using Jet information only More...
 
 L1JPTOffsetCorrectorImpl (std::shared_ptr< FactorizedJetCorrectorCalculator const > corrector, const reco::JetCorrector *offsetService)
 
virtual bool refRequired () const override
 if correction needs the jet reference More...
 
- Public Member Functions inherited from reco::JetCorrectorImpl
virtual double correction (const reco::Jet &fJet, const edm::RefToBase< reco::Jet > &fJetRef) const
 apply correction using Ref More...
 
virtual double correction (const reco::Jet &fJet, const edm::RefToBase< reco::Jet > &fJetRef, LorentzVector &corrected) const
 Apply vectorial correction. More...
 
 JetCorrectorImpl ()
 
virtual bool vectorialCorrection () const
 if vectorial correction is provided More...
 
virtual ~JetCorrectorImpl ()
 

Private Attributes

std::shared_ptr
< FactorizedJetCorrectorCalculator
const > 
corrector_
 
const reco::JetCorrectoroffsetService_
 

Detailed Description

Definition at line 37 of file L1JPTOffsetCorrectorImpl.h.

Member Typedef Documentation

Definition at line 40 of file L1JPTOffsetCorrectorImpl.h.

Constructor & Destructor Documentation

L1JPTOffsetCorrectorImpl::L1JPTOffsetCorrectorImpl ( std::shared_ptr< FactorizedJetCorrectorCalculator const >  corrector,
const reco::JetCorrector offsetService 
)

Definition at line 65 of file L1JPTOffsetCorrectorImpl.cc.

66  :
67  offsetService_(offsetService),
69 {
70 }
const reco::JetCorrector * offsetService_
std::shared_ptr< FactorizedJetCorrectorCalculator const > corrector_
tuple corrector
Definition: mvaPFMET_cff.py:86

Member Function Documentation

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

get correction using Jet information only

Implements reco::JetCorrectorImpl.

Definition at line 78 of file L1JPTOffsetCorrectorImpl.cc.

References Exception.

79 {
80  throw cms::Exception("EventRequired")
81  <<"Wrong interface correction(LorentzVector), event required!";
82  return 1.0;
83 }
double L1JPTOffsetCorrectorImpl::correction ( const reco::Jet fJet) const
overridevirtual

apply correction using Jet information only

Implements reco::JetCorrectorImpl.

Definition at line 87 of file L1JPTOffsetCorrectorImpl.cc.

References reco::JetCorrector::correction(), corrector_, reco::LeafCandidate::energy(), reco::JPTJet::getCaloJetRef(), hltrates_dqm_sourceclient-live_cfg::offset, offsetService_, reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::LeafCandidate::pz(), mps_fire::result, and makeHLTPrescaleTable::values.

88 {
89  double result = 1.;
90  const reco::JPTJet& jptjet = dynamic_cast <const reco::JPTJet&> (fJet);
91  edm::RefToBase<reco::Jet> jptjetRef = jptjet.getCaloJetRef();
92  reco::CaloJet const * rawcalojet = dynamic_cast<reco::CaloJet const *>( &* jptjetRef);
93  //------ access the offset correction service ----------------
94  double offset = 1.0;
95  if (offsetService_) {
96  offset = offsetService_->correction(*rawcalojet);
97  }
98  //------ calculate the correction for the JPT jet ------------
99  TLorentzVector JPTrawP4(rawcalojet->px(),rawcalojet->py(),rawcalojet->pz(),rawcalojet->energy());
101  values.setJPTrawP4(JPTrawP4);
102  values.setJPTrawOff(offset);
103  values.setJetE(fJet.energy());
104  result = corrector_->getCorrection(values);
105  return result;
106 }
Jets made from CaloTowers.
Definition: CaloJet.h:29
virtual double energy() const final
energy
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:47
const edm::RefToBase< reco::Jet > & getCaloJetRef() const
Definition: JPTJet.h:130
const reco::JetCorrector * offsetService_
tuple result
Definition: mps_fire.py:84
std::shared_ptr< FactorizedJetCorrectorCalculator const > corrector_
Jets made from CaloJets corrected for ZSP and tracks.
Definition: JPTJet.h:29
virtual double py() const final
y coordinate of momentum vector
virtual double pz() const final
z coordinate of momentum vector
virtual double px() const final
x coordinate of momentum vector
virtual bool L1JPTOffsetCorrectorImpl::refRequired ( ) const
inlineoverridevirtual

if correction needs the jet reference

Implements reco::JetCorrectorImpl.

Definition at line 53 of file L1JPTOffsetCorrectorImpl.h.

53 {return false;}

Member Data Documentation

std::shared_ptr<FactorizedJetCorrectorCalculator const> L1JPTOffsetCorrectorImpl::corrector_
private

Definition at line 58 of file L1JPTOffsetCorrectorImpl.h.

Referenced by correction().

const reco::JetCorrector* L1JPTOffsetCorrectorImpl::offsetService_
private

Definition at line 57 of file L1JPTOffsetCorrectorImpl.h.

Referenced by correction().