R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandCalToHitTask.cxx
Go to the documentation of this file.
1/******************************************************************************
2 * Copyright (C) 2019 GSI Helmholtzzentrum für Schwerionenforschung GmbH *
3 * Copyright (C) 2019-2023 Members of R3B Collaboration *
4 * *
5 * This software is distributed under the terms of the *
6 * GNU General Public Licence (GPL) version 3, *
7 * copied verbatim in the file "LICENSE". *
8 * *
9 * In applying this license GSI does not waive the privileges and immunities *
10 * granted to it by virtue of its status as an Intergovernmental Organization *
11 * or submit itself to any jurisdiction. *
12 ******************************************************************************/
13
15#include <R3BLogger.h>
17#include <R3BShared.h>
18
19namespace
20{
21 using R3B::ValueErrorD;
24
25 inline auto get_hit_energy(double first_e, double second_e, const HitModulePar& par)
26 {
27 // const auto attenuation_value = std::exp(R3B::Neuland::TotalBarLength / par.lightAttenuationLength.value);
28 return par.light_attenuation_factor.value * std::sqrt(first_e * second_e);
29 }
30
31 inline auto get_hit_position(double tdc_left, double tdc_right, const HitModulePar& par)
32 {
33 const auto plane_id = ::R3B::Neuland::ModuleID2PlaneID(static_cast<int>(par.module_num - 1));
34 const auto bar_num = par.module_num % ::R3B::Neuland::BarsPerPlane;
35 R3BLOG(debug2,
36 fmt::format("Calculating position with left tdc: {}, right tdc {}, effective speed: {}, tdc_diff: {}",
37 tdc_left,
38 tdc_right,
39 par.effective_speed,
40 tdc_right - tdc_left));
41 const auto pos_along_bar = par.effective_speed.value * (tdc_right - tdc_left);
42 const auto pos_perp_bar = (bar_num - 0.5 - ::R3B::Neuland::BarsPerPlane * 0.5) * ::R3B::Neuland::BarSize_XY;
43 const auto pos_z = (plane_id + 0.5) * ::R3B::Neuland::BarSize_Z;
44 R3BLOG(debug2,
45 fmt::format(
46 "pos along the bar: {} cm, pos perp to bar {} cm, z: {} cm", pos_along_bar, pos_perp_bar, pos_z));
47
48 auto is_horizontal = ::R3B::Neuland::IsPlaneIDHorizontal(plane_id);
49 return is_horizontal ? TVector3{ pos_along_bar, pos_perp_bar, pos_z }
50 : TVector3{ pos_perp_bar, pos_along_bar, pos_z };
51 }
52
53 inline auto get_hit_pixel(const TVector3& position)
54 {
55 const auto pixel_x =
56 std::floor(position.X() / ::R3B::Neuland::BarSize_XY) + (::R3B::Neuland::BarsPerPlane / 2.);
57 const auto pixel_y =
58 std::floor(position.Y() / ::R3B::Neuland::BarSize_XY) + (::R3B::Neuland::BarsPerPlane / 2.);
59 const auto pixel_z = std::floor(position.Z() / ::R3B::Neuland::BarSize_Z);
60
61 return TVector3{ pixel_x, pixel_y, pixel_z };
62 }
63} // namespace
64
65namespace R3B::Neuland
66{
67
68 Cal2HitTask::Cal2HitTask(std::string_view name, int iVerbose)
69 : CalibrationTask(name, iVerbose)
70 {
71 }
72
74
75 void Cal2HitTask::ExtraInit(FairRootManager* /*rootMan*/)
76 {
77 cal_data_.init();
78 hit_data_.init(not IsHistDisabled());
79 }
80
81 void Cal2HitTask::SetExtraPar(FairRuntimeDb* rtdb) {}
82
83 void Cal2HitTask::TriggeredExec() { calibrate(); }
84
85 void Cal2HitTask::calibrate()
86 {
87 for (const auto& calBar : cal_data_)
88 {
89 temp_left_signals_.clear();
90 temp_right_signals_.clear();
91
92 R3BLOG(debug1, fmt::format("Input calBar: {}", calBar));
93 if (calBar.module_num == 0)
94 {
95 throw R3B::runtime_error("cal-level bar signal has invalid moudule number 0!");
96 }
97
98 calculate_calibrated_signals(calBar, temp_left_signals_, Side::left);
99 calculate_calibrated_signals(calBar, temp_right_signals_, Side::right);
100
101 const auto& module_par = cal_to_hit_par_->GetModuleParAt(calBar.module_num);
102 construct_hits(temp_left_signals_, temp_right_signals_, module_par, hit_data_.get());
103 }
104 }
105
106 void Cal2HitTask::calculate_calibrated_signals(const BarCalData& calBar,
107 /* inout */ std::vector<CalibratedSignal>& signals,
108 Side side)
109 {
110 const auto calBar_signals = (side == Side::left) ? calBar.left : calBar.right;
111
112 const auto& module_par = cal_to_hit_par_->GetModuleParAt(calBar.module_num);
113
114 for (const auto& cal_signal : calBar_signals)
115 {
116 signals.push_back(to_calibrated_signal(cal_signal, module_par, side));
117 }
118 }
119
120 auto Cal2HitTask::construct_hit(const LRPair<CalibratedSignal>& signalPair,
121 const HitModulePar& par) const -> R3BNeulandHit
122 {
123 auto hit = R3BNeulandHit{};
124
125 R3BLOG(debug2,
126 fmt::format("Input left calibrated signal: {} and right calibrated signal: {}",
127 signalPair.left(),
128 signalPair.right()));
129 // hit module id is 0-based
130 hit.module_id = static_cast<int>(par.module_num - 1);
131 hit.tdc_left = signalPair.left().time.value;
132 hit.tdc_right = signalPair.right().time.value;
133 hit.time = get_hit_time(hit.tdc_left, hit.tdc_right);
134 hit.qdc_left = signalPair.left().energy.value;
135 hit.qdc_right = signalPair.right().energy.value;
136 hit.energy = get_hit_energy(hit.qdc_left, hit.qdc_right, par);
137 hit.position = get_hit_position(hit.tdc_left, hit.tdc_right, par);
138 hit.pixel = get_hit_pixel(hit.position);
139 R3BLOG(debug, fmt::format("Adding a new NeulandHit: {}\n", hit));
140 return hit;
141 }
142
143 void Cal2HitTask::construct_hits(const std::vector<CalibratedSignal>& left_signals,
144 const std::vector<CalibratedSignal>& right_signals,
145 const HitModulePar& par,
146 /* inout */ std::vector<R3BNeulandHit>& hits)
147 {
148 // TODO: Multi-hits needs to be implemented here
149 if (left_signals.size() == 1 and right_signals.size() == 1)
150 {
151 const auto& left_signal = left_signals.front();
152 const auto& right_signal = right_signals.front();
153 // signal_match_checking(left_signal, right_signal, par);
154 hits.push_back(construct_hit(R3B::LRPair<CalibratedSignal>{ left_signal, right_signal }, par));
155 }
156 }
157
158 auto Cal2HitTask::signal_match_checking(const CalibratedSignal& first_signal,
159 const CalibratedSignal& second_signal,
160 const HitModulePar& par) -> bool
161 {
162 const auto first_input = SignalMatcher::Input{ first_signal.time.value, first_signal.energy.value };
163 const auto second_input = SignalMatcher::Input{ second_signal.time.value, second_signal.energy.value };
164
165 const auto match_par =
166 SignalMatcher::Par{ BarLength / par.light_attenuation_length.value, par.effective_speed.value };
167 const auto match_goodness = SignalMatcher::GetGoodnessOfMatch(first_input, second_input, match_par);
168 const auto match_result = std::log10(match_goodness);
169 // FIXME: BAD comparison values. They should be 0. Why? Need fixing
170 GetHistMonitor().get("match_values")->Fill(match_result);
171 return true;
172 }
173
175
176 // auto Cal2Hit::CheckConditions() const -> bool
177 // {
178 // const auto t_start = GetEventHeader()->GetTStart();
179 // const auto is_beam_on = not std::isnan(t_start);
180 // return is_beam_on;
181 // }
182
183 auto Cal2HitTask::get_hit_time(double first_t, double second_t) const -> double
184 {
185 auto larger_t = (first_t > second_t) ? first_t : second_t;
186 auto smaller_t = (first_t < second_t) ? first_t : second_t;
187
188 // check if the smaller value is overflowed
189 // TODO: overflow should be checked in cal level for better time calibration. But the cal level calibration
190 // doesn't have singal paring. Solution: Do it in cal level only for one hit bar signal and calibration should
191 // only use those with one hit.
192 if (larger_t - smaller_t > 0.5 * R3B::Neuland::MaxCalTime)
193 {
194 larger_t -= R3B::Neuland::MaxCalTime;
195 }
196 return std::remainder(((larger_t + smaller_t) / 2.) - global_time_offset_ - GetEventHeader()->GetTStart(),
198 }
199
200 auto Cal2HitTask::get_calibrated_energy(const CalDataSignal& calSignal,
201 const HitModulePar& par,
202 R3B::Side side) -> ValueErrorD
203 {
204 const auto tot_no_offset = calSignal.time_over_threshold - par.pedestal.get(side);
205
206 // apply minimum 1 ns:
207 return (tot_no_offset.value < 1)
208 ? ValueErrorD{}
209 : tot_no_offset / (par.energy_gain.get(side) - par.pmt_saturation.get(side) * tot_no_offset);
210 }
211
212 auto Cal2HitTask::get_calibrated_time(const CalDataSignal& calSignal,
213 const HitModulePar& par,
214 R3B::Side side) -> ValueErrorD
215 {
216 // TODO: why positive for left?
217 const auto time_offset =
218 (side == R3B::Side::left) ? (par.t_sync - par.t_diff / 2) : (par.t_sync + par.t_diff / 2);
219 return calSignal.leading_time - calSignal.trigger_time - time_offset;
220 }
221
222 auto Cal2HitTask::to_calibrated_signal(const CalDataSignal& calSignal,
223 const HitModulePar& par,
225 {
226 const auto energy = get_calibrated_energy(calSignal, par, side);
227 const auto time = get_calibrated_time(calSignal, par, side);
228 return CalibratedSignal{ energy, time };
229 }
230} // namespace R3B::Neuland
#define R3BLOG(severity, x)
Definition R3BLogger.h:35
auto GetModuleParAt(unsigned int module_num) const -> const ::R3B::Neuland::HitModulePar &
void SetExtraPar(FairRuntimeDb *rtdb) override
void ExtraInit(FairRootManager *rootMan) override
Cal2HitTask(std::string_view name="R3BNeulandCal2Hit", int iVerbose=1)
void HistogramInit(DataMonitor &histograms) override
Simulation of NeuLAND Bar/Paddle.
constexpr auto BarsPerPlane
constexpr auto ModuleID2PlaneID(int moduleID) -> int
constexpr auto BarLength
constexpr auto BarSize_Z
constexpr auto BarSize_XY
constexpr auto IsPlaneIDHorizontal(int plane_id) -> bool
constexpr auto MaxCalTime
ValueError< double > ValueErrorD
ValueError< double > light_attenuation_factor
static auto GetGoodnessOfMatch(const Input &firstSignal, const Input &secondSignal, const Par &par) -> float