File indexing completed on 2024-04-06 12:19:24
0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0003
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/Framework/interface/ESHandle.h"
0009 #include "JetMETCorrections/JetCorrector/interface/JetCorrector.h"
0010 #include "DataFormats/VertexReco/interface/Vertex.h"
0011 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0012 #include "DataFormats/JetReco/interface/CaloJet.h"
0013 #include "DataFormats/JetReco/interface/PFJet.h"
0014 #include "DataFormats/JetReco/interface/JPTJet.h"
0015 #include "FWCore/ServiceRegistry/interface/Service.h"
0016 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0017
0018 #include "TH1F.h"
0019
0020 using namespace edm;
0021 using namespace reco;
0022 using namespace std;
0023
0024
0025
0026
0027 template <class Jet>
0028 class JetCorrectorOnTheFly : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0029 public:
0030 explicit JetCorrectorOnTheFly(const edm::ParameterSet&);
0031 ~JetCorrectorOnTheFly() override;
0032
0033 private:
0034 typedef std::vector<Jet> JetCollection;
0035 void beginJob() override;
0036 void analyze(const edm::Event&, const edm::EventSetup&) override;
0037 void endJob() override;
0038
0039 edm::Service<TFileService> fs;
0040 edm::EDGetTokenT<reco::JetCorrector> mJetCorrector;
0041 edm::EDGetTokenT<JetCollection> mJetToken;
0042 edm::EDGetTokenT<reco::VertexCollection> mVertexToken;
0043 double mMinRawJetPt;
0044 bool mDebug;
0045 TH1F *mRawPt, *mCorPt;
0046 };
0047
0048
0049
0050
0051 template <class Jet>
0052 JetCorrectorOnTheFly<Jet>::JetCorrectorOnTheFly(const edm::ParameterSet& iConfig) {
0053 mJetCorrector = consumes<reco::JetCorrector>(iConfig.getParameter<edm::InputTag>("JetCorrector"));
0054 mJetToken = consumes<JetCollection>(edm::InputTag(iConfig.getParameter<edm::InputTag>("JetCollectionName")));
0055 mVertexToken = consumes<reco::VertexCollection>(edm::InputTag("offlinePrimaryVertices"));
0056 mMinRawJetPt = iConfig.getParameter<double>("MinRawJetPt");
0057 mDebug = iConfig.getParameter<bool>("Debug");
0058
0059 usesResource(TFileService::kSharedResource);
0060 }
0061
0062 template <class Jet>
0063 JetCorrectorOnTheFly<Jet>::~JetCorrectorOnTheFly() {}
0064
0065 template <class Jet>
0066 void JetCorrectorOnTheFly<Jet>::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0067 edm::Handle<reco::JetCorrector> corrector;
0068 iEvent.getByToken(mJetCorrector, corrector);
0069 Handle<JetCollection> jets;
0070 iEvent.getByToken(mJetToken, jets);
0071
0072 edm::Handle<reco::VertexCollection> recVtxs;
0073 iEvent.getByToken(mVertexToken, recVtxs);
0074 int NPV(0);
0075 for (unsigned int ind = 0; ind < recVtxs->size(); ind++) {
0076 if (!((*recVtxs)[ind].isFake())) {
0077 NPV++;
0078 }
0079 }
0080 typename JetCollection::const_iterator i_jet;
0081
0082 for (i_jet = jets->begin(); i_jet != jets->end(); i_jet++) {
0083 if (i_jet->pt() < mMinRawJetPt)
0084 continue;
0085
0086 double scale = corrector->correction(*i_jet);
0087 if (mDebug) {
0088 std::cout << "energy = " << i_jet->energy() << ", "
0089 << "eta = " << i_jet->eta() << ", "
0090 << "raw pt = " << i_jet->pt() << ", "
0091 << "NPV = " << NPV << ", "
0092 << "correction = " << scale << ", "
0093 << "cor pt = " << scale * i_jet->pt() << endl;
0094 }
0095 mRawPt->Fill(i_jet->pt());
0096 mCorPt->Fill(scale * i_jet->pt());
0097 }
0098 }
0099
0100 template <class Jet>
0101 void JetCorrectorOnTheFly<Jet>::beginJob() {
0102 mRawPt = fs->make<TH1F>("RawJetPt", "RawJetPt", 1000, 0, 1000);
0103 mCorPt = fs->make<TH1F>("CorJetPt", "CorJetPt", 1000, 0, 1000);
0104 }
0105
0106 template <class Jet>
0107 void JetCorrectorOnTheFly<Jet>::endJob() {}
0108
0109 typedef JetCorrectorOnTheFly<CaloJet> CaloJetCorrectorOnTheFly;
0110 DEFINE_FWK_MODULE(CaloJetCorrectorOnTheFly);
0111
0112 typedef JetCorrectorOnTheFly<PFJet> PFJetCorrectorOnTheFly;
0113 DEFINE_FWK_MODULE(PFJetCorrectorOnTheFly);
0114
0115 typedef JetCorrectorOnTheFly<JPTJet> JPTJetCorrectorOnTheFly;
0116 DEFINE_FWK_MODULE(JPTJetCorrectorOnTheFly);