REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestVolumeHits.cxx
1
18
19#include "TRestVolumeHits.h"
20
21using namespace std;
22
23ClassImp(TRestVolumeHits);
24
25TRestVolumeHits::TRestVolumeHits() {
26 // TRestVolumeHits default constructor
27}
28
29TRestVolumeHits::~TRestVolumeHits() {
30 // TRestVolumeHits destructor
31}
32
33void TRestVolumeHits::AddHit(Double_t x, Double_t y, Double_t z, Double_t en, Double_t time,
34 REST_HitType type, Double_t sigmaX, Double_t sigmaY, Double_t sigmaZ) {
35 if (fType.size() > 0 && type != fType[0]) {
36 cout << "Error! Cannot add different typed hits into TRestVolumeHits!" << endl;
37 return;
38 }
39
40 TRestHits::AddHit({x, y, z}, en, time, type);
41 fSigmaX.push_back((Float_t)sigmaX);
42 fSigmaY.push_back((Float_t)sigmaY);
43 fSigmaZ.push_back((Float_t)sigmaZ);
44}
45
46void TRestVolumeHits::AddHit(const TVector3& pos, Double_t en, Double_t time, REST_HitType type,
47 const TVector3& sigma) {
48 if (fType.size() > 0 && type != fType[0]) {
49 cout << "Error! Cannot add different typed hits into TRestVolumeHits!" << endl;
50 return;
51 }
52
53 TRestHits::AddHit(pos, en, time, type);
54 fSigmaX.push_back((Float_t)sigma.X());
55 fSigmaY.push_back((Float_t)sigma.Y());
56 fSigmaZ.push_back((Float_t)sigma.Z());
57}
58
59void TRestVolumeHits::AddHit(const TRestVolumeHits& hits, Int_t n) {
60 Double_t x = hits.GetX(n);
61 Double_t y = hits.GetY(n);
62 Double_t z = hits.GetZ(n);
63 Double_t t = hits.GetTime(n);
64 Double_t en = hits.GetEnergy(n);
65 REST_HitType type = hits.GetType(n);
66
67 Double_t sx = hits.GetSigmaX(n);
68 Double_t sy = hits.GetSigmaY(n);
69 Double_t sz = hits.GetSigmaZ(n);
70
71 if (fType.size() > 0 && type != fType[0]) {
72 cout << "Error! Cannot add different typed hits into TRestVolumeHits!" << endl;
73 return;
74 }
75
76 AddHit(x, y, z, en, t, type, sx, sy, sz);
77}
78
80 // cout << "Removing hits" << endl;
82 fSigmaX.clear();
83 fSigmaY.clear();
84 fSigmaZ.clear();
85}
86
87Bool_t TRestVolumeHits::areXY() const {
88 if (fType.size() == 0)
89 return TMath::IsNaN(fZ[0]) && !TMath::IsNaN(fX[0]) && !TMath::IsNaN(fY[0]);
90 else
91 return fType[0] == XY;
92}
93Bool_t TRestVolumeHits::areXZ() const {
94 if (fType.size() == 0)
95 return !TMath::IsNaN(fZ[0]) && !TMath::IsNaN(fX[0]) && TMath::IsNaN(fY[0]);
96 else
97 return fType[0] == XZ;
98}
99Bool_t TRestVolumeHits::areYZ() const {
100 if (fType.size() == 0)
101 return !TMath::IsNaN(fZ[0]) && TMath::IsNaN(fX[0]) && !TMath::IsNaN(fY[0]);
102 else
103 return fType[0] == YZ;
104}
106 if (fType.size() == 0)
107 return !TMath::IsNaN(fZ[0]) && !TMath::IsNaN(fX[0]) && !TMath::IsNaN(fY[0]);
108 else
109 return fType[0] == XYZ;
110}
111
112void TRestVolumeHits::MergeHits(Int_t n, Int_t m) {
113 Double_t totalEnergy = fEnergy[n] + fEnergy[m];
114
115 // TODO : This is wrong but not very important for the moment
116 fSigmaX[n] = (fSigmaX[n] * fEnergy[n] + fSigmaX[m] * fEnergy[m]) / totalEnergy;
117 fSigmaY[n] = (fSigmaY[n] * fEnergy[n] + fSigmaY[m] * fEnergy[m]) / totalEnergy;
118 fSigmaZ[n] = (fSigmaZ[n] * fEnergy[n] + fSigmaZ[m] * fEnergy[m]) / totalEnergy;
119
120 fSigmaX.erase(fSigmaX.begin() + m);
121 fSigmaY.erase(fSigmaY.begin() + m);
122 fSigmaZ.erase(fSigmaZ.begin() + m);
123
125}
126
129
130 fSigmaX.erase(fSigmaX.begin() + n);
131 fSigmaY.erase(fSigmaY.begin() + n);
132 fSigmaZ.erase(fSigmaZ.begin() + n);
133}
134
135void TRestVolumeHits::SortByEnergy() {
136 while (!isSortedByEnergy()) {
137 for (unsigned int i = 0; i < GetNumberOfHits(); i++) {
138 for (unsigned int j = i + 1; j < GetNumberOfHits(); j++) {
139 if (GetEnergy(i) < GetEnergy(j)) SwapHits(i, j);
140 }
141 }
142 }
143}
144
145void TRestVolumeHits::SwapHits(Int_t i, Int_t j) {
146 iter_swap(fSigmaX.begin() + i, fSigmaX.begin() + j);
147 iter_swap(fSigmaY.begin() + i, fSigmaY.begin() + j);
148 iter_swap(fSigmaZ.begin() + i, fSigmaZ.begin() + j);
149
151}
152
153TVector3 TRestVolumeHits::GetSigma(int n) const {
154 return TVector3(((Double_t)fSigmaX[n]), ((Double_t)fSigmaY[n]), ((Double_t)fSigmaZ[n]));
155}
156
157void TRestVolumeHits::PrintHits() const {
158 for (unsigned int n = 0; n < GetNumberOfHits(); n++) {
159 cout << "Hit " << n << " X: " << GetX(n) << " Y: " << GetY(n) << " Z: " << GetZ(n)
160 << " sX: " << GetSigmaX(n) << " sY: " << GetSigmaY(n) << " sZ: " << GetSigmaZ(n)
161 << " Energy: " << GetEnergy(n) << endl;
162 }
163}
164
165void TRestVolumeHits::kMeansClustering(TRestVolumeHits* hits, TRestVolumeHits& vHits, int maxIt) {
166 const int nodes = vHits.GetNumberOfHits();
167 vector<TRestVolumeHits> volHits(nodes);
168 // std::cout<<"Nhits "<<hits->GetNumberOfHits()<<" Nodes "<<nodes<<std::endl;
169 TVector3 nullVector = TVector3(0, 0, 0);
170 std::vector<TVector3> centroid(nodes);
171 std::vector<TVector3> centroidOld(nodes, nullVector);
172
173 for (int h = 0; h < nodes; h++) centroid[h] = vHits.GetPosition(h);
174
175 for (int it = 0; it < maxIt; it++) {
176 for (int n = 0; n < nodes; n++) volHits[n].RemoveHits();
177 for (unsigned int i = 0; i < hits->GetNumberOfHits(); i++) {
178 double minDist = 1E9;
179 int clIndex = -1;
180 for (int n = 0; n < nodes; n++) {
181 TVector3 hitPos = hits->GetPosition(i);
182 double dist = (centroid[n] - hitPos).Mag();
183 if (dist < minDist) {
184 minDist = dist;
185 clIndex = n;
186 }
187 }
188 // cout<<minDist<<" "<<clIndex<<endl;
189 volHits[clIndex].AddHit(*hits, i);
190 }
191 bool converge = true;
192 for (int n = 0; n < nodes; n++) {
193 centroid[n] = volHits[n].GetMeanPosition();
194 converge &= (centroid[n] == centroidOld[n]);
195 centroidOld[n] = centroid[n];
196 }
197 if (converge) {
198 break;
199 }
200 }
201
202 vHits.RemoveHits();
203 const TVector3 sigma(0., 0., 0.);
204 for (int n = 0; n < nodes; n++) {
205 if (volHits[n].GetNumberOfHits() > 0)
206 vHits.AddHit(volHits[n].GetMeanPosition(), volHits[n].GetTotalEnergy(), 0, volHits[n].GetType(0),
207 sigma);
208 }
209}
TVector3 GetMeanPosition() const
It calculates the mean position weighting with the energy of the hits. Each coordinate is calculated ...
Definition: TRestHits.cxx:658
virtual void RemoveHits()
It removes all hits inside the class.
Definition: TRestHits.cxx:371
virtual void MergeHits(int n, int m)
It merges hits n and m being the resulting hit placed at the weighted center and being its final ener...
Definition: TRestHits.cxx:458
Bool_t isSortedByEnergy() const
It returns true if the hits are ordered in increasing energies.
Definition: TRestHits.cxx:490
Double_t GetSigmaX() const
It calculates the hits standard deviation in the X-coordinate.
Definition: TRestHits.cxx:685
virtual void RemoveHit(int n)
It removes the hit at position n from the list.
Definition: TRestHits.cxx:503
std::vector< REST_HitType > fType
The type of hit X,Y,XY,XYZ, ...
Definition: TRestHits.h:62
std::vector< Float_t > fZ
Position on Z axis for each punctual deposition (units mm)
Definition: TRestHits.h:53
std::vector< Float_t > fY
Position on Y axis for each punctual deposition (units mm)
Definition: TRestHits.h:50
std::vector< Float_t > fX
Position on X axis for each punctual deposition (units mm)
Definition: TRestHits.h:47
std::vector< Float_t > fEnergy
Energy deposited at each 3-coordinate position (units keV)
Definition: TRestHits.h:59
void AddHit(const TVector3 &position, Double_t energy, Double_t time=0, REST_HitType type=XYZ)
Adds a new hit to the list of hits using a TVector3.
Definition: TRestHits.cxx:345
Double_t GetSigmaY() const
It calculates the hits standard deviation in the Y-coordinate.
Definition: TRestHits.cxx:703
virtual void SwapHits(Int_t i, Int_t j)
It exchanges hits n and m affecting to the ordering of the hits inside the list of hits.
Definition: TRestHits.cxx:478
TVector3 GetPosition(int n) const
It returns the position of hit number n.
Definition: TRestHits.cxx:515
void RemoveHit(int n)
It removes the hit at position n from the list.
Bool_t areXZ() const
It will return true only if all the hits inside are of type XZ.
Bool_t areXYZ() const
It will return true only if all the hits inside are of type XYZ.
void SwapHits(Int_t i, Int_t j)
It exchanges hits n and m affecting to the ordering of the hits inside the list of hits.
Bool_t areYZ() const
It will return true only if all the hits inside are of type YZ.
Bool_t areXY() const
It will return true only if all the hits inside are of type XY.
void RemoveHits()
It removes all hits inside the class.