00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020 #include <cmath>
00021
00022
00023 #include "fftjet/FrequencyKernelConvolver.hh"
00024 #include "fftjet/DiscreteGauss2d.hh"
00025 #include "fftjet/EquidistantSequence.hh"
00026 #include "fftjet/interpolate.hh"
00027 #include "fftjet/FrequencyKernelConvolver.hh"
00028 #include "fftjet/FrequencySequentialConvolver.hh"
00029 #include "fftjet/DiscreteGauss1d.hh"
00030 #include "fftjet/DiscreteGauss2d.hh"
00031
00032
00033 #include "FWCore/Framework/interface/Frameworkfwd.h"
00034 #include "FWCore/Framework/interface/EDProducer.h"
00035 #include "FWCore/Framework/interface/MakerMacros.h"
00036 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00037
00038 #include "DataFormats/Common/interface/View.h"
00039 #include "DataFormats/Common/interface/Handle.h"
00040 #include "DataFormats/Histograms/interface/MEtoEDMFormat.h"
00041
00042
00043 #include "RecoJets/FFTJetProducers/interface/FFTJetParameterParser.h"
00044
00045
00046 #include "RecoJets/FFTJetProducers/interface/FFTJetInterface.h"
00047
00048 using namespace fftjetcms;
00049
00050
00051
00052
00053 class FFTJetEFlowSmoother : public edm::EDProducer, public FFTJetInterface
00054 {
00055 public:
00056 explicit FFTJetEFlowSmoother(const edm::ParameterSet&);
00057 ~FFTJetEFlowSmoother();
00058
00059 protected:
00060
00061 void beginJob() ;
00062 void produce(edm::Event&, const edm::EventSetup&);
00063 void endJob() ;
00064
00065 private:
00066 FFTJetEFlowSmoother();
00067 FFTJetEFlowSmoother(const FFTJetEFlowSmoother&);
00068 FFTJetEFlowSmoother& operator=(const FFTJetEFlowSmoother&);
00069
00070 void buildKernelConvolver(const edm::ParameterSet&);
00071
00072
00073 std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> > convolvedFlow;
00074
00075
00076 std::auto_ptr<std::vector<double> > iniScales;
00077
00078
00079 std::auto_ptr<MyFFTEngine> engine;
00080 std::auto_ptr<MyFFTEngine> anotherEngine;
00081
00082
00083 std::auto_ptr<fftjet::AbsFrequencyKernel> kernel2d;
00084 std::auto_ptr<fftjet::AbsFrequencyKernel1d> etaKernel;
00085 std::auto_ptr<fftjet::AbsFrequencyKernel1d> phiKernel;
00086
00087
00088 std::auto_ptr<fftjet::AbsConvolverBase<Real> > convolver;
00089
00090
00091 double scalePower;
00092
00093
00094 double etConversionFactor;
00095 };
00096
00097
00098
00099
00100 FFTJetEFlowSmoother::FFTJetEFlowSmoother(const edm::ParameterSet& ps)
00101 : FFTJetInterface(ps),
00102 scalePower(ps.getParameter<double>("scalePower")),
00103 etConversionFactor(ps.getParameter<double>("etConversionFactor"))
00104 {
00105
00106 energyFlow = fftjet_Grid2d_parser(
00107 ps.getParameter<edm::ParameterSet>("GridConfiguration"));
00108 checkConfig(energyFlow, "invalid discretization grid");
00109
00110
00111 convolvedFlow = std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >(
00112 new fftjet::Grid2d<fftjetcms::Real>(*energyFlow));
00113
00114
00115 iniScales = fftjet_ScaleSet_parser(
00116 ps.getParameter<edm::ParameterSet>("InitialScales"));
00117 checkConfig(iniScales, "invalid set of scales");
00118
00119
00120
00121 buildKernelConvolver(ps);
00122
00123
00124 produces<TH3F>(outputLabel);
00125 }
00126
00127
00128 void FFTJetEFlowSmoother::buildKernelConvolver(const edm::ParameterSet& ps)
00129 {
00130
00131
00132 const std::vector<double> etaDependentScaleFactors(
00133 ps.getParameter<std::vector<double> >("etaDependentScaleFactors"));
00134
00135
00136 const bool use2dKernel = etaDependentScaleFactors.empty();
00137 if (!use2dKernel)
00138 if (etaDependentScaleFactors.size() != energyFlow->nEta())
00139 throw cms::Exception("FFTJetBadConfig")
00140 << "invalid number of eta-dependent scale factors"
00141 << std::endl;
00142
00143
00144 double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
00145 const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
00146 if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
00147 throw cms::Exception("FFTJetBadConfig")
00148 << "invalid kernel scale" << std::endl;
00149
00150
00151
00152 kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
00153
00154
00155 const bool fixEfficiency = ps.getParameter<bool>("fixEfficiency");
00156
00157
00158 unsigned convolverMinBin = 0, convolverMaxBin = 0;
00159 if (fixEfficiency || !use2dKernel)
00160 {
00161 convolverMinBin = ps.getParameter<unsigned>("convolverMinBin");
00162 convolverMaxBin = ps.getParameter<unsigned>("convolverMaxBin");
00163 }
00164
00165 if (use2dKernel)
00166 {
00167
00168 engine = std::auto_ptr<MyFFTEngine>(
00169 new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
00170
00171
00172 kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
00173 new fftjet::DiscreteGauss2d(
00174 kernelEtaScale, kernelPhiScale,
00175 energyFlow->nEta(), energyFlow->nPhi()));
00176
00177
00178 convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
00179 new fftjet::FrequencyKernelConvolver<Real,Complex>(
00180 engine.get(), kernel2d.get(),
00181 convolverMinBin, convolverMaxBin));
00182 }
00183 else
00184 {
00185
00186 engine = std::auto_ptr<MyFFTEngine>(
00187 new MyFFTEngine(1, energyFlow->nEta()));
00188 anotherEngine = std::auto_ptr<MyFFTEngine>(
00189 new MyFFTEngine(1, energyFlow->nPhi()));
00190
00191
00192 etaKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
00193 new fftjet::DiscreteGauss1d(kernelEtaScale, energyFlow->nEta()));
00194
00195 phiKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
00196 new fftjet::DiscreteGauss1d(kernelPhiScale, energyFlow->nPhi()));
00197
00198
00199 convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
00200 new fftjet::FrequencySequentialConvolver<Real,Complex>(
00201 engine.get(), anotherEngine.get(),
00202 etaKernel.get(), phiKernel.get(),
00203 etaDependentScaleFactors, convolverMinBin,
00204 convolverMaxBin, fixEfficiency));
00205 }
00206 }
00207
00208
00209 FFTJetEFlowSmoother::~FFTJetEFlowSmoother()
00210 {
00211 }
00212
00213
00214
00215 void FFTJetEFlowSmoother::produce(
00216 edm::Event& iEvent, const edm::EventSetup& iSetup)
00217 {
00218 loadInputCollection(iEvent);
00219 discretizeEnergyFlow();
00220
00221
00222 const fftjet::Grid2d<Real>& g(*energyFlow);
00223 const unsigned nScales = iniScales->size();
00224 const double* scales = &(*iniScales)[0];
00225 Real* convData = const_cast<Real*>(convolvedFlow->data());
00226 const unsigned nEta = g.nEta();
00227 const unsigned nPhi = g.nPhi();
00228 const double bin0edge = g.phiBin0Edge();
00229
00230
00231 std::auto_ptr<TH3F> pTable(
00232 new TH3F("FFTJetEFlowSmoother", "FFTJetEFlowSmoother",
00233 nScales+1U, -1.5, nScales-0.5,
00234 nEta, g.etaMin(), g.etaMax(),
00235 nPhi, bin0edge, bin0edge+2.0*M_PI));
00236 TH3F* h = pTable.get();
00237 h->SetDirectory(0);
00238 h->GetXaxis()->SetTitle("Scale");
00239 h->GetYaxis()->SetTitle("Eta");
00240 h->GetZaxis()->SetTitle("Phi");
00241
00242
00243 double factor = etConversionFactor*pow(getEventScale(), scalePower);
00244 for (unsigned ieta=0; ieta<nEta; ++ieta)
00245 for (unsigned iphi=0; iphi<nPhi; ++iphi)
00246 h->SetBinContent(1U, ieta+1U, iphi+1U,
00247 factor*g.binValue(ieta, iphi));
00248
00249
00250 convolver->setEventData(g.data(), nEta, nPhi);
00251 for (unsigned iscale=0; iscale<nScales; ++iscale)
00252 {
00253 factor = etConversionFactor*pow(scales[iscale], scalePower);
00254
00255
00256 convolver->convolveWithKernel(
00257 scales[iscale], convData,
00258 convolvedFlow->nEta(), convolvedFlow->nPhi());
00259
00260
00261 for (unsigned ieta=0; ieta<nEta; ++ieta)
00262 {
00263 const Real* etaData = convData + ieta*nPhi;
00264 for (unsigned iphi=0; iphi<nPhi; ++iphi)
00265 h->SetBinContent(iscale+2U, ieta+1U, iphi+1U,
00266 factor*etaData[iphi]);
00267 }
00268 }
00269
00270 iEvent.put(pTable, outputLabel);
00271 }
00272
00273
00274
00275 void FFTJetEFlowSmoother::beginJob()
00276 {
00277 }
00278
00279
00280
00281 void FFTJetEFlowSmoother::endJob()
00282 {
00283 }
00284
00285
00286
00287 DEFINE_FWK_MODULE(FFTJetEFlowSmoother);