BmnRoot
Loading...
Searching...
No Matches
BmnNdetCluster.cxx
Go to the documentation of this file.
1
8#include "BmnNdetCluster.h"
9
10#include "TMath.h"
11
13 : TObject()
14{
15 fNCells = 0;
16}
17
19
21 : TObject()
22{
23 fCell.reserve(size);
24 fColumn.reserve(size);
25 fRow.reserve(size);
26 fLayer.reserve(size);
27 fXcoord.reserve(size);
28 fYcoord.reserve(size);
29 fZcoord.reserve(size);
30 fTime.reserve(size);
31 fEdep.reserve(size);
32 fBeta.reserve(size);
33 fNCells = 0;
34}
35
37 Int_t column,
38 Int_t row,
39 Int_t layer,
40 Double_t x,
41 Double_t y,
42 Double_t z,
43 Double_t t,
44 Double_t e_dep,
45 Double_t beta)
46{
47 fCell.push_back(cell);
48 fColumn.push_back(column);
49 fRow.push_back(row);
50 fLayer.push_back(layer);
51 fXcoord.push_back(x);
52 fYcoord.push_back(y);
53 fZcoord.push_back(z);
54 fTime.push_back(t);
55 fEdep.push_back(e_dep);
56
57 if (beta < 0.) {
58 Double_t distance = TMath::Sqrt(TMath::Sq(x) + TMath::Sq(y) + TMath::Sq(z));
59 Double_t _beta = t >= 0 ? distance * 1.E-2 / (t * 1E-9) / TMath::C() : -1.;
60 fBeta.push_back(_beta);
61 } else {
62 fBeta.push_back(beta);
63 }
64 ++fNCells;
65}
66
68{
69 return fNCells;
70}
71
72Int_t BmnNdetCluster::GetNCells(Double_t threshold)
73{
74 Int_t nCells = 0.;
75 for (Int_t iCell = 0; iCell < fNCells; ++iCell) {
76 if (fEdep[iCell] > threshold)
77 ++nCells;
78 }
79 return nCells;
80}
81
83{
84 return fCell;
85}
87{
88 return fColumn;
89}
91{
92 return fRow;
93}
95{
96 return fLayer;
97}
98vector<Double_t> BmnNdetCluster::GetX()
99{
100 return fXcoord;
101}
102vector<Double_t> BmnNdetCluster::GetY()
103{
104 return fYcoord;
105}
106vector<Double_t> BmnNdetCluster::GetZ()
107{
108 return fZcoord;
109}
110vector<Double_t> BmnNdetCluster::GetTime()
111{
112 return fTime;
113}
114vector<Double_t> BmnNdetCluster::GetEdep()
115{
116 return fEdep;
117}
118vector<Double_t> BmnNdetCluster::GetBeta()
119{
120 return fBeta;
121}
122
123Bool_t BmnNdetCluster::ContainsVetoCells(Int_t veto_layer)
124{
125 for (Int_t iCell = 0; iCell < fNCells; ++iCell) {
126 if (fLayer[iCell] == veto_layer)
127 return kTRUE;
128 }
129 return kFALSE;
130}
131
132Bool_t BmnNdetCluster::StartsOnBorder(Int_t size_x, Int_t size_y, Int_t size_z)
133{
134 Int_t iCell = FirstCell();
135 if (fLayer[iCell] == 0)
136 return kTRUE;
137 if (fLayer[iCell] == size_z - 1)
138 return kTRUE;
139 if (fColumn[iCell] == 0)
140 return kTRUE;
141 if (fColumn[iCell] == size_x - 1)
142 return kTRUE;
143 if (fRow[iCell] == 0)
144 return kTRUE;
145 if (fRow[iCell] == size_y - 1)
146 return kTRUE;
147
148 return kFALSE;
149}
150
151Double_t BmnNdetCluster::EnergyDepositedAverage(Double_t threshold)
152{
153 Double_t edep = 0.;
154 Double_t nCells = 0.;
155 for (Int_t iCell = 0; iCell < fNCells; ++iCell) {
156 if (fEdep[iCell] > threshold) {
157 edep += fEdep[iCell];
158 nCells += 1.;
159 }
160 }
161 if (nCells > 0.)
162 return edep / nCells;
163 return 0.;
164}
165
166Double_t BmnNdetCluster::EnergyDeposited(Double_t threshold)
167{
168 Double_t edep = 0.;
169 for (Int_t iCell = 0; iCell < fNCells; ++iCell) {
170 if (fEdep[iCell] > threshold) {
171 edep += fEdep[iCell];
172 }
173 }
174 return edep;
175}
176
178{
179 Double_t t = 100.;
180 Int_t iFirst = 0;
181 for (Int_t iCell = 0; iCell < fNCells; ++iCell) {
182 if (fTime[iCell] < t) {
183 t = fTime[iCell];
184 iFirst = iCell;
185 }
186 }
187 return iFirst;
188}
189
191{
192 Double_t beta = 0.;
193 Int_t iFastest = 0;
194 for (Int_t iCell = 0; iCell < fNCells; ++iCell) {
195 if (fBeta[iCell] > beta) {
196 beta = fBeta[iCell];
197 iFastest = iCell;
198 }
199 }
200 return iFastest;
201}
202
204{
205 Int_t iFirst = FirstCell();
206 return fBeta[iFirst];
207}
208
210{
211 Int_t iFastest = FastestCell();
212 return fBeta[iFastest];
213}
214
215Double_t BmnNdetCluster::EnergyTofFirst(Double_t mass)
216{
217 Double_t beta = BetaFirst();
218 if (beta >= 1.)
219 return -1.;
220 Double_t gamma = 1. / TMath::Sqrt(1. - TMath::Sq(beta));
221 return (gamma - 1.) * mass;
222}
223
224Double_t BmnNdetCluster::EnergyTofFastest(Double_t mass)
225{
226 Double_t beta = BetaFastest();
227 if (beta >= 1.)
228 return -1.;
229 Double_t gamma = 1. / TMath::Sqrt(1. - TMath::Sq(beta));
230 return (gamma - 1.) * mass;
231}
232
233Double_t BmnNdetCluster::EnergyTofAverage(Double_t mass)
234{
235 Double_t t = 0.;
236 Double_t x = 0.;
237 Double_t y = 0.;
238 Double_t z = 0.;
239 for (Int_t iCell = 0; iCell < fNCells; ++iCell) {
240 t += fTime[iCell];
241 x += fXcoord[iCell];
242 y += fYcoord[iCell];
243 z += fZcoord[iCell];
244 }
245 Double_t distance = TMath::Sqrt(TMath::Sq(x) + TMath::Sq(y) + TMath::Sq(z));
246 if (t <= 0.)
247 return -1.;
248 Double_t beta = (distance * 1.E-2 / (t * 1.E-9)) / TMath::C();
249 if (beta >= 1.)
250 return -1.;
251 Double_t gamma = 1. / TMath::Sqrt(1. - TMath::Sq(beta));
252 return (gamma - 1.) * mass;
253}
254
256{
257 Int_t iFirst = FirstCell();
258 Double_t xTarget = fXcoord[iFirst];
259 Double_t yTarget = fYcoord[iFirst];
260 Double_t zTarget = fZcoord[iFirst];
261 Double_t xCluster = 0;
262 Double_t yCluster = 0;
263 Double_t zCluster = 0;
264 for (Int_t iCell = 0; iCell < fNCells; ++iCell) {
265 Double_t x = fXcoord[iCell] - fXcoord[iFirst];
266 Double_t y = fYcoord[iCell] - fYcoord[iFirst];
267 Double_t z = fZcoord[iCell] - fZcoord[iFirst];
268
269 x *= fEdep[iCell];
270 y *= fEdep[iCell];
271 z *= fEdep[iCell];
272
273 xCluster += x;
274 yCluster += y;
275 zCluster += z;
276 }
277 Double_t scalarProduct = xCluster * xTarget + yCluster * yTarget + zCluster * zTarget;
278 Double_t clusterLength = TMath::Sqrt(TMath::Sq(xCluster) + TMath::Sq(yCluster) + TMath::Sq(zCluster));
279 Double_t targetLength = TMath::Sqrt(TMath::Sq(xTarget) + TMath::Sq(yTarget) + TMath::Sq(zTarget));
280 scalarProduct /= (clusterLength * targetLength);
281 return TMath::ACos(scalarProduct) * TMath::RadToDeg();
282}
283
284Double_t BmnNdetCluster::AngleToDirection(Double_t x, Double_t y, Double_t z)
285{
286 Int_t iFirst = FirstCell();
287 Double_t xTarget = x;
288 Double_t yTarget = y;
289 Double_t zTarget = z;
290 Double_t xCluster = 0;
291 Double_t yCluster = 0;
292 Double_t zCluster = 0;
293 for (Int_t iCell = 0; iCell < fNCells; ++iCell) {
294 Double_t dx = fXcoord[iCell] - fXcoord[iFirst];
295 Double_t dy = fYcoord[iCell] - fYcoord[iFirst];
296 Double_t dz = fZcoord[iCell] - fZcoord[iFirst];
297
298 dx *= fEdep[iCell];
299 dy *= fEdep[iCell];
300 dz *= fEdep[iCell];
301
302 xCluster += dx;
303 yCluster += dy;
304 zCluster += dz;
305 }
306 Double_t scalarProduct = xCluster * xTarget + yCluster * yTarget + zCluster * zTarget;
307 Double_t clusterLength = TMath::Sqrt(TMath::Sq(xCluster) + TMath::Sq(yCluster) + TMath::Sq(zCluster));
308 Double_t targetLength = TMath::Sqrt(TMath::Sq(xTarget) + TMath::Sq(yTarget) + TMath::Sq(zTarget));
309 scalarProduct /= (clusterLength * targetLength);
310 return TMath::ACos(scalarProduct) * TMath::RadToDeg();
311}
312
314{
315
316 for (Int_t iCell = 0; iCell < other->GetNCells(); ++iCell) {
317 this->AddCell(other->GetCells()[iCell], other->GetColumns()[iCell], other->GetRows()[iCell],
318 other->GetLayers()[iCell], other->GetX()[iCell], other->GetY()[iCell], other->GetZ()[iCell],
319 other->GetTime()[iCell], other->GetEdep()[iCell], other->GetBeta()[iCell]);
320 }
321 return this;
322}
323
325{
326 fCell.clear();
327 fColumn.clear();
328 fRow.clear();
329 fLayer.clear();
330 fXcoord.clear();
331 fYcoord.clear();
332 fZcoord.clear();
333 fTime.clear();
334 fEdep.clear();
335 fBeta.clear();
336
337 fNCells = 0;
338}
339
341 : TObject(other)
342{}
343
345 : TObject(move(other))
346{}
347
348void BmnNdetCluster::Print(const Option_t* option)
349{
350 cout << "----------------------------------- BmnNdetCluster -----------------------------" << endl;
351 for (size_t i = 0; i < fEdep.size(); ++i) {
352 cout << fCell[i] << " cell (" << fColumn[i] << ", " << fRow[i] << ", " << fLayer[i] << "): ";
353 cout << "\tX: " << fXcoord[i] << "[cm], Y: " << fYcoord[i] << "[cm], Z: " << fZcoord[i] << "[cm], ";
354 cout << "\tt: " << fTime[i] << "[ns], Edep: " << 1.E3 * fEdep[i] << "[MeV], beta: " << fBeta[i] << endl;
355 }
356 cout << "------------------------------------------------------------------------------------" << endl;
357}
int i
Definition P4_F32vec4.h:22
Class for description of recognized cluster in Bmn Ndet detector.
Double_t EnergyTofFastest(Double_t mass=0.939565420)
vector< Double_t > GetZ()
Double_t AngleToDirection(Double_t x=0, Double_t y=0, Double_t z=1)
Double_t EnergyDeposited(Double_t threshold=0.003)
vector< Double_t > GetBeta()
Double_t BetaFastest()
void Clear()
Clear all cell contents, set fNCells=0.
virtual void Print(const Option_t *option="")
Print.
vector< Int_t > GetLayers()
void AddCell(Int_t cell, Int_t column, Int_t row, Int_t layer, Double_t x, Double_t y, Double_t z, Double_t t, Double_t e_dep, Double_t beta=-1.)
BmnNdetCluster * AddCluster(BmnNdetCluster *other)
Merging of clusters.
vector< Double_t > GetTime()
Double_t EnergyTofFirst(Double_t mass=0.939565420)
vector< Double_t > GetY()
vector< Double_t > GetEdep()
vector< Int_t > GetColumns()
BmnNdetCluster()
Default constructor.
Bool_t ContainsVetoCells(Int_t veto_layer=0)
vector< Int_t > GetCells()
Double_t EnergyTofAverage(Double_t mass=0.939565420)
~BmnNdetCluster()
Default destructor.
Double_t AngleToTarget()
vector< Double_t > GetX()
Double_t EnergyDepositedAverage(Double_t threshold=0.003)
Bool_t StartsOnBorder(Int_t size_x=11, Int_t size_y=11, Int_t size_z=8)
vector< Int_t > GetRows()