Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 ////////////////////////////////////////////////////////////////////////////////
0002 //
0003 // L6SLBCorrectorImpl
0004 // --------------
0005 //
0006 //           25/10/2009   Hauke Held             <hauke.held@cern.ch>
0007 //                        Philipp Schieferdecker <philipp.schieferdecker@cern.ch
0008 ////////////////////////////////////////////////////////////////////////////////
0009 
0010 #include "JetMETCorrections/Algorithms/interface/L6SLBCorrectorImpl.h"
0011 
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0015 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0016 #include "FWCore/ParameterSet/interface/FileInPath.h"
0017 #include "FWCore/Framework/interface/Event.h"
0018 #include "FWCore/Framework/interface/ConsumesCollector.h"
0019 
0020 #include "DataFormats/TrackReco/interface/Track.h"
0021 
0022 #include <string>
0023 
0024 using namespace std;
0025 
0026 L6SLBCorrectorImplMaker::L6SLBCorrectorImplMaker(edm::ParameterSet const& fConfig, edm::ConsumesCollector fCollector)
0027     : JetCorrectorImplMakerBase(fConfig, fCollector),
0028       elecToken_(fCollector.consumes<std::vector<reco::SoftLeptonTagInfo>>(
0029           fConfig.getParameter<edm::InputTag>("srcBTagInfoElectron"))),
0030       muonToken_(fCollector.consumes<std::vector<reco::SoftLeptonTagInfo>>(
0031           fConfig.getParameter<edm::InputTag>("srcBTagInfoMuon"))),
0032       addMuonToJet_(fConfig.getParameter<bool>("addMuonToJet")) {}
0033 
0034 std::unique_ptr<reco::JetCorrectorImpl> L6SLBCorrectorImplMaker::make(edm::Event const& fEvent,
0035                                                                       edm::EventSetup const& fSetup) {
0036   edm::Handle<std::vector<reco::SoftLeptonTagInfo>> muoninfos;
0037   fEvent.getByToken(muonToken_, muoninfos);
0038   edm::RefProd<std::vector<reco::SoftLeptonTagInfo>> muonProd(muoninfos);
0039 
0040   edm::Handle<std::vector<reco::SoftLeptonTagInfo>> elecinfos;
0041   fEvent.getByToken(elecToken_, elecinfos);
0042   edm::RefProd<std::vector<reco::SoftLeptonTagInfo>> elecProd(elecinfos);
0043 
0044   auto calculator = getCalculator(fSetup, [](std::string const&) {});
0045   return std::unique_ptr<reco::JetCorrectorImpl>(new L6SLBCorrectorImpl(calculator, muonProd, elecProd, addMuonToJet_));
0046 }
0047 
0048 void L6SLBCorrectorImplMaker::fillDescriptions(edm::ConfigurationDescriptions& iDescriptions) {
0049   edm::ParameterSetDescription desc;
0050   addToDescription(desc);
0051   desc.add<edm::InputTag>("srcBTagInfoElectron");
0052   desc.add<edm::InputTag>("srcBTagInfoMuon");
0053   desc.add<bool>("addMuonToJet");
0054   iDescriptions.addDefault(desc);
0055 }
0056 
0057 ////////////////////////////////////////////////////////////////////////////////
0058 // construction / destruction
0059 ////////////////////////////////////////////////////////////////////////////////
0060 
0061 //______________________________________________________________________________
0062 L6SLBCorrectorImpl::L6SLBCorrectorImpl(std::shared_ptr<FactorizedJetCorrectorCalculator const> corrector,
0063                                        edm::RefProd<std::vector<reco::SoftLeptonTagInfo>> const& bTagInfoMuon,
0064                                        edm::RefProd<std::vector<reco::SoftLeptonTagInfo>> const& bTagInfoElec,
0065                                        bool addMuonToJet)
0066     : corrector_(corrector), bTagInfoMuon_(bTagInfoMuon), bTagInfoElec_(bTagInfoElec), addMuonToJet_(addMuonToJet) {}
0067 
0068 ////////////////////////////////////////////////////////////////////////////////
0069 // implementation of member functions
0070 ////////////////////////////////////////////////////////////////////////////////
0071 
0072 //______________________________________________________________________________
0073 double L6SLBCorrectorImpl::correction(const LorentzVector& fJet) const {
0074   throw cms::Exception("EventRequired") << "Wrong interface correction(LorentzVector), event required!";
0075   return 1.0;
0076 }
0077 
0078 //______________________________________________________________________________
0079 double L6SLBCorrectorImpl::correction(const reco::Jet& fJet) const {
0080   throw cms::Exception("EventRequired") << "Wrong interface correction(reco::Jet), event required!";
0081   return 1.0;
0082 }
0083 
0084 //______________________________________________________________________________
0085 double L6SLBCorrectorImpl::correction(const reco::Jet& fJet, const edm::RefToBase<reco::Jet>& refToRawJet) const {
0086   FactorizedJetCorrectorCalculator::VariableValues values;
0087   values.setJetPt(fJet.pt());
0088   values.setJetEta(fJet.eta());
0089   values.setJetPhi(fJet.phi());
0090   values.setJetE(fJet.energy());
0091 
0092   const reco::SoftLeptonTagInfo& sltMuon = (*bTagInfoMuon_)[getBTagInfoIndex(refToRawJet, *bTagInfoMuon_)];
0093   if (sltMuon.leptons() > 0) {
0094     const edm::RefToBase<reco::Track>& trackRef = sltMuon.lepton(0);
0095     values.setLepPx(trackRef->px());
0096     values.setLepPy(trackRef->py());
0097     values.setLepPz(trackRef->pz());
0098     values.setAddLepToJet(addMuonToJet_);
0099     return corrector_->getCorrection(values);
0100   } else {
0101     const reco::SoftLeptonTagInfo& sltElec = (*bTagInfoElec_)[getBTagInfoIndex(refToRawJet, *bTagInfoElec_)];
0102     if (sltElec.leptons() > 0) {
0103       const edm::RefToBase<reco::Track>& trackRef = sltElec.lepton(0);
0104       values.setLepPx(trackRef->px());
0105       values.setLepPy(trackRef->py());
0106       values.setLepPz(trackRef->pz());
0107       values.setAddLepToJet(false);
0108       return corrector_->getCorrection(values);
0109     }
0110   }
0111   return 1.0;
0112 }
0113 
0114 ////////////////////////////////////////////////////////////////////////////////
0115 // implementation of private member functions
0116 ////////////////////////////////////////////////////////////////////////////////
0117 
0118 //______________________________________________________________________________
0119 int L6SLBCorrectorImpl::getBTagInfoIndex(const edm::RefToBase<reco::Jet>& refToRawJet,
0120                                          const vector<reco::SoftLeptonTagInfo>& tags) const {
0121   for (unsigned int i = 0; i < tags.size(); i++)
0122     if (tags[i].jet().get() == refToRawJet.get())
0123       return i;
0124   return -1;
0125 }