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.1 2010/12/06 17:33:19 igv Exp $
17 //
18 //
19 
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 // 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() ;
78  void produce(edm::Event&, const edm::EventSetup&);
79  void endJob() ;
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:
142 };
143 
144 //
145 // constructors and destructor
146 //
148  : FFTJetInterface(ps),
149  clusteringTree(0),
150  completeEventDataCutoff(ps.getParameter<double>("completeEventDataCutoff")),
151  makeClusteringTree(ps.getParameter<bool>("makeClusteringTree")),
152  verifyDataConversion(ps.getUntrackedParameter<bool>("verifyDataConversion",false)),
153  storeDiscretizationGrid(ps.getParameter<bool>("storeDiscretizationGrid")),
154  sparsify(ps.getParameter<bool>("sparsify"))
155 {
156  // register your products
157  if (makeClusteringTree)
158  {
160  produces<reco::PattRecoTree<float,reco::PattRecoPeak<float> > >(outputLabel);
161  else
162  produces<reco::PattRecoTree<double,reco::PattRecoPeak<double> > >(outputLabel);
163  }
165  produces<DiscretizedEnergyFlow>(outputLabel);
166 
168  {
169  throw cms::Exception("FFTJetBadConfig")
170  << "FFTJetPatRecoProducer is not configured to produce anything"
171  << std::endl;
172  }
173 
174  // now do whatever other initialization is needed
175 
176  // Build the discretization grid
178  ps.getParameter<edm::ParameterSet>("GridConfiguration"));
179  checkConfig(energyFlow, "invalid discretization grid");
180 
181  // Build the FFT engine(s), pattern recognition kernel(s),
182  // and the kernel convolver
184 
185  // Build the peak selector
187  ps.getParameter<edm::ParameterSet>("PeakSelectorConfiguration"));
188  checkConfig(peakSelector, "invalid peak selector");
189 
190  // Build the initial set of pattern recognition scales
191  std::auto_ptr<std::vector<double> > iniScales = fftjet_ScaleSet_parser(
192  ps.getParameter<edm::ParameterSet>("InitialScales"));
193  checkConfig(iniScales, "invalid set of scales");
194 
195  // Do we want to use the adaptive clustering tree algorithm?
196  const unsigned maxAdaptiveScales =
197  ps.getParameter<unsigned>("maxAdaptiveScales");
198  const double minAdaptiveRatioLog =
199  ps.getParameter<double>("minAdaptiveRatioLog");
200  if (minAdaptiveRatioLog <= 0.0)
201  throw cms::Exception("FFTJetBadConfig")
202  << "bad adaptive ratio logarithm limit" << std::endl;
203 
204  // Make sure that all standard scales are larger than the
205  // complete event scale
206  if (getEventScale() > 0.0)
207  {
208  const double cs = getEventScale();
209  const unsigned nscales = iniScales->size();
210  for (unsigned i=0; i<nscales; ++i)
211  if (cs >= (*iniScales)[i])
212  throw cms::Exception("FFTJetBadConfig")
213  << "incompatible scale for complete event" << std::endl;
214  }
215 
216  // At this point we are ready to build the clustering sequencer
217  sequencer = std::auto_ptr<Sequencer>(new Sequencer(
218  convolver.get(), peakSelector.get(), buildPeakFinder(ps),
219  *iniScales, maxAdaptiveScales, minAdaptiveRatioLog));
220 
221  // Build the clustering tree sparsifier
222  const edm::ParameterSet& SparsifierConfiguration(
223  ps.getParameter<edm::ParameterSet>("SparsifierConfiguration"));
225  SparsifierConfiguration);
226  checkConfig(sparsifier, "invalid sparsifier parameters");
227 
228  // Build the distance calculator for the clustering tree
229  const edm::ParameterSet& TreeDistanceCalculator(
230  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 
238 
240 {
241  // Check the parameter named "etaDependentScaleFactors". If the vector
242  // of scales is empty we will use 2d kernel, otherwise use 1d kernels
243  const std::vector<double> etaDependentScaleFactors(
244  ps.getParameter<std::vector<double> >("etaDependentScaleFactors"));
245 
246  // Make sure that the number of scale factors provided is correct
247  const bool use2dKernel = etaDependentScaleFactors.empty();
248  if (!use2dKernel)
249  if (etaDependentScaleFactors.size() != energyFlow->nEta())
250  throw cms::Exception("FFTJetBadConfig")
251  << "invalid number of eta-dependent scale factors"
252  << std::endl;
253 
254  // Get the eta and phi scales for the kernel(s)
255  double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
256  const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
257  if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
258  throw cms::Exception("FFTJetBadConfig")
259  << "invalid kernel scale" << std::endl;
260 
261  // FFT assumes that the grid extent in eta is 2*Pi. Adjust the
262  // kernel scale in eta to compensate.
263  kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
264 
265  // Are we going to try to fix the efficiency near detector edges?
266  const bool fixEfficiency = ps.getParameter<bool>("fixEfficiency");
267 
268  // Minimum and maximum eta bin for the convolver
269  unsigned convolverMinBin = 0, convolverMaxBin = 0;
270  if (fixEfficiency || !use2dKernel)
271  {
272  convolverMinBin = ps.getParameter<unsigned>("convolverMinBin");
273  convolverMaxBin = ps.getParameter<unsigned>("convolverMaxBin");
274  }
275 
276  if (use2dKernel)
277  {
278  // Build the FFT engine
279  engine = std::auto_ptr<MyFFTEngine>(
280  new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
281 
282  // 2d kernel
283  kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
284  new fftjet::DiscreteGauss2d(
285  kernelEtaScale, kernelPhiScale,
286  energyFlow->nEta(), energyFlow->nPhi()));
287 
288  // 2d convolver
289  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
290  new fftjet::FrequencyKernelConvolver<Real,Complex>(
291  engine.get(), kernel2d.get(),
292  convolverMinBin, convolverMaxBin));
293  }
294  else
295  {
296  // Need separate FFT engines for eta and phi
297  engine = std::auto_ptr<MyFFTEngine>(
298  new MyFFTEngine(1, energyFlow->nEta()));
299  anotherEngine = std::auto_ptr<MyFFTEngine>(
300  new MyFFTEngine(1, energyFlow->nPhi()));
301 
302  // 1d kernels
303  etaKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
304  new fftjet::DiscreteGauss1d(kernelEtaScale, energyFlow->nEta()));
305 
306  phiKernel = std::auto_ptr<fftjet::AbsFrequencyKernel1d>(
307  new fftjet::DiscreteGauss1d(kernelPhiScale, energyFlow->nPhi()));
308 
309  // Sequential convolver
310  convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
311  new fftjet::FrequencySequentialConvolver<Real,Complex>(
312  engine.get(), anotherEngine.get(),
313  etaKernel.get(), phiKernel.get(),
314  etaDependentScaleFactors, convolverMinBin,
315  convolverMaxBin, fixEfficiency));
316  }
317 }
318 
319 
321 {
322  const double peakFinderMaxEta = ps.getParameter<double>("peakFinderMaxEta");
323  if (peakFinderMaxEta <= 0.0)
324  throw cms::Exception("FFTJetBadConfig")
325  << "invalid peak finder eta cut" << std::endl;
326  const double maxMagnitude = ps.getParameter<double>("peakFinderMaxMagnitude");
327  int minBin = energyFlow->getEtaBin(-peakFinderMaxEta);
328  if (minBin < 0)
329  minBin = 0;
330  int maxBin = energyFlow->getEtaBin(peakFinderMaxEta) + 1;
331  if (maxBin < 0)
332  maxBin = 0;
333  return fftjet::PeakFinder(maxMagnitude, true, minBin, maxBin);
334 }
335 
336 
338 {
339  // do anything here that needs to be done at desctruction time
340  // (e.g. close files, deallocate resources etc.)
341  delete clusteringTree;
342 }
343 
344 
345 //
346 // member functions
347 //
348 template<class Real>
350 {
352 
353  std::auto_ptr<StoredTree> tree(new StoredTree());
354 
356  sequencer->maxAdaptiveScales(),
357  tree.get());
358 
359  // Check that we can restore the tree
361  {
363  const std::vector<double>& scalesUsed(sequencer->getInitialScales());
364  sparsePeakTreeFromStorable(*tree, &scalesUsed, getEventScale(), &check);
365  if (sparseTree != check)
366  throw cms::Exception("FFTJetInterface")
367  << "Data conversion failed for sparse clustering tree"
368  << std::endl;
369  }
370 
371  ev.put(tree, outputLabel);
372 }
373 
374 
375 template<class Real>
377 {
379 
380  std::auto_ptr<StoredTree> tree(new StoredTree());
381 
383  sequencer->maxAdaptiveScales(),
384  tree.get());
385 
386  // Check that we can restore the tree
388  {
390  const std::vector<double>& scalesUsed(sequencer->getInitialScales());
391  densePeakTreeFromStorable(*tree, &scalesUsed, getEventScale(), &check);
392  if (*clusteringTree != check)
393  throw cms::Exception("FFTJetInterface")
394  << "Data conversion failed for dense clustering tree"
395  << std::endl;
396  }
397 
398  ev.put(tree, outputLabel);
399 }
400 
401 
402 // ------------ method called to produce the data ------------
404  edm::Event& iEvent, const edm::EventSetup& iSetup)
405 {
406  loadInputCollection(iEvent);
408 
409  if (makeClusteringTree)
410  {
412  if (getEventScale() > 0.0)
413  sequencer->insertCompleteEvent(getEventScale(), *energyFlow,
415 
416  if (sparsify)
417  {
418  sparsifier->sparsify(*clusteringTree, &sparseTree);
419 
420  // Do not call the "sortNodes" method of the sparse tree here.
421  // Currently, the nodes are sorted by daughter number.
422  // This is the way we want it in storage because the stored
423  // tree does not include daughter ordering info explicitly.
424 
426  buildSparseProduct<float>(iEvent);
427  else
428  buildSparseProduct<double>(iEvent);
429  }
430  else
431  {
433  buildDenseProduct<float>(iEvent);
434  else
435  buildDenseProduct<double>(iEvent);
436  }
437  }
438 
440  {
441  const fftjet::Grid2d<Real>& g(*energyFlow);
442 
443  std::auto_ptr<DiscretizedEnergyFlow> flow(
445  g.data(), g.title(), g.etaMin(), g.etaMax(),
446  g.phiBin0Edge(), g.nEta(), g.nPhi()));
447 
449  {
450  fftjet::Grid2d<Real> check(
451  flow->nEtaBins(), flow->etaMin(), flow->etaMax(),
452  flow->nPhiBins(), flow->phiBin0Edge(), flow->title());
453  check.blockSet(flow->data(), flow->nEtaBins(), flow->nPhiBins());
454  assert(g == check);
455  }
456 
457  iEvent.put(flow, outputLabel);
458  }
459 }
460 
461 
462 // ------------ method called once each job just before starting event loop
464 {
465 }
466 
467 
468 // ------------ method called once each job just after ending the event loop
470 {
471 }
472 
473 
474 //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
std::auto_ptr< fftjet::Functor1< bool, fftjet::Peak > > peakSelector
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
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:84
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)
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
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
const std::string outputLabel
std::auto_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)