CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups 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 #include <memory>
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 #include "fftjet/FrequencyKernelConvolver.hh"
28 #include "fftjet/FrequencySequentialConvolver.hh"
29 #include "fftjet/DiscreteGauss1d.hh"
30 #include "fftjet/DiscreteGauss2d.hh"
31 
32 // framework include files
36 
38 #include <TH3F.h>
39 
40 // parameter parser header
42 
43 // useful utilities collected in the second base
45 
46 using namespace fftjetcms;
47 
48 //
49 // class declaration
50 //
52 public:
53  explicit FFTJetEFlowSmoother(const edm::ParameterSet&);
54  FFTJetEFlowSmoother() = delete;
56  FFTJetEFlowSmoother& operator=(const FFTJetEFlowSmoother&) = delete;
57  ~FFTJetEFlowSmoother() override;
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:
66  void buildKernelConvolver(const edm::ParameterSet&);
67 
68  // Storage for convolved energy flow
69  std::unique_ptr<fftjet::Grid2d<fftjetcms::Real> > convolvedFlow;
70 
71  // Filtering scales
72  std::unique_ptr<std::vector<double> > iniScales;
73 
74  // The FFT engine(s)
75  std::unique_ptr<MyFFTEngine> engine;
76  std::unique_ptr<MyFFTEngine> anotherEngine;
77 
78  // The pattern recognition kernel(s)
79  std::unique_ptr<fftjet::AbsFrequencyKernel> kernel2d;
80  std::unique_ptr<fftjet::AbsFrequencyKernel1d> etaKernel;
81  std::unique_ptr<fftjet::AbsFrequencyKernel1d> phiKernel;
82 
83  // The kernel convolver
84  std::unique_ptr<fftjet::AbsConvolverBase<Real> > convolver;
85 
86  // The scale power to use for scaling Et
87  double scalePower;
88 
89  // Overall factor to multiply Et with
91 };
92 
93 //
94 // constructors and destructor
95 //
97  : FFTJetInterface(ps),
98  scalePower(ps.getParameter<double>("scalePower")),
99  etConversionFactor(ps.getParameter<double>("etConversionFactor")) {
100  // Build the discretization grid
102  checkConfig(energyFlow, "invalid discretization grid");
103 
104  // Copy of the grid which will be used for convolutions
105  convolvedFlow = std::make_unique<fftjet::Grid2d<fftjetcms::Real> >(*energyFlow);
106 
107  // Build the initial set of pattern recognition scales
109  checkConfig(iniScales, "invalid set of scales");
110 
111  // Build the FFT engine(s), pattern recognition kernel(s),
112  // and the kernel convolver
114 
115  // Build the set of pattern recognition scales
116  produces<TH3F>(outputLabel);
117 }
118 
120  // Check the parameter named "etaDependentScaleFactors". If the vector
121  // of scales is empty we will use 2d kernel, otherwise use 1d kernels
122  const std::vector<double> etaDependentScaleFactors(ps.getParameter<std::vector<double> >("etaDependentScaleFactors"));
123 
124  // Make sure that the number of scale factors provided is correct
125  const bool use2dKernel = etaDependentScaleFactors.empty();
126  if (!use2dKernel)
127  if (etaDependentScaleFactors.size() != energyFlow->nEta())
128  throw cms::Exception("FFTJetBadConfig") << "invalid number of eta-dependent scale factors" << std::endl;
129 
130  // Get the eta and phi scales for the kernel(s)
131  double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
132  const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
133  if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
134  throw cms::Exception("FFTJetBadConfig") << "invalid kernel scale" << std::endl;
135 
136  // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
137  // kernel scale in eta to compensate.
138  kernelEtaScale *= (2.0 * M_PI / (energyFlow->etaMax() - energyFlow->etaMin()));
139 
140  // Are we going to try to fix the efficiency near detector edges?
141  const bool fixEfficiency = ps.getParameter<bool>("fixEfficiency");
142 
143  // Minimum and maximum eta bin for the convolver
144  unsigned convolverMinBin = 0, convolverMaxBin = 0;
145  if (fixEfficiency || !use2dKernel) {
146  convolverMinBin = ps.getParameter<unsigned>("convolverMinBin");
147  convolverMaxBin = ps.getParameter<unsigned>("convolverMaxBin");
148  }
149 
150  if (use2dKernel) {
151  // Build the FFT engine
152  engine = std::make_unique<MyFFTEngine>(energyFlow->nEta(), energyFlow->nPhi());
153 
154  // 2d kernel
155  kernel2d = std::unique_ptr<fftjet::AbsFrequencyKernel>(
156  new fftjet::DiscreteGauss2d(kernelEtaScale, kernelPhiScale, energyFlow->nEta(), energyFlow->nPhi()));
157 
158  // 2d convolver
159  convolver = std::unique_ptr<fftjet::AbsConvolverBase<Real> >(new fftjet::FrequencyKernelConvolver<Real, Complex>(
160  engine.get(), kernel2d.get(), convolverMinBin, convolverMaxBin));
161  } else {
162  // Need separate FFT engines for eta and phi
163  engine = std::make_unique<MyFFTEngine>(1, energyFlow->nEta());
164  anotherEngine = std::make_unique<MyFFTEngine>(1, energyFlow->nPhi());
165 
166  // 1d kernels
167  etaKernel =
168  std::unique_ptr<fftjet::AbsFrequencyKernel1d>(new fftjet::DiscreteGauss1d(kernelEtaScale, energyFlow->nEta()));
169 
170  phiKernel =
171  std::unique_ptr<fftjet::AbsFrequencyKernel1d>(new fftjet::DiscreteGauss1d(kernelPhiScale, energyFlow->nPhi()));
172 
173  // Sequential convolver
174  convolver = std::unique_ptr<fftjet::AbsConvolverBase<Real> >(
175  new fftjet::FrequencySequentialConvolver<Real, Complex>(engine.get(),
176  anotherEngine.get(),
177  etaKernel.get(),
178  phiKernel.get(),
179  etaDependentScaleFactors,
180  convolverMinBin,
181  convolverMaxBin,
182  fixEfficiency));
183  }
184 }
185 
187 
188 // ------------ method called to produce the data ------------
190  loadInputCollection(iEvent);
192 
193  // Various useful variables
194  const fftjet::Grid2d<Real>& g(*energyFlow);
195  const unsigned nScales = iniScales->size();
196  const double* scales = &(*iniScales)[0];
197  Real* convData = const_cast<Real*>(convolvedFlow->data());
198  const unsigned nEta = g.nEta();
199  const unsigned nPhi = g.nPhi();
200  const double bin0edge = g.phiBin0Edge();
201 
202  // We will fill the following histo
203  auto pTable = std::make_unique<TH3F>("FFTJetEFlowSmoother",
204  "FFTJetEFlowSmoother",
205  nScales + 1U,
206  -1.5,
207  nScales - 0.5,
208  nEta,
209  g.etaMin(),
210  g.etaMax(),
211  nPhi,
212  bin0edge,
213  bin0edge + 2.0 * M_PI);
214  TH3F* h = pTable.get();
215  h->SetDirectory(nullptr);
216  h->GetXaxis()->SetTitle("Scale");
217  h->GetYaxis()->SetTitle("Eta");
218  h->GetZaxis()->SetTitle("Phi");
219 
220  // Fill the original thing
221  double factor = etConversionFactor * pow(getEventScale(), scalePower);
222  for (unsigned ieta = 0; ieta < nEta; ++ieta)
223  for (unsigned iphi = 0; iphi < nPhi; ++iphi)
224  h->SetBinContent(1U, ieta + 1U, iphi + 1U, factor * g.binValue(ieta, iphi));
225 
226  // Go over all scales and perform the convolutions
227  convolver->setEventData(g.data(), nEta, nPhi);
228  for (unsigned iscale = 0; iscale < nScales; ++iscale) {
229  factor = etConversionFactor * pow(scales[iscale], scalePower);
230 
231  // Perform the convolution
232  convolver->convolveWithKernel(scales[iscale], convData, convolvedFlow->nEta(), convolvedFlow->nPhi());
233 
234  // Fill the output histo
235  for (unsigned ieta = 0; ieta < nEta; ++ieta) {
236  const Real* etaData = convData + ieta * nPhi;
237  for (unsigned iphi = 0; iphi < nPhi; ++iphi)
238  h->SetBinContent(iscale + 2U, ieta + 1U, iphi + 1U, factor * etaData[iphi]);
239  }
240  }
241 
242  iEvent.put(std::move(pTable), outputLabel);
243 }
244 
245 // ------------ method called once each job just before starting event loop
247 
248 // ------------ method called once each job just after ending the event loop
250 
251 //define this as a plug-in
std::unique_ptr< MyFFTEngine > engine
std::unique_ptr< fftjet::AbsFrequencyKernel1d > etaKernel
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
std::unique_ptr< fftjet::AbsConvolverBase< Real > > convolver
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void loadInputCollection(const edm::Event &)
std::unique_ptr< MyFFTEngine > anotherEngine
FFTJetEFlowSmoother()=delete
void checkConfig(const Ptr &ptr, const char *message)
std::unique_ptr< fftjet::Grid2d< Real > > fftjet_Grid2d_parser(const edm::ParameterSet &ps)
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:14
int iEvent
Definition: GenABIO.cc:224
std::unique_ptr< std::vector< double > > iniScales
std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > convolvedFlow
std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
def move
Definition: eostools.py:511
void produce(edm::Event &, const edm::EventSetup &) override
#define M_PI
std::unique_ptr< fftjet::AbsFrequencyKernel > kernel2d
double Real
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void buildKernelConvolver(const edm::ParameterSet &)
std::unique_ptr< fftjet::AbsFrequencyKernel1d > phiKernel
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
std::unique_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)
const std::string outputLabel
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29