R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandMillepede.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 "R3BNeulandMillepede.h"
15#include "ParResultReader.h"
16#include "R3BDataMonitor.h"
17#include "R3BLogger.h"
18#include "R3BNeulandCalData2.h"
21#include <R3BException.h>
23#include <R3BNeulandCommon.h>
24#include <SteerWriter.h>
25
26#include <TGraphErrors.h>
27#include <TH1.h>
28#include <algorithm>
29#include <array>
30#include <cstdlib>
31#include <fmt/core.h>
32#include <iterator>
33#include <memory>
34#include <numeric>
35#include <optional>
36#include <range/v3/algorithm/all_of.hpp>
37#include <range/v3/iterator/operations.hpp>
38#include <range/v3/numeric/accumulate.hpp>
39#include <range/v3/view/all.hpp>
40#include <range/v3/view/filter.hpp>
41#include <range/v3/view/sliding.hpp>
42#include <range/v3/view/transform.hpp>
43#include <stdexcept>
44#include <string>
45#include <utility>
46#include <vector>
47
48namespace rng = ranges;
49
50constexpr auto DEFAULT_RES_FILENAME = "millepede.res";
51constexpr auto SCALE_FACTOR = 10.F;
52// constexpr auto REFERENCE_BAR_NUM = 25;
53constexpr auto MILLE_BUFFER_SIZE = std::size_t{ 100000 };
54constexpr auto DEFAULT_T_ERROR = 2; // ns
55
56namespace
57{
58 void calculate_time_offset(R3B::Neuland::Cal2HitPar& cal_to_hit_par)
59 {
60 auto& module_pars = cal_to_hit_par.GetListOfModuleParRef();
61 for (auto& [module_num, module_par] : module_pars)
62 {
63 if (module_par.effective_speed.value != 0)
64 {
65 module_par.t_diff = module_par.t_diff / module_par.effective_speed;
66 }
67 }
68 }
69
70 void change_time_offset(R3B::Neuland::Cal2HitPar& cal_to_hit_par)
71 {
72 auto& module_pars = cal_to_hit_par.GetListOfModuleParRef();
73 for (auto& [module_num, module_par] : module_pars)
74 {
75 module_par.t_diff = module_par.t_diff * module_par.effective_speed;
76 }
77 }
78} // namespace
79
81{
83 {
84 cal_to_hit_par_ = GetTask()->GetCal2HitPar();
85
87 pede_launcher_.set_steer_filename(pede_steer_filename_);
88 pede_launcher_.set_parameter_filename(parameter_filename_);
89 if (const auto* r3b_dir = std::getenv("R3BROOTPATH"); r3b_dir != nullptr)
90 {
91 pede_launcher_.set_binary_dir(fmt::format("{}/bin", r3b_dir));
92 }
93 else
94 {
96 "Environment variable R3BROOTPATH is not defined! Did you forget to source the \"config.sh\" file?");
97 }
99 data_preprocessor_ = std::make_unique<MilleDataProcessor>(GetModuleSize());
100 data_preprocessor_->set_p_value_cut(p_value_cut_);
101
104 }
105
106 // output: module_num & global label
107 inline auto MillepedeEngine::to_module_num_label(int par_num) -> std::pair<int, GlobalLabel>
108 {
109 const auto num_of_module = GetModuleSize();
110 auto res = std::pair<int, GlobalLabel>{};
111 const auto factor = (par_num - 1) / num_of_module;
112 res.first = (par_num - 1) % num_of_module + 1;
113
114 switch (factor)
115 {
116 // case 0:
117 // res.second = GlobalLabel::tsync;
118 // break;
119 // case 1:
120 // res.second = GlobalLabel::offset_effective_c;
121 // break;
122 // case 2:
123 // res.second = GlobalLabel::effective_c;
124 // break;
125 case 0:
127 break;
128 case 1:
129 res.second = GlobalLabel::effective_c;
130 break;
131 default:
132 throw R3B::logic_error(fmt::format("An error occured with unrecognized global par id: {}", par_num));
133 }
134
135 return res;
136 }
137
138 inline auto MillepedeEngine::get_global_label_id(int module_num, GlobalLabel label) -> int
139 {
140 const auto num_of_module = GetModuleSize();
141 switch (label)
142 {
143 // case GlobalLabel::tsync:
144 // return module_num;
145 // case GlobalLabel::offset_effective_c:
146 // return module_num + num_of_module;
147 // case GlobalLabel::effective_c:
148 // return module_num + (2 * num_of_module);
150 return module_num;
152 return module_num + num_of_module;
153 default:
154 throw R3B::logic_error("An error occured with unrecognized global tag");
155 }
156 }
157
159 Neuland::Cal2HitPar& cal_to_hit_par)
160 {
161 change_time_offset(cal_to_hit_par);
162 const auto& pars = result.get_pars();
163 for (const auto& [par_id, par] : pars)
164 {
165 const auto [module_num, global_label] = to_module_num_label(par_id);
166 auto& module_pars = cal_to_hit_par.GetListOfModuleParRef();
167
168 auto& par_ref = module_pars.emplace(module_num, HitModulePar{}).first->second;
169 switch (global_label)
170 {
172 par_ref.t_sync.value = par.value * SCALE_FACTOR;
173 par_ref.t_sync.error = par.error * SCALE_FACTOR;
174 break;
176 // The value here is the product of tDiff and effectiveSped. Real tDiff will be calculated later
177 par_ref.t_diff.value += par.value * SCALE_FACTOR;
178 par_ref.t_diff.error += par.error * SCALE_FACTOR;
179 break;
181 par_ref.effective_speed.value += par.value;
182 par_ref.effective_speed.error += par.error;
183 break;
184 default:
185 throw std::runtime_error("An error occured with unrecognized global tag");
186 }
187 }
188
189 calculate_time_offset(cal_to_hit_par);
190 }
191
192 auto MillepedeEngine::set_minimum_values(const std::vector<R3B::Neuland::BarCalData>& signals) -> bool
193 {
194 // make sure only one hit exists in one bar
195 auto filtered_signals = rng::filter_view(
196 signals | rng::views::all,
197 [](const auto& bar_signal) { return bar_signal.left.size() == 1 and bar_signal.right.size() == 1; });
198 if (filtered_signals.empty())
199 {
200 return false;
201 }
202
203 if (not average_t_sum_.has_value())
204 {
205 auto t_sum_view = filtered_signals | rng::views::transform(
206 [](const auto& bar_signal)
207 {
208 const auto& left_signal = bar_signal.left.front();
209 const auto& right_signal = bar_signal.right.front();
210 return (left_signal.leading_time - left_signal.trigger_time +
211 right_signal.leading_time - right_signal.trigger_time)
212 .value;
213 });
214 auto sum = rng::accumulate(t_sum_view, 0.F);
215 average_t_sum_ = sum / static_cast<float>(rng::distance(t_sum_view.begin(), t_sum_view.end()));
216 R3BLOG(info, fmt::format("Average t_sum is calculated to be {}", average_t_sum_.value()));
217 }
218 return true;
219 }
220
221 auto MillepedeEngine::SignalFilter(const std::vector<BarCalData>& signals) -> bool
222 {
223 // select out rays with few hits
224 if (signals.size() < minimum_hit_)
225 {
226 return false;
227 }
228
229 // select out vertical cosmic rays
230 if (rng::all_of(signals |
231 rng::views::transform([](const auto& bar_signal)
232 { return ModuleID2PlaneID(bar_signal.module_num - 1); }) |
233 rng::views::sliding(2),
234 [](const auto& pair) { return pair.front() == pair.back(); }))
235 {
236 return false;
237 }
238
239 if (not set_minimum_values(signals))
240 {
241 return false;
242 }
243
244 return true;
245 }
246
248 {
249 buffer_clear();
250 const auto module_num = static_cast<int>(signal.module_num);
251 const auto pos_z = ModuleNum2ZPos<float>(static_cast<int>(module_num));
252
253 auto init_effective_c = cal_to_hit_par_->GetModuleParAt(module_num).effective_speed.value;
254
255 const auto& left_signal = signal.left;
256 const auto& right_signal = signal.right;
257 const auto t_sum = (left_signal.leading_time - left_signal.trigger_time) +
258 (right_signal.leading_time - right_signal.trigger_time) - average_t_sum_.value_or(0.F);
259
260 input_data_buffer_.measurement =
261 static_cast<float>((t_sum.value / SCALE_FACTOR / 2.F) - (BarLength / SCALE_FACTOR / init_effective_c));
262 input_data_buffer_.sigma = static_cast<float>(t_sum.error / SCALE_FACTOR / 2. * error_scale_factor_);
263 // input_data_buffer_.sigma = static_cast<float>(DEFAULT_MEAS_ERROR);
264 const auto local_derivs_t = std::array{ 0.F, 0.F, pos_z / SCALE_FACTOR, 0.F, 0.F, 1.F };
265 std::copy(local_derivs_t.begin(), local_derivs_t.end(), std::back_inserter(input_data_buffer_.locals));
266 input_data_buffer_.globals.emplace_back(get_global_label_id(module_num, GlobalLabel::tsync), 1.F);
267 input_data_buffer_.globals.emplace_back(get_global_label_id(module_num, GlobalLabel::effective_c),
268 -BarLength / SCALE_FACTOR / 2.F / init_effective_c / init_effective_c);
269
271 }
272
273 void MillepedeEngine::add_spacial_local_constraint(int plane_id, const std::vector<MilleCalData>& plane_signals)
274 {
275 buffer_clear();
276 const auto pos_z = PlaneID2ZPos<float>(plane_id);
277 const auto is_horizontal = IsPlaneIDHorizontal(plane_id);
278 // const auto pos_bar_vert_disp = GetBarVerticalDisplacement(module_num);
279 const auto pos_bar_vert_disp = std::accumulate(plane_signals.begin(),
280 plane_signals.end(),
281 0.,
282 [](double sum, const auto& signal) {
283 return sum + GetBarVerticalDisplacement(signal.module_num);
284 }) /
285 static_cast<double>(plane_signals.size());
286 const auto local_derivs = is_horizontal ? std::array{ 0.F, pos_z / SCALE_FACTOR, 0.F, 1.F }
287 : std::array{ pos_z / SCALE_FACTOR, 0.F, 1.F, 0.F };
288 // const auto local_derivs = is_horizontal ? std::array{ 0.F, pos_z / SCALE_FACTOR, 0.F, 0.F, 1.F, 0.F }
289 // : std::array{ pos_z / SCALE_FACTOR, 0.F, 0.F, 1.F, 0.F, 0.F };
290
291 input_data_buffer_.measurement = static_cast<float>(pos_bar_vert_disp / SCALE_FACTOR);
292 input_data_buffer_.sigma = static_cast<float>(BarSize_XY / SQRT_12 / SCALE_FACTOR * error_scale_factor_);
293
294 // if (not is_horizontal)
295 // {
296 // fmt::println("c_value_1.append({})\na_value_1.append({})\nb_value_1.append({})",
297 // pos_bar_vert_disp / SCALE_FACTOR,
298 // pos_z / SCALE_FACTOR,
299 // 1.);
300 // }
301 std::copy(local_derivs.begin(), local_derivs.end(), std::back_inserter(input_data_buffer_.locals));
303 }
304
306 {
307 buffer_clear();
308 const auto module_num = static_cast<int>(signal.module_num);
309 const auto plane_id = ModuleID2PlaneID(static_cast<int>(module_num) - 1);
310 const auto is_horizontal = IsPlaneIDHorizontal(plane_id);
311 const auto& module_par = cal_to_hit_par_->GetModuleParAt(module_num);
312 auto init_effective_c = module_par.effective_speed.value;
313 auto init_t_offset = module_par.t_diff.value;
314 const auto pos_z = static_cast<float>(PlaneID2ZPos(plane_id));
315
316 const auto& left_signal = signal.left;
317 const auto& right_signal = signal.right;
318 const auto t_diff = (right_signal.leading_time - right_signal.trigger_time) -
319 (left_signal.leading_time - left_signal.trigger_time);
320
321 const auto t_error = t_diff.error == 0 ? DEFAULT_T_ERROR : t_diff.error;
322 input_data_buffer_.measurement = static_cast<float>((init_effective_c * init_t_offset / 2. / SCALE_FACTOR) -
323 (init_effective_c * t_diff.value / SCALE_FACTOR / 2.));
324 input_data_buffer_.sigma =
325 static_cast<float>(t_error / SCALE_FACTOR / 2. * std::abs(init_effective_c) * error_scale_factor_);
326 const auto local_derivs = is_horizontal ? std::array{ pos_z / SCALE_FACTOR, 0.F, 1.F, 0.F }
327 : std::array{ 0.F, pos_z / SCALE_FACTOR, 0.F, 1.F };
328 // const auto local_derivs = is_horizontal ? std::array{ pos_z / SCALE_FACTOR, 0.F, 0.F, 1.F, 0.F, 0.F }
329 // : std::array{ 0.F, pos_z / SCALE_FACTOR, 0.F, 0.F, 1.F, 0.F };
330 std::copy(local_derivs.begin(), local_derivs.end(), std::back_inserter(input_data_buffer_.locals));
331 // fmt::println("Adding global: {}", get_global_label_id(module_num, GlobalLabel::offset_effective_c));
333 -0.5F);
334 input_data_buffer_.globals.emplace_back(get_global_label_id(module_num, GlobalLabel::effective_c),
335 static_cast<float>(t_diff.value / SCALE_FACTOR / 2.));
336
337 // if (is_horizontal)
338 // {
339 // const auto c_val = (0.5 * (module_par.t_diff * module_par.effective_speed).value / SCALE_FACTOR) -
340 // (module_par.effective_speed.value * t_diff.value / SCALE_FACTOR / 2.);
341 // fmt::println(
342 // "c_value_2.append({})\na_value_2.append({})\nb_value_2.append({})", c_val, pos_z / SCALE_FACTOR, 1.);
343 // }
345 R3BLOG(
346 debug,
347 fmt::format(
348 "Writting Mille data to binary file with meas = {} and z = {}", input_data_buffer_.measurement, pos_z));
349 }
350
351 auto MillepedeEngine::select_t_diff_signal(const std::vector<MilleCalData>& plane_data)
352 {
353 if (plane_data.empty())
354 {
355 return plane_data.end();
356 }
357 auto calculate_residual = [this](const MilleCalData& signal) -> double
358 {
359 const auto t_diff = (signal.right.leading_time - signal.right.trigger_time) -
360 (signal.left.leading_time - signal.left.trigger_time);
361 const auto& module_par = cal_to_hit_par_->GetModuleParAt(signal.module_num);
362 const auto position = (-t_diff + module_par.t_diff) / 2 * module_par.effective_speed;
363 const auto res =
364 data_preprocessor_->calculate_residual(position.value, static_cast<int>(signal.module_num));
365 return res;
366 };
367 auto iter = std::min_element(plane_data.begin(),
368 plane_data.end(),
369 [calculate_residual](const auto& first, const auto& second)
370 { return calculate_residual(first) < calculate_residual(second); });
371 if (iter != plane_data.end())
372 {
373 auto residual = calculate_residual(*iter);
374 hist_t_offset_residual_->Fill(residual);
375 if (residual > t_diff_residual_cut_)
376 {
377 return plane_data.end();
378 }
379 // fmt::println("selected signal residual: {}", residual);
380 }
381 return iter;
382 }
383
384 void MillepedeEngine::AddSignals(const std::vector<BarCalData>& signals)
385 {
386
387 // fmt::println("==================new event===============\n");
388 if (not data_preprocessor_->filter(signals))
389 {
390 return;
391 }
392 const auto& processed_data = data_preprocessor_->get_data();
393 for (const auto& [plane_id, plane_signals] :
394 processed_data |
395 rng::views::filter([](const auto& planeid_signals) { return not planeid_signals.second.empty(); }))
396 {
397 // fmt::println("\n--------------------\n");
398 auto iter = select_t_diff_signal(plane_signals);
399 if (iter == plane_signals.end())
400 {
401 continue;
402 }
403 // fmt::println("selected signal: {}", *iter);
404 // for (const auto& signal : plane_signals)
405 // {
406 // // add_signal_t_sum(signal);
407 add_signal_t_diff(*iter);
408 add_spacial_local_constraint(plane_id, plane_signals);
409 // }
410 }
411 }
412
414 {
415 binary_data_writer_.close();
416 R3BLOG(info, "Launching pede algorithm..");
417 pede_launcher_.launch();
418 pede_launcher_.end();
419
420 par_result_.read();
422 fill_data_to_figure(hit_par);
423 }
424
426 {
427 const auto& pars = hit_par.GetListOfModulePar();
428 for (const auto& [module_num, par] : pars)
429 {
430 graph_time_offset_->SetPoint(static_cast<int>(module_num), module_num, par.t_diff.value);
431 graph_time_offset_->SetPointError(static_cast<int>(module_num), 0., par.t_diff.error);
432 graph_time_sync_->SetPoint(static_cast<int>(module_num), module_num, par.t_sync.value);
433 graph_time_sync_->SetPointError(static_cast<int>(module_num), 0., par.t_sync.error);
434 graph_effective_c_->SetPoint(static_cast<int>(module_num), module_num, par.effective_speed.value);
435 graph_effective_c_->SetPointError(static_cast<int>(module_num), 0., par.effective_speed.error);
436 }
437 }
438
439 void MillepedeEngine::EndOfEvent(unsigned int /*event_num*/)
440 {
441 // TODO: could be an empty event
443 data_preprocessor_->reset();
444 }
445
447
449 {
450 const auto module_size = GetModuleSize();
451
452 graph_time_offset_ = histograms.add_graph("time_offset", std::make_unique<TGraphErrors>(module_size));
453 graph_time_offset_->SetTitle("Time offset vs BarNum");
454
455 graph_time_sync_ = histograms.add_graph("time_sync", std::make_unique<TGraphErrors>(module_size));
456 graph_time_sync_->SetTitle("Time sync vs BarNum");
457
458 graph_effective_c_ = histograms.add_graph("effective_c", std::make_unique<TGraphErrors>(module_size));
459 graph_effective_c_->SetTitle("Effective c vs BarNum");
460
461 static constexpr auto RESIDUAL_BIN_NUM = 500;
462 hist_t_offset_residual_ = histograms.add_hist<TH1D>(
463 "t_diff_residual", "Residual values of the positios calculated from t_dff", RESIDUAL_BIN_NUM, 0., 1000.);
464 }
465
467 {
468 input_data_buffer_.locals.clear();
469 input_data_buffer_.globals.clear();
470 input_data_buffer_.measurement = 0.F;
471 input_data_buffer_.sigma = 0.F;
472 }
473
475
477 {
478 average_t_sum_.reset();
479 buffer_clear();
480 }
481
483 {
484 if (cal_to_hit_par_ == nullptr)
485 {
486 throw R3B::runtime_error("Pointer to cal_to_hit_par is nullptr!");
487 }
488 }
489
491 {
492 auto steer_writer = SteerWriter{};
494 steer_writer.set_parameter_file(parameter_filename_);
495 steer_writer.set_data_filepath(input_data_filename_);
496 static constexpr auto NUMBER_OF_ITERARTION = 3.F;
497 static constexpr auto CONVERGENCE_RECOGNITION = 0.001F;
498 steer_writer.add_method(SteerWriter::Method::inversion,
499 std::make_pair(NUMBER_OF_ITERARTION, CONVERGENCE_RECOGNITION));
500 steer_writer.add_other_options(std::vector<std::string>{ "hugecut", "50000" });
501 steer_writer.add_other_options(std::vector<std::string>{ "outlierdownweighting", "4" });
502
503 // const auto module_size = GetModuleSize();
504 // for (int module_num{ 1 }; module_num <= module_size; ++module_num)
505 // {
506 // const auto& module_par = cal_to_hit_par_->GetModuleParAt(module_num);
507 // steer_writer.add_parameter_default(
508 // get_global_label_id(module_num, GlobalLabel::effective_c),
509 // // std::make_pair(module_par.effective_speed.value, module_par.effective_speed.error));
510 // std::make_pair(DEFAULT_EFFECTIVE_C, module_par.effective_speed.error));
511
512 // const auto offset_effective_c = module_par.t_diff * module_par.effective_speed;
513 // steer_writer.add_parameter_default(
514 // get_global_label_id(module_num, GlobalLabel::offset_effective_c),
515 // std::make_pair(offset_effective_c.value / SCALE_FACTOR, offset_effective_c.error / SCALE_FACTOR));
516 // }
517 // steer_writer.add_parameter_default(get_global_label_id(REFERENCE_BAR_NUM, GlobalLabel::tsync),
518 // std::make_pair(0.F, -1.F));
519 steer_writer.write();
520 }
521} // namespace R3B::Neuland::Calibration
#define R3BLOG(severity, x)
Definition R3BLogger.h:32
constexpr auto SCALE_FACTOR
constexpr auto DEFAULT_T_ERROR
constexpr auto DEFAULT_RES_FILENAME
constexpr auto MILLE_BUFFER_SIZE
auto add_hist(std::unique_ptr< TH1 > hist) -> TH1 *
auto add_graph(std::string_view graph_name, std::unique_ptr< GraphType > graph) -> GraphType *
auto get_pars() const -> const auto &
auto GetListOfModuleParRef() -> auto &
auto GetListOfModulePar() const -> const std::unordered_map< unsigned int, ::R3B::Neuland::HitModulePar > &
auto to_module_num_label(int par_num) -> std::pair< int, GlobalLabel >
void Calibrate(Cal2HitPar &hit_par) override
void add_spacial_local_constraint(int plane_id, const std::vector< MilleCalData > &plane_signals)
auto set_minimum_values(const std::vector< R3B::Neuland::BarCalData > &signals) -> bool
void EndOfEvent(unsigned int event_num=0) override
auto select_t_diff_signal(const std::vector< MilleCalData > &plane_data)
void fill_module_parameters(const Millepede::ResultReader &result, Neuland::Cal2HitPar &cal_to_hit_par)
std::unique_ptr< MilleDataProcessor > data_preprocessor_
void add_signal_t_sum(const MilleCalData &signal)
void add_signal_t_diff(const MilleCalData &signal)
auto SignalFilter(const std::vector< BarCalData > &signals) -> bool override
void AddSignals(const std::vector< BarCalData > &signals) override
auto get_global_label_id(int module_num, GlobalLabel label) -> int
void HistInit(DataMonitor &histograms) override
void set_filepath(std::string_view filepath)
Definition SteerWriter.h:48
STL class.
constexpr auto ModuleID2PlaneID(int moduleID) -> int
constexpr auto SQRT_12
constexpr auto BarLength
constexpr auto BarSize_XY
constexpr auto ModuleNum2ZPos(int module_num) -> T
constexpr auto IsPlaneIDHorizontal(int plane_id) -> bool
constexpr auto PlaneID2ZPos(int plane_id) -> T