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
32 
35 
37 
38 // parameter parser header
40 
41 // useful utilities collected in the second base
43 
44 // Loader for the lookup tables
46 
47 
48 using namespace fftjetcms;
49 
50 //
51 // class declaration
52 //
54 {
55 public:
58 
59 protected:
60  // methods
61  void beginJob() override ;
62  void produce(edm::Event&, const edm::EventSetup&) override;
63  void endJob() override ;
64 
65 private:
69 
70  void buildKernelConvolver(const edm::ParameterSet&);
71  void mixExtraGrid();
72  void loadFlatteningFactors(const edm::EventSetup& iSetup);
73 
74  // The FFT engine
75  std::auto_ptr<MyFFTEngine> engine;
76 
77  // The pattern recognition kernel(s)
78  std::auto_ptr<fftjet::AbsFrequencyKernel> kernel2d;
79 
80  // The kernel convolver
81  std::auto_ptr<fftjet::AbsConvolverBase<Real> > convolver;
82 
83  // Storage for convolved energy flow
84  std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> > convolvedFlow;
85 
86  // Filtering scales
87  std::auto_ptr<fftjet::EquidistantInLogSpace> filterScales;
88 
89  // Eta-dependent factors to use for flattening the distribution
90  // _after_ the filtering
91  std::vector<double> etaFlatteningFactors;
92 
93  // Number of percentile points to use
94  unsigned nPercentiles;
95 
96  // Bin range. Both values of 0 means use the whole range.
97  unsigned convolverMinBin;
98  unsigned convolverMaxBin;
99 
100  // Density conversion factor
102 
103  // Variable related to mixing additional grids
104  std::vector<std::string> externalGridFiles;
105  std::ifstream gridStream;
107  unsigned currentFileNum;
108 
109  // Some memory to hold the percentiles found
110  std::vector<double> percentileData;
111 
112  // Variables to load the flattening factors from
113  // the calibration database (this has to be done
114  // in sync with configuring the appropriate event
115  // setup source)
120 };
121 
122 //
123 // constructors and destructor
124 //
126  : FFTJetInterface(ps),
127  etaFlatteningFactors(
128  ps.getParameter<std::vector<double> >("etaFlatteningFactors")),
129  nPercentiles(ps.getParameter<unsigned>("nPercentiles")),
130  convolverMinBin(ps.getParameter<unsigned>("convolverMinBin")),
131  convolverMaxBin(ps.getParameter<unsigned>("convolverMaxBin")),
132  pileupEtaPhiArea(ps.getParameter<double>("pileupEtaPhiArea")),
133  externalGridFiles(ps.getParameter<std::vector<std::string> >("externalGridFiles")),
134  externalGridMaxEnergy(ps.getParameter<double>("externalGridMaxEnergy")),
135  currentFileNum(externalGridFiles.size() + 1U),
136  flatteningTableRecord(ps.getParameter<std::string>("flatteningTableRecord")),
137  flatteningTableName(ps.getParameter<std::string>("flatteningTableName")),
138  flatteningTableCategory(ps.getParameter<std::string>("flatteningTableCategory")),
139  loadFlatteningFactorsFromDB(ps.getParameter<bool>("loadFlatteningFactorsFromDB"))
140 {
141  // Build the discretization grid
143  ps.getParameter<edm::ParameterSet>("GridConfiguration"));
144  checkConfig(energyFlow, "invalid discretization grid");
145 
146  // Copy of the grid which will be used for convolutions
147  convolvedFlow = std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >(
148  new fftjet::Grid2d<fftjetcms::Real>(*energyFlow));
149 
150  // Make sure the size of flattening factors is appropriate
151  if (!etaFlatteningFactors.empty())
152  {
153  if (etaFlatteningFactors.size() != convolvedFlow->nEta())
154  throw cms::Exception("FFTJetBadConfig")
155  << "ERROR in FFTJetPileupProcessor constructor:"
156  " number of elements in the \"etaFlatteningFactors\""
157  " vector is inconsistent with the discretization grid binning"
158  << std::endl;
159  }
160 
161  // Build the FFT engine(s), pattern recognition kernel(s),
162  // and the kernel convolver
164 
165  // Build the set of pattern recognition scales
166  const unsigned nScales = ps.getParameter<unsigned>("nScales");
167  const double minScale = ps.getParameter<double>("minScale");
168  const double maxScale = ps.getParameter<double>("maxScale");
169  if (minScale <= 0.0 || maxScale < minScale || nScales == 0U)
170  throw cms::Exception("FFTJetBadConfig")
171  << "invalid filter scales" << std::endl;
172 
173  filterScales = std::auto_ptr<fftjet::EquidistantInLogSpace>(
174  new fftjet::EquidistantInLogSpace(minScale, maxScale, nScales));
175 
176  percentileData.resize(nScales*nPercentiles);
177 
178  produces<reco::DiscretizedEnergyFlow>(outputLabel);
179  produces<std::pair<double,double> >(outputLabel);
180 }
181 
182 
184 {
185  // Get the eta and phi scales for the kernel(s)
186  double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
187  const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
188  if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
189  throw cms::Exception("FFTJetBadConfig")
190  << "invalid kernel scales" << std::endl;
191 
194  convolverMaxBin > energyFlow->nEta())
195  throw cms::Exception("FFTJetBadConfig")
196  << "invalid convolver bin range" << std::endl;
197 
198  // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
199  // kernel scale in eta to compensate.
200  kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
201 
202  // Build the FFT engine
203  engine = std::auto_ptr<MyFFTEngine>(
204  new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
205 
206  // 2d kernel
207  kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
208  new fftjet::DiscreteGauss2d(
209  kernelEtaScale, kernelPhiScale,
210  energyFlow->nEta(), energyFlow->nPhi()));
211 
212  // 2d convolver
213  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
214  new fftjet::FrequencyKernelConvolver<Real,Complex>(
215  engine.get(), kernel2d.get(),
217 }
218 
219 
221 {
222 }
223 
224 
225 // ------------ method called to produce the data ------------
227  edm::Event& iEvent, const edm::EventSetup& iSetup)
228 {
229  loadInputCollection(iEvent);
231 
232  // Determine the average Et density for this event.
233  // Needs to be done here, before mixing in another grid.
234  const fftjet::Grid2d<Real>& g(*energyFlow);
235  const double densityBeforeMixing = g.sum()/pileupEtaPhiArea;
236 
237  // Mix an extra grid (if requested)
238  double densityAfterMixing = -1.0;
239  if (!externalGridFiles.empty())
240  {
241  mixExtraGrid();
242  densityAfterMixing = g.sum()/pileupEtaPhiArea;
243  }
244 
245  // Various useful variables
246  const unsigned nScales = filterScales->size();
247  const double* scales = &(*filterScales)[0];
248  Real* convData = const_cast<Real*>(convolvedFlow->data());
249  Real* sortData = convData + convolverMinBin*convolvedFlow->nPhi();
250  const unsigned dataLen = convolverMaxBin ?
252  convolvedFlow->nEta()*convolvedFlow->nPhi();
253 
254  // Load flattenning factors from DB (if requested)
256  loadFlatteningFactors(iSetup);
257 
258  // Go over all scales and perform the convolutions
259  convolver->setEventData(g.data(), g.nEta(), g.nPhi());
260  for (unsigned iscale=0; iscale<nScales; ++iscale)
261  {
262  // Perform the convolution
263  convolver->convolveWithKernel(
264  scales[iscale], convData,
265  convolvedFlow->nEta(), convolvedFlow->nPhi());
266 
267  // Apply the flattening factors
268  if (!etaFlatteningFactors.empty())
269  convolvedFlow->scaleData(&etaFlatteningFactors[0],
270  etaFlatteningFactors.size());
271 
272  // Sort the convolved data
273  std::sort(sortData, sortData+dataLen);
274 
275  // Determine the percentile points
276  for (unsigned iper=0; iper<nPercentiles; ++iper)
277  {
278  // Map percentile 0 into point number 0,
279  // 1 into point number dataLen-1
280  const double q = (iper + 0.5)/nPercentiles;
281  const double dindex = q*(dataLen-1U);
282  const unsigned ilow = static_cast<unsigned>(std::floor(dindex));
283  const double percentile = fftjet::lin_interpolate_1d(
284  ilow, ilow+1U, sortData[ilow], sortData[ilow+1U], dindex);
285 
286  // Store the calculated percentile
287  percentileData[iscale*nPercentiles + iper] = percentile;
288  }
289  }
290 
291  // Convert percentile data into a more convenient storable object
292  // and put it into the event record
293  std::auto_ptr<reco::DiscretizedEnergyFlow> pTable(
295  &percentileData[0], "FFTJetPileupProcessor",
296  -0.5, nScales-0.5, 0.0, nScales, nPercentiles));
297  iEvent.put(pTable, outputLabel);
298 
299  std::auto_ptr<std::pair<double,double> > etSum(
300  new std::pair<double,double>(densityBeforeMixing, densityAfterMixing));
301  iEvent.put(etSum, outputLabel);
302 }
303 
304 
306 {
307  const unsigned nFiles = externalGridFiles.size();
308  if (currentFileNum > nFiles)
309  {
310  // This is the first time this function is called
311  currentFileNum = 0;
313  std::ios_base::in | std::ios_base::binary);
314  if (!gridStream.is_open())
315  throw cms::Exception("FFTJetBadConfig")
316  << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
317  " failed to open external grid file "
318  << externalGridFiles[currentFileNum] << std::endl;
319  }
320 
321  const fftjet::Grid2d<float>* g = 0;
322  const unsigned maxFail = 100U;
323  unsigned nEnergyRejected = 0;
324 
325  while(!g)
326  {
328 
329  // If we can't read the grid, we need to switch to another file
330  for (unsigned ntries=0; ntries<nFiles && g == 0; ++ntries)
331  {
332  gridStream.close();
333  currentFileNum = (currentFileNum + 1U) % nFiles;
335  std::ios_base::in | std::ios_base::binary);
336  if (!gridStream.is_open())
337  throw cms::Exception("FFTJetBadConfig")
338  << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
339  " failed to open external grid file "
340  << externalGridFiles[currentFileNum] << std::endl;
342  }
343 
344  if (g)
345  if (g->sum() > externalGridMaxEnergy)
346  {
347  delete g;
348  g = 0;
349  if (++nEnergyRejected >= maxFail)
350  throw cms::Exception("FFTJetBadConfig")
351  << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
352  " too many grids in a row (" << nEnergyRejected
353  << ") failed the maximum energy cut" << std::endl;
354  }
355  }
356 
357  if (g)
358  {
359  add_Grid2d_data(energyFlow.get(), *g);
360  delete g;
361  }
362  else
363  {
364  // Too bad, no useful file found
365  throw cms::Exception("FFTJetBadConfig")
366  << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
367  " no valid grid records found" << std::endl;
368  }
369 }
370 
371 
372 // ------------ method called once each job just before starting event loop
374 {
375 }
376 
377 
378 // ------------ method called once each job just after ending the event loop
380 {
381 }
382 
383 
385 {
388  iSetup, flatteningTableRecord, h);
389  boost::shared_ptr<npstat::StorableMultivariateFunctor> f =
391 
392  // Fill out the table of flattening factors as a function of eta
393  const unsigned nEta = energyFlow->nEta();
394  etaFlatteningFactors.clear();
395  etaFlatteningFactors.reserve(nEta);
396  for (unsigned ieta=0; ieta<nEta; ++ieta)
397  {
398  const double eta = energyFlow->etaBinCenter(ieta);
399  etaFlatteningFactors.push_back((*f)(&eta, 1U));
400  }
401 }
402 
403 
404 //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:230
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:113
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
std::auto_ptr< fftjet::AbsFrequencyKernel > kernel2d
void buildKernelConvolver(const edm::ParameterSet &)
#define M_PI
std::vector< double > etaFlatteningFactors
void loadFlatteningFactors(const edm::EventSetup &iSetup)
double Real
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