#include <RecoJets/FFTJetProducers/plugins/FFTJetEFlowSmoother.cc>
Public Member Functions | |
FFTJetEFlowSmoother (const edm::ParameterSet &) | |
~FFTJetEFlowSmoother () | |
Protected Member Functions | |
void | beginJob () |
void | endJob () |
void | produce (edm::Event &, const edm::EventSetup &) |
Private Member Functions | |
void | buildKernelConvolver (const edm::ParameterSet &) |
FFTJetEFlowSmoother (const FFTJetEFlowSmoother &) | |
FFTJetEFlowSmoother () | |
FFTJetEFlowSmoother & | operator= (const FFTJetEFlowSmoother &) |
Private Attributes | |
std::auto_ptr< MyFFTEngine > | anotherEngine |
std::auto_ptr< fftjet::Grid2d < fftjetcms::Real > > | convolvedFlow |
std::auto_ptr < fftjet::AbsConvolverBase < Real > > | convolver |
std::auto_ptr< MyFFTEngine > | engine |
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 |
Description: Runs FFTJet filtering code for multiple scales and saves the results
Implementation: [Notes on implementation]
Definition at line 53 of file FFTJetEFlowSmoother.cc.
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.
: FFTJetInterface(ps), scalePower(ps.getParameter<double>("scalePower")), etConversionFactor(ps.getParameter<double>("etConversionFactor")) { // Build the discretization grid energyFlow = fftjet_Grid2d_parser( ps.getParameter<edm::ParameterSet>("GridConfiguration")); checkConfig(energyFlow, "invalid discretization grid"); // Copy of the grid which will be used for convolutions convolvedFlow = std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >( new fftjet::Grid2d<fftjetcms::Real>(*energyFlow)); // Build the initial set of pattern recognition scales iniScales = fftjet_ScaleSet_parser( ps.getParameter<edm::ParameterSet>("InitialScales")); checkConfig(iniScales, "invalid set of scales"); // Build the FFT engine(s), pattern recognition kernel(s), // and the kernel convolver buildKernelConvolver(ps); // Build the set of pattern recognition scales produces<TH3F>(outputLabel); }
FFTJetEFlowSmoother::~FFTJetEFlowSmoother | ( | ) |
Definition at line 209 of file FFTJetEFlowSmoother.cc.
{ }
FFTJetEFlowSmoother::FFTJetEFlowSmoother | ( | ) | [private] |
FFTJetEFlowSmoother::FFTJetEFlowSmoother | ( | const FFTJetEFlowSmoother & | ) | [private] |
void FFTJetEFlowSmoother::beginJob | ( | void | ) | [protected, virtual] |
void FFTJetEFlowSmoother::buildKernelConvolver | ( | const edm::ParameterSet & | ps | ) | [private] |
Definition at line 128 of file FFTJetEFlowSmoother.cc.
References anotherEngine, convolver, fftjetcms::FFTJetInterface::energyFlow, engine, etaKernel, Exception, edm::ParameterSet::getParameter(), kernel2d, M_PI, and phiKernel.
Referenced by FFTJetEFlowSmoother().
{ // Check the parameter named "etaDependentScaleFactors". If the vector // of scales is empty we will use 2d kernel, otherwise use 1d kernels const std::vector<double> etaDependentScaleFactors( ps.getParameter<std::vector<double> >("etaDependentScaleFactors")); // Make sure that the number of scale factors provided is correct const bool use2dKernel = etaDependentScaleFactors.empty(); if (!use2dKernel) if (etaDependentScaleFactors.size() != energyFlow->nEta()) throw cms::Exception("FFTJetBadConfig") << "invalid number of eta-dependent scale factors" << std::endl; // Get the eta and phi scales for the kernel(s) double kernelEtaScale = ps.getParameter<double>("kernelEtaScale"); const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale"); if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0) throw cms::Exception("FFTJetBadConfig") << "invalid kernel scale" << std::endl; // FFT assumes that the grid extent in eta is 2*Pi. Adjust the // kernel scale in eta to compensate. kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin())); // Are we going to try to fix the efficiency near detector edges? const bool fixEfficiency = ps.getParameter<bool>("fixEfficiency"); // Minimum and maximum eta bin for the convolver unsigned convolverMinBin = 0, convolverMaxBin = 0; if (fixEfficiency || !use2dKernel) { convolverMinBin = ps.getParameter<unsigned>("convolverMinBin"); convolverMaxBin = ps.getParameter<unsigned>("convolverMaxBin"); } if (use2dKernel) { // Build the FFT engine engine = std::auto_ptr<MyFFTEngine>( new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi())); // 2d kernel kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>( new fftjet::DiscreteGauss2d( kernelEtaScale, kernelPhiScale, energyFlow->nEta(), energyFlow->nPhi())); // 2d convolver convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >( new fftjet::FrequencyKernelConvolver<Real,Complex>( engine.get(), kernel2d.get(), convolverMinBin, convolverMaxBin)); } else { // Need separate FFT engines for eta and phi engine = std::auto_ptr<MyFFTEngine>( new MyFFTEngine(1, energyFlow->nEta())); anotherEngine = std::auto_ptr<MyFFTEngine>( new MyFFTEngine(1, energyFlow->nPhi())); // 1d kernels etaKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>( new fftjet::DiscreteGauss1d(kernelEtaScale, energyFlow->nEta())); phiKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>( new fftjet::DiscreteGauss1d(kernelPhiScale, energyFlow->nPhi())); // Sequential convolver convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >( new fftjet::FrequencySequentialConvolver<Real,Complex>( engine.get(), anotherEngine.get(), etaKernel.get(), phiKernel.get(), etaDependentScaleFactors, convolverMinBin, convolverMaxBin, fixEfficiency)); } }
void FFTJetEFlowSmoother::endJob | ( | void | ) | [protected, virtual] |
FFTJetEFlowSmoother& FFTJetEFlowSmoother::operator= | ( | const FFTJetEFlowSmoother & | ) | [private] |
void FFTJetEFlowSmoother::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [protected, virtual] |
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.
{ loadInputCollection(iEvent); discretizeEnergyFlow(); // Various useful variables const fftjet::Grid2d<Real>& g(*energyFlow); const unsigned nScales = iniScales->size(); const double* scales = &(*iniScales)[0]; Real* convData = const_cast<Real*>(convolvedFlow->data()); const unsigned nEta = g.nEta(); const unsigned nPhi = g.nPhi(); const double bin0edge = g.phiBin0Edge(); // We will fill the following histo std::auto_ptr<TH3F> pTable( new TH3F("FFTJetEFlowSmoother", "FFTJetEFlowSmoother", nScales+1U, -1.5, nScales-0.5, nEta, g.etaMin(), g.etaMax(), nPhi, bin0edge, bin0edge+2.0*M_PI)); TH3F* h = pTable.get(); h->SetDirectory(0); h->GetXaxis()->SetTitle("Scale"); h->GetYaxis()->SetTitle("Eta"); h->GetZaxis()->SetTitle("Phi"); // Fill the original thing double factor = etConversionFactor*pow(getEventScale(), scalePower); for (unsigned ieta=0; ieta<nEta; ++ieta) for (unsigned iphi=0; iphi<nPhi; ++iphi) h->SetBinContent(1U, ieta+1U, iphi+1U, factor*g.binValue(ieta, iphi)); // Go over all scales and perform the convolutions convolver->setEventData(g.data(), nEta, nPhi); for (unsigned iscale=0; iscale<nScales; ++iscale) { factor = etConversionFactor*pow(scales[iscale], scalePower); // Perform the convolution convolver->convolveWithKernel( scales[iscale], convData, convolvedFlow->nEta(), convolvedFlow->nPhi()); // Fill the output histo for (unsigned ieta=0; ieta<nEta; ++ieta) { const Real* etaData = convData + ieta*nPhi; for (unsigned iphi=0; iphi<nPhi; ++iphi) h->SetBinContent(iscale+2U, ieta+1U, iphi+1U, factor*etaData[iphi]); } } iEvent.put(pTable, outputLabel); }
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().