#include <CaloJetMETcorrInputProducerT.h>
Public Member Functions | |
CaloJetMETcorrInputProducerT (const edm::ParameterSet &cfg) | |
~CaloJetMETcorrInputProducerT () | |
Private Member Functions | |
void | produce (edm::Event &evt, const edm::EventSetup &es) |
Private Attributes | |
double | jetCorrEtaMax_ |
Textractor | jetCorrExtractor_ |
std::string | jetCorrLabel_ |
std::string | moduleLabel_ |
std::string | offsetCorrLabel_ |
bool | skipEM_ |
double | skipEMfractionThreshold_ |
edm::InputTag | src_ |
edm::InputTag | srcMET_ |
double | type1JetPtThreshold_ |
Definition at line 60 of file CaloJetMETcorrInputProducerT.h.
CaloJetMETcorrInputProducerT< T, Textractor >::CaloJetMETcorrInputProducerT | ( | const edm::ParameterSet & | cfg | ) | [inline, explicit] |
Definition at line 64 of file CaloJetMETcorrInputProducerT.h.
References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), CaloJetMETcorrInputProducerT< T, Textractor >::jetCorrEtaMax_, CaloJetMETcorrInputProducerT< T, Textractor >::jetCorrLabel_, CaloJetMETcorrInputProducerT< T, Textractor >::offsetCorrLabel_, CaloJetMETcorrInputProducerT< T, Textractor >::skipEM_, CaloJetMETcorrInputProducerT< T, Textractor >::skipEMfractionThreshold_, CaloJetMETcorrInputProducerT< T, Textractor >::src_, CaloJetMETcorrInputProducerT< T, Textractor >::srcMET_, and CaloJetMETcorrInputProducerT< T, Textractor >::type1JetPtThreshold_.
: moduleLabel_(cfg.getParameter<std::string>("@module_label")) { src_ = cfg.getParameter<edm::InputTag>("src"); offsetCorrLabel_ = ( cfg.exists("offsetCorrLabel") ) ? cfg.getParameter<std::string>("offsetCorrLabel") : ""; jetCorrLabel_ = cfg.getParameter<std::string>("jetCorrLabel"); jetCorrEtaMax_ = ( cfg.exists("jetCorrEtaMax") ) ? cfg.getParameter<double>("jetCorrEtaMax") : 9.9; type1JetPtThreshold_ = cfg.getParameter<double>("type1JetPtThreshold"); skipEM_ = cfg.getParameter<bool>("skipEM"); if ( skipEM_ ) { skipEMfractionThreshold_ = cfg.getParameter<double>("skipEMfractionThreshold"); } if ( cfg.exists("srcMET") ) { srcMET_ = cfg.getParameter<edm::InputTag>("srcMET"); } produces<CorrMETData>("type1"); produces<CorrMETData>("type2"); produces<CorrMETData>("offset"); }
CaloJetMETcorrInputProducerT< T, Textractor >::~CaloJetMETcorrInputProducerT | ( | ) | [inline] |
Definition at line 91 of file CaloJetMETcorrInputProducerT.h.
{}
void CaloJetMETcorrInputProducerT< T, Textractor >::produce | ( | edm::Event & | evt, |
const edm::EventSetup & | es | ||
) | [inline, private, virtual] |
Implements edm::EDProducer.
Definition at line 95 of file CaloJetMETcorrInputProducerT.h.
References edm::Event::getByLabel(), CaloJetMETcorrInputProducerT< T, Textractor >::jetCorrEtaMax_, CaloJetMETcorrInputProducerT< T, Textractor >::jetCorrExtractor_, CaloJetMETcorrInputProducerT< T, Textractor >::jetCorrLabel_, analyzePatCleaning_cfg::jets, edm::InputTag::label(), CaloMET_cfi::met, CaloJetMETcorrInputProducerT< T, Textractor >::offsetCorrLabel_, edm::Event::put(), CaloJetMETcorrInputProducerT< T, Textractor >::skipEM_, CaloJetMETcorrInputProducerT< T, Textractor >::skipEMfractionThreshold_, CaloJetMETcorrInputProducerT< T, Textractor >::src_, CaloJetMETcorrInputProducerT< T, Textractor >::srcMET_, and CaloJetMETcorrInputProducerT< T, Textractor >::type1JetPtThreshold_.
{ std::auto_ptr<CorrMETData> type1Correction(new CorrMETData()); std::auto_ptr<CorrMETData> unclEnergySum(new CorrMETData()); std::auto_ptr<CorrMETData> offsetEnergySum(new CorrMETData()); typedef std::vector<T> JetCollection; edm::Handle<JetCollection> jets; evt.getByLabel(src_, jets); typedef edm::View<reco::MET> METView; edm::Handle<METView> met; if ( srcMET_.label() != "" ) { evt.getByLabel(srcMET_, met); if ( met->size() != 1 ) throw cms::Exception("CaloJetMETcorrInputProducer::produce") << "Failed to find unique MET in the event, src = " << srcMET_.label() << " !!\n"; //--- compute "unclustered energy" by sutracting from the reconstructed MET // (i.e. from the negative vectorial sum of all particles reconstructed in the event) // the momenta of (high Pt) jets which enter Type 1 MET corrections // // NOTE: MET = -(jets + muons + "unclustered energy"), // so "unclustered energy" = -(MET + jets + muons), // i.e. (high Pt) jets enter the sum of "unclustered energy" with negative sign. // unclEnergySum->mex = -met->front().px(); unclEnergySum->mey = -met->front().py(); unclEnergySum->sumet = met->front().sumEt(); } int numJets = jets->size(); for ( int jetIndex = 0; jetIndex < numJets; ++jetIndex ) { const T& rawJet = jets->at(jetIndex); static CaloJetMETcorrInputProducer_namespace::InputTypeCheckerT<T> checkInputType; checkInputType(rawJet); static CaloJetMETcorrInputProducer_namespace::RawJetExtractorT<T> rawJetExtractor; reco::Candidate::LorentzVector rawJetP4 = rawJetExtractor(rawJet); reco::Candidate::LorentzVector corrJetP4 = jetCorrExtractor_(rawJet, jetCorrLabel_, &evt, &es, jetCorrEtaMax_); if ( corrJetP4.pt() > type1JetPtThreshold_ ) { unclEnergySum->mex -= rawJetP4.px(); unclEnergySum->mey -= rawJetP4.py(); unclEnergySum->sumet -= rawJetP4.Et(); if ( skipEM_ && rawJet.emEnergyFraction() > skipEMfractionThreshold_ ) continue; reco::Candidate::LorentzVector rawJetP4offsetCorr = rawJetP4; if ( offsetCorrLabel_ != "" ) { rawJetP4offsetCorr = jetCorrExtractor_(rawJet, offsetCorrLabel_, &evt, &es, jetCorrEtaMax_); offsetEnergySum->mex += (rawJetP4.px() - rawJetP4offsetCorr.px()); offsetEnergySum->mey += (rawJetP4.py() - rawJetP4offsetCorr.py()); offsetEnergySum->sumet += (rawJetP4.Et() - rawJetP4offsetCorr.Et()); } //--- MET balances momentum of reconstructed particles, // hence correction to jets and corresponding Type 1 MET correction are of opposite sign type1Correction->mex -= (corrJetP4.px() - rawJetP4offsetCorr.px()); type1Correction->mey -= (corrJetP4.py() - rawJetP4offsetCorr.py()); type1Correction->sumet += (corrJetP4.Et() - rawJetP4offsetCorr.Et()); } } //--- add // o Type 1 MET correction (difference corrected-uncorrected jet energy for jets of (corrected) Pt > 20 GeV) // o momentum sum of "unclustered energy" (jets of (corrected) Pt < 20 GeV) // o momentum sum of "offset energy" (sum of energy attributed to pile-up/underlying event) // to the event evt.put(type1Correction, "type1"); evt.put(unclEnergySum, "type2"); evt.put(offsetEnergySum, "offset"); }
double CaloJetMETcorrInputProducerT< T, Textractor >::jetCorrEtaMax_ [private] |
Definition at line 181 of file CaloJetMETcorrInputProducerT.h.
Referenced by CaloJetMETcorrInputProducerT< T, Textractor >::CaloJetMETcorrInputProducerT(), and CaloJetMETcorrInputProducerT< T, Textractor >::produce().
Textractor CaloJetMETcorrInputProducerT< T, Textractor >::jetCorrExtractor_ [private] |
Definition at line 179 of file CaloJetMETcorrInputProducerT.h.
Referenced by CaloJetMETcorrInputProducerT< T, Textractor >::produce().
std::string CaloJetMETcorrInputProducerT< T, Textractor >::jetCorrLabel_ [private] |
Definition at line 178 of file CaloJetMETcorrInputProducerT.h.
Referenced by CaloJetMETcorrInputProducerT< T, Textractor >::CaloJetMETcorrInputProducerT(), and CaloJetMETcorrInputProducerT< T, Textractor >::produce().
std::string CaloJetMETcorrInputProducerT< T, Textractor >::moduleLabel_ [private] |
Definition at line 173 of file CaloJetMETcorrInputProducerT.h.
std::string CaloJetMETcorrInputProducerT< T, Textractor >::offsetCorrLabel_ [private] |
Definition at line 177 of file CaloJetMETcorrInputProducerT.h.
Referenced by CaloJetMETcorrInputProducerT< T, Textractor >::CaloJetMETcorrInputProducerT(), and CaloJetMETcorrInputProducerT< T, Textractor >::produce().
bool CaloJetMETcorrInputProducerT< T, Textractor >::skipEM_ [private] |
Definition at line 191 of file CaloJetMETcorrInputProducerT.h.
Referenced by CaloJetMETcorrInputProducerT< T, Textractor >::CaloJetMETcorrInputProducerT(), and CaloJetMETcorrInputProducerT< T, Textractor >::produce().
double CaloJetMETcorrInputProducerT< T, Textractor >::skipEMfractionThreshold_ [private] |
Definition at line 193 of file CaloJetMETcorrInputProducerT.h.
Referenced by CaloJetMETcorrInputProducerT< T, Textractor >::CaloJetMETcorrInputProducerT(), and CaloJetMETcorrInputProducerT< T, Textractor >::produce().
edm::InputTag CaloJetMETcorrInputProducerT< T, Textractor >::src_ [private] |
Definition at line 175 of file CaloJetMETcorrInputProducerT.h.
Referenced by CaloJetMETcorrInputProducerT< T, Textractor >::CaloJetMETcorrInputProducerT(), and CaloJetMETcorrInputProducerT< T, Textractor >::produce().
edm::InputTag CaloJetMETcorrInputProducerT< T, Textractor >::srcMET_ [private] |
Definition at line 195 of file CaloJetMETcorrInputProducerT.h.
Referenced by CaloJetMETcorrInputProducerT< T, Textractor >::CaloJetMETcorrInputProducerT(), and CaloJetMETcorrInputProducerT< T, Textractor >::produce().
double CaloJetMETcorrInputProducerT< T, Textractor >::type1JetPtThreshold_ [private] |
Definition at line 187 of file CaloJetMETcorrInputProducerT.h.
Referenced by CaloJetMETcorrInputProducerT< T, Textractor >::CaloJetMETcorrInputProducerT(), and CaloJetMETcorrInputProducerT< T, Textractor >::produce().