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>
25#include <R3BLogger.h>
27#include <R3BShared.h>
28#include <TVector3.h>
29#include <cmath>
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 R3BLOG(debug2,
51 fmt::format("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 R3BLOG(debug2,
60 fmt::format(
61 "pos along the bar: {} cm, pos perp to bar {} cm, z: {} cm", pos_along_bar, pos_perp_bar, pos_z));
62
63 auto is_horizontal = ::R3B::Neuland::IsPlaneIDHorizontal(plane_id);
64 return is_horizontal ? TVector3{ pos_along_bar, pos_perp_bar, pos_z }
65 : TVector3{ pos_perp_bar, pos_along_bar, pos_z };
66 }
67
68 inline auto get_hit_pixel(const TVector3& position)
69 {
70 const auto pixel_x =
71 std::floor(position.X() / ::R3B::Neuland::BarSize_XY) + (::R3B::Neuland::BarsPerPlane / 2.);
72 const auto pixel_y =
73 std::floor(position.Y() / ::R3B::Neuland::BarSize_XY) + (::R3B::Neuland::BarsPerPlane / 2.);
74 const auto pixel_z = std::floor(position.Z() / ::R3B::Neuland::BarSize_Z);
75
76 return TVector3{ pixel_x, pixel_y, pixel_z };
77 }
78} // namespace
79
80namespace R3B::Neuland
81{
82
83 Cal2HitTask::Cal2HitTask(std::string_view name, int iVerbose)
84 : CalibrationTask(name, iVerbose)
85 {
86 }
87
89
90 void Cal2HitTask::ExtraInit(FairRootManager* /*rootMan*/)
91 {
92 cal_data_.init();
93 hit_data_.init(not IsHistDisabled());
94 }
95
96 void Cal2HitTask::SetExtraPar(FairRuntimeDb* rtdb) {}
97
99
101 {
102 for (const auto& calBar : cal_data_)
103 {
104 temp_left_signals_.clear();
105 temp_right_signals_.clear();
106
107 R3BLOG(debug1, fmt::format("Input calBar: {}", calBar));
108 if (calBar.module_num == 0)
109 {
110 throw R3B::runtime_error("cal-level bar signal has invalid moudule number 0!");
111 }
112
115
116 const auto& module_par = cal_to_hit_par_->GetModuleParAt(calBar.module_num);
118 }
119 }
120
122 /* inout */ std::vector<CalibratedSignal>& signals,
123 Side side)
124 {
125 const auto calBar_signals = (side == Side::left) ? calBar.left : calBar.right;
126
127 const auto& module_par = cal_to_hit_par_->GetModuleParAt(calBar.module_num);
128
129 for (const auto& cal_signal : calBar_signals)
130 {
131 signals.push_back(to_calibrated_signal(cal_signal, module_par, side));
132 }
133 }
134
136 const HitModulePar& par) const -> R3BNeulandHit
137 {
138 auto hit = R3BNeulandHit{};
139
140 R3BLOG(debug2,
141 fmt::format("Input left calibrated signal: {} and right calibrated signal: {}",
142 signalPair.left(),
143 signalPair.right()));
144 // hit module id is 0-based
145 hit.module_id = static_cast<int>(par.module_num - 1);
146 hit.tdc_left = signalPair.left().time.value;
147 hit.tdc_right = signalPair.right().time.value;
148 hit.time = get_hit_time(hit.tdc_left, hit.tdc_right);
149 hit.qdc_left = signalPair.left().energy.value;
150 hit.qdc_right = signalPair.right().energy.value;
151 hit.energy = get_hit_energy(hit.qdc_left, hit.qdc_right, par);
152 hit.position = get_hit_position(hit.tdc_left, hit.tdc_right, par);
153 hit.pixel = get_hit_pixel(hit.position);
154 R3BLOG(debug, fmt::format("Adding a new NeulandHit: {}\n", hit));
155 return hit;
156 }
157
158 void Cal2HitTask::construct_hits(const std::vector<CalibratedSignal>& left_signals,
159 const std::vector<CalibratedSignal>& right_signals,
160 const HitModulePar& par,
161 /* inout */ std::vector<R3BNeulandHit>& hits)
162 {
163 // TODO: Multi-hits needs to be implemented here
164 if (left_signals.size() == 1 and right_signals.size() == 1)
165 {
166 const auto& left_signal = left_signals.front();
167 const auto& right_signal = right_signals.front();
168 // signal_match_checking(left_signal, right_signal, par);
169 hits.push_back(construct_hit(R3B::LRPair<CalibratedSignal>{ left_signal, right_signal }, par));
170 }
171 }
172
174 const CalibratedSignal& second_signal,
175 const HitModulePar& par) -> bool
176 {
177 const auto first_input = SignalMatcher::Input{ first_signal.time.value, first_signal.energy.value };
178 const auto second_input = SignalMatcher::Input{ second_signal.time.value, second_signal.energy.value };
179
180 const auto match_par =
181 SignalMatcher::Par{ BarLength / par.light_attenuation_length.value, par.effective_speed.value };
182 const auto match_goodness = SignalMatcher::GetGoodnessOfMatch(first_input, second_input, match_par);
183 const auto match_result = std::log10(match_goodness);
184 // FIXME: BAD comparison values. They should be 0. Why? Need fixing
185 GetHistMonitor().get("match_values")->Fill(match_result);
186 return true;
187 }
188
190
191 // auto Cal2Hit::CheckConditions() const -> bool
192 // {
193 // const auto t_start = GetEventHeader()->GetTStart();
194 // const auto is_beam_on = not std::isnan(t_start);
195 // return is_beam_on;
196 // }
197
198 auto Cal2HitTask::get_hit_time(double first_t, double second_t) const -> double
199 {
200 auto larger_t = (first_t > second_t) ? first_t : second_t;
201 auto smaller_t = (first_t < second_t) ? first_t : second_t;
202
203 // check if the smaller value is overflowed
204 // TODO: overflow should be checked in cal level for better time calibration. But the cal level calibration
205 // doesn't have singal paring. Solution: Do it in cal level only for one hit bar signal and calibration should
206 // only use those with one hit.
207 if (larger_t - smaller_t > 0.5 * R3B::Neuland::MaxCalTime)
208 {
209 larger_t -= R3B::Neuland::MaxCalTime;
210 }
211 return std::remainder(((larger_t + smaller_t) / 2.) - global_time_offset_ - GetEventHeader()->GetTStart(),
213 }
214
216 const HitModulePar& par,
217 R3B::Side side) -> ValueErrorD
218 {
219 const auto tot_no_offset = calSignal.time_over_threshold - par.pedestal.get(side);
220
221 // apply minimum 1 ns:
222 return (tot_no_offset.value < 1)
223 ? ValueErrorD{}
224 : tot_no_offset / (par.energy_gain.get(side) - par.pmt_saturation.get(side) * tot_no_offset);
225 }
226
228 const HitModulePar& par,
229 R3B::Side side) -> ValueErrorD
230 {
231 // TODO: why positive for left?
232 const auto time_offset =
233 (side == R3B::Side::left) ? (par.t_sync - par.t_diff / 2) : (par.t_sync + par.t_diff / 2);
234 return calSignal.leading_time - calSignal.trigger_time - time_offset;
235 }
236
238 const HitModulePar& par,
240 {
241 const auto energy = get_calibrated_energy(calSignal, par, side);
242 const auto time = get_calibrated_time(calSignal, par, side);
243 return CalibratedSignal{ energy, time };
244 }
245} // namespace R3B::Neuland
#define R3BLOG(severity, x)
Definition R3BLogger.h:32
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
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
Cal2HitTask(std::string_view name="R3BNeulandCal2Hit", int iVerbose=1)
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_
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
static auto GetGoodnessOfMatch(const Input &firstSignal, const Input &secondSignal, const Par &par) -> float