CMS 3D CMS Logo

PFJetAnalyzerDQM.cc
Go to the documentation of this file.
13 
14 #include <algorithm>
15 #include <numeric>
16 #include <regex>
17 #include <sstream>
18 #include <vector>
19 #include <memory>
20 
22 public:
24  void analyze(const edm::Event&, const edm::EventSetup&) override;
25 
26 protected:
27  //Book histograms
28  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
29 
30 private:
31  class Plot1DInBin {
32  public:
34  const uint32_t nbins;
35  const float min, max;
38 
39  Plot1DInBin(const std::string _name,
40  const std::string _title,
41  const uint32_t _nbins,
42  const float _min,
43  const float _max,
44  float _ptbin_low,
45  float _ptbin_high,
46  float _etabin_low,
47  float _etabin_high)
48  : name(_name),
49  title(_title),
50  nbins(_nbins),
51  min(_min),
52  max(_max),
53  ptbin_low(_ptbin_low),
54  ptbin_high(_ptbin_high),
55  etabin_low(_etabin_low),
56  etabin_high(_etabin_high) {}
57 
58  void book(DQMStore::IBooker& booker) { plot_ = booker.book1D(name, title, nbins, min, max); }
59 
60  void fill(float value) {
61  assert(plot_ != nullptr);
62  plot_->Fill(value);
63  }
64 
65  //Check if a jet with a value v would be in the bin that applies to this plot
66  bool isInBin(float v, float low, float high) { return v >= low && v < high; }
67 
68  bool isInPtBin(float pt) { return isInBin(pt, ptbin_low, ptbin_high); }
69 
70  bool isInEtaBin(float eta) { return isInBin(eta, etabin_low, etabin_high); }
71 
72  bool isInPtEtaBin(float pt, float eta) { return isInPtBin(pt) && isInEtaBin(eta); }
73  };
74 
76  public:
78  std::unique_ptr<TH1F> base_hist;
81 
83  const std::string _title,
84  std::unique_ptr<TH1F> _base_hist,
85  float _ptbin_low,
86  float _ptbin_high,
87  float _etabin_low,
88  float _etabin_high)
89  : name(_name),
90  title(_title),
91  base_hist(std::move(_base_hist)),
92  ptbin_low(_ptbin_low),
93  ptbin_high(_ptbin_high),
94  etabin_low(_etabin_low),
95  etabin_high(_etabin_high) {}
96 
97  void book(DQMStore::IBooker& booker) { plot_ = booker.book1D(name.c_str(), base_hist.get()); }
98 
99  void fill(float value) {
100  assert(plot_ != nullptr);
101  plot_->Fill(value);
102  }
103 
104  //Check if a jet with a value v would be in the bin that applies to this plot
105  bool isInBin(float v, float low, float high) { return v >= low && v < high; }
106 
107  bool isInPtBin(float pt) { return isInBin(pt, ptbin_low, ptbin_high); }
108 
109  bool isInEtaBin(float eta) { return isInBin(eta, etabin_low, etabin_high); }
110 
111  bool isInPtEtaBin(float pt, float eta) { return isInPtBin(pt) && isInEtaBin(eta); }
112  };
113 
114  std::vector<Plot1DInBin> jetResponsePlots;
115  std::vector<Plot1DInBin> jetResponsePlots_noJEC;
116  std::vector<Plot1DInBinVariable> genJetPlots;
117 
118  // Is this data or MC?
119  bool isMC;
120 
121  float jetDeltaR;
122 
123  bool genJetsOn;
124 
126 
132 
134  void prepareJetResponsePlots(const std::vector<edm::ParameterSet>& genjet_plots_pset);
135  void prepareGenJetPlots(const std::vector<edm::ParameterSet>& genjet_plots_pset);
136 };
137 
138 void PFJetAnalyzerDQM::prepareJetResponsePlots(const std::vector<edm::ParameterSet>& response_plots) {
139  for (auto& pset : response_plots) {
140  //Low and high edges of the pt and eta bins for jets to pass to be filled into this histogram
141  const auto ptbin_low = pset.getParameter<double>("ptBinLow");
142  const auto ptbin_high = pset.getParameter<double>("ptBinHigh");
143  const auto etabin_low = pset.getParameter<double>("etaBinLow");
144  const auto etabin_high = pset.getParameter<double>("etaBinHigh");
145 
146  const auto response_nbins = pset.getParameter<uint32_t>("responseNbins");
147  const auto response_low = pset.getParameter<double>("responseLow");
148  const auto response_high = pset.getParameter<double>("responseHigh");
149 
150  const auto name = pset.getParameter<std::string>("name");
151  const auto title = pset.getParameter<std::string>("title");
152 
153  // for title of raw jet response histograms
154  auto rawTitle = title;
155  rawTitle = rawTitle.replace(rawTitle.begin(), rawTitle.begin(), "Raw ");
156 
157  jetResponsePlots.push_back(Plot1DInBin(
158  name, title, response_nbins, response_low, response_high, ptbin_low, ptbin_high, etabin_low, etabin_high));
159 
161  name, rawTitle, response_nbins, response_low, response_high, ptbin_low, ptbin_high, etabin_low, etabin_high));
162  }
163  if (jetResponsePlots.size() > 200) {
164  throw std::runtime_error("Requested too many jet response plots, aborting as this seems unusual.");
165  }
166 }
167 
168 void PFJetAnalyzerDQM::prepareGenJetPlots(const std::vector<edm::ParameterSet>& genjet_plots_pset) {
169  for (auto& pset : genjet_plots_pset) {
170  const auto name = pset.getParameter<std::string>("name");
171  const auto title = pset.getParameter<std::string>("title");
172 
173  //Low and high edges of the eta bins for jets to pass to be filled into this histogram
174  const auto ptbins_d = pset.getParameter<std::vector<double>>("ptBins");
175  std::vector<float> ptbins(ptbins_d.begin(), ptbins_d.end());
176 
177  const auto etabin_low = pset.getParameter<double>("etaBinLow");
178  const auto etabin_high = pset.getParameter<double>("etaBinHigh");
179 
181  name,
182  title,
183  std::make_unique<TH1F>(name.c_str(), title.c_str(), static_cast<int>(ptbins.size()) - 1, ptbins.data()),
184  0.0,
185  0.0,
186  etabin_low,
187  etabin_high));
188  }
189 }
190 
192  recoJetsLabel = iConfig.getParameter<edm::InputTag>("recoJetCollection");
193  genJetsLabel = iConfig.getParameter<edm::InputTag>("genJetCollection");
194 
195  //label for making new folder
197 
198  //DeltaR for reco to gen jet matching
199  jetDeltaR = iConfig.getParameter<double>("jetDeltaR");
200 
201  //for turn genJet on/off
202  genJetsOn = iConfig.getParameter<bool>("genJetsOn");
203 
204  //Create all jet response plots in bins of genjet pt and eta
205  const auto& response_plots = iConfig.getParameter<std::vector<edm::ParameterSet>>("responsePlots");
206  prepareJetResponsePlots(response_plots);
207 
208  const auto& genjet_plots = iConfig.getParameter<std::vector<edm::ParameterSet>>("genJetPlots");
209  prepareGenJetPlots(genjet_plots);
210 
211  recoJetsToken = consumes<edm::View<pat::Jet>>(recoJetsLabel);
212  genJetsToken = consumes<edm::View<reco::Jet>>(genJetsLabel);
213 }
214 
216  //match gen jets to reco jets, require minimum jetDeltaR, choose closest, do not try to match charge
217  std::vector<int> matchIndices;
218  PFB::match(genJetCollection, recoJetCollection, matchIndices, false, jetDeltaR);
219 
220  for (unsigned int i = 0; i < genJetCollection.size(); i++) {
221  const auto& genJet = genJetCollection.at(i);
222  const auto pt_gen = genJet.pt();
223  const auto eta_gen = abs(genJet.eta());
224  const int iMatch = matchIndices[i];
225 
226  //Fill genjet pt if genJetOn
227  if (genJetsOn) {
228  for (auto& plot : genJetPlots) {
229  if (plot.isInEtaBin(eta_gen)) {
230  plot.fill(pt_gen);
231  }
232  }
233  }
234 
235  //If gen jet had a matched reco jet
236  if (iMatch != -1) {
237  const auto& recoJet = recoJetCollection[iMatch];
238  auto pt_reco = recoJet.pt();
239 
240  const auto response = pt_reco / pt_gen;
241  const auto response_raw = pt_reco * recoJet.jecFactor("Uncorrected") / pt_gen;
242 
243  //Loop linearly through all plots and check if they match the pt and eta bin
244  //this is not algorithmically optimal but we don't expect to more than a few hundred plots
245  //If this turns out to be a problem, can easily make a 2D-map for indices
246  for (auto& plot : jetResponsePlots) {
247  if (plot.isInPtEtaBin(pt_gen, eta_gen)) {
248  plot.fill(response);
249  }
250  }
251  // this loop should be for NoJEC plots
252  for (auto& plot : jetResponsePlots_noJEC) {
253  if (plot.isInPtEtaBin(pt_gen, eta_gen)) {
254  plot.fill(response_raw);
255  }
256  }
257  }
258  }
259 }
260 
262  booker.setCurrentFolder("ParticleFlow/JetResponse/" + jetCollectionName + "/JEC/");
263  for (auto& plot : jetResponsePlots) {
264  plot.book(booker);
265  }
266  //Book plots for noJEC
267  booker.setCurrentFolder("ParticleFlow/JetResponse/" + jetCollectionName + "/noJEC/");
268  for (auto& plot : jetResponsePlots_noJEC) {
269  plot.book(booker);
270  }
271 
272  //Book plots for gen-jet pt spectra
273  if (genJetsOn) {
274  booker.setCurrentFolder("ParticleFlow/GenJets/");
275  for (auto& plot : genJetPlots) {
276  plot.book(booker);
277  }
278  }
279 }
280 
282  edm::Handle<edm::View<pat::Jet>> recoJetCollectionHandle;
283  iEvent.getByToken(recoJetsToken, recoJetCollectionHandle);
284  auto recoJetCollection = *recoJetCollectionHandle;
285 
286  isMC = !iEvent.isRealData();
287 
288  if (isMC) {
289  edm::Handle<edm::View<reco::Jet>> genJetCollectionHandle;
290  iEvent.getByToken(genJetsToken, genJetCollectionHandle);
291  auto genJetCollection = *genJetCollectionHandle;
292 
294  }
295 }
296 
edm::EDGetTokenT< edm::View< pat::Jet > > recoJetsToken
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Plot1DInBin(const std::string _name, const std::string _title, const uint32_t _nbins, const float _min, const float _max, float _ptbin_low, float _ptbin_high, float _etabin_low, float _etabin_high)
void book(DQMStore::IBooker &booker)
void prepareJetResponsePlots(const std::vector< edm::ParameterSet > &genjet_plots_pset)
void match(const C &candCollection, const M &matchedCandCollection, std::vector< int > &matchIndices, bool matchCharge=false, float dRMax=-1)
Definition: Matchers.h:17
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< Plot1DInBin > jetResponsePlots
bool isInBin(float v, float low, float high)
void fillJetResponse(edm::View< pat::Jet > &recoJetCollection, edm::View< reco::Jet > &genJetCollection)
std::string const & label() const
Definition: InputTag.h:36
assert(be >=bs)
edm::EDGetTokenT< edm::View< reco::Jet > > genJetsToken
PFJetAnalyzerDQM(const edm::ParameterSet &)
bool isInBin(float v, float low, float high)
bool isInPtEtaBin(float pt, float eta)
edm::InputTag genJetsLabel
void Fill(long long x)
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
int iEvent
Definition: GenABIO.cc:224
void book(DQMStore::IBooker &booker)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::string jetCollectionName
Definition: value.py:1
edm::InputTag recoJetsLabel
edm::EDGetTokenT< reco::CandViewMatchMap > srcRefToJetMap
std::vector< Plot1DInBinVariable > genJetPlots
void analyze(const edm::Event &, const edm::EventSetup &) override
void prepareGenJetPlots(const std::vector< edm::ParameterSet > &genjet_plots_pset)
Plot1DInBinVariable(const std::string _name, const std::string _title, std::unique_ptr< TH1F > _base_hist, float _ptbin_low, float _ptbin_high, float _etabin_low, float _etabin_high)
bool isInPtEtaBin(float pt, float eta)
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
std::vector< Plot1DInBin > jetResponsePlots_noJEC
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45