CMS 3D CMS Logo

FFTJetPileupAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: FFTJetPileupAnalyzer
4 // Class: FFTJetPileupAnalyzer
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Thu Apr 21 15:52:11 CDT 2011
16 //
17 //
18 
19 #include <cassert>
20 #include <sstream>
21 #include <numeric>
22 
23 #include "TNtuple.h"
24 
25 // user include files
32 
34 #include <TH2D.h>
38 
40 
42 
44 
46 
47 #define init_param(type, varname) varname(ps.getParameter<type>(#varname))
48 
49 //
50 // class declaration
51 //
52 class FFTJetPileupAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
53 public:
54  explicit FFTJetPileupAnalyzer(const edm::ParameterSet&);
55  FFTJetPileupAnalyzer() = delete;
58  ~FFTJetPileupAnalyzer() override;
59 
60 private:
61  // The following method should take all necessary info from
62  // PileupSummaryInfo and fill out the ntuple
63  void analyzePileup(const std::vector<PileupSummaryInfo>& pInfo);
64 
65  void beginJob() override;
66  void analyze(const edm::Event&, const edm::EventSetup&) override;
67  void endJob() override;
68 
76 
85 
97 
98  double vertexNdofCut;
100 
101  std::vector<float> ntupleData;
102  TNtuple* nt;
103  int totalNpu;
104  int totalNPV;
105  unsigned long counter;
106 };
107 
108 //
109 // constructors and destructor
110 //
130  init_param(double, vertexNdofCut),
131  init_param(double, crazyEnergyCut),
132  nt(nullptr),
133  totalNpu(-1),
134  totalNPV(-1),
135  counter(0) {
136  usesResource(TFileService::kSharedResource);
137 
139  pileupToken = consumes<std::vector<PileupSummaryInfo> >(pileupLabel);
140 
141  if (collectHistos)
142  histoToken = consumes<TH2D>(histoLabel);
143 
144  if (collectSummaries)
145  summaryToken = consumes<reco::FFTJetPileupSummary>(summaryLabel);
146 
147  if (collectFastJetRho) {
148  fastJetRhoToken = consumes<double>(fastJetRhoLabel);
149  fastJetSigmaToken = consumes<double>(fastJetSigmaLabel);
150  }
151 
152  if (collectGrids)
153  gridToken = consumes<reco::DiscretizedEnergyFlow>(gridLabel);
154 
155  if (collectGridDensity)
156  etSumToken = consumes<std::pair<double, double> >(histoLabel);
157 
158  if (collectVertexInfo)
159  srcPVsToken = consumes<reco::VertexCollection>(srcPVs);
160 }
161 
163 
164 //
165 // member functions
166 //
167 void FFTJetPileupAnalyzer::analyzePileup(const std::vector<PileupSummaryInfo>& info) {
168  const unsigned nBx = info.size();
169  if (collectPileup)
170  ntupleData.push_back(static_cast<float>(nBx));
171 
172  double sumpt_Lo = 0.0, sumpt_Hi = 0.0;
173  totalNpu = 0;
174 
175  int npu_by_Bx[3] = {
176  0,
177  };
178  double sumpt_Lo_by_Bx[3] =
179  {
180  0.0,
181  },
182  sumpt_Hi_by_Bx[3] = {
183  0.0,
184  };
185 
186  if (verbosePileupInfo)
187  std::cout << "\n**** Pileup info begin" << std::endl;
188 
189  bool isCrazy = false;
190  for (unsigned ibx = 0; ibx < nBx; ++ibx) {
191  const PileupSummaryInfo& puInfo(info[ibx]);
192 
193  const int bx = puInfo.getBunchCrossing();
194  const int npu = puInfo.getPU_NumInteractions();
195  const std::vector<float>& lopt(puInfo.getPU_sumpT_lowpT());
196  const std::vector<float>& hipt(puInfo.getPU_sumpT_highpT());
197  const double losum = std::accumulate(lopt.begin(), lopt.end(), 0.0);
198  const double hisum = std::accumulate(hipt.begin(), hipt.end(), 0.0);
199 
200  if (losum >= crazyEnergyCut)
201  isCrazy = true;
202  if (hisum >= crazyEnergyCut)
203  isCrazy = true;
204 
205  totalNpu += npu;
206  sumpt_Lo += losum;
207  sumpt_Hi += hisum;
208 
209  const unsigned idx = bx < 0 ? 0U : (bx == 0 ? 1U : 2U);
210  npu_by_Bx[idx] += npu;
211  sumpt_Lo_by_Bx[idx] += losum;
212  sumpt_Hi_by_Bx[idx] += hisum;
213 
214  if (verbosePileupInfo)
215  std::cout << "ibx " << ibx << " bx " << bx << " npu " << npu << " losum " << losum << " hisum " << hisum
216  << std::endl;
217  }
218 
219  if (verbosePileupInfo)
220  std::cout << "**** Pileup info end\n" << std::endl;
221 
222  if (isCrazy) {
223  totalNpu = -1;
224  sumpt_Lo = 0.0;
225  sumpt_Hi = 0.0;
226  for (unsigned ibx = 0; ibx < 3; ++ibx) {
227  npu_by_Bx[ibx] = -1;
228  sumpt_Lo_by_Bx[ibx] = 0.0;
229  sumpt_Hi_by_Bx[ibx] = 0.0;
230  }
231  }
232 
233  if (collectPileup) {
234  ntupleData.push_back(totalNpu);
235  ntupleData.push_back(sumpt_Lo);
236  ntupleData.push_back(sumpt_Hi);
237  }
238 
239  if (collectOOTPileup)
240  for (unsigned ibx = 0; ibx < 3; ++ibx) {
241  ntupleData.push_back(npu_by_Bx[ibx]);
242  ntupleData.push_back(sumpt_Lo_by_Bx[ibx]);
243  ntupleData.push_back(sumpt_Hi_by_Bx[ibx]);
244  }
245 }
246 
247 // ------------ method called once each job just before starting event loop
249  // Come up with the list of variables
250  std::string vars = "cnt:run:event";
251  if (collectPileup)
252  vars += ":nbx:npu:sumptLowCut:sumptHiCut";
253  if (collectOOTPileup) {
254  vars += ":npu_negbx:sumptLowCut_negbx:sumptHiCut_negbx";
255  vars += ":npu_0bx:sumptLowCut_0bx:sumptHiCut_0bx";
256  vars += ":npu_posbx:sumptLowCut_posbx:sumptHiCut_posbx";
257  }
258  if (collectSummaries)
259  vars += ":estimate:pileup:uncert:uncertCode";
260  if (collectFastJetRho)
261  vars += ":fjrho:fjsigma";
262  if (collectGridDensity)
263  vars += ":gridEtDensity:gridEtDensityMixed";
264  if (collectVertexInfo)
265  vars += ":nPV";
266 
267  // Book the ntuple
269  nt = fs->make<TNtuple>(ntupleName.c_str(), ntupleTitle.c_str(), vars.c_str());
270  ntupleData.reserve(nt->GetNvar());
271 }
272 
273 // ------------ method called to for each event ------------
275  ntupleData.clear();
276  ntupleData.push_back(counter);
277  totalNpu = -1;
278  totalNPV = -1;
279 
280  edm::RunNumber_t const runnumber = iEvent.id().run();
281  edm::EventNumber_t const eventnumber = iEvent.id().event();
282  ntupleData.push_back(runnumber);
283  ntupleData.push_back(eventnumber);
284 
285  // Get pileup information from the pile-up information module
288  if (iEvent.getByToken(pileupToken, puInfo))
289  analyzePileup(*puInfo);
290  else {
291  if (collectPileup) {
292  ntupleData.push_back(-1);
293  ntupleData.push_back(-1);
294  ntupleData.push_back(0.f);
295  ntupleData.push_back(0.f);
296  }
297  if (collectOOTPileup)
298  for (unsigned ibx = 0; ibx < 3; ++ibx) {
299  ntupleData.push_back(-1);
300  ntupleData.push_back(0.f);
301  ntupleData.push_back(0.f);
302  }
303  }
304  }
305 
306  if (collectHistos) {
308  iEvent.getByToken(histoToken, input);
309 
311  TH2D* copy = new TH2D(*input);
312 
313  std::ostringstream os;
314  os << copy->GetName() << '_' << counter << '_' << totalNpu << '_' << runnumber << '_' << eventnumber;
315  const std::string& newname(os.str());
316  copy->SetNameTitle(newname.c_str(), newname.c_str());
317 
318  copy->SetDirectory(fs->getBareDirectory());
319  }
320 
321  if (collectSummaries) {
323  iEvent.getByToken(summaryToken, summary);
324 
325  ntupleData.push_back(summary->uncalibratedQuantile());
326  ntupleData.push_back(summary->pileupRho());
327  ntupleData.push_back(summary->pileupRhoUncertainty());
328  ntupleData.push_back(summary->uncertaintyCode());
329  }
330 
331  if (collectFastJetRho) {
332  edm::Handle<double> fjrho, fjsigma;
333  iEvent.getByToken(fastJetRhoToken, fjrho);
334  iEvent.getByToken(fastJetSigmaToken, fjsigma);
335 
336  ntupleData.push_back(*fjrho);
337  ntupleData.push_back(*fjsigma);
338  }
339 
340  if (collectGrids) {
342  iEvent.getByToken(gridToken, input);
343 
344  // Make sure the input grid is reasonable
345  const double* data = input->data();
346  assert(data);
347  assert(input->phiBin0Edge() == 0.0);
348  const unsigned nEta = input->nEtaBins();
349  const unsigned nPhi = input->nPhiBins();
350 
351  // Generate a name for the output histogram
352  std::ostringstream os;
353  os << "FFTJetGrid_" << counter << '_' << totalNpu << '_' << runnumber << '_' << eventnumber;
354  const std::string& newname(os.str());
355 
356  // Make a histogram and copy the grid data into it
358  TH2F* h =
359  fs->make<TH2F>(newname.c_str(), newname.c_str(), nEta, input->etaMin(), input->etaMax(), nPhi, 0.0, 2.0 * M_PI);
360  h->GetXaxis()->SetTitle("Eta");
361  h->GetYaxis()->SetTitle("Phi");
362  h->GetZaxis()->SetTitle("Transverse Energy");
363 
364  for (unsigned ieta = 0; ieta < nEta; ++ieta)
365  for (unsigned iphi = 0; iphi < nPhi; ++iphi)
366  h->SetBinContent(ieta + 1U, iphi + 1U, data[ieta * nPhi + iphi]);
367  }
368 
369  if (collectGridDensity) {
371  iEvent.getByToken(etSumToken, etSum);
372 
373  ntupleData.push_back(etSum->first);
374  ntupleData.push_back(etSum->second);
375  }
376 
377  if (collectVertexInfo) {
379  iEvent.getByToken(srcPVsToken, pvCollection);
380  totalNPV = 0;
381  if (!pvCollection->empty())
382  for (reco::VertexCollection::const_iterator pv = pvCollection->begin(); pv != pvCollection->end(); ++pv) {
383  const double ndof = pv->ndof();
384  if (!pv->isFake() && ndof > vertexNdofCut)
385  ++totalNPV;
386  }
387  ntupleData.push_back(totalNPV);
388  }
389 
390  assert(ntupleData.size() == static_cast<unsigned>(nt->GetNvar()));
391  nt->Fill(&ntupleData[0]);
392 
393  ++counter;
394 }
395 
396 // ------------ method called once each job just after ending the event loop
398 
399 //define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
static const TGPicture * info(bool iBackgroundIsBlack)
#define init_param(type, varname)
edm::EDGetTokenT< reco::VertexCollection > srcPVsToken
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > pileupToken
unsigned long long EventNumber_t
FFTJetPileupAnalyzer()=delete
assert(be >=bs)
void analyzePileup(const std::vector< PileupSummaryInfo > &pInfo)
static std::string const input
Definition: EdmProvDump.cc:50
FFTJetPileupAnalyzer & operator=(const FFTJetPileupAnalyzer &)=delete
edm::EDGetTokenT< TH2D > histoToken
const std::vector< float > & getPU_sumpT_lowpT() const
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< reco::FFTJetPileupSummary > summaryToken
def pv(vc)
Definition: MetAnalyzer.py:7
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< double > fastJetSigmaToken
int nt
Definition: AMPTWrapper.h:42
#define M_PI
const int getBunchCrossing() const
HLT enums.
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
edm::EDGetTokenT< double > fastJetRhoToken
const std::vector< float > & getPU_sumpT_highpT() const
unsigned int RunNumber_t
vars
Definition: DeepTauId.cc:30
edm::EDGetTokenT< reco::DiscretizedEnergyFlow > gridToken
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
std::vector< float > ntupleData
edm::EDGetTokenT< std::pair< double, double > > etSumToken
void analyze(const edm::Event &, const edm::EventSetup &) override
const int getPU_NumInteractions() const