R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandCalToHitPar.h
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#pragma once
15
16#include "R3BShared.h"
17#include <R3BNeulandParSet.h>
18#include <unordered_map>
19
20namespace R3B::Neuland
21{
36
37 class Cal2HitPar : public ParSet
38 {
39 public:
40 explicit Cal2HitPar(std::string_view name = "NeulandHitPar",
41 std::string_view title = "Neuland Cal2Hit calibration",
42 std::string_view context = "TestDefaultContext",
43 Bool_t own = kTRUE);
44
45 Cal2HitPar(const Cal2HitPar&) = default;
46 Cal2HitPar(Cal2HitPar&&) = default;
47 auto operator=(const Cal2HitPar&) -> Cal2HitPar& = default;
48 auto operator=(Cal2HitPar&&) -> Cal2HitPar& = default;
49 ~Cal2HitPar() override = default;
50
52 void Reset() { clear(); }
53
54 // setters:
55 void SetDistanceToTarget(double distance) { distance_to_target_ = distance; }
56 void SetEnergyCutoff(double cutoff) { energy_cut_ = cutoff; }
57 void SetGlobalTimeOffset(double offset) { global_time_offset_ = offset; }
58 void SetNumOfModules(int num) { num_of_modules = num; }
59 void AddModulePar(const HitModulePar& module_par)
60 {
61 const auto mNum = module_par.module_num;
62 module_pars_.insert_or_assign(mNum, module_par);
63 }
64
65 // getter:
66 auto GetDistanceToTarget() const { return distance_to_target_; }
67 auto GetEnergyCutoff() const { return energy_cut_; }
68 auto GetGlobalTimeOffset() const { return global_time_offset_; }
69 [[deprecated("Use GetDistancesToFirstPlane instead")]] auto GetDistanceToFirstPlane(
70 unsigned int plane_num) const
71 {
72 return distances_to_first_plane_.at(plane_num);
73 }
74 auto GetDistancesToFirstPlane() const -> const auto& { return distances_to_first_plane_; }
75 auto GetNumModulePar() const { return module_pars_.size(); }
76 auto GetModuleParAt(unsigned int module_num) const -> const ::R3B::Neuland::HitModulePar&
77 {
78 return module_pars_.at(module_num);
79 }
80 auto HasModuleParAt(int module_num) const -> bool
81 {
82 return module_pars_.find(module_num) != module_pars_.end();
83 }
84 auto GetModulePars() const -> const std::unordered_map<unsigned int, ::R3B::Neuland::HitModulePar>&
85 {
86 return module_pars_;
87 }
88 // no auto because pybind from pyROOT
89 auto GetListOfModulePar() const -> const std::unordered_map<unsigned int, ::R3B::Neuland::HitModulePar>&
90 {
91 return module_pars_;
92 }
93 auto GetListOfModuleParRef() -> auto& { return module_pars_; }
94 auto GetNumOfModules() const -> int { return num_of_modules; }
95
96 private:
97 int num_of_modules = 0;
98 double global_time_offset_ = 0.; // in ns
99 double distance_to_target_ = 0.; // in cm
100 double energy_cut_ = 0.; // in MeV
101 std::vector<double> distances_to_first_plane_;
102 std::unordered_map<unsigned int, ::R3B::Neuland::HitModulePar> module_pars_;
103 void clear() override
104 {
105 global_time_offset_ = 0.;
106 distance_to_target_ = 0.;
107 energy_cut_ = 0.;
108 distances_to_first_plane_.clear();
109 module_pars_.clear();
110 }
111
112 public:
114 };
115
116} // namespace R3B::Neuland
117
R3B::Neuland::Cal2HitPar R3BNeulandHitPar2
void AddModulePar(const HitModulePar &module_par)
auto GetDistancesToFirstPlane() const -> const auto &
Cal2HitPar(std::string_view name="NeulandHitPar", std::string_view title="Neuland Cal2Hit calibration", std::string_view context="TestDefaultContext", Bool_t own=kTRUE)
void SetDistanceToTarget(double distance)
void SetGlobalTimeOffset(double offset)
ClassDefOverride(Cal2HitPar, 2)
auto GetListOfModuleParRef() -> auto &
auto operator=(const Cal2HitPar &) -> Cal2HitPar &=default
auto GetModulePars() const -> const std::unordered_map< unsigned int, ::R3B::Neuland::HitModulePar > &
auto GetNumOfModules() const -> int
auto GetModuleParAt(unsigned int module_num) const -> const ::R3B::Neuland::HitModulePar &
void SetEnergyCutoff(double cutoff)
auto operator=(Cal2HitPar &&) -> Cal2HitPar &=default
auto HasModuleParAt(int module_num) const -> bool
Cal2HitPar(const Cal2HitPar &)=default
~Cal2HitPar() override=default
Cal2HitPar(Cal2HitPar &&)=default
auto GetListOfModulePar() const -> const std::unordered_map< unsigned int, ::R3B::Neuland::HitModulePar > &
auto GetDistanceToFirstPlane(unsigned int plane_num) const
ParSet(std::string_view name="parSet", std::string_view title="Neuland parameter", std::string_view context="TestDefaultContext", bool own=true)
Simulation of NeuLAND Bar/Paddle.
ValueError< double > light_attenuation_factor
LRPair< ValueError< double > > pedestal
ValueError< double > light_attenuation_length
ValueError< double > effective_speed
LRPair< ValueError< double > > pmt_threshold
ClassDefNV(HitModulePar, 2)
LRPair< ValueError< double > > energy_gain
LRPair< ValueError< double > > pmt_saturation