CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 
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
35 
45 
46 // Energy flow object
48 
49 // parameter parser header
51 
52 // functions which manipulate storable trees
54 
55 // functions which manipulate energy discretization grids
57 
58 // useful utilities collected in the second base
60 
61 using namespace fftjetcms;
62 
63 //
64 // class declaration
65 //
67 {
68 public:
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::auto_ptr<Sequencer> sequencer;
100  std::auto_ptr<Sparsifier> sparsifier;
101 
102  // The FFT engine(s)
103  std::auto_ptr<MyFFTEngine> engine;
104  std::auto_ptr<MyFFTEngine> anotherEngine;
105 
106  // The pattern recognition kernel(s)
107  std::auto_ptr<fftjet::AbsFrequencyKernel> kernel2d;
108  std::auto_ptr<fftjet::AbsFrequencyKernel1d> etaKernel;
109  std::auto_ptr<fftjet::AbsFrequencyKernel1d> phiKernel;
110 
111  // The kernel convolver
112  std::auto_ptr<fftjet::AbsConvolverBase<Real> > convolver;
113 
114  // The peak selector for the clustering tree
115  std::auto_ptr<fftjet::Functor1<bool,fftjet::Peak> > peakSelector;
116 
117  // Distance calculator for the clustering tree
118  std::auto_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:
144  FFTJetPatRecoProducer& operator=(const FFTJetPatRecoProducer&);
145 
146  // Members needed for storing grids externally
147  std::ofstream externalGridStream;
149  fftjet::Grid2d<float>* extGrid;
150 };
151 
152 //
153 // constructors and destructor
154 //
156  : FFTJetInterface(ps),
157  clusteringTree(0),
158  completeEventDataCutoff(ps.getParameter<double>("completeEventDataCutoff")),
159  makeClusteringTree(ps.getParameter<bool>("makeClusteringTree")),
160  verifyDataConversion(ps.getUntrackedParameter<bool>("verifyDataConversion",false)),
161  storeDiscretizationGrid(ps.getParameter<bool>("storeDiscretizationGrid")),
162  sparsify(ps.getParameter<bool>("sparsify")),
163  extGrid(0)
164 {
165  // register your products
166  if (makeClusteringTree)
167  {
169  produces<reco::PattRecoTree<float,reco::PattRecoPeak<float> > >(outputLabel);
170  else
171  produces<reco::PattRecoTree<double,reco::PattRecoPeak<double> > >(outputLabel);
172  }
174  produces<reco::DiscretizedEnergyFlow>(outputLabel);
175 
176  // Check if we want to write the grids into an external file
177  const std::string externalGridFile(ps.getParameter<std::string>("externalGridFile"));
178  storeGridsExternally = externalGridFile.size() > 0;
180  {
181  externalGridStream.open(externalGridFile.c_str(), std::ios_base::out |
182  std::ios_base::binary);
183  if (!externalGridStream.is_open())
184  throw cms::Exception("FFTJetBadConfig")
185  << "FFTJetPatRecoProducer failed to open file "
186  << externalGridFile << std::endl;
187  }
188 
190  {
191  throw cms::Exception("FFTJetBadConfig")
192  << "FFTJetPatRecoProducer is not configured to produce anything"
193  << std::endl;
194  }
195 
196  // now do whatever other initialization is needed
197 
198  // Build the discretization grid
200  ps.getParameter<edm::ParameterSet>("GridConfiguration"));
201  checkConfig(energyFlow, "invalid discretization grid");
202 
203  // Build the FFT engine(s), pattern recognition kernel(s),
204  // and the kernel convolver
206 
207  // Build the peak selector
209  ps.getParameter<edm::ParameterSet>("PeakSelectorConfiguration"));
210  checkConfig(peakSelector, "invalid peak selector");
211 
212  // Build the initial set of pattern recognition scales
213  std::auto_ptr<std::vector<double> > iniScales = fftjet_ScaleSet_parser(
214  ps.getParameter<edm::ParameterSet>("InitialScales"));
215  checkConfig(iniScales, "invalid set of scales");
216 
217  // Do we want to use the adaptive clustering tree algorithm?
218  const unsigned maxAdaptiveScales =
219  ps.getParameter<unsigned>("maxAdaptiveScales");
220  const double minAdaptiveRatioLog =
221  ps.getParameter<double>("minAdaptiveRatioLog");
222  if (minAdaptiveRatioLog <= 0.0)
223  throw cms::Exception("FFTJetBadConfig")
224  << "bad adaptive ratio logarithm limit" << std::endl;
225 
226  // Make sure that all standard scales are larger than the
227  // complete event scale
228  if (getEventScale() > 0.0)
229  {
230  const double cs = getEventScale();
231  const unsigned nscales = iniScales->size();
232  for (unsigned i=0; i<nscales; ++i)
233  if (cs >= (*iniScales)[i])
234  throw cms::Exception("FFTJetBadConfig")
235  << "incompatible scale for complete event" << std::endl;
236  }
237 
238  // At this point we are ready to build the clustering sequencer
239  sequencer = std::auto_ptr<Sequencer>(new Sequencer(
240  convolver.get(), peakSelector.get(), buildPeakFinder(ps),
241  *iniScales, maxAdaptiveScales, minAdaptiveRatioLog));
242 
243  // Build the clustering tree sparsifier
244  const edm::ParameterSet& SparsifierConfiguration(
245  ps.getParameter<edm::ParameterSet>("SparsifierConfiguration"));
247  SparsifierConfiguration);
248  checkConfig(sparsifier, "invalid sparsifier parameters");
249 
250  // Build the distance calculator for the clustering tree
251  const edm::ParameterSet& TreeDistanceCalculator(
252  ps.getParameter<edm::ParameterSet>("TreeDistanceCalculator"));
253  distanceCalc = fftjet_DistanceCalculator_parser(TreeDistanceCalculator);
254  checkConfig(distanceCalc, "invalid tree distance calculator");
255 
256  // Build the clustering tree itself
258 }
259 
260 
262 {
263  // Check the parameter named "etaDependentScaleFactors". If the vector
264  // of scales is empty we will use 2d kernel, otherwise use 1d kernels
265  const std::vector<double> etaDependentScaleFactors(
266  ps.getParameter<std::vector<double> >("etaDependentScaleFactors"));
267 
268  // Make sure that the number of scale factors provided is correct
269  const bool use2dKernel = etaDependentScaleFactors.empty();
270  if (!use2dKernel)
271  if (etaDependentScaleFactors.size() != energyFlow->nEta())
272  throw cms::Exception("FFTJetBadConfig")
273  << "invalid number of eta-dependent scale factors"
274  << std::endl;
275 
276  // Get the eta and phi scales for the kernel(s)
277  double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
278  const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
279  if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
280  throw cms::Exception("FFTJetBadConfig")
281  << "invalid kernel scale" << std::endl;
282 
283  // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
284  // kernel scale in eta to compensate.
285  kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
286 
287  // Are we going to try to fix the efficiency near detector edges?
288  const bool fixEfficiency = ps.getParameter<bool>("fixEfficiency");
289 
290  // Minimum and maximum eta bin for the convolver
291  unsigned convolverMinBin = 0, convolverMaxBin = 0;
292  if (fixEfficiency || !use2dKernel)
293  {
294  convolverMinBin = ps.getParameter<unsigned>("convolverMinBin");
295  convolverMaxBin = ps.getParameter<unsigned>("convolverMaxBin");
296  }
297 
298  if (use2dKernel)
299  {
300  // Build the FFT engine
301  engine = std::auto_ptr<MyFFTEngine>(
302  new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
303 
304  // 2d kernel
305  kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
306  new fftjet::DiscreteGauss2d(
307  kernelEtaScale, kernelPhiScale,
308  energyFlow->nEta(), energyFlow->nPhi()));
309 
310  // 2d convolver
311  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
312  new fftjet::FrequencyKernelConvolver<Real,Complex>(
313  engine.get(), kernel2d.get(),
314  convolverMinBin, convolverMaxBin));
315  }
316  else
317  {
318  // Need separate FFT engines for eta and phi
319  engine = std::auto_ptr<MyFFTEngine>(
320  new MyFFTEngine(1, energyFlow->nEta()));
321  anotherEngine = std::auto_ptr<MyFFTEngine>(
322  new MyFFTEngine(1, energyFlow->nPhi()));
323 
324  // 1d kernels
325  etaKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
326  new fftjet::DiscreteGauss1d(kernelEtaScale, energyFlow->nEta()));
327 
328  phiKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
329  new fftjet::DiscreteGauss1d(kernelPhiScale, energyFlow->nPhi()));
330 
331  // Sequential convolver
332  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
333  new fftjet::FrequencySequentialConvolver<Real,Complex>(
334  engine.get(), anotherEngine.get(),
335  etaKernel.get(), phiKernel.get(),
336  etaDependentScaleFactors, convolverMinBin,
337  convolverMaxBin, fixEfficiency));
338  }
339 }
340 
341 
343 {
344  const double peakFinderMaxEta = ps.getParameter<double>("peakFinderMaxEta");
345  if (peakFinderMaxEta <= 0.0)
346  throw cms::Exception("FFTJetBadConfig")
347  << "invalid peak finder eta cut" << std::endl;
348  const double maxMagnitude = ps.getParameter<double>("peakFinderMaxMagnitude");
349  int minBin = energyFlow->getEtaBin(-peakFinderMaxEta);
350  if (minBin < 0)
351  minBin = 0;
352  int maxBin = energyFlow->getEtaBin(peakFinderMaxEta) + 1;
353  if (maxBin < 0)
354  maxBin = 0;
355  return fftjet::PeakFinder(maxMagnitude, true, minBin, maxBin);
356 }
357 
358 
360 {
361  // do anything here that needs to be done at desctruction time
362  // (e.g. close files, deallocate resources etc.)
363  delete clusteringTree;
364  delete extGrid;
365 }
366 
367 
368 //
369 // member functions
370 //
371 template<class Real>
373 {
375 
376  std::auto_ptr<StoredTree> tree(new StoredTree());
377 
379  sequencer->maxAdaptiveScales(),
380  tree.get());
381 
382  // Check that we can restore the tree
384  {
386  const std::vector<double>& scalesUsed(sequencer->getInitialScales());
387  sparsePeakTreeFromStorable(*tree, &scalesUsed, getEventScale(), &check);
388  if (sparseTree != check)
389  throw cms::Exception("FFTJetInterface")
390  << "Data conversion failed for sparse clustering tree"
391  << std::endl;
392  }
393 
394  ev.put(tree, outputLabel);
395 }
396 
397 
398 template<class Real>
400 {
402 
403  std::auto_ptr<StoredTree> tree(new StoredTree());
404 
406  sequencer->maxAdaptiveScales(),
407  tree.get());
408 
409  // Check that we can restore the tree
411  {
413  const std::vector<double>& scalesUsed(sequencer->getInitialScales());
414  densePeakTreeFromStorable(*tree, &scalesUsed, getEventScale(), &check);
415  if (*clusteringTree != check)
416  throw cms::Exception("FFTJetInterface")
417  << "Data conversion failed for dense clustering tree"
418  << std::endl;
419  }
420 
421  ev.put(tree, outputLabel);
422 }
423 
424 
425 // ------------ method called to produce the data ------------
427  edm::Event& iEvent, const edm::EventSetup& iSetup)
428 {
429  loadInputCollection(iEvent);
431 
432  if (makeClusteringTree)
433  {
435  if (getEventScale() > 0.0)
436  sequencer->insertCompleteEvent(getEventScale(), *energyFlow,
438 
439  if (sparsify)
440  {
441  sparsifier->sparsify(*clusteringTree, &sparseTree);
442 
443  // Do not call the "sortNodes" method of the sparse tree here.
444  // Currently, the nodes are sorted by daughter number.
445  // This is the way we want it in storage because the stored
446  // tree does not include daughter ordering info explicitly.
447 
449  buildSparseProduct<float>(iEvent);
450  else
451  buildSparseProduct<double>(iEvent);
452  }
453  else
454  {
456  buildDenseProduct<float>(iEvent);
457  else
458  buildDenseProduct<double>(iEvent);
459  }
460  }
461 
463  {
464  const fftjet::Grid2d<Real>& g(*energyFlow);
465 
466  std::auto_ptr<reco::DiscretizedEnergyFlow> flow(
468  g.data(), g.title(), g.etaMin(), g.etaMax(),
469  g.phiBin0Edge(), g.nEta(), g.nPhi()));
470 
472  {
473  fftjet::Grid2d<Real> check(
474  flow->nEtaBins(), flow->etaMin(), flow->etaMax(),
475  flow->nPhiBins(), flow->phiBin0Edge(), flow->title());
476  check.blockSet(flow->data(), flow->nEtaBins(), flow->nPhiBins());
477  assert(g == check);
478  }
479 
480  iEvent.put(flow, outputLabel);
481  }
482 
484  {
485  if (extGrid)
487  else
489  if (!extGrid->write(externalGridStream))
490  {
491  throw cms::Exception("FFTJetPatRecoProducer::produce")
492  << "Failed to write grid data into an external file"
493  << std::endl;
494  }
495  }
496 }
497 
498 
499 // ------------ method called once each job just before starting event loop
501 {
502 }
503 
504 
505 // ------------ method called once each job just after ending the event loop
507 {
509  externalGridStream.close();
510 }
511 
512 
513 //define this as a plug-in
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
fftjet::Grid2d< float > * convert_Grid2d_to_float(const fftjet::Grid2d< Numeric > &grid)
std::auto_ptr< fftjet::Functor1< bool, fftjet::Peak > > peakSelector
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
auto_ptr< ClusterSequence > cs
Class for storing FFTJet sparse clustering trees.
Definition: PattRecoTree.h:18
std::auto_ptr< fftjet::AbsFrequencyKernel1d > phiKernel
std::auto_ptr< fftjet::Grid2d< Real > > fftjet_Grid2d_parser(const edm::ParameterSet &ps)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void loadInputCollection(const edm::Event &)
void buildSparseProduct(edm::Event &) const
void buildKernelConvolver(const edm::ParameterSet &)
void checkConfig(const Ptr &ptr, const char *message)
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
void beginJob()
Definition: Breakpoints.cc:15
std::auto_ptr< MyFFTEngine > anotherEngine
std::auto_ptr< fftjet::Functor1< bool, fftjet::Peak > > fftjet_PeakSelector_parser(const edm::ParameterSet &ps)
int iEvent
Definition: GenABIO.cc:243
std::auto_ptr< fftjet::AbsFrequencyKernel > kernel2d
void produce(edm::Event &, const edm::EventSetup &) override
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
void sparsePeakTreeToStorable(const fftjet::SparseClusteringTree< fftjet::Peak, long > &in, bool writeOutScaleInfo, reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > *out)
bool check(const DataFrame &df, bool capcheck, bool dvercheck)
ClusteringTree * clusteringTree
fftjet::ClusteringSequencer< Real > Sequencer
std::auto_ptr< MyFFTEngine > engine
std::auto_ptr< Sparsifier > sparsifier
std::auto_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > distanceCalc
bool storeInSinglePrecision() const
void densePeakTreeToStorable(const fftjet::AbsClusteringTree< fftjet::Peak, long > &in, bool writeOutScaleInfo, reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > *out)
tuple out
Definition: dbtoconf.py:99
void sparsePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::SparseClusteringTree< fftjet::Peak, long > *out)
#define M_PI
Definition: BFit3D.cc:3
void densePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::AbsClusteringTree< fftjet::Peak, long > *out)
fftjet::ClusteringTreeSparsifier< fftjet::Peak, long > Sparsifier
std::auto_ptr< fftjet::AbsFrequencyKernel1d > etaKernel
void copy_Grid2d_data(fftjet::Grid2d< F2 > *to, const fftjet::Grid2d< F1 > &from)
std::auto_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > fftjet_DistanceCalculator_parser(const edm::ParameterSet &ps)
fftjet::PeakFinder buildPeakFinder(const edm::ParameterSet &)
std::auto_ptr< Sequencer > sequencer
std::auto_ptr< fftjet::AbsConvolverBase< Real > > convolver
std::auto_ptr< fftjet::ClusteringTreeSparsifier< fftjet::Peak, long > > fftjet_ClusteringTreeSparsifier_parser(const edm::ParameterSet &ps)
fftjet::FFTWDoubleEngine MyFFTEngine
volatile std::atomic< bool > shutdown_flag false
fftjet::ProximityClusteringTree< fftjet::Peak, long > ClusteringTree
fftjet::SparseClusteringTree< fftjet::Peak, long > SparseTree
fftjet::Grid2d< float > * extGrid
const std::string outputLabel
std::auto_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)