Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // Implementation of class LXXXCorrectorImpl.
0002 // Generic LX jet corrector class.
0003 
0004 #include "JetMETCorrections/Algorithms/interface/LXXXCorrectorImpl.h"
0005 #include "CondFormats/JetMETObjects/interface/FactorizedJetCorrectorCalculator.h"
0006 
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 #include "FWCore/Framework/interface/ConsumesCollector.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0012 #include "DataFormats/JetReco/interface/Jet.h"
0013 #include "DataFormats/JetReco/interface/CaloJet.h"
0014 #include "DataFormats/JetReco/interface/JPTJet.h"
0015 
0016 using namespace std;
0017 
0018 LXXXCorrectorImplMaker::LXXXCorrectorImplMaker(edm::ParameterSet const& fConfig, edm::ConsumesCollector iC)
0019     : JetCorrectorImplMakerBase(fConfig, std::move(iC)) {}
0020 
0021 std::unique_ptr<reco::JetCorrectorImpl> LXXXCorrectorImplMaker::make(edm::Event const&, edm::EventSetup const& fSetup) {
0022   unsigned int level = 0;
0023   auto calculator = getCalculator(fSetup, [&level](std::string const& levelName) {
0024     if (levelName == "L2Relative")
0025       level = 2;
0026     else if (levelName == "L3Absolute")
0027       level = 3;
0028     else if (levelName == "L4EMF")
0029       level = 4;
0030     else if (levelName == "L5Flavor")
0031       level = 5;
0032     else if (levelName == "L7Parton")
0033       level = 7;
0034     else if (levelName == "L2L3Residual")
0035       level = 8;
0036     else
0037       throw cms::Exception("LXXXCorrectorImpl") << " unknown correction level " << levelName;
0038   });
0039   return std::unique_ptr<reco::JetCorrectorImpl>(new LXXXCorrectorImpl(calculator, level));
0040 }
0041 
0042 void LXXXCorrectorImplMaker::fillDescriptions(edm::ConfigurationDescriptions& iDescriptions) {
0043   edm::ParameterSetDescription desc;
0044   addToDescription(desc);
0045 
0046   iDescriptions.addDefault(desc);
0047 }
0048 
0049 //------------------------------------------------------------------------
0050 //--- LXXXCorrectorImpl constructor ------------------------------------------
0051 //------------------------------------------------------------------------
0052 LXXXCorrectorImpl::LXXXCorrectorImpl(std::shared_ptr<FactorizedJetCorrectorCalculator const> calculator,
0053                                      unsigned int level)
0054     : mLevel(level), mCorrector(calculator) {}
0055 //------------------------------------------------------------------------
0056 //--- Returns correction for a given 4-vector ----------------------------
0057 //------------------------------------------------------------------------
0058 double LXXXCorrectorImpl::correction(const LorentzVector& fJet) const {
0059   // L4 correction requires more information than a simple 4-vector
0060   if (mLevel == 4) {
0061     throw cms::Exception("Invalid jet type") << "L4EMFCorrection is applicable to CaloJets only";
0062     return 1;
0063   }
0064 
0065   FactorizedJetCorrectorCalculator::VariableValues values;
0066   values.setJetEta(fJet.eta());
0067   values.setJetE(fJet.energy());
0068   values.setJetPt(fJet.pt());
0069   values.setJetPhi(fJet.phi());
0070 
0071   return mCorrector->getCorrection(values);
0072 }
0073 //------------------------------------------------------------------------
0074 //--- Returns correction for a given jet ---------------------------------
0075 //------------------------------------------------------------------------
0076 double LXXXCorrectorImpl::correction(const reco::Jet& fJet) const {
0077   double result = 1.;
0078   // L4 correction applies to Calojets only
0079   if (mLevel == 4) {
0080     const reco::CaloJet& caloJet = dynamic_cast<const reco::CaloJet&>(fJet);
0081     FactorizedJetCorrectorCalculator::VariableValues values;
0082     values.setJetEta(fJet.eta());
0083     values.setJetPt(fJet.pt());
0084     values.setJetEMF(caloJet.emEnergyFraction());
0085     result = mCorrector->getCorrection(values);
0086   } else
0087     result = correction(fJet.p4());
0088   return result;
0089 }