CMS 3D CMS Logo

FFTJetProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: FFTJetProducers
4 // Class: FFTJetProducer
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Sun Jun 20 14:32:36 CDT 2010
16 //
17 //
18 
19 #include <iostream>
20 #include <fstream>
21 #include <functional>
22 #include <algorithm>
23 
24 #include "fftjet/peakEtLifetime.hh"
25 
26 // Header for this class
28 
29 // Framework include files
32 
33 // Data formats
38 
40 
41 // Loader for the lookup tables
43 
44 
45 #define make_param(type, varname) const \
46  type & varname (ps.getParameter< type >( #varname ))
47 
48 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
49 
50 // A generic switch statement based on jet type.
51 // Defining it in a single place avoids potential errors
52 // in case another jet type is introduced in the future.
53 //
54 // JPTJet is omitted for now: there is no reco::writeSpecific method
55 // for it (see header JetSpecific.h in the JetProducers package)
56 //
57 #define jet_type_switch(method, arg1, arg2) do {\
58  switch (jetType)\
59  {\
60  case CALOJET:\
61  method <reco::CaloJet> ( arg1 , arg2 );\
62  break;\
63  case PFJET:\
64  method <reco::PFJet> ( arg1 , arg2 );\
65  break;\
66  case GENJET:\
67  method <reco::GenJet> ( arg1 , arg2 );\
68  break;\
69  case TRACKJET:\
70  method <reco::TrackJet> ( arg1 , arg2 );\
71  break;\
72  case BASICJET:\
73  method <reco::BasicJet> ( arg1 , arg2 );\
74  break;\
75  default:\
76  assert(!"ERROR in FFTJetProducer : invalid jet type."\
77  " This is a bug. Please report.");\
78  }\
79 } while(0);
80 
81 namespace {
82  struct LocalSortByPt
83  {
84  template<class Jet>
85  inline bool operator()(const Jet& l, const Jet& r) const
86  {return l.pt() > r.pt();}
87  };
88 }
89 
90 using namespace fftjetcms;
91 
93  const std::string& name)
94 {
95  if (!name.compare("fixed"))
96  return FIXED;
97  else if (!name.compare("maximallyStable"))
98  return MAXIMALLY_STABLE;
99  else if (!name.compare("globallyAdaptive"))
100  return GLOBALLY_ADAPTIVE;
101  else if (!name.compare("locallyAdaptive"))
102  return LOCALLY_ADAPTIVE;
103  else if (!name.compare("fromGenJets"))
104  return FROM_GENJETS;
105  else
106  throw cms::Exception("FFTJetBadConfig")
107  << "Invalid resolution specification \""
108  << name << "\"" << std::endl;
109 }
110 
111 
112 template <typename T>
114  const std::string& alias, const std::string& tag)
115 {
116  produces<std::vector<reco::FFTAnyJet<T> > >(tag).setBranchAlias(alias);
117 }
118 
119 //
120 // constructors and destructor
121 //
123  : FFTJetInterface(ps),
124  myConfiguration(ps),
125  init_param(edm::InputTag, treeLabel),
128  init_param(unsigned, maxIterations),
136  init_param(edm::InputTag, pileupLabel),
137  init_param(double, fixedScale),
138  init_param(double, minStableScale),
139  init_param(double, maxStableScale),
140  init_param(double, stabilityAlpha),
141  init_param(double, noiseLevel),
142  init_param(unsigned, nClustersRequested),
143  init_param(double, gridScanMaxEta),
146  init_param(double, unlikelyBgWeight),
148  init_param(edm::InputTag, genJetsLabel),
150  resolution(parse_resolution(ps.getParameter<std::string>("resolution"))),
155 
156  minLevel(0),
157  maxLevel(0),
158  usedLevel(0),
159  unused(0.0),
160  iterationsPerformed(1U),
161  constituents(200)
162 {
163  // Check that the settings make sense
165  throw cms::Exception("FFTJetBadConfig")
166  << "Can't resum constituents if they are not assigned"
167  << std::endl;
168 
169  produces<reco::FFTJetProducerSummary>(outputLabel);
171  "alias", outputLabel));
173 
174  // Build the set of pattern recognition scales.
175  // This is needed in order to read the clustering tree
176  // from the event record.
178  ps.getParameter<edm::ParameterSet>("InitialScales"));
179  checkConfig(iniScales, "invalid set of scales");
180  std::sort(iniScales->begin(), iniScales->end(), std::greater<double>());
181 
182  input_recotree_token_ = consumes<reco::PattRecoTree<fftjetcms::Real,reco::PattRecoPeak<fftjetcms::Real> > >(treeLabel);
183  input_genjet_token_ = consumes<std::vector<reco::FFTAnyJet<reco::GenJet> > >(genJetsLabel);
184  input_energyflow_token_ = consumes<reco::DiscretizedEnergyFlow>(treeLabel);
185  input_pusummary_token_ = consumes<reco::FFTJetPileupSummary>(pileupLabel);
186 
187 
188  // Most of the configuration has to be performed inside
189  // the "beginJob" method. This is because chaining of the
190  // parsers between this base class and the derived classes
191  // can not work from the constructor of the base class.
192 }
193 
194 
196 {
197 }
198 
199 //
200 // member functions
201 //
202 template<class Real>
204 {
206 
207  // Get the input
209  iEvent.getByToken(input_recotree_token_, input);
210 
211  if (!input->isSparse())
212  throw cms::Exception("FFTJetBadConfig")
213  << "The stored clustering tree is not sparse" << std::endl;
214 
217  sparseTree.sortNodes();
218  fftjet::updateSplitMergeTimes(sparseTree, sparseTree.minScale(),
219  sparseTree.maxScale());
220 }
221 
222 
224  const SparseTree& /* tree */,
225  edm::Event& iEvent, const edm::EventSetup& /* iSetup */,
226  const fftjet::Functor1<bool,fftjet::Peak>& peakSelect,
227  std::vector<fftjet::Peak>* preclusters)
228 {
229  typedef reco::FFTAnyJet<reco::GenJet> InputJet;
230  typedef std::vector<InputJet> InputCollection;
231 
233  iEvent.getByToken(input_genjet_token_, input);
234 
235  const unsigned sz = input->size();
236  preclusters->reserve(sz);
237  for (unsigned i=0; i<sz; ++i)
238  {
239  const RecoFFTJet& jet(jetFromStorable((*input)[i].getFFTSpecific()));
240  fftjet::Peak p(jet.precluster());
241  const double scale(p.scale());
242  p.setEtaPhi(jet.vec().Eta(), jet.vec().Phi());
243  p.setMagnitude(jet.vec().Pt()/scale/scale);
244  p.setStatus(resolution);
245  if (peakSelect(p))
246  preclusters->push_back(p);
247  }
248 }
249 
250 
252  const SparseTree& tree,
253  const fftjet::Functor1<bool,fftjet::Peak>& peakSelect,
254  std::vector<fftjet::Peak>* preclusters)
255 {
256  nodes.clear();
257  selectTreeNodes(tree, peakSelect, &nodes);
258 
259  // Fill out the vector of preclusters using the tree node ids
260  const unsigned nNodes = nodes.size();
261  const SparseTree::NodeId* pnodes = nNodes ? &nodes[0] : nullptr;
262  preclusters->reserve(nNodes);
263  for (unsigned i=0; i<nNodes; ++i)
264  preclusters->push_back(
265  sparseTree.uncheckedNode(pnodes[i]).getCluster());
266 
267  // Remember the node id in the precluster and set
268  // the status word to indicate the resolution scheme used
269  fftjet::Peak* clusters = nNodes ? &(*preclusters)[0] : nullptr;
270  for (unsigned i=0; i<nNodes; ++i)
271  {
272  clusters[i].setCode(pnodes[i]);
273  clusters[i].setStatus(resolution);
274  }
275 }
276 
277 
279  const SparseTree& tree,
280  const fftjet::Functor1<bool,fftjet::Peak>& peakSelect,
281  std::vector<SparseTree::NodeId>* mynodes)
282 {
283  minLevel = maxLevel = usedLevel = 0;
284 
285  // Get the tree nodes which pass the cuts
286  // (according to the selected resolution strategy)
287  switch (resolution)
288  {
289  case FIXED:
290  {
291  usedLevel = tree.getLevel(fixedScale);
292  tree.getPassingNodes(usedLevel, peakSelect, mynodes);
293  }
294  break;
295 
296  case MAXIMALLY_STABLE:
297  {
298  const unsigned minStartingLevel = maxStableScale > 0.0 ?
299  tree.getLevel(maxStableScale) : 0;
300  const unsigned maxStartingLevel = minStableScale > 0.0 ?
301  tree.getLevel(minStableScale) : UINT_MAX;
302 
303  if (tree.stableClusterCount(
304  peakSelect, &minLevel, &maxLevel, stabilityAlpha,
305  minStartingLevel, maxStartingLevel))
306  {
307  usedLevel = (minLevel + maxLevel)/2;
308  tree.getPassingNodes(usedLevel, peakSelect, mynodes);
309  }
310  }
311  break;
312 
313  case GLOBALLY_ADAPTIVE:
314  {
315  const bool stable = tree.clusterCountLevels(
316  nClustersRequested, peakSelect, &minLevel, &maxLevel);
317  if (minLevel || maxLevel)
318  {
319  usedLevel = (minLevel + maxLevel)/2;
320  if (!stable)
321  {
322  const int maxlev = tree.maxStoredLevel();
323  bool levelFound = false;
324  for (int delta=0; delta<=maxlev && !levelFound; ++delta)
325  for (int ifac=1; ifac>-2 && !levelFound; ifac-=2)
326  {
327  const int level = usedLevel + ifac*delta;
328  if (level > 0 && level <= maxlev)
329  if (occupancy[level] == nClustersRequested)
330  {
331  usedLevel = level;
332  levelFound = true;
333  }
334  }
335  assert(levelFound);
336  }
337  }
338  else
339  {
340  // Can't find that exact number of preclusters.
341  // Try to get the next best thing.
342  usedLevel = 1;
343  const unsigned occ1 = occupancy[1];
344  if (nClustersRequested >= occ1)
345  {
346  const unsigned maxlev = tree.maxStoredLevel();
347  if (nClustersRequested > occupancy[maxlev])
348  usedLevel = maxlev;
349  else
350  {
351  // It would be nice to use "lower_bound" here,
352  // but the occupancy is not necessarily monotonous.
353  unsigned bestDelta = nClustersRequested > occ1 ?
354  nClustersRequested - occ1 : occ1 - nClustersRequested;
355  for (unsigned level=2; level<=maxlev; ++level)
356  {
357  const unsigned n = occupancy[level];
358  const unsigned d = nClustersRequested > n ?
359  nClustersRequested - n : n - nClustersRequested;
360  if (d < bestDelta)
361  {
362  bestDelta = d;
363  usedLevel = level;
364  }
365  }
366  }
367  }
368  }
369  tree.getPassingNodes(usedLevel, peakSelect, mynodes);
370  }
371  break;
372 
373  case LOCALLY_ADAPTIVE:
374  {
375  usedLevel = tree.getLevel(fixedScale);
376  tree.getMagS2OptimalNodes(peakSelect, nClustersRequested,
377  usedLevel, mynodes, &thresholds);
378  }
379  break;
380 
381  default:
382  assert(!"ERROR in FFTJetProducer::selectTreeNodes : "
383  "should never get here! This is a bug. Please report.");
384  }
385 }
386 
387 
389 {
390  const unsigned nClus = preclusters.size();
391  if (nClus)
392  {
393  fftjet::Peak* clus = &preclusters[0];
394  fftjet::Functor1<double,fftjet::Peak>&
395  scaleCalc(*recoScaleCalcPeak);
396  fftjet::Functor1<double,fftjet::Peak>&
397  ratioCalc(*recoScaleRatioCalcPeak);
398  fftjet::Functor1<double,fftjet::Peak>&
399  factorCalc(*memberFactorCalcPeak);
400 
401  for (unsigned i=0; i<nClus; ++i)
402  {
403  clus[i].setRecoScale(scaleCalc(clus[i]));
404  clus[i].setRecoScaleRatio(ratioCalc(clus[i]));
405  clus[i].setMembershipFactor(factorCalc(clus[i]));
406  }
407  }
408 }
409 
410 
412 {
413  int minBin = energyFlow->getEtaBin(-gridScanMaxEta);
414  if (minBin < 0)
415  minBin = 0;
416  int maxBin = energyFlow->getEtaBin(gridScanMaxEta) + 1;
417  if (maxBin < 0)
418  maxBin = 0;
419 
420  fftjet::DefaultRecombinationAlgFactory<
421  Real,VectorLike,BgData,VBuilder> factory;
422  if (factory[recombinationAlgorithm] == nullptr)
423  throw cms::Exception("FFTJetBadConfig")
424  << "Invalid grid recombination algorithm \""
425  << recombinationAlgorithm << "\"" << std::endl;
426  gridAlg = std::auto_ptr<GridAlg>(
427  factory[recombinationAlgorithm]->create(
428  jetMembershipFunction.get(),
429  bgMembershipFunction.get(),
431  isCrisp, false, assignConstituents, minBin, maxBin));
432 }
433 
434 
436  const edm::Event& iEvent,
437  std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >& flow)
438 {
440  iEvent.getByToken(input_energyflow_token_, input);
441 
442  // Make sure that the grid is compatible with the stored one
443  bool rebuildGrid = flow.get() == nullptr;
444  if (!rebuildGrid)
445  rebuildGrid =
446  !(flow->nEta() == input->nEtaBins() &&
447  flow->nPhi() == input->nPhiBins() &&
448  flow->etaMin() == input->etaMin() &&
449  flow->etaMax() == input->etaMax() &&
450  flow->phiBin0Edge() == input->phiBin0Edge());
451  if (rebuildGrid)
452  {
453  // We should not get here very often...
454  flow = std::auto_ptr<fftjet::Grid2d<Real> >(
455  new fftjet::Grid2d<Real>(
456  input->nEtaBins(), input->etaMin(), input->etaMax(),
457  input->nPhiBins(), input->phiBin0Edge(), input->title()));
458  }
459  flow->blockSet(input->data(), input->nEtaBins(), input->nPhiBins());
460  return rebuildGrid;
461 }
462 
463 
464 bool FFTJetProducer::checkConvergence(const std::vector<RecoFFTJet>& previous,
465  std::vector<RecoFFTJet>& nextSet)
466 {
467  fftjet::Functor2<double,RecoFFTJet,RecoFFTJet>&
468  distanceCalc(*jetDistanceCalc);
469 
470  const unsigned nJets = previous.size();
471  if (nJets != nextSet.size())
472  return false;
473 
474  const RecoFFTJet* prev = &previous[0];
475  RecoFFTJet* next = &nextSet[0];
476 
477  // Calculate convergence distances for all jets
478  bool converged = true;
479  for (unsigned i=0; i<nJets; ++i)
480  {
481  const double d = distanceCalc(prev[i], next[i]);
482  next[i].setConvergenceDistance(d);
483  if (i < nJetsRequiredToConverge && d > convergenceDistance)
484  converged = false;
485  }
486 
487  return converged;
488 }
489 
490 
492 {
493  fftjet::Functor1<double,RecoFFTJet>& scaleCalc(*recoScaleCalcJet);
494  fftjet::Functor1<double,RecoFFTJet>& ratioCalc(*recoScaleRatioCalcJet);
495  fftjet::Functor1<double,RecoFFTJet>& factorCalc(*memberFactorCalcJet);
496 
497  unsigned nJets = recoJets.size();
498  unsigned iterNum = 1U;
499  bool converged = false;
500  for (; iterNum<maxIterations && !converged; ++iterNum)
501  {
502  // Recreate the vector of preclusters using the jets
503  const RecoFFTJet* jets = &recoJets[0];
504  iterPreclusters.clear();
505  iterPreclusters.reserve(nJets);
506  for (unsigned i=0; i<nJets; ++i)
507  {
508  const RecoFFTJet& jet(jets[i]);
509  fftjet::Peak p(jet.precluster());
510  p.setEtaPhi(jet.vec().Eta(), jet.vec().Phi());
511  p.setRecoScale(scaleCalc(jet));
512  p.setRecoScaleRatio(ratioCalc(jet));
513  p.setMembershipFactor(factorCalc(jet));
514  iterPreclusters.push_back(p);
515  }
516 
517  // Run the algorithm
518  int status = 0;
520  status = gridAlg->run(iterPreclusters, *energyFlow,
521  &noiseLevel, 1U, 1U,
523  else
524  status = recoAlg->run(iterPreclusters, eventData, &noiseLevel, 1U,
526  if (status)
527  throw cms::Exception("FFTJetInterface")
528  << "FFTJet algorithm failed" << std::endl;
529 
530  // As it turns out, it is possible, in very rare cases,
531  // to have iterJets.size() != nJets at this point
532 
533  // Figure out if the iterations have converged
534  converged = checkConvergence(recoJets, iterJets);
535 
536  // Prepare for the next cycle
537  iterJets.swap(recoJets);
538  nJets = recoJets.size();
539  }
540 
541  // Check that we have the correct number of preclusters
542  if (preclusters.size() != nJets)
543  {
544  assert(nJets < preclusters.size());
546  assert(preclusters.size() == nJets);
547  }
548 
549  // Plug in the original precluster coordinates into the result
550  RecoFFTJet* jets = &recoJets[0];
551  for (unsigned i=0; i<nJets; ++i)
552  {
553  const fftjet::Peak& oldp(preclusters[i]);
554  jets[i].setPeakEtaPhi(oldp.eta(), oldp.phi());
555  }
556 
557  // If we have converged on the last cycle, the result
558  // would be indistinguishable from no convergence.
559  // Because of this, raise the counter by one to indicate
560  // the case when the convergence is not achieved.
561  if (!converged)
562  ++iterNum;
563 
564  return iterNum;
565 }
566 
567 
569 {
570  const unsigned nJets = recoJets.size();
571  const unsigned* clusterMask = gridAlg->getClusterMask();
572  const int nEta = gridAlg->getLastNEta();
573  const int nPhi = gridAlg->getLastNPhi();
574  const fftjet::Grid2d<Real>& g(*energyFlow);
575 
576  const unsigned nInputs = eventData.size();
577  const VectorLike* inp = nInputs ? &eventData[0] : nullptr;
578  const unsigned* candIdx = nInputs ? &candidateIndex[0] : nullptr;
579  for (unsigned i=0; i<nInputs; ++i)
580  {
581  const VectorLike& item(inp[i]);
582  const int iPhi = g.getPhiBin(item.Phi());
583  const int iEta = g.getEtaBin(item.Eta());
584  const unsigned mask = iEta >= 0 && iEta < nEta ?
585  clusterMask[iEta*nPhi + iPhi] : 0;
586  assert(mask <= nJets);
587  constituents[mask].push_back(inputCollection->ptrAt(candIdx[i]));
588  }
589 }
590 
591 
593 {
594  const unsigned nJets = recoJets.size();
595  const unsigned* clusterMask = recoAlg->getClusterMask();
596  const unsigned maskLength = recoAlg->getLastNData();
597  assert(maskLength == eventData.size());
598 
599  const unsigned* candIdx = maskLength ? &candidateIndex[0] : nullptr;
600  for (unsigned i=0; i<maskLength; ++i)
601  {
602  // In FFTJet, the mask value of 0 corresponds to unclustered
603  // energy. We will do the same here. Jet numbers are therefore
604  // shifted by 1 wrt constituents vector, and constituents[1]
605  // corresponds to jet number 0.
606  const unsigned mask = clusterMask[i];
607  assert(mask <= nJets);
608  constituents[mask].push_back(inputCollection->ptrAt(candIdx[i]));
609  }
610 }
611 
612 
613 // The following code more or less coincides with the similar method
614 // implemented in VirtualJetProducer
615 template <typename T>
617  const edm::EventSetup& iSetup)
618 {
619  using namespace reco;
620 
621  typedef FFTAnyJet<T> OutputJet;
622  typedef std::vector<OutputJet> OutputCollection;
623 
624  // Area of a single eta-phi cell for jet area calculations.
625  // Set it to 0 in case the module configuration does not allow
626  // us to calculate jet areas reliably.
627  double cellArea = useGriddedAlgorithm &&
629  energyFlow->etaBinWidth() * energyFlow->phiBinWidth() : 0.0;
630 
631  if (calculatePileup)
632  cellArea = pileupEnergyFlow->etaBinWidth() *
633  pileupEnergyFlow->phiBinWidth();
634 
635  // allocate output jet collection
636  auto jets = std::make_unique<OutputCollection>();
637  const unsigned nJets = recoJets.size();
638  jets->reserve(nJets);
639 
640  bool sorted = true;
641  double previousPt = DBL_MAX;
642  for (unsigned ijet=0; ijet<nJets; ++ijet)
643  {
644  RecoFFTJet& myjet(recoJets[ijet]);
645 
646  // Check if we should resum jet constituents
647  VectorLike jet4vec(myjet.vec());
648  if (resumConstituents)
649  {
650  VectorLike sum(0.0, 0.0, 0.0, 0.0);
651  const unsigned nCon = constituents[ijet+1].size();
652  const reco::CandidatePtr* cn = nCon ? &constituents[ijet+1][0] : nullptr;
653  for (unsigned i=0; i<nCon; ++i)
654  sum += cn[i]->p4();
655  jet4vec = sum;
656  setJetStatusBit(&myjet, CONSTITUENTS_RESUMMED, true);
657  }
658 
659  // Subtract the pile-up
661  {
662  jet4vec = adjustForPileup(jet4vec, pileup[ijet],
666  else
667  setJetStatusBit(&myjet, PILEUP_SUBTRACTED_PT, true);
668  }
669 
670  // Write the specifics to the jet (simultaneously sets 4-vector,
671  // vertex, constituents). These are overridden functions that will
672  // call the appropriate specific code.
673  T jet;
674  writeSpecific(jet, jet4vec, vertexUsed(),
675  constituents[ijet+1], iSetup);
676 
677  // calcuate the jet area
678  double ncells = myjet.ncells();
679  if (calculatePileup)
680  {
681  ncells = cellCountsVec[ijet];
682  setJetStatusBit(&myjet, PILEUP_CALCULATED, true);
683  }
684  jet.setJetArea(cellArea*ncells);
685 
686  // add jet to the list
687  FFTJet<float> fj(jetToStorable<float>(myjet));
688  fj.setFourVec(jet4vec);
689  if (calculatePileup)
690  {
691  fj.setPileup(pileup[ijet]);
692  fj.setNCells(ncells);
693  }
694  jets->push_back(OutputJet(jet, fj));
695 
696  // Check whether the sequence remains sorted by pt
697  const double pt = jet.pt();
698  if (pt > previousPt)
699  sorted = false;
700  previousPt = pt;
701  }
702 
703  // Sort the collection
704  if (!sorted)
705  std::sort(jets->begin(), jets->end(), LocalSortByPt());
706 
707  // put the collection into the event
708  iEvent.put(std::move(jets), outputLabel);
709 }
710 
711 
713  const unsigned nPreclustersFound)
714 {
715  // Write recombined jets
716  jet_type_switch(writeJets, ev, iSetup);
717 
718  // Check if we should resum unclustered energy constituents
719  VectorLike unclusE(unclustered);
720  if (resumConstituents)
721  {
722  VectorLike sum(0.0, 0.0, 0.0, 0.0);
723  const unsigned nCon = constituents[0].size();
724  const reco::CandidatePtr* cn = nCon ? &constituents[0][0] : nullptr;
725  for (unsigned i=0; i<nCon; ++i)
726  sum += cn[i]->p4();
727  unclusE = sum;
728  }
729 
730  // Write the jet reconstruction summary
731  const double minScale = minLevel ? sparseTree.getScale(minLevel) : 0.0;
732  const double maxScale = maxLevel ? sparseTree.getScale(maxLevel) : 0.0;
733  const double scaleUsed = usedLevel ? sparseTree.getScale(usedLevel) : 0.0;
734 
735  ev.put(std::make_unique<reco::FFTJetProducerSummary>(
736  thresholds, occupancy, unclusE,
737  constituents[0], unused,
738  minScale, maxScale, scaleUsed,
739  nPreclustersFound, iterationsPerformed,
740  iterationsPerformed == 1U ||
742 }
743 
744 
745 // ------------ method called to for each event ------------
747  const edm::EventSetup& iSetup)
748 {
749  // Load the clustering tree made by FFTJetPatRecoProducer
751  loadSparseTreeData<float>(iEvent);
752  else
753  loadSparseTreeData<double>(iEvent);
754 
755  // Do we need to load the candidate collection?
757  loadInputCollection(iEvent);
758 
759  // Do we need to have discretized energy flow?
761  {
762  if (reuseExistingGrid)
763  {
764  if (loadEnergyFlow(iEvent, energyFlow))
765  buildGridAlg();
766  }
767  else
769  }
770 
771  // Calculate cluster occupancy as a function of level number
772  sparseTree.occupancyInScaleSpace(*peakSelector, &occupancy);
773 
774  // Select the preclusters using the requested resolution scheme
775  preclusters.clear();
776  if (resolution == FROM_GENJETS)
777  genJetPreclusters(sparseTree, iEvent, iSetup,
779  else
781  if (preclusters.size() > maxInitialPreclusters)
782  {
783  std::sort(preclusters.begin(), preclusters.end(), std::greater<fftjet::Peak>());
785  }
786 
787  // Prepare to run the jet recombination procedure
789 
790  // Assign membership functions to preclusters. If this function
791  // is not overriden in a derived class, default algorithm membership
792  // function will be used for every cluster.
794 
795  // Count the preclusters going in
796  unsigned nPreclustersFound = 0U;
797  const unsigned npre = preclusters.size();
798  for (unsigned i=0; i<npre; ++i)
799  if (preclusters[i].membershipFactor() > 0.0)
800  ++nPreclustersFound;
801 
802  // Run the recombination algorithm once
803  int status = 0;
805  status = gridAlg->run(preclusters, *energyFlow,
806  &noiseLevel, 1U, 1U,
808  else
809  status = recoAlg->run(preclusters, eventData, &noiseLevel, 1U,
811  if (status)
812  throw cms::Exception("FFTJetInterface")
813  << "FFTJet algorithm failed (first iteration)" << std::endl;
814 
815  // If requested, iterate the jet recombination procedure
816  if (maxIterations > 1U && !recoJets.empty())
817  {
818  // It is possible to have a smaller number of jets than we had
819  // preclusters. Fake preclusters are possible, but for a good
820  // choice of pattern recognition kernel their presence should
821  // be infrequent. However, any fake preclusters will throw the
822  // iterative reconstruction off balance. Deal with the problem now.
823  const unsigned nJets = recoJets.size();
824  if (preclusters.size() != nJets)
825  {
826  assert(nJets < preclusters.size());
828  }
830  }
831  else
833 
834  // Determine jet constituents. FFTJet returns a map
835  // of constituents which is inverse to what we need here.
836  const unsigned nJets = recoJets.size();
837  if (constituents.size() <= nJets)
838  constituents.resize(nJets + 1U);
839  if (assignConstituents)
840  {
841  for (unsigned i=0; i<=nJets; ++i)
842  constituents[i].clear();
845  else
847  }
848 
849  // Figure out the pile-up
850  if (calculatePileup)
851  {
852  if (loadPileupFromDB)
853  determinePileupDensityFromDB(iEvent, iSetup,
855  else
858  determinePileup();
859  assert(pileup.size() == recoJets.size());
860  }
861 
862  // Write out the results
863  saveResults(iEvent, iSetup, nPreclustersFound);
864 }
865 
866 
867 std::auto_ptr<fftjet::Functor1<bool,fftjet::Peak> >
869 {
871  ps.getParameter<edm::ParameterSet>("PeakSelectorConfiguration"));
872 }
873 
874 
875 // Parser for the jet membership function
876 std::auto_ptr<fftjet::ScaleSpaceKernel>
878 {
880  ps.getParameter<edm::ParameterSet>("jetMembershipFunction"));
881 }
882 
883 
884 // Parser for the background membership function
885 std::auto_ptr<AbsBgFunctor>
887 {
889  ps.getParameter<edm::ParameterSet>("bgMembershipFunction"));
890 }
891 
892 
893 // Calculator for the recombination scale
894 std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
896 {
898  ps.getParameter<edm::ParameterSet>("recoScaleCalcPeak"));
899 }
900 
901 
902 // Pile-up density calculator
903 std::auto_ptr<fftjetcms::AbsPileupCalculator>
905 {
907  ps.getParameter<edm::ParameterSet>("pileupDensityCalc"));
908 }
909 
910 
911 // Calculator for the recombination scale ratio
912 std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
914 {
916  ps.getParameter<edm::ParameterSet>("recoScaleRatioCalcPeak"));
917 }
918 
919 
920 // Calculator for the membership function factor
921 std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
923 {
925  ps.getParameter<edm::ParameterSet>("memberFactorCalcPeak"));
926 }
927 
928 
929 std::auto_ptr<fftjet::Functor1<double,FFTJetProducer::RecoFFTJet> >
931 {
933  ps.getParameter<edm::ParameterSet>("recoScaleCalcJet"));
934 }
935 
936 
937 std::auto_ptr<fftjet::Functor1<double,FFTJetProducer::RecoFFTJet> >
939 {
941  ps.getParameter<edm::ParameterSet>("recoScaleRatioCalcJet"));
942 }
943 
944 
945 std::auto_ptr<fftjet::Functor1<double,FFTJetProducer::RecoFFTJet> >
947 {
949  ps.getParameter<edm::ParameterSet>("memberFactorCalcJet"));
950 }
951 
952 
953 std::auto_ptr<fftjet::Functor2<
956 {
958  ps.getParameter<edm::ParameterSet>("jetDistanceCalc"));
959 }
960 
961 
962 void FFTJetProducer::assignMembershipFunctions(std::vector<fftjet::Peak>*)
963 {
964 }
965 
966 
967 // ------------ method called once each job just before starting event loop
969 {
971 
972  // Parse the peak selector definition
974  checkConfig(peakSelector, "invalid peak selector");
975 
977  checkConfig(jetMembershipFunction, "invalid jet membership function");
978 
980  checkConfig(bgMembershipFunction, "invalid noise membership function");
981 
982  // Build the energy recombination algorithm
983  if (!useGriddedAlgorithm)
984  {
985  fftjet::DefaultVectorRecombinationAlgFactory<
986  VectorLike,BgData,VBuilder> factory;
987  if (factory[recombinationAlgorithm] == nullptr)
988  throw cms::Exception("FFTJetBadConfig")
989  << "Invalid vector recombination algorithm \""
990  << recombinationAlgorithm << "\"" << std::endl;
991  recoAlg = std::auto_ptr<RecoAlg>(
992  factory[recombinationAlgorithm]->create(
993  jetMembershipFunction.get(),
994  &VectorLike::Et, &VectorLike::Eta, &VectorLike::Phi,
995  bgMembershipFunction.get(),
997  }
998  else if (!reuseExistingGrid)
999  {
1001  ps.getParameter<edm::ParameterSet>("GridConfiguration"));
1002  checkConfig(energyFlow, "invalid discretization grid");
1003  buildGridAlg();
1004  }
1005 
1006  // Create the grid subsequently used for pile-up subtraction
1007  if (calculatePileup)
1008  {
1010  ps.getParameter<edm::ParameterSet>("PileupGridConfiguration"));
1011  checkConfig(pileupEnergyFlow, "invalid pileup density grid");
1012 
1013  if (!loadPileupFromDB)
1014  {
1016  checkConfig(pileupDensityCalc, "invalid pile-up density calculator");
1017  }
1018  }
1019 
1020  // Parse the calculator of the recombination scale
1022  checkConfig(recoScaleCalcPeak, "invalid spec for the "
1023  "reconstruction scale calculator from peaks");
1024 
1025  // Parse the calculator of the recombination scale ratio
1027  checkConfig(recoScaleRatioCalcPeak, "invalid spec for the "
1028  "reconstruction scale ratio calculator from peaks");
1029 
1030  // Calculator for the membership function factor
1032  checkConfig(memberFactorCalcPeak, "invalid spec for the "
1033  "membership function factor calculator from peaks");
1034 
1035  if (maxIterations > 1)
1036  {
1037  // We are going to run iteratively. Make required objects.
1039  checkConfig(recoScaleCalcJet, "invalid spec for the "
1040  "reconstruction scale calculator from jets");
1041 
1043  checkConfig(recoScaleRatioCalcJet, "invalid spec for the "
1044  "reconstruction scale ratio calculator from jets");
1045 
1047  checkConfig(memberFactorCalcJet, "invalid spec for the "
1048  "membership function factor calculator from jets");
1049 
1051  checkConfig(memberFactorCalcJet, "invalid spec for the "
1052  "jet distance calculator");
1053  }
1054 }
1055 
1056 
1058 {
1059  // There are two possible reasons for fake preclusters:
1060  // 1. Membership factor was set to 0
1061  // 2. Genuine problem with pattern recognition
1062  //
1063  // Anyway, we need to match jets to preclusters and keep
1064  // only those preclusters that have been matched
1065  //
1066  std::vector<int> matchTable;
1067  const unsigned nmatched = matchOneToOne(
1068  recoJets, preclusters, JetToPeakDistance(), &matchTable);
1069 
1070  // Ensure that all jets have been matched.
1071  // If not, we must have a bug somewhere.
1072  assert(nmatched == recoJets.size());
1073 
1074  // Collect all matched preclusters
1075  iterPreclusters.clear();
1076  iterPreclusters.reserve(nmatched);
1077  for (unsigned i=0; i<nmatched; ++i)
1078  iterPreclusters.push_back(preclusters[matchTable[i]]);
1080 }
1081 
1082 
1084  const int mask, const bool value)
1085 {
1086  int status = jet->status();
1087  if (value)
1088  status |= mask;
1089  else
1090  status &= ~mask;
1091  jet->setStatus(status);
1092 }
1093 
1094 
1096  const edm::Event& iEvent,
1097  std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >& density)
1098 {
1100  iEvent.getByToken(input_pusummary_token_, summary);
1101 
1102  const reco::FFTJetPileupSummary& s(*summary);
1103  const AbsPileupCalculator& calc(*pileupDensityCalc);
1104  const bool phiDependent = calc.isPhiDependent();
1105 
1106  fftjet::Grid2d<Real>& g(*density);
1107  const unsigned nEta = g.nEta();
1108  const unsigned nPhi = g.nPhi();
1109 
1110  for (unsigned ieta=0; ieta<nEta; ++ieta)
1111  {
1112  const double eta(g.etaBinCenter(ieta));
1113 
1114  if (phiDependent)
1115  {
1116  for (unsigned iphi=0; iphi<nPhi; ++iphi)
1117  {
1118  const double phi(g.phiBinCenter(iphi));
1119  g.uncheckedSetBin(ieta, iphi, calc(eta, phi, s));
1120  }
1121  }
1122  else
1123  {
1124  const double pil = calc(eta, 0.0, s);
1125  for (unsigned iphi=0; iphi<nPhi; ++iphi)
1126  g.uncheckedSetBin(ieta, iphi, pil);
1127  }
1128  }
1129 }
1130 
1131 
1133  const edm::Event& iEvent, const edm::EventSetup& iSetup,
1134  std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >& density)
1135 {
1138  iSetup, pileupTableRecord, h);
1139  boost::shared_ptr<npstat::StorableMultivariateFunctor> f =
1141 
1143  iEvent.getByToken(input_pusummary_token_, summary);
1144 
1145  const float rho = summary->pileupRho();
1146  const bool phiDependent = f->minDim() == 3U;
1147 
1148  fftjet::Grid2d<Real>& g(*density);
1149  const unsigned nEta = g.nEta();
1150  const unsigned nPhi = g.nPhi();
1151 
1152  double functorArg[3] = {0.0, 0.0, 0.0};
1153  if (phiDependent)
1154  functorArg[2] = rho;
1155  else
1156  functorArg[1] = rho;
1157 
1158  for (unsigned ieta=0; ieta<nEta; ++ieta)
1159  {
1160  const double eta(g.etaBinCenter(ieta));
1161  functorArg[0] = eta;
1162 
1163  if (phiDependent)
1164  {
1165  for (unsigned iphi=0; iphi<nPhi; ++iphi)
1166  {
1167  functorArg[1] = g.phiBinCenter(iphi);
1168  g.uncheckedSetBin(ieta, iphi, (*f)(functorArg, 3U));
1169  }
1170  }
1171  else
1172  {
1173  const double pil = (*f)(functorArg, 2U);
1174  for (unsigned iphi=0; iphi<nPhi; ++iphi)
1175  g.uncheckedSetBin(ieta, iphi, pil);
1176  }
1177  }
1178 }
1179 
1180 
1182 {
1183  // This function works with crisp clustering only
1184  if (!isCrisp)
1185  assert(!"Pile-up subtraction for fuzzy clustering "
1186  "is not implemented yet");
1187 
1188  // Clear the pileup vector
1189  const unsigned nJets = recoJets.size();
1190  pileup.resize(nJets);
1191  if (nJets == 0)
1192  return;
1193  const VectorLike zero;
1194  for (unsigned i=0; i<nJets; ++i)
1195  pileup[i] = zero;
1196 
1197  // Pileup energy flow grid
1198  const fftjet::Grid2d<Real>& g(*pileupEnergyFlow);
1199  const unsigned nEta = g.nEta();
1200  const unsigned nPhi = g.nPhi();
1201  const double cellArea = g.etaBinWidth() * g.phiBinWidth();
1202 
1203  // Various calculators
1204  fftjet::Functor1<double,RecoFFTJet>& scaleCalc(*recoScaleCalcJet);
1205  fftjet::Functor1<double,RecoFFTJet>& ratioCalc(*recoScaleRatioCalcJet);
1206  fftjet::Functor1<double,RecoFFTJet>& factorCalc(*memberFactorCalcJet);
1207 
1208  // Make sure we have enough memory
1209  memFcns2dVec.resize(nJets);
1210  fftjet::AbsKernel2d** memFcns2d = &memFcns2dVec[0];
1211 
1212  doubleBuf.resize(nJets*4U + nJets*nPhi);
1213  double* recoScales = &doubleBuf[0];
1214  double* recoScaleRatios = recoScales + nJets;
1215  double* memFactors = recoScaleRatios + nJets;
1216  double* dEta = memFactors + nJets;
1217  double* dPhiBuf = dEta + nJets;
1218 
1219  cellCountsVec.resize(nJets);
1220  unsigned* cellCounts = &cellCountsVec[0];
1221 
1222  // Go over jets and collect the necessary info
1223  for (unsigned ijet=0; ijet<nJets; ++ijet)
1224  {
1225  const RecoFFTJet& jet(recoJets[ijet]);
1226  const fftjet::Peak& peak(jet.precluster());
1227 
1228  // Make sure we are using 2-d membership functions.
1229  // Pile-up subtraction scheme for 3-d functions should be different.
1230  fftjet::AbsMembershipFunction* m3d =
1231  dynamic_cast<fftjet::AbsMembershipFunction*>(
1232  peak.membershipFunction());
1233  if (m3d == nullptr)
1234  m3d = dynamic_cast<fftjet::AbsMembershipFunction*>(
1235  jetMembershipFunction.get());
1236  if (m3d)
1237  {
1238  assert(!"Pile-up subtraction for 3-d membership functions "
1239  "is not implemented yet");
1240  }
1241  else
1242  {
1243  fftjet::AbsKernel2d* m2d =
1244  dynamic_cast<fftjet::AbsKernel2d*>(
1245  peak.membershipFunction());
1246  if (m2d == nullptr)
1247  m2d = dynamic_cast<fftjet::AbsKernel2d*>(
1248  jetMembershipFunction.get());
1249  assert(m2d);
1250  memFcns2d[ijet] = m2d;
1251  }
1252  recoScales[ijet] = scaleCalc(jet);
1253  recoScaleRatios[ijet] = ratioCalc(jet);
1254  memFactors[ijet] = factorCalc(jet);
1255  cellCounts[ijet] = 0U;
1256 
1257  const double jetPhi = jet.vec().Phi();
1258  for (unsigned iphi=0; iphi<nPhi; ++iphi)
1259  {
1260  double dphi = g.phiBinCenter(iphi) - jetPhi;
1261  while (dphi > M_PI)
1262  dphi -= (2.0*M_PI);
1263  while (dphi < -M_PI)
1264  dphi += (2.0*M_PI);
1265  dPhiBuf[iphi*nJets+ijet] = dphi;
1266  }
1267  }
1268 
1269  // Go over all grid points and integrate
1270  // the pile-up energy density
1271  VBuilder vMaker;
1272  for (unsigned ieta=0; ieta<nEta; ++ieta)
1273  {
1274  const double eta(g.etaBinCenter(ieta));
1275  const Real* databuf = g.data() + ieta*nPhi;
1276 
1277  // Figure out dEta for each jet
1278  for (unsigned ijet=0; ijet<nJets; ++ijet)
1279  dEta[ijet] = eta - recoJets[ijet].vec().Eta();
1280 
1281  for (unsigned iphi=0; iphi<nPhi; ++iphi)
1282  {
1283  double maxW(0.0);
1284  unsigned maxWJet(nJets);
1285  const double* dPhi = dPhiBuf + iphi*nJets;
1286 
1287  for (unsigned ijet=0; ijet<nJets; ++ijet)
1288  {
1289  if (recoScaleRatios[ijet] > 0.0)
1290  memFcns2d[ijet]->setScaleRatio(recoScaleRatios[ijet]);
1291  const double f = memFactors[ijet]*
1292  (*memFcns2d[ijet])(dEta[ijet], dPhi[ijet],
1293  recoScales[ijet]);
1294  if (f > maxW)
1295  {
1296  maxW = f;
1297  maxWJet = ijet;
1298  }
1299  }
1300 
1301  if (maxWJet < nJets)
1302  {
1303  pileup[maxWJet] += vMaker(cellArea*databuf[iphi],
1304  eta, g.phiBinCenter(iphi));
1305  cellCounts[maxWJet]++;
1306  }
1307  }
1308  }
1309 }
1310 
1311 
1312 // ------------ method called once each job just after ending the event loop
1314 {
1315 }
1316 
1317 
1318 //define this as a plug-in
std::vector< std::vector< reco::CandidatePtr > > constituents
dbl * delta
Definition: mlp_gen.cc:36
void produce(edm::Event &, const edm::EventSetup &) override
unsigned matchOneToOne(const std::vector< T1 > &v1, const std::vector< T2 > &v2, const DistanceCalculator &calc, std::vector< int > *matchFrom1To2, const double maxMatchingDistance=1.0e300)
Definition: matchOneToOne.h:41
T getParameter(std::string const &) const
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > pileupEnergyFlow
T getUntrackedParameter(std::string const &, T const &) const
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
std::auto_ptr< RecoAlg > recoAlg
std::auto_ptr< fftjet::Functor1< double, fftjet::Peak > > memberFactorCalcPeak
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
unsigned iterationsPerformed
const bool useGriddedAlgorithm
std::vector< fftjetcms::VectorLike > pileup
edm::Handle< reco::CandidateView > inputCollection
static Resolution parse_resolution(const std::string &name)
const bool resumConstituents
void selectTreeNodes(const SparseTree &tree, const fftjet::Functor1< bool, fftjet::Peak > &peakSelect, std::vector< SparseTree::NodeId > *nodes)
Class for storing FFTJet sparse clustering trees.
Definition: PattRecoTree.h:20
const double gridScanMaxEta
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
void saveResults(edm::Event &iEvent, const edm::EventSetup &, unsigned nPreclustersFound)
std::auto_ptr< fftjet::Grid2d< Real > > fftjet_Grid2d_parser(const edm::ParameterSet &ps)
virtual std::auto_ptr< fftjet::Functor1< double, fftjet::Peak > > parse_recoScaleRatioCalcPeak(const edm::ParameterSet &)
std::auto_ptr< GridAlg > gridAlg
std::auto_ptr< fftjetcms::AbsPileupCalculator > pileupDensityCalc
virtual void assignMembershipFunctions(std::vector< fftjet::Peak > *preclusters)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
Ptr< value_type > ptrAt(size_type i) const
const reco::Particle::Point & vertexUsed() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void loadInputCollection(const edm::Event &)
virtual std::auto_ptr< fftjet::Functor1< bool, fftjet::Peak > > parse_peakSelector(const edm::ParameterSet &)
virtual std::auto_ptr< fftjetcms::AbsPileupCalculator > parse_pileupDensityCalc(const edm::ParameterSet &ps)
double BgData
const std::string recombinationAlgorithm
static const int stable
Definition: TopGenEvent.h:11
const double fixedScale
virtual std::auto_ptr< fftjetcms::AbsBgFunctor > parse_bgMembershipFunction(const edm::ParameterSet &)
const bool subtractPileupAs4Vec
fftjet::RecombinedJet< VectorLike > jetFromStorable(const reco::FFTJet< Real > &jet)
Definition: jetConverters.h:63
std::auto_ptr< fftjet::ScaleSpaceKernel > fftjet_MembershipFunction_parser(const edm::ParameterSet &ps)
void checkConfig(const Ptr &ptr, const char *message)
const double unlikelyBgWeight
bool ev
void makeProduces(const std::string &alias, const std::string &tag)
const bool subtractPileup
Resolution resolution
std::auto_ptr< AbsPileupCalculator > fftjet_PileupCalculator_parser(const edm::ParameterSet &ps)
edm::EDGetTokenT< reco::PattRecoTree< fftjetcms::Real, reco::PattRecoPeak< fftjetcms::Real > > > input_recotree_token_
std::vector< double > doubleBuf
std::auto_ptr< fftjet::Functor1< double, fftjet::Peak > > recoScaleCalcPeak
std::vector< fftjet::AbsKernel2d * > memFcns2dVec
std::vector< RecoFFTJet > recoJets
const bool isCrisp
const double recombinationDataCutoff
Summary info for pile-up determined by Gaussian filtering.
std::vector< unsigned > occupancy
const double minStableScale
const double maxStableScale
#define init_param(type, varname)
std::vector< fftjet::Peak > preclusters
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
void determineVectorConstituents()
FFTJetProducer()=delete
static std::string const input
Definition: EdmProvDump.cc:44
virtual std::auto_ptr< fftjet::Functor1< double, fftjet::Peak > > parse_memberFactorCalcPeak(const edm::ParameterSet &)
const double * data() const
void loadSparseTreeData(const edm::Event &)
std::auto_ptr< fftjetcms::AbsBgFunctor > bgMembershipFunction
const bool reuseExistingGrid
std::auto_ptr< fftjet::Functor1< bool, fftjet::Peak > > fftjet_PeakSelector_parser(const edm::ParameterSet &ps)
const unsigned nClustersRequested
int iEvent
Definition: GenABIO.cc:230
std::vector< unsigned > candidateIndex
std::auto_ptr< fftjet::Functor2< double, fftjet::RecombinedJet< VectorLike >, fftjet::RecombinedJet< VectorLike > > > fftjet_JetDistance_parser(const edm::ParameterSet &ps)
Definition: Jet.py:1
std::auto_ptr< std::vector< double > > iniScales
math::XYZTLorentzVector VectorLike
unsigned usedLevel
double p4[4]
Definition: TauolaWrapper.h:92
std::auto_ptr< fftjet::Functor1< double, RecoFFTJet > > memberFactorCalcJet
void clear(CLHEP::HepGenMatrix &m)
Helper function: Reset all elements of a matrix to 0.
Definition: matutil.cc:167
virtual std::auto_ptr< fftjet::Functor1< double, RecoFFTJet > > parse_recoScaleCalcJet(const edm::ParameterSet &)
vector< PseudoJet > jets
void endJob() override
Implements inheritance relationships for FFTJet jets.
Definition: FFTAnyJet.h:16
const double noiseLevel
static void setJetStatusBit(RecoFFTJet *jet, int mask, bool value)
std::vector< RecoFFTJet > iterJets
double f[11][100]
fftjet::SparseClusteringTree< fftjet::Peak, long > SparseTree
std::vector< unsigned > cellCountsVec
Definition: value.py:1
std::string pileupTableRecord
std::string pileupTableName
std::vector< double > thresholds
virtual std::auto_ptr< fftjet::Functor2< double, RecoFFTJet, RecoFFTJet > > parse_jetDistanceCalc(const edm::ParameterSet &)
unsigned iterateJetReconstruction()
void setNCells(const double nc)
Definition: FFTJet.h:84
const bool calculatePileup
const bool assignConstituents
bool storeInSinglePrecision() const
void setPileup(const math::XYZTLorentzVector &p)
Definition: FFTJet.h:80
#define M_PI
std::auto_ptr< fftjet::Functor1< double, fftjet::Peak > > fftjet_PeakFunctor_parser(const edm::ParameterSet &ps)
bool checkConvergence(const std::vector< RecoFFTJet > &previousIterResult, std::vector< RecoFFTJet > &thisIterResult)
const edm::InputTag treeLabel
virtual void selectPreclusters(const SparseTree &tree, const fftjet::Functor1< bool, fftjet::Peak > &peakSelector, std::vector< fftjet::Peak > *preclusters)
virtual void determinePileupDensityFromDB(const edm::Event &iEvent, const edm::EventSetup &iSetup, std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > &density)
double Real
void sparsePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::SparseClusteringTree< fftjet::Peak, long > *out)
const edm::InputTag genJetsLabel
const edm::ParameterSet myConfiguration
edm::EDGetTokenT< reco::DiscretizedEnergyFlow > input_energyflow_token_
fftjet::RecombinedJet< fftjetcms::VectorLike > RecoFFTJet
edm::EDGetTokenT< reco::FFTJetPileupSummary > input_pusummary_token_
virtual bool isPhiDependent() const =0
const double stabilityAlpha
std::auto_ptr< fftjet::Functor1< double, fftjet::Peak > > recoScaleRatioCalcPeak
virtual std::auto_ptr< fftjet::Functor1< double, RecoFFTJet > > parse_memberFactorCalcJet(const edm::ParameterSet &)
const unsigned maxInitialPreclusters
std::string pileupTableCategory
SparseTree sparseTree
fixed size matrix
HLT enums.
virtual std::auto_ptr< fftjet::Functor1< double, fftjet::Peak > > parse_recoScaleCalcPeak(const edm::ParameterSet &)
void setFourVec(const math::XYZTLorentzVector &p)
Definition: FFTJet.h:81
std::vector< SparseTree::NodeId > nodes
void determineGriddedConstituents()
const double convergenceDistance
edm::EDGetTokenT< std::vector< reco::FFTAnyJet< reco::GenJet > > > input_genjet_token_
void writeJets(edm::Event &iEvent, const edm::EventSetup &)
virtual std::auto_ptr< fftjet::ScaleSpaceKernel > parse_jetMembershipFunction(const edm::ParameterSet &)
~FFTJetProducer() override
virtual std::auto_ptr< fftjet::Functor1< double, RecoFFTJet > > parse_recoScaleRatioCalcJet(const edm::ParameterSet &)
std::auto_ptr< fftjet::Functor1< double, RecoFFTJet > > recoScaleRatioCalcJet
const edm::InputTag pileupLabel
void beginJob() override
#define jet_type_switch(method, arg1, arg2)
Definition: tree.py:1
bool loadEnergyFlow(const edm::Event &iEvent, std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > &flow)
virtual void genJetPreclusters(const SparseTree &tree, edm::Event &, const edm::EventSetup &, const fftjet::Functor1< bool, fftjet::Peak > &peakSelector, std::vector< fftjet::Peak > *preclusters)
std::auto_ptr< AbsBgFunctor > fftjet_BgFunctor_parser(const edm::ParameterSet &ps)
long double T
static const Mapper & instance()
std::vector< fftjetcms::VectorLike > eventData
const std::string outputLabel
math::XYZTLorentzVector adjustForPileup(const math::XYZTLorentzVector &jet, const math::XYZTLorentzVector &pileup, bool subtractPileupAs4Vec)
fftjetcms::VectorLike unclustered
def move(src, dest)
Definition: eostools.py:510
std::auto_ptr< fftjet::Functor1< bool, fftjet::Peak > > peakSelector
std::vector< fftjet::Peak > iterPreclusters
std::auto_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)
virtual void determinePileupDensityFromConfig(const edm::Event &iEvent, std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > &density)
std::auto_ptr< fftjet::Functor2< double, RecoFFTJet, RecoFFTJet > > jetDistanceCalc
std::auto_ptr< fftjet::ScaleSpaceKernel > jetMembershipFunction
std::auto_ptr< fftjet::Functor1< double, fftjet::RecombinedJet< VectorLike > > > fftjet_JetFunctor_parser(const edm::ParameterSet &ps)
void writeSpecific(reco::CaloJet &jet, reco::Particle::LorentzVector const &p4, reco::Particle::Point const &point, std::vector< reco::CandidatePtr > const &constituents, edm::EventSetup const &c)
Definition: JetSpecific.cc:41
void removeFakePreclusters()
const unsigned maxIterations
std::auto_ptr< fftjet::Functor1< double, RecoFFTJet > > recoScaleCalcJet
void prepareRecombinationScales()