CMS 3D CMS Logo

PatJetAnalyzer.cc
Go to the documentation of this file.
1 #include <map>
2 #include <string>
3 #include <iomanip>
4 #include <sstream>
5 #include <iostream>
6 
7 #include "TH1.h"
8 
15 
18 
21 static const unsigned int MAXBIN = 8;
25 static const float BINS[] = {30., 40., 50., 60., 70., 80., 100., 125., 150.};
26 
43 public:
45  explicit PatJetAnalyzer(const edm::ParameterSet& cfg);
47  ~PatJetAnalyzer() override{};
48 
49 private:
51  void analyze(const edm::Event& event, const edm::EventSetup& setup) override;
52 
54  bool booked(const std::string histName) const { return hists_.find(histName) != hists_.end(); };
56  void fill(const std::string histName, double value) const {
57  if (booked(histName))
58  hists_.find(histName)->second->Fill(value);
59  };
60  // print jet pt for each level of energy corrections
61  void print(edm::View<pat::Jet>::const_iterator& jet, unsigned int idx);
62 
63 private:
69  std::map<std::string, TH1F*> hists_;
70 };
71 
73  : corrLevel_(cfg.getParameter<std::string>("corrLevel")),
74  jetsToken_(consumes<edm::View<pat::Jet> >(cfg.getParameter<edm::InputTag>("src"))) {
75  // register TFileService
77 
78  // jet multiplicity
79  hists_["mult"] = fs->make<TH1F>("mult", "N_{Jet}", 15, 0., 15.);
80  // jet pt (for all jets)
81  hists_["pt"] = fs->make<TH1F>("pt", "p_{T}(Jet) [GeV]", 60, 0., 300.);
82  // jet eta (for all jets)
83  hists_["eta"] = fs->make<TH1F>("eta", "#eta (Jet)", 60, -3., 3.);
84  // jet phi (for all jets)
85  hists_["phi"] = fs->make<TH1F>("phi", "#phi (Jet)", 60, 3.2, 3.2);
86  // dijet mass (if available)
87  hists_["mass"] = fs->make<TH1F>("mass", "M_{jj} [GeV]", 50, 0., 500.);
88  // basic histograms for jet energy response
89  for (unsigned int idx = 0; idx < MAXBIN; ++idx) {
90  char buffer[10];
91  sprintf(buffer, "jes_%i", idx);
92  char title[50];
93  sprintf(title, "p_{T}^{rec}/p_{T}^{gen} [%i GeV - %i GeV]", (int)BINS[idx], (int)BINS[idx + 1]);
94  hists_[buffer] = fs->make<TH1F>(buffer, title, 100, 0., 2.);
95  }
96 }
97 
99  // recieve jet collection label
101  event.getByToken(jetsToken_, jets);
102 
103  // loop jets
104  for (edm::View<pat::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); ++jet) {
105  // print jec factors
106  // print(jet, jet-jets->begin());
107 
108  // fill basic kinematics
109  fill("pt", jet->correctedJet(corrLevel_).pt());
110  fill("eta", jet->eta());
111  fill("phi", jet->phi());
112  // basic plots for jet responds plot as a function of pt
113  if (jet->genJet()) {
114  double resp = jet->correctedJet(corrLevel_).pt() / jet->genJet()->pt();
115  for (unsigned int idx = 0; idx < MAXBIN; ++idx) {
116  if (BINS[idx] <= jet->genJet()->pt() && jet->genJet()->pt() < BINS[idx + 1]) {
117  char buffer[10];
118  sprintf(buffer, "jes_%i", idx);
119  fill(buffer, resp);
120  }
121  }
122  }
123  }
124  // jet multiplicity
125  fill("mult", jets->size());
126  // invariant dijet mass for first two leading jets
127  if (jets->size() > 1) {
128  fill("mass", ((*jets)[0].p4() + (*jets)[1].p4()).mass());
129  }
130 }
131 
133  //edm::LogInfo log("JEC");
134  std::cout << "[" << idx << "] :: eta=" << std::setw(10) << jet->eta() << " phi=" << std::setw(10) << jet->phi()
135  << " size: " << jet->availableJECLevels().size() << std::endl;
136  for (unsigned int idx = 0; idx < jet->availableJECLevels().size(); ++idx) {
137  pat::Jet correctedJet;
138  if (jet->availableJECLevels()[idx].find("L5Flavor") != std::string::npos ||
139  jet->availableJECLevels()[idx].find("L7Parton") != std::string::npos) {
140  correctedJet = jet->correctedJet(idx, pat::JetCorrFactors::UDS);
141  } else {
142  correctedJet = jet->correctedJet(idx, pat::JetCorrFactors::NONE);
143  }
144  std::cout << std::setw(10) << correctedJet.currentJECLevel() << " pt=" << std::setw(10) << correctedJet.pt()
145  << std::endl;
146  }
147 }
148 
std::string corrLevel_
correction level for pat jet
PatJetAnalyzer(const edm::ParameterSet &cfg)
default contructor
size_type size() const
double pt() const final
transverse momentum
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
void fill(const std::string histName, double value) const
fill histogram if it had been booked before
void analyze(const edm::Event &event, const edm::EventSetup &setup) override
everything that needs to be done during the even loop
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Definition: HeavyIon.h:7
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< edm::View< pat::Jet > > jetsToken_
pat jets
Definition: Jet.py:1
double p4[4]
Definition: TauolaWrapper.h:92
static const float BINS[]
Definition: value.py:1
std::string currentJECLevel() const
return the name of the current step of jet energy corrections
Definition: Jet.h:141
Analysis-level calorimeter jet class.
Definition: Jet.h:77
Module to analyze pat::Jets in the context of a more complex exercise.
bool booked(const std::string histName) const
check if histogram was booked
static const unsigned int MAXBIN
std::map< std::string, TH1F * > hists_
management of 1d histograms
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
~PatJetAnalyzer() override
default destructor
Definition: event.py:1
void print(edm::View< pat::Jet >::const_iterator &jet, unsigned int idx)