CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FFTJetEFlowSmoother.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: FFTJetProducers
4 // Class: FFTJetEFlowSmoother
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Thu Jun 2 18:49:49 CDT 2011
16 //
17 //
18 
19 #include <cmath>
20 
21 // FFTJet headers
22 #include "fftjet/FrequencyKernelConvolver.hh"
23 #include "fftjet/DiscreteGauss2d.hh"
24 #include "fftjet/EquidistantSequence.hh"
25 #include "fftjet/interpolate.hh"
26 #include "fftjet/FrequencyKernelConvolver.hh"
27 #include "fftjet/FrequencySequentialConvolver.hh"
28 #include "fftjet/DiscreteGauss1d.hh"
29 #include "fftjet/DiscreteGauss2d.hh"
30 
31 // framework include files
35 
38 
39 // parameter parser header
41 
42 // useful utilities collected in the second base
44 
45 using namespace fftjetcms;
46 
47 //
48 // class declaration
49 //
51 {
52 public:
53  explicit FFTJetEFlowSmoother(const edm::ParameterSet&);
55 
56 protected:
57  // methods
58  void beginJob() override ;
59  void produce(edm::Event&, const edm::EventSetup&) override;
60  void endJob() override ;
61 
62 private:
65  FFTJetEFlowSmoother& operator=(const FFTJetEFlowSmoother&);
66 
67  void buildKernelConvolver(const edm::ParameterSet&);
68 
69  // Storage for convolved energy flow
70  std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> > convolvedFlow;
71 
72  // Filtering scales
73  std::auto_ptr<std::vector<double> > iniScales;
74 
75  // The FFT engine(s)
76  std::auto_ptr<MyFFTEngine> engine;
77  std::auto_ptr<MyFFTEngine> anotherEngine;
78 
79  // The pattern recognition kernel(s)
80  std::auto_ptr<fftjet::AbsFrequencyKernel> kernel2d;
81  std::auto_ptr<fftjet::AbsFrequencyKernel1d> etaKernel;
82  std::auto_ptr<fftjet::AbsFrequencyKernel1d> phiKernel;
83 
84  // The kernel convolver
85  std::auto_ptr<fftjet::AbsConvolverBase<Real> > convolver;
86 
87  // The scale power to use for scaling Et
88  double scalePower;
89 
90  // Overall factor to multiply Et with
92 };
93 
94 //
95 // constructors and destructor
96 //
98  : FFTJetInterface(ps),
99  scalePower(ps.getParameter<double>("scalePower")),
100  etConversionFactor(ps.getParameter<double>("etConversionFactor"))
101 {
102  // Build the discretization grid
104  ps.getParameter<edm::ParameterSet>("GridConfiguration"));
105  checkConfig(energyFlow, "invalid discretization grid");
106 
107  // Copy of the grid which will be used for convolutions
108  convolvedFlow = std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >(
109  new fftjet::Grid2d<fftjetcms::Real>(*energyFlow));
110 
111  // Build the initial set of pattern recognition scales
113  ps.getParameter<edm::ParameterSet>("InitialScales"));
114  checkConfig(iniScales, "invalid set of scales");
115 
116  // Build the FFT engine(s), pattern recognition kernel(s),
117  // and the kernel convolver
119 
120  // Build the set of pattern recognition scales
121  produces<TH3F>(outputLabel);
122 }
123 
124 
126 {
127  // Check the parameter named "etaDependentScaleFactors". If the vector
128  // of scales is empty we will use 2d kernel, otherwise use 1d kernels
129  const std::vector<double> etaDependentScaleFactors(
130  ps.getParameter<std::vector<double> >("etaDependentScaleFactors"));
131 
132  // Make sure that the number of scale factors provided is correct
133  const bool use2dKernel = etaDependentScaleFactors.empty();
134  if (!use2dKernel)
135  if (etaDependentScaleFactors.size() != energyFlow->nEta())
136  throw cms::Exception("FFTJetBadConfig")
137  << "invalid number of eta-dependent scale factors"
138  << std::endl;
139 
140  // Get the eta and phi scales for the kernel(s)
141  double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
142  const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
143  if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
144  throw cms::Exception("FFTJetBadConfig")
145  << "invalid kernel scale" << std::endl;
146 
147  // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
148  // kernel scale in eta to compensate.
149  kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
150 
151  // Are we going to try to fix the efficiency near detector edges?
152  const bool fixEfficiency = ps.getParameter<bool>("fixEfficiency");
153 
154  // Minimum and maximum eta bin for the convolver
155  unsigned convolverMinBin = 0, convolverMaxBin = 0;
156  if (fixEfficiency || !use2dKernel)
157  {
158  convolverMinBin = ps.getParameter<unsigned>("convolverMinBin");
159  convolverMaxBin = ps.getParameter<unsigned>("convolverMaxBin");
160  }
161 
162  if (use2dKernel)
163  {
164  // Build the FFT engine
165  engine = std::auto_ptr<MyFFTEngine>(
166  new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
167 
168  // 2d kernel
169  kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
170  new fftjet::DiscreteGauss2d(
171  kernelEtaScale, kernelPhiScale,
172  energyFlow->nEta(), energyFlow->nPhi()));
173 
174  // 2d convolver
175  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
176  new fftjet::FrequencyKernelConvolver<Real,Complex>(
177  engine.get(), kernel2d.get(),
178  convolverMinBin, convolverMaxBin));
179  }
180  else
181  {
182  // Need separate FFT engines for eta and phi
183  engine = std::auto_ptr<MyFFTEngine>(
184  new MyFFTEngine(1, energyFlow->nEta()));
185  anotherEngine = std::auto_ptr<MyFFTEngine>(
186  new MyFFTEngine(1, energyFlow->nPhi()));
187 
188  // 1d kernels
189  etaKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
190  new fftjet::DiscreteGauss1d(kernelEtaScale, energyFlow->nEta()));
191 
192  phiKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
193  new fftjet::DiscreteGauss1d(kernelPhiScale, energyFlow->nPhi()));
194 
195  // Sequential convolver
196  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
197  new fftjet::FrequencySequentialConvolver<Real,Complex>(
198  engine.get(), anotherEngine.get(),
199  etaKernel.get(), phiKernel.get(),
200  etaDependentScaleFactors, convolverMinBin,
201  convolverMaxBin, fixEfficiency));
202  }
203 }
204 
205 
207 {
208 }
209 
210 
211 // ------------ method called to produce the data ------------
213  edm::Event& iEvent, const edm::EventSetup& iSetup)
214 {
215  loadInputCollection(iEvent);
217 
218  // Various useful variables
219  const fftjet::Grid2d<Real>& g(*energyFlow);
220  const unsigned nScales = iniScales->size();
221  const double* scales = &(*iniScales)[0];
222  Real* convData = const_cast<Real*>(convolvedFlow->data());
223  const unsigned nEta = g.nEta();
224  const unsigned nPhi = g.nPhi();
225  const double bin0edge = g.phiBin0Edge();
226 
227  // We will fill the following histo
228  std::auto_ptr<TH3F> pTable(
229  new TH3F("FFTJetEFlowSmoother", "FFTJetEFlowSmoother",
230  nScales+1U, -1.5, nScales-0.5,
231  nEta, g.etaMin(), g.etaMax(),
232  nPhi, bin0edge, bin0edge+2.0*M_PI));
233  TH3F* h = pTable.get();
234  h->SetDirectory(0);
235  h->GetXaxis()->SetTitle("Scale");
236  h->GetYaxis()->SetTitle("Eta");
237  h->GetZaxis()->SetTitle("Phi");
238 
239  // Fill the original thing
240  double factor = etConversionFactor*pow(getEventScale(), scalePower);
241  for (unsigned ieta=0; ieta<nEta; ++ieta)
242  for (unsigned iphi=0; iphi<nPhi; ++iphi)
243  h->SetBinContent(1U, ieta+1U, iphi+1U,
244  factor*g.binValue(ieta, iphi));
245 
246  // Go over all scales and perform the convolutions
247  convolver->setEventData(g.data(), nEta, nPhi);
248  for (unsigned iscale=0; iscale<nScales; ++iscale)
249  {
250  factor = etConversionFactor*pow(scales[iscale], scalePower);
251 
252  // Perform the convolution
253  convolver->convolveWithKernel(
254  scales[iscale], convData,
255  convolvedFlow->nEta(), convolvedFlow->nPhi());
256 
257  // Fill the output histo
258  for (unsigned ieta=0; ieta<nEta; ++ieta)
259  {
260  const Real* etaData = convData + ieta*nPhi;
261  for (unsigned iphi=0; iphi<nPhi; ++iphi)
262  h->SetBinContent(iscale+2U, ieta+1U, iphi+1U,
263  factor*etaData[iphi]);
264  }
265  }
266 
267  iEvent.put(pTable, outputLabel);
268 }
269 
270 
271 // ------------ method called once each job just before starting event loop
273 {
274 }
275 
276 
277 // ------------ method called once each job just after ending the event loop
279 {
280 }
281 
282 
283 //define this as a plug-in
T getParameter(std::string const &) const
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
std::auto_ptr< MyFFTEngine > anotherEngine
std::auto_ptr< fftjet::AbsFrequencyKernel1d > etaKernel
std::auto_ptr< fftjet::Grid2d< Real > > fftjet_Grid2d_parser(const edm::ParameterSet &ps)
std::auto_ptr< fftjet::AbsFrequencyKernel1d > phiKernel
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void loadInputCollection(const edm::Event &)
std::auto_ptr< MyFFTEngine > engine
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
std::auto_ptr< fftjet::AbsFrequencyKernel > kernel2d
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
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< fftjet::Grid2d< fftjetcms::Real > > convolvedFlow
std::auto_ptr< std::vector< double > > iniScales
#define M_PI
double Real
fftjet::FFTWDoubleEngine MyFFTEngine
void buildKernelConvolver(const edm::ParameterSet &)
std::auto_ptr< fftjet::AbsConvolverBase< Real > > convolver
const std::string outputLabel
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
std::auto_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)