CMS 3D CMS Logo

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