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:
fftjetcms::FFTJetInterface edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 FFTJetEFlowSmoother (const edm::ParameterSet &)
 
 ~FFTJetEFlowSmoother ()
 
- Public Member Functions inherited from fftjetcms::FFTJetInterface
virtual ~FFTJetInterface ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 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
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Protected Member Functions

void beginJob () override
 
void endJob () override
 
void produce (edm::Event &, const edm::EventSetup &) override
 
- 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
 
- 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)
 

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
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- 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 50 of file FFTJetEFlowSmoother.cc.

Constructor & Destructor Documentation

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

Definition at line 97 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.

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 }
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 206 of file FFTJetEFlowSmoother.cc.

207 {
208 }
FFTJetEFlowSmoother::FFTJetEFlowSmoother ( )
private
FFTJetEFlowSmoother::FFTJetEFlowSmoother ( const FFTJetEFlowSmoother )
private

Member Function Documentation

void FFTJetEFlowSmoother::beginJob ( void  )
overrideprotectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 272 of file FFTJetEFlowSmoother.cc.

273 {
274 }
void FFTJetEFlowSmoother::buildKernelConvolver ( const edm::ParameterSet ps)
private

Definition at line 125 of file FFTJetEFlowSmoother.cc.

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

Referenced by FFTJetEFlowSmoother().

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 }
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
fftjet::FFTWDoubleEngine MyFFTEngine
std::auto_ptr< fftjet::AbsConvolverBase< Real > > convolver
void FFTJetEFlowSmoother::endJob ( void  )
overrideprotectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 278 of file FFTJetEFlowSmoother.cc.

279 {
280 }
FFTJetEFlowSmoother& FFTJetEFlowSmoother::operator= ( const FFTJetEFlowSmoother )
private
void FFTJetEFlowSmoother::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprotectedvirtual

Implements edm::EDProducer.

Definition at line 212 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, HLT_25ns14e33_v1_cff::nEta, HLT_25ns14e33_v1_cff::nPhi, fftjetcms::FFTJetInterface::outputLabel, funct::pow(), edm::Event::put(), and scalePower.

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 }
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:120
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::Grid2d< fftjetcms::Real > > convolvedFlow
std::auto_ptr< std::vector< double > > iniScales
#define M_PI
double Real
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 77 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

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

Definition at line 70 of file FFTJetEFlowSmoother.cc.

Referenced by FFTJetEFlowSmoother(), and produce().

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

Definition at line 85 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver(), and produce().

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

Definition at line 76 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

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

Definition at line 81 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

double FFTJetEFlowSmoother::etConversionFactor
private

Definition at line 91 of file FFTJetEFlowSmoother.cc.

Referenced by produce().

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

Definition at line 73 of file FFTJetEFlowSmoother.cc.

Referenced by FFTJetEFlowSmoother(), and produce().

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

Definition at line 80 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

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

Definition at line 82 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

double FFTJetEFlowSmoother::scalePower
private

Definition at line 88 of file FFTJetEFlowSmoother.cc.

Referenced by produce().