CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
L1JPTOffsetCorrector.cc
Go to the documentation of this file.
1 // Implementation of class L1JPTOffsetCorrector.
2 // L1JPTOffset jet corrector class.
3 
6 
17 
18 using namespace std;
19 
20 
21 //------------------------------------------------------------------------
22 //--- L1OffsetCorrector constructor ------------------------------------------
23 //------------------------------------------------------------------------
25 {
26  mOffsetService = fConfig.getParameter<std::string>("offsetService");
27  mIsOffsetSet = false;
28  if (mOffsetService != "")
29  mIsOffsetSet = true;
30  if (fParam.definitions().level() != "L1JPTOffset")
31  throw cms::Exception("L1OffsetCorrector")<<" correction level: "<<fParam.definitions().level()<<" is not L1JPTOffset";
32  vector<JetCorrectorParameters> vParam;
33  vParam.push_back(fParam);
34  mCorrector = new FactorizedJetCorrector(vParam);
35 }
36 //------------------------------------------------------------------------
37 //--- L1OffsetCorrector destructor -------------------------------------------
38 //------------------------------------------------------------------------
40 {
41  delete mCorrector;
42 }
43 //------------------------------------------------------------------------
44 //--- Returns correction for a given 4-vector ----------------------------
45 //------------------------------------------------------------------------
47 {
48  throw cms::Exception("EventRequired")
49  <<"Wrong interface correction(LorentzVector), event required!";
50  return 1.0;
51 }
52 //------------------------------------------------------------------------
53 //--- Returns correction for a given jet ---------------------------------
54 //------------------------------------------------------------------------
56 {
57  throw cms::Exception("EventRequired")
58  <<"Wrong interface correction(reco::Jet), event required!";
59  return 1.0;
60 }
61 //------------------------------------------------------------------------
62 //--- Returns correction for a given jet using event indormation ---------
63 //------------------------------------------------------------------------
65  const edm::Event& fEvent,
66  const edm::EventSetup& fSetup) const
67 {
68  double result = 1.;
69  const reco::JPTJet& jptjet = dynamic_cast <const reco::JPTJet&> (fJet);
70  edm::RefToBase<reco::Jet> jptjetRef = jptjet.getCaloJetRef();
71  reco::CaloJet const * rawcalojet = dynamic_cast<reco::CaloJet const *>( &* jptjetRef);
72  //------ access the offset correction service ----------------
73  double offset = 1.0;
74  if (mIsOffsetSet) {
75  const JetCorrector* OffsetCorrector = JetCorrector::getJetCorrector(mOffsetService,fSetup);
76  offset = OffsetCorrector->correction(*rawcalojet,fEvent,fSetup);
77  }
78  //------ calculate the correction for the JPT jet ------------
79  TLorentzVector JPTrawP4(rawcalojet->px(),rawcalojet->py(),rawcalojet->pz(),rawcalojet->energy());
80  mCorrector->setJPTrawP4(JPTrawP4);
81  mCorrector->setJPTrawOff(offset);
82  mCorrector->setJetE(fJet.energy());
83  result = mCorrector->getCorrection();
84  return result;
85 }
86 
T getParameter(std::string const &) const
Jets made from CaloTowers.
Definition: CaloJet.h:30
virtual double correction(const LorentzVector &fJet) const =0
get correction using Jet information only
Base class for all types of Jets.
Definition: Jet.h:21
const Definitions & definitions() const
const edm::RefToBase< reco::Jet > & getCaloJetRef() const
Definition: JPTJet.h:130
virtual double correction(const LorentzVector &fJet) const
get correction using Jet information only
L1JPTOffsetCorrector(const JetCorrectorParameters &fConfig, const edm::ParameterSet &fParameters)
virtual double energy() const
energy
Jets made from CaloJets corrected for ZSP and tracks.
Definition: JPTJet.h:29
tuple result
Definition: query.py:137
unsigned int offset(bool)
virtual double px() const
x 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:51
virtual double pz() const
z coordinate of momentum vector
virtual double py() const
y coordinate of momentum vector
reco::Particle::LorentzVector LorentzVector
Definition: JetCorrector.h:24