CMS 3D CMS Logo

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 ()=delete
 
 FFTJetEFlowSmoother (const FFTJetEFlowSmoother &)=delete
 
FFTJetEFlowSmootheroperator= (const FFTJetEFlowSmoother &)=delete
 
 ~FFTJetEFlowSmoother () override
 
- Public Member Functions inherited from fftjetcms::FFTJetInterface
 FFTJetInterface ()=delete
 
 FFTJetInterface (const FFTJetInterface &)=delete
 
FFTJetInterfaceoperator= (const FFTJetInterface &)=delete
 
 ~FFTJetInterface () override
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDProducer () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
std::vector< bool > const & recordProvenanceList () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
TypeLabelList const & typeLabelList () const
 used by the fwk to register the list of products of this module More...
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESProxyIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const *> const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

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::ProducerBase
template<Transition Tr = Transition::Event>
auto produces (std::string instanceName) noexcept
 declare what type of product will make and with which optional label More...
 
template<Transition B>
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
template<BranchType B>
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
template<typename ProductType , Transition B>
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<class ProductType >
BranchAliasSetterT< ProductType > produces ()
 
template<typename ProductType , BranchType B>
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<typename ProductType , BranchType B>
BranchAliasSetterT< ProductType > produces ()
 
template<class ProductType >
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<typename ProductType , Transition B>
BranchAliasSetterT< ProductType > produces ()
 
template<Transition Tr = Transition::Event>
auto produces () noexcept
 
ProducesCollector producesCollector ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

Private Member Functions

void buildKernelConvolver (const edm::ParameterSet &)
 

Private Attributes

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

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex > >
 
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)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsInputProcessBlocks ()
 
static bool wantsProcessBlocks ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Attributes inherited from fftjetcms::FFTJetInterface
const AnomalousTower anomalous
 
std::vector< unsigned > candidateIndex
 
const bool doPVCorrection
 
std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
 
const std::vector< double > etaDependentMagnutideFactors
 
std::vector< fftjetcms::VectorLikeeventData
 
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 51 of file FFTJetEFlowSmoother.cc.

Constructor & Destructor Documentation

◆ FFTJetEFlowSmoother() [1/3]

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

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

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 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void checkConfig(const Ptr &ptr, const char *message)
std::unique_ptr< fftjet::Grid2d< Real > > fftjet_Grid2d_parser(const edm::ParameterSet &ps)
std::unique_ptr< std::vector< double > > iniScales
std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > convolvedFlow
std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
void buildKernelConvolver(const edm::ParameterSet &)
std::unique_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)
const std::string outputLabel

◆ FFTJetEFlowSmoother() [2/3]

FFTJetEFlowSmoother::FFTJetEFlowSmoother ( )
delete

◆ FFTJetEFlowSmoother() [3/3]

FFTJetEFlowSmoother::FFTJetEFlowSmoother ( const FFTJetEFlowSmoother )
delete

◆ ~FFTJetEFlowSmoother()

FFTJetEFlowSmoother::~FFTJetEFlowSmoother ( )
override

Definition at line 186 of file FFTJetEFlowSmoother.cc.

186 {}

Member Function Documentation

◆ beginJob()

void FFTJetEFlowSmoother::beginJob ( void  )
overrideprotectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 246 of file FFTJetEFlowSmoother.cc.

246 {}

◆ buildKernelConvolver()

void FFTJetEFlowSmoother::buildKernelConvolver ( const edm::ParameterSet ps)
private

Definition at line 119 of file FFTJetEFlowSmoother.cc.

References anotherEngine, convolver, fftjeteflowsmoother_cfi::convolverMaxBin, fftjeteflowsmoother_cfi::convolverMinBin, fftjetcms::FFTJetInterface::energyFlow, engine, fftjeteflowsmoother_cfi::etaDependentScaleFactors, etaKernel, Exception, fftjeteflowsmoother_cfi::fixEfficiency, edm::ParameterSet::getParameter(), kernel2d, fftjeteflowsmoother_cfi::kernelEtaScale, fftjeteflowsmoother_cfi::kernelPhiScale, M_PI, and phiKernel.

Referenced by FFTJetEFlowSmoother().

119  {
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>(
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(),
182  fixEfficiency));
183  }
184 }
std::unique_ptr< MyFFTEngine > engine
std::unique_ptr< fftjet::AbsFrequencyKernel1d > etaKernel
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::unique_ptr< fftjet::AbsConvolverBase< Real > > convolver
std::unique_ptr< MyFFTEngine > anotherEngine
std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
#define M_PI
std::unique_ptr< fftjet::AbsFrequencyKernel > kernel2d
std::unique_ptr< fftjet::AbsFrequencyKernel1d > phiKernel

◆ endJob()

void FFTJetEFlowSmoother::endJob ( void  )
overrideprotectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 249 of file FFTJetEFlowSmoother.cc.

249 {}

◆ operator=()

FFTJetEFlowSmoother& FFTJetEFlowSmoother::operator= ( const FFTJetEFlowSmoother )
delete

◆ produce()

void FFTJetEFlowSmoother::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprotectedvirtual

Implements edm::EDProducer.

Definition at line 189 of file FFTJetEFlowSmoother.cc.

References convolvedFlow, convolver, fftjetcms::FFTJetInterface::discretizeEnergyFlow(), fftjetcms::FFTJetInterface::energyFlow, etConversionFactor, DQMScaleToClient_cfi::factor, g, fftjetcms::FFTJetInterface::getEventScale(), h, LEDCalibrationChannels::ieta, iEvent, iniScales, LEDCalibrationChannels::iphi, fftjetcms::FFTJetInterface::loadInputCollection(), M_PI, eostools::move(), HLT_2022v12_cff::nEta, HLT_2022v12_cff::nPhi, fftjetcommon_cfi::nScales, fftjetcms::FFTJetInterface::outputLabel, funct::pow(), scalePower, and mitigatedMETSequence_cff::U.

189  {
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
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 }
std::unique_ptr< fftjet::AbsConvolverBase< Real > > convolver
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
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
#define M_PI
double Real
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
const std::string outputLabel
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ anotherEngine

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

Definition at line 76 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

◆ convolvedFlow

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

Definition at line 69 of file FFTJetEFlowSmoother.cc.

Referenced by FFTJetEFlowSmoother(), and produce().

◆ convolver

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

Definition at line 84 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver(), and produce().

◆ engine

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

Definition at line 75 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

◆ etaKernel

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

Definition at line 80 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

◆ etConversionFactor

double FFTJetEFlowSmoother::etConversionFactor
private

Definition at line 90 of file FFTJetEFlowSmoother.cc.

Referenced by produce().

◆ iniScales

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

Definition at line 72 of file FFTJetEFlowSmoother.cc.

Referenced by FFTJetEFlowSmoother(), and produce().

◆ kernel2d

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

Definition at line 79 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

◆ phiKernel

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

Definition at line 81 of file FFTJetEFlowSmoother.cc.

Referenced by buildKernelConvolver().

◆ scalePower

double FFTJetEFlowSmoother::scalePower
private

Definition at line 87 of file FFTJetEFlowSmoother.cc.

Referenced by produce().