CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FFTJetPileupProcessor.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoJets/FFTJetProducers
4 // Class: FFTJetPileupProcessor
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Wed Apr 20 13:52:23 CDT 2011
16 //
17 //
18 
19 #include <cmath>
20 #include <fstream>
21 
22 // FFTJet headers
23 #include "fftjet/FrequencyKernelConvolver.hh"
24 #include "fftjet/DiscreteGauss2d.hh"
25 #include "fftjet/EquidistantSequence.hh"
26 #include "fftjet/interpolate.hh"
27 
28 // framework include files
33 
36 // #include "DataFormats/Histograms/interface/MEtoEDMFormat.h"
38 
40 
41 // parameter parser header
43 
44 // useful utilities collected in the second base
46 
47 // Loader for the lookup tables
49 
50 
51 using namespace fftjetcms;
52 
53 //
54 // class declaration
55 //
57 {
58 public:
61 
62 protected:
63  // methods
64  void beginJob() override ;
65  void produce(edm::Event&, const edm::EventSetup&) override;
66  void endJob() override ;
67 
68 private:
72 
73  void buildKernelConvolver(const edm::ParameterSet&);
74  void mixExtraGrid();
75  void loadFlatteningFactors(const edm::EventSetup& iSetup);
76 
77  // The FFT engine
78  std::auto_ptr<MyFFTEngine> engine;
79 
80  // The pattern recognition kernel(s)
81  std::auto_ptr<fftjet::AbsFrequencyKernel> kernel2d;
82 
83  // The kernel convolver
84  std::auto_ptr<fftjet::AbsConvolverBase<Real> > convolver;
85 
86  // Storage for convolved energy flow
87  std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> > convolvedFlow;
88 
89  // Filtering scales
90  std::auto_ptr<fftjet::EquidistantInLogSpace> filterScales;
91 
92  // Eta-dependent factors to use for flattening the distribution
93  // _after_ the filtering
94  std::vector<double> etaFlatteningFactors;
95 
96  // Number of percentile points to use
97  unsigned nPercentiles;
98 
99  // Bin range. Both values of 0 means use the whole range.
100  unsigned convolverMinBin;
101  unsigned convolverMaxBin;
102 
103  // Density conversion factor
105 
106  // Variable related to mixing additional grids
107  std::vector<std::string> externalGridFiles;
108  std::ifstream gridStream;
110  unsigned currentFileNum;
111 
112  // Some memory to hold the percentiles found
113  std::vector<double> percentileData;
114 
115  // Variables to load the flattening factors from
116  // the calibration database (this has to be done
117  // in sync with configuring the appropriate event
118  // setup source)
123 };
124 
125 //
126 // constructors and destructor
127 //
129  : FFTJetInterface(ps),
130  etaFlatteningFactors(
131  ps.getParameter<std::vector<double> >("etaFlatteningFactors")),
132  nPercentiles(ps.getParameter<unsigned>("nPercentiles")),
133  convolverMinBin(ps.getParameter<unsigned>("convolverMinBin")),
134  convolverMaxBin(ps.getParameter<unsigned>("convolverMaxBin")),
135  pileupEtaPhiArea(ps.getParameter<double>("pileupEtaPhiArea")),
136  externalGridFiles(ps.getParameter<std::vector<std::string> >("externalGridFiles")),
137  externalGridMaxEnergy(ps.getParameter<double>("externalGridMaxEnergy")),
138  currentFileNum(externalGridFiles.size() + 1U),
139  flatteningTableRecord(ps.getParameter<std::string>("flatteningTableRecord")),
140  flatteningTableName(ps.getParameter<std::string>("flatteningTableName")),
141  flatteningTableCategory(ps.getParameter<std::string>("flatteningTableCategory")),
142  loadFlatteningFactorsFromDB(ps.getParameter<bool>("loadFlatteningFactorsFromDB"))
143 {
144  // Build the discretization grid
146  ps.getParameter<edm::ParameterSet>("GridConfiguration"));
147  checkConfig(energyFlow, "invalid discretization grid");
148 
149  // Copy of the grid which will be used for convolutions
150  convolvedFlow = std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >(
151  new fftjet::Grid2d<fftjetcms::Real>(*energyFlow));
152 
153  // Make sure the size of flattening factors is appropriate
154  if (!etaFlatteningFactors.empty())
155  {
156  if (etaFlatteningFactors.size() != convolvedFlow->nEta())
157  throw cms::Exception("FFTJetBadConfig")
158  << "ERROR in FFTJetPileupProcessor constructor:"
159  " number of elements in the \"etaFlatteningFactors\""
160  " vector is inconsistent with the discretization grid binning"
161  << std::endl;
162  }
163 
164  // Build the FFT engine(s), pattern recognition kernel(s),
165  // and the kernel convolver
167 
168  // Build the set of pattern recognition scales
169  const unsigned nScales = ps.getParameter<unsigned>("nScales");
170  const double minScale = ps.getParameter<double>("minScale");
171  const double maxScale = ps.getParameter<double>("maxScale");
172  if (minScale <= 0.0 || maxScale < minScale || nScales == 0U)
173  throw cms::Exception("FFTJetBadConfig")
174  << "invalid filter scales" << std::endl;
175 
176  filterScales = std::auto_ptr<fftjet::EquidistantInLogSpace>(
177  new fftjet::EquidistantInLogSpace(minScale, maxScale, nScales));
178 
179  percentileData.resize(nScales*nPercentiles);
180 
181  produces<reco::DiscretizedEnergyFlow>(outputLabel);
182  produces<std::pair<double,double> >(outputLabel);
183 }
184 
185 
187 {
188  // Get the eta and phi scales for the kernel(s)
189  double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
190  const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
191  if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
192  throw cms::Exception("FFTJetBadConfig")
193  << "invalid kernel scales" << std::endl;
194 
197  convolverMaxBin > energyFlow->nEta())
198  throw cms::Exception("FFTJetBadConfig")
199  << "invalid convolver bin range" << std::endl;
200 
201  // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
202  // kernel scale in eta to compensate.
203  kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
204 
205  // Build the FFT engine
206  engine = std::auto_ptr<MyFFTEngine>(
207  new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
208 
209  // 2d kernel
210  kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
211  new fftjet::DiscreteGauss2d(
212  kernelEtaScale, kernelPhiScale,
213  energyFlow->nEta(), energyFlow->nPhi()));
214 
215  // 2d convolver
216  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
217  new fftjet::FrequencyKernelConvolver<Real,Complex>(
218  engine.get(), kernel2d.get(),
220 }
221 
222 
224 {
225 }
226 
227 
228 // ------------ method called to produce the data ------------
230  edm::Event& iEvent, const edm::EventSetup& iSetup)
231 {
232  loadInputCollection(iEvent);
234 
235  // Determine the average Et density for this event.
236  // Needs to be done here, before mixing in another grid.
237  const fftjet::Grid2d<Real>& g(*energyFlow);
238  const double densityBeforeMixing = g.sum()/pileupEtaPhiArea;
239 
240  // Mix an extra grid (if requested)
241  double densityAfterMixing = -1.0;
242  if (!externalGridFiles.empty())
243  {
244  mixExtraGrid();
245  densityAfterMixing = g.sum()/pileupEtaPhiArea;
246  }
247 
248  // Various useful variables
249  const unsigned nScales = filterScales->size();
250  const double* scales = &(*filterScales)[0];
251  Real* convData = const_cast<Real*>(convolvedFlow->data());
252  Real* sortData = convData + convolverMinBin*convolvedFlow->nPhi();
253  const unsigned dataLen = convolverMaxBin ?
255  convolvedFlow->nEta()*convolvedFlow->nPhi();
256 
257  // Load flattenning factors from DB (if requested)
259  loadFlatteningFactors(iSetup);
260 
261  // Go over all scales and perform the convolutions
262  convolver->setEventData(g.data(), g.nEta(), g.nPhi());
263  for (unsigned iscale=0; iscale<nScales; ++iscale)
264  {
265  // Perform the convolution
266  convolver->convolveWithKernel(
267  scales[iscale], convData,
268  convolvedFlow->nEta(), convolvedFlow->nPhi());
269 
270  // Apply the flattening factors
271  if (!etaFlatteningFactors.empty())
272  convolvedFlow->scaleData(&etaFlatteningFactors[0],
273  etaFlatteningFactors.size());
274 
275  // Sort the convolved data
276  std::sort(sortData, sortData+dataLen);
277 
278  // Determine the percentile points
279  for (unsigned iper=0; iper<nPercentiles; ++iper)
280  {
281  // Map percentile 0 into point number 0,
282  // 1 into point number dataLen-1
283  const double q = (iper + 0.5)/nPercentiles;
284  const double dindex = q*(dataLen-1U);
285  const unsigned ilow = static_cast<unsigned>(std::floor(dindex));
286  const double percentile = fftjet::lin_interpolate_1d(
287  ilow, ilow+1U, sortData[ilow], sortData[ilow+1U], dindex);
288 
289  // Store the calculated percentile
290  percentileData[iscale*nPercentiles + iper] = percentile;
291  }
292  }
293 
294  // Convert percentile data into a more convenient storable object
295  // and put it into the event record
296  std::auto_ptr<reco::DiscretizedEnergyFlow> pTable(
298  &percentileData[0], "FFTJetPileupProcessor",
299  -0.5, nScales-0.5, 0.0, nScales, nPercentiles));
300  iEvent.put(pTable, outputLabel);
301 
302  std::auto_ptr<std::pair<double,double> > etSum(
303  new std::pair<double,double>(densityBeforeMixing, densityAfterMixing));
304  iEvent.put(etSum, outputLabel);
305 }
306 
307 
309 {
310  const unsigned nFiles = externalGridFiles.size();
311  if (currentFileNum > nFiles)
312  {
313  // This is the first time this function is called
314  currentFileNum = 0;
316  std::ios_base::in | std::ios_base::binary);
317  if (!gridStream.is_open())
318  throw cms::Exception("FFTJetBadConfig")
319  << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
320  " failed to open external grid file "
321  << externalGridFiles[currentFileNum] << std::endl;
322  }
323 
324  const fftjet::Grid2d<float>* g = 0;
325  const unsigned maxFail = 100U;
326  unsigned nEnergyRejected = 0;
327 
328  while(!g)
329  {
331 
332  // If we can't read the grid, we need to switch to another file
333  for (unsigned ntries=0; ntries<nFiles && g == 0; ++ntries)
334  {
335  gridStream.close();
336  currentFileNum = (currentFileNum + 1U) % nFiles;
338  std::ios_base::in | std::ios_base::binary);
339  if (!gridStream.is_open())
340  throw cms::Exception("FFTJetBadConfig")
341  << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
342  " failed to open external grid file "
343  << externalGridFiles[currentFileNum] << std::endl;
345  }
346 
347  if (g)
348  if (g->sum() > externalGridMaxEnergy)
349  {
350  delete g;
351  g = 0;
352  if (++nEnergyRejected >= maxFail)
353  throw cms::Exception("FFTJetBadConfig")
354  << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
355  " too many grids in a row (" << nEnergyRejected
356  << ") failed the maximum energy cut" << std::endl;
357  }
358  }
359 
360  if (g)
361  {
362  add_Grid2d_data(energyFlow.get(), *g);
363  delete g;
364  }
365  else
366  {
367  // Too bad, no useful file found
368  throw cms::Exception("FFTJetBadConfig")
369  << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
370  " no valid grid records found" << std::endl;
371  }
372 }
373 
374 
375 // ------------ method called once each job just before starting event loop
377 {
378 }
379 
380 
381 // ------------ method called once each job just after ending the event loop
383 {
384 }
385 
386 
388 {
391  iSetup, flatteningTableRecord, h);
392  boost::shared_ptr<npstat::StorableMultivariateFunctor> f =
394 
395  // Fill out the table of flattening factors as a function of eta
396  const unsigned nEta = energyFlow->nEta();
397  etaFlatteningFactors.clear();
398  etaFlatteningFactors.reserve(nEta);
399  for (unsigned ieta=0; ieta<nEta; ++ieta)
400  {
401  const double eta = energyFlow->etaBinCenter(ieta);
402  etaFlatteningFactors.push_back((*f)(&eta, 1U));
403  }
404 }
405 
406 
407 //define this as a plug-in
T getParameter(std::string const &) const
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
std::vector< std::string > externalGridFiles
std::auto_ptr< fftjet::EquidistantInLogSpace > filterScales
std::auto_ptr< fftjet::Grid2d< Real > > fftjet_Grid2d_parser(const edm::ParameterSet &ps)
std::auto_ptr< fftjet::AbsConvolverBase< Real > > convolver
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void loadInputCollection(const edm::Event &)
void checkConfig(const Ptr &ptr, const char *message)
T eta() const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
void beginJob()
Definition: Breakpoints.cc:15
int iEvent
Definition: GenABIO.cc:243
std::vector< double > percentileData
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > convolvedFlow
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
double f[11][100]
std::auto_ptr< fftjet::AbsFrequencyKernel > kernel2d
void buildKernelConvolver(const edm::ParameterSet &)
std::vector< double > etaFlatteningFactors
void loadFlatteningFactors(const edm::EventSetup &iSetup)
double Real
#define M_PI
Definition: BFit3D.cc:3
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
void produce(edm::Event &, const edm::EventSetup &) override
std::auto_ptr< MyFFTEngine > engine
fftjet::FFTWDoubleEngine MyFFTEngine
static const Mapper & instance()
void add_Grid2d_data(fftjet::Grid2d< F2 > *to, const fftjet::Grid2d< F1 > &from)
tuple size
Write out results.
const std::string outputLabel