CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
L1OffsetCorrector.cc
Go to the documentation of this file.
1 // Implementation of class L1OffsetCorrector.
2 // L1Offset jet corrector class.
3 
6 
15 
16 using namespace std;
17 
18 
19 //------------------------------------------------------------------------
20 //--- L1OffsetCorrector constructor ------------------------------------------
21 //------------------------------------------------------------------------
23 {
24  mVertexCollName = fConfig.getParameter<std::string>("vertexCollection");
25  mMinVtxNdof = fConfig.getParameter<int>("minVtxNdof");
26  if (fParam.definitions().level() != "L1Offset")
27  throw cms::Exception("L1OffsetCorrector")<<" correction level: "<<fParam.definitions().level()<<" is not L1Offset";
28  vector<JetCorrectorParameters> vParam;
29  vParam.push_back(fParam);
30  mCorrector = new FactorizedJetCorrectorCalculator(vParam);
31 }
32 //------------------------------------------------------------------------
33 //--- L1OffsetCorrector destructor -------------------------------------------
34 //------------------------------------------------------------------------
36 {
37  delete mCorrector;
38 }
39 //------------------------------------------------------------------------
40 //--- Returns correction for a given 4-vector ----------------------------
41 //------------------------------------------------------------------------
43 {
44  throw cms::Exception("EventRequired")
45  <<"Wrong interface correction(LorentzVector), event required!";
46  return 1.0;
47 }
48 //------------------------------------------------------------------------
49 //--- Returns correction for a given jet ---------------------------------
50 //------------------------------------------------------------------------
51 double L1OffsetCorrector::correction(const reco::Jet& fJet) const
52 {
53  throw cms::Exception("EventRequired")
54  <<"Wrong interface correction(reco::Jet), event required!";
55  return 1.0;
56 }
57 //------------------------------------------------------------------------
58 //--- Returns correction for a given jet using event indormation ---------
59 //------------------------------------------------------------------------
61  const edm::Event& fEvent,
62  const edm::EventSetup& fSetup) const
63 {
64  double result = 1.;
66  fEvent.getByLabel(mVertexCollName,recVtxs);
67  int NPV(0);
68  for(unsigned int ind=0;ind<recVtxs->size();ind++) {
69  if (!((*recVtxs)[ind].isFake()) && (*recVtxs)[ind].ndof() > mMinVtxNdof) {
70  NPV++;
71  }
72  }
73  if (NPV > 0) {
75  values.setJetEta(fJet.eta());
76  values.setJetPt(fJet.pt());
77  values.setJetE(fJet.energy());
78  values.setNPV(NPV);
79  result = mCorrector->getCorrection(values);
80  }
81  return result;
82 }
83 
T getParameter(std::string const &) const
Base class for all types of Jets.
Definition: Jet.h:20
const Definitions & definitions() const
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
virtual ~L1OffsetCorrector()
virtual double energy() const
energy
virtual double correction(const LorentzVector &fJet) const
get correction using Jet information only
tuple result
Definition: query.py:137
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
L1OffsetCorrector(const JetCorrectorParameters &fConfig, const edm::ParameterSet &fParameters)
reco::Particle::LorentzVector LorentzVector
Definition: JetCorrector.h:23