ClusterShapes.h
Go to the documentation of this file.
1 #ifndef ClusterShapes_h
2 #define ClusterShapes_h
3 
4 #include <stdio.h>
5 #include <iostream>
6 #include <iomanip>
7 #include <string>
8 #include <sstream>
9 #include <cstdlib>
10 #include <math.h>
11 #include <vector>
12 
13 namespace util {
14 
15  /**
16  * Utility class to derive properties of clusters, such as centre of gravity,
17  * axes of inertia, fits of the cluster shape and so on. All the details are
18  * explained in the documentation of the methods. Several classes of the GSL
19  * (GNU Scientific Library) are needed in this class.
20  *
21  * @authors V. Morgunov (ITEP/DESY), A. Raspereza (DESY), O. Wendt (DESY)
22  * @version $Id: ClusterShapes.h,v 1.14 2007-04-27 13:56:53 owendt Exp $
23  *
24  */
25  class ClusterShapes {
26 
27  public:
28 
29  /**
30  * Constructor
31  * @param nhits : number of hits in the cluster
32  * @param a : amplitudes of elements ('cells') of the cluster. Stored in
33  * an array, with one entry for each element ('cell'). Each entry
34  * is depending on coordinates x,y,z (Cartesian), which are stored
35  * in the arrays x,y,z.
36  * @param x,y,z : array of coordinates corresponding to the array of amplitudes a.
37  *
38  *
39  */
40  ClusterShapes(int nhits, float* a, float* x, float* y, float* z);
41 
42  ClusterShapes(int nhits, float* a, float* t, float* x, float* y, float* z);
43 
44  ClusterShapes(int nhits, std::vector<float> a, std::vector<float> x, std::vector<float> y, std::vector<float> z);
45 
46  ClusterShapes(int nhits, std::vector<float> a, std::vector<float> t, std::vector<float> x, std::vector<float> y, std::vector<float> z);
47 
48  /**
49  * Destructor
50  */
52 
53  /**
54  * returns the number of elements of the cluster
55  */
56  int getNumberOfHits();
57 
58  /**
59  * returns the accumulated amplitude for the whole cluster (just the sum of the
60  * energy in all the entries of the cluster)
61  */
62  float getTotalAmplitude();
63 
64  /**
65  * returns the average time over the whole cluster
66  */
67  float getAverageTime();
68 
69  /**
70  * returns an array, which represents a vector from the origin of the
71  * coordiante system, i.\ e.\ IP, to the centre of gravity of the cluster. The centre
72  * of gravity is calculated with the energy of the entries of the cluster.
73  */
74  float* getCenterOfGravity();
75 
76  /**
77  * array of the inertias of mass (i.\ e.\ energy) corresponding to the three main axes
78  * of inertia. The array is sorted in ascending order.
79  */
80  float* getEigenValInertia();
81 
82  /**
83  * array of the three main axes of inertia (9 entries) starting
84  * with the axis corresponding to the smallest inertia of mass
85  * (main principal axis). All axes are normalised to a length
86  * of 1.
87  */
88  float* getEigenVecInertia();
89 
90  /**
91  * 'mean' width of the cluster perpendicular to the main
92  * principal axis, defined as:
93  * width := sqrt( 1/TotalAmplitude * Sum(a[i]*d[i]*d[i]) ),
94  * where d[i] is the distance of the i-th point to the main
95  * principal axis.
96  */
97  float getWidth();
98 
99  /**
100  * distance to the centre of gravity measured from IP
101  * (absolut value of the vector to the centre of gravity)
102  */
103  inline float radius() { return _radius; }
104 
105  /**
106  * largest spatial axis length of the ellipsoid derived
107  * by the inertia tensor (by their eigenvalues and eigen-
108  * vectors)
109  */
110  inline float getElipsoid_r1() { return _r1; }
111 
112  /**
113  * medium spatial axis length of the ellipsoid derived
114  * by the inertia tensor (by their eigenvalues and eigen-
115  * vectors)
116  */
117  inline float getElipsoid_r2() { return _r2; }
118 
119  /**
120  * smallest spatial axis length of the ellipsoid derived
121  * by the inertia tensor (by their eigenvalues and eigen-
122  * vectors)
123  */
124  inline float getElipsoid_r3() { return _r3; }
125 
126  /**
127  * volume of the ellipsoid
128  */
129  inline float getElipsoid_vol() { return _vol; }
130 
131  /**
132  * average radius of the ellipsoid (qubic root of volume)
133  */
134  inline float getElipsoid_r_ave() { return _r_ave; }
135 
136  /**
137  * density of the ellipsoid defined by: totAmpl/vol
138  */
139  inline float getElipsoid_density() { return _density; }
140 
141  /**
142  * eccentricity of the ellipsoid defined by:
143  * Width/r1
144  */
145  inline float getElipsoid_eccentricity() { return _eccentricity; }
146 
147  /**
148  * distance from centre of gravity to the point most far
149  * away from IP projected on the main principal axis
150  */
151  inline float getElipsoid_r_forw() { return _r1_forw; }
152 
153  /**
154  * distance from centre of gravity to the point nearest
155  * to IP projected on the main principal axis
156  */
157  inline float getElipsoid_r_back() { return _r1_back; }
158 
159  private:
160 
161  int _nHits;
162 
163  std::vector<float> _aHit;
164  std::vector<float> _tHit;
165  std::vector<float> _xHit;
166  std::vector<float> _yHit;
167  std::vector<float> _zHit;
168 
169  int _ifNotGravity = 1;
170  float _totAmpl = 0.0;
171  float _totTime = 0.0;
172  float _radius = 0.0;
173  float _xgr = 0.0;
174  float _ygr = 0.0;
175  float _zgr = 0.0;
176  float _analogGravity[3] = {0.0, 0.0, 0.0};
177 
178  int _ifNotWidth = 1;
179  float _analogWidth = 0.0;
180 
181  int _ifNotInertia = 1;
182  float _ValAnalogInertia[3] = {0., 0., 0.};
183  float _VecAnalogInertia[9] = {0., 0., 0., 0., 0., 0., 0., 0., 0.};
184 
185  int _ifNotElipsoid = 1;
186  float _r1 = 0.0; // Cluster spatial axis length -- the largest
187  float _r2 = 0.0; // Cluster spatial axis length -- less
188  float _r3 = 0.0; // Cluster spatial axis length -- less
189  float _vol = 0.0; // Cluster ellipsoid volume
190  float _r_ave = 0.0; // Cluster average radius (qubic root)
191  float _density = 0.0; // Cluster density
192  float _eccentricity = 0.0; // Cluster Eccentricity
193  float _r1_forw = 0.0;
194  float _r1_back = 0.0;
195 
196  void findElipsoid();
197  void findGravity();
198  void findInertia();
199  void findWidth();
200  float findDistance(int i);
201 
202  };
203 }
204 
205 #endif
Namespace for general, non-LArSoft-specific utilities.
std::vector< float > _aHit
std::vector< float > _zHit
std::vector< float > _tHit
std::vector< float > _yHit
float * getEigenValInertia()
float _ValAnalogInertia[3]
float findDistance(int i)
const double a
float _VecAnalogInertia[9]
ClusterShapes(int nhits, float *a, float *x, float *y, float *z)
float getElipsoid_eccentricity()
float * getEigenVecInertia()
std::vector< float > _xHit
list x
Definition: train.py:276
float * getCenterOfGravity()