R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandCluster.cxx
Go to the documentation of this file.
1/******************************************************************************
2 * Copyright (C) 2019 GSI Helmholtzzentrum für Schwerionenforschung GmbH *
3 * Copyright (C) 2019-2025 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 "R3BNeulandCluster.h"
15#include <algorithm>
16#include <numeric>
17#include <stdexcept>
18
20{
21 auto min = std::min_element(fHits.cbegin(),
22 fHits.cend(),
23 [](const R3BNeulandHit& a, const R3BNeulandHit& b)
24 { return a.GetPosition().Z() < b.GetPosition().Z(); });
25 if (min == fHits.end())
26 {
27 throw std::logic_error("R3BNeulandCluster::GetFirstHit(): Cluster has no Hits!");
28 }
29 return *min;
30}
31
33{
34 auto min = std::min_element(fHits.cbegin(),
35 fHits.cend(),
36 [](const R3BNeulandHit& a, const R3BNeulandHit& b) { return a.GetT() < b.GetT(); });
37 if (min == fHits.end())
38 {
39 throw std::logic_error("R3BNeulandCluster::GetFirstHit(): Cluster has no Hits!");
40 }
41 return *min;
42}
43
45{
46 auto max = std::max_element(fHits.cbegin(),
47 fHits.cend(),
48 [](const R3BNeulandHit& a, const R3BNeulandHit& b) { return a.GetT() < b.GetT(); });
49 if (max == fHits.end())
50 {
51 throw std::logic_error("R3BNeulandCluster::GetLastHit(): Cluster has no Hits!");
52 }
53 return *max;
54}
55
57{
58 auto max = std::max_element(fHits.cbegin(),
59 fHits.cend(),
60 [](const R3BNeulandHit& a, const R3BNeulandHit& b) { return a.GetE() < b.GetE(); });
61 if (max == fHits.end())
62 {
63 throw std::logic_error("R3BNeulandCluster::GetLastHit(): Cluster has no Hits!");
64 }
65 return *max;
66}
67
69{
70 return std::accumulate(
71 fHits.begin(), fHits.end(), 0., [](const Double_t sum, const R3BNeulandHit& hit) { return sum + hit.GetE(); });
72}
73
74Double_t R3BNeulandCluster::GetT() const { return GetFirstHit().GetT(); }
75
77
79{
80 // analog to Geometrical Centroid \vec{c} = \frac{\sum_i (\vec{r}_{i} \cdot V_i)}{\sum_i V_i}
81 TVector3 centroid = std::accumulate(fHits.cbegin(),
82 fHits.cend(),
83 TVector3(),
84 [](const TVector3& c, const R3BNeulandHit& hit)
85 { return c + (hit.GetPosition() * hit.GetE()); });
86 return centroid * (1. / GetE());
87}
88
90{
91 const TVector3 centroid = GetEnergyCentroid();
92 Double_t mom = std::accumulate(fHits.cbegin(),
93 fHits.cend(),
94 0.,
95 [&](const Double_t c, const R3BNeulandHit& hit)
96 { return c + (hit.GetPosition() - centroid).Mag() * hit.GetE(); });
97 return mom / GetE();
98}
99
100Double_t R3BNeulandCluster::GetRCluster(Double_t beta) const
101{
102 // Equation 4.2 in TDR (Page 55).
103 return std::abs(beta - GetBeta()) / GetE();
104}
105
106Double_t R3BNeulandCluster::GetRECluster(Double_t ekin) const { return std::abs(ekin - GetEToF()) / GetE(); }
107
108std::ostream& operator<<(std::ostream& os, const R3BNeulandCluster& cluster)
109{
110 os << "R3BNeulandCluster: NeuLAND Cluster with size " << cluster.GetSize() << std::endl;
111 if (cluster.GetSize() > 0)
112 {
113 os << " Sum Energy: " << cluster.GetE() << std::endl
114 << " FirstHit: " << cluster.GetFirstHit() << " LastHit: " << cluster.GetLastHit();
115 }
116 return os;
117}
118
119void R3BNeulandCluster::Print(const Option_t*) const { std::cout << *this; }
120
ClassImp(R3B::Neuland::Cal2HitPar)
std::ostream & operator<<(std::ostream &os, const R3BNeulandCluster &cluster)
static const Double_t c
Double_t GetE() const
R3BNeulandHit GetFirstHit() const
TVector3 GetEnergyCentroid() const
void Print(const Option_t *) const override
Double_t GetRCluster(Double_t beta) const
Size_t GetSize() const
R3BNeulandHit GetMaxEnergyHit() const
Double_t GetT() const
Double_t GetBeta() const
Double_t GetEToF() const
Double_t GetRECluster(Double_t ekin) const
TVector3 GetPosition() const
Double_t GetEnergyMoment() const
R3BNeulandHit GetLastHit() const
R3BNeulandHit GetForemostHit() const
auto GetT() const -> double
auto GetPosition() const -> TVector3
auto GetE() const -> double