stablejet is hosted by Hepforge, IPPP Durham
StableJet
ClusterAnalysis.h
1 #ifndef FastJetStability_ClusterAnalysis_h
2 #define FastJetStability_ClusterAnalysis_h
3 
4 #include "JetStabilityDefs.h"
5 #include "ClusteringComparisonIndex.h"
6 #include "AbsFuzzinessCalculator.h"
7 
8 #include "fastjet/ClusterSequence.hh"
9 
10 namespace stab {
11 
15 enum StabilityCalcMethod {
16  PAIRWISE=0,
17  EXHAUSTIVE,
18  UNSMEARED,
19  RESAMPLED
20 };
21 
28 
29 public:
31  inline ClusterAnalysis (const double dCut)
32  :m_dCut(dCut)
33  ,m_nEvent(0U),m_nEventSize(100U),m_event(m_nEventSize)
34  ,m_nJetSize(50U),m_jets(m_nJetSize)
35  ,m_nJetsInCycle(0U),m_nJets(0U)
36  ,m_jetCycleMap(m_nJetSize),m_partMap(m_nEventSize,0U)
37  ,m_nCycle(0U)
38  { }
39 
41  inline ClusterAnalysis ()
42  :m_dCut(-1.)
43  ,m_nEvent(0U),m_nEventSize(100U),m_event(m_nEventSize)
44  ,m_nJetSize(50U),m_jets(m_nJetSize)
45  ,m_nJetsInCycle(0U),m_nJets(0U)
46  ,m_jetCycleMap(m_nJetSize),m_partMap(m_nEventSize,0U)
47  ,m_nCycle(0U)
48  { }
49 
50  inline virtual ~ClusterAnalysis () { }
51 
54  void initialize(const fastjet::ClusterSequence &);
55 
57  void clear();
58 
61  void fill(const fastjet::ClusterSequence &);
62 
64  double doStability(const ClusteringComparisonIndex &,const StabilityCalcMethod method=PAIRWISE);
65 
67  double doFuzziness(const AbsFuzzinessCalculator &);
68 
70  inline double dCut() const { return m_dCut; }
71 
72 private:
73  //---------------------------- private methods -------------------
75  double pairwiseStability(const ClusteringComparisonIndex &);
76  double exhaustiveStability(const ClusteringComparisonIndex &);
77  double unsmearedStability(const ClusteringComparisonIndex &);
78  double resampledStability(const ClusteringComparisonIndex &);
79 
80 
81  //---------------------------- private data ----------------------
82  // distance cut-off for this analysis
83  double m_dCut;
84 
85  unsigned m_nEvent; // number of particles in the event
86  unsigned m_nEventSize; // size of m_event vector
87  // the initial event
88  Event m_event;
89 
90  // list of clustered jets to keep track of each event
91  // This grows as it encounters more jets
92  // should perform one-to-one match with each fill
93  unsigned m_nJetSize; // size of capacity of m_jets
94  Event m_jets;
95 
96  unsigned m_nJetsInCycle;
97  unsigned m_nJets;
98 
99  // this vector points for each element in the cycle's jets to
100  // an index in m_jets
101  std::vector<unsigned> m_jetCycleMap;
102 
103  // map the initial event into the m_jets vector
104  std::vector<unsigned> m_partMap;
105  // a collection of m_partMap for each cycle of smearing
106  std::vector<std::vector<unsigned> > m_partMapVec;
107 
108  // the number of cycles encountered (i.e. the number of time fill has been called)
109  unsigned m_nCycle;
110 
111 
112 }; // end ClusterAnalysis declaration
113 
114 } // end stab namespace
115 
116 #endif
double doStability(const ClusteringComparisonIndex &, const StabilityCalcMethod method=PAIRWISE)
perform the stability calculation
ClusterAnalysis()
default constructor
double doFuzziness(const AbsFuzzinessCalculator &)
perform the fuzziness calculation
ClusterAnalysis(const double dCut)
Constructor. Initialize the distance paramter (scale cut-off) by dCut.
void fill(const fastjet::ClusterSequence &)
Add results of latest clustering to stability consideration.
double dCut() const
return the distance cut-off value
void initialize(const fastjet::ClusterSequence &)
Initialize internal data from first clustering results.
void clear()
clear the contents