CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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::JetCorrector
jetCorrResToken_
 
edm::EDGetTokenT
< reco::JetCorrector
jetCorrToken_
 
std::string moduleLabel_
 
edm::InputTag offsetCorrLabel_
 
edm::EDGetTokenT
< reco::JetCorrector
offsetCorrToken_
 
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 67 of file PFJetMETcorrInputProducerT.h.

Constructor & Destructor Documentation

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

Definition at line 69 of file PFJetMETcorrInputProducerT.h.

References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), 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_.

70  : moduleLabel_(cfg.getParameter<std::string>("@module_label")),
71  offsetCorrLabel_(""),
72  skipMuonSelection_(nullptr) {
73  token_ = consumes<std::vector<T> >(cfg.getParameter<edm::InputTag>("src"));
74 
75  if (cfg.exists("offsetCorrLabel")) {
76  offsetCorrLabel_ = cfg.getParameter<edm::InputTag>("offsetCorrLabel");
77  offsetCorrToken_ = consumes<reco::JetCorrector>(offsetCorrLabel_);
78  }
79  jetCorrLabel_ = cfg.getParameter<edm::InputTag>("jetCorrLabel"); //for MC
80  jetCorrLabelRes_ = cfg.getParameter<edm::InputTag>("jetCorrLabelRes"); //for data
81  jetCorrToken_ = mayConsume<reco::JetCorrector>(jetCorrLabel_);
82  jetCorrResToken_ = mayConsume<reco::JetCorrector>(jetCorrLabelRes_);
83 
84  jetCorrEtaMax_ = (cfg.exists("jetCorrEtaMax")) ? cfg.getParameter<double>("jetCorrEtaMax") : 9.9;
85 
86  type1JetPtThreshold_ = cfg.getParameter<double>("type1JetPtThreshold");
87 
88  skipEM_ = cfg.getParameter<bool>("skipEM");
89  if (skipEM_) {
90  skipEMfractionThreshold_ = cfg.getParameter<double>("skipEMfractionThreshold");
91  }
92 
93  skipMuons_ = cfg.getParameter<bool>("skipMuons");
94  if (skipMuons_) {
95  std::string skipMuonSelection_string = cfg.getParameter<std::string>("skipMuonSelection");
96  skipMuonSelection_ = new StringCutObjectSelector<reco::Candidate>(skipMuonSelection_string, true);
97  }
98 
99  if (cfg.exists("type2Binning")) {
100  typedef std::vector<edm::ParameterSet> vParameterSet;
101  vParameterSet cfgType2Binning = cfg.getParameter<vParameterSet>("type2Binning");
102  for (vParameterSet::const_iterator cfgType2BinningEntry = cfgType2Binning.begin();
103  cfgType2BinningEntry != cfgType2Binning.end();
104  ++cfgType2BinningEntry) {
105  type2Binning_.push_back(new type2BinningEntryType(*cfgType2BinningEntry));
106  }
107  } else {
108  type2Binning_.push_back(new type2BinningEntryType());
109  }
110 
111  produces<CorrMETData>("type1");
112  for (typename std::vector<type2BinningEntryType*>::const_iterator type2BinningEntry = type2Binning_.begin();
113  type2BinningEntry != type2Binning_.end();
114  ++type2BinningEntry) {
115  produces<CorrMETData>((*type2BinningEntry)->getInstanceLabel_full("type2"));
116  produces<CorrMETData>((*type2BinningEntry)->getInstanceLabel_full("offset"));
117  }
118  }
edm::EDGetTokenT< reco::JetCorrector > jetCorrResToken_
edm::EDGetTokenT< reco::JetCorrector > jetCorrToken_
StringCutObjectSelector< reco::Candidate > * skipMuonSelection_
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< reco::JetCorrector > offsetCorrToken_
std::vector< type2BinningEntryType * > type2Binning_
edm::EDGetTokenT< std::vector< T > > token_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
template<typename T , typename Textractor >
PFJetMETcorrInputProducerT< T, Textractor >::~PFJetMETcorrInputProducerT ( )
inlineoverride

Definition at line 119 of file PFJetMETcorrInputProducerT.h.

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

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

Member Function Documentation

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

Definition at line 129 of file PFJetMETcorrInputProducerT.h.

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

129  {
131  desc.add<edm::InputTag>("src", edm::InputTag("ak4PFJetsCHS"));
132  desc.add<edm::InputTag>("offsetCorrLabel", edm::InputTag("ak4PFCHSL1FastjetCorrector"));
133  desc.add<edm::InputTag>("jetCorrLabel", edm::InputTag("ak4PFCHSL1FastL2L3Corrector"));
134  desc.add<edm::InputTag>("jetCorrLabelRes", edm::InputTag("ak4PFCHSL1FastL2L3ResidualCorrector"));
135  desc.add<double>("jetCorrEtaMax", 9.9);
136  desc.add<double>("type1JetPtThreshold", 15);
137  desc.add<bool>("skipEM", true);
138  desc.add<double>("skipEMfractionThreshold", 0.90);
139  desc.add<bool>("skipMuons", true);
140  desc.add<std::string>("skipMuonSelection", "isGlobalMuon | isStandAloneMuon");
142  }
std::string defaultModuleLabel()
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
template<typename T , typename Textractor >
void PFJetMETcorrInputProducerT< T, Textractor >::produce ( edm::Event evt,
const edm::EventSetup es 
)
inlineoverrideprivate

Definition at line 145 of file PFJetMETcorrInputProducerT.h.

References HLT_FULL_cff::cands, edm::Ref< C, T, F >::get(), reco::LeafCandidate::get(), edm::Event::getByToken(), edm::Ref< C, T, F >::isNonnull(), 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_, fwrapper::jets, edm::InputTag::label(), eostools::move(), RPCpg::mu, reco::PFCandidate::muonRef(), PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrLabel_, PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrToken_, 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_.

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

Member Data Documentation

template<typename T , typename Textractor >
double PFJetMETcorrInputProducerT< T, Textractor >::jetCorrEtaMax_
private
template<typename T , typename Textractor >
Textractor PFJetMETcorrInputProducerT< T, Textractor >::jetCorrExtractor_
private
template<typename T , typename Textractor >
edm::InputTag PFJetMETcorrInputProducerT< T, Textractor >::jetCorrLabel_
private
template<typename T , typename Textractor >
edm::InputTag PFJetMETcorrInputProducerT< T, Textractor >::jetCorrLabelRes_
private
template<typename T , typename Textractor >
edm::EDGetTokenT<reco::JetCorrector> PFJetMETcorrInputProducerT< T, Textractor >::jetCorrResToken_
private
template<typename T , typename Textractor >
edm::EDGetTokenT<reco::JetCorrector> PFJetMETcorrInputProducerT< T, Textractor >::jetCorrToken_
private
template<typename T , typename Textractor >
std::string PFJetMETcorrInputProducerT< T, Textractor >::moduleLabel_
private
template<typename T , typename Textractor >
edm::InputTag PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrLabel_
private
template<typename T , typename Textractor >
edm::EDGetTokenT<reco::JetCorrector> PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrToken_
private
template<typename T , typename Textractor >
bool PFJetMETcorrInputProducerT< T, Textractor >::skipEM_
private
template<typename T , typename Textractor >
double PFJetMETcorrInputProducerT< T, Textractor >::skipEMfractionThreshold_
private
template<typename T , typename Textractor >
bool PFJetMETcorrInputProducerT< T, Textractor >::skipMuons_
private
template<typename T , typename Textractor >
StringCutObjectSelector<reco::Candidate>* PFJetMETcorrInputProducerT< T, Textractor >::skipMuonSelection_
private
template<typename T , typename Textractor >
edm::EDGetTokenT<std::vector<T> > PFJetMETcorrInputProducerT< T, Textractor >::token_
private
template<typename T , typename Textractor >
double PFJetMETcorrInputProducerT< T, Textractor >::type1JetPtThreshold_
private
template<typename T , typename Textractor >
std::vector<type2BinningEntryType*> PFJetMETcorrInputProducerT< T, Textractor >::type2Binning_
private