CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/RecoJets/FFTJetProducers/plugins/FFTJetProducer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    FFTJetProducers
00004 // Class:      FFTJetProducer
00005 // 
00013 //
00014 // Original Author:  Igor Volobouev
00015 //         Created:  Sun Jun 20 14:32:36 CDT 2010
00016 // $Id: FFTJetProducer.cc,v 1.11 2011/10/26 22:57:17 igv Exp $
00017 //
00018 //
00019 
00020 #include <iostream>
00021 #include <fstream>
00022 #include <algorithm>
00023 
00024 // Header for this class
00025 #include "RecoJets/FFTJetProducers/plugins/FFTJetProducer.h"
00026 
00027 // Additional FFTJet headers
00028 #include "fftjet/VectorRecombinationAlgFactory.hh"
00029 #include "fftjet/RecombinationAlgFactory.hh"
00030 
00031 // Framework include files
00032 #include "FWCore/Framework/interface/MakerMacros.h"
00033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00034 
00035 // Data formats
00036 #include "DataFormats/Common/interface/View.h"
00037 #include "DataFormats/Common/interface/Handle.h"
00038 #include "DataFormats/VertexReco/interface/Vertex.h"
00039 #include "DataFormats/JetReco/interface/FFTCaloJetCollection.h"
00040 #include "DataFormats/JetReco/interface/FFTGenJetCollection.h"
00041 #include "DataFormats/JetReco/interface/FFTPFJetCollection.h"
00042 #include "DataFormats/JetReco/interface/FFTJPTJetCollection.h"
00043 #include "DataFormats/JetReco/interface/FFTBasicJetCollection.h"
00044 #include "DataFormats/JetReco/interface/FFTTrackJetCollection.h"
00045 #include "DataFormats/Candidate/interface/Candidate.h"
00046 #include "DataFormats/JetReco/interface/FFTJetProducerSummary.h"
00047 
00048 #include "RecoJets/FFTJetProducers/interface/FFTJetParameterParser.h"
00049 
00050 #include "RecoJets/FFTJetAlgorithms/interface/clusteringTreeConverters.h"
00051 #include "RecoJets/FFTJetAlgorithms/interface/jetConverters.h"
00052 #include "RecoJets/FFTJetAlgorithms/interface/matchOneToOne.h"
00053 #include "RecoJets/FFTJetAlgorithms/interface/JetToPeakDistance.h"
00054 
00055 #include "DataFormats/JetReco/interface/DiscretizedEnergyFlow.h"
00056 
00057 #include "RecoJets/JetProducers/interface/JetSpecific.h"
00058 
00059 
00060 #define make_param(type, varname) const \
00061     type & varname (ps.getParameter< type >( #varname ))
00062 
00063 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
00064 
00065 // A generic switch statement based on jet type.
00066 // Defining it in a single place avoids potential errors
00067 // in case another jet type is introduced in the future.
00068 //
00069 // JPTJet is omitted for now: there is no reco::writeSpecific method
00070 // for it (see header JetSpecific.h in the JetProducers package) 
00071 //
00072 #define jet_type_switch(method, arg1, arg2) do {\
00073     switch (jetType)\
00074     {\
00075     case CALOJET:\
00076         method <reco::CaloJet> ( arg1 , arg2 );\
00077         break;\
00078     case PFJET:\
00079         method <reco::PFJet> ( arg1 , arg2 );\
00080         break;\
00081     case GENJET:\
00082         method <reco::GenJet> ( arg1 , arg2 );\
00083         break;\
00084     case TRACKJET:\
00085         method <reco::TrackJet> ( arg1 , arg2 );\
00086         break;\
00087     case BASICJET:\
00088         method <reco::BasicJet> ( arg1 , arg2 );\
00089         break;\
00090     default:\
00091         assert(!"ERROR in FFTJetProducer : invalid jet type."\
00092                " This is a bug. Please report.");\
00093     }\
00094 } while(0);
00095 
00096 using namespace fftjetcms;
00097 
00098 FFTJetProducer::Resolution FFTJetProducer::parse_resolution(
00099     const std::string& name)
00100 {
00101     if (!name.compare("fixed"))
00102         return FIXED;
00103     else if (!name.compare("maximallyStable"))
00104         return MAXIMALLY_STABLE;
00105     else if (!name.compare("globallyAdaptive"))
00106         return GLOBALLY_ADAPTIVE;
00107     else if (!name.compare("locallyAdaptive"))
00108         return LOCALLY_ADAPTIVE;
00109     else
00110         throw cms::Exception("FFTJetBadConfig")
00111             << "Invalid resolution specification \""
00112             << name << "\"" << std::endl;
00113 }
00114 
00115 
00116 template <typename T>
00117 void FFTJetProducer::makeProduces(
00118     const std::string& alias, const std::string& tag)
00119 {
00120     produces<std::vector<reco::FFTAnyJet<T> > >(tag).setBranchAlias(alias);
00121 }
00122 
00123 //
00124 // constructors and destructor
00125 //
00126 FFTJetProducer::FFTJetProducer(const edm::ParameterSet& ps)
00127     : FFTJetInterface(ps),
00128       myConfiguration(ps),
00129       init_param(edm::InputTag, treeLabel),
00130       init_param(bool, useGriddedAlgorithm),
00131       init_param(bool, reuseExistingGrid),
00132       init_param(unsigned, maxIterations),
00133       init_param(unsigned, nJetsRequiredToConverge),
00134       init_param(double, convergenceDistance),
00135       init_param(bool, assignConstituents),
00136       init_param(bool, resumConstituents),
00137       init_param(bool, calculatePileup),
00138       init_param(bool, subtractPileup),
00139       init_param(bool, subtractPileupAs4Vec),
00140       init_param(edm::InputTag, pileupLabel),
00141       init_param(double, fixedScale),
00142       init_param(double, minStableScale),
00143       init_param(double, maxStableScale),
00144       init_param(double, stabilityAlpha),
00145       init_param(double, noiseLevel),
00146       init_param(unsigned, nClustersRequested),
00147       init_param(double, gridScanMaxEta),
00148       init_param(std::string, recombinationAlgorithm),
00149       init_param(bool, isCrisp),
00150       init_param(double, unlikelyBgWeight),
00151       init_param(double, recombinationDataCutoff),
00152       resolution(parse_resolution(ps.getParameter<std::string>("resolution"))),
00153 
00154       minLevel(0),
00155       maxLevel(0),
00156       usedLevel(0),
00157       unused(0.0),
00158       iterationsPerformed(1U),
00159       constituents(200)
00160 {
00161     // Check that the settings make sense
00162     if (resumConstituents && !assignConstituents)
00163         throw cms::Exception("FFTJetBadConfig") 
00164             << "Can't resum constituents if they are not assigned"
00165             << std::endl;
00166 
00167     produces<reco::FFTJetProducerSummary>(outputLabel);
00168     const std::string alias(ps.getUntrackedParameter<std::string>(
00169                                 "alias", outputLabel));
00170     jet_type_switch(makeProduces, alias, outputLabel);
00171 
00172     // Build the set of pattern recognition scales.
00173     // This is needed in order to read the clustering tree
00174     // from the event record.
00175     iniScales = fftjet_ScaleSet_parser(
00176         ps.getParameter<edm::ParameterSet>("InitialScales"));
00177     checkConfig(iniScales, "invalid set of scales");
00178     std::sort(iniScales->begin(), iniScales->end(), std::greater<double>());
00179 
00180     // Most of the configuration has to be performed inside
00181     // the "beginJob" method. This is because chaining of the
00182     // parsers between this base class and the derived classes
00183     // can not work from the constructor of the base class.
00184 }
00185 
00186 
00187 FFTJetProducer::~FFTJetProducer()
00188 {
00189 }
00190 
00191 //
00192 // member functions
00193 //
00194 template<class Real>
00195 void FFTJetProducer::loadSparseTreeData(const edm::Event& iEvent)
00196 {
00197     typedef reco::PattRecoTree<Real,reco::PattRecoPeak<Real> > StoredTree;
00198 
00199     // Get the input
00200     edm::Handle<StoredTree> input;
00201     iEvent.getByLabel(treeLabel, input);
00202 
00203     if (!input->isSparse())
00204         throw cms::Exception("FFTJetBadConfig") 
00205             << "The stored clustering tree is not sparse" << std::endl;
00206 
00207     sparsePeakTreeFromStorable(*input, iniScales.get(),
00208                                getEventScale(), &sparseTree);
00209     sparseTree.sortNodes();
00210 }
00211 
00212 
00213 void FFTJetProducer::selectPreclusters(
00214     const SparseTree& tree,
00215     const fftjet::Functor1<bool,fftjet::Peak>& peakSelector,
00216     std::vector<fftjet::Peak>* preclusters)
00217 {
00218     nodes.clear();
00219     selectTreeNodes(tree, peakSelector, &nodes);
00220 
00221     // Fill out the vector of preclusters using the tree node ids
00222     const unsigned nNodes = nodes.size();
00223     const SparseTree::NodeId* pnodes = nNodes ? &nodes[0] : 0;
00224     preclusters->reserve(nNodes);
00225     for (unsigned i=0; i<nNodes; ++i)
00226         preclusters->push_back(
00227             sparseTree.uncheckedNode(pnodes[i]).getCluster());
00228 
00229     // Remember the node id in the precluster and set
00230     // the status word to indicate the resolution scheme used
00231     fftjet::Peak* clusters = nNodes ? &(*preclusters)[0] : 0;
00232     for (unsigned i=0; i<nNodes; ++i)
00233     {
00234         clusters[i].setCode(pnodes[i]);
00235         clusters[i].setStatus(resolution);
00236     }
00237 }
00238 
00239 
00240 void FFTJetProducer::selectTreeNodes(
00241     const SparseTree& tree,
00242     const fftjet::Functor1<bool,fftjet::Peak>& peakSelect,
00243     std::vector<SparseTree::NodeId>* mynodes)
00244 {
00245     minLevel = maxLevel = usedLevel = 0;
00246 
00247     // Get the tree nodes which pass the cuts
00248     // (according to the selected resolution strategy)
00249     switch (resolution)
00250     {
00251     case FIXED:
00252     {
00253         usedLevel = tree.getLevel(fixedScale);
00254         tree.getPassingNodes(usedLevel, peakSelect, mynodes);
00255     }
00256     break;
00257 
00258     case MAXIMALLY_STABLE:
00259     {
00260         const unsigned minStartingLevel = maxStableScale > 0.0 ? 
00261             tree.getLevel(maxStableScale) : 0;
00262         const unsigned maxStartingLevel = minStableScale > 0.0 ?
00263             tree.getLevel(minStableScale) : UINT_MAX;
00264 
00265         if (tree.stableClusterCount(
00266                 peakSelect, &minLevel, &maxLevel, stabilityAlpha,
00267                 minStartingLevel, maxStartingLevel))
00268         {
00269             usedLevel = (minLevel + maxLevel)/2;
00270             tree.getPassingNodes(usedLevel, peakSelect, mynodes);
00271         }
00272     }
00273     break;
00274 
00275     case GLOBALLY_ADAPTIVE:
00276     {
00277         const bool stable = tree.clusterCountLevels(
00278             nClustersRequested, peakSelect, &minLevel, &maxLevel);
00279         if (minLevel || maxLevel)
00280         {
00281             usedLevel = (minLevel + maxLevel)/2;
00282             if (!stable)
00283             {
00284                 const int maxlev = tree.maxStoredLevel();
00285                 bool levelFound = false;
00286                 for (int delta=0; delta<=maxlev && !levelFound; ++delta)
00287                     for (int ifac=1; ifac>-2 && !levelFound; ifac-=2)
00288                     {
00289                         const int level = usedLevel + ifac*delta;
00290                         if (level > 0 && level <= maxlev)
00291                             if (occupancy[level] == nClustersRequested)
00292                             {
00293                                 usedLevel = level;
00294                                 levelFound = true;
00295                             }
00296                     }
00297                 assert(levelFound);
00298             }
00299         }
00300         else
00301         {
00302             // Can't find that exact number of preclusters.
00303             // Try to get the next best thing.
00304             usedLevel = 1;
00305             const unsigned occ1 = occupancy[1];
00306             if (nClustersRequested >= occ1)
00307             {
00308                 const unsigned maxlev = tree.maxStoredLevel();
00309                 if (nClustersRequested > occupancy[maxlev])
00310                     usedLevel = maxlev;
00311                 else
00312                 {
00313                     // It would be nice to use "lower_bound" here,
00314                     // but the occupancy is not necessarily monotonous.
00315                     unsigned bestDelta = nClustersRequested > occ1 ?
00316                         nClustersRequested - occ1 : occ1 - nClustersRequested;
00317                     for (unsigned level=2; level<=maxlev; ++level)
00318                     {
00319                         const unsigned n = occupancy[level];
00320                         const unsigned d = nClustersRequested > n ? 
00321                             nClustersRequested - n : n - nClustersRequested;
00322                         if (d < bestDelta)
00323                         {
00324                             bestDelta = d;
00325                             usedLevel = level;
00326                         }
00327                     }
00328                 }
00329             }
00330         }
00331         tree.getPassingNodes(usedLevel, peakSelect, mynodes);
00332     }
00333     break;
00334 
00335     case LOCALLY_ADAPTIVE:
00336     {
00337         usedLevel = tree.getLevel(fixedScale);
00338         tree.getMagS2OptimalNodes(peakSelect, nClustersRequested,
00339                                   usedLevel, mynodes, &thresholds);
00340     }
00341     break;
00342 
00343     default:
00344         assert(!"ERROR in FFTJetProducer::selectTreeNodes : "
00345                "should never get here! This is a bug. Please report.");
00346     }
00347 }
00348 
00349 
00350 void FFTJetProducer::prepareRecombinationScales()
00351 {
00352     const unsigned nClus = preclusters.size();
00353     if (nClus)
00354     {
00355         fftjet::Peak* clus = &preclusters[0];
00356         fftjet::Functor1<double,fftjet::Peak>& 
00357             scaleCalc(*recoScaleCalcPeak);
00358         fftjet::Functor1<double,fftjet::Peak>& 
00359             ratioCalc(*recoScaleRatioCalcPeak);
00360         fftjet::Functor1<double,fftjet::Peak>& 
00361             factorCalc(*memberFactorCalcPeak);
00362 
00363         for (unsigned i=0; i<nClus; ++i)
00364         {
00365             clus[i].setRecoScale(scaleCalc(clus[i]));
00366             clus[i].setRecoScaleRatio(ratioCalc(clus[i]));
00367             clus[i].setMembershipFactor(factorCalc(clus[i]));
00368         }
00369     }
00370 }
00371 
00372 
00373 void FFTJetProducer::buildGridAlg()
00374 {
00375     int minBin = energyFlow->getEtaBin(-gridScanMaxEta);
00376     if (minBin < 0)
00377         minBin = 0;
00378     int maxBin = energyFlow->getEtaBin(gridScanMaxEta) + 1;
00379     if (maxBin < 0)
00380         maxBin = 0;
00381     
00382     fftjet::DefaultRecombinationAlgFactory<
00383         Real,VectorLike,BgData,VBuilder> factory;
00384     if (factory[recombinationAlgorithm] == NULL)
00385         throw cms::Exception("FFTJetBadConfig")
00386             << "Invalid grid recombination algorithm \""
00387             << recombinationAlgorithm << "\"" << std::endl;
00388     gridAlg = std::auto_ptr<GridAlg>(
00389         factory[recombinationAlgorithm]->create(
00390             jetMembershipFunction.get(),
00391             bgMembershipFunction.get(),
00392             unlikelyBgWeight, recombinationDataCutoff,
00393             isCrisp, false, assignConstituents, minBin, maxBin));
00394 }
00395 
00396 
00397 bool FFTJetProducer::loadEnergyFlow(
00398     const edm::Event& iEvent, const edm::InputTag& label,
00399     std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >& flow)
00400 {
00401     edm::Handle<reco::DiscretizedEnergyFlow> input;
00402     iEvent.getByLabel(label, input);
00403 
00404     // Make sure that the grid is compatible with the stored one
00405     bool rebuildGrid = flow.get() == NULL;
00406     if (!rebuildGrid)
00407         rebuildGrid = 
00408             !(flow->nEta() == input->nEtaBins() &&
00409               flow->nPhi() == input->nPhiBins() &&
00410               flow->etaMin() == input->etaMin() &&
00411               flow->etaMax() == input->etaMax() &&
00412               flow->phiBin0Edge() == input->phiBin0Edge());
00413     if (rebuildGrid)
00414     {
00415         // We should not get here very often...
00416         flow = std::auto_ptr<fftjet::Grid2d<Real> >(
00417             new fftjet::Grid2d<Real>(
00418                 input->nEtaBins(), input->etaMin(), input->etaMax(),
00419                 input->nPhiBins(), input->phiBin0Edge(), input->title()));
00420     }
00421     flow->blockSet(input->data(), input->nEtaBins(), input->nPhiBins());
00422     return rebuildGrid;
00423 }
00424 
00425 
00426 bool FFTJetProducer::checkConvergence(const std::vector<RecoFFTJet>& previous,
00427                                       std::vector<RecoFFTJet>& nextSet)
00428 {
00429     fftjet::Functor2<double,RecoFFTJet,RecoFFTJet>&
00430         distanceCalc(*jetDistanceCalc);
00431 
00432     const unsigned nJets = previous.size();
00433     if (nJets != nextSet.size())
00434         return false;
00435 
00436     const RecoFFTJet* prev = &previous[0];
00437     RecoFFTJet* next = &nextSet[0];
00438 
00439     // Calculate convergence distances for all jets
00440     bool converged = true;
00441     for (unsigned i=0; i<nJets; ++i)
00442     {
00443         const double d = distanceCalc(prev[i], next[i]);
00444         next[i].setConvergenceDistance(d);
00445         if (i < nJetsRequiredToConverge && d > convergenceDistance)
00446             converged = false;
00447     }
00448 
00449     return converged;
00450 }
00451 
00452 
00453 unsigned FFTJetProducer::iterateJetReconstruction()
00454 {
00455     fftjet::Functor1<double,RecoFFTJet>& scaleCalc(*recoScaleCalcJet);
00456     fftjet::Functor1<double,RecoFFTJet>& ratioCalc(*recoScaleRatioCalcJet);
00457     fftjet::Functor1<double,RecoFFTJet>& factorCalc(*memberFactorCalcJet);
00458 
00459     unsigned nJets = recoJets.size();
00460     unsigned iterNum = 1U;
00461     bool converged = false;
00462     for (; iterNum<maxIterations && !converged; ++iterNum)
00463     {
00464         // Recreate the vector of preclusters using the jets
00465         const RecoFFTJet* jets = &recoJets[0];
00466         iterPreclusters.clear();
00467         iterPreclusters.reserve(nJets);
00468         for (unsigned i=0; i<nJets; ++i)
00469         {
00470             const RecoFFTJet& jet(jets[i]);
00471             fftjet::Peak p(jet.precluster());
00472             p.setEtaPhi(jet.vec().Eta(), jet.vec().Phi());
00473             p.setRecoScale(scaleCalc(jet));
00474             p.setRecoScaleRatio(ratioCalc(jet));
00475             p.setMembershipFactor(factorCalc(jet));
00476             iterPreclusters.push_back(p);
00477         }
00478 
00479         // Run the algorithm
00480         int status = 0;
00481         if (useGriddedAlgorithm)
00482             status = gridAlg->run(iterPreclusters, *energyFlow,
00483                                   &noiseLevel, 1U, 1U,
00484                                   &iterJets, &unclustered, &unused);
00485         else
00486             status = recoAlg->run(iterPreclusters, eventData, &noiseLevel, 1U,
00487                                   &iterJets, &unclustered, &unused);
00488         if (status)
00489             throw cms::Exception("FFTJetInterface")
00490                 << "FFTJet algorithm failed" << std::endl;
00491 
00492         // As it turns out, it is possible, in very rare cases,
00493         // to have iterJets.size() != nJets at this point
00494 
00495         // Figure out if the iterations have converged
00496         converged = checkConvergence(recoJets, iterJets);
00497 
00498         // Prepare for the next cycle
00499         iterJets.swap(recoJets);
00500         nJets = recoJets.size();
00501     }
00502 
00503     // Check that we have the correct number of preclusters
00504     if (preclusters.size() != nJets)
00505     {
00506         assert(nJets < preclusters.size());
00507         removeFakePreclusters();
00508         assert(preclusters.size() == nJets);
00509     }
00510 
00511     // Plug in the original precluster coordinates into the result
00512     RecoFFTJet* jets = &recoJets[0];
00513     for (unsigned i=0; i<nJets; ++i)
00514     {
00515         const fftjet::Peak& oldp(preclusters[i]);
00516         jets[i].setPeakEtaPhi(oldp.eta(), oldp.phi());
00517     }
00518 
00519     // If we have converged on the last cycle, the result
00520     // would be indistinguishable from no convergence.
00521     // Because of this, raise the counter by one to indicate
00522     // the case when the convergence is not achieved.
00523     if (!converged)
00524         ++iterNum;
00525 
00526     return iterNum;
00527 }
00528 
00529 
00530 void FFTJetProducer::determineGriddedConstituents()
00531 {
00532     const unsigned nJets = recoJets.size();
00533     const unsigned* clusterMask = gridAlg->getClusterMask();
00534     const int nEta = gridAlg->getLastNEta();
00535     const int nPhi = gridAlg->getLastNPhi();
00536     const fftjet::Grid2d<Real>& g(*energyFlow);
00537 
00538     const unsigned nInputs = eventData.size();
00539     const VectorLike* inp = nInputs ? &eventData[0] : 0;
00540     const unsigned* candIdx = nInputs ? &candidateIndex[0] : 0;
00541     for (unsigned i=0; i<nInputs; ++i)
00542     {
00543         const VectorLike& item(inp[i]);
00544         const int iPhi = g.getPhiBin(item.Phi());
00545         const int iEta = g.getEtaBin(item.Eta());
00546         const unsigned mask = iEta >= 0 && iEta < nEta ?
00547             clusterMask[iEta*nPhi + iPhi] : 0;
00548         assert(mask <= nJets);
00549         constituents[mask].push_back(inputCollection->ptrAt(candIdx[i]));
00550     }
00551 }
00552 
00553 
00554 void FFTJetProducer::determineVectorConstituents()
00555 {
00556     const unsigned nJets = recoJets.size();
00557     const unsigned* clusterMask = recoAlg->getClusterMask();
00558     const unsigned maskLength = recoAlg->getLastNData();
00559     assert(maskLength == eventData.size());
00560 
00561     const unsigned* candIdx = maskLength ? &candidateIndex[0] : 0;
00562     for (unsigned i=0; i<maskLength; ++i)
00563     {
00564         // In FFTJet, the mask value of 0 corresponds to unclustered
00565         // energy. We will do the same here. Jet numbers are therefore
00566         // shifted by 1 wrt constituents vector, and constituents[1]
00567         // corresponds to jet number 0.
00568         const unsigned mask = clusterMask[i];
00569         assert(mask <= nJets);
00570         constituents[mask].push_back(inputCollection->ptrAt(candIdx[i]));
00571     }
00572 }
00573 
00574 
00575 // The following code more or less coincides with the similar method
00576 // implemented in VirtualJetProducer
00577 template <typename T>
00578 void FFTJetProducer::writeJets(edm::Event& iEvent,
00579                                const edm::EventSetup& iSetup)
00580 {
00581     using namespace reco;
00582 
00583     typedef FFTAnyJet<T> OutputJet;
00584     typedef std::vector<OutputJet> OutputCollection;
00585 
00586     // Area of a single eta-phi cell for jet area calculations.
00587     // Set it to 0 in case the module configuration does not allow
00588     // us to calculate jet areas reliably.
00589     double cellArea = useGriddedAlgorithm && 
00590                       recombinationDataCutoff < 0.0 ?
00591         energyFlow->etaBinWidth() * energyFlow->phiBinWidth() : 0.0;
00592 
00593     if (calculatePileup)
00594         cellArea = pileupEnergyFlow->etaBinWidth() *
00595                    pileupEnergyFlow->phiBinWidth();
00596 
00597     // allocate output jet collection
00598     std::auto_ptr<OutputCollection> jets(new OutputCollection());
00599     const unsigned nJets = recoJets.size();
00600     jets->reserve(nJets);
00601 
00602     for (unsigned ijet=0; ijet<nJets; ++ijet)
00603     {
00604         RecoFFTJet& myjet(recoJets[ijet]);
00605 
00606         // Check if we should resum jet constituents
00607         VectorLike jet4vec(myjet.vec());
00608         if (resumConstituents)
00609         {
00610             VectorLike sum(0.0, 0.0, 0.0, 0.0);
00611             const unsigned nCon = constituents[ijet+1].size();
00612             const reco::CandidatePtr* cn = nCon ? &constituents[ijet+1][0] : 0;
00613             for (unsigned i=0; i<nCon; ++i)
00614                 sum += cn[i]->p4();
00615             jet4vec = sum;
00616             setJetStatusBit(&myjet, CONSTITUENTS_RESUMMED, true);
00617         }
00618 
00619         // Subtract the pile-up
00620         if (calculatePileup && subtractPileup)
00621         {
00622             if (subtractPileupAs4Vec)
00623             {
00624                 jet4vec -= pileup[ijet];
00625                 setJetStatusBit(&myjet, PILEUP_SUBTRACTED_4VEC, true);
00626             }
00627             else
00628             {
00629                 const double pt = jet4vec.Pt();
00630                 if (pt > 0.0)
00631                 {
00632                     const double pileupPt = pileup[ijet].Pt();
00633                     jet4vec *= ((pt - pileupPt)/pt);
00634                 }
00635                 setJetStatusBit(&myjet, PILEUP_SUBTRACTED_PT, true);
00636             }
00637         }
00638 
00639         // Write the specifics to the jet (simultaneously sets 4-vector,
00640         // vertex, constituents). These are overridden functions that will
00641         // call the appropriate specific code.
00642         T jet;
00643         writeSpecific(jet, jet4vec, vertexUsed(),
00644                       constituents[ijet+1], iSetup);
00645 
00646         // calcuate the jet area
00647         double ncells = myjet.ncells();
00648         if (calculatePileup)
00649         {
00650             ncells = cellCountsVec[ijet];
00651             setJetStatusBit(&myjet, PILEUP_CALCULATED, true);
00652         }
00653         jet.setJetArea(cellArea*ncells);
00654 
00655         // add jet to the list
00656         FFTJet<float> fj(jetToStorable<float>(myjet));
00657         fj.setFourVec(jet4vec);
00658         if (calculatePileup)
00659         {
00660             fj.setPileup(pileup[ijet]);
00661             fj.setNCells(ncells);
00662         }
00663         jets->push_back(OutputJet(jet, fj));
00664     }
00665 
00666     // put the collection into the event
00667     iEvent.put(jets, outputLabel);
00668 }
00669 
00670 
00671 void FFTJetProducer::saveResults(edm::Event& ev, const edm::EventSetup& iSetup,
00672                                  const unsigned nPreclustersFound)
00673 {
00674     // Write recombined jets
00675     jet_type_switch(writeJets, ev, iSetup);
00676 
00677     // Check if we should resum unclustered energy constituents
00678     VectorLike unclusE(unclustered);
00679     if (resumConstituents)
00680     {
00681         VectorLike sum(0.0, 0.0, 0.0, 0.0);
00682         const unsigned nCon = constituents[0].size();
00683         const reco::CandidatePtr* cn = nCon ? &constituents[0][0] : 0;
00684         for (unsigned i=0; i<nCon; ++i)
00685             sum += cn[i]->p4();
00686         unclusE = sum;
00687     }
00688 
00689     // Write the jet reconstruction summary
00690     const double minScale = minLevel ? sparseTree.getScale(minLevel) : 0.0;
00691     const double maxScale = maxLevel ? sparseTree.getScale(maxLevel) : 0.0;
00692     const double scaleUsed = usedLevel ? sparseTree.getScale(usedLevel) : 0.0;
00693 
00694     std::auto_ptr<reco::FFTJetProducerSummary> summary(
00695         new reco::FFTJetProducerSummary(
00696             thresholds, occupancy, unclusE,
00697             constituents[0], unused,
00698             minScale, maxScale, scaleUsed,
00699             nPreclustersFound, iterationsPerformed,
00700             iterationsPerformed == 1U ||
00701             iterationsPerformed <= maxIterations));
00702     ev.put(summary, outputLabel);
00703 }
00704 
00705 
00706 // ------------ method called to for each event  ------------
00707 void FFTJetProducer::produce(edm::Event& iEvent,
00708                              const edm::EventSetup& iSetup)
00709 {
00710     // Load the clustering tree made by FFTJetPatRecoProducer
00711     if (storeInSinglePrecision())
00712         loadSparseTreeData<float>(iEvent);
00713     else
00714         loadSparseTreeData<double>(iEvent);
00715 
00716     // Do we need to load the candidate collection?
00717     if (assignConstituents || !(useGriddedAlgorithm && reuseExistingGrid))
00718         loadInputCollection(iEvent);
00719 
00720     // Do we need to have discretized energy flow?
00721     if (useGriddedAlgorithm)
00722     {
00723         if (reuseExistingGrid)
00724         {
00725             if (loadEnergyFlow(iEvent, treeLabel, energyFlow))
00726                 buildGridAlg();
00727         }
00728         else
00729             discretizeEnergyFlow();
00730     }
00731 
00732     // Calculate cluster occupancy as a function of level number
00733     sparseTree.occupancyInScaleSpace(*peakSelector, &occupancy);
00734 
00735     // Select the preclusters using the requested resolution scheme
00736     preclusters.clear();
00737     selectPreclusters(sparseTree, *peakSelector, &preclusters);
00738 
00739     // Prepare to run the jet recombination procedure
00740     prepareRecombinationScales();
00741 
00742     // Assign membership functions to preclusters. If this function
00743     // is not overriden in a derived class, default algorithm membership
00744     // function will be used for every cluster.
00745     assignMembershipFunctions(&preclusters);
00746 
00747     // Count the preclusters going in
00748     unsigned nPreclustersFound = 0U;
00749     const unsigned npre = preclusters.size();
00750     for (unsigned i=0; i<npre; ++i)
00751         if (preclusters[i].membershipFactor() > 0.0)
00752             ++nPreclustersFound;
00753 
00754     // Run the recombination algorithm once
00755     int status = 0;
00756     if (useGriddedAlgorithm)
00757         status = gridAlg->run(preclusters, *energyFlow,
00758                               &noiseLevel, 1U, 1U,
00759                               &recoJets, &unclustered, &unused);
00760     else
00761         status = recoAlg->run(preclusters, eventData, &noiseLevel, 1U,
00762                               &recoJets, &unclustered, &unused);
00763     if (status)
00764         throw cms::Exception("FFTJetInterface")
00765             << "FFTJet algorithm failed (first iteration)" << std::endl;
00766 
00767     // If requested, iterate the jet recombination procedure
00768     if (maxIterations > 1U && !recoJets.empty())
00769     {
00770         // It is possible to have a smaller number of jets than we had
00771         // preclusters. Fake preclusters are possible, but for a good
00772         // choice of pattern recognition kernel their presence should
00773         // be infrequent. However, any fake preclusters will throw the
00774         // iterative reconstruction off balance. Deal with the problem now.
00775         const unsigned nJets = recoJets.size();
00776         if (preclusters.size() != nJets)
00777         {
00778             assert(nJets < preclusters.size());
00779             removeFakePreclusters();
00780         }
00781         iterationsPerformed = iterateJetReconstruction();
00782     }
00783     else
00784         iterationsPerformed = 1U;
00785 
00786     // Determine jet constituents. FFTJet returns a map
00787     // of constituents which is inverse to what we need here.
00788     const unsigned nJets = recoJets.size();
00789     if (constituents.size() <= nJets)
00790         constituents.resize(nJets + 1U);
00791     if (assignConstituents)
00792     {
00793         for (unsigned i=0; i<=nJets; ++i)
00794             constituents[i].clear();
00795         if (useGriddedAlgorithm)
00796             determineGriddedConstituents();
00797         else
00798             determineVectorConstituents();
00799     }
00800 
00801     // Figure out the pile-up
00802     if (calculatePileup)
00803     {
00804         determinePileupDensity(iEvent, pileupLabel, pileupEnergyFlow);
00805         determinePileup();
00806         assert(pileup.size() == recoJets.size());
00807     }
00808 
00809     // Write out the results
00810     saveResults(iEvent, iSetup, nPreclustersFound);
00811 }
00812 
00813 
00814 std::auto_ptr<fftjet::Functor1<bool,fftjet::Peak> >
00815 FFTJetProducer::parse_peakSelector(const edm::ParameterSet& ps)
00816 {
00817     return fftjet_PeakSelector_parser(
00818         ps.getParameter<edm::ParameterSet>("PeakSelectorConfiguration"));
00819 }
00820 
00821 
00822 // Parser for the jet membership function
00823 std::auto_ptr<fftjet::ScaleSpaceKernel>
00824 FFTJetProducer::parse_jetMembershipFunction(const edm::ParameterSet& ps)
00825 {
00826     return fftjet_MembershipFunction_parser(
00827         ps.getParameter<edm::ParameterSet>("jetMembershipFunction"));
00828 }
00829 
00830 
00831 // Parser for the background membership function
00832 std::auto_ptr<AbsBgFunctor>
00833 FFTJetProducer::parse_bgMembershipFunction(const edm::ParameterSet& ps)
00834 {
00835     return fftjet_BgFunctor_parser(
00836         ps.getParameter<edm::ParameterSet>("bgMembershipFunction"));
00837 }
00838 
00839 
00840 // Calculator for the recombination scale
00841 std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
00842 FFTJetProducer::parse_recoScaleCalcPeak(const edm::ParameterSet& ps)
00843 {
00844     return fftjet_PeakFunctor_parser(
00845         ps.getParameter<edm::ParameterSet>("recoScaleCalcPeak"));
00846 }
00847 
00848 
00849 // Pile-up density calculator
00850 std::auto_ptr<fftjetcms::AbsPileupCalculator>
00851 FFTJetProducer::parse_pileupDensityCalc(const edm::ParameterSet& ps)
00852 {
00853     return fftjet_PileupCalculator_parser(
00854         ps.getParameter<edm::ParameterSet>("pileupDensityCalc"));
00855 }
00856 
00857 
00858 // Calculator for the recombination scale ratio
00859 std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
00860 FFTJetProducer::parse_recoScaleRatioCalcPeak(const edm::ParameterSet& ps)
00861 {
00862     return fftjet_PeakFunctor_parser(
00863         ps.getParameter<edm::ParameterSet>("recoScaleRatioCalcPeak"));
00864 }
00865 
00866 
00867 // Calculator for the membership function factor
00868 std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
00869 FFTJetProducer::parse_memberFactorCalcPeak(const edm::ParameterSet& ps)
00870 {
00871     return fftjet_PeakFunctor_parser(
00872         ps.getParameter<edm::ParameterSet>("memberFactorCalcPeak"));
00873 }
00874 
00875 
00876 std::auto_ptr<fftjet::Functor1<double,FFTJetProducer::RecoFFTJet> >
00877 FFTJetProducer::parse_recoScaleCalcJet(const edm::ParameterSet& ps)
00878 {
00879     return fftjet_JetFunctor_parser(
00880         ps.getParameter<edm::ParameterSet>("recoScaleCalcJet"));
00881 }
00882 
00883 
00884 std::auto_ptr<fftjet::Functor1<double,FFTJetProducer::RecoFFTJet> >
00885 FFTJetProducer::parse_recoScaleRatioCalcJet(const edm::ParameterSet& ps)
00886 {
00887     return fftjet_JetFunctor_parser(
00888         ps.getParameter<edm::ParameterSet>("recoScaleRatioCalcJet"));
00889 }
00890 
00891 
00892 std::auto_ptr<fftjet::Functor1<double,FFTJetProducer::RecoFFTJet> >
00893 FFTJetProducer::parse_memberFactorCalcJet(const edm::ParameterSet& ps)
00894 {
00895     return fftjet_JetFunctor_parser(
00896         ps.getParameter<edm::ParameterSet>("memberFactorCalcJet"));
00897 }
00898 
00899 
00900 std::auto_ptr<fftjet::Functor2<
00901     double,FFTJetProducer::RecoFFTJet,FFTJetProducer::RecoFFTJet> >
00902 FFTJetProducer::parse_jetDistanceCalc(const edm::ParameterSet& ps)
00903 {
00904     return fftjet_JetDistance_parser(
00905         ps.getParameter<edm::ParameterSet>("jetDistanceCalc"));
00906 }
00907 
00908 
00909 void FFTJetProducer::assignMembershipFunctions(std::vector<fftjet::Peak>*)
00910 {
00911 }
00912 
00913 
00914 // ------------ method called once each job just before starting event loop
00915 void FFTJetProducer::beginJob()
00916 {
00917     const edm::ParameterSet& ps(myConfiguration);
00918 
00919     // Parse the peak selector definition
00920     peakSelector = parse_peakSelector(ps);
00921     checkConfig(peakSelector, "invalid peak selector");
00922 
00923     jetMembershipFunction = parse_jetMembershipFunction(ps);
00924     checkConfig(jetMembershipFunction, "invalid jet membership function");
00925 
00926     bgMembershipFunction = parse_bgMembershipFunction(ps);
00927     checkConfig(bgMembershipFunction, "invalid noise membership function");
00928 
00929     // Build the energy recombination algorithm
00930     if (!useGriddedAlgorithm)
00931     {
00932         fftjet::DefaultVectorRecombinationAlgFactory<
00933             VectorLike,BgData,VBuilder> factory;
00934         if (factory[recombinationAlgorithm] == NULL)
00935             throw cms::Exception("FFTJetBadConfig")
00936                 << "Invalid vector recombination algorithm \""
00937                 << recombinationAlgorithm << "\"" << std::endl;
00938         recoAlg = std::auto_ptr<RecoAlg>(
00939             factory[recombinationAlgorithm]->create(
00940                 jetMembershipFunction.get(),
00941                 &VectorLike::Et, &VectorLike::Eta, &VectorLike::Phi,
00942                 bgMembershipFunction.get(),
00943                 unlikelyBgWeight, isCrisp, false, assignConstituents));
00944     }
00945     else if (!reuseExistingGrid)
00946     {
00947         energyFlow = fftjet_Grid2d_parser(
00948             ps.getParameter<edm::ParameterSet>("GridConfiguration"));
00949         checkConfig(energyFlow, "invalid discretization grid");
00950         buildGridAlg();
00951     }
00952 
00953     // Create the grid subsequently used for pile-up subtraction
00954     if (calculatePileup)
00955     {
00956         pileupEnergyFlow = fftjet_Grid2d_parser(
00957             ps.getParameter<edm::ParameterSet>("PileupGridConfiguration"));
00958         checkConfig(pileupEnergyFlow, "invalid pileup density grid");
00959 
00960         pileupDensityCalc = parse_pileupDensityCalc(ps);
00961         checkConfig(pileupDensityCalc, "invalid pile-up density calculator");
00962     }
00963 
00964     // Parse the calculator of the recombination scale
00965     recoScaleCalcPeak = parse_recoScaleCalcPeak(ps);
00966     checkConfig(recoScaleCalcPeak, "invalid spec for the "
00967                 "reconstruction scale calculator from peaks");
00968 
00969     // Parse the calculator of the recombination scale ratio
00970     recoScaleRatioCalcPeak = parse_recoScaleRatioCalcPeak(ps);
00971     checkConfig(recoScaleRatioCalcPeak, "invalid spec for the "
00972                 "reconstruction scale ratio calculator from peaks");
00973 
00974     // Calculator for the membership function factor
00975     memberFactorCalcPeak = parse_memberFactorCalcPeak(ps);
00976     checkConfig(memberFactorCalcPeak, "invalid spec for the "
00977                 "membership function factor calculator from peaks");
00978 
00979     if (maxIterations > 1)
00980     {
00981         // We are going to run iteratively. Make required objects.
00982         recoScaleCalcJet = parse_recoScaleCalcJet(ps);
00983         checkConfig(recoScaleCalcJet, "invalid spec for the "
00984                     "reconstruction scale calculator from jets");
00985 
00986         recoScaleRatioCalcJet = parse_recoScaleRatioCalcJet(ps);
00987         checkConfig(recoScaleRatioCalcJet, "invalid spec for the "
00988                     "reconstruction scale ratio calculator from jets");
00989 
00990         memberFactorCalcJet = parse_memberFactorCalcJet(ps);
00991         checkConfig(memberFactorCalcJet, "invalid spec for the "
00992                     "membership function factor calculator from jets");
00993 
00994         jetDistanceCalc = parse_jetDistanceCalc(ps);
00995         checkConfig(memberFactorCalcJet, "invalid spec for the "
00996                     "jet distance calculator");
00997     }
00998 }
00999 
01000 
01001 void FFTJetProducer::removeFakePreclusters()
01002 {
01003     // There are two possible reasons for fake preclusters:
01004     // 1. Membership factor was set to 0
01005     // 2. Genuine problem with pattern recognition
01006     //
01007     // Anyway, we need to match jets to preclusters and keep
01008     // only those preclusters that have been matched
01009     //
01010     std::vector<int> matchTable;
01011     const unsigned nmatched = matchOneToOne(
01012         recoJets, preclusters, JetToPeakDistance(), &matchTable);
01013 
01014     // Ensure that all jets have been matched.
01015     // If not, we must have a bug somewhere.
01016     assert(nmatched == recoJets.size());
01017 
01018     // Collect all matched preclusters
01019     iterPreclusters.clear();
01020     iterPreclusters.reserve(nmatched);
01021     for (unsigned i=0; i<nmatched; ++i)
01022         iterPreclusters.push_back(preclusters[matchTable[i]]);
01023     iterPreclusters.swap(preclusters);
01024 }
01025 
01026 
01027 void FFTJetProducer::setJetStatusBit(RecoFFTJet* jet,
01028                                      const int mask, const bool value)
01029 {
01030     int status = jet->status();
01031     if (value)
01032         status |= mask;
01033     else
01034         status &= ~mask;
01035     jet->setStatus(status);
01036 }
01037 
01038 
01039 void FFTJetProducer::determinePileupDensity(
01040     const edm::Event& iEvent, const edm::InputTag& label,
01041     std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >& density)
01042 {
01043     edm::Handle<reco::FFTJetPileupSummary> summary;
01044     iEvent.getByLabel(label, summary);
01045 
01046     const reco::FFTJetPileupSummary& s(*summary);
01047     const AbsPileupCalculator& calc(*pileupDensityCalc);
01048     const bool phiDependent = calc.isPhiDependent();
01049 
01050     fftjet::Grid2d<Real>& g(*density);
01051     const unsigned nEta = g.nEta();
01052     const unsigned nPhi = g.nPhi();
01053 
01054     for (unsigned ieta=0; ieta<nEta; ++ieta)
01055     {
01056         const double eta(g.etaBinCenter(ieta));
01057 
01058         if (phiDependent)
01059         {
01060             for (unsigned iphi=0; iphi<nPhi; ++iphi)
01061             {
01062                 const double phi(g.phiBinCenter(iphi));
01063                 g.uncheckedSetBin(ieta, iphi, calc(eta, phi, s));
01064             }
01065         }
01066         else
01067         {
01068             const double pil = calc(eta, 0.0, s);
01069             for (unsigned iphi=0; iphi<nPhi; ++iphi)
01070                 g.uncheckedSetBin(ieta, iphi, pil);
01071         }
01072     }
01073 }
01074 
01075 
01076 void FFTJetProducer::determinePileup()
01077 {
01078     // This function works with crisp clustering only
01079     if (!isCrisp)
01080         assert(!"Pile-up subtraction for fuzzy clustering "
01081                "is not implemented yet");
01082 
01083     // Clear the pileup vector
01084     const unsigned nJets = recoJets.size();
01085     pileup.resize(nJets);
01086     if (nJets == 0)
01087         return;
01088     const VectorLike zero;
01089     for (unsigned i=0; i<nJets; ++i)
01090         pileup[i] = zero;
01091 
01092     // Pileup energy flow grid
01093     const fftjet::Grid2d<Real>& g(*pileupEnergyFlow);
01094     const unsigned nEta = g.nEta();
01095     const unsigned nPhi = g.nPhi();
01096     const double cellArea = g.etaBinWidth() * g.phiBinWidth();
01097 
01098     // Various calculators
01099     fftjet::Functor1<double,RecoFFTJet>& scaleCalc(*recoScaleCalcJet);
01100     fftjet::Functor1<double,RecoFFTJet>& ratioCalc(*recoScaleRatioCalcJet);
01101     fftjet::Functor1<double,RecoFFTJet>& factorCalc(*memberFactorCalcJet);
01102 
01103     // Make sure we have enough memory
01104     memFcns2dVec.resize(nJets);
01105     fftjet::AbsKernel2d** memFcns2d = &memFcns2dVec[0];
01106 
01107     doubleBuf.resize(nJets*4U + nJets*nPhi);
01108     double* recoScales = &doubleBuf[0];
01109     double* recoScaleRatios = recoScales + nJets;
01110     double* memFactors = recoScaleRatios + nJets;
01111     double* dEta = memFactors + nJets;
01112     double* dPhiBuf = dEta + nJets;
01113 
01114     cellCountsVec.resize(nJets);
01115     unsigned* cellCounts = &cellCountsVec[0];
01116 
01117     // Go over jets and collect the necessary info
01118     for (unsigned ijet=0; ijet<nJets; ++ijet)
01119     {
01120         const RecoFFTJet& jet(recoJets[ijet]);
01121         const fftjet::Peak& peak(jet.precluster());
01122 
01123         // Make sure we are using 2-d membership functions.
01124         // Pile-up subtraction scheme for 3-d functions should be different.
01125         fftjet::AbsMembershipFunction* m3d = 
01126             dynamic_cast<fftjet::AbsMembershipFunction*>(
01127                 peak.membershipFunction());
01128         if (m3d == 0)
01129             m3d = dynamic_cast<fftjet::AbsMembershipFunction*>(
01130                 jetMembershipFunction.get());
01131         if (m3d)
01132         {
01133             assert(!"Pile-up subtraction for 3-d membership functions "
01134                    "is not implemented yet");
01135         }
01136         else
01137         {
01138             fftjet::AbsKernel2d* m2d =
01139                 dynamic_cast<fftjet::AbsKernel2d*>(
01140                     peak.membershipFunction());
01141             if (m2d == 0)
01142                 m2d = dynamic_cast<fftjet::AbsKernel2d*>(
01143                     jetMembershipFunction.get());
01144             assert(m2d);
01145             memFcns2d[ijet] = m2d;
01146         }
01147         recoScales[ijet] = scaleCalc(jet);
01148         recoScaleRatios[ijet] = ratioCalc(jet);
01149         memFactors[ijet] = factorCalc(jet);
01150         cellCounts[ijet] = 0U;
01151 
01152         const double jetPhi = jet.vec().Phi();
01153         for (unsigned iphi=0; iphi<nPhi; ++iphi)
01154         {
01155             double dphi = g.phiBinCenter(iphi) - jetPhi;
01156             while (dphi > M_PI)
01157                 dphi -= (2.0*M_PI);
01158             while (dphi < -M_PI)
01159                 dphi += (2.0*M_PI);
01160             dPhiBuf[iphi*nJets+ijet] = dphi;
01161         }
01162     }
01163 
01164     // Go over all grid points and integrate
01165     // the pile-up energy density
01166     VBuilder vMaker;
01167     for (unsigned ieta=0; ieta<nEta; ++ieta)
01168     {
01169         const double eta(g.etaBinCenter(ieta));
01170         const Real* databuf = g.data() + ieta*nPhi;
01171 
01172         // Figure out dEta for each jet
01173         for (unsigned ijet=0; ijet<nJets; ++ijet)
01174             dEta[ijet] = eta - recoJets[ijet].vec().Eta();
01175 
01176         for (unsigned iphi=0; iphi<nPhi; ++iphi)
01177         {
01178             double maxW(0.0);
01179             unsigned maxWJet(nJets);
01180             const double* dPhi = dPhiBuf + iphi*nJets;
01181 
01182             for (unsigned ijet=0; ijet<nJets; ++ijet)
01183             {
01184                 if (recoScaleRatios[ijet] > 0.0)
01185                     memFcns2d[ijet]->setScaleRatio(recoScaleRatios[ijet]);
01186                 const double f = memFactors[ijet]*
01187                     (*memFcns2d[ijet])(dEta[ijet], dPhi[ijet],
01188                                        recoScales[ijet]);
01189                 if (f > maxW)
01190                 {
01191                     maxW = f;
01192                     maxWJet = ijet;
01193                 }
01194             }
01195 
01196             if (maxWJet < nJets)
01197             {
01198                 pileup[maxWJet] += vMaker(cellArea*databuf[iphi],
01199                                           eta, g.phiBinCenter(iphi));
01200                 cellCounts[maxWJet]++;
01201             }
01202         }
01203     }
01204 }
01205 
01206 
01207 // ------------ method called once each job just after ending the event loop
01208 void FFTJetProducer::endJob()
01209 {
01210 }
01211 
01212 
01213 //define this as a plug-in
01214 DEFINE_FWK_MODULE(FFTJetProducer);