CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
L6SLBCorrectorImpl.cc
Go to the documentation of this file.
1 //
3 // L6SLBCorrectorImpl
4 // --------------
5 //
6 // 25/10/2009 Hauke Held <hauke.held@cern.ch>
7 // Philipp Schieferdecker <philipp.schieferdecker@cern.ch
9 
11 
19 
21 
22 
23 #include <string>
24 
25 
26 using namespace std;
27 
29  edm::ConsumesCollector fCollector)
30  : JetCorrectorImplMakerBase(fConfig)
31  , elecToken_(fCollector.consumes<std::vector<reco::SoftLeptonTagInfo>>(fConfig.getParameter<edm::InputTag>("srcBTagInfoElectron")))
32  , muonToken_(fCollector.consumes<std::vector<reco::SoftLeptonTagInfo>>(fConfig.getParameter<edm::InputTag>("srcBTagInfoMuon")))
33  , addMuonToJet_(fConfig.getParameter<bool>("addMuonToJet"))
34 {
35 }
36 
37 std::unique_ptr<reco::JetCorrectorImpl>
39 {
41  fEvent.getByToken(muonToken_,muoninfos);
43 
45  fEvent.getByToken(elecToken_,elecinfos);
47 
48  auto calculator = getCalculator(fSetup,
49  [](std::string const& ) {});
50  return std::unique_ptr<reco::JetCorrectorImpl>(new L6SLBCorrectorImpl(calculator,
51  muonProd,
52  elecProd,
53  addMuonToJet_));
54 }
55 
56 void
58 {
60  addToDescription(desc);
61  desc.add<edm::InputTag>("srcBTagInfoElectron");
62  desc.add<edm::InputTag>("srcBTagInfoMuon");
63  desc.add<bool>("addMuonToJet");
64  iDescriptions.addDefault(desc);
65 }
66 
68 // construction / destruction
70 
71 //______________________________________________________________________________
72 L6SLBCorrectorImpl::L6SLBCorrectorImpl(std::shared_ptr<FactorizedJetCorrectorCalculator const> corrector,
73  edm::RefProd<std::vector<reco::SoftLeptonTagInfo>> const& bTagInfoMuon,
74  edm::RefProd<std::vector<reco::SoftLeptonTagInfo>> const& bTagInfoElec,
75  bool addMuonToJet):
76  corrector_(corrector),
77  bTagInfoMuon_(bTagInfoMuon),
78  bTagInfoElec_(bTagInfoElec),
79  addMuonToJet_(addMuonToJet)
80 {
81 }
82 
84 // implementation of member functions
86 
87 //______________________________________________________________________________
89 {
90  throw cms::Exception("EventRequired")
91  <<"Wrong interface correction(LorentzVector), event required!";
92  return 1.0;
93 }
94 
95 
96 //______________________________________________________________________________
97 double L6SLBCorrectorImpl::correction(const reco::Jet& fJet) const
98 {
99  throw cms::Exception("EventRequired")
100  <<"Wrong interface correction(reco::Jet), event required!";
101  return 1.0;
102 }
103 
104 
105 //______________________________________________________________________________
107  const edm::RefToBase<reco::Jet>& refToRawJet) const
108 {
110  values.setJetPt(fJet.pt());
111  values.setJetEta(fJet.eta());
112  values.setJetPhi(fJet.phi());
113  values.setJetE(fJet.energy());
114 
115  const reco::SoftLeptonTagInfo& sltMuon =
116  (*bTagInfoMuon_)[getBTagInfoIndex(refToRawJet,*bTagInfoMuon_)];
117  if (sltMuon.leptons()>0) {
118  edm::RefToBase<reco::Track> trackRef = sltMuon.lepton(0);
119  values.setLepPx(trackRef->px());
120  values.setLepPy(trackRef->py());
121  values.setLepPz(trackRef->pz());
123  return corrector_->getCorrection(values);
124  }
125  else {
126  const reco::SoftLeptonTagInfo& sltElec =
127  (*bTagInfoElec_)[getBTagInfoIndex(refToRawJet,*bTagInfoElec_)];
128  if (sltElec.leptons()>0) {
129  edm::RefToBase<reco::Track> trackRef = sltElec.lepton(0);
130  values.setLepPx(trackRef->px());
131  values.setLepPy(trackRef->py());
132  values.setLepPz(trackRef->pz());
133  values.setAddLepToJet(false);
134  return corrector_->getCorrection(values);
135  }
136  }
137  return 1.0;
138 }
139 
140 
142 // implementation of private member functions
144 
145 //______________________________________________________________________________
147  const vector<reco::SoftLeptonTagInfo>& tags)
148  const
149 {
150  for (unsigned int i=0;i<tags.size();i++)
151  if (tags[i].jet().get()==refToRawJet.get()) return i;
152  return -1;
153 }
value_type const * get() const
Definition: RefToBase.h:225
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
TemplatedSoftLeptonTagInfo< TrackBaseRef > SoftLeptonTagInfo
Base class for all types of Jets.
Definition: Jet.h:20
std::shared_ptr< FactorizedJetCorrectorCalculator const > corrector_
reco::Particle::LorentzVector LorentzVector
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:614
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
L6SLBCorrectorImplMaker(edm::ParameterSet const &, edm::ConsumesCollector)
virtual double energy() const
energy
void addDefault(ParameterSetDescription const &psetDescription)
L6SLBCorrectorImpl(std::shared_ptr< FactorizedJetCorrectorCalculator const > corrector, edm::RefProd< std::vector< reco::SoftLeptonTagInfo >> const &bTagInfoMuon, edm::RefProd< std::vector< reco::SoftLeptonTagInfo >> const &bTagInfoElec, bool addMuonToJet)
tuple corrector
Definition: mvaPFMET_cff.py:86
edm::EDGetTokenT< std::vector< reco::SoftLeptonTagInfo > > muonToken_
int getBTagInfoIndex(const edm::RefToBase< reco::Jet > &refToRawJet, const std::vector< reco::SoftLeptonTagInfo > &tags) const
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:626
tuple tags
Definition: o2o.py:248
virtual double correction(const LorentzVector &fJet) const override
apply correction using Jet information only
edm::RefProd< std::vector< reco::SoftLeptonTagInfo > > bTagInfoMuon_
static void fillDescriptions(edm::ConfigurationDescriptions &iDescriptions)
static void addToDescription(edm::ParameterSetDescription &iDescription)
edm::RefProd< std::vector< reco::SoftLeptonTagInfo > > bTagInfoElec_
virtual double phi() const
momentum azimuthal angle
std::shared_ptr< FactorizedJetCorrectorCalculator const > getCalculator(edm::EventSetup const &, std::function< void(std::string const &)> levelCheck)
edm::EDGetTokenT< std::vector< reco::SoftLeptonTagInfo > > elecToken_
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:620
std::unique_ptr< reco::JetCorrectorImpl > make(edm::Event const &, edm::EventSetup const &)