Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 14:21:24

0001 // Implementation of class L1OffsetCorrector.
0002 // L1Offset jet corrector class.
0003 
0004 #include "JetMETCorrections/Algorithms/interface/L1OffsetCorrector.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/FileInPath.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "DataFormats/JetReco/interface/Jet.h"
0013 #include "DataFormats/VertexReco/interface/Vertex.h"
0014 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0015 
0016 using namespace std;
0017 
0018 //------------------------------------------------------------------------
0019 //--- L1OffsetCorrector constructor ------------------------------------------
0020 //------------------------------------------------------------------------
0021 L1OffsetCorrector::L1OffsetCorrector(const JetCorrectorParameters& fParam, const edm::ParameterSet& fConfig) {
0022   mVertexCollName = fConfig.getParameter<std::string>("vertexCollection");
0023   mMinVtxNdof = fConfig.getParameter<int>("minVtxNdof");
0024   if (fParam.definitions().level() != "L1Offset")
0025     throw cms::Exception("L1OffsetCorrector")
0026         << " correction level: " << fParam.definitions().level() << " is not L1Offset";
0027   vector<JetCorrectorParameters> vParam;
0028   vParam.push_back(fParam);
0029   mCorrector = new FactorizedJetCorrectorCalculator(vParam);
0030 }
0031 //------------------------------------------------------------------------
0032 //--- L1OffsetCorrector destructor -------------------------------------------
0033 //------------------------------------------------------------------------
0034 L1OffsetCorrector::~L1OffsetCorrector() { delete mCorrector; }
0035 //------------------------------------------------------------------------
0036 //--- Returns correction for a given 4-vector ----------------------------
0037 //------------------------------------------------------------------------
0038 double L1OffsetCorrector::correction(const LorentzVector& fJet) const {
0039   throw cms::Exception("EventRequired") << "Wrong interface correction(LorentzVector), event required!";
0040   return 1.0;
0041 }
0042 //------------------------------------------------------------------------
0043 //--- Returns correction for a given jet ---------------------------------
0044 //------------------------------------------------------------------------
0045 double L1OffsetCorrector::correction(const reco::Jet& fJet) const {
0046   throw cms::Exception("EventRequired") << "Wrong interface correction(reco::Jet), event required!";
0047   return 1.0;
0048 }
0049 //------------------------------------------------------------------------
0050 //--- Returns correction for a given jet using event indormation ---------
0051 //------------------------------------------------------------------------
0052 double L1OffsetCorrector::correction(const reco::Jet& fJet,
0053                                      const edm::Event& fEvent,
0054                                      const edm::EventSetup& fSetup) const {
0055   double result = 1.;
0056   edm::Handle<reco::VertexCollection> recVtxs;
0057   fEvent.getByLabel(mVertexCollName, recVtxs);
0058   int NPV(0);
0059   for (unsigned int ind = 0; ind < recVtxs->size(); ind++) {
0060     if (!((*recVtxs)[ind].isFake()) && (*recVtxs)[ind].ndof() > mMinVtxNdof) {
0061       NPV++;
0062     }
0063   }
0064   if (NPV > 0) {
0065     FactorizedJetCorrectorCalculator::VariableValues values;
0066     values.setJetEta(fJet.eta());
0067     values.setJetPt(fJet.pt());
0068     values.setJetE(fJet.energy());
0069     values.setNPV(NPV);
0070     result = mCorrector->getCorrection(values);
0071   }
0072   return result;
0073 }