CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/RecoJets/FFTJetProducers/plugins/FFTJetPileupProcessor.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    FFTJetProducers
00004 // Class:      FFTJetPileupProcessor
00005 // 
00013 //
00014 // Original Author:  Igor Volobouev
00015 //         Created:  Wed Apr 20 13:52:23 CDT 2011
00016 // $Id: FFTJetPileupProcessor.cc,v 1.8 2011/07/05 08:24:15 igv Exp $
00017 //
00018 //
00019 
00020 #include <cmath>
00021 #include <fstream>
00022 
00023 // FFTJet headers
00024 #include "fftjet/FrequencyKernelConvolver.hh"
00025 #include "fftjet/DiscreteGauss2d.hh"
00026 #include "fftjet/EquidistantSequence.hh"
00027 #include "fftjet/interpolate.hh"
00028 
00029 // framework include files
00030 #include "FWCore/Framework/interface/Frameworkfwd.h"
00031 #include "FWCore/Framework/interface/EDProducer.h"
00032 #include "FWCore/Framework/interface/MakerMacros.h"
00033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00034 
00035 #include "DataFormats/Common/interface/View.h"
00036 #include "DataFormats/Common/interface/Handle.h"
00037 #include "DataFormats/Histograms/interface/MEtoEDMFormat.h"
00038 
00039 #include "RecoJets/FFTJetAlgorithms/interface/gridConverters.h"
00040 
00041 // parameter parser header
00042 #include "RecoJets/FFTJetProducers/interface/FFTJetParameterParser.h"
00043 
00044 // useful utilities collected in the second base
00045 #include "RecoJets/FFTJetProducers/interface/FFTJetInterface.h"
00046 
00047 using namespace fftjetcms;
00048 
00049 //
00050 // class declaration
00051 //
00052 class FFTJetPileupProcessor : public edm::EDProducer, public FFTJetInterface
00053 {
00054 public:
00055     explicit FFTJetPileupProcessor(const edm::ParameterSet&);
00056     ~FFTJetPileupProcessor();
00057 
00058 protected:
00059     // methods
00060     void beginJob() ;
00061     void produce(edm::Event&, const edm::EventSetup&);
00062     void endJob() ;
00063 
00064 private:
00065     FFTJetPileupProcessor();
00066     FFTJetPileupProcessor(const FFTJetPileupProcessor&);
00067     FFTJetPileupProcessor& operator=(const FFTJetPileupProcessor&);
00068 
00069     void buildKernelConvolver(const edm::ParameterSet&);
00070     void mixExtraGrid();
00071 
00072     // The FFT engine
00073     std::auto_ptr<MyFFTEngine> engine;
00074 
00075     // The pattern recognition kernel(s)
00076     std::auto_ptr<fftjet::AbsFrequencyKernel> kernel2d;
00077 
00078     // The kernel convolver
00079     std::auto_ptr<fftjet::AbsConvolverBase<Real> > convolver;
00080 
00081     // Storage for convolved energy flow
00082     std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> > convolvedFlow;
00083 
00084     // Filtering scales
00085     std::auto_ptr<fftjet::EquidistantInLogSpace> filterScales;
00086 
00087     // Eta-dependent factors to use for flattening the distribution
00088     // _after_ the filtering
00089     std::vector<double> etaFlatteningFactors;
00090 
00091     // Number of percentile points to use
00092     unsigned nPercentiles;
00093 
00094     // Bin range. Both values of 0 means use the whole range.
00095     unsigned convolverMinBin;
00096     unsigned convolverMaxBin;
00097 
00098     // Density conversion factor
00099     double pileupEtaPhiArea;
00100 
00101     // Variable related to mixing additional grids
00102     std::vector<std::string> externalGridFiles;
00103     std::ifstream gridStream;
00104     unsigned currentFileNum;
00105 };
00106 
00107 //
00108 // constructors and destructor
00109 //
00110 FFTJetPileupProcessor::FFTJetPileupProcessor(const edm::ParameterSet& ps)
00111     : FFTJetInterface(ps),
00112       etaFlatteningFactors(
00113           ps.getParameter<std::vector<double> >("etaFlatteningFactors")),
00114       nPercentiles(ps.getParameter<unsigned>("nPercentiles")),
00115       convolverMinBin(ps.getParameter<unsigned>("convolverMinBin")),
00116       convolverMaxBin(ps.getParameter<unsigned>("convolverMaxBin")),
00117       pileupEtaPhiArea(ps.getParameter<double>("pileupEtaPhiArea")),
00118       externalGridFiles(ps.getParameter<std::vector<std::string> >("externalGridFiles")),
00119       currentFileNum(externalGridFiles.size() + 1U)
00120 {
00121     // Build the discretization grid
00122     energyFlow = fftjet_Grid2d_parser(
00123         ps.getParameter<edm::ParameterSet>("GridConfiguration"));
00124     checkConfig(energyFlow, "invalid discretization grid");
00125 
00126     // Copy of the grid which will be used for convolutions
00127     convolvedFlow = std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >(
00128         new fftjet::Grid2d<fftjetcms::Real>(*energyFlow));
00129 
00130     // Make sure the size of flattening factors is appropriate
00131     if (!etaFlatteningFactors.empty())
00132     {
00133         if (etaFlatteningFactors.size() != convolvedFlow->nEta())
00134             throw cms::Exception("FFTJetBadConfig")
00135                 << "ERROR in FFTJetPileupProcessor constructor:"
00136                 " number of elements in the \"etaFlatteningFactors\""
00137                 " vector is inconsistent with the iscretization grid binning"
00138                 << std::endl;
00139     }
00140 
00141     // Build the FFT engine(s), pattern recognition kernel(s),
00142     // and the kernel convolver
00143     buildKernelConvolver(ps);
00144 
00145     // Build the set of pattern recognition scales
00146     const unsigned nScales = ps.getParameter<unsigned>("nScales");
00147     const double minScale = ps.getParameter<double>("minScale");
00148     const double maxScale = ps.getParameter<double>("maxScale");
00149     if (minScale <= 0.0 || maxScale < minScale || nScales == 0U)
00150         throw cms::Exception("FFTJetBadConfig")
00151             << "invalid filter scales" << std::endl;
00152 
00153     filterScales = std::auto_ptr<fftjet::EquidistantInLogSpace>(
00154         new fftjet::EquidistantInLogSpace(minScale, maxScale, nScales));
00155 
00156     produces<TH2D>(outputLabel);
00157     produces<std::pair<double,double> >(outputLabel);
00158 }
00159 
00160 
00161 void FFTJetPileupProcessor::buildKernelConvolver(const edm::ParameterSet& ps)
00162 {
00163     // Get the eta and phi scales for the kernel(s)
00164     double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
00165     const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
00166     if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
00167         throw cms::Exception("FFTJetBadConfig")
00168             << "invalid kernel scales" << std::endl;
00169 
00170     if (convolverMinBin || convolverMaxBin)
00171         if (convolverMinBin >= convolverMaxBin || 
00172             convolverMaxBin > energyFlow->nEta())
00173             throw cms::Exception("FFTJetBadConfig")
00174                 << "invalid convolver bin range" << std::endl;
00175 
00176     // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
00177     // kernel scale in eta to compensate.
00178     kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
00179 
00180     // Build the FFT engine
00181     engine = std::auto_ptr<MyFFTEngine>(
00182         new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
00183 
00184     // 2d kernel
00185     kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
00186         new fftjet::DiscreteGauss2d(
00187             kernelEtaScale, kernelPhiScale,
00188             energyFlow->nEta(), energyFlow->nPhi()));
00189 
00190     // 2d convolver
00191     convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
00192         new fftjet::FrequencyKernelConvolver<Real,Complex>(
00193             engine.get(), kernel2d.get(),
00194             convolverMinBin, convolverMaxBin));
00195 }
00196 
00197 
00198 FFTJetPileupProcessor::~FFTJetPileupProcessor()
00199 {
00200 }
00201 
00202 
00203 // ------------ method called to produce the data  ------------
00204 void FFTJetPileupProcessor::produce(
00205     edm::Event& iEvent, const edm::EventSetup& iSetup)
00206 {
00207     loadInputCollection(iEvent);
00208     discretizeEnergyFlow();
00209 
00210     // Determine the average Et density for this event.
00211     // Needs to be done here, before mixing in another grid.
00212     const fftjet::Grid2d<Real>& g(*energyFlow);
00213     const double densityBeforeMixing = g.sum()/pileupEtaPhiArea;
00214 
00215     // Mix an extra grid (if requested)
00216     double densityAfterMixing = -1.0;
00217     if (!externalGridFiles.empty())
00218     {
00219         mixExtraGrid();
00220         densityAfterMixing = g.sum()/pileupEtaPhiArea;
00221     }
00222 
00223     // Various useful variables
00224     const unsigned nScales = filterScales->size();
00225     const double* scales = &(*filterScales)[0];
00226     Real* convData = const_cast<Real*>(convolvedFlow->data());
00227     Real* sortData = convData + convolverMinBin*convolvedFlow->nPhi();
00228     const unsigned dataLen = convolverMaxBin ? 
00229         (convolverMaxBin - convolverMinBin)*convolvedFlow->nPhi() :
00230         convolvedFlow->nEta()*convolvedFlow->nPhi();
00231 
00232     // We will fill the following histo
00233     std::auto_ptr<TH2D> pTable(new TH2D("FFTJetPileupProcessor",
00234                                         "FFTJetPileupProcessor",
00235                                         nScales, -0.5, nScales-0.5,
00236                                         nPercentiles, 0.0, 1.0));
00237     pTable->SetDirectory(0);
00238     pTable->GetXaxis()->SetTitle("Filter Number");
00239     pTable->GetYaxis()->SetTitle("Et CDF");
00240     pTable->GetZaxis()->SetTitle("Et Density");
00241 
00242     // Go over all scales and perform the convolutions
00243     convolver->setEventData(g.data(), g.nEta(), g.nPhi());
00244     for (unsigned iscale=0; iscale<nScales; ++iscale)
00245     {
00246         // Perform the convolution
00247         convolver->convolveWithKernel(
00248             scales[iscale], convData,
00249             convolvedFlow->nEta(), convolvedFlow->nPhi());
00250 
00251         // Apply the flattening factors
00252         if (!etaFlatteningFactors.empty())
00253             convolvedFlow->scaleData(&etaFlatteningFactors[0],
00254                                      etaFlatteningFactors.size());
00255 
00256         // Sort the convolved data
00257         std::sort(sortData, sortData+dataLen);
00258 
00259         // Determine the percentile points
00260         for (unsigned iper=0; iper<nPercentiles; ++iper)
00261         {
00262             // Map percentile 0 into point number 0, 
00263             // 1 into point number dataLen-1
00264             const double q = (iper + 0.5)/nPercentiles;
00265             const double dindex = q*(dataLen-1U);
00266             const unsigned ilow = static_cast<unsigned>(std::floor(dindex));
00267             const double percentile = fftjet::lin_interpolate_1d(
00268                 ilow, ilow+1U, sortData[ilow], sortData[ilow+1U], dindex);
00269 
00270             // Store the calculated percentile
00271             pTable->SetBinContent(iscale+1U, iper+1U, percentile);
00272         }
00273     }
00274 
00275     iEvent.put(pTable, outputLabel);
00276 
00277     std::auto_ptr<std::pair<double,double> > etSum(
00278         new std::pair<double,double>(densityBeforeMixing, densityAfterMixing));
00279     iEvent.put(etSum, outputLabel);
00280 }
00281 
00282 
00283 void FFTJetPileupProcessor::mixExtraGrid()
00284 {
00285     const unsigned nFiles = externalGridFiles.size();
00286     if (currentFileNum > nFiles)
00287     {
00288         // This is the first time this function is called
00289         currentFileNum = 0;
00290         gridStream.open(externalGridFiles[currentFileNum].c_str(),
00291                         std::ios_base::in | std::ios_base::binary);
00292         if (!gridStream.is_open())
00293             throw cms::Exception("FFTJetBadConfig")
00294                 << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
00295                 " failed to open external grid file "
00296                 << externalGridFiles[currentFileNum] << std::endl;
00297     }
00298 
00299     const fftjet::Grid2d<float>* g = fftjet::Grid2d<float>::read(gridStream);
00300 
00301     // If we can't read the grid, we need to switch to another file
00302     for (unsigned ntries=0; ntries<nFiles && g == 0; ++ntries)
00303     {
00304         gridStream.close();
00305         currentFileNum = (currentFileNum + 1U) % nFiles;
00306         gridStream.open(externalGridFiles[currentFileNum].c_str(),
00307                         std::ios_base::in | std::ios_base::binary);
00308         if (!gridStream.is_open())
00309             throw cms::Exception("FFTJetBadConfig")
00310                 << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
00311                 " failed to open external grid file "
00312                 << externalGridFiles[currentFileNum] << std::endl;
00313         g = fftjet::Grid2d<float>::read(gridStream);
00314     }
00315 
00316     if (g)
00317     {
00318         add_Grid2d_data(energyFlow.get(), *g);
00319         delete g;
00320     }
00321     else
00322     {
00323         // Too bad, no useful file found
00324         throw cms::Exception("FFTJetBadConfig")
00325             << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
00326             " no valid grid records found" << std::endl;
00327     }
00328 }
00329 
00330 
00331 // ------------ method called once each job just before starting event loop
00332 void FFTJetPileupProcessor::beginJob()
00333 {
00334 }
00335 
00336 
00337 // ------------ method called once each job just after ending the event loop
00338 void FFTJetPileupProcessor::endJob()
00339 {
00340 }
00341 
00342 
00343 //define this as a plug-in
00344 DEFINE_FWK_MODULE(FFTJetPileupProcessor);