CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
FFTJetEFlowSmoother Class Reference

#include <RecoJets/FFTJetProducers/plugins/FFTJetEFlowSmoother.cc>

Inheritance diagram for FFTJetEFlowSmoother:
edm::EDProducer fftjetcms::FFTJetInterface edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 FFTJetEFlowSmoother (const edm::ParameterSet &)
 
 ~FFTJetEFlowSmoother ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 
- Public Member Functions inherited from fftjetcms::FFTJetInterface
virtual ~FFTJetInterface ()
 

Protected Member Functions

void beginJob ()
 
void endJob ()
 
void produce (edm::Event &, const edm::EventSetup &)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
- Protected Member Functions inherited from fftjetcms::FFTJetInterface
template<class Ptr >
void checkConfig (const Ptr &ptr, const char *message)
 
void discretizeEnergyFlow ()
 
 FFTJetInterface (const edm::ParameterSet &)
 
double getEventScale () const
 
void loadInputCollection (const edm::Event &)
 
bool storeInSinglePrecision () const
 
const reco::Particle::PointvertexUsed () const
 

Private Member Functions

void buildKernelConvolver (const edm::ParameterSet &)
 
 FFTJetEFlowSmoother ()
 
 FFTJetEFlowSmoother (const FFTJetEFlowSmoother &)
 
FFTJetEFlowSmootheroperator= (const FFTJetEFlowSmoother &)
 

Private Attributes

std::auto_ptr< MyFFTEngineanotherEngine
 
std::auto_ptr< fftjet::Grid2d
< fftjetcms::Real > > 
convolvedFlow
 
std::auto_ptr
< fftjet::AbsConvolverBase
< Real > > 
convolver
 
std::auto_ptr< MyFFTEngineengine
 
std::auto_ptr
< fftjet::AbsFrequencyKernel1d > 
etaKernel
 
double etConversionFactor
 
std::auto_ptr< std::vector
< double > > 
iniScales
 
std::auto_ptr
< fftjet::AbsFrequencyKernel > 
kernel2d
 
std::auto_ptr
< fftjet::AbsFrequencyKernel1d > 
phiKernel
 
double scalePower
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Attributes inherited from fftjetcms::FFTJetInterface
const AnomalousTower anomalous
 
std::vector< unsigned > candidateIndex
 
const bool doPVCorrection
 
std::auto_ptr< fftjet::Grid2d
< fftjetcms::Real > > 
energyFlow
 
const std::vector< double > etaDependentMagnutideFactors
 
std::vector
< fftjetcms::VectorLike
eventData
 
edm::Handle< reco::CandidateViewinputCollection
 
const edm::InputTag inputLabel
 
const JetType jetType
 
const std::string outputLabel
 
const edm::InputTag srcPVs
 

Detailed Description

Description: Runs FFTJet filtering code for multiple scales and saves the results

Implementation: [Notes on implementation]

Definition at line 53 of file FFTJetEFlowSmoother.cc.

Constructor & Destructor Documentation

FFTJetEFlowSmoother::FFTJetEFlowSmoother ( const edm::ParameterSet ps)
explicit

Definition at line 100 of file FFTJetEFlowSmoother.cc.

References buildKernelConvolver(), fftjetcms::FFTJetInterface::checkConfig(), convolvedFlow, fftjetcms::FFTJetInterface::energyFlow, fftjetcms::fftjet_Grid2d_parser(), fftjetcms::fftjet_ScaleSet_parser(), edm::ParameterSet::getParameter(), iniScales, and fftjetcms::FFTJetInterface::outputLabel.

101  : FFTJetInterface(ps),
102  scalePower(ps.getParameter<double>("scalePower")),
103  etConversionFactor(ps.getParameter<double>("etConversionFactor"))
104 {
105  // Build the discretization grid
107  ps.getParameter<edm::ParameterSet>("GridConfiguration"));
108  checkConfig(energyFlow, "invalid discretization grid");
109 
110  // Copy of the grid which will be used for convolutions
111  convolvedFlow = std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >(
112  new fftjet::Grid2d<fftjetcms::Real>(*energyFlow));
113 
114  // Build the initial set of pattern recognition scales
116  ps.getParameter<edm::ParameterSet>("InitialScales"));
117  checkConfig(iniScales, "invalid set of scales");
118 
119  // Build the FFT engine(s), pattern recognition kernel(s),
120  // and the kernel convolver
122 
123  // Build the set of pattern recognition scales
124  produces<TH3F>(outputLabel);
125 }
T getParameter(std::string const &) const
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
std::auto_ptr< fftjet::Grid2d< Real > > fftjet_Grid2d_parser(const edm::ParameterSet &ps)
void checkConfig(const Ptr &ptr, const char *message)
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > convolvedFlow
std::auto_ptr< std::vector< double > > iniScales
void buildKernelConvolver(const edm::ParameterSet &)
const std::string outputLabel
std::auto_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)
FFTJetEFlowSmoother::~FFTJetEFlowSmoother ( )

Definition at line 209 of file FFTJetEFlowSmoother.cc.

210 {
211 }
FFTJetEFlowSmoother::FFTJetEFlowSmoother ( )
private
FFTJetEFlowSmoother::FFTJetEFlowSmoother ( const FFTJetEFlowSmoother )
private

Member Function Documentation

void FFTJetEFlowSmoother::beginJob ( void  )
protectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 275 of file FFTJetEFlowSmoother.cc.

276 {
277 }
void FFTJetEFlowSmoother::buildKernelConvolver ( const edm::ParameterSet ps)
private

Definition at line 128 of file FFTJetEFlowSmoother.cc.

References anotherEngine, convolver, fftjetcms::FFTJetInterface::energyFlow, engine, etaKernel, edm::hlt::Exception, edm::ParameterSet::getParameter(), kernel2d, M_PI, and phiKernel.

Referenced by FFTJetEFlowSmoother().

129 {
130  // Check the parameter named "etaDependentScaleFactors". If the vector
131  // of scales is empty we will use 2d kernel, otherwise use 1d kernels
132  const std::vector<double> etaDependentScaleFactors(
133  ps.getParameter<std::vector<double> >("etaDependentScaleFactors"));
134 
135  // Make sure that the number of scale factors provided is correct
136  const bool use2dKernel = etaDependentScaleFactors.empty();
137  if (!use2dKernel)
138  if (etaDependentScaleFactors.size() != energyFlow->nEta())
139  throw cms::Exception("FFTJetBadConfig")
140  << "invalid number of eta-dependent scale factors"
141  << std::endl;
142 
143  // Get the eta and phi scales for the kernel(s)
144  double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
145  const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
146  if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
147  throw cms::Exception("FFTJetBadConfig")
148  << "invalid kernel scale" << std::endl;
149 
150  // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
151  // kernel scale in eta to compensate.
152  kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
153 
154  // Are we going to try to fix the efficiency near detector edges?
155  const bool fixEfficiency = ps.getParameter<bool>("fixEfficiency");
156 
157  // Minimum and maximum eta bin for the convolver
158  unsigned convolverMinBin = 0, convolverMaxBin = 0;
159  if (fixEfficiency || !use2dKernel)
160  {
161  convolverMinBin = ps.getParameter<unsigned>("convolverMinBin");
162  convolverMaxBin = ps.getParameter<unsigned>("convolverMaxBin");
163  }
164 
165  if (use2dKernel)
166  {
167  // Build the FFT engine
168  engine = std::auto_ptr<MyFFTEngine>(
169  new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
170 
171  // 2d kernel
172  kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
173  new fftjet::DiscreteGauss2d(
174  kernelEtaScale, kernelPhiScale,
175  energyFlow->nEta(), energyFlow->nPhi()));
176 
177  // 2d convolver
178  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
179  new fftjet::FrequencyKernelConvolver<Real,Complex>(
180  engine.get(), kernel2d.get(),
181  convolverMinBin, convolverMaxBin));
182  }
183  else
184  {
185  // Need separate FFT engines for eta and phi
186  engine = std::auto_ptr<MyFFTEngine>(
187  new MyFFTEngine(1, energyFlow->nEta()));
188  anotherEngine = std::auto_ptr<MyFFTEngine>(
189  new MyFFTEngine(1, energyFlow->nPhi()));
190 
191  // 1d kernels
192  etaKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
193  new fftjet::DiscreteGauss1d(kernelEtaScale, energyFlow->nEta()));
194 
195  phiKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
196  new fftjet::DiscreteGauss1d(kernelPhiScale, energyFlow->nPhi()));
197 
198  // Sequential convolver
199  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
200  new fftjet::FrequencySequentialConvolver<Real,Complex>(
201  engine.get(), anotherEngine.get(),
202  etaKernel.get(), phiKernel.get(),
203  etaDependentScaleFactors, convolverMinBin,
204  convolverMaxBin, fixEfficiency));
205  }
206 }
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::AbsFrequencyKernel1d > phiKernel
std::auto_ptr< MyFFTEngine > engine
std::auto_ptr< fftjet::AbsFrequencyKernel > kernel2d
#define M_PI
Definition: BFit3D.cc:3
fftjet::FFTWDoubleEngine MyFFTEngine
std::auto_ptr< fftjet::AbsConvolverBase< Real > > convolver
void FFTJetEFlowSmoother::endJob ( void  )
protectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 281 of file FFTJetEFlowSmoother.cc.

282 {
283 }
FFTJetEFlowSmoother& FFTJetEFlowSmoother::operator= ( const FFTJetEFlowSmoother )
private
void FFTJetEFlowSmoother::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
protectedvirtual

Implements edm::EDProducer.

Definition at line 215 of file FFTJetEFlowSmoother.cc.

References convolvedFlow, convolver, fftjetcms::FFTJetInterface::discretizeEnergyFlow(), fftjetcms::FFTJetInterface::energyFlow, etConversionFactor, g, fftjetcms::FFTJetInterface::getEventScale(), h, iniScales, fftjetcms::FFTJetInterface::loadInputCollection(), M_PI, fftjetcms::FFTJetInterface::outputLabel, funct::pow(), edm::Event::put(), and scalePower.

217 {
218  loadInputCollection(iEvent);
220 
221  // Various useful variables
222  const fftjet::Grid2d<Real>& g(*energyFlow);
223  const unsigned nScales = iniScales->size();
224  const double* scales = &(*iniScales)[0];
225  Real* convData = const_cast<Real*>(convolvedFlow->data());
226  const unsigned nEta = g.nEta();
227  const unsigned nPhi = g.nPhi();
228  const double bin0edge = g.phiBin0Edge();
229 
230  // We will fill the following histo
231  std::auto_ptr<TH3F> pTable(
232  new TH3F("FFTJetEFlowSmoother", "FFTJetEFlowSmoother",
233  nScales+1U, -1.5, nScales-0.5,
234  nEta, g.etaMin(), g.etaMax(),
235  nPhi, bin0edge, bin0edge+2.0*M_PI));
236  TH3F* h = pTable.get();
237  h->SetDirectory(0);
238  h->GetXaxis()->SetTitle("Scale");
239  h->GetYaxis()->SetTitle("Eta");
240  h->GetZaxis()->SetTitle("Phi");
241 
242  // Fill the original thing
243  double factor = etConversionFactor*pow(getEventScale(), scalePower);
244  for (unsigned ieta=0; ieta<nEta; ++ieta)
245  for (unsigned iphi=0; iphi<nPhi; ++iphi)
246  h->SetBinContent(1U, ieta+1U, iphi+1U,
247  factor*g.binValue(ieta, iphi));
248 
249  // Go over all scales and perform the convolutions
250  convolver->setEventData(g.data(), nEta, nPhi);
251  for (unsigned iscale=0; iscale<nScales; ++iscale)
252  {
253  factor = etConversionFactor*pow(scales[iscale], scalePower);
254 
255  // Perform the convolution
256  convolver->convolveWithKernel(
257  scales[iscale], convData,
258  convolvedFlow->nEta(), convolvedFlow->nPhi());
259 
260  // Fill the output histo
261  for (unsigned ieta=0; ieta<nEta; ++ieta)
262  {
263  const Real* etaData = convData + ieta*nPhi;
264  for (unsigned iphi=0; iphi<nPhi; ++iphi)
265  h->SetBinContent(iscale+2U, ieta+1U, iphi+1U,
266  factor*etaData[iphi]);
267  }
268  }
269 
270  iEvent.put(pTable, outputLabel);
271 }
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
void loadInputCollection(const edm::Event &)
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
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
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
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

Member Data Documentation

std::auto_ptr<MyFFTEngine> FFTJetEFlowSmoother::anotherEngine
private

Definition at line 80 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> > FFTJetEFlowSmoother::convolvedFlow
private

Definition at line 73 of file FFTJetEFlowSmoother.cc.

Referenced by FFTJetEFlowSmoother(), and produce().

std::auto_ptr<fftjet::AbsConvolverBase<Real> > FFTJetEFlowSmoother::convolver
private

Definition at line 88 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver(), and produce().

std::auto_ptr<MyFFTEngine> FFTJetEFlowSmoother::engine
private

Definition at line 79 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

std::auto_ptr<fftjet::AbsFrequencyKernel1d> FFTJetEFlowSmoother::etaKernel
private

Definition at line 84 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

double FFTJetEFlowSmoother::etConversionFactor
private

Definition at line 94 of file FFTJetEFlowSmoother.cc.

Referenced by produce().

std::auto_ptr<std::vector<double> > FFTJetEFlowSmoother::iniScales
private

Definition at line 76 of file FFTJetEFlowSmoother.cc.

Referenced by FFTJetEFlowSmoother(), and produce().

std::auto_ptr<fftjet::AbsFrequencyKernel> FFTJetEFlowSmoother::kernel2d
private

Definition at line 83 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

std::auto_ptr<fftjet::AbsFrequencyKernel1d> FFTJetEFlowSmoother::phiKernel
private

Definition at line 85 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

double FFTJetEFlowSmoother::scalePower
private

Definition at line 91 of file FFTJetEFlowSmoother.cc.

Referenced by produce().