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