CMS 3D CMS Logo

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