test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SmearedJetProducerT.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_PatUtils_SmearedJetProducerT_h
2 #define PhysicsTools_PatUtils_SmearedJetProducerT_h
3 
24 
27 
31 
35 
38 
39 #include <TFile.h>
40 #include <TFormula.h>
41 #include <TH2.h>
42 #include <TMath.h>
43 #include <TRandom3.h>
44 #include <TString.h>
45 
46 #include <iostream>
47 #include <iomanip>
48 
49 namespace SmearedJetProducer_namespace
50 {
51  template <typename T>
53  {
54  public:
55 
57  : srcGenJetsToken_(iC.consumes<reco::GenJetCollection>(cfg.getParameter<edm::InputTag>("srcGenJets"))),
59  {
60  TString dRmaxGenJetMatch_formula = cfg.getParameter<std::string>("dRmaxGenJetMatch").data();
61  dRmaxGenJetMatch_formula.ReplaceAll("genJetPt", "x");
62  dRmaxGenJetMatch_ = new TFormula("dRmaxGenJetMatch", dRmaxGenJetMatch_formula.Data());
63  }
65  {
66  delete dRmaxGenJetMatch_;
67  }
68 
69  const reco::GenJet* operator()(const T& jet, edm::Event* evt = 0) const
70  {
71  assert(evt);
72 
74  evt->getByToken(srcGenJetsToken_, genJets);
75 
76  const reco::GenJet* retVal = 0;
77 
78  double dRbestMatch = 1.e+6;
79  for ( reco::GenJetCollection::const_iterator genJet = genJets->begin();
80  genJet != genJets->end(); ++genJet ) {
81  double dRmax = dRmaxGenJetMatch_->Eval(genJet->pt());
82  //std::cout << "genJetPt = " << genJet->pt() << ": dRmax = " << dRmax << std::endl;
83  double dR = deltaR(jet.p4(), genJet->p4());
84  if ( dR < dRbestMatch && dR < dRmax ) {
85  retVal = &(*genJet);
86  dRbestMatch = dR;
87  }
88  }
89 
90  return retVal;
91  }
92 
93  private:
94 
95 //--- configuration parameter
98 
99  TFormula* dRmaxGenJetMatch_;
100  };
101 
102  template <typename T>
104  {
105  public:
106 
109 
110  double operator()(const T&) const
111  {
112  throw cms::Exception("SmearedJetProducer::produce")
113  << " Jets of type other than PF not supported yet !!\n";
114  }
115  };
116 
117  template <typename T>
118  class RawJetExtractorT // this template is neccessary to support pat::Jets
119  // (because pat::Jet->p4() returns the JES corrected, not the raw, jet momentum)
120  {
121  public:
122 
125  {
126  return jet.p4();
127  }
128  };
129 
130  template <>
131  class RawJetExtractorT<pat::Jet>
132  {
133  public:
134 
137  {
138  if ( jet.jecSetsAvailable() ) return jet.correctedP4("Uncorrected");
139  else return jet.p4();
140  }
141  };
142 }
143 
144 template <typename T, typename Textractor>
146 {
147  typedef std::vector<T> JetCollection;
148 
149  public:
150 
152  : moduleLabel_(cfg.getParameter<std::string>("@module_label")),
154  jetResolutionExtractor_(cfg.getParameter<edm::ParameterSet>("jetResolutions")),
156  {
157  //std::cout << "<SmearedJetProducer::SmearedJetProducer>:" << std::endl;
158  //std::cout << " moduleLabel = " << moduleLabel_ << std::endl;
159 
160  src_ = cfg.getParameter<edm::InputTag>("src");
161  srcToken_ = consumes<JetCollection>(src_);
162 
164  std::string lutName = cfg.getParameter<std::string>("lutName");
165  if (inputFileName.location() == edm::FileInPath::Unknown)
166  throw cms::Exception("JetMETsmearInputProducer")
167  << " Failed to find File = " << inputFileName << " !!\n";
168 
169  inputFile_ = new TFile(inputFileName.fullPath().data());
170  lut_ = dynamic_cast<TH2*>(inputFile_->Get(lutName.data()));
171  if ( !lut_ )
172  throw cms::Exception("SmearedJetProducer")
173  << " Failed to load LUT = " << lutName.data() << " from file = " << inputFileName.fullPath().data() << " !!\n";
174 
175  jetCorrLabel_ = ( cfg.exists("jetCorrLabel") ) ?
176  cfg.getParameter<std::string>("jetCorrLabel") : "";
177  jetCorrEtaMax_ = ( cfg.exists("jetCorrEtaMax") ) ?
178  cfg.getParameter<double>("jetCorrEtaMax") : 9.9;
179 
180  sigmaMaxGenJetMatch_ = cfg.getParameter<double>("sigmaMaxGenJetMatch");
181 
182  smearBy_ = ( cfg.exists("smearBy") ) ? cfg.getParameter<double>("smearBy") : 1.0;
183  //std::cout << "smearBy = " << smearBy_ << std::endl;
184 
185  shiftBy_ = ( cfg.exists("shiftBy") ) ? cfg.getParameter<double>("shiftBy") : 0.;
186  //std::cout << "shiftBy = " << shiftBy_ << std::endl;
187 
188  if ( cfg.exists("skipJetSelection") ) {
189  std::string skipJetSelection_string = cfg.getParameter<std::string>("skipJetSelection");
190  skipJetSelection_ = new StringCutObjectSelector<T>(skipJetSelection_string);
191  }
192 
193  skipRawJetPtThreshold_ = ( cfg.exists("skipRawJetPtThreshold") ) ?
194  cfg.getParameter<double>("skipRawJetPtThreshold") : 1.e-2;
195  skipCorrJetPtThreshold_ = ( cfg.exists("skipCorrJetPtThreshold") ) ?
196  cfg.getParameter<double>("skipCorrJetPtThreshold") : 1.e-2;
197 
198  verbosity_ = ( cfg.exists("verbosity") ) ?
199  cfg.getParameter<int>("verbosity") : 0;
200 
201  produces<JetCollection>();
202  }
204  {
205  delete skipJetSelection_;
206  delete inputFile_;
207  delete lut_;
208  }
209 
210  private:
211 
212  virtual void produce(edm::Event& evt, const edm::EventSetup& es)
213  {
214  if ( verbosity_ ) {
215  std::cout << "<SmearedJetProducerT::produce>:" << std::endl;
216  std::cout << " moduleLabel = " << moduleLabel_ << std::endl;
217  std::cout << " src = " << src_.label() << std::endl;
218  }
219 
220  std::auto_ptr<JetCollection> smearedJets(new JetCollection);
221 
223  evt.getByToken(srcToken_, jets);
224 
225  int numJets = jets->size();
226  for ( int jetIndex = 0; jetIndex < numJets; ++jetIndex ) {
227  const T& jet = jets->at(jetIndex);
228 
229  static const SmearedJetProducer_namespace::RawJetExtractorT<T> rawJetExtractor;
230  reco::Candidate::LorentzVector rawJetP4 = rawJetExtractor(jet);
231  if ( verbosity_ ) {
232  std::cout << "rawJet: Pt = " << rawJetP4.pt() << ", eta = " << rawJetP4.eta() << ", phi = " << rawJetP4.phi() << std::endl;
233  }
234 
235  reco::Candidate::LorentzVector corrJetP4 = jet.p4();
236  if ( jetCorrLabel_ != "" ) corrJetP4 = jetCorrExtractor_(jet, jetCorrLabel_, &evt, &es, jetCorrEtaMax_, &rawJetP4);
237  if ( verbosity_ ) {
238  std::cout << "corrJet: Pt = " << corrJetP4.pt() << ", eta = " << corrJetP4.eta() << ", phi = " << corrJetP4.phi() << std::endl;
239  }
240 
241  double smearFactor = 1.;
242  double x = TMath::Abs(corrJetP4.eta());
243  double y = corrJetP4.pt();
244  if ( x > lut_->GetXaxis()->GetXmin() && x < lut_->GetXaxis()->GetXmax() &&
245  y > lut_->GetYaxis()->GetXmin() && y < lut_->GetYaxis()->GetXmax() ) {
246  int binIndex = lut_->FindBin(x, y);
247 
248  if ( smearBy_ > 0. ) smearFactor += smearBy_*(lut_->GetBinContent(binIndex) - 1.);
249  double smearFactorErr = lut_->GetBinError(binIndex);
250  if ( verbosity_ ) std::cout << "smearFactor = " << smearFactor << " +/- " << smearFactorErr << std::endl;
251 
252  if ( shiftBy_ != 0. ) {
253  smearFactor += (shiftBy_*smearFactorErr);
254  if ( verbosity_ ) std::cout << "smearFactor(shifted) = " << smearFactor << std::endl;
255  }
256  }
257 
258  double smearedJetEn = jet.energy();
259  double sigmaEn = jetResolutionExtractor_(jet)*TMath::Sqrt(smearFactor*smearFactor - 1.);
260  const reco::GenJet* genJet = genJetMatcher_(jet, &evt);
261  bool isGenMatched = false;
262  if ( genJet ) {
263  if ( verbosity_ ) {
264  std::cout << "genJet: Pt = " << genJet->pt() << ", eta = " << genJet->eta() << ", phi = " << genJet->phi() << std::endl;
265  }
266  double dEn = corrJetP4.E() - genJet->energy();
267  if ( dEn < (sigmaMaxGenJetMatch_*sigmaEn) ) {
268 //--- case 1: reconstructed jet matched to generator level jet,
269 // smear difference between reconstructed and "true" jet energy
270 
271  if ( verbosity_ ) {
272  std::cout << " successfully matched to genJet" << std::endl;
273  std::cout << "corrJetEn = " << corrJetP4.E() << ", genJetEn = " << genJet->energy() << " --> dEn = " << dEn << std::endl;
274  }
275 
276  smearedJetEn = jet.energy()*(1. + (smearFactor - 1.)*dEn/TMath::Max(rawJetP4.E(), corrJetP4.E()));
277  isGenMatched = true;
278  }
279  }
280  if ( !isGenMatched ) {
281 //--- case 2: reconstructed jet **not** matched to generator level jet,
282 // smear jet energy using MC resolution functions implemented in PFMEt significance algorithm (CMS AN-10/400)
283 
284  if ( verbosity_ ) {
285  std::cout << " not matched to genJet" << std::endl;
286  std::cout << "corrJetEn = " << corrJetP4.E() << ", sigmaEn = " << sigmaEn << std::endl;
287  }
288 
289  if ( smearFactor > 1. ) {
290  // CV: MC resolution already accounted for in reconstructed jet,
291  // add additional Gaussian smearing of width = sqrt(smearFactor^2 - 1)
292  // to account for Data/MC **difference** in jet resolutions.
293  // Take maximum(rawJetEn, corrJetEn) to avoid pathological cases
294  // (e.g. corrJetEn << rawJetEn, due to L1Fastjet corrections)
295 
296  smearedJetEn = jet.energy()*(1. + rnd_.Gaus(0., sigmaEn)/TMath::Max(rawJetP4.E(), corrJetP4.E()));
297  }
298  }
299 
300  // CV: keep minimum jet energy, in order not to loose direction information
301  const double minJetEn = 1.e-2;
302  if ( smearedJetEn < minJetEn ) smearedJetEn = minJetEn;
303 
304  // CV: skip smearing in case either "raw" or "corrected" jet energy is very low
305  // or jet passes selection configurable via python
306  // (allows for protection against "pathological cases",
307  // cf. PhysicsTools/PatUtils/python/tools/metUncertaintyTools.py)
308  reco::Candidate::LorentzVector smearedJetP4 = jet.p4();
309  if ( !((skipJetSelection_ && (*skipJetSelection_)(jet)) ||
310  rawJetP4.pt() < skipRawJetPtThreshold_ ||
311  corrJetP4.pt() < skipCorrJetPtThreshold_ ) ) {
312  if ( verbosity_ ) {
313  std::cout << " smearing jetP4 by factor = " << (smearedJetEn/jet.energy()) << " --> smearedJetEn = " << smearedJetEn << std::endl;
314  }
315  smearedJetP4 *= (smearedJetEn/jet.energy());
316  }
317 
318  if ( verbosity_ ) {
319  std::cout << "smearedJet: Pt = " << smearedJetP4.pt() << ", eta = " << smearedJetP4.eta() << ", phi = " << smearedJetP4.phi() << std::endl;
320  std::cout << " dPt = " << (smearedJetP4.pt() - jet.pt())
321  << " (Px = " << (smearedJetP4.px() - jet.px()) << ", Py = " << (smearedJetP4.py() - jet.py()) << ")" << std::endl;
322  }
323 
324  T smearedJet = (jet);
325  smearedJet.setP4(smearedJetP4);
326 
327  smearedJets->push_back(smearedJet);
328  }
329 
330 //--- add collection of "smeared" jets to the event
331  evt.put(smearedJets);
332  }
333 
335 
337 
338 //--- configuration parameters
339 
340  // collection of pat::Jets (with L2L3/L2L3Residual corrections applied)
343 
344  TFile* inputFile_;
345  TH2* lut_;
346 
348  TRandom3 rnd_;
349 
350  std::string jetCorrLabel_; // e.g. 'ak5PFJetL1FastL2L3' (reco::PFJets) / '' (pat::Jets)
351  double jetCorrEtaMax_; // do not use JEC factors for |eta| above this threshold (recommended default = 4.7),
352  // in order to work around problem with CMSSW_4_2_x JEC factors at high eta,
353  // reported in
354  // https://hypernews.cern.ch/HyperNews/CMS/get/jes/270.html
355  // https://hypernews.cern.ch/HyperNews/CMS/get/JetMET/1259/1.html
356  Textractor jetCorrExtractor_;
357 
358  double sigmaMaxGenJetMatch_; // maximum difference between energy of reconstructed jet and matched generator level jet
359  // (if the difference between reconstructed and generated jet energy exceeds this threshold,
360  // the jet is considered to have substantial pile-up contributions are is considered to be unmatched)
361 
362  double smearBy_; // option to "smear" jet energy by N standard-deviations, useful for template morphing
363 
364  double shiftBy_; // option to increase/decrease within uncertainties the jet energy resolution used for smearing
365 
366  StringCutObjectSelector<T>* skipJetSelection_; // jets passing this cut are **not** smeared
367  double skipRawJetPtThreshold_; // jets with transverse momenta below this value (either on "raw" or "corrected" level)
368  double skipCorrJetPtThreshold_; // are **not** smeared
369 
370  int verbosity_; // flag to enabled/disable debug output
371 };
372 
373 #endif
T getParameter(std::string const &) const
reco::Candidate::LorentzVector operator()(const pat::Jet &jet) const
edm::EDGetTokenT< JetCollection > srcToken_
StringCutObjectSelector< T > * skipJetSelection_
bool jecSetsAvailable() const
Definition: Jet.h:119
virtual float pt() const
transverse momentum
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual float phi() const
momentum azimuthal angle
std::vector< GenJet > GenJetCollection
collection of GenJet objects
bool exists(std::string const &parameterName) const
checks if a parameter exists
const LorentzVector & correctedP4(const std::string &level, const std::string &flavor="none", const std::string &set="") const
Definition: Jet.h:146
SmearedJetProducer_namespace::GenJetMatcherT< T > genJetMatcher_
SmearedJetProducerT(const edm::ParameterSet &cfg)
virtual double energy() const
energy
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
virtual float eta() const
momentum pseudorapidity
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
vector< PseudoJet > jets
T Abs(T a)
Definition: MathUtil.h:49
Jets made from MC generator particles.
Definition: GenJet.h:24
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
T Max(T a, T b)
Definition: MathUtil.h:44
GenJetMatcherT(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
SmearedJetProducer_namespace::JetResolutionExtractorT< T > jetResolutionExtractor_
Analysis-level calorimeter jet class.
Definition: Jet.h:73
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
std::string const & label() const
Definition: InputTag.h:42
const reco::GenJet * operator()(const T &jet, edm::Event *evt=0) const
edm::EDGetTokenT< reco::GenJetCollection > srcGenJetsToken_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
reco::Candidate::LorentzVector operator()(const T &jet) const
tuple cout
Definition: gather_cfg.py:121
Definition: DDAxes.h:10
long double T
std::vector< T > JetCollection
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
virtual void produce(edm::Event &evt, const edm::EventSetup &es)