CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FFTJetProducerSummary.h
Go to the documentation of this file.
1 
9 #ifndef DataFormats_JetReco_FFTJetProducerSummary_h
10 #define DataFormats_JetReco_FFTJetProducerSummary_h
11 
12 #include <vector>
13 
16 
17 namespace reco {
19  {
20  public:
21  // Must have a default constructor
23  : unused_(0.f), minScale_(0.f), maxScale_(0.f), scaleUsed_(0.f),
25 
26  // Meaningful constructor
27  FFTJetProducerSummary(const std::vector<double>& thresholds,
28  const std::vector<unsigned>& levelOccupancy,
30  const std::vector<CandidatePtr>& constituents,
31  double unused, double minScale, double maxScale,
32  double scaleUsed, unsigned preclustersFound,
33  unsigned iterationsPerformed, bool converged);
34 
35  // Vector of occupancy thresholds generated by
36  // the locally adaptive resolution scheme
37  inline const std::vector<float>& thresholds() const
38  {return thresholds_;}
39 
40  // Clustering tree occupancy as a function of the level number
41  inline const std::vector<unsigned>& levelOccupancy() const
42  {return levelOccupancy_;}
43 
44  // Raw unclustered energy 4-vector
45  inline const math::XYZTLorentzVector& unclustered() const
46  {return unclustered_;}
47 
48  // Pointers to unclustered energy constituents
49  inline const std::vector<CandidatePtr>& unclusteredConstituents() const
50  {return unclusConstituents_;}
51 
52  // Scalar sum of unclustered transverse energies
53  inline float unusedEt() const {return unused_;}
54 
55  // Minimum and maximum resolution scales to be used for
56  // configuration stability calculations
57  inline float minScale() const {return minScale_;}
58  inline float maxScale() const {return maxScale_;}
59 
60  // The scale actually used to lookup clusters.
61  // This info is useful mostly for maximally stable
62  // and globally adaptive resolution schemes.
63  inline float scaleUsed() const {return scaleUsed_;}
64 
65  // Number of preclusters found by the peak selection code.
66  // Should normally coincide with the number of jets.
67  inline unsigned preclustersFound() const {return preclustersFound_;}
68 
69  // Actual number of iterations made by the jet reconstruction
70  // algorithm
71  inline unsigned iterationsPerformed() const
72  {return iterationsPerformed_;}
73 
74  // Did the jet reconstruction algorithm converge?
75  inline bool iterationsConverged() const {return converged_;}
76 
77  private:
78  std::vector<float> thresholds_;
79  std::vector<unsigned> levelOccupancy_;
81  std::vector<CandidatePtr> unclusConstituents_;
82  float unused_;
83  float minScale_;
84  float maxScale_;
85  float scaleUsed_;
88  bool converged_;
89  };
90 }
91 
92 #endif // DataFormats_JetReco_FFTJetProducerSummary_h
const std::vector< float > & thresholds() const
const std::vector< CandidatePtr > & unclusteredConstituents() const
std::vector< unsigned > levelOccupancy_
Data processing summary generated by FFTJetProducer.
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
const std::vector< unsigned > & levelOccupancy() const
std::vector< CandidatePtr > unclusConstituents_
double f[11][100]
std::vector< float > thresholds_
math::XYZTLorentzVector unclustered_
const math::XYZTLorentzVector & unclustered() const