CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 #include <memory>
21 
22 // FFTJet headers
23 #include "fftjet/ProximityClusteringTree.hh"
24 #include "fftjet/ClusteringSequencer.hh"
25 #include "fftjet/ClusteringTreeSparsifier.hh"
26 #include "fftjet/FrequencyKernelConvolver.hh"
27 #include "fftjet/FrequencySequentialConvolver.hh"
28 #include "fftjet/DiscreteGauss1d.hh"
29 #include "fftjet/DiscreteGauss2d.hh"
30 
31 // framework include files
35 
43 
44 // Energy flow object
46 
47 // parameter parser header
49 
50 // functions which manipulate storable trees
52 
53 // functions which manipulate energy discretization grids
55 
56 // useful utilities collected in the second base
58 
59 using namespace fftjetcms;
60 
61 //
62 // class declaration
63 //
65 public:
67  FFTJetPatRecoProducer() = delete;
69  FFTJetPatRecoProducer& operator=(const FFTJetPatRecoProducer&) = delete;
70  ~FFTJetPatRecoProducer() override;
71 
72 protected:
73  // Useful local typedefs
74  typedef fftjet::ProximityClusteringTree<fftjet::Peak, long> ClusteringTree;
75  typedef fftjet::SparseClusteringTree<fftjet::Peak, long> SparseTree;
76  typedef fftjet::ClusteringSequencer<Real> Sequencer;
77  typedef fftjet::ClusteringTreeSparsifier<fftjet::Peak, long> Sparsifier;
78 
79  // methods
80  void beginJob() override;
81  void produce(edm::Event&, const edm::EventSetup&) override;
82  void endJob() override;
83 
84  void buildKernelConvolver(const edm::ParameterSet&);
85  fftjet::PeakFinder buildPeakFinder(const edm::ParameterSet&);
86 
87  template <class Real>
88  void buildSparseProduct(edm::Event&) const;
89 
90  template <class Real>
91  void buildDenseProduct(edm::Event&) const;
92 
93  // The complete clustering tree
95 
96  // Basically, we need to create FFTJet objects
97  // ClusteringSequencer and ClusteringTreeSparsifier
98  // which will subsequently perform most of the work
99  std::unique_ptr<Sequencer> sequencer;
100  std::unique_ptr<Sparsifier> sparsifier;
101 
102  // The FFT engine(s)
103  std::unique_ptr<MyFFTEngine> engine;
104  std::unique_ptr<MyFFTEngine> anotherEngine;
105 
106  // The pattern recognition kernel(s)
107  std::unique_ptr<fftjet::AbsFrequencyKernel> kernel2d;
108  std::unique_ptr<fftjet::AbsFrequencyKernel1d> etaKernel;
109  std::unique_ptr<fftjet::AbsFrequencyKernel1d> phiKernel;
110 
111  // The kernel convolver
112  std::unique_ptr<fftjet::AbsConvolverBase<Real> > convolver;
113 
114  // The peak selector for the clustering tree
115  std::unique_ptr<fftjet::Functor1<bool, fftjet::Peak> > peakSelector;
116 
117  // Distance calculator for the clustering tree
118  std::unique_ptr<fftjet::AbsDistanceCalculator<fftjet::Peak> > distanceCalc;
119 
120  // The sparse clustering tree
122 
123  // The following parameters will define the behavior
124  // of the algorithm wrt insertion of the complete event
125  // into the clustering tree
127 
128  // Are we going to make clustering trees?
129  const bool makeClusteringTree;
130 
131  // Are we going to verify the data conversion for double precision
132  // storage?
134 
135  // Are we going to store the discretization grid?
137 
138  // Sparsify the clustering tree?
139  const bool sparsify;
140 
141 private:
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"));
172  storeGridsExternally = !externalGridFile.empty();
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::make_unique<Sequencer>(
222  convolver.get(), peakSelector.get(), buildPeakFinder(ps), *iniScales, maxAdaptiveScales, minAdaptiveRatioLog);
223 
224  // Build the clustering tree sparsifier
225  const edm::ParameterSet& SparsifierConfiguration(ps.getParameter<edm::ParameterSet>("SparsifierConfiguration"));
226  sparsifier = fftjet_ClusteringTreeSparsifier_parser(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"));
231  distanceCalc = fftjet_DistanceCalculator_parser(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::make_unique<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>(
279  engine.get(), kernel2d.get(), convolverMinBin, convolverMaxBin));
280  } else {
281  // Need separate FFT engines for eta and phi
282  engine = std::make_unique<MyFFTEngine>(1, energyFlow->nEta());
283  anotherEngine = std::make_unique<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(),
298  etaDependentScaleFactors,
299  convolverMinBin,
300  convolverMaxBin,
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());
341  sparsePeakTreeFromStorable(*tree, &scalesUsed, getEventScale(), &check);
342  if (sparseTree != check)
343  throw cms::Exception("FFTJetInterface") << "Data conversion failed for sparse clustering tree" << std::endl;
344  }
345 
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 
367 }
368 
369 // ------------ method called to produce the data ------------
371  loadInputCollection(iEvent);
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 
412  iEvent.put(std::move(flow), outputLabel);
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
fftjet::Grid2d< float > * convert_Grid2d_to_float(const fftjet::Grid2d< Numeric > &grid)
std::unique_ptr< Sparsifier > sparsifier
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
std::unique_ptr< fftjet::AbsFrequencyKernel1d > etaKernel
Class for storing FFTJet sparse clustering trees.
Definition: PattRecoTree.h:20
constexpr unsigned int maxBin
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void loadInputCollection(const edm::Event &)
unique_ptr< ClusterSequence > cs
void buildSparseProduct(edm::Event &) const
void buildKernelConvolver(const edm::ParameterSet &)
void checkConfig(const Ptr &ptr, const char *message)
bool ev
std::unique_ptr< fftjet::Grid2d< Real > > fftjet_Grid2d_parser(const edm::ParameterSet &ps)
fftjet::SparseClusteringTree< fftjet::Peak, long > SparseTree
std::unique_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > distanceCalc
assert(be >=bs)
fftjet::ClusteringTreeSparsifier< fftjet::Peak, long > Sparsifier
void buildDenseProduct(edm::Event &) const
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
std::unique_ptr< fftjet::AbsFrequencyKernel1d > phiKernel
FFTJetPatRecoProducer()=delete
void beginJob()
Definition: Breakpoints.cc:14
std::unique_ptr< Sequencer > sequencer
std::unique_ptr< fftjet::Functor1< bool, fftjet::Peak > > peakSelector
int iEvent
Definition: GenABIO.cc:224
void produce(edm::Event &, const edm::EventSetup &) override
std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
void sparsePeakTreeToStorable(const fftjet::SparseClusteringTree< fftjet::Peak, long > &in, bool writeOutScaleInfo, reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > *out)
std::unique_ptr< MyFFTEngine > anotherEngine
ClusteringTree * clusteringTree
def move
Definition: eostools.py:511
fftjet::ClusteringSequencer< Real > Sequencer
std::unique_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > fftjet_DistanceCalculator_parser(const edm::ParameterSet &ps)
bool storeInSinglePrecision() const
#define M_PI
void densePeakTreeToStorable(const fftjet::AbsClusteringTree< fftjet::Peak, long > &in, bool writeOutScaleInfo, reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > *out)
void sparsePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::SparseClusteringTree< fftjet::Peak, long > *out)
void densePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::AbsClusteringTree< fftjet::Peak, long > *out)
std::unique_ptr< fftjet::AbsConvolverBase< Real > > convolver
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void copy_Grid2d_data(fftjet::Grid2d< F2 > *to, const fftjet::Grid2d< F1 > &from)
std::unique_ptr< fftjet::AbsFrequencyKernel > kernel2d
fftjet::PeakFinder buildPeakFinder(const edm::ParameterSet &)
std::unique_ptr< fftjet::Functor1< bool, fftjet::Peak > > fftjet_PeakSelector_parser(const edm::ParameterSet &ps)
std::unique_ptr< fftjet::ClusteringTreeSparsifier< fftjet::Peak, long > > fftjet_ClusteringTreeSparsifier_parser(const edm::ParameterSet &ps)
fftjet::ProximityClusteringTree< fftjet::Peak, long > ClusteringTree
fftjet::Grid2d< float > * extGrid
std::unique_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)
const std::string outputLabel
std::unique_ptr< MyFFTEngine > engine
constexpr unsigned int minBin