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 
45 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
46 
47 //
48 // class declaration
49 //
51 {
52 public:
53  explicit FFTJetPileupAnalyzer(const edm::ParameterSet&);
55 
56 private:
60 
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  virtual void beginJob() override ;
66  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
67  virtual 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 //
112  : init_param(edm::InputTag, histoLabel),
113  init_param(edm::InputTag, summaryLabel),
114  init_param(edm::InputTag, fastJetRhoLabel),
115  init_param(edm::InputTag, fastJetSigmaLabel),
116  init_param(edm::InputTag, gridLabel),
117  init_param(edm::InputTag, srcPVs),
118  init_param(std::string, pileupLabel),
119  init_param(std::string, ntupleName),
120  init_param(std::string, ntupleTitle),
121  init_param(bool, collectHistos),
122  init_param(bool, collectSummaries),
123  init_param(bool, collectFastJetRho),
124  init_param(bool, collectPileup),
125  init_param(bool, collectOOTPileup),
126  init_param(bool, collectGrids),
127  init_param(bool, collectGridDensity),
128  init_param(bool, collectVertexInfo),
129  init_param(bool, verbosePileupInfo),
130  init_param(double, vertexNdofCut),
131  init_param(double, crazyEnergyCut),
132  nt(0),
133  totalNpu(-1),
134  totalNPV(-1),
135  counter(0)
136 {
138  pileupToken = consumes<std::vector<PileupSummaryInfo> >(pileupLabel);
139 
140  if (collectHistos)
141  histoToken = consumes<TH2D>(histoLabel);
142 
143  if (collectSummaries)
144  summaryToken = consumes<reco::FFTJetPileupSummary>(summaryLabel);
145 
146  if (collectFastJetRho)
147  {
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 
162 
164 {
165 }
166 
167 
168 //
169 // member functions
170 //
172  const std::vector<PileupSummaryInfo>& info)
173 {
174  const unsigned nBx = info.size();
175  if (collectPileup)
176  ntupleData.push_back(static_cast<float>(nBx));
177 
178  double sumpt_Lo = 0.0, sumpt_Hi = 0.0;
179  totalNpu = 0;
180 
181  int npu_by_Bx[3] = {0,};
182  double sumpt_Lo_by_Bx[3] = {0.0,}, sumpt_Hi_by_Bx[3] = {0.0,};
183 
184  if (verbosePileupInfo)
185  std::cout << "\n**** Pileup info begin" << std::endl;
186 
187  bool isCrazy = false;
188  for (unsigned ibx = 0; ibx < nBx; ++ibx)
189  {
190  const PileupSummaryInfo& puInfo(info[ibx]);
191 
192  const int bx = puInfo.getBunchCrossing();
193  const int npu = puInfo.getPU_NumInteractions();
194  const std::vector<float>& lopt(puInfo.getPU_sumpT_lowpT());
195  const std::vector<float>& hipt(puInfo.getPU_sumpT_highpT());
196  const double losum = std::accumulate(lopt.begin(), lopt.end(), 0.0);
197  const double hisum = std::accumulate(hipt.begin(), hipt.end(), 0.0);
198 
199  if (losum >= crazyEnergyCut)
200  isCrazy = true;
201  if (hisum >= crazyEnergyCut)
202  isCrazy = true;
203 
204  totalNpu += npu;
205  sumpt_Lo += losum;
206  sumpt_Hi += hisum;
207 
208  const unsigned idx = bx < 0 ? 0U : (bx == 0 ? 1U : 2U);
209  npu_by_Bx[idx] += npu;
210  sumpt_Lo_by_Bx[idx] += losum;
211  sumpt_Hi_by_Bx[idx] += hisum;
212 
213  if (verbosePileupInfo)
214  std::cout << "ibx " << ibx << " bx " << bx
215  << " npu " << npu << " losum " << losum
216  << " hisum " << hisum
217  << std::endl;
218  }
219 
220  if (verbosePileupInfo)
221  std::cout << "**** Pileup info end\n" << std::endl;
222 
223  if (isCrazy)
224  {
225  totalNpu = -1;
226  sumpt_Lo = 0.0;
227  sumpt_Hi = 0.0;
228  for (unsigned ibx = 0; ibx < 3; ++ibx)
229  {
230  npu_by_Bx[ibx] = -1;
231  sumpt_Lo_by_Bx[ibx] = 0.0;
232  sumpt_Hi_by_Bx[ibx] = 0.0;
233  }
234  }
235 
236  if (collectPileup)
237  {
238  ntupleData.push_back(totalNpu);
239  ntupleData.push_back(sumpt_Lo);
240  ntupleData.push_back(sumpt_Hi);
241  }
242 
243  if (collectOOTPileup)
244  for (unsigned ibx = 0; ibx < 3; ++ibx)
245  {
246  ntupleData.push_back(npu_by_Bx[ibx]);
247  ntupleData.push_back(sumpt_Lo_by_Bx[ibx]);
248  ntupleData.push_back(sumpt_Hi_by_Bx[ibx]);
249  }
250 }
251 
252 
253 // ------------ method called once each job just before starting event loop
255 {
256  // Come up with the list of variables
257  std::string vars = "cnt:run:event";
258  if (collectPileup)
259  vars += ":nbx:npu:sumptLowCut:sumptHiCut";
260  if (collectOOTPileup)
261  {
262  vars += ":npu_negbx:sumptLowCut_negbx:sumptHiCut_negbx";
263  vars += ":npu_0bx:sumptLowCut_0bx:sumptHiCut_0bx";
264  vars += ":npu_posbx:sumptLowCut_posbx:sumptHiCut_posbx";
265  }
266  if (collectSummaries)
267  vars += ":estimate:pileup:uncert:uncertCode";
268  if (collectFastJetRho)
269  vars += ":fjrho:fjsigma";
270  if (collectGridDensity)
271  vars += ":gridEtDensity:gridEtDensityMixed";
272  if (collectVertexInfo)
273  vars += ":nPV";
274 
275  // Book the ntuple
277  nt = fs->make<TNtuple>(ntupleName.c_str(), ntupleTitle.c_str(),
278  vars.c_str());
279  ntupleData.reserve(nt->GetNvar());
280 }
281 
282 
283 // ------------ method called to for each event ------------
285  const edm::EventSetup& iSetup)
286 {
287  ntupleData.clear();
288  ntupleData.push_back(counter);
289  totalNpu = -1;
290  totalNPV = -1;
291 
292  const long runnumber = iEvent.id().run();
293  const long eventnumber = iEvent.id().event();
294  ntupleData.push_back(runnumber);
295  ntupleData.push_back(eventnumber);
296 
297  // Get pileup information from the pile-up information module
299  {
301  if (iEvent.getByToken(pileupToken, puInfo))
302  analyzePileup(*puInfo);
303  else
304  {
305  if (collectPileup)
306  {
307  ntupleData.push_back(-1);
308  ntupleData.push_back(-1);
309  ntupleData.push_back(0.f);
310  ntupleData.push_back(0.f);
311  }
312  if (collectOOTPileup)
313  for (unsigned ibx = 0; ibx < 3; ++ibx)
314  {
315  ntupleData.push_back(-1);
316  ntupleData.push_back(0.f);
317  ntupleData.push_back(0.f);
318  }
319  }
320  }
321 
322  if (collectHistos)
323  {
325  iEvent.getByToken(histoToken, input);
326 
328  TH2D* copy = new TH2D(*input);
329 
330  std::ostringstream os;
331  os << copy->GetName() << '_' << counter << '_'
332  << totalNpu << '_' << runnumber << '_' << eventnumber;
333  const std::string& newname(os.str());
334  copy->SetNameTitle(newname.c_str(), newname.c_str());
335 
336  copy->SetDirectory(fs->getBareDirectory());
337  }
338 
339  if (collectSummaries)
340  {
342  iEvent.getByToken(summaryToken, summary);
343 
344  ntupleData.push_back(summary->uncalibratedQuantile());
345  ntupleData.push_back(summary->pileupRho());
346  ntupleData.push_back(summary->pileupRhoUncertainty());
347  ntupleData.push_back(summary->uncertaintyCode());
348  }
349 
350  if (collectFastJetRho)
351  {
352  edm::Handle<double> fjrho, fjsigma;
353  iEvent.getByToken(fastJetRhoToken, fjrho);
354  iEvent.getByToken(fastJetSigmaToken, fjsigma);
355 
356  ntupleData.push_back(*fjrho);
357  ntupleData.push_back(*fjsigma);
358  }
359 
360  if (collectGrids)
361  {
363  iEvent.getByToken(gridToken, input);
364 
365  // Make sure the input grid is reasonable
366  const double* data = input->data();
367  assert(data);
368  assert(input->phiBin0Edge() == 0.0);
369  const unsigned nEta = input->nEtaBins();
370  const unsigned nPhi = input->nPhiBins();
371 
372  // Generate a name for the output histogram
373  std::ostringstream os;
374  os << "FFTJetGrid_" << counter << '_'
375  << totalNpu << '_' << runnumber << '_' << eventnumber;
376  const std::string& newname(os.str());
377 
378  // Make a histogram and copy the grid data into it
380  TH2F* h = fs->make<TH2F>(newname.c_str(), newname.c_str(),
381  nEta, input->etaMin(), input->etaMax(),
382  nPhi, 0.0, 2.0*M_PI);
383  h->GetXaxis()->SetTitle("Eta");
384  h->GetYaxis()->SetTitle("Phi");
385  h->GetZaxis()->SetTitle("Transverse Energy");
386 
387  for (unsigned ieta=0; ieta<nEta; ++ieta)
388  for (unsigned iphi=0; iphi<nPhi; ++iphi)
389  h->SetBinContent(ieta+1U, iphi+1U, data[ieta*nPhi + iphi]);
390  }
391 
392  if (collectGridDensity)
393  {
395  iEvent.getByToken(etSumToken, etSum);
396 
397  ntupleData.push_back(etSum->first);
398  ntupleData.push_back(etSum->second);
399  }
400 
401  if (collectVertexInfo)
402  {
404  iEvent.getByToken(srcPVsToken, pvCollection);
405  totalNPV = 0;
406  if (!pvCollection->empty())
407  for (reco::VertexCollection::const_iterator pv = pvCollection->begin();
408  pv != pvCollection->end(); ++pv)
409  {
410  const double ndof = pv->ndof();
411  if (!pv->isFake() && ndof > vertexNdofCut)
412  ++totalNPV;
413  }
414  ntupleData.push_back(totalNPV);
415  }
416 
417  assert(ntupleData.size() == static_cast<unsigned>(nt->GetNvar()));
418  nt->Fill(&ntupleData[0]);
419 
420  ++counter;
421 }
422 
423 
424 // ------------ method called once each job just after ending the event loop
426 {
427 }
428 
429 
430 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:42
EventNumber_t event() const
Definition: EventID.h:44
static const TGPicture * info(bool iBackgroundIsBlack)
#define init_param(type, varname)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#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
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]
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
#define M_PI
edm::EDGetTokenT< double > fastJetSigmaToken
int nt
Definition: AMPTWrapper.h:32
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:56
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:121
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