CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/RecoJets/FFTJetProducers/plugins/FFTJetPatRecoProducer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    FFTJetProducers
00004 // Class:      FFTJetPatRecoProducer
00005 // 
00013 //
00014 // Original Author:  Igor Volobouev
00015 //         Created:  Tue Jun 15 12:45:45 CDT 2010
00016 // $Id: FFTJetPatRecoProducer.cc,v 1.5 2011/07/18 17:08:24 igv Exp $
00017 //
00018 //
00019 
00020 #include <fstream>
00021 
00022 // FFTJet headers
00023 #include "fftjet/ProximityClusteringTree.hh"
00024 #include "fftjet/ClusteringSequencer.hh"
00025 #include "fftjet/ClusteringTreeSparsifier.hh"
00026 #include "fftjet/FrequencyKernelConvolver.hh"
00027 #include "fftjet/FrequencySequentialConvolver.hh"
00028 #include "fftjet/DiscreteGauss1d.hh"
00029 #include "fftjet/DiscreteGauss2d.hh"
00030 
00031 // framework include files
00032 #include "FWCore/Framework/interface/Frameworkfwd.h"
00033 #include "FWCore/Framework/interface/EDProducer.h"
00034 #include "FWCore/Framework/interface/MakerMacros.h"
00035 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00036 
00037 #include "DataFormats/Common/interface/View.h"
00038 #include "DataFormats/Common/interface/Handle.h"
00039 #include "DataFormats/VertexReco/interface/Vertex.h"
00040 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00041 #include "DataFormats/JetReco/interface/GenJetCollection.h"
00042 #include "DataFormats/JetReco/interface/PFJetCollection.h"
00043 #include "DataFormats/JetReco/interface/BasicJetCollection.h"
00044 #include "DataFormats/JetReco/interface/TrackJetCollection.h"
00045 #include "DataFormats/Candidate/interface/Candidate.h"
00046 
00047 // Energy flow object
00048 #include "DataFormats/JetReco/interface/DiscretizedEnergyFlow.h"
00049 
00050 // parameter parser header
00051 #include "RecoJets/FFTJetProducers/interface/FFTJetParameterParser.h"
00052 
00053 // functions which manipulate storable trees
00054 #include "RecoJets/FFTJetAlgorithms/interface/clusteringTreeConverters.h"
00055 
00056 // functions which manipulate energy discretization grids
00057 #include "RecoJets/FFTJetAlgorithms/interface/gridConverters.h"
00058 
00059 // useful utilities collected in the second base
00060 #include "RecoJets/FFTJetProducers/interface/FFTJetInterface.h"
00061 
00062 using namespace fftjetcms;
00063 
00064 //
00065 // class declaration
00066 //
00067 class FFTJetPatRecoProducer : public edm::EDProducer, public FFTJetInterface
00068 {
00069 public:
00070     explicit FFTJetPatRecoProducer(const edm::ParameterSet&);
00071     ~FFTJetPatRecoProducer();
00072 
00073 protected:
00074     // Useful local typedefs
00075     typedef fftjet::ProximityClusteringTree<fftjet::Peak,long> ClusteringTree;
00076     typedef fftjet::SparseClusteringTree<fftjet::Peak,long> SparseTree;
00077     typedef fftjet::ClusteringSequencer<Real> Sequencer;
00078     typedef fftjet::ClusteringTreeSparsifier<fftjet::Peak,long> Sparsifier;
00079 
00080     // methods
00081     void beginJob() ;
00082     void produce(edm::Event&, const edm::EventSetup&);
00083     void endJob() ;
00084 
00085     void buildKernelConvolver(const edm::ParameterSet&);
00086     fftjet::PeakFinder buildPeakFinder(const edm::ParameterSet&);
00087 
00088     template<class Real>
00089     void buildSparseProduct(edm::Event&) const;
00090 
00091     template<class Real>
00092     void buildDenseProduct(edm::Event&) const;
00093 
00094     // The complete clustering tree
00095     ClusteringTree* clusteringTree;
00096 
00097     // Basically, we need to create FFTJet objects
00098     // ClusteringSequencer and ClusteringTreeSparsifier
00099     // which will subsequently perform most of the work
00100     std::auto_ptr<Sequencer> sequencer;
00101     std::auto_ptr<Sparsifier> sparsifier;
00102 
00103     // The FFT engine(s)
00104     std::auto_ptr<MyFFTEngine> engine;
00105     std::auto_ptr<MyFFTEngine> anotherEngine;
00106 
00107     // The pattern recognition kernel(s)
00108     std::auto_ptr<fftjet::AbsFrequencyKernel> kernel2d;
00109     std::auto_ptr<fftjet::AbsFrequencyKernel1d> etaKernel;
00110     std::auto_ptr<fftjet::AbsFrequencyKernel1d> phiKernel;
00111 
00112     // The kernel convolver
00113     std::auto_ptr<fftjet::AbsConvolverBase<Real> > convolver;
00114 
00115     // The peak selector for the clustering tree
00116     std::auto_ptr<fftjet::Functor1<bool,fftjet::Peak> > peakSelector;
00117 
00118     // Distance calculator for the clustering tree
00119     std::auto_ptr<fftjet::AbsDistanceCalculator<fftjet::Peak> > distanceCalc;
00120 
00121     // The sparse clustering tree
00122     SparseTree sparseTree;
00123 
00124     // The following parameters will define the behavior
00125     // of the algorithm wrt insertion of the complete event
00126     // into the clustering tree
00127     const double completeEventDataCutoff;
00128 
00129     // Are we going to make clustering trees?
00130     const bool makeClusteringTree;
00131 
00132     // Are we going to verify the data conversion for double precision
00133     // storage?
00134     const bool verifyDataConversion;
00135 
00136     // Are we going to store the discretization grid?
00137     const bool storeDiscretizationGrid;
00138 
00139     // Sparsify the clustering tree?
00140     const bool sparsify;
00141 
00142 private:
00143     FFTJetPatRecoProducer();
00144     FFTJetPatRecoProducer(const FFTJetPatRecoProducer&);
00145     FFTJetPatRecoProducer& operator=(const FFTJetPatRecoProducer&);
00146 
00147     // Members needed for storing grids externally
00148     std::ofstream externalGridStream;
00149     bool storeGridsExternally;
00150     fftjet::Grid2d<float>* extGrid;
00151 };
00152 
00153 //
00154 // constructors and destructor
00155 //
00156 FFTJetPatRecoProducer::FFTJetPatRecoProducer(const edm::ParameterSet& ps)
00157     : FFTJetInterface(ps),
00158       clusteringTree(0),
00159       completeEventDataCutoff(ps.getParameter<double>("completeEventDataCutoff")),
00160       makeClusteringTree(ps.getParameter<bool>("makeClusteringTree")),
00161       verifyDataConversion(ps.getUntrackedParameter<bool>("verifyDataConversion",false)),
00162       storeDiscretizationGrid(ps.getParameter<bool>("storeDiscretizationGrid")),
00163       sparsify(ps.getParameter<bool>("sparsify")),
00164       extGrid(0)
00165 {
00166     // register your products
00167     if (makeClusteringTree)
00168     {
00169         if (storeInSinglePrecision())
00170             produces<reco::PattRecoTree<float,reco::PattRecoPeak<float> > >(outputLabel);
00171         else
00172             produces<reco::PattRecoTree<double,reco::PattRecoPeak<double> > >(outputLabel);
00173     }
00174     if (storeDiscretizationGrid)
00175         produces<reco::DiscretizedEnergyFlow>(outputLabel);
00176 
00177     // Check if we want to write the grids into an external file
00178     const std::string externalGridFile(ps.getParameter<std::string>("externalGridFile"));
00179     storeGridsExternally = externalGridFile.size() > 0;
00180     if (storeGridsExternally)
00181     {
00182         externalGridStream.open(externalGridFile.c_str(), std::ios_base::out | 
00183                                                           std::ios_base::binary);
00184         if (!externalGridStream.is_open())
00185             throw cms::Exception("FFTJetBadConfig")
00186                 << "FFTJetPatRecoProducer failed to open file "
00187                 << externalGridFile << std::endl;
00188     }
00189 
00190     if (!makeClusteringTree && !storeDiscretizationGrid && !storeGridsExternally)
00191     {
00192         throw cms::Exception("FFTJetBadConfig")
00193             << "FFTJetPatRecoProducer is not configured to produce anything"
00194             << std::endl;
00195     }
00196 
00197     // now do whatever other initialization is needed
00198 
00199     // Build the discretization grid
00200     energyFlow = fftjet_Grid2d_parser(
00201         ps.getParameter<edm::ParameterSet>("GridConfiguration"));
00202     checkConfig(energyFlow, "invalid discretization grid");
00203 
00204     // Build the FFT engine(s), pattern recognition kernel(s),
00205     // and the kernel convolver
00206     buildKernelConvolver(ps);
00207 
00208     // Build the peak selector
00209     peakSelector = fftjet_PeakSelector_parser(
00210         ps.getParameter<edm::ParameterSet>("PeakSelectorConfiguration"));
00211     checkConfig(peakSelector, "invalid peak selector");
00212 
00213     // Build the initial set of pattern recognition scales
00214     std::auto_ptr<std::vector<double> > iniScales = fftjet_ScaleSet_parser(
00215         ps.getParameter<edm::ParameterSet>("InitialScales"));
00216     checkConfig(iniScales, "invalid set of scales");
00217 
00218     // Do we want to use the adaptive clustering tree algorithm?
00219     const unsigned maxAdaptiveScales = 
00220         ps.getParameter<unsigned>("maxAdaptiveScales");
00221     const double minAdaptiveRatioLog = 
00222         ps.getParameter<double>("minAdaptiveRatioLog");
00223     if (minAdaptiveRatioLog <= 0.0)
00224         throw cms::Exception("FFTJetBadConfig")
00225             << "bad adaptive ratio logarithm limit" << std::endl;
00226 
00227     // Make sure that all standard scales are larger than the
00228     // complete event scale
00229     if (getEventScale() > 0.0)
00230     {
00231       const double cs = getEventScale();
00232       const unsigned nscales = iniScales->size();
00233       for (unsigned i=0; i<nscales; ++i)
00234         if (cs >= (*iniScales)[i])
00235           throw cms::Exception("FFTJetBadConfig")
00236             << "incompatible scale for complete event" << std::endl;
00237     }
00238 
00239     // At this point we are ready to build the clustering sequencer
00240     sequencer = std::auto_ptr<Sequencer>(new Sequencer(
00241         convolver.get(), peakSelector.get(), buildPeakFinder(ps),
00242         *iniScales, maxAdaptiveScales, minAdaptiveRatioLog));
00243 
00244     // Build the clustering tree sparsifier
00245     const edm::ParameterSet& SparsifierConfiguration(
00246         ps.getParameter<edm::ParameterSet>("SparsifierConfiguration"));
00247     sparsifier = fftjet_ClusteringTreeSparsifier_parser(
00248         SparsifierConfiguration);
00249     checkConfig(sparsifier, "invalid sparsifier parameters");
00250 
00251     // Build the distance calculator for the clustering tree
00252     const edm::ParameterSet& TreeDistanceCalculator(
00253         ps.getParameter<edm::ParameterSet>("TreeDistanceCalculator"));
00254     distanceCalc = fftjet_DistanceCalculator_parser(TreeDistanceCalculator);
00255     checkConfig(distanceCalc, "invalid tree distance calculator");
00256 
00257     // Build the clustering tree itself
00258     clusteringTree = new ClusteringTree(distanceCalc.get());
00259 }
00260 
00261 
00262 void FFTJetPatRecoProducer::buildKernelConvolver(const edm::ParameterSet& ps)
00263 {
00264     // Check the parameter named "etaDependentScaleFactors". If the vector
00265     // of scales is empty we will use 2d kernel, otherwise use 1d kernels
00266     const std::vector<double> etaDependentScaleFactors(
00267         ps.getParameter<std::vector<double> >("etaDependentScaleFactors"));
00268 
00269     // Make sure that the number of scale factors provided is correct
00270     const bool use2dKernel = etaDependentScaleFactors.empty();
00271     if (!use2dKernel)
00272         if (etaDependentScaleFactors.size() != energyFlow->nEta())
00273             throw cms::Exception("FFTJetBadConfig")
00274                 << "invalid number of eta-dependent scale factors"
00275                 << std::endl;
00276 
00277     // Get the eta and phi scales for the kernel(s)
00278     double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
00279     const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
00280     if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
00281         throw cms::Exception("FFTJetBadConfig")
00282             << "invalid kernel scale" << std::endl;
00283 
00284     // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
00285     // kernel scale in eta to compensate.
00286     kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
00287 
00288     // Are we going to try to fix the efficiency near detector edges?
00289     const bool fixEfficiency = ps.getParameter<bool>("fixEfficiency");
00290 
00291     // Minimum and maximum eta bin for the convolver
00292     unsigned convolverMinBin = 0, convolverMaxBin = 0;
00293     if (fixEfficiency || !use2dKernel)
00294     {
00295         convolverMinBin = ps.getParameter<unsigned>("convolverMinBin");
00296         convolverMaxBin = ps.getParameter<unsigned>("convolverMaxBin");
00297     }
00298 
00299     if (use2dKernel)
00300     {
00301         // Build the FFT engine
00302         engine = std::auto_ptr<MyFFTEngine>(
00303             new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
00304 
00305         // 2d kernel
00306         kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
00307             new fftjet::DiscreteGauss2d(
00308                 kernelEtaScale, kernelPhiScale,
00309                 energyFlow->nEta(), energyFlow->nPhi()));
00310 
00311         // 2d convolver
00312         convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
00313             new fftjet::FrequencyKernelConvolver<Real,Complex>(
00314                 engine.get(), kernel2d.get(),
00315                 convolverMinBin, convolverMaxBin));
00316     }
00317     else
00318     {
00319         // Need separate FFT engines for eta and phi
00320         engine = std::auto_ptr<MyFFTEngine>(
00321             new MyFFTEngine(1, energyFlow->nEta()));
00322         anotherEngine = std::auto_ptr<MyFFTEngine>(
00323             new MyFFTEngine(1, energyFlow->nPhi()));
00324 
00325         // 1d kernels
00326         etaKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
00327             new fftjet::DiscreteGauss1d(kernelEtaScale, energyFlow->nEta()));
00328 
00329         phiKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
00330             new fftjet::DiscreteGauss1d(kernelPhiScale, energyFlow->nPhi()));
00331 
00332         // Sequential convolver
00333         convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
00334             new fftjet::FrequencySequentialConvolver<Real,Complex>(
00335                 engine.get(), anotherEngine.get(),
00336                 etaKernel.get(), phiKernel.get(),
00337                 etaDependentScaleFactors, convolverMinBin,
00338                 convolverMaxBin, fixEfficiency));
00339     }
00340 }
00341 
00342 
00343 fftjet::PeakFinder FFTJetPatRecoProducer::buildPeakFinder(const edm::ParameterSet& ps)
00344 {
00345     const double peakFinderMaxEta = ps.getParameter<double>("peakFinderMaxEta");
00346     if (peakFinderMaxEta <= 0.0)
00347         throw cms::Exception("FFTJetBadConfig")
00348             << "invalid peak finder eta cut" << std::endl;
00349     const double maxMagnitude = ps.getParameter<double>("peakFinderMaxMagnitude");
00350     int minBin = energyFlow->getEtaBin(-peakFinderMaxEta);
00351     if (minBin < 0)
00352         minBin = 0;
00353     int maxBin = energyFlow->getEtaBin(peakFinderMaxEta) + 1;
00354     if (maxBin < 0)
00355         maxBin = 0;
00356     return fftjet::PeakFinder(maxMagnitude, true, minBin, maxBin);
00357 }
00358 
00359 
00360 FFTJetPatRecoProducer::~FFTJetPatRecoProducer()
00361 {
00362     // do anything here that needs to be done at desctruction time
00363     // (e.g. close files, deallocate resources etc.)
00364     delete clusteringTree;
00365     delete extGrid;
00366 }
00367 
00368 
00369 //
00370 // member functions
00371 //
00372 template<class Real>
00373 void FFTJetPatRecoProducer::buildSparseProduct(edm::Event& ev) const
00374 {
00375     typedef reco::PattRecoTree<Real,reco::PattRecoPeak<Real> > StoredTree;
00376 
00377     std::auto_ptr<StoredTree> tree(new StoredTree());
00378 
00379     sparsePeakTreeToStorable(sparseTree,
00380                              sequencer->maxAdaptiveScales(),
00381                              tree.get());
00382 
00383     // Check that we can restore the tree
00384     if (verifyDataConversion && !storeInSinglePrecision())
00385     {
00386         SparseTree check;
00387         const std::vector<double>& scalesUsed(sequencer->getInitialScales());
00388         sparsePeakTreeFromStorable(*tree, &scalesUsed, getEventScale(), &check);
00389         if (sparseTree != check)
00390             throw cms::Exception("FFTJetInterface")
00391                 << "Data conversion failed for sparse clustering tree"
00392                 << std::endl;
00393     }
00394 
00395     ev.put(tree, outputLabel);
00396 }
00397 
00398 
00399 template<class Real>
00400 void FFTJetPatRecoProducer::buildDenseProduct(edm::Event& ev) const
00401 {
00402     typedef reco::PattRecoTree<Real,reco::PattRecoPeak<Real> > StoredTree;
00403 
00404     std::auto_ptr<StoredTree> tree(new StoredTree());
00405 
00406     densePeakTreeToStorable(*clusteringTree,
00407                             sequencer->maxAdaptiveScales(),
00408                             tree.get());
00409 
00410     // Check that we can restore the tree
00411     if (verifyDataConversion && !storeInSinglePrecision())
00412     {
00413         ClusteringTree check(distanceCalc.get());
00414         const std::vector<double>& scalesUsed(sequencer->getInitialScales());
00415         densePeakTreeFromStorable(*tree, &scalesUsed, getEventScale(), &check);
00416         if (*clusteringTree != check)
00417             throw cms::Exception("FFTJetInterface")
00418                 << "Data conversion failed for dense clustering tree"
00419                 << std::endl;
00420     }
00421 
00422     ev.put(tree, outputLabel);
00423 }
00424 
00425 
00426 // ------------ method called to produce the data  ------------
00427 void FFTJetPatRecoProducer::produce(
00428     edm::Event& iEvent, const edm::EventSetup& iSetup)
00429 {
00430     loadInputCollection(iEvent);
00431     discretizeEnergyFlow();
00432 
00433     if (makeClusteringTree)
00434     {
00435         sequencer->run(*energyFlow, clusteringTree);
00436         if (getEventScale() > 0.0)
00437             sequencer->insertCompleteEvent(getEventScale(), *energyFlow,
00438                                            clusteringTree, completeEventDataCutoff);
00439 
00440         if (sparsify)
00441         {
00442             sparsifier->sparsify(*clusteringTree, &sparseTree);
00443 
00444             // Do not call the "sortNodes" method of the sparse tree here.
00445             // Currently, the nodes are sorted by daughter number.
00446             // This is the way we want it in storage because the stored
00447             // tree does not include daughter ordering info explicitly.
00448 
00449             if (storeInSinglePrecision())
00450                 buildSparseProduct<float>(iEvent);
00451             else
00452                 buildSparseProduct<double>(iEvent);
00453         }
00454         else
00455         {
00456             if (storeInSinglePrecision())
00457                 buildDenseProduct<float>(iEvent);
00458             else
00459                 buildDenseProduct<double>(iEvent);
00460         }
00461     }
00462 
00463     if (storeDiscretizationGrid)
00464     {
00465         const fftjet::Grid2d<Real>& g(*energyFlow);
00466 
00467         std::auto_ptr<reco::DiscretizedEnergyFlow> flow(
00468             new reco::DiscretizedEnergyFlow(
00469                 g.data(), g.title(), g.etaMin(), g.etaMax(),
00470                 g.phiBin0Edge(), g.nEta(), g.nPhi()));
00471 
00472         if (verifyDataConversion)
00473         {
00474             fftjet::Grid2d<Real> check(
00475                 flow->nEtaBins(), flow->etaMin(), flow->etaMax(),
00476                 flow->nPhiBins(), flow->phiBin0Edge(), flow->title());
00477             check.blockSet(flow->data(), flow->nEtaBins(), flow->nPhiBins());
00478             assert(g == check);
00479         }
00480 
00481         iEvent.put(flow, outputLabel);
00482     }
00483 
00484     if (storeGridsExternally)
00485     {
00486         if (extGrid)
00487             copy_Grid2d_data(extGrid, *energyFlow);
00488         else
00489             extGrid = convert_Grid2d_to_float(*energyFlow);
00490         if (!extGrid->write(externalGridStream))
00491         {
00492             throw cms::Exception("FFTJetPatRecoProducer::produce")
00493                 << "Failed to write grid data into an external file"
00494                 << std::endl;
00495         }
00496     }
00497 }
00498 
00499 
00500 // ------------ method called once each job just before starting event loop
00501 void FFTJetPatRecoProducer::beginJob()
00502 {
00503 }
00504 
00505 
00506 // ------------ method called once each job just after ending the event loop
00507 void FFTJetPatRecoProducer::endJob()
00508 {
00509     if (storeGridsExternally)
00510         externalGridStream.close();
00511 }
00512 
00513 
00514 //define this as a plug-in
00515 DEFINE_FWK_MODULE(FFTJetPatRecoProducer);