CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FFTJetInterface.h
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: FFTJetProducers
4 // Class: FFTJetInterface
5 //
11 //
12 // Original Author: Igor Volobouev
13 // Created: June 29 2010
14 //
15 //
16 
17 #ifndef RecoJets_FFTJetProducers_FFTJetInterface_h
18 #define RecoJets_FFTJetProducers_FFTJetInterface_h
19 
20 #include <memory>
21 #include <vector>
22 #include <cassert>
23 
24 // FFTJet headers
25 #include "fftjet/Grid2d.hh"
26 
27 // framework include files
31 
33 
35 
39 
40 // local FFTJet-related definitions
43 
44 //
45 // class declaration
46 //
47 namespace fftjetcms {
49  {
50  public:
51  virtual ~FFTJetInterface() {}
52 
53  protected:
54  explicit FFTJetInterface(const edm::ParameterSet&);
55 
56  template<class Ptr>
57  void checkConfig(const Ptr& ptr, const char* message)
58  {
59  if (ptr.get() == NULL)
60  throw cms::Exception("FFTJetBadConfig") << message << std::endl;
61  }
62 
63  void loadInputCollection(const edm::Event&);
64  void discretizeEnergyFlow();
65  double getEventScale() const;
66  bool storeInSinglePrecision() const;
67 
68  const reco::Particle::Point& vertexUsed() const {return vertex_;}
69 
70  // Label for the input collection
72 
73  // Label for the output objects
75 
76  // Jet type to produce
78 
79  // Vertex correction-related stuff
80  const bool doPVCorrection;
81 
82  // Label for the vertices
84 
85  // Try to equalize magnitudes in the energy flow grid?
86  const std::vector<double> etaDependentMagnutideFactors;
87 
88  // Functor for finding anomalous towers
90 
91  // Event data 4-vectors
92  std::vector<fftjetcms::VectorLike> eventData;
93 
94  // Candidate number which corresponds to the item in the "eventData"
95  std::vector<unsigned> candidateIndex;
96 
97  // The energy discretization grid
98  std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> > energyFlow;
99 
100  // The input handle for the collection of candidates
102 
103  private:
104  // Explicitly disable other ways to construct this object
105  FFTJetInterface();
108 
110  const double completeEventScale;
112  };
113 }
114 
115 #endif // RecoJets_FFTJetProducers_FFTJetInterface_h
std::auto_ptr< fftjet::Grid2d< fftjetcms::Real > > energyFlow
const edm::InputTag inputLabel
edm::Handle< reco::CandidateView > inputCollection
const reco::Particle::Point & vertexUsed() const
void loadInputCollection(const edm::Event &)
const edm::InputTag srcPVs
#define NULL
Definition: scimark2.h:8
void checkConfig(const Ptr &ptr, const char *message)
std::vector< unsigned > candidateIndex
reco::Particle::Point vertex_
const std::vector< double > etaDependentMagnutideFactors
math::XYZPoint Point
point in the space
Definition: Particle.h:31
const AnomalousTower anomalous
bool storeInSinglePrecision() const
FFTJetInterface & operator=(const FFTJetInterface &)
std::vector< fftjetcms::VectorLike > eventData
const std::string outputLabel