stablejet is hosted by Hepforge, IPPP Durham
StableJet
ClusterAnalysisHighPt.h
1 #ifndef FastJetStability_ClusterAnalysisHighPt_h
2 #define FastJetStability_ClusterAnalysisHighPt_h
3 
4 #include "JetStabilityDefs.h"
5 #include "ClusteringComparisonIndexCP.h"
6 #include "AbsFuzzinessCalculator.h"
7 #include "ClusterAnalysis.h"
8 
9 #include "fastjet/ClusterSequence.hh"
10 
11 namespace stab {
12 
23 
24 public:
26  inline ClusterAnalysisHighPt (const double dCut, const bool useJunkJet = true)
27  :m_dCut(dCut)
28  ,m_nEvent(0U),m_nEventSize(100U),m_event(m_nEventSize)
29  ,m_nJetSize(50U),m_jets(m_nJetSize)
30  ,m_nJetsInCycle(0U),m_nJets(0U)
31  ,m_jetCycleMap(m_nJetSize),m_partMap(m_nEventSize,0U)
32  ,m_nCycle(0U),m_useJunkJet(useJunkJet)
33  { }
34 
37  :m_dCut(-1.)
38  ,m_nEvent(0U),m_nEventSize(100U),m_event(m_nEventSize)
39  ,m_nJetSize(50U),m_jets(m_nJetSize)
40  ,m_nJetsInCycle(0U),m_nJets(0U)
41  ,m_jetCycleMap(m_nJetSize),m_partMap(m_nEventSize,0U)
42  ,m_nCycle(0U),m_useJunkJet(true)
43  { }
44 
45  inline virtual ~ClusterAnalysisHighPt () { }
46 
49  void initialize(const fastjet::ClusterSequence &);
50 
52  void clear();
53 
56  void fill(const fastjet::ClusterSequence &);
57 
59  double doStability(const ClusteringComparisonIndexCP &,const StabilityCalcMethod method=PAIRWISE);
60 
62  double doFuzziness(const AbsFuzzinessCalculator &);
63 
65  inline double dCut() const { return m_dCut; }
66 
68  inline const std::vector<double> & cpStab() const { return m_cpStab; }
69 
70 private:
71  //---------------------------- private methods -------------------
73  double pairwiseStability(const ClusteringComparisonIndexCP &);
74  double exhaustiveStability(const ClusteringComparisonIndexCP &);
75  double unsmearedStability(const ClusteringComparisonIndexCP &);
76  double resampledStability(const ClusteringComparisonIndexCP &);
77 
78 
79  //---------------------------- private data ----------------------
80  // distance cut-off for this analysis
81  double m_dCut;
82 
83  unsigned m_nEvent; // number of particles in the event
84  unsigned m_nEventSize; // size of m_event vector
85  // the initial event
86  Event m_event;
87 
88  // list of clustered jets to keep track of each event
89  // This grows as it encounters more jets
90  // should perform one-to-one match with each fill
91  unsigned m_nJetSize; // size of capacity of m_jets
92  Event m_jets;
93 
94  unsigned m_nJetsInCycle;
95  unsigned m_nJets;
96 
97  // this vector points for each element in the cycle's jets to
98  // an index in m_jets
99  std::vector<unsigned> m_jetCycleMap;
100 
101  // map the initial event into the m_jets vector
102  std::vector<unsigned> m_partMap;
103  // a collection of m_partMap for each cycle of smearing
104  std::vector<std::vector<unsigned> > m_partMapVec;
105 
106  // constituent perspective clustering stability map
107  std::vector<double> m_cpStab;
108 
109  // the number of cycles encountered (i.e. the number of time fill has been called)
110  unsigned m_nCycle;
111 
112  // use a third jet for everything outside of the first in pt
113  bool m_useJunkJet;
114 
115 
116 }; // end ClusterAnalysisHighPt declaration
117 
118 } // end stab namespace
119 
120 #endif
void clear()
clear the contents
ClusterAnalysisHighPt()
default constructor
void fill(const fastjet::ClusterSequence &)
Add results of latest clustering to stability consideration.
double doStability(const ClusteringComparisonIndexCP &, const StabilityCalcMethod method=PAIRWISE)
perform the stability calculation
double doFuzziness(const AbsFuzzinessCalculator &)
perform the fuzziness calculation
double dCut() const
return the distance cut-off value
ClusterAnalysisHighPt(const double dCut, const bool useJunkJet=true)
Constructor. Initialize the distance paramter (scale cut-off) by dCut.
void initialize(const fastjet::ClusterSequence &)
Initialize internal data from first clustering results.
const std::vector< double > & cpStab() const
return a reference to the constituent perspective stability vector