R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandNeutronsRValue.cxx
Go to the documentation of this file.
2#include <FairRootManager.h>
3#include <IsElastic.h>
4
6 std::string_view inputMult,
7 std::string_view inputCluster,
8 std::string_view output)
9 : FairTask("R3BNeulandNeutronsRValue")
10 , fEkinRefMeV(EkinRefMeV)
11 , fInputMultName(inputMult)
12 , fMultiplicity(nullptr)
13 , fClusters(inputCluster)
14 , fNeutrons(output)
15{
16}
17
19{
20 auto* ioman = FairRootManager::Instance();
21 if (ioman == nullptr)
22 {
23 throw R3B::runtime_error("TCAInputConnector: No FairRootManager");
24 }
25 fMultiplicity = ioman->InitObjectAs<const R3BNeulandMultiplicity*>(fInputMultName.c_str());
26
27 fClusters.init();
28 fNeutrons.init();
29 return kSUCCESS;
30}
31
32void R3BNeulandNeutronsRValue::Exec(Option_t* /*option*/)
33{
34 fNeutrons.clear();
35 cluster_buffer_.clear();
36
37 cluster_buffer_ = fClusters.get();
38
39 // Recreate R3BNeutronTracker2D Advanced Method
40 // FilterClustersByElasticScattering(clusters); // Check all pairs of clusters. Remove clusters from elastic
41 // scattering FilterClustersByEnergyDeposit(clusters); FilterClustersByKineticEnergy(clusters);
42 SortClustersByRValue(cluster_buffer_);
43 PrioritizeTimeWiseFirstCluster(cluster_buffer_);
44
45 const auto mult = fMultiplicity->GetMultiplicity();
46 for (size_t index = 0; index < cluster_buffer_.size() && index < mult; index++)
47 {
48 fNeutrons.get().emplace_back(cluster_buffer_.at(index));
49 }
50}
51
52void R3BNeulandNeutronsRValue::SortClustersByRValue(std::vector<R3BNeulandCluster>& clusters) const
53{
54 std::sort(clusters.begin(),
55 clusters.end(),
56 [this](const R3BNeulandCluster& one, const R3BNeulandCluster& other)
57 { return one.GetRECluster(fEkinRefMeV) < other.GetRECluster(fEkinRefMeV); });
58}
59
60void R3BNeulandNeutronsRValue::PrioritizeTimeWiseFirstCluster(std::vector<R3BNeulandCluster>& clusters)
61{
62 auto timewiseFirstCluster = std::min_element(clusters.begin(),
63 clusters.end(),
64 [](const R3BNeulandCluster& one, const R3BNeulandCluster& other)
65 { return one.GetT() < other.GetT(); });
66 // Put first cluster in front
67 std::rotate(clusters.begin(), timewiseFirstCluster, timewiseFirstCluster + 1);
68}
69
70void R3BNeulandNeutronsRValue::FilterClustersByEnergyDeposit(std::vector<R3BNeulandCluster>& clusters)
71{
72 clusters.erase(
73 std::remove_if(
74 clusters.begin(), clusters.end(), [&](const R3BNeulandCluster& cluster) { return cluster.GetE() < 2.5; }),
75 clusters.end());
76}
77
78void R3BNeulandNeutronsRValue::FilterClustersByKineticEnergy(std::vector<R3BNeulandCluster>& clusters) const
79{
80 clusters.erase(std::remove_if(clusters.begin(),
81 clusters.end(),
82 [this](const R3BNeulandCluster& cluster)
83 { return std::abs(cluster.GetEToF() - fEkinRefMeV) / fEkinRefMeV > 0.05; }),
84 clusters.end());
85}
86
87void R3BNeulandNeutronsRValue::FilterClustersByElasticScattering(std::vector<R3BNeulandCluster>& clusters)
88{
89 std::map<const R3BNeulandCluster*, bool> marked;
90
91 for (const auto& cluster : clusters)
92 {
93 marked[&cluster] = false;
94 }
95
96 for (const auto& one_cluster : clusters)
97 {
98 for (const auto& other_cluster : clusters)
99 {
100 if (&one_cluster != &other_cluster && one_cluster.GetT() < other_cluster.GetT() &&
101 Neuland::IsElastic(&one_cluster, &other_cluster))
102 {
103 marked[&other_cluster] = true;
104 }
105 }
106 }
107
108 clusters.erase(
109 std::remove_if(
110 clusters.begin(), clusters.end(), [&](const R3BNeulandCluster& cluster) { return marked.at(&cluster); }),
111 clusters.end());
112}
113
ClassImp(R3B::Neuland::Cal2HitPar)
void Exec(Option_t *) override
auto Init() -> InitStatus override
R3BNeulandNeutronsRValue(double EkinRefMeV, std::string_view inputMult="NeulandMultiplicity", std::string_view inputCluster="NeulandClusters", std::string_view output="NeulandNeutrons")
bool IsElastic(const R3BNeulandCluster *, const R3BNeulandCluster *)
Definition IsElastic.cxx:19