File indexing completed on 2023-03-17 11:10:47
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015 #include <memory>
0016
0017
0018 #include "FWCore/Framework/interface/Frameworkfwd.h"
0019 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0020
0021 #include "FWCore/Framework/interface/Event.h"
0022 #include "FWCore/Framework/interface/MakerMacros.h"
0023
0024 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0025
0026
0027 #include "FWCore/Framework/interface/ESHandle.h"
0028 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0029 #include "DataFormats/JetReco/interface/CaloJet.h"
0030
0031 #include "JetMETCorrections/JetCorrector/interface/JetCorrector.h"
0032 #include "DataFormats/TrackReco/interface/Track.h"
0033 #include "DataFormats/BTauReco/interface/IsolatedTauTagInfo.h"
0034 #include "DataFormats/JetReco/interface/Jet.h"
0035 #include <DataFormats/Common/interface/RefVector.h>
0036 #include <DataFormats/Common/interface/RefVectorIterator.h>
0037 #include <DataFormats/Common/interface/Ref.h>
0038
0039
0040
0041 #include <iostream>
0042 #include <string>
0043
0044 #include <TFile.h>
0045 #include <TH1D.h>
0046 #include <TMath.h>
0047 #include <TTree.h>
0048 #include <TLorentzVector.h>
0049 #include <TVector3.h>
0050
0051
0052
0053
0054 class TauJetCorrectorExample : public edm::one::EDAnalyzer<> {
0055 public:
0056 explicit TauJetCorrectorExample(const edm::ParameterSet&);
0057 ~TauJetCorrectorExample() override;
0058
0059 private:
0060 void beginJob() override;
0061 void analyze(const edm::Event&, const edm::EventSetup&) override;
0062 void endJob() override;
0063
0064
0065 std::string jetname;
0066 edm::EDGetTokenT<reco::IsolatedTauTagInfoCollection> tautoken;
0067 edm::EDGetTokenT<reco::JetCorrector> tauCorrectortoken;
0068
0069 int nEvt;
0070 int njets;
0071 };
0072
0073
0074
0075
0076
0077
0078
0079
0080
0081
0082
0083
0084
0085 using namespace reco;
0086
0087 TauJetCorrectorExample::TauJetCorrectorExample(const edm::ParameterSet& iConfig)
0088 : jetname(iConfig.getUntrackedParameter<std::string>("JetHandle", "iterativeCone5CaloJets")),
0089 tautoken(consumes<reco::IsolatedTauTagInfoCollection>(
0090 edm::InputTag(iConfig.getUntrackedParameter<std::string>("TauHandle", "coneIsolation")))),
0091 tauCorrectortoken(consumes<reco::JetCorrector>(
0092 iConfig.getUntrackedParameter<std::string>("tauCorrHandle", "TauJetCorrectorIcone5"))),
0093 nEvt(0),
0094 njets(0) {
0095
0096 }
0097
0098 TauJetCorrectorExample::~TauJetCorrectorExample() {
0099
0100
0101 }
0102
0103
0104
0105
0106
0107
0108 void TauJetCorrectorExample::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0109 using namespace edm;
0110
0111 Handle<reco::JetCorrector> taucorrector;
0112 iEvent.getByToken(tauCorrectortoken, taucorrector);
0113
0114
0115
0116 ++nEvt;
0117 if ((nEvt % 10 == 0 && nEvt <= 100) || (nEvt % 100 == 0 && nEvt > 100))
0118 std::cout << "reading event " << nEvt << std::endl;
0119
0120
0121 Handle<reco::IsolatedTauTagInfoCollection> tauTagInfoHandle;
0122 iEvent.getByToken(tautoken, tauTagInfoHandle);
0123 reco::IsolatedTauTagInfoCollection::const_iterator tau = tauTagInfoHandle->begin();
0124
0125
0126 njets = 0;
0127
0128 std::cout << "starting tau loop" << std::endl;
0129 for (tau = tauTagInfoHandle->begin(); tau != tauTagInfoHandle->end() && njets < 10; ++tau) {
0130
0131
0132 double pt = tau->jet().get()->et();
0133
0134
0135 double scale = taucorrector->correction(tau->jet().get()->p4());
0136 double ptcorr = tau->jet().get()->et() * scale;
0137
0138 std::cout << "Tau jet: Original Et = " << pt << " Corrected Et = " << ptcorr << std::endl;
0139
0140 ++njets;
0141 }
0142 }
0143
0144
0145 void
0146
0147 TauJetCorrectorExample::beginJob() {}
0148
0149
0150 void TauJetCorrectorExample::endJob() {}
0151
0152
0153 DEFINE_FWK_MODULE(TauJetCorrectorExample);