Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:19:15

0001 // Implementation of class L1OffsetCorrectorImpl.
0002 // L1Offset jet corrector class.
0003 
0004 #include "JetMETCorrections/Algorithms/interface/L1OffsetCorrectorImpl.h"
0005 #include "CondFormats/JetMETObjects/interface/FactorizedJetCorrector.h"
0006 
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0010 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/EventSetup.h"
0013 #include "FWCore/Framework/interface/ConsumesCollector.h"
0014 #include "DataFormats/JetReco/interface/Jet.h"
0015 #include "DataFormats/VertexReco/interface/Vertex.h"
0016 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0017 #include "CondFormats/JetMETObjects/interface/JetCorrectorParameters.h"
0018 
0019 using namespace std;
0020 
0021 L1OffsetCorrectorImplMaker::L1OffsetCorrectorImplMaker(edm::ParameterSet const& fConfig,
0022                                                        edm::ConsumesCollector fCollector)
0023     : JetCorrectorImplMakerBase(fConfig, fCollector),
0024       verticesToken_(
0025           fCollector.consumes<reco::VertexCollection>(fConfig.getParameter<edm::InputTag>("vertexCollection"))),
0026       minVtxNdof_(fConfig.getParameter<int>("minVtxNdof")) {}
0027 std::unique_ptr<reco::JetCorrectorImpl> L1OffsetCorrectorImplMaker::make(edm::Event const& fEvent,
0028                                                                          edm::EventSetup const& fSetup) {
0029   edm::Handle<reco::VertexCollection> recVtxs;
0030   fEvent.getByToken(verticesToken_, recVtxs);
0031   int NPV(0);
0032   for (auto const& vertex : *recVtxs) {
0033     if ((not vertex.isFake()) and (vertex.ndof() > minVtxNdof_)) {
0034       NPV++;
0035     }
0036   }
0037 
0038   auto calculator = getCalculator(fSetup, [](std::string const& level) {
0039     if (level != "L1Offset") {
0040       throw cms::Exception("L1OffsetCorrectorImpl") << " correction level: " << level << " is not L1Offset";
0041     }
0042   });
0043   return std::unique_ptr<reco::JetCorrectorImpl>(new L1OffsetCorrectorImpl(calculator, NPV));
0044 }
0045 
0046 void L1OffsetCorrectorImplMaker::fillDescriptions(edm::ConfigurationDescriptions& iDescriptions) {
0047   edm::ParameterSetDescription desc;
0048   addToDescription(desc);
0049   desc.add<edm::InputTag>("vertexCollection");
0050   desc.add<int>("minVtxNdof");
0051   iDescriptions.addDefault(desc);
0052 }
0053 
0054 //------------------------------------------------------------------------
0055 //--- L1OffsetCorrectorImpl constructor ------------------------------------------
0056 //------------------------------------------------------------------------
0057 L1OffsetCorrectorImpl::L1OffsetCorrectorImpl(std::shared_ptr<FactorizedJetCorrectorCalculator const> calculator,
0058                                              int npv)
0059     : corrector_(calculator), npv_(npv) {}
0060 
0061 //------------------------------------------------------------------------
0062 //--- Returns correction for a given 4-vector ----------------------------
0063 //------------------------------------------------------------------------
0064 double L1OffsetCorrectorImpl::correction(const LorentzVector& fJet) const {
0065   throw cms::Exception("EventRequired") << "Wrong interface correction(LorentzVector), event required!";
0066   return 1.0;
0067 }
0068 //------------------------------------------------------------------------
0069 //--- Returns correction for a given jet ---------------------------------
0070 //------------------------------------------------------------------------
0071 double L1OffsetCorrectorImpl::correction(const reco::Jet& fJet) const {
0072   double result = 1.;
0073   if (npv_ > 0) {
0074     FactorizedJetCorrectorCalculator::VariableValues values;
0075     values.setJetEta(fJet.eta());
0076     values.setJetPt(fJet.pt());
0077     values.setJetE(fJet.energy());
0078     values.setNPV(npv_);
0079     result = corrector_->getCorrection(values);
0080   }
0081   return result;
0082 }