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