Public Member Functions | Private Member Functions | Private Attributes | List of all members
util::ClusterShapes Class Reference

#include <ClusterShapes.h>

Public Member Functions

 ClusterShapes (int nhits, float *a, float *x, float *y, float *z)
 
 ClusterShapes (int nhits, float *a, float *t, float *x, float *y, float *z)
 
 ClusterShapes (int nhits, std::vector< float > a, std::vector< float > x, std::vector< float > y, std::vector< float > z)
 
 ClusterShapes (int nhits, std::vector< float > a, std::vector< float > t, std::vector< float > x, std::vector< float > y, std::vector< float > z)
 
 ~ClusterShapes ()
 
int getNumberOfHits ()
 
float getTotalAmplitude ()
 
float getAverageTime ()
 
float * getCenterOfGravity ()
 
float * getEigenValInertia ()
 
float * getEigenVecInertia ()
 
float getWidth ()
 
float radius ()
 
float getElipsoid_r1 ()
 
float getElipsoid_r2 ()
 
float getElipsoid_r3 ()
 
float getElipsoid_vol ()
 
float getElipsoid_r_ave ()
 
float getElipsoid_density ()
 
float getElipsoid_eccentricity ()
 
float getElipsoid_r_forw ()
 
float getElipsoid_r_back ()
 

Private Member Functions

void findElipsoid ()
 
void findGravity ()
 
void findInertia ()
 
void findWidth ()
 
float findDistance (int i)
 

Private Attributes

int _nHits
 
std::vector< float > _aHit
 
std::vector< float > _tHit
 
std::vector< float > _xHit
 
std::vector< float > _yHit
 
std::vector< float > _zHit
 
int _ifNotGravity = 1
 
float _totAmpl = 0.0
 
float _totTime = 0.0
 
float _radius = 0.0
 
float _xgr = 0.0
 
float _ygr = 0.0
 
float _zgr = 0.0
 
float _analogGravity [3] = {0.0, 0.0, 0.0}
 
int _ifNotWidth = 1
 
float _analogWidth = 0.0
 
int _ifNotInertia = 1
 
float _ValAnalogInertia [3] = {0., 0., 0.}
 
float _VecAnalogInertia [9] = {0., 0., 0., 0., 0., 0., 0., 0., 0.}
 
int _ifNotElipsoid = 1
 
float _r1 = 0.0
 
float _r2 = 0.0
 
float _r3 = 0.0
 
float _vol = 0.0
 
float _r_ave = 0.0
 
float _density = 0.0
 
float _eccentricity = 0.0
 
float _r1_forw = 0.0
 
float _r1_back = 0.0
 

Detailed Description

Utility class to derive properties of clusters, such as centre of gravity, axes of inertia, fits of the cluster shape and so on. All the details are explained in the documentation of the methods. Several classes of the GSL (GNU Scientific Library) are needed in this class.

Authors
V. Morgunov (ITEP/DESY), A. Raspereza (DESY), O. Wendt (DESY)
Version
Id
ClusterShapes.h,v 1.14 2007-04-27 13:56:53 owendt Exp

Definition at line 25 of file ClusterShapes.h.

Constructor & Destructor Documentation

util::ClusterShapes::ClusterShapes ( int  nhits,
float *  a,
float *  x,
float *  y,
float *  z 
)

Constructor

Parameters
nhits: number of hits in the cluster
a: amplitudes of elements ('cells') of the cluster. Stored in an array, with one entry for each element ('cell'). Each entry is depending on coordinates x,y,z (Cartesian), which are stored in the arrays x,y,z.
x,y,z: array of coordinates corresponding to the array of amplitudes a.

Definition at line 13 of file ClusterShapes.cxx.

14  : _nHits(nhits), _aHit(nhits, 0.0), _tHit(nhits, 0.0), _xHit(nhits, 0.0), _yHit(nhits, 0.0), _zHit(nhits, 0.0), _ifNotGravity(1), _ifNotWidth(1), _ifNotInertia(1), _ifNotElipsoid(1)
15  {
16  for (int i(0); i < nhits; ++i)
17  {
18  _aHit[i] = a[i];
19  _tHit[i] = 0.;
20  _xHit[i] = x[i];
21  _yHit[i] = y[i];
22  _zHit[i] = z[i];
23  }
24  }
std::vector< float > _aHit
std::vector< float > _zHit
std::vector< float > _tHit
std::vector< float > _yHit
const double a
std::vector< float > _xHit
list x
Definition: train.py:276
util::ClusterShapes::ClusterShapes ( int  nhits,
float *  a,
float *  t,
float *  x,
float *  y,
float *  z 
)

Definition at line 26 of file ClusterShapes.cxx.

27  : _nHits(nhits), _aHit(nhits, 0.0), _tHit(nhits, 0.0), _xHit(nhits, 0.0), _yHit(nhits, 0.0), _zHit(nhits, 0.0), _ifNotGravity(1), _ifNotWidth(1), _ifNotInertia(1), _ifNotElipsoid(1)
28  {
29  for (int i(0); i < nhits; ++i)
30  {
31  _aHit[i] = a[i];
32  _tHit[i] = t[i];
33  _xHit[i] = x[i];
34  _yHit[i] = y[i];
35  _zHit[i] = z[i];
36  }
37  }
std::vector< float > _aHit
std::vector< float > _zHit
std::vector< float > _tHit
std::vector< float > _yHit
const double a
std::vector< float > _xHit
list x
Definition: train.py:276
util::ClusterShapes::ClusterShapes ( int  nhits,
std::vector< float >  a,
std::vector< float >  x,
std::vector< float >  y,
std::vector< float >  z 
)

Definition at line 41 of file ClusterShapes.cxx.

42  : _nHits(nhits), _aHit(nhits, 0.0), _xHit(nhits, 0.0), _yHit(nhits, 0.0), _zHit(nhits, 0.0), _ifNotGravity(1), _ifNotWidth(1), _ifNotInertia(1), _ifNotElipsoid(1)
43  {
44  for (int i(0); i < nhits; ++i)
45  {
46  _aHit[i] = a[i];
47  _tHit[i] = 0.;
48  _xHit[i] = x[i];
49  _yHit[i] = y[i];
50  _zHit[i] = z[i];
51  }
52  }
std::vector< float > _aHit
std::vector< float > _zHit
std::vector< float > _tHit
std::vector< float > _yHit
std::vector< float > _xHit
util::ClusterShapes::ClusterShapes ( int  nhits,
std::vector< float >  a,
std::vector< float >  t,
std::vector< float >  x,
std::vector< float >  y,
std::vector< float >  z 
)

Definition at line 54 of file ClusterShapes.cxx.

55  : _nHits(nhits), _aHit(nhits, 0.0), _tHit(nhits, 0.0), _xHit(nhits, 0.0), _yHit(nhits, 0.0), _zHit(nhits, 0.0), _ifNotGravity(1), _ifNotWidth(1), _ifNotInertia(1), _ifNotElipsoid(1)
56  {
57  for (int i(0); i < nhits; ++i)
58  {
59  _aHit[i] = a[i];
60  _tHit[i] = t[i];
61  _xHit[i] = x[i];
62  _yHit[i] = y[i];
63  _zHit[i] = z[i];
64  }
65  }
std::vector< float > _aHit
std::vector< float > _zHit
std::vector< float > _tHit
std::vector< float > _yHit
std::vector< float > _xHit
util::ClusterShapes::~ClusterShapes ( )

Destructor

Definition at line 70 of file ClusterShapes.cxx.

70  {
71 
72  }

Member Function Documentation

float util::ClusterShapes::findDistance ( int  i)
private

Definition at line 304 of file ClusterShapes.cxx.

304  {
305 
306  float cx = 0.0;
307  float cy = 0.0;
308  float cz = 0.0;
309  float dx = 0.0;
310  float dy = 0.0;
311  float dz = 0.0;
312 
313  cx = _VecAnalogInertia[0] ;
314  cy = _VecAnalogInertia[1] ;
315  cz = _VecAnalogInertia[2] ;
316  dx = _analogGravity[0] - _xHit[i] ;
317  dy = _analogGravity[1] - _yHit[i] ;
318  dz = _analogGravity[2] - _zHit[i] ;
319 
320  float tx = cy*dz - cz*dy ;
321  float ty = cz*dx - cx*dz ;
322  float tz = cx*dy - cy*dx ;
323  float tt = sqrt(tx*tx+ty*ty+tz*tz) ;
324  float ti = sqrt(cx*cx+cy*cy+cz*cz) ;
325  float f = tt / ti ;
326 
327  return f ;
328  }
std::vector< float > _zHit
std::vector< float > _yHit
Definition: type_traits.h:61
float _VecAnalogInertia[9]
std::vector< float > _xHit
void util::ClusterShapes::findElipsoid ( )
private

Elipsoid parameter calculations see cluster_proper.f

Definition at line 122 of file ClusterShapes.cxx.

122  {
123 
124  /** Elipsoid parameter calculations see cluster_proper.f */
125  float cx,cy,cz ;
126  float dx,dy,dz ;
127  float r_hit_max, d_begn, d_last, r_max, proj;
128 
129  if (_ifNotInertia == 1) findInertia() ;
130 
131  // Normalize the eigen values of inertia tensor
132  float wr1 = sqrt(_ValAnalogInertia[0]/_totAmpl);
133  float wr2 = sqrt(_ValAnalogInertia[1]/_totAmpl);
134  float wr3 = sqrt(_ValAnalogInertia[2]/_totAmpl);
135 
136  _r1 = sqrt(wr2*wr3); // spatial axis length -- the largest
137  _r2 = sqrt(wr1*wr3); // spatial axis length -- less
138  _r3 = sqrt(wr1*wr2); // spatial axis length -- even more less
139  _vol = 4.*M_PI*_r1*_r2*_r3/3.; // ellipsoid volume
140  _r_ave = pow(_vol,1/3); // average radius (quibc root)
141  _density = _totAmpl/_vol; // density
142  _eccentricity =_analogWidth/_r1; // Cluster Eccentricity
143 
144  // Find Minumal and Maximal Lenght for Principal axis
145  r_hit_max = -100000.;
146  d_begn = 100000.;
147  d_last = -100000.;
148  cx = _VecAnalogInertia[0] ;
149  cy = _VecAnalogInertia[1] ;
150  cz = _VecAnalogInertia[2] ;
151 
152  for (int i(0); i < _nHits; ++i)
153  {
154  dx = _xHit[i] - _xgr;
155  dy = _yHit[i] - _ygr;
156  dz = _zHit[i] - _zgr;
157  r_max = sqrt(dx*dx + dy*dy + dz*dz);
158 
159  if(r_max > r_hit_max) r_hit_max = r_max;
160  proj = dx*cx + dy*cy + dz*cz;
161 
162  if(proj < d_begn)
163  d_begn = proj;
164 
165  if(proj > d_last)
166  d_last = proj;
167  }
168 
169  _r1_forw = fabs(d_last);
170  _r1_back = fabs(d_begn);
171 
172  _ifNotElipsoid = 0;
173  }
constexpr T pow(T x)
Definition: pow.h:72
std::vector< float > _zHit
std::vector< float > _yHit
float _ValAnalogInertia[3]
float _VecAnalogInertia[9]
#define M_PI
Definition: includeROOT.h:54
std::vector< float > _xHit
void util::ClusterShapes::findGravity ( )
private

Definition at line 177 of file ClusterShapes.cxx.

177  {
178 
179  _totAmpl = 0. ;
180 
181  for (int i(0); i < 3; ++i)
182  _analogGravity[i] = 0.0 ;
183 
184  for (int i(0); i < _nHits; ++i)
185  {
186  _totAmpl+=_aHit[i] ;
187  _totTime+=_tHit[i] ;
188  _analogGravity[0]+=_aHit[i]*_xHit[i] ;
189  _analogGravity[1]+=_aHit[i]*_yHit[i] ;
190  _analogGravity[2]+=_aHit[i]*_zHit[i] ;
191  }
192 
193  for (int i(0); i < 3; ++i)
195 
196  _xgr = _analogGravity[0];
197  _ygr = _analogGravity[1];
198  _zgr = _analogGravity[2];
199 
200  _ifNotGravity = 0;
201  }
std::vector< float > _aHit
std::vector< float > _zHit
std::vector< float > _tHit
std::vector< float > _yHit
std::vector< float > _xHit
void util::ClusterShapes::findInertia ( )
private

Definition at line 205 of file ClusterShapes.cxx.

205  {
206 
207  double aIne[3][3];
208  float radius2 = 0.0;
209 
210  findGravity();
211 
212  for (int i(0); i < 3; ++i) {
213  for (int j(0); j < 3; ++j) {
214  aIne[i][j] = 0.0;
215  }
216  }
217 
218  for (int i(0); i < _nHits; ++i)
219  {
220  float dX = _xHit[i] - _analogGravity[0];
221  float dY = _yHit[i] - _analogGravity[1];
222  float dZ = _zHit[i] - _analogGravity[2];
223  aIne[0][0] += _aHit[i]*(dY*dY+dZ*dZ);
224  aIne[1][1] += _aHit[i]*(dX*dX+dZ*dZ);
225  aIne[2][2] += _aHit[i]*(dX*dX+dY*dY);
226  aIne[0][1] -= _aHit[i]*dX*dY;
227  aIne[0][2] -= _aHit[i]*dX*dZ;
228  aIne[1][2] -= _aHit[i]*dY*dZ;
229  }
230 
231  for (int i(0); i < 2; ++i) {
232  for (int j = i+1; j < 3; ++j) {
233  aIne[j][i] = aIne[i][j];
234  }
235  }
236  //****************************************
237  // analog Inertia
238  //****************************************
239 
240  gsl_matrix_view aMatrix = gsl_matrix_view_array((double*)aIne,3,3);
241  gsl_vector* aVector = gsl_vector_alloc(3);
242  gsl_matrix* aEigenVec = gsl_matrix_alloc(3,3);
243  gsl_eigen_symmv_workspace* wa = gsl_eigen_symmv_alloc(3);
244  gsl_eigen_symmv(&aMatrix.matrix,aVector,aEigenVec,wa);
245  gsl_eigen_symmv_free(wa);
246  gsl_eigen_symmv_sort(aVector,aEigenVec,GSL_EIGEN_SORT_ABS_ASC);
247 
248  for (int i(0); i < 3; i++) {
249  _ValAnalogInertia[i] = gsl_vector_get(aVector,i);
250  for (int j(0); j < 3; j++) {
251  _VecAnalogInertia[i+3*j] = gsl_matrix_get(aEigenVec,i,j);
252  }
253  }
254 
255  // Main principal points away from IP
256  _radius = 0.;
257  radius2 = 0.;
258 
259  for (int i(0); i < 3; ++i) {
261  radius2 += (_analogGravity[i]+_VecAnalogInertia[i])*(_analogGravity[i]+_VecAnalogInertia[i]);
262  }
263 
264  if ( radius2 < _radius) {
265  for (int i(0); i < 3; ++i)
267  }
268 
269  _radius = sqrt(_radius);
270  _ifNotInertia = 0;
271 
272  // The final job
273  findWidth();
274  findElipsoid();
275 
276  gsl_vector_free(aVector);
277  gsl_matrix_free(aEigenVec);
278 
279  }
std::vector< float > _aHit
std::vector< float > _zHit
std::vector< float > _yHit
float _ValAnalogInertia[3]
float _VecAnalogInertia[9]
std::vector< float > _xHit
void util::ClusterShapes::findWidth ( )
private

Definition at line 283 of file ClusterShapes.cxx.

283  {
284 
285  float dist = 0.0;
286 
287  if (_ifNotInertia == 1) findInertia() ;
288 
289  _analogWidth = 0.0 ;
290 
291  for (int i(0); i < _nHits; ++i)
292  {
293  dist = findDistance(i) ;
294  _analogWidth+=_aHit[i]*dist*dist ;
295  }
296 
298 
299  _ifNotWidth = 0 ;
300  }
std::vector< float > _aHit
float findDistance(int i)
constexpr double dist(const TReal *x, const TReal *y, const unsigned int dimension)
float util::ClusterShapes::getAverageTime ( )

returns the average time over the whole cluster

Definition at line 89 of file ClusterShapes.cxx.

89  {
90  if (_ifNotGravity == 1) findGravity();
91  return _totTime/_nHits;
92  }
float * util::ClusterShapes::getCenterOfGravity ( )

returns an array, which represents a vector from the origin of the coordiante system, i. e. IP, to the centre of gravity of the cluster. The centre of gravity is calculated with the energy of the entries of the cluster.

Definition at line 96 of file ClusterShapes.cxx.

96  {
97  if (_ifNotGravity == 1) findGravity() ;
98  return &_analogGravity[0] ;
99  }
float * util::ClusterShapes::getEigenValInertia ( )

array of the inertias of mass (i. e. energy) corresponding to the three main axes of inertia. The array is sorted in ascending order.

Definition at line 103 of file ClusterShapes.cxx.

103  {
104  if (_ifNotInertia == 1) findInertia();
105  return &_ValAnalogInertia[0] ;
106  }
float _ValAnalogInertia[3]
float * util::ClusterShapes::getEigenVecInertia ( )

array of the three main axes of inertia (9 entries) starting with the axis corresponding to the smallest inertia of mass (main principal axis). All axes are normalised to a length of 1.

Definition at line 110 of file ClusterShapes.cxx.

110  {
111  if (_ifNotInertia == 1) findInertia();
112  return &_VecAnalogInertia[0] ;
113  }
float _VecAnalogInertia[9]
float util::ClusterShapes::getElipsoid_density ( )
inline

density of the ellipsoid defined by: totAmpl/vol

Definition at line 139 of file ClusterShapes.h.

139 { return _density; }
float util::ClusterShapes::getElipsoid_eccentricity ( )
inline

eccentricity of the ellipsoid defined by: Width/r1

Definition at line 145 of file ClusterShapes.h.

145 { return _eccentricity; }
float util::ClusterShapes::getElipsoid_r1 ( )
inline

largest spatial axis length of the ellipsoid derived by the inertia tensor (by their eigenvalues and eigen- vectors)

Definition at line 110 of file ClusterShapes.h.

110 { return _r1; }
float util::ClusterShapes::getElipsoid_r2 ( )
inline

medium spatial axis length of the ellipsoid derived by the inertia tensor (by their eigenvalues and eigen- vectors)

Definition at line 117 of file ClusterShapes.h.

117 { return _r2; }
float util::ClusterShapes::getElipsoid_r3 ( )
inline

smallest spatial axis length of the ellipsoid derived by the inertia tensor (by their eigenvalues and eigen- vectors)

Definition at line 124 of file ClusterShapes.h.

124 { return _r3; }
float util::ClusterShapes::getElipsoid_r_ave ( )
inline

average radius of the ellipsoid (qubic root of volume)

Definition at line 134 of file ClusterShapes.h.

134 { return _r_ave; }
float util::ClusterShapes::getElipsoid_r_back ( )
inline

distance from centre of gravity to the point nearest to IP projected on the main principal axis

Definition at line 157 of file ClusterShapes.h.

157 { return _r1_back; }
float util::ClusterShapes::getElipsoid_r_forw ( )
inline

distance from centre of gravity to the point most far away from IP projected on the main principal axis

Definition at line 151 of file ClusterShapes.h.

151 { return _r1_forw; }
float util::ClusterShapes::getElipsoid_vol ( )
inline

volume of the ellipsoid

Definition at line 129 of file ClusterShapes.h.

129 { return _vol; }
int util::ClusterShapes::getNumberOfHits ( )

returns the number of elements of the cluster

Definition at line 76 of file ClusterShapes.cxx.

76  {
77  return _nHits;
78  }
float util::ClusterShapes::getTotalAmplitude ( )

returns the accumulated amplitude for the whole cluster (just the sum of the energy in all the entries of the cluster)

Definition at line 82 of file ClusterShapes.cxx.

82  {
83  if (_ifNotGravity == 1) findGravity();
84  return _totAmpl;
85  }
float util::ClusterShapes::getWidth ( )

'mean' width of the cluster perpendicular to the main principal axis, defined as: width := sqrt( 1/TotalAmplitude * Sum(a[i]*d[i]*d[i]) ), where d[i] is the distance of the i-th point to the main principal axis.

Definition at line 117 of file ClusterShapes.cxx.

117  {
118  if (_ifNotWidth == 1) findWidth();
119  return _analogWidth;
120  }
float util::ClusterShapes::radius ( )
inline

distance to the centre of gravity measured from IP (absolut value of the vector to the centre of gravity)

Definition at line 103 of file ClusterShapes.h.

103 { return _radius; }

Member Data Documentation

std::vector<float> util::ClusterShapes::_aHit
private

Definition at line 163 of file ClusterShapes.h.

float util::ClusterShapes::_analogGravity[3] = {0.0, 0.0, 0.0}
private

Definition at line 176 of file ClusterShapes.h.

float util::ClusterShapes::_analogWidth = 0.0
private

Definition at line 179 of file ClusterShapes.h.

float util::ClusterShapes::_density = 0.0
private

Definition at line 191 of file ClusterShapes.h.

float util::ClusterShapes::_eccentricity = 0.0
private

Definition at line 192 of file ClusterShapes.h.

int util::ClusterShapes::_ifNotElipsoid = 1
private

Definition at line 185 of file ClusterShapes.h.

int util::ClusterShapes::_ifNotGravity = 1
private

Definition at line 169 of file ClusterShapes.h.

int util::ClusterShapes::_ifNotInertia = 1
private

Definition at line 181 of file ClusterShapes.h.

int util::ClusterShapes::_ifNotWidth = 1
private

Definition at line 178 of file ClusterShapes.h.

int util::ClusterShapes::_nHits
private

Definition at line 161 of file ClusterShapes.h.

float util::ClusterShapes::_r1 = 0.0
private

Definition at line 186 of file ClusterShapes.h.

float util::ClusterShapes::_r1_back = 0.0
private

Definition at line 194 of file ClusterShapes.h.

float util::ClusterShapes::_r1_forw = 0.0
private

Definition at line 193 of file ClusterShapes.h.

float util::ClusterShapes::_r2 = 0.0
private

Definition at line 187 of file ClusterShapes.h.

float util::ClusterShapes::_r3 = 0.0
private

Definition at line 188 of file ClusterShapes.h.

float util::ClusterShapes::_r_ave = 0.0
private

Definition at line 190 of file ClusterShapes.h.

float util::ClusterShapes::_radius = 0.0
private

Definition at line 172 of file ClusterShapes.h.

std::vector<float> util::ClusterShapes::_tHit
private

Definition at line 164 of file ClusterShapes.h.

float util::ClusterShapes::_totAmpl = 0.0
private

Definition at line 170 of file ClusterShapes.h.

float util::ClusterShapes::_totTime = 0.0
private

Definition at line 171 of file ClusterShapes.h.

float util::ClusterShapes::_ValAnalogInertia[3] = {0., 0., 0.}
private

Definition at line 182 of file ClusterShapes.h.

float util::ClusterShapes::_VecAnalogInertia[9] = {0., 0., 0., 0., 0., 0., 0., 0., 0.}
private

Definition at line 183 of file ClusterShapes.h.

float util::ClusterShapes::_vol = 0.0
private

Definition at line 189 of file ClusterShapes.h.

float util::ClusterShapes::_xgr = 0.0
private

Definition at line 173 of file ClusterShapes.h.

std::vector<float> util::ClusterShapes::_xHit
private

Definition at line 165 of file ClusterShapes.h.

float util::ClusterShapes::_ygr = 0.0
private

Definition at line 174 of file ClusterShapes.h.

std::vector<float> util::ClusterShapes::_yHit
private

Definition at line 166 of file ClusterShapes.h.

float util::ClusterShapes::_zgr = 0.0
private

Definition at line 175 of file ClusterShapes.h.

std::vector<float> util::ClusterShapes::_zHit
private

Definition at line 167 of file ClusterShapes.h.


The documentation for this class was generated from the following files: