CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
PFJetMETcorrInputProducerT< T, Textractor > Class Template Reference

#include <PFJetMETcorrInputProducerT.h>

Inheritance diagram for PFJetMETcorrInputProducerT< T, Textractor >:
edm::stream::EDProducer<>

Classes

struct  type2BinningEntryType
 

Public Member Functions

 PFJetMETcorrInputProducerT (const edm::ParameterSet &cfg)
 
 ~PFJetMETcorrInputProducerT () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void produce (edm::Event &evt, const edm::EventSetup &es) override
 

Private Attributes

double jetCorrEtaMax_
 
Textractor jetCorrExtractor_
 
edm::InputTag jetCorrLabel_
 
edm::InputTag jetCorrLabelRes_
 
edm::EDGetTokenT< reco::JetCorrectorjetCorrResToken_
 
edm::EDGetTokenT< reco::JetCorrectorjetCorrToken_
 
edm::InputTag offsetCorrLabel_
 
edm::EDGetTokenT< reco::JetCorrectoroffsetCorrToken_
 
bool skipEM_
 
double skipEMfractionThreshold_
 
bool skipMuons_
 
StringCutObjectSelector< reco::Candidate > * skipMuonSelection_
 
edm::EDGetTokenT< std::vector< T > > token_
 
double type1JetPtThreshold_
 
std::vector< type2BinningEntryType * > type2Binning_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

template<typename T, typename Textractor>
class PFJetMETcorrInputProducerT< T, Textractor >

Produce Type 1 + 2 MET corrections corresponding to differences between raw PFJets and PFJets with jet energy corrections (JECs) applied

NOTE: class is templated to that it works with reco::PFJets as well as with pat::Jets of PF-type as input

Authors
Michael Schmitt, Richard Cavanaugh, The University of Florida Florent Lacroix, University of Illinois at Chicago Christian Veelken, LLR

Definition at line 77 of file PFJetMETcorrInputProducerT.h.

Constructor & Destructor Documentation

◆ PFJetMETcorrInputProducerT()

template<typename T , typename Textractor >
PFJetMETcorrInputProducerT< T, Textractor >::PFJetMETcorrInputProducerT ( const edm::ParameterSet cfg)
inlineexplicit

Definition at line 79 of file PFJetMETcorrInputProducerT.h.

References looper::cfg, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrEtaMax_, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrLabel_, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrLabelRes_, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrResToken_, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrToken_, PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrLabel_, PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrToken_, PFJetMETcorrInputProducerT< T, Textractor >::skipEM_, PFJetMETcorrInputProducerT< T, Textractor >::skipEMfractionThreshold_, PFJetMETcorrInputProducerT< T, Textractor >::skipMuons_, PFJetMETcorrInputProducerT< T, Textractor >::skipMuonSelection_, AlCaHLTBitMon_QueryRunRegistry::string, PFJetMETcorrInputProducerT< T, Textractor >::token_, PFJetMETcorrInputProducerT< T, Textractor >::type1JetPtThreshold_, and PFJetMETcorrInputProducerT< T, Textractor >::type2Binning_.

79  : skipMuonSelection_(nullptr) {
80  token_ = consumes<std::vector<T> >(cfg.getParameter<edm::InputTag>("src"));
81  offsetCorrLabel_ = cfg.getParameter<edm::InputTag>("offsetCorrLabel");
82  offsetCorrToken_ = consumes<reco::JetCorrector>(offsetCorrLabel_);
83  jetCorrLabel_ = cfg.getParameter<edm::InputTag>("jetCorrLabel"); //for MC
84  jetCorrLabelRes_ = cfg.getParameter<edm::InputTag>("jetCorrLabelRes"); //for data
85  jetCorrToken_ = mayConsume<reco::JetCorrector>(jetCorrLabel_);
86  jetCorrResToken_ = mayConsume<reco::JetCorrector>(jetCorrLabelRes_);
87 
88  jetCorrEtaMax_ = cfg.getParameter<double>("jetCorrEtaMax");
89 
90  type1JetPtThreshold_ = cfg.getParameter<double>("type1JetPtThreshold");
91 
92  skipEM_ = cfg.getParameter<bool>("skipEM");
93  if (skipEM_) {
94  skipEMfractionThreshold_ = cfg.getParameter<double>("skipEMfractionThreshold");
95  }
96 
97  skipMuons_ = cfg.getParameter<bool>("skipMuons");
98  if (skipMuons_) {
99  std::string skipMuonSelection_string = cfg.getParameter<std::string>("skipMuonSelection");
100  skipMuonSelection_ = new StringCutObjectSelector<reco::Candidate>(skipMuonSelection_string, true);
101  }
102 
103  if (cfg.exists("type2Binning")) {
104  typedef std::vector<edm::ParameterSet> vParameterSet;
105  vParameterSet cfgType2Binning = cfg.getParameter<vParameterSet>("type2Binning");
106  for (vParameterSet::const_iterator cfgType2BinningEntry = cfgType2Binning.begin();
107  cfgType2BinningEntry != cfgType2Binning.end();
108  ++cfgType2BinningEntry) {
109  type2Binning_.push_back(new type2BinningEntryType(*cfgType2BinningEntry));
110  }
111  } else {
112  type2Binning_.push_back(new type2BinningEntryType());
113  }
114 
115  produces<CorrMETData>("type1");
116  for (typename std::vector<type2BinningEntryType*>::const_iterator type2BinningEntry = type2Binning_.begin();
117  type2BinningEntry != type2Binning_.end();
118  ++type2BinningEntry) {
119  produces<CorrMETData>((*type2BinningEntry)->getInstanceLabel_full("type2"));
120  produces<CorrMETData>((*type2BinningEntry)->getInstanceLabel_full("offset"));
121  }
122  }
edm::EDGetTokenT< reco::JetCorrector > jetCorrResToken_
edm::EDGetTokenT< reco::JetCorrector > jetCorrToken_
StringCutObjectSelector< reco::Candidate > * skipMuonSelection_
edm::EDGetTokenT< reco::JetCorrector > offsetCorrToken_
std::vector< type2BinningEntryType * > type2Binning_
edm::EDGetTokenT< std::vector< T > > token_

◆ ~PFJetMETcorrInputProducerT()

template<typename T , typename Textractor >
PFJetMETcorrInputProducerT< T, Textractor >::~PFJetMETcorrInputProducerT ( )
inlineoverride

Definition at line 123 of file PFJetMETcorrInputProducerT.h.

References PFJetMETcorrInputProducerT< T, Textractor >::skipMuonSelection_, and PFJetMETcorrInputProducerT< T, Textractor >::type2Binning_.

123  {
124  delete skipMuonSelection_;
125 
126  for (typename std::vector<type2BinningEntryType*>::const_iterator it = type2Binning_.begin();
127  it != type2Binning_.end();
128  ++it) {
129  delete (*it);
130  }
131  }
StringCutObjectSelector< reco::Candidate > * skipMuonSelection_
std::vector< type2BinningEntryType * > type2Binning_

Member Function Documentation

◆ fillDescriptions()

template<typename T , typename Textractor >
static void PFJetMETcorrInputProducerT< T, Textractor >::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
inlinestatic

Definition at line 133 of file PFJetMETcorrInputProducerT.h.

References edm::ConfigurationDescriptions::add(), defaultModuleLabel(), submitPVResolutionJobs::desc, HLT_2022v15_cff::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

133  {
135  desc.add<edm::InputTag>("src", edm::InputTag("ak4PFJetsCHS"));
136  desc.add<edm::InputTag>("offsetCorrLabel", edm::InputTag("ak4PFCHSL1FastjetCorrector"));
137  desc.add<edm::InputTag>("jetCorrLabel", edm::InputTag("ak4PFCHSL1FastL2L3Corrector"));
138  desc.add<edm::InputTag>("jetCorrLabelRes", edm::InputTag("ak4PFCHSL1FastL2L3ResidualCorrector"));
139  desc.add<double>("jetCorrEtaMax", 9.9);
140  desc.add<double>("type1JetPtThreshold", 15);
141  desc.add<bool>("skipEM", true);
142  desc.add<double>("skipEMfractionThreshold", 0.90);
143  desc.add<bool>("skipMuons", true);
144  desc.add<std::string>("skipMuonSelection", "isGlobalMuon | isStandAloneMuon");
146  }
std::string defaultModuleLabel()
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ produce()

template<typename T , typename Textractor >
void PFJetMETcorrInputProducerT< T, Textractor >::produce ( edm::Event evt,
const edm::EventSetup es 
)
inlineoverrideprivate

Definition at line 149 of file PFJetMETcorrInputProducerT.h.

References HLT_2022v15_cff::cands, edm::Event::getByToken(), edm::EventBase::isRealData(), metsig::jet, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrEtaMax_, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrExtractor_, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrLabel_, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrLabelRes_, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrResToken_, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrToken_, PDWG_EXODelayedJetMET_cff::jets, edm::InputTag::label(), eostools::move(), amptDefaultParameters_cff::mu, PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrLabel_, PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrToken_, pfDeepBoostedJetPreprocessParams_cfi::pfcand, edm::Handle< T >::product(), edm::Event::put(), PFJetMETcorrInputProducerT< T, Textractor >::skipEM_, PFJetMETcorrInputProducerT< T, Textractor >::skipEMfractionThreshold_, PFJetMETcorrInputProducerT< T, Textractor >::skipMuons_, PFJetMETcorrInputProducerT< T, Textractor >::skipMuonSelection_, PFJetMETcorrInputProducerT< T, Textractor >::token_, PFJetMETcorrInputProducerT< T, Textractor >::type1JetPtThreshold_, and PFJetMETcorrInputProducerT< T, Textractor >::type2Binning_.

149  {
150  std::unique_ptr<CorrMETData> type1Correction(new CorrMETData());
151  for (typename std::vector<type2BinningEntryType*>::iterator type2BinningEntry = type2Binning_.begin();
152  type2BinningEntry != type2Binning_.end();
153  ++type2BinningEntry) {
154  (*type2BinningEntry)->binUnclEnergySum_ = CorrMETData();
155  (*type2BinningEntry)->binOffsetEnergySum_ = CorrMETData();
156  }
157 
159  //automatic switch for residual corrections
160  if (evt.isRealData()) {
162  evt.getByToken(jetCorrResToken_, jetCorr);
163  } else {
164  evt.getByToken(jetCorrToken_, jetCorr);
165  }
166 
167  typedef std::vector<T> JetCollection;
169  evt.getByToken(token_, jets);
170 
171  int numJets = jets->size();
172  for (int jetIndex = 0; jetIndex < numJets; ++jetIndex) {
173  const T& jet = jets->at(jetIndex);
174 
176  checkInputType(jet);
177 
178  double emEnergyFraction = jet.chargedEmEnergyFraction() + jet.neutralEmEnergyFraction();
179  if (skipEM_ && emEnergyFraction > skipEMfractionThreshold_)
180  continue;
181 
182  const static PFJetMETcorrInputProducer_namespace::RawJetExtractorT<T> rawJetExtractor{};
183  reco::Candidate::LorentzVector rawJetP4 = rawJetExtractor(jet);
184 
185  if (skipMuons_) {
186  const std::vector<reco::CandidatePtr>& cands = jet.daughterPtrVector();
187  for (std::vector<reco::CandidatePtr>::const_iterator cand = cands.begin(); cand != cands.end(); ++cand) {
188  const reco::PFCandidate* pfcand = dynamic_cast<const reco::PFCandidate*>(cand->get());
189  const reco::Candidate* mu =
190  (pfcand != nullptr ? (pfcand->muonRef().isNonnull() ? pfcand->muonRef().get() : nullptr) : cand->get());
191  if (mu != nullptr && (*skipMuonSelection_)(*mu)) {
192  reco::Candidate::LorentzVector muonP4 = (*cand)->p4();
193  rawJetP4 -= muonP4;
194  }
195  }
196  }
197 
199  if (checkInputType.isPatJet(jet))
200  corrJetP4 = jetCorrExtractor_(jet, jetCorrLabel_.label(), jetCorrEtaMax_, &rawJetP4);
201  else
202  corrJetP4 = jetCorrExtractor_(jet, jetCorr.product(), jetCorrEtaMax_, &rawJetP4);
203  // retrieve additional jet energy scales in case of pat::Jets (done via specialized template defined for pat::Jets) and apply them
204  const static PFJetMETcorrInputProducer_namespace::AdditionalScalesT<T> additionalScales{};
205  corrJetP4 *= additionalScales(jet);
206 
207  if (corrJetP4.pt() > type1JetPtThreshold_) {
208  reco::Candidate::LorentzVector rawJetP4offsetCorr = rawJetP4;
209  if (!offsetCorrLabel_.label().empty()) {
211  evt.getByToken(offsetCorrToken_, offsetCorr);
212  if (checkInputType.isPatJet(jet))
213  rawJetP4offsetCorr = jetCorrExtractor_(jet, offsetCorrLabel_.label(), jetCorrEtaMax_, &rawJetP4);
214  else
215  rawJetP4offsetCorr = jetCorrExtractor_(jet, offsetCorr.product(), jetCorrEtaMax_, &rawJetP4);
216 
217  for (typename std::vector<type2BinningEntryType*>::iterator type2BinningEntry = type2Binning_.begin();
218  type2BinningEntry != type2Binning_.end();
219  ++type2BinningEntry) {
220  if (!(*type2BinningEntry)->binSelection_ || (*(*type2BinningEntry)->binSelection_)(corrJetP4)) {
221  (*type2BinningEntry)->binOffsetEnergySum_.mex += (rawJetP4.px() - rawJetP4offsetCorr.px());
222  (*type2BinningEntry)->binOffsetEnergySum_.mey += (rawJetP4.py() - rawJetP4offsetCorr.py());
223  (*type2BinningEntry)->binOffsetEnergySum_.sumet += (rawJetP4.Et() - rawJetP4offsetCorr.Et());
224  }
225  }
226  }
227 
228  //--- MET balances momentum of reconstructed particles,
229  // hence correction to jets and corresponding Type 1 MET correction are of opposite sign
230  type1Correction->mex -= (corrJetP4.px() - rawJetP4offsetCorr.px());
231  type1Correction->mey -= (corrJetP4.py() - rawJetP4offsetCorr.py());
232  type1Correction->sumet += (corrJetP4.Et() - rawJetP4offsetCorr.Et());
233  } else {
234  for (typename std::vector<type2BinningEntryType*>::iterator type2BinningEntry = type2Binning_.begin();
235  type2BinningEntry != type2Binning_.end();
236  ++type2BinningEntry) {
237  if (!(*type2BinningEntry)->binSelection_ || (*(*type2BinningEntry)->binSelection_)(corrJetP4)) {
238  (*type2BinningEntry)->binUnclEnergySum_.mex += rawJetP4.px();
239  (*type2BinningEntry)->binUnclEnergySum_.mey += rawJetP4.py();
240  (*type2BinningEntry)->binUnclEnergySum_.sumet += rawJetP4.Et();
241  }
242  }
243  }
244  }
245 
246  //--- add
247  // o Type 1 MET correction (difference corrected-uncorrected jet energy for jets of (corrected) Pt > 10 GeV)
248  // o momentum sum of "unclustered energy" (jets of (corrected) Pt < 10 GeV)
249  // o momentum sum of "offset energy" (sum of energy attributed to pile-up/underlying event)
250  // to the event
251  evt.put(std::move(type1Correction), "type1");
252  for (typename std::vector<type2BinningEntryType*>::const_iterator type2BinningEntry = type2Binning_.begin();
253  type2BinningEntry != type2Binning_.end();
254  ++type2BinningEntry) {
255  evt.put(std::unique_ptr<CorrMETData>(new CorrMETData((*type2BinningEntry)->binUnclEnergySum_)),
256  (*type2BinningEntry)->getInstanceLabel_full("type2"));
257  evt.put(std::unique_ptr<CorrMETData>(new CorrMETData((*type2BinningEntry)->binOffsetEnergySum_)),
258  (*type2BinningEntry)->getInstanceLabel_full("offset"));
259  }
260  }
edm::EDGetTokenT< reco::JetCorrector > jetCorrResToken_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
std::vector< Jet > JetCollection
Definition: Jet.h:53
edm::EDGetTokenT< reco::JetCorrector > jetCorrToken_
T const * product() const
Definition: Handle.h:70
StringCutObjectSelector< reco::Candidate > * skipMuonSelection_
edm::EDGetTokenT< reco::JetCorrector > offsetCorrToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:540
std::string const & label() const
Definition: InputTag.h:36
std::vector< type2BinningEntryType * > type2Binning_
edm::EDGetTokenT< std::vector< T > > token_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
a MET correction term
Definition: CorrMETData.h:14
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
bool isRealData() const
Definition: EventBase.h:66
long double T
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ jetCorrEtaMax_

template<typename T , typename Textractor >
double PFJetMETcorrInputProducerT< T, Textractor >::jetCorrEtaMax_
private

◆ jetCorrExtractor_

template<typename T , typename Textractor >
Textractor PFJetMETcorrInputProducerT< T, Textractor >::jetCorrExtractor_
private

◆ jetCorrLabel_

template<typename T , typename Textractor >
edm::InputTag PFJetMETcorrInputProducerT< T, Textractor >::jetCorrLabel_
private

◆ jetCorrLabelRes_

template<typename T , typename Textractor >
edm::InputTag PFJetMETcorrInputProducerT< T, Textractor >::jetCorrLabelRes_
private

◆ jetCorrResToken_

template<typename T , typename Textractor >
edm::EDGetTokenT<reco::JetCorrector> PFJetMETcorrInputProducerT< T, Textractor >::jetCorrResToken_
private

◆ jetCorrToken_

template<typename T , typename Textractor >
edm::EDGetTokenT<reco::JetCorrector> PFJetMETcorrInputProducerT< T, Textractor >::jetCorrToken_
private

◆ offsetCorrLabel_

template<typename T , typename Textractor >
edm::InputTag PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrLabel_
private

◆ offsetCorrToken_

template<typename T , typename Textractor >
edm::EDGetTokenT<reco::JetCorrector> PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrToken_
private

◆ skipEM_

template<typename T , typename Textractor >
bool PFJetMETcorrInputProducerT< T, Textractor >::skipEM_
private

◆ skipEMfractionThreshold_

template<typename T , typename Textractor >
double PFJetMETcorrInputProducerT< T, Textractor >::skipEMfractionThreshold_
private

◆ skipMuons_

template<typename T , typename Textractor >
bool PFJetMETcorrInputProducerT< T, Textractor >::skipMuons_
private

◆ skipMuonSelection_

template<typename T , typename Textractor >
StringCutObjectSelector<reco::Candidate>* PFJetMETcorrInputProducerT< T, Textractor >::skipMuonSelection_
private

◆ token_

template<typename T , typename Textractor >
edm::EDGetTokenT<std::vector<T> > PFJetMETcorrInputProducerT< T, Textractor >::token_
private

◆ type1JetPtThreshold_

template<typename T , typename Textractor >
double PFJetMETcorrInputProducerT< T, Textractor >::type1JetPtThreshold_
private

◆ type2Binning_

template<typename T , typename Textractor >
std::vector<type2BinningEntryType*> PFJetMETcorrInputProducerT< T, Textractor >::type2Binning_
private