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 "R3BDataMonitor.h"
16#include "R3BException.h"
17#include "R3BNeulandCalData2.h"
20#include "R3BNeulandCommon.h"
21#include "R3BNeulandHit.h"
22#include "R3BValueError.h"
23#include <FairRootManager.h>
24#include <FairRuntimeDb.h>
26#include <R3BShared.h>
27#include <TVector3.h>
28#include <cmath>
29#include <fairlogger/Logger.h>
30#include <fmt/core.h>
31#include <string_view>
32#include <vector>
33
34namespace
35{
36 using R3B::ValueErrorD;
39
40 inline auto get_hit_energy(double first_e, double second_e, const HitModulePar& par)
41 {
42 // const auto attenuation_value = std::exp(R3B::Neuland::TotalBarLength / par.lightAttenuationLength.value);
43 return par.light_attenuation_factor.value * std::sqrt(first_e * second_e);
44 }
45
46 inline auto get_hit_position(double tdc_left, double tdc_right, const HitModulePar& par)
47 {
48 const auto plane_id = ::R3B::Neuland::ModuleID2PlaneID(static_cast<int>(par.module_num - 1));
49 const auto bar_num = par.module_num % ::R3B::Neuland::BarsPerPlane;
50 LOGP(debug2,
51 "Calculating position with left tdc: {}, right tdc {}, effective speed: {}, tdc_diff: {}",
52 tdc_left,
53 tdc_right,
54 par.effective_speed,
55 tdc_right - tdc_left);
56 const auto pos_along_bar = par.effective_speed.value * (tdc_right - tdc_left);
57 const auto pos_perp_bar = (bar_num - 0.5 - ::R3B::Neuland::BarsPerPlane * 0.5) * ::R3B::Neuland::BarSize_XY;
58 const auto pos_z = (plane_id + 0.5) * ::R3B::Neuland::BarSize_Z;
59 LOGP(debug2, "pos along the bar: {} cm, pos perp to bar {} cm, z: {} cm", pos_along_bar, pos_perp_bar, pos_z);
60
61 auto is_horizontal = ::R3B::Neuland::IsPlaneIDHorizontal(plane_id);
62 return is_horizontal ? TVector3{ pos_along_bar, pos_perp_bar, pos_z }
63 : TVector3{ pos_perp_bar, pos_along_bar, pos_z };
64 }
65
66 inline auto get_hit_pixel(const TVector3& position)
67 {
68 const auto pixel_x =
69 std::floor(position.X() / ::R3B::Neuland::BarSize_XY) + (::R3B::Neuland::BarsPerPlane / 2.);
70 const auto pixel_y =
71 std::floor(position.Y() / ::R3B::Neuland::BarSize_XY) + (::R3B::Neuland::BarsPerPlane / 2.);
72 const auto pixel_z = std::floor(position.Z() / ::R3B::Neuland::BarSize_Z);
73
74 return TVector3{ pixel_x, pixel_y, pixel_z };
75 }
76} // namespace
77
78namespace R3B::Neuland
79{
80
81 Cal2HitTask::Cal2HitTask(std::string_view input_cal_data_name,
82 std::string_view output_hit_data_name,
83 std::string_view input_cal_2_hit_par_name)
84 : CalibrationTask("R3BNeulandCal2Hit", 1)
85 , cal_data_{ input_cal_data_name }
86 , hit_data_{ output_hit_data_name }
87 , cal_to_hit_par_{ input_cal_2_hit_par_name }
88 {
89 }
90
92
93 void Cal2HitTask::ExtraInit(FairRootManager* /*rootMan*/)
94 {
95 cal_to_hit_par_.init(this);
96 cal_data_.init();
97 hit_data_.init(not IsHistDisabled());
98 }
99
100 void Cal2HitTask::SetExtraPar(FairRuntimeDb* rtdb) {}
101
103
105 {
106 for (const auto& calBar : cal_data_)
107 {
108 temp_left_signals_.clear();
109 temp_right_signals_.clear();
110
111 LOGP(debug1, "Input calBar: {}", calBar);
112 if (calBar.module_num == 0)
113 {
114 throw R3B::runtime_error("cal-level bar signal has invalid moudule number 0!");
115 }
116
119
120 const auto& module_par = cal_to_hit_par_->GetModuleParAt(calBar.module_num);
122 }
123 }
124
126 /* inout */ std::vector<CalibratedSignal>& signals,
127 Side side)
128 {
129 const auto calBar_signals = (side == Side::left) ? calBar.left : calBar.right;
130
131 const auto& module_par = cal_to_hit_par_->GetModuleParAt(calBar.module_num);
132
133 for (const auto& cal_signal : calBar_signals)
134 {
135 signals.push_back(to_calibrated_signal(cal_signal, module_par, side));
136 }
137 }
138
140 const HitModulePar& par) const -> R3BNeulandHit
141 {
142 auto hit = R3BNeulandHit{};
143
144 LOGP(debug2,
145 "Input left calibrated signal: {} and right calibrated signal: {}",
146 signalPair.left(),
147 signalPair.right());
148 // hit module id is 0-based
149 hit.module_id = static_cast<int>(par.module_num - 1);
150 hit.tdc_left = signalPair.left().time.value;
151 hit.tdc_right = signalPair.right().time.value;
152 hit.time = get_hit_time(hit.tdc_left, hit.tdc_right);
153 hit.qdc_left = signalPair.left().energy.value;
154 hit.qdc_right = signalPair.right().energy.value;
155 hit.energy = get_hit_energy(hit.qdc_left, hit.qdc_right, par);
156 hit.position = get_hit_position(hit.tdc_left, hit.tdc_right, par);
157 hit.pixel = get_hit_pixel(hit.position);
158 LOGP(debug, "Adding a new NeulandHit: {}\n", hit);
159 return hit;
160 }
161
162 void Cal2HitTask::construct_hits(const std::vector<CalibratedSignal>& left_signals,
163 const std::vector<CalibratedSignal>& right_signals,
164 const HitModulePar& par,
165 /* inout */ std::vector<R3BNeulandHit>& hits)
166 {
167 // TODO: Multi-hits needs to be implemented here
168 if (left_signals.size() == 1 and right_signals.size() == 1)
169 {
170 const auto& left_signal = left_signals.front();
171 const auto& right_signal = right_signals.front();
172 // signal_match_checking(left_signal, right_signal, par);
173 hits.push_back(construct_hit(R3B::LRPair<CalibratedSignal>{ left_signal, right_signal }, par));
174 }
175 }
176
178 const CalibratedSignal& second_signal,
179 const HitModulePar& par) -> bool
180 {
181 const auto first_input = SignalMatcher::Input{ first_signal.time.value, first_signal.energy.value };
182 const auto second_input = SignalMatcher::Input{ second_signal.time.value, second_signal.energy.value };
183
184 const auto match_par =
185 SignalMatcher::Par{ BarLength / par.light_attenuation_length.value, par.effective_speed.value };
186 const auto match_goodness = SignalMatcher::GetGoodnessOfMatch(first_input, second_input, match_par);
187 const auto match_result = std::log10(match_goodness);
188 // FIXME: BAD comparison values. They should be 0. Why? Need fixing
189 GetHistMonitor().get("match_values")->Fill(match_result);
190 return true;
191 }
192
194
195 // auto Cal2Hit::CheckConditions() const -> bool
196 // {
197 // const auto t_start = GetEventHeader()->GetTStart();
198 // const auto is_beam_on = not std::isnan(t_start);
199 // return is_beam_on;
200 // }
201
202 auto Cal2HitTask::get_hit_time(double first_t, double second_t) const -> double
203 {
204 auto larger_t = (first_t > second_t) ? first_t : second_t;
205 auto smaller_t = (first_t < second_t) ? first_t : second_t;
206
207 // check if the smaller value is overflowed
208 // TODO: overflow should be checked in cal level for better time calibration. But the cal level calibration
209 // doesn't have singal paring. Solution: Do it in cal level only for one hit bar signal and calibration should
210 // only use those with one hit.
211 if (larger_t - smaller_t > 0.5 * R3B::Neuland::MaxCalTime)
212 {
213 larger_t -= R3B::Neuland::MaxCalTime;
214 }
215 return std::remainder(((larger_t + smaller_t) / 2.) - global_time_offset_ - GetEventHeader()->GetTStart(),
217 }
218
220 const HitModulePar& par,
221 R3B::Side side) -> ValueErrorD
222 {
223 const auto tot_no_offset = calSignal.time_over_threshold - par.pedestal.get(side);
224
225 // apply minimum 1 ns:
226 return (tot_no_offset.value < 1)
227 ? ValueErrorD{}
228 : tot_no_offset / (par.energy_gain.get(side) - par.pmt_saturation.get(side) * tot_no_offset);
229 }
230
232 const HitModulePar& par,
233 R3B::Side side) -> ValueErrorD
234 {
235 // TODO: why positive for left?
236 const auto time_offset =
237 (side == R3B::Side::left) ? (par.t_sync - par.t_diff / 2) : (par.t_sync + par.t_diff / 2);
238 return calSignal.leading_time - calSignal.trigger_time - time_offset;
239 }
240
242 const HitModulePar& par,
244 {
245 const auto energy = get_calibrated_energy(calSignal, par, side);
246 const auto time = get_calibrated_time(calSignal, par, side);
247 return CalibratedSignal{ energy, time };
248 }
249} // namespace R3B::Neuland
static auto get_calibrated_energy(const CalDataSignal &calSignal, const HitModulePar &par, R3B::Side side) -> ValueErrorD
auto construct_hit(const LRPair< CalibratedSignal > &signalPair, const HitModulePar &par) const -> R3BNeulandHit
void SetExtraPar(FairRuntimeDb *rtdb) override
void ExtraInit(FairRootManager *rootMan) override
InputParView< Cal2HitPar > cal_to_hit_par_
static auto get_calibrated_time(const CalDataSignal &calSignal, const HitModulePar &par, R3B::Side side) -> ValueErrorD
static auto to_calibrated_signal(const CalDataSignal &calSignal, const HitModulePar &par, R3B::Side side) -> CalibratedSignal
std::vector< CalibratedSignal > temp_left_signals_
void calculate_calibrated_signals(const BarCalData &calBar, std::vector< CalibratedSignal > &signals, Side side)
OutputVectorConnector< R3BNeulandHit > hit_data_
auto signal_match_checking(const CalibratedSignal &first_signal, const CalibratedSignal &second_signal, const HitModulePar &par) -> bool
void HistogramInit(DataMonitor &histograms) override
void construct_hits(const std::vector< CalibratedSignal > &left_signals, const std::vector< CalibratedSignal > &right_signals, const HitModulePar &par, std::vector< R3BNeulandHit > &hits)
InputVectorConnector< BarCalData > cal_data_
auto get_hit_time(double first_t, double second_t) const -> double
std::vector< CalibratedSignal > temp_right_signals_
Cal2HitTask(std::string_view input_cal_data_name="NeulandCalData", std::string_view output_hit_data_name="NeulandHits", std::string_view input_cal_2_hit_par_name="NeulandHitPar")
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
std::vector< CalDataSignal > left
std::vector< CalDataSignal > right
ValueError< double > light_attenuation_factor
exp(alpha*L/2)
static auto GetGoodnessOfMatch(const Input &firstSignal, const Input &secondSignal, const Par &par) -> float