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 (std::vector< sndScifiHit * > snd_hits, Configuration configuration, Scifi *scifi_geometry, int station)
 
const int GetStation () const
 
const std::vector< ScifiHitGetHits () const
 
const ROOT::Math::XYZPoint GetCentroid () const
 
const ROOT::Math::XYZPoint 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 ROOT::Math::XYZPoint GetCluster (int max_gap) const
 
void FindCentroid ()
 
bool HasShower () const
 
void TimeFilter (double min_timestamp, double max_timestamp)
 
xy_pair< double > GetPointQdc (const ROOT::Math::XYZPoint &point, double radius) const
 

Private Attributes

std::vector< ScifiHithits_
 
Configuration configuration_
 
ROOT::Math::XYZPoint centroid_
 
ROOT::Math::XYZPoint centroid_error_
 
int station_
 

Detailed Description

Definition at line 13 of file sndScifiPlane.h.

Constructor & Destructor Documentation

◆ ScifiPlane()

snd::analysis_tools::ScifiPlane::ScifiPlane ( std::vector< sndScifiHit * >  snd_hits,
Configuration  configuration,
Scifi scifi_geometry,
int  station 
)

Definition at line 15 of file sndScifiPlane.cxx.

15 : configuration_(configuration), centroid_(std::nan(""), std::nan(""), std::nan("")), centroid_error_(std::nan(""), std::nan(""), std::nan("")), station_(station)
16{
17 for ( auto snd_hit : snd_hits)
18 {
19 ScifiHit hit;
20 hit.channel_index = 512 * snd_hit->GetMat() + 128 * snd_hit->GetSiPM() + snd_hit->GetSiPMChan();
21 hit.timestamp = (scifi_geometry->GetCorrectedTime(snd_hit->GetDetectorID(), snd_hit->GetTime(0)*ShipUnit::snd_TDC2ns, 0) / ShipUnit::snd_TDC2ns); // timestamp is in clock cycles
22 hit.qdc = snd_hit->GetSignal(0);
23 hit.is_x = snd_hit->isVertical();
24
25 TVector3 A, B;
26 int detectorID = snd_hit->GetDetectorID();
27 scifi_geometry->GetSiPMPosition(detectorID, A, B);
28 hit.z = A.Z();
29 if (hit.is_x)
30 {
31 hit.x = A.X();
32 hit.y = std::nan("");
33 }
34 else
35 {
36 hit.x = std::nan("");
37 hit.y = A.Y();
38 }
39 hits_.push_back(hit);
40 }
41}
void GetSiPMPosition(Int_t SiPMChan, TVector3 &A, TVector3 &B)
Definition Scifi.cxx:727
Double_t GetCorrectedTime(Int_t fDetectorID, Double_t rawTime, Double_t L)
Definition Scifi.cxx:614
ROOT::Math::XYZPoint centroid_error_
std::vector< ScifiHit > hits_
ROOT::Math::XYZPoint centroid_

Member Function Documentation

◆ FindCentroid()

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

Definition at line 222 of file sndScifiPlane.cxx.

223{
224 centroid_.SetXYZ(0, 0, 0);
225 double tot_qdc_x{0};
226 double tot_qdc_y{0};
227 std::vector<ScifiHit> cleaned_hits = hits_;
228
229 cleaned_hits.erase(std::remove_if(cleaned_hits.begin(), cleaned_hits.end(),
230 [&](auto &hit)
231 { return hit.qdc <= 0; }),
232 cleaned_hits.end());
233 int counts_x = std::count_if(cleaned_hits.begin(), cleaned_hits.end(), [](auto &hit)
234 { return hit.is_x; });
235 int counts_y = static_cast<int>(cleaned_hits.size())-counts_x;
237 centroid_ = ROOT::Math::XYZPoint(std::nan(""), std::nan(""), std::nan(""));
238 return;
239 }
240
241 for (auto &hit : cleaned_hits)
242 {
243 if (hit.qdc > 0)
244 {
245 if (hit.is_x)
246 {
247 centroid_.SetX(centroid_.X() + hit.x * hit.qdc);
248 tot_qdc_x += hit.qdc;
249 }
250 else
251 {
252 centroid_.SetY(centroid_.Y() + hit.y * hit.qdc);
253 tot_qdc_y += hit.qdc;
254 }
255
256 centroid_.SetZ(centroid_.Z() + hit.z * hit.qdc);
257 }
258 }
259 centroid_ = ROOT::Math::XYZPoint((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(""));
261}

◆ GetCentroid()

const ROOT::Math::XYZPoint snd::analysis_tools::ScifiPlane::GetCentroid ( ) const
inline

Definition at line 39 of file sndScifiPlane.h.

39{ return centroid_; };

◆ GetCentroidError()

const ROOT::Math::XYZPoint snd::analysis_tools::ScifiPlane::GetCentroidError ( ) const
inline

Definition at line 40 of file sndScifiPlane.h.

40{ return centroid_error_; }

◆ GetCluster()

const ROOT::Math::XYZPoint snd::analysis_tools::ScifiPlane::GetCluster ( int  max_gap) const

Definition at line 106 of file sndScifiPlane.cxx.

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

◆ 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 43 of file sndScifiPlane.cxx.

44{
45 xy_pair<int> counts{0, 0};
46 counts.x = std::count_if(hits_.begin(), hits_.end(), [](auto &hit)
47 { return hit.is_x; });
48 counts.y = static_cast<int>(hits_.size()) - counts.x;
49
50 return counts;
51}

◆ GetPointQdc()

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

Definition at line 208 of file sndScifiPlane.cxx.

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

◆ 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 274 of file sndScifiPlane.cxx.

275{
276 xy_pair<double> energy{0.0, 0.0};
277 auto qdc = GetTotQdc(only_positive);
278
281
282 return energy;
283}
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 263 of file sndScifiPlane.cxx.

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

◆ HasShower()

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

Definition at line 53 of file sndScifiPlane.cxx.

54{
56 {
57 throw std::runtime_error{"min_hits > window_width"};
58 }
59
61 {
62 throw std::runtime_error{"window_width > n_channels_per_plane"};
63 }
64
65 xy_pair<std::vector<int>> is_hit;
68
69 for (auto &hit : hits_)
70 {
71 if (hit.channel_index < 0 ||
72 hit.channel_index >= configuration_.scifi_n_channels_per_plane)
73 {
74 throw std::out_of_range{"hit.channel_index out of range"};
75 }
76 (hit.is_x ? is_hit.x : is_hit.y)[hit.channel_index] = 1;
77 }
78
79 auto density = [&](const std::vector<int> &hit_arr)
80 {
81 int count{0};
82
83 // Initial count for the first window
84 for (int i{0}; i < configuration_.scifi_shower_window_width; ++i)
85 {
86 count += hit_arr[i];
87 }
88
90 return true;
91
92 // Slide the window across the array
94 {
95 count += hit_arr[i] - hit_arr[i - configuration_.scifi_shower_window_width];
97 return true;
98 }
99
100 return false;
101 };
102
103 return density(is_hit.x) && density(is_hit.y);
104}

◆ TimeFilter()

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

Definition at line 200 of file sndScifiPlane.cxx.

201{
202 hits_.erase(std::remove_if(hits_.begin(), hits_.end(),
203 [&](auto &hit)
204 { return hit.timestamp < min_timestamp || hit.timestamp > max_timestamp; }),
205 hits_.end());
206}

Member Data Documentation

◆ centroid_

ROOT::Math::XYZPoint snd::analysis_tools::ScifiPlane::centroid_
private

Definition at line 56 of file sndScifiPlane.h.

◆ centroid_error_

ROOT::Math::XYZPoint 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: