SND@LHC Software
Loading...
Searching...
No Matches
snd::analysis_tools::ScifiPlane Class Reference

#include <sndScifiPlane.h>

Collaboration diagram for snd::analysis_tools::ScifiPlane:

Classes

struct  ScifiHit
 
struct  xy_pair
 

Public Member Functions

 ScifiPlane (TClonesArray *snd_hits, Configuration configuration, Scifi *scifi_geometry, int index_begin, int index_end, int station)
 
const int GetStation () const
 
const std::vector< ScifiHitGetHits () const
 
const TVector3 GetCentroid () const
 
const TVector3 GetCentroidError () const
 
const xy_pair< double > GetTotQdc (bool only_positive=false) const
 
const xy_pair< double > GetTotEnergy (bool only_positive=false) const
 
const xy_pair< int > GetNHits () const
 
const TVector3 GetCluster (int max_gap) const
 
void FindCentroid ()
 
bool IsShower () const
 
void TimeFilter (double min_timestamp, double max_timestamp)
 
xy_pair< double > GetPointQdc (const TVector3 &point, double radius) const
 

Private Attributes

std::vector< ScifiHithits_
 
Configuration configuration_
 
TVector3 centroid_
 
TVector3 centroid_error_
 
int station_
 

Detailed Description

Definition at line 13 of file sndScifiPlane.h.

Constructor & Destructor Documentation

◆ ScifiPlane()

snd::analysis_tools::ScifiPlane::ScifiPlane ( TClonesArray *  snd_hits,
snd::Configuration  configuration,
Scifi scifi_geometry,
int  index_begin,
int  index_end,
int  station 
)

Definition at line 14 of file sndScifiPlane.cxx.

14 : configuration_(configuration), centroid_(std::nan(""), std::nan(""), std::nan("")), centroid_error_(std::nan(""), std::nan(""), std::nan("")), station_(station)
15{
16 if (index_begin > index_end)
17 {
18 throw std::runtime_error{"Begin index > end index"};
19 }
20
21 for (int i{index_begin}; i < index_end; ++i)
22 {
23 auto snd_hit = static_cast<sndScifiHit *>(snd_hits->At(i));
24 ScifiHit hit;
25 hit.channel_index = 512 * snd_hit->GetMat() + 64 * snd_hit->GetTofpetID(0) + 63 - snd_hit->Getchannel(0);
26 hit.timestamp = snd_hit->GetTime(0);
27 hit.qdc = snd_hit->GetSignal(0);
28 hit.is_x = snd_hit->isVertical();
29
30 TVector3 A, B;
31 int detectorID = snd_hit->GetDetectorID();
32 scifi_geometry->GetSiPMPosition(detectorID, A, B);
33 hit.z = A.Z();
34 if (hit.is_x)
35 {
36 hit.x = A.X();
37 hit.y = std::nan("");
38 }
39 else
40 {
41 hit.x = std::nan("");
42 hit.y = A.Y();
43 }
44 hits_.push_back(hit);
45 }
46}
void GetSiPMPosition(Int_t SiPMChan, TVector3 &A, TVector3 &B)
Definition Scifi.cxx:727
Int_t GetMat()
Definition sndScifiHit.h:36
std::vector< ScifiHit > hits_
int i
Definition ShipAna.py:86

Member Function Documentation

◆ FindCentroid()

void snd::analysis_tools::ScifiPlane::FindCentroid ( )

Definition at line 217 of file sndScifiPlane.cxx.

218{
219 centroid_.SetXYZ(0, 0, 0);
220 double tot_qdc_x{0};
221 double tot_qdc_y{0};
222 std::vector<ScifiHit> cleaned_hits = hits_;
223
224 cleaned_hits.erase(std::remove_if(cleaned_hits.begin(), cleaned_hits.end(),
225 [&](auto &hit)
226 { return hit.qdc <= 0; }),
227 cleaned_hits.end());
228 int counts_x = std::count_if(cleaned_hits.begin(), cleaned_hits.end(), [](auto &hit)
229 { return hit.is_x; });
230 int counts_y = cleaned_hits.size()-counts_x;
232 centroid_.SetXYZ(std::nan(""), std::nan(""), std::nan(""));
233 return;
234 }
235
236 for (auto &hit : cleaned_hits)
237 {
238 if (hit.qdc > 0)
239 {
240 if (hit.is_x)
241 {
242 centroid_.SetX(centroid_.X() + hit.x * hit.qdc);
243 tot_qdc_x += hit.qdc;
244 }
245 else
246 {
247 centroid_.SetY(centroid_.Y() + hit.y * hit.qdc);
248 tot_qdc_y += hit.qdc;
249 }
250
251 centroid_.SetZ(centroid_.Z() + hit.z * hit.qdc);
252 }
253 }
254 centroid_.SetXYZ((tot_qdc_x > 0) ? centroid_.X() / tot_qdc_x : std::nan(""), (tot_qdc_y > 0) ? centroid_.Y() / tot_qdc_y : std::nan(""), (tot_qdc_x+tot_qdc_y > 0) ? centroid_.Z() / (tot_qdc_x+tot_qdc_y) : std::nan(""));
256}

◆ GetCentroid()

const TVector3 snd::analysis_tools::ScifiPlane::GetCentroid ( ) const
inline

Definition at line 39 of file sndScifiPlane.h.

39{ return centroid_; };

◆ GetCentroidError()

const TVector3 snd::analysis_tools::ScifiPlane::GetCentroidError ( ) const
inline

Definition at line 40 of file sndScifiPlane.h.

40{ return centroid_error_; }

◆ GetCluster()

const TVector3 snd::analysis_tools::ScifiPlane::GetCluster ( int  max_gap) const

Definition at line 101 of file sndScifiPlane.cxx.

102{
103 std::vector<double> pos_x(configuration_.scifi_n_channels_per_plane, std::nan(""));
104 std::vector<double> pos_y(configuration_.scifi_n_channels_per_plane, std::nan(""));
105
106 for (auto &hit : hits_)
107 {
108 if (hit.qdc > 0)
109 {
110 if (hit.is_x)
111 {
112 pos_x[hit.channel_index] = hit.x;
113 }
114 else
115 {
116 pos_y[hit.channel_index] = hit.y;
117 }
118 }
119 }
120
121 auto largest_cluster = [&](const std::vector<double> &positions)
122 {
123 int n = positions.size();
124
125 int best_start = -1, best_end = -1, best_size = 0;
126 int start = -1, gap_count = 0, size = 0;
127
128 // Find the largest cluster
129 for (int i = 0; i < n; ++i)
130 {
131 if (!std::isnan(positions[i]))
132 {
133 if (start == -1)
134 start = i; // Start a new cluster
135 size++;
136 gap_count = 0; // Reset consecutive gap counter
137 }
138 else
139 {
140 gap_count++;
141 if (gap_count > max_gap)
142 { // Too many consecutive gaps, finalize previous cluster
143 if (size > best_size)
144 {
145 best_start = start;
146 best_end = i - gap_count; // End before the excessive gaps
147 best_size = size;
148 }
149 // Reset for a new potential cluster
150 start = -1;
151 gap_count = 0;
152 size = 0;
153 }
154 else
155 {
156 size++; // Gaps within max_gap still count in cluster size
157 }
158 }
159 }
160
161 // Check last cluster
162 if (size > best_size)
163 {
164 best_start = start;
165 best_end = n - 1;
166 best_size = size;
167 }
168
169 // Compute the average of non-gap values in the best cluster
170 if (best_start == -1 || best_end == -1)
171 return std::nan(""); // No valid cluster found
172
173 double sum = 0.0;
174 int count = 0;
175 for (int i = best_start; i <= best_end; ++i)
176 {
177 if (!std::isnan(positions[i]))
178 {
179 sum += positions[i];
180 count++;
181 }
182 }
183
184 return (count > 0) ? sum / count : std::nan("");
185 };
186
187 double cluster_x = largest_cluster(pos_x);
188 double cluster_y = largest_cluster(pos_y);
189 if (!(std::isnan(cluster_x) || std::isnan(cluster_y))) {
190 return TVector3(cluster_x, cluster_y, hits_[0].z);
191 }
192 return TVector3(std::nan(""), std::nan(""), std::nan(""));
193}

◆ GetHits()

const std::vector< ScifiHit > snd::analysis_tools::ScifiPlane::GetHits ( ) const
inline

Definition at line 38 of file sndScifiPlane.h.

38{ return hits_; };

◆ GetNHits()

const snd::analysis_tools::ScifiPlane::xy_pair< int > snd::analysis_tools::ScifiPlane::GetNHits ( ) const

Definition at line 48 of file sndScifiPlane.cxx.

49{
50 xy_pair<int> counts{0, 0};
51 counts.x = std::count_if(hits_.begin(), hits_.end(), [](auto &hit)
52 { return hit.is_x; });
53 counts.y = hits_.size() - counts.x;
54
55 return counts;
56}

◆ GetPointQdc()

snd::analysis_tools::ScifiPlane::xy_pair< double > snd::analysis_tools::ScifiPlane::GetPointQdc ( const TVector3 &  point,
double  radius 
) const

Definition at line 203 of file sndScifiPlane.cxx.

204{
205 xy_pair<double> qdc{0.0, 0.0};
206 for (const auto &hit : hits_) {
207 if (hit.is_x && std::abs(hit.x - point.X()) <= radius) {
208 qdc.x += hit.qdc;
209 }
210 else if (!hit.is_x && std::abs(hit.y - point.Y()) <= radius) {
211 qdc.y += hit.qdc;
212 }
213 }
214 return qdc;
215}

◆ GetStation()

const int snd::analysis_tools::ScifiPlane::GetStation ( ) const
inline

Definition at line 37 of file sndScifiPlane.h.

37{ return station_; };

◆ GetTotEnergy()

const snd::analysis_tools::ScifiPlane::xy_pair< double > snd::analysis_tools::ScifiPlane::GetTotEnergy ( bool  only_positive = false) const

Definition at line 269 of file sndScifiPlane.cxx.

270{
271 xy_pair<double> energy{0.0, 0.0};
272 auto qdc = GetTotQdc(only_positive);
273
276
277 return energy;
278}
const xy_pair< double > GetTotQdc(bool only_positive=false) const

◆ GetTotQdc()

const snd::analysis_tools::ScifiPlane::xy_pair< double > snd::analysis_tools::ScifiPlane::GetTotQdc ( bool  only_positive = false) const

Definition at line 258 of file sndScifiPlane.cxx.

259{
260 xy_pair<double> qdc_sum{0.0, 0.0};
261 qdc_sum.x = std::accumulate(hits_.begin(), hits_.end(), 0.0, [&](double current_sum, auto &hit)
262 { return (hit.is_x && (hit.qdc > 0 || !only_positive)) ? (current_sum + hit.qdc) : current_sum; });
263 qdc_sum.y = std::accumulate(hits_.begin(), hits_.end(), 0.0, [&](double current_sum, auto &hit)
264 { return (!hit.is_x && (hit.qdc > 0 || !only_positive)) ? (current_sum + hit.qdc) : current_sum; });
265
266 return qdc_sum;
267}

◆ IsShower()

bool snd::analysis_tools::ScifiPlane::IsShower ( ) const

Definition at line 58 of file sndScifiPlane.cxx.

59{
61 {
62 throw std::runtime_error{"min_hits > window_width"};
63 }
64
65 xy_pair<std::vector<int>> is_hit;
68
69 for (auto &hit : hits_)
70 {
71 (hit.is_x ? is_hit.x : is_hit.y)[hit.channel_index] = 1;
72 }
73
74 auto density = [&](std::vector<int> &hit_arr)
75 {
76 int count{0};
77
78 // Initial count for the first window
79 for (int i{0}; i < configuration_.scifi_shower_window_width; ++i)
80 {
81 count += hit_arr[i];
82 }
83
85 return true;
86
87 // Slide the window across the array
89 {
90 count += hit_arr[i] - hit_arr[i - configuration_.scifi_shower_window_width];
92 return true;
93 }
94
95 return false;
96 };
97
98 return density(is_hit.x) && density(is_hit.y);
99}

◆ TimeFilter()

void snd::analysis_tools::ScifiPlane::TimeFilter ( double  min_timestamp,
double  max_timestamp 
)

Definition at line 195 of file sndScifiPlane.cxx.

196{
197 hits_.erase(std::remove_if(hits_.begin(), hits_.end(),
198 [&](auto &hit)
199 { return hit.timestamp < min_timestamp || hit.timestamp > max_timestamp; }),
200 hits_.end());
201}

Member Data Documentation

◆ centroid_

TVector3 snd::analysis_tools::ScifiPlane::centroid_
private

Definition at line 56 of file sndScifiPlane.h.

◆ centroid_error_

TVector3 snd::analysis_tools::ScifiPlane::centroid_error_
private

Definition at line 57 of file sndScifiPlane.h.

◆ configuration_

Configuration snd::analysis_tools::ScifiPlane::configuration_
private

Definition at line 55 of file sndScifiPlane.h.

◆ hits_

std::vector<ScifiHit> snd::analysis_tools::ScifiPlane::hits_
private

Definition at line 54 of file sndScifiPlane.h.

◆ station_

int snd::analysis_tools::ScifiPlane::station_
private

Definition at line 59 of file sndScifiPlane.h.


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