Go to the documentation of this file.00001 #ifndef PhysicsTools_PatUtils_SmearedJetProducerT_h
00002 #define PhysicsTools_PatUtils_SmearedJetProducerT_h
00003
00020 #include "FWCore/Framework/interface/Frameworkfwd.h"
00021 #include "FWCore/Framework/interface/EDProducer.h"
00022 #include "FWCore/Framework/interface/Event.h"
00023 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00024 #include "FWCore/Utilities/interface/InputTag.h"
00025
00026 #include "FWCore/ParameterSet/interface/FileInPath.h"
00027 #include "FWCore/Utilities/interface/Exception.h"
00028
00029 #include "DataFormats/PatCandidates/interface/Jet.h"
00030 #include "DataFormats/METReco/interface/CorrMETData.h"
00031 #include "DataFormats/Common/interface/Handle.h"
00032
00033 #include "DataFormats/JetReco/interface/GenJet.h"
00034 #include "DataFormats/JetReco/interface/GenJetCollection.h"
00035 #include "DataFormats/Math/interface/deltaR.h"
00036
00037 #include "JetMETCorrections/Type1MET/interface/JetCorrExtractorT.h"
00038 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
00039
00040 #include <TFile.h>
00041 #include <TFormula.h>
00042 #include <TH2.h>
00043 #include <TMath.h>
00044 #include <TRandom3.h>
00045 #include <TString.h>
00046
00047 #include <iostream>
00048 #include <iomanip>
00049
00050 namespace SmearedJetProducer_namespace
00051 {
00052 template <typename T>
00053 class GenJetMatcherT
00054 {
00055 public:
00056
00057 GenJetMatcherT(const edm::ParameterSet& cfg)
00058 : srcGenJets_(cfg.getParameter<edm::InputTag>("srcGenJets")),
00059 dRmaxGenJetMatch_(0)
00060 {
00061 TString dRmaxGenJetMatch_formula = cfg.getParameter<std::string>("dRmaxGenJetMatch").data();
00062 dRmaxGenJetMatch_formula.ReplaceAll("genJetPt", "x");
00063 dRmaxGenJetMatch_ = new TFormula("dRmaxGenJetMatch", dRmaxGenJetMatch_formula.Data());
00064 }
00065 ~GenJetMatcherT()
00066 {
00067 delete dRmaxGenJetMatch_;
00068 }
00069
00070 const reco::GenJet* operator()(const T& jet, edm::Event* evt = 0) const
00071 {
00072 assert(evt);
00073
00074 edm::Handle<reco::GenJetCollection> genJets;
00075 evt->getByLabel(srcGenJets_, genJets);
00076
00077 const reco::GenJet* retVal = 0;
00078
00079 double dRbestMatch = 1.e+6;
00080 for ( reco::GenJetCollection::const_iterator genJet = genJets->begin();
00081 genJet != genJets->end(); ++genJet ) {
00082 double dRmax = dRmaxGenJetMatch_->Eval(genJet->pt());
00083
00084 double dR = deltaR(jet.p4(), genJet->p4());
00085 if ( dR < dRbestMatch && dR < dRmax ) {
00086 retVal = &(*genJet);
00087 dRbestMatch = dR;
00088 }
00089 }
00090
00091 return retVal;
00092 }
00093
00094 private:
00095
00096
00097 edm::InputTag srcGenJets_;
00098
00099 TFormula* dRmaxGenJetMatch_;
00100 };
00101
00102 template <typename T>
00103 class JetResolutionExtractorT
00104 {
00105 public:
00106
00107 JetResolutionExtractorT(const edm::ParameterSet&) {}
00108 ~JetResolutionExtractorT() {}
00109
00110 double operator()(const T&) const
00111 {
00112 throw cms::Exception("SmearedJetProducer::produce")
00113 << " Jets of type other than PF not supported yet !!\n";
00114 }
00115 };
00116
00117 template <typename T>
00118 class RawJetExtractorT
00119
00120 {
00121 public:
00122
00123 reco::Candidate::LorentzVector operator()(const T& jet) const
00124 {
00125 return jet.p4();
00126 }
00127 };
00128
00129 template <>
00130 class RawJetExtractorT<pat::Jet>
00131 {
00132 public:
00133
00134 reco::Candidate::LorentzVector operator()(const pat::Jet& jet) const
00135 {
00136 if ( jet.jecSetsAvailable() ) return jet.correctedP4("Uncorrected");
00137 else return jet.p4();
00138 }
00139 };
00140 }
00141
00142 template <typename T, typename Textractor>
00143 class SmearedJetProducerT : public edm::EDProducer
00144 {
00145 typedef std::vector<T> JetCollection;
00146
00147 public:
00148
00149 explicit SmearedJetProducerT(const edm::ParameterSet& cfg)
00150 : moduleLabel_(cfg.getParameter<std::string>("@module_label")),
00151 genJetMatcher_(cfg),
00152 jetResolutionExtractor_(cfg.getParameter<edm::ParameterSet>("jetResolutions")),
00153 skipJetSelection_(0)
00154 {
00155
00156
00157
00158 src_ = cfg.getParameter<edm::InputTag>("src");
00159
00160 edm::FileInPath inputFileName = cfg.getParameter<edm::FileInPath>("inputFileName");
00161 std::string lutName = cfg.getParameter<std::string>("lutName");
00162 if ( !inputFileName.isLocal() )
00163 throw cms::Exception("JetMETsmearInputProducer")
00164 << " Failed to find File = " << inputFileName << " !!\n";
00165
00166 inputFile_ = new TFile(inputFileName.fullPath().data());
00167 lut_ = dynamic_cast<TH2*>(inputFile_->Get(lutName.data()));
00168 if ( !lut_ )
00169 throw cms::Exception("SmearedJetProducer")
00170 << " Failed to load LUT = " << lutName.data() << " from file = " << inputFileName.fullPath().data() << " !!\n";
00171
00172 jetCorrLabel_ = ( cfg.exists("jetCorrLabel") ) ?
00173 cfg.getParameter<std::string>("jetCorrLabel") : "";
00174 jetCorrEtaMax_ = ( cfg.exists("jetCorrEtaMax") ) ?
00175 cfg.getParameter<double>("jetCorrEtaMax") : 9.9;
00176
00177 sigmaMaxGenJetMatch_ = cfg.getParameter<double>("sigmaMaxGenJetMatch");
00178
00179 smearBy_ = ( cfg.exists("smearBy") ) ? cfg.getParameter<double>("smearBy") : 1.0;
00180
00181
00182 shiftBy_ = ( cfg.exists("shiftBy") ) ? cfg.getParameter<double>("shiftBy") : 0.;
00183
00184
00185 if ( cfg.exists("skipJetSelection") ) {
00186 std::string skipJetSelection_string = cfg.getParameter<std::string>("skipJetSelection");
00187 skipJetSelection_ = new StringCutObjectSelector<T>(skipJetSelection_string);
00188 }
00189
00190 skipRawJetPtThreshold_ = ( cfg.exists("skipRawJetPtThreshold") ) ?
00191 cfg.getParameter<double>("skipRawJetPtThreshold") : 1.e-2;
00192 skipCorrJetPtThreshold_ = ( cfg.exists("skipCorrJetPtThreshold") ) ?
00193 cfg.getParameter<double>("skipCorrJetPtThreshold") : 1.e-2;
00194
00195 verbosity_ = ( cfg.exists("verbosity") ) ?
00196 cfg.getParameter<int>("verbosity") : 0;
00197
00198 produces<JetCollection>();
00199 }
00200 ~SmearedJetProducerT()
00201 {
00202 delete skipJetSelection_;
00203 delete inputFile_;
00204 delete lut_;
00205 }
00206
00207 private:
00208
00209 virtual void produce(edm::Event& evt, const edm::EventSetup& es)
00210 {
00211 if ( verbosity_ ) {
00212 std::cout << "<SmearedJetProducerT::produce>:" << std::endl;
00213 std::cout << " moduleLabel = " << moduleLabel_ << std::endl;
00214 std::cout << " src = " << src_.label() << std::endl;
00215 }
00216
00217 std::auto_ptr<JetCollection> smearedJets(new JetCollection);
00218
00219 edm::Handle<JetCollection> jets;
00220 evt.getByLabel(src_, jets);
00221
00222 int numJets = jets->size();
00223 for ( int jetIndex = 0; jetIndex < numJets; ++jetIndex ) {
00224 const T& jet = jets->at(jetIndex);
00225
00226 static SmearedJetProducer_namespace::RawJetExtractorT<T> rawJetExtractor;
00227 reco::Candidate::LorentzVector rawJetP4 = rawJetExtractor(jet);
00228 if ( verbosity_ ) {
00229 std::cout << "rawJet: Pt = " << rawJetP4.pt() << ", eta = " << rawJetP4.eta() << ", phi = " << rawJetP4.phi() << std::endl;
00230 }
00231
00232 reco::Candidate::LorentzVector corrJetP4 = jet.p4();
00233 if ( jetCorrLabel_ != "" ) corrJetP4 = jetCorrExtractor_(jet, jetCorrLabel_, &evt, &es, jetCorrEtaMax_, &rawJetP4);
00234 if ( verbosity_ ) {
00235 std::cout << "corrJet: Pt = " << corrJetP4.pt() << ", eta = " << corrJetP4.eta() << ", phi = " << corrJetP4.phi() << std::endl;
00236 }
00237
00238 double smearFactor = 1.;
00239 double x = TMath::Abs(corrJetP4.eta());
00240 double y = corrJetP4.pt();
00241 if ( x > lut_->GetXaxis()->GetXmin() && x < lut_->GetXaxis()->GetXmax() &&
00242 y > lut_->GetYaxis()->GetXmin() && y < lut_->GetYaxis()->GetXmax() ) {
00243 int binIndex = lut_->FindBin(x, y);
00244
00245 if ( smearBy_ > 0. ) smearFactor += smearBy_*(lut_->GetBinContent(binIndex) - 1.);
00246 double smearFactorErr = lut_->GetBinError(binIndex);
00247 if ( verbosity_ ) std::cout << "smearFactor = " << smearFactor << " +/- " << smearFactorErr << std::endl;
00248
00249 if ( shiftBy_ != 0. ) {
00250 smearFactor += (shiftBy_*smearFactorErr);
00251 if ( verbosity_ ) std::cout << "smearFactor(shifted) = " << smearFactor << std::endl;
00252 }
00253 }
00254
00255 double smearedJetEn = jet.energy();
00256 double sigmaEn = jetResolutionExtractor_(jet)*TMath::Sqrt(smearFactor*smearFactor - 1.);
00257 const reco::GenJet* genJet = genJetMatcher_(jet, &evt);
00258 bool isGenMatched = false;
00259 if ( genJet ) {
00260 if ( verbosity_ ) {
00261 std::cout << "genJet: Pt = " << genJet->pt() << ", eta = " << genJet->eta() << ", phi = " << genJet->phi() << std::endl;
00262 }
00263 double dEn = corrJetP4.E() - genJet->energy();
00264 if ( dEn < (sigmaMaxGenJetMatch_*sigmaEn) ) {
00265
00266
00267
00268 if ( verbosity_ ) {
00269 std::cout << " successfully matched to genJet" << std::endl;
00270 std::cout << "corrJetEn = " << corrJetP4.E() << ", genJetEn = " << genJet->energy() << " --> dEn = " << dEn << std::endl;
00271 }
00272
00273 smearedJetEn = jet.energy()*(1. + (smearFactor - 1.)*dEn/TMath::Max(rawJetP4.E(), corrJetP4.E()));
00274 isGenMatched = true;
00275 }
00276 }
00277 if ( !isGenMatched ) {
00278
00279
00280
00281 if ( verbosity_ ) {
00282 std::cout << " not matched to genJet" << std::endl;
00283 std::cout << "corrJetEn = " << corrJetP4.E() << ", sigmaEn = " << sigmaEn << std::endl;
00284 }
00285
00286 if ( smearFactor > 1. ) {
00287
00288
00289
00290
00291
00292
00293 smearedJetEn = jet.energy()*(1. + rnd_.Gaus(0., sigmaEn)/TMath::Max(rawJetP4.E(), corrJetP4.E()));
00294 }
00295 }
00296
00297
00298 const double minJetEn = 1.e-2;
00299 if ( smearedJetEn < minJetEn ) smearedJetEn = minJetEn;
00300
00301
00302
00303
00304
00305 reco::Candidate::LorentzVector smearedJetP4 = jet.p4();
00306 if ( !((skipJetSelection_ && (*skipJetSelection_)(jet)) ||
00307 rawJetP4.pt() < skipRawJetPtThreshold_ ||
00308 corrJetP4.pt() < skipCorrJetPtThreshold_ ) ) {
00309 if ( verbosity_ ) {
00310 std::cout << " smearing jetP4 by factor = " << (smearedJetEn/jet.energy()) << " --> smearedJetEn = " << smearedJetEn << std::endl;
00311 }
00312 smearedJetP4 *= (smearedJetEn/jet.energy());
00313 }
00314
00315 if ( verbosity_ ) {
00316 std::cout << "smearedJet: Pt = " << smearedJetP4.pt() << ", eta = " << smearedJetP4.eta() << ", phi = " << smearedJetP4.phi() << std::endl;
00317 std::cout << " dPt = " << (smearedJetP4.pt() - jet.pt())
00318 << " (Px = " << (smearedJetP4.px() - jet.px()) << ", Py = " << (smearedJetP4.py() - jet.py()) << ")" << std::endl;
00319 }
00320
00321 T smearedJet = (jet);
00322 smearedJet.setP4(smearedJetP4);
00323
00324 smearedJets->push_back(smearedJet);
00325 }
00326
00327
00328 evt.put(smearedJets);
00329 }
00330
00331 std::string moduleLabel_;
00332
00333 SmearedJetProducer_namespace::GenJetMatcherT<T> genJetMatcher_;
00334
00335
00336
00337
00338 edm::InputTag src_;
00339
00340 TFile* inputFile_;
00341 TH2* lut_;
00342
00343 SmearedJetProducer_namespace::JetResolutionExtractorT<T> jetResolutionExtractor_;
00344 TRandom3 rnd_;
00345
00346 std::string jetCorrLabel_;
00347 double jetCorrEtaMax_;
00348
00349
00350
00351
00352 Textractor jetCorrExtractor_;
00353
00354 double sigmaMaxGenJetMatch_;
00355
00356
00357
00358 double smearBy_;
00359
00360 double shiftBy_;
00361
00362 StringCutObjectSelector<T>* skipJetSelection_;
00363 double skipRawJetPtThreshold_;
00364 double skipCorrJetPtThreshold_;
00365
00366 int verbosity_;
00367 };
00368
00369 #endif