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
36 
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:
55  explicit FFTJetEFlowSmoother(const edm::ParameterSet&);
57 
58 protected:
59  // methods
60  void beginJob() override ;
61  void produce(edm::Event&, const edm::EventSetup&) override;
62  void endJob() override ;
63 
64 private:
67  FFTJetEFlowSmoother& operator=(const FFTJetEFlowSmoother&);
68 
69  void buildKernelConvolver(const edm::ParameterSet&);
70 
71  // Storage for convolved energy flow
72  std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> > convolvedFlow;
73 
74  // Filtering scales
75  std::auto_ptr<std::vector<double> > iniScales;
76 
77  // The FFT engine(s)
78  std::auto_ptr<MyFFTEngine> engine;
79  std::auto_ptr<MyFFTEngine> anotherEngine;
80 
81  // The pattern recognition kernel(s)
82  std::auto_ptr<fftjet::AbsFrequencyKernel> kernel2d;
83  std::auto_ptr<fftjet::AbsFrequencyKernel1d> etaKernel;
84  std::auto_ptr<fftjet::AbsFrequencyKernel1d> phiKernel;
85 
86  // The kernel convolver
87  std::auto_ptr<fftjet::AbsConvolverBase<Real> > convolver;
88 
89  // The scale power to use for scaling Et
90  double scalePower;
91 
92  // Overall factor to multiply Et with
94 };
95 
96 //
97 // constructors and destructor
98 //
100  : FFTJetInterface(ps),
101  scalePower(ps.getParameter<double>("scalePower")),
102  etConversionFactor(ps.getParameter<double>("etConversionFactor"))
103 {
104  // Build the discretization grid
106  ps.getParameter<edm::ParameterSet>("GridConfiguration"));
107  checkConfig(energyFlow, "invalid discretization grid");
108 
109  // Copy of the grid which will be used for convolutions
110  convolvedFlow = std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >(
111  new fftjet::Grid2d<fftjetcms::Real>(*energyFlow));
112 
113  // Build the initial set of pattern recognition scales
115  ps.getParameter<edm::ParameterSet>("InitialScales"));
116  checkConfig(iniScales, "invalid set of scales");
117 
118  // Build the FFT engine(s), pattern recognition kernel(s),
119  // and the kernel convolver
121 
122  // Build the set of pattern recognition scales
123  produces<TH3F>(outputLabel);
124 }
125 
126 
128 {
129  // Check the parameter named "etaDependentScaleFactors". If the vector
130  // of scales is empty we will use 2d kernel, otherwise use 1d kernels
131  const std::vector<double> etaDependentScaleFactors(
132  ps.getParameter<std::vector<double> >("etaDependentScaleFactors"));
133 
134  // Make sure that the number of scale factors provided is correct
135  const bool use2dKernel = etaDependentScaleFactors.empty();
136  if (!use2dKernel)
137  if (etaDependentScaleFactors.size() != energyFlow->nEta())
138  throw cms::Exception("FFTJetBadConfig")
139  << "invalid number of eta-dependent scale factors"
140  << std::endl;
141 
142  // Get the eta and phi scales for the kernel(s)
143  double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
144  const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
145  if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
146  throw cms::Exception("FFTJetBadConfig")
147  << "invalid kernel scale" << std::endl;
148 
149  // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
150  // kernel scale in eta to compensate.
151  kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
152 
153  // Are we going to try to fix the efficiency near detector edges?
154  const bool fixEfficiency = ps.getParameter<bool>("fixEfficiency");
155 
156  // Minimum and maximum eta bin for the convolver
157  unsigned convolverMinBin = 0, convolverMaxBin = 0;
158  if (fixEfficiency || !use2dKernel)
159  {
160  convolverMinBin = ps.getParameter<unsigned>("convolverMinBin");
161  convolverMaxBin = ps.getParameter<unsigned>("convolverMaxBin");
162  }
163 
164  if (use2dKernel)
165  {
166  // Build the FFT engine
167  engine = std::auto_ptr<MyFFTEngine>(
168  new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
169 
170  // 2d kernel
171  kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
172  new fftjet::DiscreteGauss2d(
173  kernelEtaScale, kernelPhiScale,
174  energyFlow->nEta(), energyFlow->nPhi()));
175 
176  // 2d convolver
177  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
178  new fftjet::FrequencyKernelConvolver<Real,Complex>(
179  engine.get(), kernel2d.get(),
180  convolverMinBin, convolverMaxBin));
181  }
182  else
183  {
184  // Need separate FFT engines for eta and phi
185  engine = std::auto_ptr<MyFFTEngine>(
186  new MyFFTEngine(1, energyFlow->nEta()));
187  anotherEngine = std::auto_ptr<MyFFTEngine>(
188  new MyFFTEngine(1, energyFlow->nPhi()));
189 
190  // 1d kernels
191  etaKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
192  new fftjet::DiscreteGauss1d(kernelEtaScale, energyFlow->nEta()));
193 
194  phiKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
195  new fftjet::DiscreteGauss1d(kernelPhiScale, energyFlow->nPhi()));
196 
197  // Sequential convolver
198  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
199  new fftjet::FrequencySequentialConvolver<Real,Complex>(
200  engine.get(), anotherEngine.get(),
201  etaKernel.get(), phiKernel.get(),
202  etaDependentScaleFactors, convolverMinBin,
203  convolverMaxBin, fixEfficiency));
204  }
205 }
206 
207 
209 {
210 }
211 
212 
213 // ------------ method called to produce the data ------------
215  edm::Event& iEvent, const edm::EventSetup& iSetup)
216 {
217  loadInputCollection(iEvent);
219 
220  // Various useful variables
221  const fftjet::Grid2d<Real>& g(*energyFlow);
222  const unsigned nScales = iniScales->size();
223  const double* scales = &(*iniScales)[0];
224  Real* convData = const_cast<Real*>(convolvedFlow->data());
225  const unsigned nEta = g.nEta();
226  const unsigned nPhi = g.nPhi();
227  const double bin0edge = g.phiBin0Edge();
228 
229  // We will fill the following histo
230  std::auto_ptr<TH3F> pTable(
231  new TH3F("FFTJetEFlowSmoother", "FFTJetEFlowSmoother",
232  nScales+1U, -1.5, nScales-0.5,
233  nEta, g.etaMin(), g.etaMax(),
234  nPhi, bin0edge, bin0edge+2.0*M_PI));
235  TH3F* h = pTable.get();
236  h->SetDirectory(0);
237  h->GetXaxis()->SetTitle("Scale");
238  h->GetYaxis()->SetTitle("Eta");
239  h->GetZaxis()->SetTitle("Phi");
240 
241  // Fill the original thing
242  double factor = etConversionFactor*pow(getEventScale(), scalePower);
243  for (unsigned ieta=0; ieta<nEta; ++ieta)
244  for (unsigned iphi=0; iphi<nPhi; ++iphi)
245  h->SetBinContent(1U, ieta+1U, iphi+1U,
246  factor*g.binValue(ieta, iphi));
247 
248  // Go over all scales and perform the convolutions
249  convolver->setEventData(g.data(), nEta, nPhi);
250  for (unsigned iscale=0; iscale<nScales; ++iscale)
251  {
252  factor = etConversionFactor*pow(scales[iscale], scalePower);
253 
254  // Perform the convolution
255  convolver->convolveWithKernel(
256  scales[iscale], convData,
257  convolvedFlow->nEta(), convolvedFlow->nPhi());
258 
259  // Fill the output histo
260  for (unsigned ieta=0; ieta<nEta; ++ieta)
261  {
262  const Real* etaData = convData + ieta*nPhi;
263  for (unsigned iphi=0; iphi<nPhi; ++iphi)
264  h->SetBinContent(iscale+2U, ieta+1U, iphi+1U,
265  factor*etaData[iphi]);
266  }
267  }
268 
269  iEvent.put(pTable, outputLabel);
270 }
271 
272 
273 // ------------ method called once each job just before starting event loop
275 {
276 }
277 
278 
279 // ------------ method called once each job just after ending the event loop
281 {
282 }
283 
284 
285 //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:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
void produce(edm::Event &, const edm::EventSetup &) override
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > convolvedFlow
std::auto_ptr< std::vector< double > > iniScales
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
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)