R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandMapToCalTask.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 "R3BNeulandCalData2.h"
18#include "R3BNeulandCommon.h"
21#include "R3BShared.h"
22#include "R3BValueError.h"
23#include <FairRootManager.h>
24#include <R3BEventHeader.h>
25#include <R3BException.h>
26#include <R3BLogger.h>
27#include <TH1.h>
28#include <TH2.h>
29#include <algorithm>
30#include <cmath>
31#include <fmt/core.h>
32#include <fmt/format.h>
33#include <iterator>
34#include <range/v3/view/map.hpp>
35#include <string_view>
36#include <vector>
37
38namespace
39{
40 // NOLINTBEGIN
41 auto GetWalkCorrection(double value) -> double
42 {
43 value = value * 17.2278 / 162.464;
44 const auto veronika3 = -4.29359 + 17.3841 / sqrt(value) + 0.073181; // veronika 3
45 const auto veronika2 = 1.59667 + 78.5274 * pow(value, -1.97051) - 4.00192 / value - 0.125473 * value +
46 0.00130958 * value * value; // veronika 2
47 const auto result = (value > 25.0000) ? veronika3 : (veronika3 + veronika2) / 2.;
48 return -1. * result;
49 }
50 // NOLINTEND
51
52} // namespace
53
54namespace R3B::Neuland
55{
56
58 : Map2CalTask{ "NeulandMapped2Cal", 1 }
59 {
60 }
61
62 Map2CalTask::Map2CalTask(std::string_view name, int iVerbose)
63 : CalibrationTask(name.data(), iVerbose)
64 {
65 }
66
68 {
69 const auto total_bar_num = plane_num_ * BarsPerPlane;
70 constexpr auto BINSIZE_TOT = 1000;
71 constexpr auto BINSIZE_LT = 1000;
72 constexpr auto BINSIZE_TIME = 8000;
73 constexpr auto HIST_MAX_TIME = 6000.;
74 constexpr auto HIST_MIN_TIME = -5000.;
75 constexpr auto HIST_MAX_TOT = 600.;
76 constexpr auto HIST_MAX_LT = 2000.;
77
78 histograms.add_hist<TH1D>("module_num", "Counts with module ids", total_bar_num, -0.5, total_bar_num + 0.5);
79 histograms.add_hist<TH1D>("ToT", "Time over threshold", BINSIZE_TOT * 2, 0., BINSIZE_TOT + 0.);
80 histograms.add_hist<TH1D>("LeadingTime", "Leading time", BINSIZE_LT, 0., HIST_MAX_LT);
81 histograms.add_hist<TH1D>("TriggerTime", "Trigger time", BINSIZE_LT, 0., HIST_MAX_LT);
82 histograms.add_hist<TH1D>("Leading_minus_trigger", "Time", BINSIZE_TIME, HIST_MIN_TIME, HIST_MAX_TIME);
83 histograms.add_hist<TH1I>(
84 "Bar_hitNum_l", "Number of hits on the left PMT per bar", total_bar_num, 0.5, 0.5 + total_bar_num);
85 histograms.add_hist<TH1I>(
86 "Bar_hitNum_r", "Number of hits on the right PMT per bar", total_bar_num, 0.5, 0.5 + total_bar_num);
87
88 histograms.add_hist<TH2D>("TimeLVsBar",
89 "Time_vs_bar_left",
90 total_bar_num,
91 0.5,
92 0.5 + total_bar_num,
93 BINSIZE_TIME,
94 HIST_MIN_TIME,
95 HIST_MAX_TIME);
96 histograms.add_hist<TH2D>("TimeRVsBar",
97 "Time_vs_bar_right",
98 total_bar_num,
99 0.5,
100 0.5 + total_bar_num,
101 BINSIZE_TIME,
102 HIST_MIN_TIME,
103 HIST_MAX_TIME);
104
105 histograms.add_hist<TH2D>(
106 "ToTLVsBar", "ToT_vs_bar_left", total_bar_num, 0.5, 0.5 + total_bar_num, BINSIZE_TOT, 0., HIST_MAX_TOT);
107 histograms.add_hist<TH2D>(
108 "ToTRVsBar", "ToT_vs_bar_right", total_bar_num, 0.5, 0.5 + total_bar_num, BINSIZE_TOT, 0., HIST_MAX_TOT);
109 }
110
111 void Map2CalTask::ExtraInit(FairRootManager* /*rootMan*/)
112 {
113 mappedData_.init();
114 trigMappedData_.init();
115 calData_.init(not IsHistDisabled());
116 set_pmt_num();
117 set_ct_freq();
118 }
119
121 {
122 if (plane_num_ = base_par_->get_num_of_planes(); plane_num_ == 0)
123 {
124 R3BLOG(warn, "plane number obtained from the calibration parameters is 0!");
125 }
126 else
127 {
129 R3BLOG(info, fmt::format("total number of PMTs set to be {}", total_pmt_nums_).c_str());
130 }
131
132 if (total_pmt_nums_ == 0)
133 {
134 throw R3B::runtime_error("The number of PMTs cannot be zero!");
135 }
136 }
137
139 {
140 if (auto ct_freq = calibrationPar_->GetSlowClockFrequency(); ct_freq == 0)
141 {
142 R3BLOG(warn, "Coarse time frequency obtained from parameteers is 0! Use default value.");
144 }
145 else
146 {
147 coarse_time_frequency_ = ct_freq;
149 R3BLOG(info, fmt::format("Coarse time frequency set to be {} MHz", ct_freq).c_str());
150 }
151 }
152
153 auto Map2CalTask::CheckConditions() const -> bool
154 {
155
156 auto signal_size = mappedData_.size();
157 R3BLOG(debug2,
158 fmt::format("Minimal signal size: {}. Current signal size: {}. Number of PMTs: {}",
160 signal_size,
162 if (signal_size < signal_min_size_)
163 {
164 R3BLOG(debug2,
165 fmt::format("condition of the minimal size is not met with current paddle signal size. Skip the "
166 "current event."));
167 return false;
168 }
169 return is_pulse_mode_ ? signal_size > total_pmt_nums_ : signal_size < total_pmt_nums_ / 2;
170 }
171
173
175 SingleEdgeSignal signal,
176 FTType ftType,
177 unsigned int module_num) const -> ValueError<double>
178 {
179 const auto& modulePar = calPar->GetParamAt(module_num);
180 const auto fineTime = modulePar.GetFineTime(ftType, signal.fine);
181
182 // signal time in ns
183 const auto coarseTime = 1000. / coarse_time_frequency_ * static_cast<float>(signal.coarse + 1);
184 return { coarseTime - fineTime.value, fineTime.error };
185 }
186
187 auto Map2CalTask::get_trigger_time(unsigned int module_num, Side side) const -> ValueError<double>
188 {
189
190 const auto& triggerMap = base_par_->get_trig_id_map();
191 auto triggerIDPair = triggerMap.find(module_num);
192 if (triggerIDPair == triggerMap.end())
193 {
194 const auto eventNum = GetEventHeader()->GetEventno();
195 R3BLOG(error,
196 fmt::format("Can't find the trigger ID for the module ID {} in the event {}", module_num, eventNum));
197 return { -1., 0 };
198 }
199 const auto triggerID = (side == Side::left) ? triggerIDPair->second.first : triggerIDPair->second.second;
200 const auto trigData = std::find_if(trigMappedData_.begin(),
201 trigMappedData_.end(),
202 [&triggerID](const auto& ele) { return ele.first == triggerID; });
203 if (trigData == trigMappedData_.end())
204 {
205 const auto eventNum = GetEventHeader()->GetEventno();
206 R3BLOG(error,
207 fmt::format("No such trigger ID {} in mappedTrigData for moduleNum {} in the event {}!",
208 triggerID,
209 module_num,
210 eventNum));
211 R3BLOG(error, fmt::format("Available trigIDs: {}", fmt::join(trigMappedData_ | ranges::views::keys, ", ")));
212 return { -1., 0 };
213 }
214 return convert_to_real_time(calibrationTrigPar_, trigData->second.signal, FTType::trigger, trigData->first);
215 }
216
218 unsigned int module_num,
219 R3B::Side module_side) const -> ValueError<double>
220 {
221 const auto leadFType = (module_side == Side::left) ? FTType::leftleading : FTType::rightleading;
222 const auto trailFType = (module_side == Side::left) ? FTType::lefttrailing : FTType::righttrailing;
223 const auto leadingT = convert_to_real_time(calibrationPar_, pmtSignal.leading, leadFType, module_num);
224 const auto trailingT = convert_to_real_time(calibrationPar_, pmtSignal.trailing, trailFType, module_num);
225 const auto time_over_thres = trailingT - leadingT;
226 R3BLOG(debug3, fmt::format("leading :{} trailing: {}", leadingT.value, trailingT.value));
227 return (time_over_thres.value > 0) ? time_over_thres : time_over_thres + max_coarse_time_;
228 }
229
231 {
232 // trigger time should come after signal time
233 if (calSignal.leading_time.value > calSignal.trigger_time.value)
234 {
235 calSignal.leading_time -= max_coarse_time_;
236 }
237 }
238
240 R3B::Side side,
241 unsigned int module_num) const -> CalDataSignal
242 {
243 auto calDataSignal = CalDataSignal{};
244 const auto ftType = (side == R3B::Side::left) ? FTType::leftleading : FTType::rightleading;
245 calDataSignal.time_over_threshold = get_tot(double_edge_signal, module_num, side);
246 const auto walk_correction =
247 (is_walk_enabled_) ? GetWalkCorrection(calDataSignal.time_over_threshold.value) : 0.;
248 calDataSignal.leading_time =
249 convert_to_real_time(calibrationPar_, double_edge_signal.leading, ftType, module_num) + walk_correction;
250 calDataSignal.trigger_time = get_trigger_time(module_num, side);
251 overflow_correct(calDataSignal);
252 R3BLOG(debug, fmt::format("Adding a new cal signal: {}", calDataSignal));
253 return calDataSignal;
254 }
255
257 unsigned int module_num,
258 R3B::Side side) const -> std::vector<CalDataSignal>
259 {
260 const auto& signals = (side == Side::left) ? map_bar_signals.left : map_bar_signals.right;
261 auto calSignals = std::vector<CalDataSignal>{};
262 calSignals.reserve(signals.size());
263 std::transform(signals.begin(),
264 signals.end(),
265 std::back_inserter(calSignals),
266 [module_num, side, this](const auto& dESignal)
267 { return doubleEdgeSignal_to_calSignal(dESignal, side, module_num); });
268 return calSignals;
269 }
270
272 {
275
276 if (not IsHistDisabled())
277 {
280 }
281 }
282
284 {
285 R3BLOG(debug2, fmt::format("mapped Data size: {}", mappedData_.size()));
286 for (const auto& planeSignals : mappedData_)
287 {
288 const auto planeNum = planeSignals.plane_num;
289 for (const auto& [barNum, barSignals] : planeSignals.bars)
290 {
291 if (not IsHistDisabled())
292 {
293 GetHistMonitor().get("module_num")->Fill(barNum);
294 }
295 const auto module_num = Neuland_PlaneBar2ModuleNum(planeNum, barNum);
296 auto& cal = calData_.get().emplace_back(module_num);
297 fill_cal_data(cal, barSignals);
298 }
299 }
300 }
301
303 {
304 auto& histograms = GetHistMonitor();
305 const auto module_num = cal.module_num;
306 const auto& calSignals = (side == Side::left) ? cal.left : cal.right;
307
308 for (const auto& signal : calSignals)
309 {
310 const auto lTime = signal.leading_time.value;
311 const auto time_over_thres = signal.time_over_threshold.value;
312 const auto triggerTime = signal.trigger_time.value;
313 histograms.get("ToT")->Fill(time_over_thres);
314 histograms.get("LeadingTime")->Fill(lTime);
315 histograms.get("TriggerTime")->Fill(triggerTime);
316 histograms.get("Leading_minus_trigger")->Fill(lTime - triggerTime);
317 if (side == Side::left)
318 {
319 histograms.get("TimeLVsBar")->Fill(module_num, lTime - triggerTime);
320 histograms.get("ToTLVsBar")->Fill(module_num, time_over_thres);
321 histograms.get("Bar_hitNum_l")->Fill(module_num);
322 }
323 else
324 {
325 histograms.get("TimeRVsBar")->Fill(module_num, lTime - triggerTime);
326 histograms.get("ToTRVsBar")->Fill(module_num, time_over_thres);
327 histograms.get("Bar_hitNum_r")->Fill(module_num);
328 }
329 }
330 }
331
333} // namespace R3B::Neuland
#define R3BLOG(severity, x)
Definition R3BLogger.h:32
R3B::Map2CalPar R3BTCalPar2
auto add_hist(std::unique_ptr< TH1 > hist) -> TH1 *
void HistogramInit(DataMonitor &histograms) override
InputMapConnector< unsigned int, PaddleTamexTrigMappedData > trigMappedData_
auto get_tot(const DoubleEdgeSignal &pmtSignal, unsigned int module_num, R3B::Side module_side) const -> ValueError< double >
OutputVectorConnector< BarCalData > calData_
auto get_trigger_time(unsigned int module_num, Side side) const -> ValueError< double >
auto mapBarSignal_to_calSignals(const MapBarSignal &map_bar_signals, unsigned int module_num, R3B::Side side) const -> std::vector< CalDataSignal >
auto convert_to_real_time(R3BTCalPar2 *calPar, SingleEdgeSignal signal, FTType ftType, unsigned int module_num) const -> ValueError< double >
auto CheckConditions() const -> bool override
InputVectorConnector< PaddleTamexMappedData > mappedData_
void overflow_correct(R3B::Neuland::CalDataSignal &calSignal) const
void ExtraInit(FairRootManager *rootMan) override
void fill_cal_data(BarCalData &cal, const MapBarSignal &signals)
auto doubleEdgeSignal_to_calSignal(const DoubleEdgeSignal &double_edge_signal, R3B::Side side, unsigned int module_num) const -> CalDataSignal
void histogram_monitor(const BarCalData &cal, Side side)
Simulation of NeuLAND Bar/Paddle.
constexpr auto COARSE_TIME_CLOCK_FREQUENCY_MHZ
constexpr auto BarsPerPlane
constexpr auto Neuland_PlaneBar2ModuleNum(unsigned int planeNum, unsigned int barNum) -> unsigned int
std::vector< CalDataSignal > left
std::vector< CalDataSignal > right
ValueError< double > trigger_time
ValueError< double > leading_time