CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
L1OffsetCorrectorImpl.cc
Go to the documentation of this file.
1 // Implementation of class L1OffsetCorrectorImpl.
2 // L1Offset jet corrector class.
3 
6 
18 
19 using namespace std;
20 
22  edm::ConsumesCollector fCollector) :
24  verticesToken_(fCollector.consumes<reco::VertexCollection>(fConfig.getParameter<edm::InputTag>("vertexCollection"))),
25  minVtxNdof_(fConfig.getParameter<int>("minVtxNdof"))
26 {
27 }
28 std::unique_ptr<reco::JetCorrectorImpl>
30 {
32  fEvent.getByToken(verticesToken_,recVtxs);
33  int NPV(0);
34  for(auto const& vertex : *recVtxs) {
35  if ((not vertex.isFake()) and (vertex.ndof() > minVtxNdof_)) {
36  NPV++;
37  }
38  }
39 
40  auto calculator = getCalculator(fSetup,
41  [](std::string const& level)
42  {
43  if ( level != "L1Offset") {
44  throw cms::Exception("L1OffsetCorrectorImpl")<<" correction level: "<<level<<" is not L1Offset";
45  }
46  });
47  return std::unique_ptr<reco::JetCorrectorImpl>(new L1OffsetCorrectorImpl(calculator,NPV));
48 }
49 
50 void
52 {
54  addToDescription(desc);
55  desc.add<edm::InputTag>("vertexCollection");
56  desc.add<int>("minVtxNdof");
57  iDescriptions.addDefault(desc);
58 }
59 
60 //------------------------------------------------------------------------
61 //--- L1OffsetCorrectorImpl constructor ------------------------------------------
62 //------------------------------------------------------------------------
63 L1OffsetCorrectorImpl::L1OffsetCorrectorImpl(std::shared_ptr<FactorizedJetCorrectorCalculator const> calculator,
64  int npv):
65  corrector_(calculator),
66  npv_(npv)
67 {
68 }
69 
70 //------------------------------------------------------------------------
71 //--- Returns correction for a given 4-vector ----------------------------
72 //------------------------------------------------------------------------
74 {
75  throw cms::Exception("EventRequired")
76  <<"Wrong interface correction(LorentzVector), event required!";
77  return 1.0;
78 }
79 //------------------------------------------------------------------------
80 //--- Returns correction for a given jet ---------------------------------
81 //------------------------------------------------------------------------
83 {
84  double result = 1.;
85  if (npv_ > 0) {
87  values.setJetEta(fJet.eta());
88  values.setJetPt(fJet.pt());
89  values.setJetE(fJet.energy());
90  values.setNPV(npv_);
91  result = corrector_->getCorrection(values);
92  }
93  return result;
94 }
95 
static void fillDescriptions(edm::ConfigurationDescriptions &iDescriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
Base class for all types of Jets.
Definition: Jet.h:20
std::unique_ptr< reco::JetCorrectorImpl > make(edm::Event const &, edm::EventSetup const &)
reco::Particle::LorentzVector LorentzVector
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
virtual double energy() const
energy
void addDefault(ParameterSetDescription const &psetDescription)
tuple result
Definition: query.py:137
edm::EDGetTokenT< reco::VertexCollection > verticesToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
virtual double correction(const LorentzVector &fJet) const override
get correction using Jet information only
std::shared_ptr< FactorizedJetCorrectorCalculator const > corrector_
static void addToDescription(edm::ParameterSetDescription &iDescription)
tuple level
Definition: testEve_cfg.py:34
L1OffsetCorrectorImpl(std::shared_ptr< FactorizedJetCorrectorCalculator const > calculator, int npv)
std::shared_ptr< FactorizedJetCorrectorCalculator const > getCalculator(edm::EventSetup const &, std::function< void(std::string const &)> levelCheck)
L1OffsetCorrectorImplMaker(edm::ParameterSet const &, edm::ConsumesCollector)