CMS 3D CMS Logo

FFTJetPatRecoProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: FFTJetProducers
4 // Class: FFTJetPatRecoProducer
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Tue Jun 15 12:45:45 CDT 2010
16 //
17 //
18 
19 #include <fstream>
20 
21 // FFTJet headers
22 #include "fftjet/ProximityClusteringTree.hh"
23 #include "fftjet/ClusteringSequencer.hh"
24 #include "fftjet/ClusteringTreeSparsifier.hh"
25 #include "fftjet/FrequencyKernelConvolver.hh"
26 #include "fftjet/FrequencySequentialConvolver.hh"
27 #include "fftjet/DiscreteGauss1d.hh"
28 #include "fftjet/DiscreteGauss2d.hh"
29 
30 // framework include files
34 
42 
43 // Energy flow object
45 
46 // parameter parser header
48 
49 // functions which manipulate storable trees
51 
52 // functions which manipulate energy discretization grids
54 
55 // useful utilities collected in the second base
57 
58 using namespace fftjetcms;
59 
60 //
61 // class declaration
62 //
64 public:
66  ~FFTJetPatRecoProducer() override;
67 
68 protected:
69  // Useful local typedefs
70  typedef fftjet::ProximityClusteringTree<fftjet::Peak, long> ClusteringTree;
71  typedef fftjet::SparseClusteringTree<fftjet::Peak, long> SparseTree;
72  typedef fftjet::ClusteringSequencer<Real> Sequencer;
73  typedef fftjet::ClusteringTreeSparsifier<fftjet::Peak, long> Sparsifier;
74 
75  // methods
76  void beginJob() override;
77  void produce(edm::Event&, const edm::EventSetup&) override;
78  void endJob() override;
79 
80  void buildKernelConvolver(const edm::ParameterSet&);
81  fftjet::PeakFinder buildPeakFinder(const edm::ParameterSet&);
82 
83  template <class Real>
84  void buildSparseProduct(edm::Event&) const;
85 
86  template <class Real>
87  void buildDenseProduct(edm::Event&) const;
88 
89  // The complete clustering tree
91 
92  // Basically, we need to create FFTJet objects
93  // ClusteringSequencer and ClusteringTreeSparsifier
94  // which will subsequently perform most of the work
95  std::unique_ptr<Sequencer> sequencer;
96  std::unique_ptr<Sparsifier> sparsifier;
97 
98  // The FFT engine(s)
99  std::unique_ptr<MyFFTEngine> engine;
100  std::unique_ptr<MyFFTEngine> anotherEngine;
101 
102  // The pattern recognition kernel(s)
103  std::unique_ptr<fftjet::AbsFrequencyKernel> kernel2d;
104  std::unique_ptr<fftjet::AbsFrequencyKernel1d> etaKernel;
105  std::unique_ptr<fftjet::AbsFrequencyKernel1d> phiKernel;
106 
107  // The kernel convolver
108  std::unique_ptr<fftjet::AbsConvolverBase<Real> > convolver;
109 
110  // The peak selector for the clustering tree
111  std::unique_ptr<fftjet::Functor1<bool, fftjet::Peak> > peakSelector;
112 
113  // Distance calculator for the clustering tree
114  std::unique_ptr<fftjet::AbsDistanceCalculator<fftjet::Peak> > distanceCalc;
115 
116  // The sparse clustering tree
118 
119  // The following parameters will define the behavior
120  // of the algorithm wrt insertion of the complete event
121  // into the clustering tree
123 
124  // Are we going to make clustering trees?
125  const bool makeClusteringTree;
126 
127  // Are we going to verify the data conversion for double precision
128  // storage?
130 
131  // Are we going to store the discretization grid?
133 
134  // Sparsify the clustering tree?
135  const bool sparsify;
136 
137 private:
138  FFTJetPatRecoProducer() = delete;
140  FFTJetPatRecoProducer& operator=(const FFTJetPatRecoProducer&) = delete;
141 
142  // Members needed for storing grids externally
143  std::ofstream externalGridStream;
145  fftjet::Grid2d<float>* extGrid;
146 };
147 
148 //
149 // constructors and destructor
150 //
152  : FFTJetInterface(ps),
153  clusteringTree(nullptr),
154  completeEventDataCutoff(ps.getParameter<double>("completeEventDataCutoff")),
155  makeClusteringTree(ps.getParameter<bool>("makeClusteringTree")),
156  verifyDataConversion(ps.getUntrackedParameter<bool>("verifyDataConversion", false)),
157  storeDiscretizationGrid(ps.getParameter<bool>("storeDiscretizationGrid")),
158  sparsify(ps.getParameter<bool>("sparsify")),
159  extGrid(nullptr) {
160  // register your products
161  if (makeClusteringTree) {
163  produces<reco::PattRecoTree<float, reco::PattRecoPeak<float> > >(outputLabel);
164  else
165  produces<reco::PattRecoTree<double, reco::PattRecoPeak<double> > >(outputLabel);
166  }
168  produces<reco::DiscretizedEnergyFlow>(outputLabel);
169 
170  // Check if we want to write the grids into an external file
171  const std::string externalGridFile(ps.getParameter<std::string>("externalGridFile"));
173  if (storeGridsExternally) {
174  externalGridStream.open(externalGridFile.c_str(), std::ios_base::out | std::ios_base::binary);
175  if (!externalGridStream.is_open())
176  throw cms::Exception("FFTJetBadConfig")
177  << "FFTJetPatRecoProducer failed to open file " << externalGridFile << std::endl;
178  }
179 
181  throw cms::Exception("FFTJetBadConfig")
182  << "FFTJetPatRecoProducer is not configured to produce anything" << std::endl;
183  }
184 
185  // now do whatever other initialization is needed
186 
187  // Build the discretization grid
189  checkConfig(energyFlow, "invalid discretization grid");
190 
191  // Build the FFT engine(s), pattern recognition kernel(s),
192  // and the kernel convolver
194 
195  // Build the peak selector
196  peakSelector = fftjet_PeakSelector_parser(ps.getParameter<edm::ParameterSet>("PeakSelectorConfiguration"));
197  checkConfig(peakSelector, "invalid peak selector");
198 
199  // Build the initial set of pattern recognition scales
200  std::unique_ptr<std::vector<double> > iniScales =
202  checkConfig(iniScales, "invalid set of scales");
203 
204  // Do we want to use the adaptive clustering tree algorithm?
205  const unsigned maxAdaptiveScales = ps.getParameter<unsigned>("maxAdaptiveScales");
206  const double minAdaptiveRatioLog = ps.getParameter<double>("minAdaptiveRatioLog");
207  if (minAdaptiveRatioLog <= 0.0)
208  throw cms::Exception("FFTJetBadConfig") << "bad adaptive ratio logarithm limit" << std::endl;
209 
210  // Make sure that all standard scales are larger than the
211  // complete event scale
212  if (getEventScale() > 0.0) {
213  const double cs = getEventScale();
214  const unsigned nscales = iniScales->size();
215  for (unsigned i = 0; i < nscales; ++i)
216  if (cs >= (*iniScales)[i])
217  throw cms::Exception("FFTJetBadConfig") << "incompatible scale for complete event" << std::endl;
218  }
219 
220  // At this point we are ready to build the clustering sequencer
221  sequencer = std::unique_ptr<Sequencer>(new Sequencer(
223 
224  // Build the clustering tree sparsifier
225  const edm::ParameterSet& SparsifierConfiguration(ps.getParameter<edm::ParameterSet>("SparsifierConfiguration"));
227  checkConfig(sparsifier, "invalid sparsifier parameters");
228 
229  // Build the distance calculator for the clustering tree
230  const edm::ParameterSet& TreeDistanceCalculator(ps.getParameter<edm::ParameterSet>("TreeDistanceCalculator"));
232  checkConfig(distanceCalc, "invalid tree distance calculator");
233 
234  // Build the clustering tree itself
236 }
237 
239  // Check the parameter named "etaDependentScaleFactors". If the vector
240  // of scales is empty we will use 2d kernel, otherwise use 1d kernels
241  const std::vector<double> etaDependentScaleFactors(ps.getParameter<std::vector<double> >("etaDependentScaleFactors"));
242 
243  // Make sure that the number of scale factors provided is correct
244  const bool use2dKernel = etaDependentScaleFactors.empty();
245  if (!use2dKernel)
246  if (etaDependentScaleFactors.size() != energyFlow->nEta())
247  throw cms::Exception("FFTJetBadConfig") << "invalid number of eta-dependent scale factors" << std::endl;
248 
249  // Get the eta and phi scales for the kernel(s)
250  double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
251  const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
252  if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
253  throw cms::Exception("FFTJetBadConfig") << "invalid kernel scale" << std::endl;
254 
255  // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
256  // kernel scale in eta to compensate.
257  kernelEtaScale *= (2.0 * M_PI / (energyFlow->etaMax() - energyFlow->etaMin()));
258 
259  // Are we going to try to fix the efficiency near detector edges?
260  const bool fixEfficiency = ps.getParameter<bool>("fixEfficiency");
261 
262  // Minimum and maximum eta bin for the convolver
263  unsigned convolverMinBin = 0, convolverMaxBin = 0;
264  if (fixEfficiency || !use2dKernel) {
265  convolverMinBin = ps.getParameter<unsigned>("convolverMinBin");
266  convolverMaxBin = ps.getParameter<unsigned>("convolverMaxBin");
267  }
268 
269  if (use2dKernel) {
270  // Build the FFT engine
271  engine = std::unique_ptr<MyFFTEngine>(new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
272 
273  // 2d kernel
274  kernel2d = std::unique_ptr<fftjet::AbsFrequencyKernel>(
275  new fftjet::DiscreteGauss2d(kernelEtaScale, kernelPhiScale, energyFlow->nEta(), energyFlow->nPhi()));
276 
277  // 2d convolver
278  convolver = std::unique_ptr<fftjet::AbsConvolverBase<Real> >(new fftjet::FrequencyKernelConvolver<Real, Complex>(
280  } else {
281  // Need separate FFT engines for eta and phi
282  engine = std::unique_ptr<MyFFTEngine>(new MyFFTEngine(1, energyFlow->nEta()));
283  anotherEngine = std::unique_ptr<MyFFTEngine>(new MyFFTEngine(1, energyFlow->nPhi()));
284 
285  // 1d kernels
286  etaKernel =
287  std::unique_ptr<fftjet::AbsFrequencyKernel1d>(new fftjet::DiscreteGauss1d(kernelEtaScale, energyFlow->nEta()));
288 
289  phiKernel =
290  std::unique_ptr<fftjet::AbsFrequencyKernel1d>(new fftjet::DiscreteGauss1d(kernelPhiScale, energyFlow->nPhi()));
291 
292  // Sequential convolver
293  convolver = std::unique_ptr<fftjet::AbsConvolverBase<Real> >(
294  new fftjet::FrequencySequentialConvolver<Real, Complex>(engine.get(),
295  anotherEngine.get(),
296  etaKernel.get(),
297  phiKernel.get(),
301  fixEfficiency));
302  }
303 }
304 
306  const double peakFinderMaxEta = ps.getParameter<double>("peakFinderMaxEta");
307  if (peakFinderMaxEta <= 0.0)
308  throw cms::Exception("FFTJetBadConfig") << "invalid peak finder eta cut" << std::endl;
309  const double maxMagnitude = ps.getParameter<double>("peakFinderMaxMagnitude");
310  int minBin = energyFlow->getEtaBin(-peakFinderMaxEta);
311  if (minBin < 0)
312  minBin = 0;
313  int maxBin = energyFlow->getEtaBin(peakFinderMaxEta) + 1;
314  if (maxBin < 0)
315  maxBin = 0;
316  return fftjet::PeakFinder(maxMagnitude, true, minBin, maxBin);
317 }
318 
320  // do anything here that needs to be done at desctruction time
321  // (e.g. close files, deallocate resources etc.)
322  delete clusteringTree;
323  delete extGrid;
324 }
325 
326 //
327 // member functions
328 //
329 template <class Real>
332 
333  auto tree = std::make_unique<StoredTree>();
334 
335  sparsePeakTreeToStorable(sparseTree, sequencer->maxAdaptiveScales(), tree.get());
336 
337  // Check that we can restore the tree
340  const std::vector<double>& scalesUsed(sequencer->getInitialScales());
342  if (sparseTree != check)
343  throw cms::Exception("FFTJetInterface") << "Data conversion failed for sparse clustering tree" << std::endl;
344  }
345 
346  ev.put(std::move(tree), outputLabel);
347 }
348 
349 template <class Real>
352 
353  auto tree = std::make_unique<StoredTree>();
354 
355  densePeakTreeToStorable(*clusteringTree, sequencer->maxAdaptiveScales(), tree.get());
356 
357  // Check that we can restore the tree
360  const std::vector<double>& scalesUsed(sequencer->getInitialScales());
362  if (*clusteringTree != check)
363  throw cms::Exception("FFTJetInterface") << "Data conversion failed for dense clustering tree" << std::endl;
364  }
365 
366  ev.put(std::move(tree), outputLabel);
367 }
368 
369 // ------------ method called to produce the data ------------
373 
374  if (makeClusteringTree) {
376  if (getEventScale() > 0.0)
378 
379  if (sparsify) {
380  sparsifier->sparsify(*clusteringTree, &sparseTree);
381 
382  // Do not call the "sortNodes" method of the sparse tree here.
383  // Currently, the nodes are sorted by daughter number.
384  // This is the way we want it in storage because the stored
385  // tree does not include daughter ordering info explicitly.
386 
388  buildSparseProduct<float>(iEvent);
389  else
390  buildSparseProduct<double>(iEvent);
391  } else {
393  buildDenseProduct<float>(iEvent);
394  else
395  buildDenseProduct<double>(iEvent);
396  }
397  }
398 
400  const fftjet::Grid2d<Real>& g(*energyFlow);
401 
402  auto flow = std::make_unique<reco::DiscretizedEnergyFlow>(
403  g.data(), g.title(), g.etaMin(), g.etaMax(), g.phiBin0Edge(), g.nEta(), g.nPhi());
404 
405  if (verifyDataConversion) {
406  fftjet::Grid2d<Real> check(
407  flow->nEtaBins(), flow->etaMin(), flow->etaMax(), flow->nPhiBins(), flow->phiBin0Edge(), flow->title());
408  check.blockSet(flow->data(), flow->nEtaBins(), flow->nPhiBins());
409  assert(g == check);
410  }
411 
413  }
414 
415  if (storeGridsExternally) {
416  if (extGrid)
418  else
420  if (!extGrid->write(externalGridStream)) {
421  throw cms::Exception("FFTJetPatRecoProducer::produce")
422  << "Failed to write grid data into an external file" << std::endl;
423  }
424  }
425 }
426 
427 // ------------ method called once each job just before starting event loop
429 
430 // ------------ method called once each job just after ending the event loop
433  externalGridStream.close();
434 }
435 
436 //define this as a plug-in
fftjetcms::densePeakTreeToStorable
void densePeakTreeToStorable(const fftjet::AbsClusteringTree< fftjet::Peak, long > &in, bool writeOutScaleInfo, reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > *out)
Definition: clusteringTreeConverters.h:210
cms::cuda::allocator::minBin
constexpr unsigned int minBin
Definition: getCachingDeviceAllocator.h:18
CaloJetCollection.h
GenJetCollection.h
bk::beginJob
void beginJob()
Definition: Breakpoints.cc:14
electrons_cff.bool
bool
Definition: electrons_cff.py:372
mps_fire.i
i
Definition: mps_fire.py:355
FFTJetPatRecoProducer::buildDenseProduct
void buildDenseProduct(edm::Event &) const
Definition: FFTJetPatRecoProducer.cc:350
FFTJetPatRecoProducer::phiKernel
std::unique_ptr< fftjet::AbsFrequencyKernel1d > phiKernel
Definition: FFTJetPatRecoProducer.cc:105
funct::false
false
Definition: Factorize.h:34
fftjetcms
Definition: AbsPileupCalculator.h:15
fftjetpatrecoproducer_cfi.sparsify
sparsify
Definition: fftjetpatrecoproducer_cfi.py:25
FFTJetPatRecoProducer
Definition: FFTJetPatRecoProducer.cc:63
fftjetcms::FFTJetInterface::checkConfig
void checkConfig(const Ptr &ptr, const char *message)
Definition: FFTJetInterface.h:60
fftjetcms::FFTJetInterface::loadInputCollection
void loadInputCollection(const edm::Event &)
Definition: FFTJetInterface.cc:40
FFTJetPatRecoProducer::engine
std::unique_ptr< MyFFTEngine > engine
Definition: FFTJetPatRecoProducer.cc:99
fftjetcms::FFTJetInterface::outputLabel
const std::string outputLabel
Definition: FFTJetInterface.h:76
fwrapper::cs
unique_ptr< ClusterSequence > cs
Definition: fastjetfortran_madfks.cc:45
FFTJetPatRecoProducer::completeEventDataCutoff
const double completeEventDataCutoff
Definition: FFTJetPatRecoProducer.cc:122
FFTJetPatRecoProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: FFTJetPatRecoProducer.cc:370
fftjetcms::fftjet_ScaleSet_parser
std::unique_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)
Definition: FFTJetParameterParser.cc:348
FFTJetPatRecoProducer::sparseTree
SparseTree sparseTree
Definition: FFTJetPatRecoProducer.cc:117
tree
Definition: tree.py:1
FFTJetPatRecoProducer::verifyDataConversion
const bool verifyDataConversion
Definition: FFTJetPatRecoProducer.cc:129
fftjetpatrecoproducer_cfi.completeEventDataCutoff
completeEventDataCutoff
Definition: fftjetpatrecoproducer_cfi.py:57
FFTJetPatRecoProducer::distanceCalc
std::unique_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > distanceCalc
Definition: FFTJetPatRecoProducer.cc:114
PFJetCollection.h
fftjeteflowsmoother_cfi.convolverMaxBin
convolverMaxBin
Definition: fftjeteflowsmoother_cfi.py:21
cms::cuda::assert
assert(be >=bs)
BasicJetCollection.h
clusteringTreeConverters.h
fftjetcms::FFTJetInterface::getEventScale
double getEventScale() const
Definition: FFTJetInterface.cc:38
fftjetpatrecoproducer_cfi.minAdaptiveRatioLog
minAdaptiveRatioLog
Definition: fftjetpatrecoproducer_cfi.py:84
FFTJetPatRecoProducer::FFTJetPatRecoProducer
FFTJetPatRecoProducer()=delete
FFTJetPatRecoProducer::buildPeakFinder
fftjet::PeakFinder buildPeakFinder(const edm::ParameterSet &)
Definition: FFTJetPatRecoProducer.cc:305
FFTJetPatRecoProducer::buildSparseProduct
void buildSparseProduct(edm::Event &) const
Definition: FFTJetPatRecoProducer.cc:330
FFTJetInterface.h
fftjetpatrecoproducer_cfi.verifyDataConversion
verifyDataConversion
Definition: fftjetpatrecoproducer_cfi.py:22
fftjetcms::densePeakTreeFromStorable
void densePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::AbsClusteringTree< fftjet::Peak, long > *out)
Definition: clusteringTreeConverters.h:271
FFTJetPatRecoProducer::~FFTJetPatRecoProducer
~FFTJetPatRecoProducer() override
Definition: FFTJetPatRecoProducer.cc:319
MakerMacros.h
RPCNoise_example.check
check
Definition: RPCNoise_example.py:71
FFTJetPatRecoProducer::kernel2d
std::unique_ptr< fftjet::AbsFrequencyKernel > kernel2d
Definition: FFTJetPatRecoProducer.cc:103
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
FFTJetPatRecoProducer::Sequencer
fftjet::ClusteringSequencer< Real > Sequencer
Definition: FFTJetPatRecoProducer.cc:72
fftjetpatrecoproducer_cfi.SparsifierConfiguration
SparsifierConfiguration
Definition: fftjetpatrecoproducer_cfi.py:115
fftjeteflowsmoother_cfi.kernelEtaScale
kernelEtaScale
Definition: fftjeteflowsmoother_cfi.py:11
FFTJetPatRecoProducer::etaKernel
std::unique_ptr< fftjet::AbsFrequencyKernel1d > etaKernel
Definition: FFTJetPatRecoProducer.cc:104
FFTJetPatRecoProducer::SparseTree
fftjet::SparseClusteringTree< fftjet::Peak, long > SparseTree
Definition: FFTJetPatRecoProducer.cc:71
FFTJetPatRecoProducer::convolver
std::unique_ptr< fftjet::AbsConvolverBase< Real > > convolver
Definition: FFTJetPatRecoProducer.cc:108
fftjetcms::FFTJetInterface
Definition: FFTJetInterface.h:52
FFTJetPatRecoProducer::ClusteringTree
fftjet::ProximityClusteringTree< fftjet::Peak, long > ClusteringTree
Definition: FFTJetPatRecoProducer.cc:70
fftjetcms::MyFFTEngine
fftjet::FFTWDoubleEngine MyFFTEngine
Definition: fftjetTypedefs.h:23
fftjetcms::sparsePeakTreeToStorable
void sparsePeakTreeToStorable(const fftjet::SparseClusteringTree< fftjet::Peak, long > &in, bool writeOutScaleInfo, reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > *out)
Definition: clusteringTreeConverters.h:75
fftjetpatrecoproducer_cfi.peakFinderMaxEta
peakFinderMaxEta
Definition: fftjetpatrecoproducer_cfi.py:36
fftjetcms::sparsePeakTreeFromStorable
void sparsePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::SparseClusteringTree< fftjet::Peak, long > *out)
Definition: clusteringTreeConverters.h:133
TrackJetCollection.h
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
FFTJetPatRecoProducer::clusteringTree
ClusteringTree * clusteringTree
Definition: FFTJetPatRecoProducer.cc:90
fftjetcms::fftjet_ClusteringTreeSparsifier_parser
std::unique_ptr< fftjet::ClusteringTreeSparsifier< fftjet::Peak, long > > fftjet_ClusteringTreeSparsifier_parser(const edm::ParameterSet &ps)
Definition: FFTJetParameterParser.cc:390
FFTJetPatRecoProducer::buildKernelConvolver
void buildKernelConvolver(const edm::ParameterSet &)
Definition: FFTJetPatRecoProducer.cc:238
reco::PattRecoTree
Class for storing FFTJet sparse clustering trees.
Definition: PattRecoTree.h:20
edm::ParameterSet
Definition: ParameterSet.h:36
FFTJetPatRecoProducer::peakSelector
std::unique_ptr< fftjet::Functor1< bool, fftjet::Peak > > peakSelector
Definition: FFTJetPatRecoProducer.cc:111
FFTJetPatRecoProducer::extGrid
fftjet::Grid2d< float > * extGrid
Definition: FFTJetPatRecoProducer.cc:145
FFTJetPatRecoProducer::storeGridsExternally
bool storeGridsExternally
Definition: FFTJetPatRecoProducer.cc:144
FFTJetPatRecoProducer::anotherEngine
std::unique_ptr< MyFFTEngine > anotherEngine
Definition: FFTJetPatRecoProducer.cc:100
cms::cuda::allocator::maxBin
constexpr unsigned int maxBin
Definition: getCachingDeviceAllocator.h:20
iEvent
int iEvent
Definition: GenABIO.cc:224
fftjetcms::FFTJetInterface::energyFlow
std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
Definition: FFTJetInterface.h:100
FFTJetPatRecoProducer::beginJob
void beginJob() override
Definition: FFTJetPatRecoProducer.cc:428
FFTJetPatRecoProducer::makeClusteringTree
const bool makeClusteringTree
Definition: FFTJetPatRecoProducer.cc:125
FFTJetPatRecoProducer::sequencer
std::unique_ptr< Sequencer > sequencer
Definition: FFTJetPatRecoProducer.cc:95
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:50
fftjetcms::FFTJetInterface::storeInSinglePrecision
bool storeInSinglePrecision() const
Definition: FFTJetInterface.cc:15
FFTJetPatRecoProducer::storeDiscretizationGrid
const bool storeDiscretizationGrid
Definition: FFTJetPatRecoProducer.cc:132
fftjetcms::fftjet_Grid2d_parser
std::unique_ptr< fftjet::Grid2d< Real > > fftjet_Grid2d_parser(const edm::ParameterSet &ps)
Definition: FFTJetParameterParser.cc:123
edm::EventSetup
Definition: EventSetup.h:57
fftjeteflowsmoother_cfi.convolverMinBin
convolverMinBin
Definition: fftjeteflowsmoother_cfi.py:20
fftjeteflowsmoother_cfi.fixEfficiency
fixEfficiency
Definition: fftjeteflowsmoother_cfi.py:15
fftjeteflowsmoother_cfi.kernelPhiScale
kernelPhiScale
Definition: fftjeteflowsmoother_cfi.py:12
fftjetpatrecoproducer_cfi.storeDiscretizationGrid
storeDiscretizationGrid
Definition: fftjetpatrecoproducer_cfi.py:28
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
FFTJetPatRecoProducer::endJob
void endJob() override
Definition: FFTJetPatRecoProducer.cc:431
fftjetcommon_cfi.flow
flow
Definition: fftjetcommon_cfi.py:191
eostools.move
def move(src, dest)
Definition: eostools.py:511
fftjetcms::copy_Grid2d_data
void copy_Grid2d_data(fftjet::Grid2d< F2 > *to, const fftjet::Grid2d< F1 > &from)
Definition: gridConverters.h:39
Frameworkfwd.h
fftjetpatrecoproducer_cfi.maxAdaptiveScales
maxAdaptiveScales
Definition: fftjetpatrecoproducer_cfi.py:79
FFTJetPatRecoProducer::sparsifier
std::unique_ptr< Sparsifier > sparsifier
Definition: FFTJetPatRecoProducer.cc:96
FFTJetPatRecoProducer::Sparsifier
fftjet::ClusteringTreeSparsifier< fftjet::Peak, long > Sparsifier
Definition: FFTJetPatRecoProducer.cc:73
fftjetdijetfilter_cfi.TreeDistanceCalculator
TreeDistanceCalculator
Definition: fftjetdijetfilter_cfi.py:25
ev
bool ev
Definition: Hydjet2Hadronizer.cc:95
Exception
Definition: hltDiff.cc:246
FFTJetParameterParser.h
FFTJetPatRecoProducer::externalGridStream
std::ofstream externalGridStream
Definition: FFTJetPatRecoProducer.cc:143
fftjetcms::convert_Grid2d_to_float
fftjet::Grid2d< float > * convert_Grid2d_to_float(const fftjet::Grid2d< Numeric > &grid)
Definition: gridConverters.h:63
fftjetpatrecoproducer_cfi.externalGridFile
externalGridFile
Definition: fftjetpatrecoproducer_cfi.py:32
fftjetcms::FFTJetInterface::discretizeEnergyFlow
void discretizeEnergyFlow()
Definition: FFTJetInterface.cc:79
MillePedeFileConverter_cfg.out
out
Definition: MillePedeFileConverter_cfg.py:31
DiscretizedEnergyFlow.h
cms::Exception
Definition: Exception.h:70
Candidate.h
fftjeteflowsmoother_cfi.etaDependentScaleFactors
etaDependentScaleFactors
Definition: fftjeteflowsmoother_cfi.py:42
View.h
ParameterSet.h
fftjetpatrecoproducer_cfi.makeClusteringTree
makeClusteringTree
Definition: fftjetpatrecoproducer_cfi.py:17
edm::Event
Definition: Event.h:73
fftjetcms::fftjet_DistanceCalculator_parser
std::unique_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > fftjet_DistanceCalculator_parser(const edm::ParameterSet &ps)
Definition: FFTJetParameterParser.cc:404
gridConverters.h
fftjetcms::fftjet_PeakSelector_parser
std::unique_ptr< fftjet::Functor1< bool, fftjet::Peak > > fftjet_PeakSelector_parser(const edm::ParameterSet &ps)
Definition: FFTJetParameterParser.cc:164
g
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
FFTJetPatRecoProducer::sparsify
const bool sparsify
Definition: FFTJetPatRecoProducer.cc:135