00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020 #include <iostream>
00021 #include <fstream>
00022 #include <functional>
00023 #include <algorithm>
00024
00025
00026 #include "RecoJets/FFTJetProducers/plugins/FFTJetProducer.h"
00027
00028
00029 #include "fftjet/VectorRecombinationAlgFactory.hh"
00030 #include "fftjet/RecombinationAlgFactory.hh"
00031
00032
00033 #include "FWCore/Framework/interface/MakerMacros.h"
00034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00035
00036
00037 #include "DataFormats/Common/interface/View.h"
00038 #include "DataFormats/Common/interface/Handle.h"
00039 #include "DataFormats/VertexReco/interface/Vertex.h"
00040 #include "DataFormats/JetReco/interface/FFTCaloJetCollection.h"
00041 #include "DataFormats/JetReco/interface/FFTGenJetCollection.h"
00042 #include "DataFormats/JetReco/interface/FFTPFJetCollection.h"
00043 #include "DataFormats/JetReco/interface/FFTJPTJetCollection.h"
00044 #include "DataFormats/JetReco/interface/FFTBasicJetCollection.h"
00045 #include "DataFormats/JetReco/interface/FFTTrackJetCollection.h"
00046 #include "DataFormats/Candidate/interface/Candidate.h"
00047 #include "DataFormats/JetReco/interface/FFTJetProducerSummary.h"
00048
00049 #include "RecoJets/FFTJetProducers/interface/FFTJetParameterParser.h"
00050
00051 #include "RecoJets/FFTJetAlgorithms/interface/clusteringTreeConverters.h"
00052 #include "RecoJets/FFTJetAlgorithms/interface/jetConverters.h"
00053 #include "RecoJets/FFTJetAlgorithms/interface/matchOneToOne.h"
00054 #include "RecoJets/FFTJetAlgorithms/interface/JetToPeakDistance.h"
00055 #include "RecoJets/FFTJetAlgorithms/interface/adjustForPileup.h"
00056
00057 #include "DataFormats/JetReco/interface/DiscretizedEnergyFlow.h"
00058
00059 #include "RecoJets/JetProducers/interface/JetSpecific.h"
00060
00061
00062 #include "JetMETCorrections/FFTJetObjects/interface/FFTJetLookupTableSequenceLoader.h"
00063
00064
00065 #define make_param(type, varname) const \
00066 type & varname (ps.getParameter< type >( #varname ))
00067
00068 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
00069
00070
00071
00072
00073
00074
00075
00076
00077 #define jet_type_switch(method, arg1, arg2) do {\
00078 switch (jetType)\
00079 {\
00080 case CALOJET:\
00081 method <reco::CaloJet> ( arg1 , arg2 );\
00082 break;\
00083 case PFJET:\
00084 method <reco::PFJet> ( arg1 , arg2 );\
00085 break;\
00086 case GENJET:\
00087 method <reco::GenJet> ( arg1 , arg2 );\
00088 break;\
00089 case TRACKJET:\
00090 method <reco::TrackJet> ( arg1 , arg2 );\
00091 break;\
00092 case BASICJET:\
00093 method <reco::BasicJet> ( arg1 , arg2 );\
00094 break;\
00095 default:\
00096 assert(!"ERROR in FFTJetProducer : invalid jet type."\
00097 " This is a bug. Please report.");\
00098 }\
00099 } while(0);
00100
00101 namespace {
00102 struct LocalSortByPt
00103 {
00104 template<class Jet>
00105 inline bool operator()(const Jet& l, const Jet& r) const
00106 {return l.pt() > r.pt();}
00107 };
00108 }
00109
00110 using namespace fftjetcms;
00111
00112 FFTJetProducer::Resolution FFTJetProducer::parse_resolution(
00113 const std::string& name)
00114 {
00115 if (!name.compare("fixed"))
00116 return FIXED;
00117 else if (!name.compare("maximallyStable"))
00118 return MAXIMALLY_STABLE;
00119 else if (!name.compare("globallyAdaptive"))
00120 return GLOBALLY_ADAPTIVE;
00121 else if (!name.compare("locallyAdaptive"))
00122 return LOCALLY_ADAPTIVE;
00123 else if (!name.compare("fromGenJets"))
00124 return FROM_GENJETS;
00125 else
00126 throw cms::Exception("FFTJetBadConfig")
00127 << "Invalid resolution specification \""
00128 << name << "\"" << std::endl;
00129 }
00130
00131
00132 template <typename T>
00133 void FFTJetProducer::makeProduces(
00134 const std::string& alias, const std::string& tag)
00135 {
00136 produces<std::vector<reco::FFTAnyJet<T> > >(tag).setBranchAlias(alias);
00137 }
00138
00139
00140
00141
00142 FFTJetProducer::FFTJetProducer(const edm::ParameterSet& ps)
00143 : FFTJetInterface(ps),
00144 myConfiguration(ps),
00145 init_param(edm::InputTag, treeLabel),
00146 init_param(bool, useGriddedAlgorithm),
00147 init_param(bool, reuseExistingGrid),
00148 init_param(unsigned, maxIterations),
00149 init_param(unsigned, nJetsRequiredToConverge),
00150 init_param(double, convergenceDistance),
00151 init_param(bool, assignConstituents),
00152 init_param(bool, resumConstituents),
00153 init_param(bool, calculatePileup),
00154 init_param(bool, subtractPileup),
00155 init_param(bool, subtractPileupAs4Vec),
00156 init_param(edm::InputTag, pileupLabel),
00157 init_param(double, fixedScale),
00158 init_param(double, minStableScale),
00159 init_param(double, maxStableScale),
00160 init_param(double, stabilityAlpha),
00161 init_param(double, noiseLevel),
00162 init_param(unsigned, nClustersRequested),
00163 init_param(double, gridScanMaxEta),
00164 init_param(std::string, recombinationAlgorithm),
00165 init_param(bool, isCrisp),
00166 init_param(double, unlikelyBgWeight),
00167 init_param(double, recombinationDataCutoff),
00168 init_param(edm::InputTag, genJetsLabel),
00169 init_param(unsigned, maxInitialPreclusters),
00170 resolution(parse_resolution(ps.getParameter<std::string>("resolution"))),
00171 init_param(std::string, pileupTableRecord),
00172 init_param(std::string, pileupTableName),
00173 init_param(std::string, pileupTableCategory),
00174 init_param(bool, loadPileupFromDB),
00175
00176 minLevel(0),
00177 maxLevel(0),
00178 usedLevel(0),
00179 unused(0.0),
00180 iterationsPerformed(1U),
00181 constituents(200)
00182 {
00183
00184 if (resumConstituents && !assignConstituents)
00185 throw cms::Exception("FFTJetBadConfig")
00186 << "Can't resum constituents if they are not assigned"
00187 << std::endl;
00188
00189 produces<reco::FFTJetProducerSummary>(outputLabel);
00190 const std::string alias(ps.getUntrackedParameter<std::string>(
00191 "alias", outputLabel));
00192 jet_type_switch(makeProduces, alias, outputLabel);
00193
00194
00195
00196
00197 iniScales = fftjet_ScaleSet_parser(
00198 ps.getParameter<edm::ParameterSet>("InitialScales"));
00199 checkConfig(iniScales, "invalid set of scales");
00200 std::sort(iniScales->begin(), iniScales->end(), std::greater<double>());
00201
00202
00203
00204
00205
00206 }
00207
00208
00209 FFTJetProducer::~FFTJetProducer()
00210 {
00211 }
00212
00213
00214
00215
00216 template<class Real>
00217 void FFTJetProducer::loadSparseTreeData(const edm::Event& iEvent)
00218 {
00219 typedef reco::PattRecoTree<Real,reco::PattRecoPeak<Real> > StoredTree;
00220
00221
00222 edm::Handle<StoredTree> input;
00223 iEvent.getByLabel(treeLabel, input);
00224
00225 if (!input->isSparse())
00226 throw cms::Exception("FFTJetBadConfig")
00227 << "The stored clustering tree is not sparse" << std::endl;
00228
00229 sparsePeakTreeFromStorable(*input, iniScales.get(),
00230 getEventScale(), &sparseTree);
00231 sparseTree.sortNodes();
00232 }
00233
00234
00235 void FFTJetProducer::genJetPreclusters(
00236 const SparseTree& ,
00237 edm::Event& iEvent, const edm::EventSetup& ,
00238 const fftjet::Functor1<bool,fftjet::Peak>& peakSelect,
00239 std::vector<fftjet::Peak>* preclusters)
00240 {
00241 typedef reco::FFTAnyJet<reco::GenJet> InputJet;
00242 typedef std::vector<InputJet> InputCollection;
00243
00244 edm::Handle<InputCollection> input;
00245 iEvent.getByLabel(genJetsLabel, input);
00246
00247 const unsigned sz = input->size();
00248 preclusters->reserve(sz);
00249 for (unsigned i=0; i<sz; ++i)
00250 {
00251 const RecoFFTJet& jet(jetFromStorable((*input)[i].getFFTSpecific()));
00252 fftjet::Peak p(jet.precluster());
00253 const double scale(p.scale());
00254 p.setEtaPhi(jet.vec().Eta(), jet.vec().Phi());
00255 p.setMagnitude(jet.vec().Pt()/scale/scale);
00256 p.setStatus(resolution);
00257 if (peakSelect(p))
00258 preclusters->push_back(p);
00259 }
00260 }
00261
00262
00263 void FFTJetProducer::selectPreclusters(
00264 const SparseTree& tree,
00265 const fftjet::Functor1<bool,fftjet::Peak>& peakSelect,
00266 std::vector<fftjet::Peak>* preclusters)
00267 {
00268 nodes.clear();
00269 selectTreeNodes(tree, peakSelect, &nodes);
00270
00271
00272 const unsigned nNodes = nodes.size();
00273 const SparseTree::NodeId* pnodes = nNodes ? &nodes[0] : 0;
00274 preclusters->reserve(nNodes);
00275 for (unsigned i=0; i<nNodes; ++i)
00276 preclusters->push_back(
00277 sparseTree.uncheckedNode(pnodes[i]).getCluster());
00278
00279
00280
00281 fftjet::Peak* clusters = nNodes ? &(*preclusters)[0] : 0;
00282 for (unsigned i=0; i<nNodes; ++i)
00283 {
00284 clusters[i].setCode(pnodes[i]);
00285 clusters[i].setStatus(resolution);
00286 }
00287 }
00288
00289
00290 void FFTJetProducer::selectTreeNodes(
00291 const SparseTree& tree,
00292 const fftjet::Functor1<bool,fftjet::Peak>& peakSelect,
00293 std::vector<SparseTree::NodeId>* mynodes)
00294 {
00295 minLevel = maxLevel = usedLevel = 0;
00296
00297
00298
00299 switch (resolution)
00300 {
00301 case FIXED:
00302 {
00303 usedLevel = tree.getLevel(fixedScale);
00304 tree.getPassingNodes(usedLevel, peakSelect, mynodes);
00305 }
00306 break;
00307
00308 case MAXIMALLY_STABLE:
00309 {
00310 const unsigned minStartingLevel = maxStableScale > 0.0 ?
00311 tree.getLevel(maxStableScale) : 0;
00312 const unsigned maxStartingLevel = minStableScale > 0.0 ?
00313 tree.getLevel(minStableScale) : UINT_MAX;
00314
00315 if (tree.stableClusterCount(
00316 peakSelect, &minLevel, &maxLevel, stabilityAlpha,
00317 minStartingLevel, maxStartingLevel))
00318 {
00319 usedLevel = (minLevel + maxLevel)/2;
00320 tree.getPassingNodes(usedLevel, peakSelect, mynodes);
00321 }
00322 }
00323 break;
00324
00325 case GLOBALLY_ADAPTIVE:
00326 {
00327 const bool stable = tree.clusterCountLevels(
00328 nClustersRequested, peakSelect, &minLevel, &maxLevel);
00329 if (minLevel || maxLevel)
00330 {
00331 usedLevel = (minLevel + maxLevel)/2;
00332 if (!stable)
00333 {
00334 const int maxlev = tree.maxStoredLevel();
00335 bool levelFound = false;
00336 for (int delta=0; delta<=maxlev && !levelFound; ++delta)
00337 for (int ifac=1; ifac>-2 && !levelFound; ifac-=2)
00338 {
00339 const int level = usedLevel + ifac*delta;
00340 if (level > 0 && level <= maxlev)
00341 if (occupancy[level] == nClustersRequested)
00342 {
00343 usedLevel = level;
00344 levelFound = true;
00345 }
00346 }
00347 assert(levelFound);
00348 }
00349 }
00350 else
00351 {
00352
00353
00354 usedLevel = 1;
00355 const unsigned occ1 = occupancy[1];
00356 if (nClustersRequested >= occ1)
00357 {
00358 const unsigned maxlev = tree.maxStoredLevel();
00359 if (nClustersRequested > occupancy[maxlev])
00360 usedLevel = maxlev;
00361 else
00362 {
00363
00364
00365 unsigned bestDelta = nClustersRequested > occ1 ?
00366 nClustersRequested - occ1 : occ1 - nClustersRequested;
00367 for (unsigned level=2; level<=maxlev; ++level)
00368 {
00369 const unsigned n = occupancy[level];
00370 const unsigned d = nClustersRequested > n ?
00371 nClustersRequested - n : n - nClustersRequested;
00372 if (d < bestDelta)
00373 {
00374 bestDelta = d;
00375 usedLevel = level;
00376 }
00377 }
00378 }
00379 }
00380 }
00381 tree.getPassingNodes(usedLevel, peakSelect, mynodes);
00382 }
00383 break;
00384
00385 case LOCALLY_ADAPTIVE:
00386 {
00387 usedLevel = tree.getLevel(fixedScale);
00388 tree.getMagS2OptimalNodes(peakSelect, nClustersRequested,
00389 usedLevel, mynodes, &thresholds);
00390 }
00391 break;
00392
00393 default:
00394 assert(!"ERROR in FFTJetProducer::selectTreeNodes : "
00395 "should never get here! This is a bug. Please report.");
00396 }
00397 }
00398
00399
00400 void FFTJetProducer::prepareRecombinationScales()
00401 {
00402 const unsigned nClus = preclusters.size();
00403 if (nClus)
00404 {
00405 fftjet::Peak* clus = &preclusters[0];
00406 fftjet::Functor1<double,fftjet::Peak>&
00407 scaleCalc(*recoScaleCalcPeak);
00408 fftjet::Functor1<double,fftjet::Peak>&
00409 ratioCalc(*recoScaleRatioCalcPeak);
00410 fftjet::Functor1<double,fftjet::Peak>&
00411 factorCalc(*memberFactorCalcPeak);
00412
00413 for (unsigned i=0; i<nClus; ++i)
00414 {
00415 clus[i].setRecoScale(scaleCalc(clus[i]));
00416 clus[i].setRecoScaleRatio(ratioCalc(clus[i]));
00417 clus[i].setMembershipFactor(factorCalc(clus[i]));
00418 }
00419 }
00420 }
00421
00422
00423 void FFTJetProducer::buildGridAlg()
00424 {
00425 int minBin = energyFlow->getEtaBin(-gridScanMaxEta);
00426 if (minBin < 0)
00427 minBin = 0;
00428 int maxBin = energyFlow->getEtaBin(gridScanMaxEta) + 1;
00429 if (maxBin < 0)
00430 maxBin = 0;
00431
00432 fftjet::DefaultRecombinationAlgFactory<
00433 Real,VectorLike,BgData,VBuilder> factory;
00434 if (factory[recombinationAlgorithm] == NULL)
00435 throw cms::Exception("FFTJetBadConfig")
00436 << "Invalid grid recombination algorithm \""
00437 << recombinationAlgorithm << "\"" << std::endl;
00438 gridAlg = std::auto_ptr<GridAlg>(
00439 factory[recombinationAlgorithm]->create(
00440 jetMembershipFunction.get(),
00441 bgMembershipFunction.get(),
00442 unlikelyBgWeight, recombinationDataCutoff,
00443 isCrisp, false, assignConstituents, minBin, maxBin));
00444 }
00445
00446
00447 bool FFTJetProducer::loadEnergyFlow(
00448 const edm::Event& iEvent, const edm::InputTag& label,
00449 std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >& flow)
00450 {
00451 edm::Handle<reco::DiscretizedEnergyFlow> input;
00452 iEvent.getByLabel(label, input);
00453
00454
00455 bool rebuildGrid = flow.get() == NULL;
00456 if (!rebuildGrid)
00457 rebuildGrid =
00458 !(flow->nEta() == input->nEtaBins() &&
00459 flow->nPhi() == input->nPhiBins() &&
00460 flow->etaMin() == input->etaMin() &&
00461 flow->etaMax() == input->etaMax() &&
00462 flow->phiBin0Edge() == input->phiBin0Edge());
00463 if (rebuildGrid)
00464 {
00465
00466 flow = std::auto_ptr<fftjet::Grid2d<Real> >(
00467 new fftjet::Grid2d<Real>(
00468 input->nEtaBins(), input->etaMin(), input->etaMax(),
00469 input->nPhiBins(), input->phiBin0Edge(), input->title()));
00470 }
00471 flow->blockSet(input->data(), input->nEtaBins(), input->nPhiBins());
00472 return rebuildGrid;
00473 }
00474
00475
00476 bool FFTJetProducer::checkConvergence(const std::vector<RecoFFTJet>& previous,
00477 std::vector<RecoFFTJet>& nextSet)
00478 {
00479 fftjet::Functor2<double,RecoFFTJet,RecoFFTJet>&
00480 distanceCalc(*jetDistanceCalc);
00481
00482 const unsigned nJets = previous.size();
00483 if (nJets != nextSet.size())
00484 return false;
00485
00486 const RecoFFTJet* prev = &previous[0];
00487 RecoFFTJet* next = &nextSet[0];
00488
00489
00490 bool converged = true;
00491 for (unsigned i=0; i<nJets; ++i)
00492 {
00493 const double d = distanceCalc(prev[i], next[i]);
00494 next[i].setConvergenceDistance(d);
00495 if (i < nJetsRequiredToConverge && d > convergenceDistance)
00496 converged = false;
00497 }
00498
00499 return converged;
00500 }
00501
00502
00503 unsigned FFTJetProducer::iterateJetReconstruction()
00504 {
00505 fftjet::Functor1<double,RecoFFTJet>& scaleCalc(*recoScaleCalcJet);
00506 fftjet::Functor1<double,RecoFFTJet>& ratioCalc(*recoScaleRatioCalcJet);
00507 fftjet::Functor1<double,RecoFFTJet>& factorCalc(*memberFactorCalcJet);
00508
00509 unsigned nJets = recoJets.size();
00510 unsigned iterNum = 1U;
00511 bool converged = false;
00512 for (; iterNum<maxIterations && !converged; ++iterNum)
00513 {
00514
00515 const RecoFFTJet* jets = &recoJets[0];
00516 iterPreclusters.clear();
00517 iterPreclusters.reserve(nJets);
00518 for (unsigned i=0; i<nJets; ++i)
00519 {
00520 const RecoFFTJet& jet(jets[i]);
00521 fftjet::Peak p(jet.precluster());
00522 p.setEtaPhi(jet.vec().Eta(), jet.vec().Phi());
00523 p.setRecoScale(scaleCalc(jet));
00524 p.setRecoScaleRatio(ratioCalc(jet));
00525 p.setMembershipFactor(factorCalc(jet));
00526 iterPreclusters.push_back(p);
00527 }
00528
00529
00530 int status = 0;
00531 if (useGriddedAlgorithm)
00532 status = gridAlg->run(iterPreclusters, *energyFlow,
00533 &noiseLevel, 1U, 1U,
00534 &iterJets, &unclustered, &unused);
00535 else
00536 status = recoAlg->run(iterPreclusters, eventData, &noiseLevel, 1U,
00537 &iterJets, &unclustered, &unused);
00538 if (status)
00539 throw cms::Exception("FFTJetInterface")
00540 << "FFTJet algorithm failed" << std::endl;
00541
00542
00543
00544
00545
00546 converged = checkConvergence(recoJets, iterJets);
00547
00548
00549 iterJets.swap(recoJets);
00550 nJets = recoJets.size();
00551 }
00552
00553
00554 if (preclusters.size() != nJets)
00555 {
00556 assert(nJets < preclusters.size());
00557 removeFakePreclusters();
00558 assert(preclusters.size() == nJets);
00559 }
00560
00561
00562 RecoFFTJet* jets = &recoJets[0];
00563 for (unsigned i=0; i<nJets; ++i)
00564 {
00565 const fftjet::Peak& oldp(preclusters[i]);
00566 jets[i].setPeakEtaPhi(oldp.eta(), oldp.phi());
00567 }
00568
00569
00570
00571
00572
00573 if (!converged)
00574 ++iterNum;
00575
00576 return iterNum;
00577 }
00578
00579
00580 void FFTJetProducer::determineGriddedConstituents()
00581 {
00582 const unsigned nJets = recoJets.size();
00583 const unsigned* clusterMask = gridAlg->getClusterMask();
00584 const int nEta = gridAlg->getLastNEta();
00585 const int nPhi = gridAlg->getLastNPhi();
00586 const fftjet::Grid2d<Real>& g(*energyFlow);
00587
00588 const unsigned nInputs = eventData.size();
00589 const VectorLike* inp = nInputs ? &eventData[0] : 0;
00590 const unsigned* candIdx = nInputs ? &candidateIndex[0] : 0;
00591 for (unsigned i=0; i<nInputs; ++i)
00592 {
00593 const VectorLike& item(inp[i]);
00594 const int iPhi = g.getPhiBin(item.Phi());
00595 const int iEta = g.getEtaBin(item.Eta());
00596 const unsigned mask = iEta >= 0 && iEta < nEta ?
00597 clusterMask[iEta*nPhi + iPhi] : 0;
00598 assert(mask <= nJets);
00599 constituents[mask].push_back(inputCollection->ptrAt(candIdx[i]));
00600 }
00601 }
00602
00603
00604 void FFTJetProducer::determineVectorConstituents()
00605 {
00606 const unsigned nJets = recoJets.size();
00607 const unsigned* clusterMask = recoAlg->getClusterMask();
00608 const unsigned maskLength = recoAlg->getLastNData();
00609 assert(maskLength == eventData.size());
00610
00611 const unsigned* candIdx = maskLength ? &candidateIndex[0] : 0;
00612 for (unsigned i=0; i<maskLength; ++i)
00613 {
00614
00615
00616
00617
00618 const unsigned mask = clusterMask[i];
00619 assert(mask <= nJets);
00620 constituents[mask].push_back(inputCollection->ptrAt(candIdx[i]));
00621 }
00622 }
00623
00624
00625
00626
00627 template <typename T>
00628 void FFTJetProducer::writeJets(edm::Event& iEvent,
00629 const edm::EventSetup& iSetup)
00630 {
00631 using namespace reco;
00632
00633 typedef FFTAnyJet<T> OutputJet;
00634 typedef std::vector<OutputJet> OutputCollection;
00635
00636
00637
00638
00639 double cellArea = useGriddedAlgorithm &&
00640 recombinationDataCutoff < 0.0 ?
00641 energyFlow->etaBinWidth() * energyFlow->phiBinWidth() : 0.0;
00642
00643 if (calculatePileup)
00644 cellArea = pileupEnergyFlow->etaBinWidth() *
00645 pileupEnergyFlow->phiBinWidth();
00646
00647
00648 std::auto_ptr<OutputCollection> jets(new OutputCollection());
00649 const unsigned nJets = recoJets.size();
00650 jets->reserve(nJets);
00651
00652 bool sorted = true;
00653 double previousPt = DBL_MAX;
00654 for (unsigned ijet=0; ijet<nJets; ++ijet)
00655 {
00656 RecoFFTJet& myjet(recoJets[ijet]);
00657
00658
00659 VectorLike jet4vec(myjet.vec());
00660 if (resumConstituents)
00661 {
00662 VectorLike sum(0.0, 0.0, 0.0, 0.0);
00663 const unsigned nCon = constituents[ijet+1].size();
00664 const reco::CandidatePtr* cn = nCon ? &constituents[ijet+1][0] : 0;
00665 for (unsigned i=0; i<nCon; ++i)
00666 sum += cn[i]->p4();
00667 jet4vec = sum;
00668 setJetStatusBit(&myjet, CONSTITUENTS_RESUMMED, true);
00669 }
00670
00671
00672 if (calculatePileup && subtractPileup)
00673 {
00674 jet4vec = adjustForPileup(jet4vec, pileup[ijet],
00675 subtractPileupAs4Vec);
00676 if (subtractPileupAs4Vec)
00677 setJetStatusBit(&myjet, PILEUP_SUBTRACTED_4VEC, true);
00678 else
00679 setJetStatusBit(&myjet, PILEUP_SUBTRACTED_PT, true);
00680 }
00681
00682
00683
00684
00685 T jet;
00686 writeSpecific(jet, jet4vec, vertexUsed(),
00687 constituents[ijet+1], iSetup);
00688
00689
00690 double ncells = myjet.ncells();
00691 if (calculatePileup)
00692 {
00693 ncells = cellCountsVec[ijet];
00694 setJetStatusBit(&myjet, PILEUP_CALCULATED, true);
00695 }
00696 jet.setJetArea(cellArea*ncells);
00697
00698
00699 FFTJet<float> fj(jetToStorable<float>(myjet));
00700 fj.setFourVec(jet4vec);
00701 if (calculatePileup)
00702 {
00703 fj.setPileup(pileup[ijet]);
00704 fj.setNCells(ncells);
00705 }
00706 jets->push_back(OutputJet(jet, fj));
00707
00708
00709 const double pt = jet.pt();
00710 if (pt > previousPt)
00711 sorted = false;
00712 previousPt = pt;
00713 }
00714
00715
00716 if (!sorted)
00717 std::sort(jets->begin(), jets->end(), LocalSortByPt());
00718
00719
00720 iEvent.put(jets, outputLabel);
00721 }
00722
00723
00724 void FFTJetProducer::saveResults(edm::Event& ev, const edm::EventSetup& iSetup,
00725 const unsigned nPreclustersFound)
00726 {
00727
00728 jet_type_switch(writeJets, ev, iSetup);
00729
00730
00731 VectorLike unclusE(unclustered);
00732 if (resumConstituents)
00733 {
00734 VectorLike sum(0.0, 0.0, 0.0, 0.0);
00735 const unsigned nCon = constituents[0].size();
00736 const reco::CandidatePtr* cn = nCon ? &constituents[0][0] : 0;
00737 for (unsigned i=0; i<nCon; ++i)
00738 sum += cn[i]->p4();
00739 unclusE = sum;
00740 }
00741
00742
00743 const double minScale = minLevel ? sparseTree.getScale(minLevel) : 0.0;
00744 const double maxScale = maxLevel ? sparseTree.getScale(maxLevel) : 0.0;
00745 const double scaleUsed = usedLevel ? sparseTree.getScale(usedLevel) : 0.0;
00746
00747 std::auto_ptr<reco::FFTJetProducerSummary> summary(
00748 new reco::FFTJetProducerSummary(
00749 thresholds, occupancy, unclusE,
00750 constituents[0], unused,
00751 minScale, maxScale, scaleUsed,
00752 nPreclustersFound, iterationsPerformed,
00753 iterationsPerformed == 1U ||
00754 iterationsPerformed <= maxIterations));
00755 ev.put(summary, outputLabel);
00756 }
00757
00758
00759
00760 void FFTJetProducer::produce(edm::Event& iEvent,
00761 const edm::EventSetup& iSetup)
00762 {
00763
00764 if (storeInSinglePrecision())
00765 loadSparseTreeData<float>(iEvent);
00766 else
00767 loadSparseTreeData<double>(iEvent);
00768
00769
00770 if (assignConstituents || !(useGriddedAlgorithm && reuseExistingGrid))
00771 loadInputCollection(iEvent);
00772
00773
00774 if (useGriddedAlgorithm)
00775 {
00776 if (reuseExistingGrid)
00777 {
00778 if (loadEnergyFlow(iEvent, treeLabel, energyFlow))
00779 buildGridAlg();
00780 }
00781 else
00782 discretizeEnergyFlow();
00783 }
00784
00785
00786 sparseTree.occupancyInScaleSpace(*peakSelector, &occupancy);
00787
00788
00789 preclusters.clear();
00790 if (resolution == FROM_GENJETS)
00791 genJetPreclusters(sparseTree, iEvent, iSetup,
00792 *peakSelector, &preclusters);
00793 else
00794 selectPreclusters(sparseTree, *peakSelector, &preclusters);
00795 if (preclusters.size() > maxInitialPreclusters)
00796 {
00797 std::sort(preclusters.begin(), preclusters.end(), std::greater<fftjet::Peak>());
00798 preclusters.erase(preclusters.begin()+maxInitialPreclusters, preclusters.end());
00799 }
00800
00801
00802 prepareRecombinationScales();
00803
00804
00805
00806
00807 assignMembershipFunctions(&preclusters);
00808
00809
00810 unsigned nPreclustersFound = 0U;
00811 const unsigned npre = preclusters.size();
00812 for (unsigned i=0; i<npre; ++i)
00813 if (preclusters[i].membershipFactor() > 0.0)
00814 ++nPreclustersFound;
00815
00816
00817 int status = 0;
00818 if (useGriddedAlgorithm)
00819 status = gridAlg->run(preclusters, *energyFlow,
00820 &noiseLevel, 1U, 1U,
00821 &recoJets, &unclustered, &unused);
00822 else
00823 status = recoAlg->run(preclusters, eventData, &noiseLevel, 1U,
00824 &recoJets, &unclustered, &unused);
00825 if (status)
00826 throw cms::Exception("FFTJetInterface")
00827 << "FFTJet algorithm failed (first iteration)" << std::endl;
00828
00829
00830 if (maxIterations > 1U && !recoJets.empty())
00831 {
00832
00833
00834
00835
00836
00837 const unsigned nJets = recoJets.size();
00838 if (preclusters.size() != nJets)
00839 {
00840 assert(nJets < preclusters.size());
00841 removeFakePreclusters();
00842 }
00843 iterationsPerformed = iterateJetReconstruction();
00844 }
00845 else
00846 iterationsPerformed = 1U;
00847
00848
00849
00850 const unsigned nJets = recoJets.size();
00851 if (constituents.size() <= nJets)
00852 constituents.resize(nJets + 1U);
00853 if (assignConstituents)
00854 {
00855 for (unsigned i=0; i<=nJets; ++i)
00856 constituents[i].clear();
00857 if (useGriddedAlgorithm)
00858 determineGriddedConstituents();
00859 else
00860 determineVectorConstituents();
00861 }
00862
00863
00864 if (calculatePileup)
00865 {
00866 if (loadPileupFromDB)
00867 determinePileupDensityFromDB(iEvent, iSetup,
00868 pileupLabel, pileupEnergyFlow);
00869 else
00870 determinePileupDensityFromConfig(iEvent, pileupLabel,
00871 pileupEnergyFlow);
00872 determinePileup();
00873 assert(pileup.size() == recoJets.size());
00874 }
00875
00876
00877 saveResults(iEvent, iSetup, nPreclustersFound);
00878 }
00879
00880
00881 std::auto_ptr<fftjet::Functor1<bool,fftjet::Peak> >
00882 FFTJetProducer::parse_peakSelector(const edm::ParameterSet& ps)
00883 {
00884 return fftjet_PeakSelector_parser(
00885 ps.getParameter<edm::ParameterSet>("PeakSelectorConfiguration"));
00886 }
00887
00888
00889
00890 std::auto_ptr<fftjet::ScaleSpaceKernel>
00891 FFTJetProducer::parse_jetMembershipFunction(const edm::ParameterSet& ps)
00892 {
00893 return fftjet_MembershipFunction_parser(
00894 ps.getParameter<edm::ParameterSet>("jetMembershipFunction"));
00895 }
00896
00897
00898
00899 std::auto_ptr<AbsBgFunctor>
00900 FFTJetProducer::parse_bgMembershipFunction(const edm::ParameterSet& ps)
00901 {
00902 return fftjet_BgFunctor_parser(
00903 ps.getParameter<edm::ParameterSet>("bgMembershipFunction"));
00904 }
00905
00906
00907
00908 std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
00909 FFTJetProducer::parse_recoScaleCalcPeak(const edm::ParameterSet& ps)
00910 {
00911 return fftjet_PeakFunctor_parser(
00912 ps.getParameter<edm::ParameterSet>("recoScaleCalcPeak"));
00913 }
00914
00915
00916
00917 std::auto_ptr<fftjetcms::AbsPileupCalculator>
00918 FFTJetProducer::parse_pileupDensityCalc(const edm::ParameterSet& ps)
00919 {
00920 return fftjet_PileupCalculator_parser(
00921 ps.getParameter<edm::ParameterSet>("pileupDensityCalc"));
00922 }
00923
00924
00925
00926 std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
00927 FFTJetProducer::parse_recoScaleRatioCalcPeak(const edm::ParameterSet& ps)
00928 {
00929 return fftjet_PeakFunctor_parser(
00930 ps.getParameter<edm::ParameterSet>("recoScaleRatioCalcPeak"));
00931 }
00932
00933
00934
00935 std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
00936 FFTJetProducer::parse_memberFactorCalcPeak(const edm::ParameterSet& ps)
00937 {
00938 return fftjet_PeakFunctor_parser(
00939 ps.getParameter<edm::ParameterSet>("memberFactorCalcPeak"));
00940 }
00941
00942
00943 std::auto_ptr<fftjet::Functor1<double,FFTJetProducer::RecoFFTJet> >
00944 FFTJetProducer::parse_recoScaleCalcJet(const edm::ParameterSet& ps)
00945 {
00946 return fftjet_JetFunctor_parser(
00947 ps.getParameter<edm::ParameterSet>("recoScaleCalcJet"));
00948 }
00949
00950
00951 std::auto_ptr<fftjet::Functor1<double,FFTJetProducer::RecoFFTJet> >
00952 FFTJetProducer::parse_recoScaleRatioCalcJet(const edm::ParameterSet& ps)
00953 {
00954 return fftjet_JetFunctor_parser(
00955 ps.getParameter<edm::ParameterSet>("recoScaleRatioCalcJet"));
00956 }
00957
00958
00959 std::auto_ptr<fftjet::Functor1<double,FFTJetProducer::RecoFFTJet> >
00960 FFTJetProducer::parse_memberFactorCalcJet(const edm::ParameterSet& ps)
00961 {
00962 return fftjet_JetFunctor_parser(
00963 ps.getParameter<edm::ParameterSet>("memberFactorCalcJet"));
00964 }
00965
00966
00967 std::auto_ptr<fftjet::Functor2<
00968 double,FFTJetProducer::RecoFFTJet,FFTJetProducer::RecoFFTJet> >
00969 FFTJetProducer::parse_jetDistanceCalc(const edm::ParameterSet& ps)
00970 {
00971 return fftjet_JetDistance_parser(
00972 ps.getParameter<edm::ParameterSet>("jetDistanceCalc"));
00973 }
00974
00975
00976 void FFTJetProducer::assignMembershipFunctions(std::vector<fftjet::Peak>*)
00977 {
00978 }
00979
00980
00981
00982 void FFTJetProducer::beginJob()
00983 {
00984 const edm::ParameterSet& ps(myConfiguration);
00985
00986
00987 peakSelector = parse_peakSelector(ps);
00988 checkConfig(peakSelector, "invalid peak selector");
00989
00990 jetMembershipFunction = parse_jetMembershipFunction(ps);
00991 checkConfig(jetMembershipFunction, "invalid jet membership function");
00992
00993 bgMembershipFunction = parse_bgMembershipFunction(ps);
00994 checkConfig(bgMembershipFunction, "invalid noise membership function");
00995
00996
00997 if (!useGriddedAlgorithm)
00998 {
00999 fftjet::DefaultVectorRecombinationAlgFactory<
01000 VectorLike,BgData,VBuilder> factory;
01001 if (factory[recombinationAlgorithm] == NULL)
01002 throw cms::Exception("FFTJetBadConfig")
01003 << "Invalid vector recombination algorithm \""
01004 << recombinationAlgorithm << "\"" << std::endl;
01005 recoAlg = std::auto_ptr<RecoAlg>(
01006 factory[recombinationAlgorithm]->create(
01007 jetMembershipFunction.get(),
01008 &VectorLike::Et, &VectorLike::Eta, &VectorLike::Phi,
01009 bgMembershipFunction.get(),
01010 unlikelyBgWeight, isCrisp, false, assignConstituents));
01011 }
01012 else if (!reuseExistingGrid)
01013 {
01014 energyFlow = fftjet_Grid2d_parser(
01015 ps.getParameter<edm::ParameterSet>("GridConfiguration"));
01016 checkConfig(energyFlow, "invalid discretization grid");
01017 buildGridAlg();
01018 }
01019
01020
01021 if (calculatePileup)
01022 {
01023 pileupEnergyFlow = fftjet_Grid2d_parser(
01024 ps.getParameter<edm::ParameterSet>("PileupGridConfiguration"));
01025 checkConfig(pileupEnergyFlow, "invalid pileup density grid");
01026
01027 if (!loadPileupFromDB)
01028 {
01029 pileupDensityCalc = parse_pileupDensityCalc(ps);
01030 checkConfig(pileupDensityCalc, "invalid pile-up density calculator");
01031 }
01032 }
01033
01034
01035 recoScaleCalcPeak = parse_recoScaleCalcPeak(ps);
01036 checkConfig(recoScaleCalcPeak, "invalid spec for the "
01037 "reconstruction scale calculator from peaks");
01038
01039
01040 recoScaleRatioCalcPeak = parse_recoScaleRatioCalcPeak(ps);
01041 checkConfig(recoScaleRatioCalcPeak, "invalid spec for the "
01042 "reconstruction scale ratio calculator from peaks");
01043
01044
01045 memberFactorCalcPeak = parse_memberFactorCalcPeak(ps);
01046 checkConfig(memberFactorCalcPeak, "invalid spec for the "
01047 "membership function factor calculator from peaks");
01048
01049 if (maxIterations > 1)
01050 {
01051
01052 recoScaleCalcJet = parse_recoScaleCalcJet(ps);
01053 checkConfig(recoScaleCalcJet, "invalid spec for the "
01054 "reconstruction scale calculator from jets");
01055
01056 recoScaleRatioCalcJet = parse_recoScaleRatioCalcJet(ps);
01057 checkConfig(recoScaleRatioCalcJet, "invalid spec for the "
01058 "reconstruction scale ratio calculator from jets");
01059
01060 memberFactorCalcJet = parse_memberFactorCalcJet(ps);
01061 checkConfig(memberFactorCalcJet, "invalid spec for the "
01062 "membership function factor calculator from jets");
01063
01064 jetDistanceCalc = parse_jetDistanceCalc(ps);
01065 checkConfig(memberFactorCalcJet, "invalid spec for the "
01066 "jet distance calculator");
01067 }
01068 }
01069
01070
01071 void FFTJetProducer::removeFakePreclusters()
01072 {
01073
01074
01075
01076
01077
01078
01079
01080 std::vector<int> matchTable;
01081 const unsigned nmatched = matchOneToOne(
01082 recoJets, preclusters, JetToPeakDistance(), &matchTable);
01083
01084
01085
01086 assert(nmatched == recoJets.size());
01087
01088
01089 iterPreclusters.clear();
01090 iterPreclusters.reserve(nmatched);
01091 for (unsigned i=0; i<nmatched; ++i)
01092 iterPreclusters.push_back(preclusters[matchTable[i]]);
01093 iterPreclusters.swap(preclusters);
01094 }
01095
01096
01097 void FFTJetProducer::setJetStatusBit(RecoFFTJet* jet,
01098 const int mask, const bool value)
01099 {
01100 int status = jet->status();
01101 if (value)
01102 status |= mask;
01103 else
01104 status &= ~mask;
01105 jet->setStatus(status);
01106 }
01107
01108
01109 void FFTJetProducer::determinePileupDensityFromConfig(
01110 const edm::Event& iEvent, const edm::InputTag& label,
01111 std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >& density)
01112 {
01113 edm::Handle<reco::FFTJetPileupSummary> summary;
01114 iEvent.getByLabel(label, summary);
01115
01116 const reco::FFTJetPileupSummary& s(*summary);
01117 const AbsPileupCalculator& calc(*pileupDensityCalc);
01118 const bool phiDependent = calc.isPhiDependent();
01119
01120 fftjet::Grid2d<Real>& g(*density);
01121 const unsigned nEta = g.nEta();
01122 const unsigned nPhi = g.nPhi();
01123
01124 for (unsigned ieta=0; ieta<nEta; ++ieta)
01125 {
01126 const double eta(g.etaBinCenter(ieta));
01127
01128 if (phiDependent)
01129 {
01130 for (unsigned iphi=0; iphi<nPhi; ++iphi)
01131 {
01132 const double phi(g.phiBinCenter(iphi));
01133 g.uncheckedSetBin(ieta, iphi, calc(eta, phi, s));
01134 }
01135 }
01136 else
01137 {
01138 const double pil = calc(eta, 0.0, s);
01139 for (unsigned iphi=0; iphi<nPhi; ++iphi)
01140 g.uncheckedSetBin(ieta, iphi, pil);
01141 }
01142 }
01143 }
01144
01145
01146 void FFTJetProducer::determinePileupDensityFromDB(
01147 const edm::Event& iEvent, const edm::EventSetup& iSetup,
01148 const edm::InputTag& label,
01149 std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >& density)
01150 {
01151 edm::ESHandle<FFTJetLookupTableSequence> h;
01152 StaticFFTJetLookupTableSequenceLoader::instance().load(
01153 iSetup, pileupTableRecord, h);
01154 boost::shared_ptr<npstat::StorableMultivariateFunctor> f =
01155 (*h)[pileupTableCategory][pileupTableName];
01156
01157 edm::Handle<reco::FFTJetPileupSummary> summary;
01158 iEvent.getByLabel(label, summary);
01159
01160 const float rho = summary->pileupRho();
01161 const bool phiDependent = f->minDim() == 3U;
01162
01163 fftjet::Grid2d<Real>& g(*density);
01164 const unsigned nEta = g.nEta();
01165 const unsigned nPhi = g.nPhi();
01166
01167 double functorArg[3] = {0.0, 0.0, 0.0};
01168 if (phiDependent)
01169 functorArg[2] = rho;
01170 else
01171 functorArg[1] = rho;
01172
01173 for (unsigned ieta=0; ieta<nEta; ++ieta)
01174 {
01175 const double eta(g.etaBinCenter(ieta));
01176 functorArg[0] = eta;
01177
01178 if (phiDependent)
01179 {
01180 for (unsigned iphi=0; iphi<nPhi; ++iphi)
01181 {
01182 functorArg[1] = g.phiBinCenter(iphi);
01183 g.uncheckedSetBin(ieta, iphi, (*f)(functorArg, 3U));
01184 }
01185 }
01186 else
01187 {
01188 const double pil = (*f)(functorArg, 2U);
01189 for (unsigned iphi=0; iphi<nPhi; ++iphi)
01190 g.uncheckedSetBin(ieta, iphi, pil);
01191 }
01192 }
01193 }
01194
01195
01196 void FFTJetProducer::determinePileup()
01197 {
01198
01199 if (!isCrisp)
01200 assert(!"Pile-up subtraction for fuzzy clustering "
01201 "is not implemented yet");
01202
01203
01204 const unsigned nJets = recoJets.size();
01205 pileup.resize(nJets);
01206 if (nJets == 0)
01207 return;
01208 const VectorLike zero;
01209 for (unsigned i=0; i<nJets; ++i)
01210 pileup[i] = zero;
01211
01212
01213 const fftjet::Grid2d<Real>& g(*pileupEnergyFlow);
01214 const unsigned nEta = g.nEta();
01215 const unsigned nPhi = g.nPhi();
01216 const double cellArea = g.etaBinWidth() * g.phiBinWidth();
01217
01218
01219 fftjet::Functor1<double,RecoFFTJet>& scaleCalc(*recoScaleCalcJet);
01220 fftjet::Functor1<double,RecoFFTJet>& ratioCalc(*recoScaleRatioCalcJet);
01221 fftjet::Functor1<double,RecoFFTJet>& factorCalc(*memberFactorCalcJet);
01222
01223
01224 memFcns2dVec.resize(nJets);
01225 fftjet::AbsKernel2d** memFcns2d = &memFcns2dVec[0];
01226
01227 doubleBuf.resize(nJets*4U + nJets*nPhi);
01228 double* recoScales = &doubleBuf[0];
01229 double* recoScaleRatios = recoScales + nJets;
01230 double* memFactors = recoScaleRatios + nJets;
01231 double* dEta = memFactors + nJets;
01232 double* dPhiBuf = dEta + nJets;
01233
01234 cellCountsVec.resize(nJets);
01235 unsigned* cellCounts = &cellCountsVec[0];
01236
01237
01238 for (unsigned ijet=0; ijet<nJets; ++ijet)
01239 {
01240 const RecoFFTJet& jet(recoJets[ijet]);
01241 const fftjet::Peak& peak(jet.precluster());
01242
01243
01244
01245 fftjet::AbsMembershipFunction* m3d =
01246 dynamic_cast<fftjet::AbsMembershipFunction*>(
01247 peak.membershipFunction());
01248 if (m3d == 0)
01249 m3d = dynamic_cast<fftjet::AbsMembershipFunction*>(
01250 jetMembershipFunction.get());
01251 if (m3d)
01252 {
01253 assert(!"Pile-up subtraction for 3-d membership functions "
01254 "is not implemented yet");
01255 }
01256 else
01257 {
01258 fftjet::AbsKernel2d* m2d =
01259 dynamic_cast<fftjet::AbsKernel2d*>(
01260 peak.membershipFunction());
01261 if (m2d == 0)
01262 m2d = dynamic_cast<fftjet::AbsKernel2d*>(
01263 jetMembershipFunction.get());
01264 assert(m2d);
01265 memFcns2d[ijet] = m2d;
01266 }
01267 recoScales[ijet] = scaleCalc(jet);
01268 recoScaleRatios[ijet] = ratioCalc(jet);
01269 memFactors[ijet] = factorCalc(jet);
01270 cellCounts[ijet] = 0U;
01271
01272 const double jetPhi = jet.vec().Phi();
01273 for (unsigned iphi=0; iphi<nPhi; ++iphi)
01274 {
01275 double dphi = g.phiBinCenter(iphi) - jetPhi;
01276 while (dphi > M_PI)
01277 dphi -= (2.0*M_PI);
01278 while (dphi < -M_PI)
01279 dphi += (2.0*M_PI);
01280 dPhiBuf[iphi*nJets+ijet] = dphi;
01281 }
01282 }
01283
01284
01285
01286 VBuilder vMaker;
01287 for (unsigned ieta=0; ieta<nEta; ++ieta)
01288 {
01289 const double eta(g.etaBinCenter(ieta));
01290 const Real* databuf = g.data() + ieta*nPhi;
01291
01292
01293 for (unsigned ijet=0; ijet<nJets; ++ijet)
01294 dEta[ijet] = eta - recoJets[ijet].vec().Eta();
01295
01296 for (unsigned iphi=0; iphi<nPhi; ++iphi)
01297 {
01298 double maxW(0.0);
01299 unsigned maxWJet(nJets);
01300 const double* dPhi = dPhiBuf + iphi*nJets;
01301
01302 for (unsigned ijet=0; ijet<nJets; ++ijet)
01303 {
01304 if (recoScaleRatios[ijet] > 0.0)
01305 memFcns2d[ijet]->setScaleRatio(recoScaleRatios[ijet]);
01306 const double f = memFactors[ijet]*
01307 (*memFcns2d[ijet])(dEta[ijet], dPhi[ijet],
01308 recoScales[ijet]);
01309 if (f > maxW)
01310 {
01311 maxW = f;
01312 maxWJet = ijet;
01313 }
01314 }
01315
01316 if (maxWJet < nJets)
01317 {
01318 pileup[maxWJet] += vMaker(cellArea*databuf[iphi],
01319 eta, g.phiBinCenter(iphi));
01320 cellCounts[maxWJet]++;
01321 }
01322 }
01323 }
01324 }
01325
01326
01327
01328 void FFTJetProducer::endJob()
01329 {
01330 }
01331
01332
01333
01334 DEFINE_FWK_MODULE(FFTJetProducer);