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
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 {
65 public:
68 
69 protected:
70  // Useful local typedefs
71  typedef fftjet::ProximityClusteringTree<fftjet::Peak,long> ClusteringTree;
72  typedef fftjet::SparseClusteringTree<fftjet::Peak,long> SparseTree;
73  typedef fftjet::ClusteringSequencer<Real> Sequencer;
74  typedef fftjet::ClusteringTreeSparsifier<fftjet::Peak,long> Sparsifier;
75 
76  // methods
77  void beginJob() override ;
78  void produce(edm::Event&, const edm::EventSetup&) override;
79  void endJob() override ;
80 
81  void buildKernelConvolver(const edm::ParameterSet&);
82  fftjet::PeakFinder buildPeakFinder(const edm::ParameterSet&);
83 
84  template<class Real>
85  void buildSparseProduct(edm::Event&) const;
86 
87  template<class Real>
88  void buildDenseProduct(edm::Event&) const;
89 
90  // The complete clustering tree
92 
93  // Basically, we need to create FFTJet objects
94  // ClusteringSequencer and ClusteringTreeSparsifier
95  // which will subsequently perform most of the work
96  std::auto_ptr<Sequencer> sequencer;
97  std::auto_ptr<Sparsifier> sparsifier;
98 
99  // The FFT engine(s)
100  std::auto_ptr<MyFFTEngine> engine;
101  std::auto_ptr<MyFFTEngine> anotherEngine;
102 
103  // The pattern recognition kernel(s)
104  std::auto_ptr<fftjet::AbsFrequencyKernel> kernel2d;
105  std::auto_ptr<fftjet::AbsFrequencyKernel1d> etaKernel;
106  std::auto_ptr<fftjet::AbsFrequencyKernel1d> phiKernel;
107 
108  // The kernel convolver
109  std::auto_ptr<fftjet::AbsConvolverBase<Real> > convolver;
110 
111  // The peak selector for the clustering tree
112  std::auto_ptr<fftjet::Functor1<bool,fftjet::Peak> > peakSelector;
113 
114  // Distance calculator for the clustering tree
115  std::auto_ptr<fftjet::AbsDistanceCalculator<fftjet::Peak> > distanceCalc;
116 
117  // The sparse clustering tree
119 
120  // The following parameters will define the behavior
121  // of the algorithm wrt insertion of the complete event
122  // into the clustering tree
124 
125  // Are we going to make clustering trees?
126  const bool makeClusteringTree;
127 
128  // Are we going to verify the data conversion for double precision
129  // storage?
131 
132  // Are we going to store the discretization grid?
134 
135  // Sparsify the clustering tree?
136  const bool sparsify;
137 
138 private:
141  FFTJetPatRecoProducer& operator=(const FFTJetPatRecoProducer&);
142 
143  // Members needed for storing grids externally
144  std::ofstream externalGridStream;
146  fftjet::Grid2d<float>* extGrid;
147 };
148 
149 //
150 // constructors and destructor
151 //
153  : FFTJetInterface(ps),
154  clusteringTree(0),
155  completeEventDataCutoff(ps.getParameter<double>("completeEventDataCutoff")),
156  makeClusteringTree(ps.getParameter<bool>("makeClusteringTree")),
157  verifyDataConversion(ps.getUntrackedParameter<bool>("verifyDataConversion",false)),
158  storeDiscretizationGrid(ps.getParameter<bool>("storeDiscretizationGrid")),
159  sparsify(ps.getParameter<bool>("sparsify")),
160  extGrid(0)
161 {
162  // register your products
163  if (makeClusteringTree)
164  {
166  produces<reco::PattRecoTree<float,reco::PattRecoPeak<float> > >(outputLabel);
167  else
168  produces<reco::PattRecoTree<double,reco::PattRecoPeak<double> > >(outputLabel);
169  }
171  produces<reco::DiscretizedEnergyFlow>(outputLabel);
172 
173  // Check if we want to write the grids into an external file
174  const std::string externalGridFile(ps.getParameter<std::string>("externalGridFile"));
175  storeGridsExternally = externalGridFile.size() > 0;
177  {
178  externalGridStream.open(externalGridFile.c_str(), std::ios_base::out |
179  std::ios_base::binary);
180  if (!externalGridStream.is_open())
181  throw cms::Exception("FFTJetBadConfig")
182  << "FFTJetPatRecoProducer failed to open file "
183  << externalGridFile << std::endl;
184  }
185 
187  {
188  throw cms::Exception("FFTJetBadConfig")
189  << "FFTJetPatRecoProducer is not configured to produce anything"
190  << std::endl;
191  }
192 
193  // now do whatever other initialization is needed
194 
195  // Build the discretization grid
197  ps.getParameter<edm::ParameterSet>("GridConfiguration"));
198  checkConfig(energyFlow, "invalid discretization grid");
199 
200  // Build the FFT engine(s), pattern recognition kernel(s),
201  // and the kernel convolver
203 
204  // Build the peak selector
206  ps.getParameter<edm::ParameterSet>("PeakSelectorConfiguration"));
207  checkConfig(peakSelector, "invalid peak selector");
208 
209  // Build the initial set of pattern recognition scales
210  std::auto_ptr<std::vector<double> > iniScales = fftjet_ScaleSet_parser(
211  ps.getParameter<edm::ParameterSet>("InitialScales"));
212  checkConfig(iniScales, "invalid set of scales");
213 
214  // Do we want to use the adaptive clustering tree algorithm?
215  const unsigned maxAdaptiveScales =
216  ps.getParameter<unsigned>("maxAdaptiveScales");
217  const double minAdaptiveRatioLog =
218  ps.getParameter<double>("minAdaptiveRatioLog");
219  if (minAdaptiveRatioLog <= 0.0)
220  throw cms::Exception("FFTJetBadConfig")
221  << "bad adaptive ratio logarithm limit" << std::endl;
222 
223  // Make sure that all standard scales are larger than the
224  // complete event scale
225  if (getEventScale() > 0.0)
226  {
227  const double cs = getEventScale();
228  const unsigned nscales = iniScales->size();
229  for (unsigned i=0; i<nscales; ++i)
230  if (cs >= (*iniScales)[i])
231  throw cms::Exception("FFTJetBadConfig")
232  << "incompatible scale for complete event" << std::endl;
233  }
234 
235  // At this point we are ready to build the clustering sequencer
236  sequencer = std::auto_ptr<Sequencer>(new Sequencer(
237  convolver.get(), peakSelector.get(), buildPeakFinder(ps),
238  *iniScales, maxAdaptiveScales, minAdaptiveRatioLog));
239 
240  // Build the clustering tree sparsifier
241  const edm::ParameterSet& SparsifierConfiguration(
242  ps.getParameter<edm::ParameterSet>("SparsifierConfiguration"));
244  SparsifierConfiguration);
245  checkConfig(sparsifier, "invalid sparsifier parameters");
246 
247  // Build the distance calculator for the clustering tree
248  const edm::ParameterSet& TreeDistanceCalculator(
249  ps.getParameter<edm::ParameterSet>("TreeDistanceCalculator"));
250  distanceCalc = fftjet_DistanceCalculator_parser(TreeDistanceCalculator);
251  checkConfig(distanceCalc, "invalid tree distance calculator");
252 
253  // Build the clustering tree itself
255 }
256 
257 
259 {
260  // Check the parameter named "etaDependentScaleFactors". If the vector
261  // of scales is empty we will use 2d kernel, otherwise use 1d kernels
262  const std::vector<double> etaDependentScaleFactors(
263  ps.getParameter<std::vector<double> >("etaDependentScaleFactors"));
264 
265  // Make sure that the number of scale factors provided is correct
266  const bool use2dKernel = etaDependentScaleFactors.empty();
267  if (!use2dKernel)
268  if (etaDependentScaleFactors.size() != energyFlow->nEta())
269  throw cms::Exception("FFTJetBadConfig")
270  << "invalid number of eta-dependent scale factors"
271  << std::endl;
272 
273  // Get the eta and phi scales for the kernel(s)
274  double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
275  const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
276  if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
277  throw cms::Exception("FFTJetBadConfig")
278  << "invalid kernel scale" << std::endl;
279 
280  // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
281  // kernel scale in eta to compensate.
282  kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
283 
284  // Are we going to try to fix the efficiency near detector edges?
285  const bool fixEfficiency = ps.getParameter<bool>("fixEfficiency");
286 
287  // Minimum and maximum eta bin for the convolver
288  unsigned convolverMinBin = 0, convolverMaxBin = 0;
289  if (fixEfficiency || !use2dKernel)
290  {
291  convolverMinBin = ps.getParameter<unsigned>("convolverMinBin");
292  convolverMaxBin = ps.getParameter<unsigned>("convolverMaxBin");
293  }
294 
295  if (use2dKernel)
296  {
297  // Build the FFT engine
298  engine = std::auto_ptr<MyFFTEngine>(
299  new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
300 
301  // 2d kernel
302  kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
303  new fftjet::DiscreteGauss2d(
304  kernelEtaScale, kernelPhiScale,
305  energyFlow->nEta(), energyFlow->nPhi()));
306 
307  // 2d convolver
308  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
309  new fftjet::FrequencyKernelConvolver<Real,Complex>(
310  engine.get(), kernel2d.get(),
311  convolverMinBin, convolverMaxBin));
312  }
313  else
314  {
315  // Need separate FFT engines for eta and phi
316  engine = std::auto_ptr<MyFFTEngine>(
317  new MyFFTEngine(1, energyFlow->nEta()));
318  anotherEngine = std::auto_ptr<MyFFTEngine>(
319  new MyFFTEngine(1, energyFlow->nPhi()));
320 
321  // 1d kernels
322  etaKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
323  new fftjet::DiscreteGauss1d(kernelEtaScale, energyFlow->nEta()));
324 
325  phiKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
326  new fftjet::DiscreteGauss1d(kernelPhiScale, energyFlow->nPhi()));
327 
328  // Sequential convolver
329  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
330  new fftjet::FrequencySequentialConvolver<Real,Complex>(
331  engine.get(), anotherEngine.get(),
332  etaKernel.get(), phiKernel.get(),
333  etaDependentScaleFactors, convolverMinBin,
334  convolverMaxBin, fixEfficiency));
335  }
336 }
337 
338 
340 {
341  const double peakFinderMaxEta = ps.getParameter<double>("peakFinderMaxEta");
342  if (peakFinderMaxEta <= 0.0)
343  throw cms::Exception("FFTJetBadConfig")
344  << "invalid peak finder eta cut" << std::endl;
345  const double maxMagnitude = ps.getParameter<double>("peakFinderMaxMagnitude");
346  int minBin = energyFlow->getEtaBin(-peakFinderMaxEta);
347  if (minBin < 0)
348  minBin = 0;
349  int maxBin = energyFlow->getEtaBin(peakFinderMaxEta) + 1;
350  if (maxBin < 0)
351  maxBin = 0;
352  return fftjet::PeakFinder(maxMagnitude, true, minBin, maxBin);
353 }
354 
355 
357 {
358  // do anything here that needs to be done at desctruction time
359  // (e.g. close files, deallocate resources etc.)
360  delete clusteringTree;
361  delete extGrid;
362 }
363 
364 
365 //
366 // member functions
367 //
368 template<class Real>
370 {
372 
373  std::auto_ptr<StoredTree> tree(new StoredTree());
374 
376  sequencer->maxAdaptiveScales(),
377  tree.get());
378 
379  // Check that we can restore the tree
381  {
383  const std::vector<double>& scalesUsed(sequencer->getInitialScales());
384  sparsePeakTreeFromStorable(*tree, &scalesUsed, getEventScale(), &check);
385  if (sparseTree != check)
386  throw cms::Exception("FFTJetInterface")
387  << "Data conversion failed for sparse clustering tree"
388  << std::endl;
389  }
390 
391  ev.put(tree, outputLabel);
392 }
393 
394 
395 template<class Real>
397 {
399 
400  std::auto_ptr<StoredTree> tree(new StoredTree());
401 
403  sequencer->maxAdaptiveScales(),
404  tree.get());
405 
406  // Check that we can restore the tree
408  {
410  const std::vector<double>& scalesUsed(sequencer->getInitialScales());
411  densePeakTreeFromStorable(*tree, &scalesUsed, getEventScale(), &check);
412  if (*clusteringTree != check)
413  throw cms::Exception("FFTJetInterface")
414  << "Data conversion failed for dense clustering tree"
415  << std::endl;
416  }
417 
418  ev.put(tree, outputLabel);
419 }
420 
421 
422 // ------------ method called to produce the data ------------
424  edm::Event& iEvent, const edm::EventSetup& iSetup)
425 {
426  loadInputCollection(iEvent);
428 
429  if (makeClusteringTree)
430  {
432  if (getEventScale() > 0.0)
433  sequencer->insertCompleteEvent(getEventScale(), *energyFlow,
435 
436  if (sparsify)
437  {
438  sparsifier->sparsify(*clusteringTree, &sparseTree);
439 
440  // Do not call the "sortNodes" method of the sparse tree here.
441  // Currently, the nodes are sorted by daughter number.
442  // This is the way we want it in storage because the stored
443  // tree does not include daughter ordering info explicitly.
444 
446  buildSparseProduct<float>(iEvent);
447  else
448  buildSparseProduct<double>(iEvent);
449  }
450  else
451  {
453  buildDenseProduct<float>(iEvent);
454  else
455  buildDenseProduct<double>(iEvent);
456  }
457  }
458 
460  {
461  const fftjet::Grid2d<Real>& g(*energyFlow);
462 
463  std::auto_ptr<reco::DiscretizedEnergyFlow> flow(
465  g.data(), g.title(), g.etaMin(), g.etaMax(),
466  g.phiBin0Edge(), g.nEta(), g.nPhi()));
467 
469  {
470  fftjet::Grid2d<Real> check(
471  flow->nEtaBins(), flow->etaMin(), flow->etaMax(),
472  flow->nPhiBins(), flow->phiBin0Edge(), flow->title());
473  check.blockSet(flow->data(), flow->nEtaBins(), flow->nPhiBins());
474  assert(g == check);
475  }
476 
477  iEvent.put(flow, outputLabel);
478  }
479 
481  {
482  if (extGrid)
484  else
486  if (!extGrid->write(externalGridStream))
487  {
488  throw cms::Exception("FFTJetPatRecoProducer::produce")
489  << "Failed to write grid data into an external file"
490  << std::endl;
491  }
492  }
493 }
494 
495 
496 // ------------ method called once each job just before starting event loop
498 {
499 }
500 
501 
502 // ------------ method called once each job just after ending the event loop
504 {
506  externalGridStream.close();
507 }
508 
509 
510 //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 &)
assert(m_qm.get())
void buildSparseProduct(edm::Event &) const
void buildKernelConvolver(const edm::ParameterSet &)
void checkConfig(const Ptr &ptr, const char *message)
bool ev
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:230
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:121
void sparsePeakTreeToStorable(const fftjet::SparseClusteringTree< fftjet::Peak, long > &in, bool writeOutScaleInfo, reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > *out)
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
#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)
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)