R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandTamexReader2.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
14#include "FairRootManager.h"
15
16#include "R3BException.h"
17#include "R3BLogger.h"
20#include "R3BReader.h"
21#include "ext_data_struct_info.hh"
22#include "ext_h101_raw_nnp_tamex.h"
23#include <R3BShared.h>
24#include <Rtypes.h>
25#include <TH1.h>
26#include <array>
27
28// TODO: C++20 std::span
29#include <boost/regex/v5/regex.hpp>
30#include <boost/regex/v5/regex_fwd.hpp>
31#include <boost/regex/v5/regex_match.hpp>
32#include <boost/regex/v5/regex_search.hpp>
33#include <cstddef>
34#include <fairlogger/Logger.h>
35#include <fmt/core.h>
36#include <gsl/span>
37// TODO: C++20 std::range
38#include <map>
39#include <range/v3/algorithm/equal.hpp>
40#include <range/v3/view/iota.hpp>
41#include <range/v3/view/take.hpp>
42#include <range/v3/view/zip.hpp>
43#include <string>
44#include <string_view>
45#include <utility>
46
51
52using gsl::span;
53namespace views = ranges::views;
54using R3B::Side;
56
57namespace
58{
59 constexpr size_t PRINT_ERROR_MAX = 100000;
60 constexpr auto BarsPerPlane = 50;
61 const auto errorStrings = std::map<Errors, std::string>{
62 { Errors::module_size, "Incompatible module size" },
63 { Errors::data_size, "Incompatible signal size" },
64 { Errors::divider, "Incompatible signal dividers" },
65 { Errors::indices, "Incompatible module IDs" },
66 };
67
68 const auto SIDES = std::array<Side, 2>{ Side::left, Side::right };
69 template <typename BinaryOperation, typename Item0, typename Item1>
70 constexpr auto CheckAmong(BinaryOperation optn, const Item0& item0, const Item1& item1) -> bool
71 {
72 return optn(item0, item1);
73 }
74
75 template <typename BinaryOperation, typename Item0, typename Item1, typename... Items>
76 constexpr auto CheckAmong(BinaryOperation optn, const Item0& item0, const Item1& item1, const Items&... items)
77 -> bool
78 {
79 return CheckAmong(optn, item0, item1) && CheckAmong(optn, item1, items...);
80 }
81
82 template <typename UnaryOperation, typename... Items>
83 inline auto CheckAny(UnaryOperation optn, const Items&... items)
84 {
85 return (optn(items) || ...);
86 }
87
88 template <typename UnaryFunction, typename... Items>
89 inline auto Zip_from(UnaryFunction optn, const Items&... items)
90 {
91 return ranges::zip_view(optn(items)...);
92 }
93
94 template <typename Signals, typename ModuleNum_divider, typename UnaryFunction>
95 inline void FillData(const Signals& signals_view,
96 const ModuleNum_divider& moduleNum_divider_view,
97 UnaryFunction optn)
98 {
99 auto range_begin = signals_view.begin();
100 auto range_iter = signals_view.begin();
101 for (const auto& [module_num, divider] : moduleNum_divider_view)
102 {
103 for (; !(range_iter == signals_view.end() || range_iter == (range_begin + divider)); ++range_iter)
104 {
105 optn(*range_iter, module_num);
106 }
107 }
108 }
109} // namespace
110
112 : R3BReader("R3BNeulandTamexReader2")
113 , max_limit_(PRINT_ERROR_MAX)
114 , fOffset(offset)
115 , numPlanes_(static_cast<int>(R3B::GetSize(inputData_->NN_P)))
116 , inputData_{ data }
119{
120 error_log_.insert_or_assign(Errors::module_size, 0);
121 error_log_.insert_or_assign(Errors::data_size, 0);
122 error_log_.insert_or_assign(Errors::indices, 0);
123 error_log_.insert_or_assign(Errors::divider, 0);
124}
125
126bool R3BNeulandTamexReader2::Init(ext_data_struct_info* a_struct_info) // NOLINT
127{
128 auto status_ok = 1;
129 R3BLOG(info, "");
130 EXT_STR_h101_raw_nnp_tamex_ITEMS_INFO(status_ok, *a_struct_info, fOffset, EXT_STR_h101_raw_nnp_tamex, 0); // NOLINT
131
132 if (status_ok == 0)
133 {
134 R3BLOG(error, "Failed to setup structure information.");
135 return false;
136 }
137
138 R3BLOG(info, "Number of planes " << numPlanes_);
139
140 // Register output arrays in tree
141 auto* rootMan = FairRootManager::Instance();
142 rootMan->RegisterAny("NeulandMappedData", mappedDataPtr_, !is_online_);
143 if (is_triggered_)
144 {
145 rootMan->RegisterAny("NeulandTrigMappedData", mappedTrigDataPtr_, !is_online_);
146 }
147 Reset();
148 if (not is_online_)
149 {
151 }
152 return true;
153}
154
156{
157 const auto BarBinSize = numPlanes_ * 50;
158 histograms_.add_hist<TH1I>(
159 "bar_FT_l", "number of hits on left PMTS for each bar", BarBinSize, 0.5, 0.5 + BarBinSize);
160 histograms_.add_hist<TH1I>(
161 "bar_FT_r", "number of hits on right PMTS for each bar", BarBinSize, 0.5, 0.5 + BarBinSize);
162}
163
164namespace
165{
166 template <typename ErrorLog, typename... Items>
167 auto CheckCondition(ErrorLog& log, const Items&... items) -> bool
168 {
169 auto ApplyCriterium = [&log, &items...](Errors error, auto criterium) -> bool
170 {
171 if (!CheckAmong(criterium, items...))
172 {
173 ++log.at(error);
174 return false;
175 }
176 return true;
177 };
178 return ApplyCriterium(Errors::module_size,
179 [](const auto& left, const auto& right) { return left.BM == right.BM; }) &&
180 ApplyCriterium(Errors::data_size,
181 [](const auto& left, const auto& right) { return left.B == right.B; }) &&
182 ApplyCriterium(
183 Errors::indices,
184 [](const auto& left, const auto& right)
185 { return ranges::equal(views::take(left.BMI, left.BM), views::take(right.BMI, right.BM)); }) &&
186 ApplyCriterium(
187 Errors::divider,
188 [](const auto& left, const auto& right)
189 { return ranges::equal(views::take(left.BME, left.BM), views::take(right.BME, right.BM)); });
190 }
191} // namespace
192
193template <typename ViewType>
194auto R3BNeulandTamexReader2::extract_plane_signals(const ViewType& signalsPlane, int planeNum)
195{
196 auto planeSignals = R3B::PaddleTamexMappedData{};
197 const auto signals_sides_view =
198 ranges::zip_view(signalsPlane.tcl_T, signalsPlane.tfl_T, signalsPlane.tct_T, signalsPlane.tft_T, SIDES);
199 for (const auto& [coarse_leading, fine_leading, coarse_trailing, fine_trailing, side] : signals_sides_view)
200 {
201 // pass if no data
202 if (CheckAny([](const auto& item) { return item.BM == 0; },
203 coarse_leading,
204 fine_leading,
205 coarse_trailing,
206 fine_trailing))
207 {
208 R3BLOG(debug3, fmt::format("No signals at the current event from plane {}.", planeNum));
209 continue;
210 }
211
212 if (!CheckCondition(error_log_, coarse_leading, fine_leading, coarse_trailing, fine_trailing))
213 {
215 }
216
217 const auto module_size = coarse_leading.BM;
218 const auto signal_size = coarse_leading.B;
219 R3BLOG(debug, fmt::format("Signals at the current event from plane {}: {}", planeNum, signal_size));
220 const auto barNum_divider_view = ranges::zip_view(span(coarse_leading.BMI), span(coarse_leading.BME));
221
222 const auto signals_view = Zip_from([](const auto& item) { return span(item.Bv); },
223 coarse_leading,
224 fine_leading,
225 coarse_trailing,
226 fine_trailing);
227 // TODO C++20: planeID and side can be directly captured by the lambda. Before C++20, capturing struture
228 // bindings is undefined behaviour
229 auto side_temp = side;
230 FillData(
231 signals_view | ranges::views::take(signal_size),
232 barNum_divider_view | ranges::views::take(module_size),
233 [&](const auto& signals, const auto& barNum)
234 {
235 const auto& [cl_time, fl_time, ct_time, ft_time] = signals;
236 auto doubleEdgeSignal = R3B::DoubleEdgeSignal{};
237 doubleEdgeSignal.leading.fine = fl_time;
238 doubleEdgeSignal.leading.coarse = cl_time;
239 doubleEdgeSignal.trailing.fine = ft_time;
240 doubleEdgeSignal.trailing.coarse = ct_time;
241
242 planeSignals.push_back(side_temp, barNum, doubleEdgeSignal);
243 R3BLOG(
244 debug3,
245 fmt::format(
246 "Writing a doube edge signal: barNum: {}, fl_time: {}, cl_time: {}, ft_time: {}, ct_time: {}",
247 barNum,
248 fl_time,
249 cl_time,
250 ft_time,
251 ct_time));
252
253 if (is_online_)
254 {
255 return;
256 }
257
258 if (side_temp == R3B::Side::left)
259 {
260 histograms_.get("bar_FT_l")->Fill(barNum + (planeNum - 1) * BarsPerPlane);
261 }
262 else
263 {
264 histograms_.get("bar_FT_r")->Fill(barNum + (planeNum - 1) * BarsPerPlane);
265 }
266 });
267 }
268 return planeSignals;
269}
270
272{
273 const auto signalsAllPlanes = span(inputData->NN_P);
274
275 mappedData_.clear();
276 for (const auto& [signalsPlane, plane_num] : ranges::zip_view(signalsAllPlanes, ranges::iota_view(1)))
277 {
278 auto planeSignals = extract_plane_signals(signalsPlane, plane_num);
279 if (!planeSignals.empty())
280 {
281 if (not is_online_)
282 {
283 histogram_action(planeSignals);
284 }
285 mappedData_.try_emplace(plane_num, std::move(planeSignals));
286 }
287 }
288
289 return true;
290}
291
293{
294 for (const auto& [hist_name, action] : hist_actions_)
295 {
296 action(mappedData, histograms_.get(hist_name));
297 }
298}
299
300namespace
301{
302 template <typename ErrorLog, typename... Items>
303 auto Checkondition(ErrorLog& log, EXT_STR_h101_raw_nnp_tamex_onion* inputData)
304 {
305 auto ApplyCriterium = [&log](Errors error, bool criterium) -> bool
306 {
307 if (!criterium)
308 {
309 ++log.at(error);
310 return false;
311 }
312 return true;
313 };
314 return ApplyCriterium(Errors::module_size, inputData->NN_TRIGCM == inputData->NN_TRIGFM) &&
315 ApplyCriterium(Errors::data_size, inputData->NN_TRIGC == inputData->NN_TRIGF) &&
316 ApplyCriterium(Errors::indices,
317 ranges::equal(views::take(inputData->NN_TRIGCMI, inputData->NN_TRIGCM),
318 views::take(inputData->NN_TRIGFMI, inputData->NN_TRIGFM))) &&
319 ApplyCriterium(Errors::divider,
320 ranges::equal(views::take(inputData->NN_TRIGCME, inputData->NN_TRIGCM),
321 views::take(inputData->NN_TRIGFME, inputData->NN_TRIGFM)));
322 }
323} // namespace
324
326{
327 if (!is_triggered_ || inputData->NN_TRIGCM == 0 || inputData->NN_TRIGFM == 0)
328 {
329 return true;
330 }
331
332 if (!Checkondition(error_log_, inputData))
333 {
334 return false;
335 }
336
337 const auto module_size = inputData->NN_TRIGCM;
338 const auto signal_size = inputData->NN_TRIGC;
339 auto const trigger_signals_view =
340 Zip_from([](const auto& item) { return span(item); }, inputData->NN_TRIGCv, inputData->NN_TRIGFv);
341 auto const moduleNum_divider_view = ranges::zip_view(span(inputData->NN_TRIGCMI), span(inputData->NN_TRIGCME));
342 FillData(ranges::views::take(trigger_signals_view, signal_size),
343 ranges::views::take(moduleNum_divider_view, module_size),
344 [&](const auto& signal_pack, const auto& module_num)
345 {
346 const auto& [coarse_time, fine_time] = signal_pack;
347 auto trigDatum = mappedTrigData_.emplace(module_num, R3BPaddleTamexTrigMappedData{}).first;
348 trigDatum->second.signal.fine = fine_time;
349 trigDatum->second.signal.coarse = coarse_time;
350 R3BLOG(debug3,
351 fmt::format("Writing a trig signal: module ID: {}, fine time: {}, coarse_time: {}",
352 module_num,
353 fine_time,
354 coarse_time));
355 // std::cout << "trig signal: " << module_num << "\t" << fine_time << "\t" << coarse_time << "\n";
356 });
357 return true;
358}
359
361{
363 {
364 return false;
365 }
366 if (fair::Logger::GetConsoleSeverity() <= fair::Severity::info && ++counter_ == max_limit_)
367 {
368 counter_ = 0;
369 print_error();
370 }
371 return true;
372}
373
374auto R3BNeulandTamexReader2::check_trigger_needed(std::string_view item_name) const -> bool
375{
376 // TODO: not the exact regex, but it suffices
377 const auto neuland_trigger_regex = boost::regex{ "^NN_TRIG(C|F)(v|M)?(I|E)?$" };
378 auto res = is_triggered_ && boost::regex_match(item_name.data(), neuland_trigger_regex); // NOLINT
379 return res;
380}
381
382auto R3BNeulandTamexReader2::check_bar_needed(std::string_view item_name) const -> bool
383{
384 const auto neuland_bar_regex = boost::regex{ R"(^NN_P(\d*)t.*$)" };
385 auto result = boost::cmatch{};
386 if (boost::regex_search(item_name.data(), result, neuland_bar_regex)) // NOLINT
387 {
388 const auto plane_num = std::stoi(result.str(1));
389 if (plane_num != 0 and plane_num <= numPlanes_)
390 {
391 return true;
392 }
393 }
394 return false;
395}
396
397auto R3BNeulandTamexReader2::MismappedItemRequired(std::string_view item_name) const -> bool
398{
399 if (numPlanes_ == 0)
400 {
401 throw R3B::logic_error(
402 "The plane number of NeuLAND is not set. Please use R3BNeulandTamexReader2::SetMaxNbPalnes to set it.");
403 }
404
405 return check_trigger_needed(item_name) || check_bar_needed(item_name);
406}
407
409{
410 mappedData_.clear();
411 mappedTrigData_.clear();
412}
413
415{
416 fmt::print("----------Tamex Reader Error Summary:--------------\n");
417 fmt::print("{0:^30}|{1:^10}\n", "Causes", "Counts");
418 for (const auto& [key, counts] : error_log_)
419 {
420 fmt::print("{0:<30}|{1:^10}\n", errorStrings.at(key), counts);
421 }
422 fmt::print("---------------------------------------------------\n");
423}
424
#define R3BLOG(severity, x)
Definition R3BLogger.h:32
ClassImp(R3BNeulandTamexReader2)
R3BNeulandTamexReader2::Errors Errors
struct EXT_STR_h101_raw_nnp_tamex_onion_t EXT_STR_h101_raw_nnp_tamex_onion
struct EXT_STR_h101_raw_nnp_tamex_t EXT_STR_h101_raw_nnp_tamex
R3B::PaddleTamexMappedData R3BPaddleTamexMappedData2
R3B::PaddleTamexTrigMappedData R3BPaddleTamexTrigMappedData
auto R3BRead() -> bool override
std::map< std::string, std::function< void(const R3B::PaddleTamexMappedData &, TH1 *)> > hist_actions_
auto ReadSignals(EXT_STR_h101_raw_nnp_tamex_onion *inputData) -> bool
auto check_bar_needed(std::string_view item_name) const -> bool
R3BNeulandTamexReader2(EXT_STR_h101_raw_nnp_tamex_onion *, size_t)
auto Init(ext_data_struct_info *) -> bool override
auto extract_plane_signals(const ViewType &signalsPlane, int planeNum)
TrigMappedDataVector * mappedTrigDataPtr_
auto ReadTriggerSignals(EXT_STR_h101_raw_nnp_tamex_onion *inputData) -> bool
EXT_STR_h101_raw_nnp_tamex_onion * inputData_
auto check_trigger_needed(std::string_view item_name) const -> bool
std::map< Errors, uint > error_log_
TrigMappedDataVector mappedTrigData_
void histogram_action(const R3B::PaddleTamexMappedData &mappedData)
MappedDataVector * mappedDataPtr_
auto MismappedItemRequired(std::string_view item_name) const -> bool override
R3BReader(TString const &)
Definition R3BReader.cxx:16