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: FFTJetProducers
4 // Class: FFTJetPileupProcessor
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Wed Apr 20 13:52:23 CDT 2011
16 // $Id: FFTJetPileupProcessor.cc,v 1.8 2011/07/05 08:24:15 igv Exp $
17 //
18 //
19 
20 #include <cmath>
21 #include <fstream>
22 
23 // FFTJet headers
24 #include "fftjet/FrequencyKernelConvolver.hh"
25 #include "fftjet/DiscreteGauss2d.hh"
26 #include "fftjet/EquidistantSequence.hh"
27 #include "fftjet/interpolate.hh"
28 
29 // framework include files
34 
38 
40 
41 // parameter parser header
43 
44 // useful utilities collected in the second base
46 
47 using namespace fftjetcms;
48 
49 //
50 // class declaration
51 //
53 {
54 public:
57 
58 protected:
59  // methods
60  void beginJob() ;
61  void produce(edm::Event&, const edm::EventSetup&);
62  void endJob() ;
63 
64 private:
68 
69  void buildKernelConvolver(const edm::ParameterSet&);
70  void mixExtraGrid();
71 
72  // The FFT engine
73  std::auto_ptr<MyFFTEngine> engine;
74 
75  // The pattern recognition kernel(s)
76  std::auto_ptr<fftjet::AbsFrequencyKernel> kernel2d;
77 
78  // The kernel convolver
79  std::auto_ptr<fftjet::AbsConvolverBase<Real> > convolver;
80 
81  // Storage for convolved energy flow
82  std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> > convolvedFlow;
83 
84  // Filtering scales
85  std::auto_ptr<fftjet::EquidistantInLogSpace> filterScales;
86 
87  // Eta-dependent factors to use for flattening the distribution
88  // _after_ the filtering
89  std::vector<double> etaFlatteningFactors;
90 
91  // Number of percentile points to use
92  unsigned nPercentiles;
93 
94  // Bin range. Both values of 0 means use the whole range.
95  unsigned convolverMinBin;
96  unsigned convolverMaxBin;
97 
98  // Density conversion factor
100 
101  // Variable related to mixing additional grids
102  std::vector<std::string> externalGridFiles;
103  std::ifstream gridStream;
104  unsigned currentFileNum;
105 };
106 
107 //
108 // constructors and destructor
109 //
111  : FFTJetInterface(ps),
112  etaFlatteningFactors(
113  ps.getParameter<std::vector<double> >("etaFlatteningFactors")),
114  nPercentiles(ps.getParameter<unsigned>("nPercentiles")),
115  convolverMinBin(ps.getParameter<unsigned>("convolverMinBin")),
116  convolverMaxBin(ps.getParameter<unsigned>("convolverMaxBin")),
117  pileupEtaPhiArea(ps.getParameter<double>("pileupEtaPhiArea")),
118  externalGridFiles(ps.getParameter<std::vector<std::string> >("externalGridFiles")),
119  currentFileNum(externalGridFiles.size() + 1U)
120 {
121  // Build the discretization grid
123  ps.getParameter<edm::ParameterSet>("GridConfiguration"));
124  checkConfig(energyFlow, "invalid discretization grid");
125 
126  // Copy of the grid which will be used for convolutions
127  convolvedFlow = std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >(
128  new fftjet::Grid2d<fftjetcms::Real>(*energyFlow));
129 
130  // Make sure the size of flattening factors is appropriate
131  if (!etaFlatteningFactors.empty())
132  {
133  if (etaFlatteningFactors.size() != convolvedFlow->nEta())
134  throw cms::Exception("FFTJetBadConfig")
135  << "ERROR in FFTJetPileupProcessor constructor:"
136  " number of elements in the \"etaFlatteningFactors\""
137  " vector is inconsistent with the iscretization grid binning"
138  << std::endl;
139  }
140 
141  // Build the FFT engine(s), pattern recognition kernel(s),
142  // and the kernel convolver
144 
145  // Build the set of pattern recognition scales
146  const unsigned nScales = ps.getParameter<unsigned>("nScales");
147  const double minScale = ps.getParameter<double>("minScale");
148  const double maxScale = ps.getParameter<double>("maxScale");
149  if (minScale <= 0.0 || maxScale < minScale || nScales == 0U)
150  throw cms::Exception("FFTJetBadConfig")
151  << "invalid filter scales" << std::endl;
152 
153  filterScales = std::auto_ptr<fftjet::EquidistantInLogSpace>(
154  new fftjet::EquidistantInLogSpace(minScale, maxScale, nScales));
155 
156  produces<TH2D>(outputLabel);
157  produces<std::pair<double,double> >(outputLabel);
158 }
159 
160 
162 {
163  // Get the eta and phi scales for the kernel(s)
164  double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
165  const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
166  if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
167  throw cms::Exception("FFTJetBadConfig")
168  << "invalid kernel scales" << std::endl;
169 
172  convolverMaxBin > energyFlow->nEta())
173  throw cms::Exception("FFTJetBadConfig")
174  << "invalid convolver bin range" << std::endl;
175 
176  // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
177  // kernel scale in eta to compensate.
178  kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
179 
180  // Build the FFT engine
181  engine = std::auto_ptr<MyFFTEngine>(
182  new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
183 
184  // 2d kernel
185  kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
186  new fftjet::DiscreteGauss2d(
187  kernelEtaScale, kernelPhiScale,
188  energyFlow->nEta(), energyFlow->nPhi()));
189 
190  // 2d convolver
191  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
192  new fftjet::FrequencyKernelConvolver<Real,Complex>(
193  engine.get(), kernel2d.get(),
195 }
196 
197 
199 {
200 }
201 
202 
203 // ------------ method called to produce the data ------------
205  edm::Event& iEvent, const edm::EventSetup& iSetup)
206 {
207  loadInputCollection(iEvent);
209 
210  // Determine the average Et density for this event.
211  // Needs to be done here, before mixing in another grid.
212  const fftjet::Grid2d<Real>& g(*energyFlow);
213  const double densityBeforeMixing = g.sum()/pileupEtaPhiArea;
214 
215  // Mix an extra grid (if requested)
216  double densityAfterMixing = -1.0;
217  if (!externalGridFiles.empty())
218  {
219  mixExtraGrid();
220  densityAfterMixing = g.sum()/pileupEtaPhiArea;
221  }
222 
223  // Various useful variables
224  const unsigned nScales = filterScales->size();
225  const double* scales = &(*filterScales)[0];
226  Real* convData = const_cast<Real*>(convolvedFlow->data());
227  Real* sortData = convData + convolverMinBin*convolvedFlow->nPhi();
228  const unsigned dataLen = convolverMaxBin ?
230  convolvedFlow->nEta()*convolvedFlow->nPhi();
231 
232  // We will fill the following histo
233  std::auto_ptr<TH2D> pTable(new TH2D("FFTJetPileupProcessor",
234  "FFTJetPileupProcessor",
235  nScales, -0.5, nScales-0.5,
236  nPercentiles, 0.0, 1.0));
237  pTable->SetDirectory(0);
238  pTable->GetXaxis()->SetTitle("Filter Number");
239  pTable->GetYaxis()->SetTitle("Et CDF");
240  pTable->GetZaxis()->SetTitle("Et Density");
241 
242  // Go over all scales and perform the convolutions
243  convolver->setEventData(g.data(), g.nEta(), g.nPhi());
244  for (unsigned iscale=0; iscale<nScales; ++iscale)
245  {
246  // Perform the convolution
247  convolver->convolveWithKernel(
248  scales[iscale], convData,
249  convolvedFlow->nEta(), convolvedFlow->nPhi());
250 
251  // Apply the flattening factors
252  if (!etaFlatteningFactors.empty())
253  convolvedFlow->scaleData(&etaFlatteningFactors[0],
254  etaFlatteningFactors.size());
255 
256  // Sort the convolved data
257  std::sort(sortData, sortData+dataLen);
258 
259  // Determine the percentile points
260  for (unsigned iper=0; iper<nPercentiles; ++iper)
261  {
262  // Map percentile 0 into point number 0,
263  // 1 into point number dataLen-1
264  const double q = (iper + 0.5)/nPercentiles;
265  const double dindex = q*(dataLen-1U);
266  const unsigned ilow = static_cast<unsigned>(std::floor(dindex));
267  const double percentile = fftjet::lin_interpolate_1d(
268  ilow, ilow+1U, sortData[ilow], sortData[ilow+1U], dindex);
269 
270  // Store the calculated percentile
271  pTable->SetBinContent(iscale+1U, iper+1U, percentile);
272  }
273  }
274 
275  iEvent.put(pTable, outputLabel);
276 
277  std::auto_ptr<std::pair<double,double> > etSum(
278  new std::pair<double,double>(densityBeforeMixing, densityAfterMixing));
279  iEvent.put(etSum, outputLabel);
280 }
281 
282 
284 {
285  const unsigned nFiles = externalGridFiles.size();
286  if (currentFileNum > nFiles)
287  {
288  // This is the first time this function is called
289  currentFileNum = 0;
291  std::ios_base::in | std::ios_base::binary);
292  if (!gridStream.is_open())
293  throw cms::Exception("FFTJetBadConfig")
294  << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
295  " failed to open external grid file "
296  << externalGridFiles[currentFileNum] << std::endl;
297  }
298 
299  const fftjet::Grid2d<float>* g = fftjet::Grid2d<float>::read(gridStream);
300 
301  // If we can't read the grid, we need to switch to another file
302  for (unsigned ntries=0; ntries<nFiles && g == 0; ++ntries)
303  {
304  gridStream.close();
305  currentFileNum = (currentFileNum + 1U) % nFiles;
307  std::ios_base::in | std::ios_base::binary);
308  if (!gridStream.is_open())
309  throw cms::Exception("FFTJetBadConfig")
310  << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
311  " failed to open external grid file "
312  << externalGridFiles[currentFileNum] << std::endl;
314  }
315 
316  if (g)
317  {
318  add_Grid2d_data(energyFlow.get(), *g);
319  delete g;
320  }
321  else
322  {
323  // Too bad, no useful file found
324  throw cms::Exception("FFTJetBadConfig")
325  << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
326  " no valid grid records found" << std::endl;
327  }
328 }
329 
330 
331 // ------------ method called once each job just before starting event loop
333 {
334 }
335 
336 
337 // ------------ method called once each job just after ending the event loop
339 {
340 }
341 
342 
343 //define this as a plug-in
nocap nocap const skelname & operator=(const skelname &)
T getParameter(std::string const &) const
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
std::vector< std::string > externalGridFiles
void produce(edm::Event &, const edm::EventSetup &)
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)
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::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > convolvedFlow
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
std::auto_ptr< fftjet::AbsFrequencyKernel > kernel2d
void buildKernelConvolver(const edm::ParameterSet &)
std::vector< double > etaFlatteningFactors
double Real
#define M_PI
Definition: BFit3D.cc:3
std::auto_ptr< MyFFTEngine > engine
fftjet::FFTWDoubleEngine MyFFTEngine
void add_Grid2d_data(fftjet::Grid2d< F2 > *to, const fftjet::Grid2d< F1 > &from)
tuple size
Write out results.
const std::string outputLabel