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 <algorithm>
20 #include <fstream>
21 #include <functional>
22 #include <iostream>
23 #include <memory>
24 
25 #include "fftjet/peakEtLifetime.hh"
26 
27 // Header for this class
29 
30 // Framework include files
33 
34 // Data formats
39 
43 
44 // Loader for the lookup tables
46 
47 #define make_param(type, varname) const type& varname(ps.getParameter<type>(#varname))
48 
49 #define init_param(type, varname) varname(ps.getParameter<type>(#varname))
50 
51 // A generic switch statement based on jet type.
52 // Defining it in a single place avoids potential errors
53 // in case another jet type is introduced in the future.
54 //
55 // JPTJet is omitted for now: there is no reco::writeSpecific method
56 // for it (see header JetSpecific.h in the JetProducers package)
57 //
58 #define jet_type_switch(method, arg1, arg2) \
59  do { \
60  switch (jetType) { \
61  case CALOJET: \
62  method<reco::CaloJet>(arg1, arg2); \
63  break; \
64  case PFJET: \
65  method<reco::PFJet>(arg1, arg2); \
66  break; \
67  case GENJET: \
68  method<reco::GenJet>(arg1, arg2); \
69  break; \
70  case TRACKJET: \
71  method<reco::TrackJet>(arg1, arg2); \
72  break; \
73  case BASICJET: \
74  method<reco::BasicJet>(arg1, arg2); \
75  break; \
76  default: \
77  assert(!"ERROR in FFTJetProducer : invalid jet type." \
78  " This is a bug. Please report."); \
79  } \
80  } while (0);
81 
82 namespace {
83  struct LocalSortByPt {
84  template <class Jet>
85  inline bool operator()(const Jet& l, const Jet& r) const {
86  return l.pt() > r.pt();
87  }
88  };
89 } // namespace
90 
91 using namespace fftjetcms;
92 
94  if (!name.compare("fixed"))
95  return FIXED;
96  else if (!name.compare("maximallyStable"))
97  return MAXIMALLY_STABLE;
98  else if (!name.compare("globallyAdaptive"))
99  return GLOBALLY_ADAPTIVE;
100  else if (!name.compare("locallyAdaptive"))
101  return LOCALLY_ADAPTIVE;
102  else if (!name.compare("fromGenJets"))
103  return FROM_GENJETS;
104  else
105  throw cms::Exception("FFTJetBadConfig") << "Invalid resolution specification \"" << name << "\"" << std::endl;
106 }
107 
108 template <typename T>
110  produces<std::vector<reco::FFTAnyJet<T> > >(tag).setBranchAlias(alias);
111 }
112 
113 //
114 // constructors and destructor
115 //
117  : FFTJetInterface(ps),
118  myConfiguration(ps),
122  init_param(unsigned, maxIterations),
131  init_param(double, fixedScale),
132  init_param(double, minStableScale),
133  init_param(double, maxStableScale),
134  init_param(double, stabilityAlpha),
135  init_param(double, noiseLevel),
136  init_param(unsigned, nClustersRequested),
137  init_param(double, gridScanMaxEta),
140  init_param(double, unlikelyBgWeight),
144  resolution(parse_resolution(ps.getParameter<std::string>("resolution"))),
149 
150  minLevel(0),
151  maxLevel(0),
152  usedLevel(0),
153  unused(0.0),
154  iterationsPerformed(1U),
155  constituents(200) {
156  // Check that the settings make sense
158  throw cms::Exception("FFTJetBadConfig") << "Can't resum constituents if they are not assigned" << std::endl;
159 
160  produces<reco::FFTJetProducerSummary>(outputLabel);
163 
164  if (jetType == CALOJET) {
167  }
168 
169  // Build the set of pattern recognition scales.
170  // This is needed in order to read the clustering tree
171  // from the event record.
173  checkConfig(iniScales, "invalid set of scales");
174  std::sort(iniScales->begin(), iniScales->end(), std::greater<double>());
175 
177  input_recotree_token_f_ = consumes<reco::PattRecoTree<float, reco::PattRecoPeak<float> > >(treeLabel);
178  else
179  input_recotree_token_d_ = consumes<reco::PattRecoTree<double, reco::PattRecoPeak<double> > >(treeLabel);
180  input_genjet_token_ = consumes<std::vector<reco::FFTAnyJet<reco::GenJet> > >(genJetsLabel);
181  input_energyflow_token_ = consumes<reco::DiscretizedEnergyFlow>(treeLabel);
182  input_pusummary_token_ = consumes<reco::FFTJetPileupSummary>(pileupLabel);
183 
184  esLoader_.acquireToken(pileupTableRecord, consumesCollector());
185 
186  // Most of the configuration has to be performed inside
187  // the "beginStream" method. This is because chaining of the
188  // parsers between this base class and the derived classes
189  // can not work from the constructor of the base class.
190 }
191 
193 
194 //
195 // member functions
196 //
197 template <class Real>
201 
202  // Get the input
204  iEvent.getByToken(tok, input);
205 
206  if (!input->isSparse())
207  throw cms::Exception("FFTJetBadConfig") << "The stored clustering tree is not sparse" << std::endl;
208 
210  sparseTree.sortNodes();
211  fftjet::updateSplitMergeTimes(sparseTree, sparseTree.minScale(), sparseTree.maxScale());
212 }
213 
216  const edm::EventSetup& /* iSetup */,
217  const fftjet::Functor1<bool, fftjet::Peak>& peakSelect,
218  std::vector<fftjet::Peak>* preclusters) {
219  typedef reco::FFTAnyJet<reco::GenJet> InputJet;
220  typedef std::vector<InputJet> InputCollection;
221 
223  iEvent.getByToken(input_genjet_token_, input);
224 
225  const unsigned sz = input->size();
226  preclusters->reserve(sz);
227  for (unsigned i = 0; i < sz; ++i) {
228  const RecoFFTJet& jet(jetFromStorable((*input)[i].getFFTSpecific()));
229  fftjet::Peak p(jet.precluster());
230  const double scale(p.scale());
231  p.setEtaPhi(jet.vec().Eta(), jet.vec().Phi());
232  p.setMagnitude(jet.vec().Pt() / scale / scale);
233  p.setStatus(resolution);
234  if (peakSelect(p))
235  preclusters->push_back(p);
236  }
237 }
238 
240  const fftjet::Functor1<bool, fftjet::Peak>& peakSelect,
241  std::vector<fftjet::Peak>* preclusters) {
242  nodes.clear();
243  selectTreeNodes(tree, peakSelect, &nodes);
244 
245  // Fill out the vector of preclusters using the tree node ids
246  const unsigned nNodes = nodes.size();
247  const SparseTree::NodeId* pnodes = nNodes ? &nodes[0] : nullptr;
248  preclusters->reserve(nNodes);
249  for (unsigned i = 0; i < nNodes; ++i)
250  preclusters->push_back(sparseTree.uncheckedNode(pnodes[i]).getCluster());
251 
252  // Remember the node id in the precluster and set
253  // the status word to indicate the resolution scheme used
254  fftjet::Peak* clusters = nNodes ? &(*preclusters)[0] : nullptr;
255  for (unsigned i = 0; i < nNodes; ++i) {
256  clusters[i].setCode(pnodes[i]);
257  clusters[i].setStatus(resolution);
258  }
259 }
260 
262  const fftjet::Functor1<bool, fftjet::Peak>& peakSelect,
263  std::vector<SparseTree::NodeId>* mynodes) {
264  minLevel = maxLevel = usedLevel = 0;
265 
266  // Get the tree nodes which pass the cuts
267  // (according to the selected resolution strategy)
268  switch (resolution) {
269  case FIXED: {
270  usedLevel = tree.getLevel(fixedScale);
271  tree.getPassingNodes(usedLevel, peakSelect, mynodes);
272  } break;
273 
274  case MAXIMALLY_STABLE: {
275  const unsigned minStartingLevel = maxStableScale > 0.0 ? tree.getLevel(maxStableScale) : 0;
276  const unsigned maxStartingLevel = minStableScale > 0.0 ? tree.getLevel(minStableScale) : UINT_MAX;
277 
278  if (tree.stableClusterCount(
279  peakSelect, &minLevel, &maxLevel, stabilityAlpha, minStartingLevel, maxStartingLevel)) {
280  usedLevel = (minLevel + maxLevel) / 2;
281  tree.getPassingNodes(usedLevel, peakSelect, mynodes);
282  }
283  } break;
284 
285  case GLOBALLY_ADAPTIVE: {
286  const bool stable = tree.clusterCountLevels(nClustersRequested, peakSelect, &minLevel, &maxLevel);
287  if (minLevel || maxLevel) {
288  usedLevel = (minLevel + maxLevel) / 2;
289  if (!stable) {
290  const int maxlev = tree.maxStoredLevel();
291  bool levelFound = false;
292  for (int delta = 0; delta <= maxlev && !levelFound; ++delta)
293  for (int ifac = 1; ifac > -2 && !levelFound; ifac -= 2) {
294  const int level = usedLevel + ifac * delta;
295  if (level > 0 && level <= maxlev)
297  usedLevel = level;
298  levelFound = true;
299  }
300  }
301  assert(levelFound);
302  }
303  } else {
304  // Can't find that exact number of preclusters.
305  // Try to get the next best thing.
306  usedLevel = 1;
307  const unsigned occ1 = occupancy[1];
308  if (nClustersRequested >= occ1) {
309  const unsigned maxlev = tree.maxStoredLevel();
310  if (nClustersRequested > occupancy[maxlev])
311  usedLevel = maxlev;
312  else {
313  // It would be nice to use "lower_bound" here,
314  // but the occupancy is not necessarily monotonous.
315  unsigned bestDelta = nClustersRequested > occ1 ? nClustersRequested - occ1 : occ1 - nClustersRequested;
316  for (unsigned level = 2; level <= maxlev; ++level) {
317  const unsigned n = occupancy[level];
318  const unsigned d = nClustersRequested > n ? nClustersRequested - n : n - nClustersRequested;
319  if (d < bestDelta) {
320  bestDelta = d;
321  usedLevel = level;
322  }
323  }
324  }
325  }
326  }
327  tree.getPassingNodes(usedLevel, peakSelect, mynodes);
328  } break;
329 
330  case LOCALLY_ADAPTIVE: {
331  usedLevel = tree.getLevel(fixedScale);
332  tree.getMagS2OptimalNodes(peakSelect, nClustersRequested, usedLevel, mynodes, &thresholds);
333  } break;
334 
335  default:
336  assert(!"ERROR in FFTJetProducer::selectTreeNodes : "
337  "should never get here! This is a bug. Please report.");
338  }
339 }
340 
342  const unsigned nClus = preclusters.size();
343  if (nClus) {
344  fftjet::Peak* clus = &preclusters[0];
345  fftjet::Functor1<double, fftjet::Peak>& scaleCalc(*recoScaleCalcPeak);
346  fftjet::Functor1<double, fftjet::Peak>& ratioCalc(*recoScaleRatioCalcPeak);
347  fftjet::Functor1<double, fftjet::Peak>& factorCalc(*memberFactorCalcPeak);
348 
349  for (unsigned i = 0; i < nClus; ++i) {
350  clus[i].setRecoScale(scaleCalc(clus[i]));
351  clus[i].setRecoScaleRatio(ratioCalc(clus[i]));
352  clus[i].setMembershipFactor(factorCalc(clus[i]));
353  }
354  }
355 }
356 
358  int minBin = energyFlow->getEtaBin(-gridScanMaxEta);
359  if (minBin < 0)
360  minBin = 0;
361  int maxBin = energyFlow->getEtaBin(gridScanMaxEta) + 1;
362  if (maxBin < 0)
363  maxBin = 0;
364 
365  fftjet::DefaultRecombinationAlgFactory<Real, VectorLike, BgData, VBuilder> factory;
366  if (factory[recombinationAlgorithm] == nullptr)
367  throw cms::Exception("FFTJetBadConfig")
368  << "Invalid grid recombination algorithm \"" << recombinationAlgorithm << "\"" << std::endl;
369  gridAlg = std::unique_ptr<GridAlg>(factory[recombinationAlgorithm]->create(jetMembershipFunction.get(),
370  bgMembershipFunction.get(),
373  isCrisp,
374  false,
376  minBin,
377  maxBin));
378 }
379 
380 bool FFTJetProducer::loadEnergyFlow(const edm::Event& iEvent, std::unique_ptr<fftjet::Grid2d<fftjetcms::Real> >& flow) {
382  iEvent.getByToken(input_energyflow_token_, input);
383 
384  // Make sure that the grid is compatible with the stored one
385  bool rebuildGrid = flow.get() == nullptr;
386  if (!rebuildGrid)
387  rebuildGrid =
388  !(flow->nEta() == input->nEtaBins() && flow->nPhi() == input->nPhiBins() && flow->etaMin() == input->etaMin() &&
389  flow->etaMax() == input->etaMax() && flow->phiBin0Edge() == input->phiBin0Edge());
390  if (rebuildGrid) {
391  // We should not get here very often...
392  flow = std::make_unique<fftjet::Grid2d<Real> >(
393  input->nEtaBins(), input->etaMin(), input->etaMax(), input->nPhiBins(), input->phiBin0Edge(), input->title());
394  }
395  flow->blockSet(input->data(), input->nEtaBins(), input->nPhiBins());
396  return rebuildGrid;
397 }
398 
399 bool FFTJetProducer::checkConvergence(const std::vector<RecoFFTJet>& previous, std::vector<RecoFFTJet>& nextSet) {
400  fftjet::Functor2<double, RecoFFTJet, RecoFFTJet>& distanceCalc(*jetDistanceCalc);
401 
402  const unsigned nJets = previous.size();
403  if (nJets != nextSet.size())
404  return false;
405 
406  const RecoFFTJet* prev = &previous[0];
407  RecoFFTJet* next = &nextSet[0];
408 
409  // Calculate convergence distances for all jets
410  bool converged = true;
411  for (unsigned i = 0; i < nJets; ++i) {
412  const double d = distanceCalc(prev[i], next[i]);
413  next[i].setConvergenceDistance(d);
414  if (i < nJetsRequiredToConverge && d > convergenceDistance)
415  converged = false;
416  }
417 
418  return converged;
419 }
420 
422  fftjet::Functor1<double, RecoFFTJet>& scaleCalc(*recoScaleCalcJet);
423  fftjet::Functor1<double, RecoFFTJet>& ratioCalc(*recoScaleRatioCalcJet);
424  fftjet::Functor1<double, RecoFFTJet>& factorCalc(*memberFactorCalcJet);
425 
426  unsigned nJets = recoJets.size();
427  unsigned iterNum = 1U;
428  bool converged = false;
429  for (; iterNum < maxIterations && !converged; ++iterNum) {
430  // Recreate the vector of preclusters using the jets
431  const RecoFFTJet* jets = &recoJets[0];
432  iterPreclusters.clear();
433  iterPreclusters.reserve(nJets);
434  for (unsigned i = 0; i < nJets; ++i) {
435  const RecoFFTJet& jet(jets[i]);
436  fftjet::Peak p(jet.precluster());
437  p.setEtaPhi(jet.vec().Eta(), jet.vec().Phi());
438  p.setRecoScale(scaleCalc(jet));
439  p.setRecoScaleRatio(ratioCalc(jet));
440  p.setMembershipFactor(factorCalc(jet));
441  iterPreclusters.push_back(p);
442  }
443 
444  // Run the algorithm
445  int status = 0;
448  else
450  if (status)
451  throw cms::Exception("FFTJetInterface") << "FFTJet algorithm failed" << std::endl;
452 
453  // As it turns out, it is possible, in very rare cases,
454  // to have iterJets.size() != nJets at this point
455 
456  // Figure out if the iterations have converged
457  converged = checkConvergence(recoJets, iterJets);
458 
459  // Prepare for the next cycle
460  iterJets.swap(recoJets);
461  nJets = recoJets.size();
462  }
463 
464  // Check that we have the correct number of preclusters
465  if (preclusters.size() != nJets) {
466  assert(nJets < preclusters.size());
468  assert(preclusters.size() == nJets);
469  }
470 
471  // Plug in the original precluster coordinates into the result
472  RecoFFTJet* jets = &recoJets[0];
473  for (unsigned i = 0; i < nJets; ++i) {
474  const fftjet::Peak& oldp(preclusters[i]);
475  jets[i].setPeakEtaPhi(oldp.eta(), oldp.phi());
476  }
477 
478  // If we have converged on the last cycle, the result
479  // would be indistinguishable from no convergence.
480  // Because of this, raise the counter by one to indicate
481  // the case when the convergence is not achieved.
482  if (!converged)
483  ++iterNum;
484 
485  return iterNum;
486 }
487 
489  const unsigned nJets = recoJets.size();
490  const unsigned* clusterMask = gridAlg->getClusterMask();
491  const int nEta = gridAlg->getLastNEta();
492  const int nPhi = gridAlg->getLastNPhi();
493  const fftjet::Grid2d<Real>& g(*energyFlow);
494 
495  const unsigned nInputs = eventData.size();
496  const VectorLike* inp = nInputs ? &eventData[0] : nullptr;
497  const unsigned* candIdx = nInputs ? &candidateIndex[0] : nullptr;
498  for (unsigned i = 0; i < nInputs; ++i) {
499  const VectorLike& item(inp[i]);
500  const int iPhi = g.getPhiBin(item.Phi());
501  const int iEta = g.getEtaBin(item.Eta());
502  const unsigned mask = iEta >= 0 && iEta < nEta ? clusterMask[iEta * nPhi + iPhi] : 0;
503  assert(mask <= nJets);
504  constituents[mask].push_back(inputCollection->ptrAt(candIdx[i]));
505  }
506 }
507 
509  const unsigned nJets = recoJets.size();
510  const unsigned* clusterMask = recoAlg->getClusterMask();
511  const unsigned maskLength = recoAlg->getLastNData();
512  assert(maskLength == eventData.size());
513 
514  const unsigned* candIdx = maskLength ? &candidateIndex[0] : nullptr;
515  for (unsigned i = 0; i < maskLength; ++i) {
516  // In FFTJet, the mask value of 0 corresponds to unclustered
517  // energy. We will do the same here. Jet numbers are therefore
518  // shifted by 1 wrt constituents vector, and constituents[1]
519  // corresponds to jet number 0.
520  const unsigned mask = clusterMask[i];
521  assert(mask <= nJets);
522  constituents[mask].push_back(inputCollection->ptrAt(candIdx[i]));
523  }
524 }
525 
526 // The following code more or less coincides with the similar method
527 // implemented in VirtualJetProducer
528 template <typename T>
530  using namespace reco;
531 
532  typedef FFTAnyJet<T> OutputJet;
533  typedef std::vector<OutputJet> OutputCollection;
534 
535  // Area of a single eta-phi cell for jet area calculations.
536  // Set it to 0 in case the module configuration does not allow
537  // us to calculate jet areas reliably.
538  double cellArea = useGriddedAlgorithm && recombinationDataCutoff < 0.0
539  ? energyFlow->etaBinWidth() * energyFlow->phiBinWidth()
540  : 0.0;
541 
542  if (calculatePileup)
543  cellArea = pileupEnergyFlow->etaBinWidth() * pileupEnergyFlow->phiBinWidth();
544 
545  // allocate output jet collection
546  auto jets = std::make_unique<OutputCollection>();
547  const unsigned nJets = recoJets.size();
548  jets->reserve(nJets);
549 
550  bool sorted = true;
551  double previousPt = DBL_MAX;
552  for (unsigned ijet = 0; ijet < nJets; ++ijet) {
553  RecoFFTJet& myjet(recoJets[ijet]);
554 
555  // Check if we should resum jet constituents
556  VectorLike jet4vec(myjet.vec());
557  if (resumConstituents) {
558  VectorLike sum(0.0, 0.0, 0.0, 0.0);
559  const unsigned nCon = constituents[ijet + 1].size();
560  const reco::CandidatePtr* cn = nCon ? &constituents[ijet + 1][0] : nullptr;
561  for (unsigned i = 0; i < nCon; ++i)
562  sum += cn[i]->p4();
563  jet4vec = sum;
564  setJetStatusBit(&myjet, CONSTITUENTS_RESUMMED, true);
565  }
566 
567  // Subtract the pile-up
569  jet4vec = adjustForPileup(jet4vec, pileup[ijet], subtractPileupAs4Vec);
572  else
573  setJetStatusBit(&myjet, PILEUP_SUBTRACTED_PT, true);
574  }
575 
576  // Write the specifics to the jet (simultaneously sets 4-vector,
577  // vertex, constituents). These are overridden functions that will
578  // call the appropriate specific code.
579  T jet;
580  if constexpr (std::is_same_v<T, reco::CaloJet>) {
581  const CaloGeometry& geometry = iSetup.getData(geometry_token_);
582  const HcalTopology& topology = iSetup.getData(topology_token_);
583  writeSpecific(jet, jet4vec, vertexUsed(), constituents[ijet + 1], geometry, topology);
584  } else {
585  writeSpecific(jet, jet4vec, vertexUsed(), constituents[ijet + 1]);
586  }
587 
588  // calcuate the jet area
589  double ncells = myjet.ncells();
590  if (calculatePileup) {
591  ncells = cellCountsVec[ijet];
592  setJetStatusBit(&myjet, PILEUP_CALCULATED, true);
593  }
594  jet.setJetArea(cellArea * ncells);
595 
596  // add jet to the list
597  FFTJet<float> fj(jetToStorable<float>(myjet));
598  fj.setFourVec(jet4vec);
599  if (calculatePileup) {
600  fj.setPileup(pileup[ijet]);
601  fj.setNCells(ncells);
602  }
603  jets->push_back(OutputJet(jet, fj));
604 
605  // Check whether the sequence remains sorted by pt
606  const double pt = jet.pt();
607  if (pt > previousPt)
608  sorted = false;
609  previousPt = pt;
610  }
611 
612  // Sort the collection
613  if (!sorted)
614  std::sort(jets->begin(), jets->end(), LocalSortByPt());
615 
616  // put the collection into the event
618 }
619 
620 void FFTJetProducer::saveResults(edm::Event& ev, const edm::EventSetup& iSetup, const unsigned nPreclustersFound) {
621  // Write recombined jets
622  jet_type_switch(writeJets, ev, iSetup);
623 
624  // Check if we should resum unclustered energy constituents
625  VectorLike unclusE(unclustered);
626  if (resumConstituents) {
627  VectorLike sum(0.0, 0.0, 0.0, 0.0);
628  const unsigned nCon = constituents[0].size();
629  const reco::CandidatePtr* cn = nCon ? &constituents[0][0] : nullptr;
630  for (unsigned i = 0; i < nCon; ++i)
631  sum += cn[i]->p4();
632  unclusE = sum;
633  }
634 
635  // Write the jet reconstruction summary
636  const double minScale = minLevel ? sparseTree.getScale(minLevel) : 0.0;
637  const double maxScale = maxLevel ? sparseTree.getScale(maxLevel) : 0.0;
638  const double scaleUsed = usedLevel ? sparseTree.getScale(usedLevel) : 0.0;
639 
640  ev.put(
641  std::make_unique<reco::FFTJetProducerSummary>(thresholds,
642  occupancy,
643  unclusE,
644  constituents[0],
645  unused,
646  minScale,
647  maxScale,
648  scaleUsed,
649  nPreclustersFound,
652  outputLabel);
653 }
654 
655 // ------------ method called to for each event ------------
657  // Load the clustering tree made by FFTJetPatRecoProducer
659  loadSparseTreeData<float>(iEvent, input_recotree_token_f_);
660  else
661  loadSparseTreeData<double>(iEvent, input_recotree_token_d_);
662 
663  // Do we need to load the candidate collection?
666 
667  // Do we need to have discretized energy flow?
668  if (useGriddedAlgorithm) {
669  if (reuseExistingGrid) {
671  buildGridAlg();
672  } else
674  }
675 
676  // Calculate cluster occupancy as a function of level number
677  sparseTree.occupancyInScaleSpace(*peakSelector, &occupancy);
678 
679  // Select the preclusters using the requested resolution scheme
680  preclusters.clear();
681  if (resolution == FROM_GENJETS)
683  else
685  if (preclusters.size() > maxInitialPreclusters) {
686  std::sort(preclusters.begin(), preclusters.end(), std::greater<fftjet::Peak>());
688  }
689 
690  // Prepare to run the jet recombination procedure
692 
693  // Assign membership functions to preclusters. If this function
694  // is not overriden in a derived class, default algorithm membership
695  // function will be used for every cluster.
697 
698  // Count the preclusters going in
699  unsigned nPreclustersFound = 0U;
700  const unsigned npre = preclusters.size();
701  for (unsigned i = 0; i < npre; ++i)
702  if (preclusters[i].membershipFactor() > 0.0)
703  ++nPreclustersFound;
704 
705  // Run the recombination algorithm once
706  int status = 0;
709  else
711  if (status)
712  throw cms::Exception("FFTJetInterface") << "FFTJet algorithm failed (first iteration)" << std::endl;
713 
714  // If requested, iterate the jet recombination procedure
715  if (maxIterations > 1U && !recoJets.empty()) {
716  // It is possible to have a smaller number of jets than we had
717  // preclusters. Fake preclusters are possible, but for a good
718  // choice of pattern recognition kernel their presence should
719  // be infrequent. However, any fake preclusters will throw the
720  // iterative reconstruction off balance. Deal with the problem now.
721  const unsigned nJets = recoJets.size();
722  if (preclusters.size() != nJets) {
723  assert(nJets < preclusters.size());
725  }
727  } else
729 
730  // Determine jet constituents. FFTJet returns a map
731  // of constituents which is inverse to what we need here.
732  const unsigned nJets = recoJets.size();
733  if (constituents.size() <= nJets)
734  constituents.resize(nJets + 1U);
735  if (assignConstituents) {
736  for (unsigned i = 0; i <= nJets; ++i)
737  constituents[i].clear();
740  else
742  }
743 
744  // Figure out the pile-up
745  if (calculatePileup) {
746  if (loadPileupFromDB)
748  else
750  determinePileup();
751  assert(pileup.size() == recoJets.size());
752  }
753 
754  // Write out the results
755  saveResults(iEvent, iSetup, nPreclustersFound);
756 }
757 
758 std::unique_ptr<fftjet::Functor1<bool, fftjet::Peak> > FFTJetProducer::parse_peakSelector(const edm::ParameterSet& ps) {
759  return fftjet_PeakSelector_parser(ps.getParameter<edm::ParameterSet>("PeakSelectorConfiguration"));
760 }
761 
762 // Parser for the jet membership function
763 std::unique_ptr<fftjet::ScaleSpaceKernel> FFTJetProducer::parse_jetMembershipFunction(const edm::ParameterSet& ps) {
764  return fftjet_MembershipFunction_parser(ps.getParameter<edm::ParameterSet>("jetMembershipFunction"));
765 }
766 
767 // Parser for the background membership function
768 std::unique_ptr<AbsBgFunctor> FFTJetProducer::parse_bgMembershipFunction(const edm::ParameterSet& ps) {
769  return fftjet_BgFunctor_parser(ps.getParameter<edm::ParameterSet>("bgMembershipFunction"));
770 }
771 
772 // Calculator for the recombination scale
773 std::unique_ptr<fftjet::Functor1<double, fftjet::Peak> > FFTJetProducer::parse_recoScaleCalcPeak(
774  const edm::ParameterSet& ps) {
775  return fftjet_PeakFunctor_parser(ps.getParameter<edm::ParameterSet>("recoScaleCalcPeak"));
776 }
777 
778 // Pile-up density calculator
779 std::unique_ptr<fftjetcms::AbsPileupCalculator> FFTJetProducer::parse_pileupDensityCalc(const edm::ParameterSet& ps) {
780  return fftjet_PileupCalculator_parser(ps.getParameter<edm::ParameterSet>("pileupDensityCalc"));
781 }
782 
783 // Calculator for the recombination scale ratio
784 std::unique_ptr<fftjet::Functor1<double, fftjet::Peak> > FFTJetProducer::parse_recoScaleRatioCalcPeak(
785  const edm::ParameterSet& ps) {
786  return fftjet_PeakFunctor_parser(ps.getParameter<edm::ParameterSet>("recoScaleRatioCalcPeak"));
787 }
788 
789 // Calculator for the membership function factor
790 std::unique_ptr<fftjet::Functor1<double, fftjet::Peak> > FFTJetProducer::parse_memberFactorCalcPeak(
791  const edm::ParameterSet& ps) {
792  return fftjet_PeakFunctor_parser(ps.getParameter<edm::ParameterSet>("memberFactorCalcPeak"));
793 }
794 
795 std::unique_ptr<fftjet::Functor1<double, FFTJetProducer::RecoFFTJet> > FFTJetProducer::parse_recoScaleCalcJet(
796  const edm::ParameterSet& ps) {
797  return fftjet_JetFunctor_parser(ps.getParameter<edm::ParameterSet>("recoScaleCalcJet"));
798 }
799 
800 std::unique_ptr<fftjet::Functor1<double, FFTJetProducer::RecoFFTJet> > FFTJetProducer::parse_recoScaleRatioCalcJet(
801  const edm::ParameterSet& ps) {
802  return fftjet_JetFunctor_parser(ps.getParameter<edm::ParameterSet>("recoScaleRatioCalcJet"));
803 }
804 
805 std::unique_ptr<fftjet::Functor1<double, FFTJetProducer::RecoFFTJet> > FFTJetProducer::parse_memberFactorCalcJet(
806  const edm::ParameterSet& ps) {
807  return fftjet_JetFunctor_parser(ps.getParameter<edm::ParameterSet>("memberFactorCalcJet"));
808 }
809 
810 std::unique_ptr<fftjet::Functor2<double, FFTJetProducer::RecoFFTJet, FFTJetProducer::RecoFFTJet> >
812  return fftjet_JetDistance_parser(ps.getParameter<edm::ParameterSet>("jetDistanceCalc"));
813 }
814 
815 void FFTJetProducer::assignMembershipFunctions(std::vector<fftjet::Peak>*) {}
816 
817 // ------------ method called once each job just before starting event loop
820 
821  // Parse the peak selector definition
823  checkConfig(peakSelector, "invalid peak selector");
824 
826  checkConfig(jetMembershipFunction, "invalid jet membership function");
827 
829  checkConfig(bgMembershipFunction, "invalid noise membership function");
830 
831  // Build the energy recombination algorithm
832  if (!useGriddedAlgorithm) {
833  fftjet::DefaultVectorRecombinationAlgFactory<VectorLike, BgData, VBuilder> factory;
834  if (factory[recombinationAlgorithm] == nullptr)
835  throw cms::Exception("FFTJetBadConfig")
836  << "Invalid vector recombination algorithm \"" << recombinationAlgorithm << "\"" << std::endl;
837  recoAlg = std::unique_ptr<RecoAlg>(factory[recombinationAlgorithm]->create(jetMembershipFunction.get(),
838  &VectorLike::Et,
839  &VectorLike::Eta,
841  bgMembershipFunction.get(),
843  isCrisp,
844  false,
846  } else if (!reuseExistingGrid) {
848  checkConfig(energyFlow, "invalid discretization grid");
849  buildGridAlg();
850  }
851 
852  // Create the grid subsequently used for pile-up subtraction
853  if (calculatePileup) {
854  pileupEnergyFlow = fftjet_Grid2d_parser(ps.getParameter<edm::ParameterSet>("PileupGridConfiguration"));
855  checkConfig(pileupEnergyFlow, "invalid pileup density grid");
856 
857  if (!loadPileupFromDB) {
859  checkConfig(pileupDensityCalc, "invalid pile-up density calculator");
860  }
861  }
862 
863  // Parse the calculator of the recombination scale
866  "invalid spec for the "
867  "reconstruction scale calculator from peaks");
868 
869  // Parse the calculator of the recombination scale ratio
872  "invalid spec for the "
873  "reconstruction scale ratio calculator from peaks");
874 
875  // Calculator for the membership function factor
878  "invalid spec for the "
879  "membership function factor calculator from peaks");
880 
881  if (maxIterations > 1) {
882  // We are going to run iteratively. Make required objects.
885  "invalid spec for the "
886  "reconstruction scale calculator from jets");
887 
890  "invalid spec for the "
891  "reconstruction scale ratio calculator from jets");
892 
895  "invalid spec for the "
896  "membership function factor calculator from jets");
897 
900  "invalid spec for the "
901  "jet distance calculator");
902  }
903 }
904 
906  // There are two possible reasons for fake preclusters:
907  // 1. Membership factor was set to 0
908  // 2. Genuine problem with pattern recognition
909  //
910  // Anyway, we need to match jets to preclusters and keep
911  // only those preclusters that have been matched
912  //
913  std::vector<int> matchTable;
914  const unsigned nmatched = matchOneToOne(recoJets, preclusters, JetToPeakDistance(), &matchTable);
915 
916  // Ensure that all jets have been matched.
917  // If not, we must have a bug somewhere.
918  assert(nmatched == recoJets.size());
919 
920  // Collect all matched preclusters
921  iterPreclusters.clear();
922  iterPreclusters.reserve(nmatched);
923  for (unsigned i = 0; i < nmatched; ++i)
924  iterPreclusters.push_back(preclusters[matchTable[i]]);
926 }
927 
928 void FFTJetProducer::setJetStatusBit(RecoFFTJet* jet, const int mask, const bool value) {
929  int status = jet->status();
930  if (value)
931  status |= mask;
932  else
933  status &= ~mask;
934  jet->setStatus(status);
935 }
936 
938  std::unique_ptr<fftjet::Grid2d<fftjetcms::Real> >& density) {
941 
944  const bool phiDependent = calc.isPhiDependent();
945 
946  fftjet::Grid2d<Real>& g(*density);
947  const unsigned nEta = g.nEta();
948  const unsigned nPhi = g.nPhi();
949 
950  for (unsigned ieta = 0; ieta < nEta; ++ieta) {
951  const double eta(g.etaBinCenter(ieta));
952 
953  if (phiDependent) {
954  for (unsigned iphi = 0; iphi < nPhi; ++iphi) {
955  const double phi(g.phiBinCenter(iphi));
956  g.uncheckedSetBin(ieta, iphi, calc(eta, phi, s));
957  }
958  } else {
959  const double pil = calc(eta, 0.0, s);
960  for (unsigned iphi = 0; iphi < nPhi; ++iphi)
961  g.uncheckedSetBin(ieta, iphi, pil);
962  }
963  }
964 }
965 
967  const edm::EventSetup& iSetup,
968  std::unique_ptr<fftjet::Grid2d<fftjetcms::Real> >& density) {
970  std::shared_ptr<npstat::StorableMultivariateFunctor> f = (*h)[pileupTableCategory][pileupTableName];
971 
974 
975  const float rho = summary->pileupRho();
976  const bool phiDependent = f->minDim() == 3U;
977 
978  fftjet::Grid2d<Real>& g(*density);
979  const unsigned nEta = g.nEta();
980  const unsigned nPhi = g.nPhi();
981 
982  double functorArg[3] = {0.0, 0.0, 0.0};
983  if (phiDependent)
984  functorArg[2] = rho;
985  else
986  functorArg[1] = rho;
987 
988  for (unsigned ieta = 0; ieta < nEta; ++ieta) {
989  const double eta(g.etaBinCenter(ieta));
990  functorArg[0] = eta;
991 
992  if (phiDependent) {
993  for (unsigned iphi = 0; iphi < nPhi; ++iphi) {
994  functorArg[1] = g.phiBinCenter(iphi);
995  g.uncheckedSetBin(ieta, iphi, (*f)(functorArg, 3U));
996  }
997  } else {
998  const double pil = (*f)(functorArg, 2U);
999  for (unsigned iphi = 0; iphi < nPhi; ++iphi)
1000  g.uncheckedSetBin(ieta, iphi, pil);
1001  }
1002  }
1003 }
1004 
1006  // This function works with crisp clustering only
1007  if (!isCrisp)
1008  assert(!"Pile-up subtraction for fuzzy clustering "
1009  "is not implemented yet");
1010 
1011  // Clear the pileup vector
1012  const unsigned nJets = recoJets.size();
1013  pileup.resize(nJets);
1014  if (nJets == 0)
1015  return;
1016  const VectorLike zero;
1017  for (unsigned i = 0; i < nJets; ++i)
1018  pileup[i] = zero;
1019 
1020  // Pileup energy flow grid
1021  const fftjet::Grid2d<Real>& g(*pileupEnergyFlow);
1022  const unsigned nEta = g.nEta();
1023  const unsigned nPhi = g.nPhi();
1024  const double cellArea = g.etaBinWidth() * g.phiBinWidth();
1025 
1026  // Various calculators
1027  fftjet::Functor1<double, RecoFFTJet>& scaleCalc(*recoScaleCalcJet);
1028  fftjet::Functor1<double, RecoFFTJet>& ratioCalc(*recoScaleRatioCalcJet);
1029  fftjet::Functor1<double, RecoFFTJet>& factorCalc(*memberFactorCalcJet);
1030 
1031  // Make sure we have enough memory
1032  memFcns2dVec.resize(nJets);
1033  fftjet::AbsKernel2d** memFcns2d = &memFcns2dVec[0];
1034 
1035  doubleBuf.resize(nJets * 4U + nJets * nPhi);
1036  double* recoScales = &doubleBuf[0];
1037  double* recoScaleRatios = recoScales + nJets;
1038  double* memFactors = recoScaleRatios + nJets;
1039  double* dEta = memFactors + nJets;
1040  double* dPhiBuf = dEta + nJets;
1041 
1042  cellCountsVec.resize(nJets);
1043  unsigned* cellCounts = &cellCountsVec[0];
1044 
1045  // Go over jets and collect the necessary info
1046  for (unsigned ijet = 0; ijet < nJets; ++ijet) {
1047  const RecoFFTJet& jet(recoJets[ijet]);
1048  const fftjet::Peak& peak(jet.precluster());
1049 
1050  // Make sure we are using 2-d membership functions.
1051  // Pile-up subtraction scheme for 3-d functions should be different.
1052  fftjet::AbsMembershipFunction* m3d = dynamic_cast<fftjet::AbsMembershipFunction*>(peak.membershipFunction());
1053  if (m3d == nullptr)
1054  m3d = dynamic_cast<fftjet::AbsMembershipFunction*>(jetMembershipFunction.get());
1055  if (m3d) {
1056  assert(!"Pile-up subtraction for 3-d membership functions "
1057  "is not implemented yet");
1058  } else {
1059  fftjet::AbsKernel2d* m2d = dynamic_cast<fftjet::AbsKernel2d*>(peak.membershipFunction());
1060  if (m2d == nullptr)
1061  m2d = dynamic_cast<fftjet::AbsKernel2d*>(jetMembershipFunction.get());
1062  assert(m2d);
1063  memFcns2d[ijet] = m2d;
1064  }
1065  recoScales[ijet] = scaleCalc(jet);
1066  recoScaleRatios[ijet] = ratioCalc(jet);
1067  memFactors[ijet] = factorCalc(jet);
1068  cellCounts[ijet] = 0U;
1069 
1070  const double jetPhi = jet.vec().Phi();
1071  for (unsigned iphi = 0; iphi < nPhi; ++iphi) {
1072  double dphi = g.phiBinCenter(iphi) - jetPhi;
1073  while (dphi > M_PI)
1074  dphi -= (2.0 * M_PI);
1075  while (dphi < -M_PI)
1076  dphi += (2.0 * M_PI);
1077  dPhiBuf[iphi * nJets + ijet] = dphi;
1078  }
1079  }
1080 
1081  // Go over all grid points and integrate
1082  // the pile-up energy density
1083  VBuilder vMaker;
1084  for (unsigned ieta = 0; ieta < nEta; ++ieta) {
1085  const double eta(g.etaBinCenter(ieta));
1086  const Real* databuf = g.data() + ieta * nPhi;
1087 
1088  // Figure out dEta for each jet
1089  for (unsigned ijet = 0; ijet < nJets; ++ijet)
1090  dEta[ijet] = eta - recoJets[ijet].vec().Eta();
1091 
1092  for (unsigned iphi = 0; iphi < nPhi; ++iphi) {
1093  double maxW(0.0);
1094  unsigned maxWJet(nJets);
1095  const double* dPhi = dPhiBuf + iphi * nJets;
1096 
1097  for (unsigned ijet = 0; ijet < nJets; ++ijet) {
1098  if (recoScaleRatios[ijet] > 0.0)
1099  memFcns2d[ijet]->setScaleRatio(recoScaleRatios[ijet]);
1100  const double f = memFactors[ijet] * (*memFcns2d[ijet])(dEta[ijet], dPhi[ijet], recoScales[ijet]);
1101  if (f > maxW) {
1102  maxW = f;
1103  maxWJet = ijet;
1104  }
1105  }
1106 
1107  if (maxWJet < nJets) {
1108  pileup[maxWJet] += vMaker(cellArea * databuf[iphi], eta, g.phiBinCenter(iphi));
1109  cellCounts[maxWJet]++;
1110  }
1111  }
1112  }
1113 }
1114 
1115 //define this as a plug-in
void loadSparseTreeData(const edm::Event &, const edm::EDGetTokenT< reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > > &tok)
std::vector< std::vector< reco::CandidatePtr > > constituents
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:38
edm::ESHandle< DataType > load(const std::string &record, const edm::EventSetup &iSetup) const
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
virtual std::unique_ptr< fftjet::ScaleSpaceKernel > parse_jetMembershipFunction(const edm::ParameterSet &)
std::unique_ptr< fftjet::Functor2< double, RecoFFTJet, RecoFFTJet > > jetDistanceCalc
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
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
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
Preclusters from FFTJet pattern recognition stage.
Definition: PattRecoPeak.h:16
const double gridScanMaxEta
void saveResults(edm::Event &iEvent, const edm::EventSetup &, unsigned nPreclustersFound)
virtual void assignMembershipFunctions(std::vector< fftjet::Peak > *preclusters)
std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > memberFactorCalcPeak
constexpr unsigned int maxBin
void loadInputCollection(const edm::Event &)
std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > recoScaleCalcJet
const std::string recombinationAlgorithm
static const int stable
Definition: TopGenEvent.h:10
const double fixedScale
const bool subtractPileupAs4Vec
std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > fftjet_PeakFunctor_parser(const edm::ParameterSet &ps)
virtual std::unique_ptr< fftjetcms::AbsBgFunctor > parse_bgMembershipFunction(const edm::ParameterSet &)
fftjet::RecombinedJet< VectorLike > jetFromStorable(const reco::FFTJet< Real > &jet)
Definition: jetConverters.h:65
const reco::Particle::Point & vertexUsed() const
void checkConfig(const Ptr &ptr, const char *message)
const double unlikelyBgWeight
void writeSpecific(reco::CaloJet &jet, reco::Particle::LorentzVector const &p4, reco::Particle::Point const &point, std::vector< reco::CandidatePtr > const &constituents, CaloGeometry const &geometry, HcalTopology const &topology)
Definition: JetSpecific.cc:32
void makeProduces(const std::string &alias, const std::string &tag)
std::unique_ptr< AbsBgFunctor > fftjet_BgFunctor_parser(const edm::ParameterSet &ps)
const bool subtractPileup
std::unique_ptr< fftjet::Grid2d< Real > > fftjet_Grid2d_parser(const edm::ParameterSet &ps)
Resolution resolution
std::vector< double > doubleBuf
std::vector< fftjet::AbsKernel2d * > memFcns2dVec
std::vector< RecoFFTJet > recoJets
std::unique_ptr< GridAlg > gridAlg
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geometry_token_
const bool isCrisp
const double recombinationDataCutoff
Summary info for pile-up determined by Gaussian filtering.
static constexpr int nJets
void acquireToken(const std::string &record, edm::ConsumesCollector iC)
std::vector< unsigned > occupancy
std::unique_ptr< T, impl::DeviceDeleter > unique_ptr
assert(be >=bs)
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:50
std::unique_ptr< fftjet::Functor2< double, fftjet::RecombinedJet< VectorLike >, fftjet::RecombinedJet< VectorLike > > > fftjet_JetDistance_parser(const edm::ParameterSet &ps)
T getUntrackedParameter(std::string const &, T const &) const
const bool reuseExistingGrid
const unsigned nClustersRequested
virtual std::unique_ptr< fftjetcms::AbsPileupCalculator > parse_pileupDensityCalc(const edm::ParameterSet &ps)
std::unique_ptr< std::vector< double > > iniScales
int iEvent
Definition: GenABIO.cc:224
std::vector< unsigned > candidateIndex
edm::EDGetTokenT< reco::PattRecoTree< double, reco::PattRecoPeak< double > > > input_recotree_token_d_
std::unique_ptr< fftjetcms::AbsBgFunctor > bgMembershipFunction
std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
Definition: Jet.py:1
virtual std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > parse_recoScaleCalcJet(const edm::ParameterSet &)
std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > memberFactorCalcJet
math::XYZTLorentzVector VectorLike
unsigned usedLevel
fftjet::SparseClusteringTree< fftjet::Peak, long > SparseTree
edm::EDGetTokenT< reco::PattRecoTree< float, reco::PattRecoPeak< float > > > input_recotree_token_f_
std::unique_ptr< fftjet::Functor1< double, fftjet::RecombinedJet< VectorLike > > > fftjet_JetFunctor_parser(const edm::ParameterSet &ps)
virtual std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > parse_recoScaleRatioCalcPeak(const edm::ParameterSet &)
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]
virtual void determinePileupDensityFromDB(const edm::Event &iEvent, const edm::EventSetup &iSetup, std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > &density)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< unsigned > cellCountsVec
std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > recoScaleRatioCalcJet
std::unique_ptr< fftjet::ScaleSpaceKernel > jetMembershipFunction
Definition: value.py:1
std::string pileupTableRecord
std::string pileupTableName
std::vector< double > thresholds
d
Definition: ztail.py:151
unsigned iterateJetReconstruction()
std::unique_ptr< AbsPileupCalculator > fftjet_PileupCalculator_parser(const edm::ParameterSet &ps)
void setNCells(const double nc)
Definition: FFTJet.h:96
const bool calculatePileup
const bool assignConstituents
void setPileup(const math::XYZTLorentzVector &p)
Definition: FFTJet.h:92
#define M_PI
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)
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
virtual std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > parse_recoScaleRatioCalcJet(const edm::ParameterSet &)
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > topology_token_
edm::EDGetTokenT< reco::FFTJetPileupSummary > input_pusummary_token_
const double stabilityAlpha
void beginStream(edm::StreamID) override
const unsigned maxInitialPreclusters
std::string pileupTableCategory
SparseTree sparseTree
fixed size matrix
HLT enums.
std::unique_ptr< fftjetcms::AbsPileupCalculator > pileupDensityCalc
void setFourVec(const math::XYZTLorentzVector &p)
Definition: FFTJet.h:93
std::vector< SparseTree::NodeId > nodes
std::unique_ptr< fftjet::Functor1< bool, fftjet::Peak > > fftjet_PeakSelector_parser(const edm::ParameterSet &ps)
void determineGriddedConstituents()
std::unique_ptr< fftjet::Functor1< bool, fftjet::Peak > > peakSelector
const double convergenceDistance
edm::EDGetTokenT< std::vector< reco::FFTAnyJet< reco::GenJet > > > input_genjet_token_
void writeJets(edm::Event &iEvent, const edm::EventSetup &)
~FFTJetProducer() override
std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > recoScaleRatioCalcPeak
const edm::InputTag pileupLabel
virtual std::unique_ptr< fftjet::Functor1< bool, fftjet::Peak > > parse_peakSelector(const edm::ParameterSet &)
#define jet_type_switch(method, arg1, arg2)
Definition: tree.py:1
std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > pileupEnergyFlow
virtual void genJetPreclusters(const SparseTree &tree, edm::Event &, const edm::EventSetup &, const fftjet::Functor1< bool, fftjet::Peak > &peakSelector, std::vector< fftjet::Peak > *preclusters)
void clear(EGIsoObj &c)
Definition: egamma.h:82
bool loadEnergyFlow(const edm::Event &iEvent, std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > &flow)
virtual bool isPhiDependent() const =0
std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > recoScaleCalcPeak
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
long double T
FFTJetLookupTableSequenceLoader esLoader_
std::unique_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)
std::unique_ptr< RecoAlg > recoAlg
std::vector< fftjetcms::VectorLike > eventData
virtual std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > parse_memberFactorCalcPeak(const edm::ParameterSet &)
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:511
virtual std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > parse_recoScaleCalcPeak(const edm::ParameterSet &)
std::vector< fftjet::Peak > iterPreclusters
virtual void determinePileupDensityFromConfig(const edm::Event &iEvent, std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > &density)
constexpr unsigned int minBin
void removeFakePreclusters()
std::unique_ptr< fftjet::ScaleSpaceKernel > fftjet_MembershipFunction_parser(const edm::ParameterSet &ps)
const unsigned maxIterations
void prepareRecombinationScales()
virtual std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > parse_memberFactorCalcJet(const edm::ParameterSet &)
virtual std::unique_ptr< fftjet::Functor2< double, RecoFFTJet, RecoFFTJet > > parse_jetDistanceCalc(const edm::ParameterSet &)