CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Private Member Functions
FFTJetPatRecoProducer Class Reference

#include <RecoJets/FFTJetProducer/plugins/FFTJetPatRecoProducer.cc>

Inheritance diagram for FFTJetPatRecoProducer:
edm::EDProducer fftjetcms::FFTJetInterface edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

 FFTJetPatRecoProducer (const edm::ParameterSet &)
 
 ~FFTJetPatRecoProducer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from fftjetcms::FFTJetInterface
virtual ~FFTJetInterface ()
 

Protected Types

typedef
fftjet::ProximityClusteringTree
< fftjet::Peak, long > 
ClusteringTree
 
typedef
fftjet::ClusteringSequencer
< Real
Sequencer
 
typedef
fftjet::SparseClusteringTree
< fftjet::Peak, long > 
SparseTree
 
typedef
fftjet::ClusteringTreeSparsifier
< fftjet::Peak, long > 
Sparsifier
 

Protected Member Functions

void beginJob ()
 
template<class Real >
void buildDenseProduct (edm::Event &) const
 
void buildKernelConvolver (const edm::ParameterSet &)
 
fftjet::PeakFinder buildPeakFinder (const edm::ParameterSet &)
 
template<class Real >
void buildSparseProduct (edm::Event &) const
 
void endJob ()
 
void produce (edm::Event &, const edm::EventSetup &)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 
- Protected Member Functions inherited from fftjetcms::FFTJetInterface
template<class Ptr >
void checkConfig (const Ptr &ptr, const char *message)
 
void discretizeEnergyFlow ()
 
 FFTJetInterface (const edm::ParameterSet &)
 
double getEventScale () const
 
void loadInputCollection (const edm::Event &)
 
bool storeInSinglePrecision () const
 
const reco::Particle::PointvertexUsed () const
 

Protected Attributes

std::auto_ptr< MyFFTEngineanotherEngine
 
ClusteringTreeclusteringTree
 
const double completeEventDataCutoff
 
std::auto_ptr
< fftjet::AbsConvolverBase
< Real > > 
convolver
 
std::auto_ptr
< fftjet::AbsDistanceCalculator
< fftjet::Peak > > 
distanceCalc
 
std::auto_ptr< MyFFTEngineengine
 
std::auto_ptr
< fftjet::AbsFrequencyKernel1d > 
etaKernel
 
std::auto_ptr
< fftjet::AbsFrequencyKernel > 
kernel2d
 
const bool makeClusteringTree
 
std::auto_ptr
< fftjet::Functor1< bool,
fftjet::Peak > > 
peakSelector
 
std::auto_ptr
< fftjet::AbsFrequencyKernel1d > 
phiKernel
 
std::auto_ptr< Sequencersequencer
 
SparseTree sparseTree
 
std::auto_ptr< Sparsifiersparsifier
 
const bool sparsify
 
const bool storeDiscretizationGrid
 
const bool verifyDataConversion
 
- Protected Attributes inherited from fftjetcms::FFTJetInterface
const AnomalousTower anomalous
 
std::vector< unsigned > candidateIndex
 
const bool doPVCorrection
 
std::auto_ptr< fftjet::Grid2d
< fftjetcms::Real > > 
energyFlow
 
const std::vector< double > etaDependentMagnutideFactors
 
std::vector
< fftjetcms::VectorLike
eventData
 
edm::Handle< reco::CandidateViewinputCollection
 
const edm::InputTag inputLabel
 
const JetType jetType
 
const std::string outputLabel
 
const edm::InputTag srcPVs
 

Private Member Functions

 FFTJetPatRecoProducer ()
 
 FFTJetPatRecoProducer (const FFTJetPatRecoProducer &)
 
FFTJetPatRecoProduceroperator= (const FFTJetPatRecoProducer &)
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from fftjetcms::FFTJetInterface
enum  JetType {
  BASICJET = 0, GENJET, CALOJET, PFJET,
  TRACKJET, JPTJET
}
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from fftjetcms::FFTJetInterface
static JetType parse_jet_type (const std::string &name)
 

Detailed Description

Description: Runs FFTJet pattern recognition stage and saves the results

Implementation: [Notes on implementation]

Definition at line 63 of file FFTJetPatRecoProducer.cc.

Member Typedef Documentation

typedef fftjet::ProximityClusteringTree<fftjet::Peak,long> FFTJetPatRecoProducer::ClusteringTree
protected

Definition at line 71 of file FFTJetPatRecoProducer.cc.

typedef fftjet::ClusteringSequencer<Real> FFTJetPatRecoProducer::Sequencer
protected

Definition at line 73 of file FFTJetPatRecoProducer.cc.

typedef fftjet::SparseClusteringTree<fftjet::Peak,long> FFTJetPatRecoProducer::SparseTree
protected

Definition at line 72 of file FFTJetPatRecoProducer.cc.

typedef fftjet::ClusteringTreeSparsifier<fftjet::Peak,long> FFTJetPatRecoProducer::Sparsifier
protected

Definition at line 74 of file FFTJetPatRecoProducer.cc.

Constructor & Destructor Documentation

FFTJetPatRecoProducer::FFTJetPatRecoProducer ( const edm::ParameterSet ps)
explicit

Definition at line 147 of file FFTJetPatRecoProducer.cc.

References buildKernelConvolver(), buildPeakFinder(), fftjetcms::FFTJetInterface::checkConfig(), clusteringTree, convolver, distanceCalc, fftjetcms::FFTJetInterface::energyFlow, edm::hlt::Exception, fftjetcms::fftjet_ClusteringTreeSparsifier_parser(), fftjetcms::fftjet_DistanceCalculator_parser(), fftjetcms::fftjet_Grid2d_parser(), fftjetcms::fftjet_PeakSelector_parser(), fftjetcms::fftjet_ScaleSet_parser(), fftjetcms::FFTJetInterface::getEventScale(), edm::ParameterSet::getParameter(), i, makeClusteringTree, fftjetcms::FFTJetInterface::outputLabel, peakSelector, sequencer, sparsifier, storeDiscretizationGrid, and fftjetcms::FFTJetInterface::storeInSinglePrecision().

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 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T 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
std::auto_ptr< fftjet::Grid2d< Real > > fftjet_Grid2d_parser(const edm::ParameterSet &ps)
void buildKernelConvolver(const edm::ParameterSet &)
void checkConfig(const Ptr &ptr, const char *message)
std::auto_ptr< fftjet::Functor1< bool, fftjet::Peak > > fftjet_PeakSelector_parser(const edm::ParameterSet &ps)
ClusteringTree * clusteringTree
fftjet::ClusteringSequencer< Real > Sequencer
std::auto_ptr< Sparsifier > sparsifier
std::auto_ptr< fftjet::AbsDistanceCalculator< fftjet::Peak > > distanceCalc
bool storeInSinglePrecision() const
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::ProximityClusteringTree< fftjet::Peak, long > ClusteringTree
const std::string outputLabel
std::auto_ptr< std::vector< double > > fftjet_ScaleSet_parser(const edm::ParameterSet &ps)
FFTJetPatRecoProducer::~FFTJetPatRecoProducer ( )

Definition at line 337 of file FFTJetPatRecoProducer.cc.

References clusteringTree.

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 }
ClusteringTree * clusteringTree
FFTJetPatRecoProducer::FFTJetPatRecoProducer ( )
private
FFTJetPatRecoProducer::FFTJetPatRecoProducer ( const FFTJetPatRecoProducer )
private

Member Function Documentation

void FFTJetPatRecoProducer::beginJob ( void  )
protectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 463 of file FFTJetPatRecoProducer.cc.

464 {
465 }
template<class Real >
void FFTJetPatRecoProducer::buildDenseProduct ( edm::Event ev) const
protected

Definition at line 376 of file FFTJetPatRecoProducer.cc.

References CastorDataFrameFilter_impl::check(), clusteringTree, fftjetcms::densePeakTreeFromStorable(), fftjetcms::densePeakTreeToStorable(), distanceCalc, edm::hlt::Exception, fftjetcms::FFTJetInterface::getEventScale(), fftjetcms::FFTJetInterface::outputLabel, edm::Event::put(), sequencer, fftjetcms::FFTJetInterface::storeInSinglePrecision(), diffTreeTool::tree, and verifyDataConversion.

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());
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 }
Class for storing FFTJet sparse clustering trees.
Definition: PattRecoTree.h:19
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
bool check(const DataFrame &df, bool capcheck, bool dvercheck)
ClusteringTree * clusteringTree
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 densePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::AbsClusteringTree< fftjet::Peak, long > *out)
std::auto_ptr< Sequencer > sequencer
fftjet::ProximityClusteringTree< fftjet::Peak, long > ClusteringTree
const std::string outputLabel
void FFTJetPatRecoProducer::buildKernelConvolver ( const edm::ParameterSet ps)
protected

Definition at line 239 of file FFTJetPatRecoProducer.cc.

References anotherEngine, convolver, fftjetcms::FFTJetInterface::energyFlow, engine, etaKernel, edm::hlt::Exception, edm::ParameterSet::getParameter(), kernel2d, M_PI, and phiKernel.

Referenced by FFTJetPatRecoProducer().

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 }
T getParameter(std::string const &) const
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
std::auto_ptr< fftjet::AbsFrequencyKernel1d > phiKernel
std::auto_ptr< MyFFTEngine > anotherEngine
std::auto_ptr< fftjet::AbsFrequencyKernel > kernel2d
std::auto_ptr< MyFFTEngine > engine
#define M_PI
Definition: BFit3D.cc:3
std::auto_ptr< fftjet::AbsFrequencyKernel1d > etaKernel
std::auto_ptr< fftjet::AbsConvolverBase< Real > > convolver
fftjet::FFTWDoubleEngine MyFFTEngine
fftjet::PeakFinder FFTJetPatRecoProducer::buildPeakFinder ( const edm::ParameterSet ps)
protected

Definition at line 320 of file FFTJetPatRecoProducer.cc.

References fftjetcms::FFTJetInterface::energyFlow, edm::hlt::Exception, and edm::ParameterSet::getParameter().

Referenced by FFTJetPatRecoProducer().

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 }
T getParameter(std::string const &) const
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
template<class Real >
void FFTJetPatRecoProducer::buildSparseProduct ( edm::Event ev) const
protected

Definition at line 349 of file FFTJetPatRecoProducer.cc.

References CastorDataFrameFilter_impl::check(), edm::hlt::Exception, fftjetcms::FFTJetInterface::getEventScale(), fftjetcms::FFTJetInterface::outputLabel, edm::Event::put(), sequencer, fftjetcms::sparsePeakTreeFromStorable(), fftjetcms::sparsePeakTreeToStorable(), sparseTree, fftjetcms::FFTJetInterface::storeInSinglePrecision(), diffTreeTool::tree, and verifyDataConversion.

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 }
Class for storing FFTJet sparse clustering trees.
Definition: PattRecoTree.h:19
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)
bool storeInSinglePrecision() const
void sparsePeakTreeFromStorable(const reco::PattRecoTree< Real, reco::PattRecoPeak< Real > > &in, const std::vector< double > *scaleSetIfNotAdaptive, double completeEventScale, fftjet::SparseClusteringTree< fftjet::Peak, long > *out)
std::auto_ptr< Sequencer > sequencer
fftjet::SparseClusteringTree< fftjet::Peak, long > SparseTree
const std::string outputLabel
void FFTJetPatRecoProducer::endJob ( void  )
protectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 469 of file FFTJetPatRecoProducer.cc.

470 {
471 }
FFTJetPatRecoProducer& FFTJetPatRecoProducer::operator= ( const FFTJetPatRecoProducer )
private
void FFTJetPatRecoProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
protectedvirtual

Implements edm::EDProducer.

Definition at line 403 of file FFTJetPatRecoProducer.cc.

References CastorDataFrameFilter_impl::check(), clusteringTree, completeEventDataCutoff, fftjetcms::FFTJetInterface::discretizeEnergyFlow(), fftjetcms::FFTJetInterface::energyFlow, g, fftjetcms::FFTJetInterface::getEventScale(), iEvent, fftjetcms::FFTJetInterface::loadInputCollection(), makeClusteringTree, fftjetcms::FFTJetInterface::outputLabel, edm::Event::put(), sequencer, sparseTree, sparsifier, sparsify, storeDiscretizationGrid, fftjetcms::FFTJetInterface::storeInSinglePrecision(), and verifyDataConversion.

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 }
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
void loadInputCollection(const edm::Event &)
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
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
bool check(const DataFrame &df, bool capcheck, bool dvercheck)
ClusteringTree * clusteringTree
std::auto_ptr< Sparsifier > sparsifier
bool storeInSinglePrecision() const
std::auto_ptr< Sequencer > sequencer
const std::string outputLabel

Member Data Documentation

std::auto_ptr<MyFFTEngine> FFTJetPatRecoProducer::anotherEngine
protected

Definition at line 101 of file FFTJetPatRecoProducer.cc.

Referenced by buildKernelConvolver().

ClusteringTree* FFTJetPatRecoProducer::clusteringTree
protected
const double FFTJetPatRecoProducer::completeEventDataCutoff
protected

Definition at line 123 of file FFTJetPatRecoProducer.cc.

Referenced by produce().

std::auto_ptr<fftjet::AbsConvolverBase<Real> > FFTJetPatRecoProducer::convolver
protected

Definition at line 109 of file FFTJetPatRecoProducer.cc.

Referenced by buildKernelConvolver(), and FFTJetPatRecoProducer().

std::auto_ptr<fftjet::AbsDistanceCalculator<fftjet::Peak> > FFTJetPatRecoProducer::distanceCalc
protected

Definition at line 115 of file FFTJetPatRecoProducer.cc.

Referenced by buildDenseProduct(), and FFTJetPatRecoProducer().

std::auto_ptr<MyFFTEngine> FFTJetPatRecoProducer::engine
protected

Definition at line 100 of file FFTJetPatRecoProducer.cc.

Referenced by buildKernelConvolver().

std::auto_ptr<fftjet::AbsFrequencyKernel1d> FFTJetPatRecoProducer::etaKernel
protected

Definition at line 105 of file FFTJetPatRecoProducer.cc.

Referenced by buildKernelConvolver().

std::auto_ptr<fftjet::AbsFrequencyKernel> FFTJetPatRecoProducer::kernel2d
protected

Definition at line 104 of file FFTJetPatRecoProducer.cc.

Referenced by buildKernelConvolver().

const bool FFTJetPatRecoProducer::makeClusteringTree
protected

Definition at line 126 of file FFTJetPatRecoProducer.cc.

Referenced by FFTJetPatRecoProducer(), and produce().

std::auto_ptr<fftjet::Functor1<bool,fftjet::Peak> > FFTJetPatRecoProducer::peakSelector
protected

Definition at line 112 of file FFTJetPatRecoProducer.cc.

Referenced by FFTJetPatRecoProducer().

std::auto_ptr<fftjet::AbsFrequencyKernel1d> FFTJetPatRecoProducer::phiKernel
protected

Definition at line 106 of file FFTJetPatRecoProducer.cc.

Referenced by buildKernelConvolver().

std::auto_ptr<Sequencer> FFTJetPatRecoProducer::sequencer
protected
SparseTree FFTJetPatRecoProducer::sparseTree
protected

Definition at line 118 of file FFTJetPatRecoProducer.cc.

Referenced by buildSparseProduct(), and produce().

std::auto_ptr<Sparsifier> FFTJetPatRecoProducer::sparsifier
protected

Definition at line 97 of file FFTJetPatRecoProducer.cc.

Referenced by FFTJetPatRecoProducer(), and produce().

const bool FFTJetPatRecoProducer::sparsify
protected

Definition at line 136 of file FFTJetPatRecoProducer.cc.

Referenced by produce().

const bool FFTJetPatRecoProducer::storeDiscretizationGrid
protected

Definition at line 133 of file FFTJetPatRecoProducer.cc.

Referenced by FFTJetPatRecoProducer(), and produce().

const bool FFTJetPatRecoProducer::verifyDataConversion
protected

Definition at line 130 of file FFTJetPatRecoProducer.cc.

Referenced by buildDenseProduct(), buildSparseProduct(), and produce().