test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
38 
40 
42 
44 
46 
47 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
48 
49 //
50 // class declaration
51 //
53 {
54 public:
55  explicit FFTJetPileupAnalyzer(const edm::ParameterSet&);
57 
58 private:
62 
63  // The following method should take all necessary info from
64  // PileupSummaryInfo and fill out the ntuple
65  void analyzePileup(const std::vector<PileupSummaryInfo>& pInfo);
66 
67  virtual void beginJob() override ;
68  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
69  virtual void endJob() override ;
70 
78 
87 
99 
102 
103  std::vector<float> ntupleData;
104  TNtuple* nt;
105  int totalNpu;
106  int totalNPV;
107  unsigned long counter;
108 };
109 
110 //
111 // constructors and destructor
112 //
114  : init_param(edm::InputTag, histoLabel),
115  init_param(edm::InputTag, summaryLabel),
116  init_param(edm::InputTag, fastJetRhoLabel),
117  init_param(edm::InputTag, fastJetSigmaLabel),
118  init_param(edm::InputTag, gridLabel),
119  init_param(edm::InputTag, srcPVs),
120  init_param(std::string, pileupLabel),
121  init_param(std::string, ntupleName),
122  init_param(std::string, ntupleTitle),
123  init_param(bool, collectHistos),
124  init_param(bool, collectSummaries),
125  init_param(bool, collectFastJetRho),
126  init_param(bool, collectPileup),
127  init_param(bool, collectOOTPileup),
128  init_param(bool, collectGrids),
129  init_param(bool, collectGridDensity),
130  init_param(bool, collectVertexInfo),
131  init_param(bool, verbosePileupInfo),
132  init_param(double, vertexNdofCut),
133  init_param(double, crazyEnergyCut),
134  nt(0),
135  totalNpu(-1),
136  totalNPV(-1),
137  counter(0)
138 {
140  pileupToken = consumes<std::vector<PileupSummaryInfo> >(pileupLabel);
141 
142  if (collectHistos)
143  histoToken = consumes<TH2D>(histoLabel);
144 
145  if (collectSummaries)
146  summaryToken = consumes<reco::FFTJetPileupSummary>(summaryLabel);
147 
148  if (collectFastJetRho)
149  {
150  fastJetRhoToken = consumes<double>(fastJetRhoLabel);
151  fastJetSigmaToken = consumes<double>(fastJetSigmaLabel);
152  }
153 
154  if (collectGrids)
155  gridToken = consumes<reco::DiscretizedEnergyFlow>(gridLabel);
156 
157  if (collectGridDensity)
158  etSumToken = consumes<std::pair<double,double> >(histoLabel);
159 
160  if (collectVertexInfo)
161  srcPVsToken = consumes<reco::VertexCollection>(srcPVs);
162 }
163 
164 
166 {
167 }
168 
169 
170 //
171 // member functions
172 //
174  const std::vector<PileupSummaryInfo>& info)
175 {
176  const unsigned nBx = info.size();
177  if (collectPileup)
178  ntupleData.push_back(static_cast<float>(nBx));
179 
180  double sumpt_Lo = 0.0, sumpt_Hi = 0.0;
181  totalNpu = 0;
182 
183  int npu_by_Bx[3] = {0,};
184  double sumpt_Lo_by_Bx[3] = {0.0,}, sumpt_Hi_by_Bx[3] = {0.0,};
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  {
192  const PileupSummaryInfo& puInfo(info[ibx]);
193 
194  const int bx = puInfo.getBunchCrossing();
195  const int npu = puInfo.getPU_NumInteractions();
196  const std::vector<float>& lopt(puInfo.getPU_sumpT_lowpT());
197  const std::vector<float>& hipt(puInfo.getPU_sumpT_highpT());
198  const double losum = std::accumulate(lopt.begin(), lopt.end(), 0.0);
199  const double hisum = std::accumulate(hipt.begin(), hipt.end(), 0.0);
200 
201  if (losum >= crazyEnergyCut)
202  isCrazy = true;
203  if (hisum >= crazyEnergyCut)
204  isCrazy = true;
205 
206  totalNpu += npu;
207  sumpt_Lo += losum;
208  sumpt_Hi += hisum;
209 
210  const unsigned idx = bx < 0 ? 0U : (bx == 0 ? 1U : 2U);
211  npu_by_Bx[idx] += npu;
212  sumpt_Lo_by_Bx[idx] += losum;
213  sumpt_Hi_by_Bx[idx] += hisum;
214 
215  if (verbosePileupInfo)
216  std::cout << "ibx " << ibx << " bx " << bx
217  << " npu " << npu << " losum " << losum
218  << " hisum " << hisum
219  << std::endl;
220  }
221 
222  if (verbosePileupInfo)
223  std::cout << "**** Pileup info end\n" << std::endl;
224 
225  if (isCrazy)
226  {
227  totalNpu = -1;
228  sumpt_Lo = 0.0;
229  sumpt_Hi = 0.0;
230  for (unsigned ibx = 0; ibx < 3; ++ibx)
231  {
232  npu_by_Bx[ibx] = -1;
233  sumpt_Lo_by_Bx[ibx] = 0.0;
234  sumpt_Hi_by_Bx[ibx] = 0.0;
235  }
236  }
237 
238  if (collectPileup)
239  {
240  ntupleData.push_back(totalNpu);
241  ntupleData.push_back(sumpt_Lo);
242  ntupleData.push_back(sumpt_Hi);
243  }
244 
245  if (collectOOTPileup)
246  for (unsigned ibx = 0; ibx < 3; ++ibx)
247  {
248  ntupleData.push_back(npu_by_Bx[ibx]);
249  ntupleData.push_back(sumpt_Lo_by_Bx[ibx]);
250  ntupleData.push_back(sumpt_Hi_by_Bx[ibx]);
251  }
252 }
253 
254 
255 // ------------ method called once each job just before starting event loop
257 {
258  // Come up with the list of variables
259  std::string vars = "cnt:run:event";
260  if (collectPileup)
261  vars += ":nbx:npu:sumptLowCut:sumptHiCut";
262  if (collectOOTPileup)
263  {
264  vars += ":npu_negbx:sumptLowCut_negbx:sumptHiCut_negbx";
265  vars += ":npu_0bx:sumptLowCut_0bx:sumptHiCut_0bx";
266  vars += ":npu_posbx:sumptLowCut_posbx:sumptHiCut_posbx";
267  }
268  if (collectSummaries)
269  vars += ":estimate:pileup:uncert:uncertCode";
270  if (collectFastJetRho)
271  vars += ":fjrho:fjsigma";
272  if (collectGridDensity)
273  vars += ":gridEtDensity:gridEtDensityMixed";
274  if (collectVertexInfo)
275  vars += ":nPV";
276 
277  // Book the ntuple
279  nt = fs->make<TNtuple>(ntupleName.c_str(), ntupleTitle.c_str(),
280  vars.c_str());
281  ntupleData.reserve(nt->GetNvar());
282 }
283 
284 
285 // ------------ method called to for each event ------------
287  const edm::EventSetup& iSetup)
288 {
289  ntupleData.clear();
290  ntupleData.push_back(counter);
291  totalNpu = -1;
292  totalNPV = -1;
293 
294  edm::RunNumber_t const runnumber = iEvent.id().run();
295  edm::EventNumber_t const eventnumber = iEvent.id().event();
296  ntupleData.push_back(runnumber);
297  ntupleData.push_back(eventnumber);
298 
299  // Get pileup information from the pile-up information module
301  {
303  if (iEvent.getByToken(pileupToken, puInfo))
304  analyzePileup(*puInfo);
305  else
306  {
307  if (collectPileup)
308  {
309  ntupleData.push_back(-1);
310  ntupleData.push_back(-1);
311  ntupleData.push_back(0.f);
312  ntupleData.push_back(0.f);
313  }
314  if (collectOOTPileup)
315  for (unsigned ibx = 0; ibx < 3; ++ibx)
316  {
317  ntupleData.push_back(-1);
318  ntupleData.push_back(0.f);
319  ntupleData.push_back(0.f);
320  }
321  }
322  }
323 
324  if (collectHistos)
325  {
327  iEvent.getByToken(histoToken, input);
328 
330  TH2D* copy = new TH2D(*input);
331 
332  std::ostringstream os;
333  os << copy->GetName() << '_' << counter << '_'
334  << totalNpu << '_' << runnumber << '_' << eventnumber;
335  const std::string& newname(os.str());
336  copy->SetNameTitle(newname.c_str(), newname.c_str());
337 
338  copy->SetDirectory(fs->getBareDirectory());
339  }
340 
341  if (collectSummaries)
342  {
344  iEvent.getByToken(summaryToken, summary);
345 
346  ntupleData.push_back(summary->uncalibratedQuantile());
347  ntupleData.push_back(summary->pileupRho());
348  ntupleData.push_back(summary->pileupRhoUncertainty());
349  ntupleData.push_back(summary->uncertaintyCode());
350  }
351 
352  if (collectFastJetRho)
353  {
354  edm::Handle<double> fjrho, fjsigma;
355  iEvent.getByToken(fastJetRhoToken, fjrho);
356  iEvent.getByToken(fastJetSigmaToken, fjsigma);
357 
358  ntupleData.push_back(*fjrho);
359  ntupleData.push_back(*fjsigma);
360  }
361 
362  if (collectGrids)
363  {
365  iEvent.getByToken(gridToken, input);
366 
367  // Make sure the input grid is reasonable
368  const double* data = input->data();
369  assert(data);
370  assert(input->phiBin0Edge() == 0.0);
371  const unsigned nEta = input->nEtaBins();
372  const unsigned nPhi = input->nPhiBins();
373 
374  // Generate a name for the output histogram
375  std::ostringstream os;
376  os << "FFTJetGrid_" << counter << '_'
377  << totalNpu << '_' << runnumber << '_' << eventnumber;
378  const std::string& newname(os.str());
379 
380  // Make a histogram and copy the grid data into it
382  TH2F* h = fs->make<TH2F>(newname.c_str(), newname.c_str(),
383  nEta, input->etaMin(), input->etaMax(),
384  nPhi, 0.0, 2.0*M_PI);
385  h->GetXaxis()->SetTitle("Eta");
386  h->GetYaxis()->SetTitle("Phi");
387  h->GetZaxis()->SetTitle("Transverse Energy");
388 
389  for (unsigned ieta=0; ieta<nEta; ++ieta)
390  for (unsigned iphi=0; iphi<nPhi; ++iphi)
391  h->SetBinContent(ieta+1U, iphi+1U, data[ieta*nPhi + iphi]);
392  }
393 
394  if (collectGridDensity)
395  {
397  iEvent.getByToken(etSumToken, etSum);
398 
399  ntupleData.push_back(etSum->first);
400  ntupleData.push_back(etSum->second);
401  }
402 
403  if (collectVertexInfo)
404  {
406  iEvent.getByToken(srcPVsToken, pvCollection);
407  totalNPV = 0;
408  if (!pvCollection->empty())
409  for (reco::VertexCollection::const_iterator pv = pvCollection->begin();
410  pv != pvCollection->end(); ++pv)
411  {
412  const double ndof = pv->ndof();
413  if (!pv->isFake() && ndof > vertexNdofCut)
414  ++totalNPV;
415  }
416  ntupleData.push_back(totalNPV);
417  }
418 
419  assert(ntupleData.size() == static_cast<unsigned>(nt->GetNvar()));
420  nt->Fill(&ntupleData[0]);
421 
422  ++counter;
423 }
424 
425 
426 // ------------ method called once each job just after ending the event loop
428 {
429 }
430 
431 
432 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:39
EventNumber_t event() const
Definition: EventID.h:41
static const TGPicture * info(bool iBackgroundIsBlack)
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
#define init_param(type, varname)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void endJob() override
edm::EDGetTokenT< reco::VertexCollection > srcPVsToken
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > pileupToken
virtual void beginJob() override
assert(m_qm.get())
unsigned long long EventNumber_t
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
const int getBunchCrossing() const
edm::EDGetTokenT< std::pair< double, double > > etSumToken
void analyzePileup(const std::vector< PileupSummaryInfo > &pInfo)
static std::string const input
Definition: EdmProvDump.cc:44
edm::EDGetTokenT< TH2D > histoToken
int iEvent
Definition: GenABIO.cc:230
const std::vector< float > & getPU_sumpT_highpT() const
TDirectory * getBareDirectory(const std::string &subdir="") const
Definition: TFileService.h:52
edm::EDGetTokenT< reco::FFTJetPileupSummary > summaryToken
double f[11][100]
edm::EDGetTokenT< double > fastJetSigmaToken
int nt
Definition: AMPTWrapper.h:32
#define M_PI
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
const int getPU_NumInteractions() const
FFTJetPileupAnalyzer & operator=(const FFTJetPileupAnalyzer &)
edm::EventID id() const
Definition: EventBase.h:59
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
static std::atomic< unsigned int > counter
edm::EDGetTokenT< double > fastJetRhoToken
tuple cout
Definition: gather_cfg.py:145
unsigned int RunNumber_t
edm::EDGetTokenT< reco::DiscretizedEnergyFlow > gridToken
std::vector< float > ntupleData
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
const std::vector< float > & getPU_sumpT_lowpT() const