stablejet is hosted by Hepforge, IPPP Durham
StableJet
ClusterAnalysis.cc
1 
2 #include <climits>
3 
4 #include "ClusterAnalysis.h"
5 
6 #include "matchOneToOne.hh"
7 
8 namespace stab {
9 
11 void ClusterAnalysis::initialize(const fastjet::ClusterSequence &cs)
12 {
13 
14  // clear the class
15  clear();
16 
17  // how many particles do we have
18  m_nEvent = cs.n_particles();
19 
20  // check if there is room in our maps and vectors
21  if ( m_nEvent > m_nEventSize ) {
22  m_nEventSize = m_nEvent;
23  m_event.resize(m_nEvent);
24  m_partMap.resize(m_nEvent);
25  }
26 
27  // fill the initial event into m_event
28  const Event & allJets = cs.jets();
29  const History &history = cs.history();
30  for ( unsigned jp=0; jp != m_nEvent; jp++ ) {
31  const int jpHistIndex = allJets[jp].cluster_sequence_history_index();
32  assert(jpHistIndex >= 0);
33  const int jpParent = history[jpHistIndex].parent1;
34  //assert(jpParent > 0);
35  m_event[jp] = allJets[jp];
36  assert( m_event[jp].pt() > 0. );
37  }
38 
39 
40  // perform initial clustering
41  const Event & cluJets = fastjet::sorted_by_pt(cs.exclusive_jets(m_dCut));
42 
43  // make sure there is room in m_jets and the like
44  m_nJets=cluJets.size();
45  m_nJetsInCycle = m_nJets;
46  if ( m_nJetSize < m_nJets ) {
47  m_jets.resize(m_nJets);
48  m_jetCycleMap.resize(m_nJets);
49  }
50 
51  const std::vector<int> & histMap = cs.particle_jet_indices(cluJets);
52  assert(m_nEvent == histMap.size());
53  for ( unsigned p=0; p != m_nEvent; p++ ) {
54  const int index = histMap[p];
55  if ( index >= 0 ) m_partMap[p] = index;
56  else {
57  m_partMap[p] = UINT_MAX;
58  }
59  }
60 
61 }
62 
63 
65 void ClusterAnalysis::fill(const fastjet::ClusterSequence &cs)
66 {
67 
68  // get the clustered jets
69  const Event & cluJets = fastjet::sorted_by_pt(cs.exclusive_jets(m_dCut));
70  m_nJetsInCycle = cluJets.size();
71 
72  // make sure we don't exceed the number of jets seen thus far
73  if ( m_nJetsInCycle > m_nJetSize ) {
74  m_nJetSize = m_nJetsInCycle;
75  m_jets.resize(m_nJetSize);
76  m_jetCycleMap.resize(m_nJetSize);
77  }
78 
79  // map cycle's jets to m_jets indices
80  for ( unsigned j=0; j != m_nJetsInCycle; j++ )
81  m_jetCycleMap[j]=0U;
82 
83  // match one-to-one
84  std::vector<int> matchIndex;
85  deltaRCalc calc;
86  matchOneToOne(cluJets,m_jets,m_nJetsInCycle,m_nJets,calc,&matchIndex);
87 
88  for ( unsigned j=0; j != m_nJetsInCycle; j++ ) {
89  const int index = matchIndex[j];
90  if ( index >=0 ) m_jetCycleMap[j] = index;
91  else { // we have encountered a new jet
92  if ( m_nJetSize == m_nJets ) m_jets.resize(++m_nJetSize);
93  m_jets[m_nJets] = cluJets[j];
94  m_jetCycleMap[j] = m_nJets++;
95  }
96  }
97 
98  // fill particle maps and histograms
99  std::vector<unsigned> partMap(m_nEvent);
100  const std::vector<int> & histMap = cs.particle_jet_indices(cluJets);
101  assert(m_nEvent == histMap.size());
102  for ( unsigned p=0; p != m_nEvent; p++ ) {
103  const int index = histMap[p];
104  if ( index >= 0 ) partMap[p] = index;
105  else partMap[p] = UINT_MAX;
106  }
107  m_partMapVec.push_back(partMap);
108 
109  // increment m_nCycle
110  m_nCycle++;
111 
112  assert( m_nJetsInCycle );
113 }
114 
116 {
117  m_nCycle=0U;
118 
119  // clear event
120  for ( unsigned p=0; p != m_nEventSize; p++ ) {
121  m_event[p] = PseudoJet();
122  m_partMap[p] = 0U;
123  }
124  m_nEvent = 0U;
125  m_partMapVec.clear();
126 
127  for ( unsigned j=0; j != m_nJetSize; j++ ) {
128  m_jets[j] = PseudoJet();
129  m_jetCycleMap[j] = 0U;
130  }
131  m_nJets = 0U;
132 
133 }
134 
135 
136 double ClusterAnalysis::doStability(const ClusteringComparisonIndex& index, const StabilityCalcMethod method)
137 {
138  double stab=0.;
139 
140  switch(method) {
141  case PAIRWISE:
142  stab = pairwiseStability(index);
143  break;
144  case EXHAUSTIVE:
145  stab = exhaustiveStability(index);
146  break;
147  case UNSMEARED:
148  stab = unsmearedStability(index);
149  break;
150  case RESAMPLED:
151  stab = resampledStability(index);
152  break;
153  }
154 
155  return stab;
156 }
157 
158 
159 double ClusterAnalysis::pairwiseStability(const ClusteringComparisonIndex &index)
160 {
161 
162  long double stab=0.L;
163 
164  for ( unsigned c=0; c != m_nCycle-1; c++ ) {
165  std::vector<unsigned> & clustMap1 = m_partMapVec[c];
166  std::vector<unsigned> & clustMap2 = m_partMapVec[c+1U];
167  const double temp = index.index(&m_event[0],&clustMap1[0],&clustMap2[0],m_nEvent);
168  //stab += index.index(&m_event[0],&clustMap1[0],&clustMap2[0],m_nEvent);
169  stab += temp;
170  }
171 
172  stab /= (m_nCycle-1);
173  //std::cout << "stab=" << stab << std::endl;
174 
175  return stab;
176 
177 }
178 
179 double ClusterAnalysis::exhaustiveStability(const ClusteringComparisonIndex &index)
180 {
181 
182  double stab=0.;
183 
184  for ( unsigned c1=0; c1 != m_nCycle-1U; c1++ ) {
185  std::vector<unsigned> & clustMap1 = m_partMapVec[c1];
186  for ( unsigned c2=c1+1; c2 != m_nCycle; c2++ ) {
187  std::vector<unsigned> & clustMap2 = m_partMapVec[c2];
188  stab += index.index(&m_event[0],&clustMap1[0],&clustMap2[0],m_nEvent);
189  }
190  }
191 
192  stab /= (double)m_nCycle*((double)m_nCycle-1.)/2.;
193  return stab;
194 }
195 
196 double ClusterAnalysis::unsmearedStability(const ClusteringComparisonIndex &index)
197 {
198  double stab=0.;
199 
200  for ( unsigned c=0; c != m_nCycle; c++ ) {
201  std::vector<unsigned> & clustMap = m_partMapVec[c];
202  stab += index.index(&m_event[0],&m_partMap[0],&clustMap[0],m_nEvent);
203  }
204  stab /= m_nCycle;
205  return stab;
206 }
207 
208 double ClusterAnalysis::resampledStability(const ClusteringComparisonIndex &index)
209 {
210  return 0.;
211 }
212 
214 {
215 
216  // construct weights
217  const unsigned weightSize = m_nJets*m_nEvent;
218  double weights[weightSize];
219 
220  for ( unsigned i=0; i != weightSize; i++ )
221  weights[i] = 0.0;
222 
223  // fill weights as probability of initial particle to go into a final jet
224  // ** this could use some optimization **
225  for ( unsigned i=0; i != m_nJets; i++ ) {
226  for ( unsigned j=0; j != m_nEvent; j++ ) {
227  unsigned count = 0;
228  for ( unsigned k=0; k != m_nCycle; k++ )
229  if ( m_partMapVec[k][j] == i ) count++;
230  weights[i*m_nEvent+j] = (double)count/(double)m_nCycle;
231  assert( weights[i*m_nEvent+j] >= 0. );
232  }
233  }
234 
235  bool lastFinalJetIsUnclustered=false;
236  const unsigned lenJetFuzziness = m_nJets;
237  double jetFuzziness[lenJetFuzziness];
238 
239  return fuzzyCalc.calculate(m_event, m_jets, weights, m_nJets, m_nEvent, lastFinalJetIsUnclustered,
240  jetFuzziness, lenJetFuzziness);
241 
242  return 0.;
243 }
244 
245 } // end stab namespace
double doStability(const ClusteringComparisonIndex &, const StabilityCalcMethod method=PAIRWISE)
perform the stability calculation
double doFuzziness(const AbsFuzzinessCalculator &)
perform the fuzziness calculation
void fill(const fastjet::ClusterSequence &)
Add results of latest clustering to stability consideration.
void initialize(const fastjet::ClusterSequence &)
Initialize internal data from first clustering results.
void clear()
clear the contents