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