Loading...
Searching...
No Matches
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
RadioBio::LETAccumulable Class Reference

Accumulable of LET-related data (that must be thread-local). More...

#include <Doxymodules_medical.h>

Inheritance diagram for RadioBio::LETAccumulable:
RadioBio::VRadiobiologicalAccumulable G4VAccumulable

Public Types

using array_type = std::valarray< G4double >
 

Public Member Functions

 LETAccumulable ()
 
 LETAccumulable (const LETAccumulable &other)=default
 
void Merge (const G4VAccumulable &rhs) override
 
void Reset () override
 
void Accumulate (Hit *hit) override
 
const array_type GetTotalLETT () const
 
const array_type GetTotalLETD () const
 
const array_type GetDTotalLETT () const
 
const array_type GetDTotalLETD () const
 
const std::vector< IonLetGetIonLetStore () const
 
G4int GetVerboseLevel () const
 
- Public Member Functions inherited from RadioBio::VRadiobiologicalAccumulable
 VRadiobiologicalAccumulable (G4String name)
 
virtual ~VRadiobiologicalAccumulable ()
 

Private Member Functions

void Initialize ()
 

Private Attributes

G4bool fInitialized {false}
 
array_type fTotalLETT = {}
 
array_type fTotalLETD = {}
 
array_type fDTotalLETT = {}
 
array_type fDTotalLETD = {}
 
std::vector< IonLetfIonLetStore = {}
 

Detailed Description

Accumulable of LET-related data (that must be thread-local).

It keeps the sum of alpha and beta numerators/denominator, as well as energy deposits. The class is closely tied with the singLETon LET that is used both to calculate alphas and betas, and also to store results.

This is implemented as a customized G4VAccumulable with non-scalar data.

Note
There are two levels of merging (accumulation): 1) From more threads in one run (G4VAccumulable merging is applied) 2) (Optional) inter-run merging of data (implemented in LET).
std::valarray is used (instead of C arrays or std::vectors) to accumulate data for its logical simplicity.

Definition at line 135 of file Doxymodules_medical.h.

Member Typedef Documentation

◆ array_type

using RadioBio::LETAccumulable::array_type = std::valarray<G4double>

Definition at line 76 of file LETAccumulable.hh.

Constructor & Destructor Documentation

◆ LETAccumulable() [1/2]

RadioBio::LETAccumulable::LETAccumulable ( )

Definition at line 48 of file LETAccumulable.cc.

◆ LETAccumulable() [2/2]

RadioBio::LETAccumulable::LETAccumulable ( const LETAccumulable other)
default

Member Function Documentation

◆ Merge()

void RadioBio::LETAccumulable::Merge ( const G4VAccumulable rhs)
override

Definition at line 52 of file LETAccumulable.cc.

53{
54 if (GetVerboseLevel() > 0) {
55 G4cout << "LETAccumulable::Merge()" << G4endl;
56 }
57 const LETAccumulable& other = dynamic_cast<const LETAccumulable&>(rhs);
58
59 // Merges standard counters
60 fTotalLETT += other.GetTotalLETT();
61 fDTotalLETT += other.GetDTotalLETT();
62 fDTotalLETD += other.GetDTotalLETD();
63 fTotalLETD += other.GetTotalLETD();
64
65 // Merges ion counters
66 auto rhsIonLetStore = other.GetIonLetStore();
67
68 // Loop over rhs ions
69 for (unsigned int l = 0; l < rhsIonLetStore.size(); l++) {
70 G4int PDGencoding = rhsIonLetStore[l].GetPDGencoding();
71 size_t q;
72 // Loop over lhs ions to find the one
73 for (q = 0; q < fIonLetStore.size(); q++) {
74 // Check if primary status is the same
75 if (fIonLetStore[q].GetPDGencoding() == PDGencoding) {
76 if (rhsIonLetStore[l].IsPrimary() == fIonLetStore[q].IsPrimary()) break;
77 }
78 }
79
80 if (q == fIonLetStore.size()) // Just copy the IonLet element in lhs store
81 fIonLetStore.push_back(rhsIonLetStore[l]);
82 else // Merges rhs data with lhs ones
83 fIonLetStore[q].Merge(&(rhsIonLetStore[l]));
84 }
85}
std::vector< IonLet > fIonLetStore

◆ Reset()

void RadioBio::LETAccumulable::Reset ( )
override

Definition at line 89 of file LETAccumulable.cc.

90{
91 if (GetVerboseLevel() > 0) {
92 G4cout << "LETAccumulable::Reset()" << G4endl;
93 }
94 if (fInitialized) {
95 fTotalLETT = 0.0;
96 fDTotalLETT = 0.0;
97 fDTotalLETD = 0.0;
98 fTotalLETD = 0.0;
99 fIonLetStore.clear();
100 }
101 else {
102 Initialize();
103 }
104}

◆ Accumulate()

void RadioBio::LETAccumulable::Accumulate ( Hit hit)
overridevirtual

Implements RadioBio::VRadiobiologicalAccumulable.

Definition at line 109 of file LETAccumulable.cc.

110{
111 if (GetVerboseLevel() > 1) {
112 G4cout << "LETAccumulable::Accumulate()" << G4endl;
113 }
114
115 // No need to accumulate if no energy deposit or track length is zero.
116 if (hit->GetDeltaE() == 0. || hit->GetTrackLength() == 0.) return;
117
118 // Atomic number
119 G4int Z = hit->GetPartType()->GetAtomicNumber();
120 if (Z < 1) return; // Calculate only protons and ions
121
122 G4int PDGencoding = hit->GetPartType()->GetPDGEncoding();
123
124 if (PDGencoding == 22 || PDGencoding == 11) return; // do not accumulate for electrons or gammas
125
126 // Get simple Particle data Group ID without excitation level
127 PDGencoding -= PDGencoding % 10;
128
129 // Get voxel number
130 G4int xIndex = hit->GetXindex();
131 G4int yIndex = hit->GetYindex();
132 ;
133 G4int zIndex = hit->GetZindex();
134 ;
135
136 G4int voxel =
138
139 // Get mean kinetic energy
140 G4double ekinMean = hit->GetEkinMean();
141
142 // Get particle definition
143 const G4ParticleDefinition* particleDef = hit->GetPartType();
144
145 // Get material
146 G4Material* mat = hit->GetPhysicalVolume()->GetLogicalVolume()->GetMaterial();
147
148 // ICRU stopping power calculation
149 G4EmCalculator emCal;
150 // Use the mean kinetic energy of ions in a step to calculate ICRU stopping power
151 G4double Lsn = emCal.ComputeElectronicDEDX(ekinMean, particleDef, mat);
152 // G4cout << "ERASE ME. Lsn = " << Lsn << G4endl;
153 // G4cout << "ERASE ME. ekinMean = " << ekinMean << G4endl;
154
155 // Get deposited energy
156 G4double DE = hit->GetDeltaE();
157
158 // Get deposited energy considering secondary electrons
159 G4double DEElectrons = hit->GetElectronEnergy();
160
161 // Get track length
162 G4double DX = hit->GetTrackLength();
163
164 // Get track ID
165 G4int trackID = hit->GetTrackID();
166
167 // Total LET calculation...
168 // Total dose LET Numerator, including secondary electrons energy deposit
169 fTotalLETD[voxel] += (DE + DEElectrons) * Lsn;
170 // Total dose LET Denominator, including secondary electrons energy deposit
171 fDTotalLETD[voxel] += DE + DEElectrons;
172 // Total track LET Numerator
173 fTotalLETT[voxel] += DX * Lsn;
174 // Total track LET Denominator
175 fDTotalLETT[voxel] += DX;
176
177 size_t l;
178 for (l = 0; l < fIonLetStore.size(); l++) {
179 // Judge species of the current particle and store it
180 if (fIonLetStore[l].GetPDGencoding() == PDGencoding)
181 if (((trackID == 1) && (fIonLetStore[l].IsPrimary()))
182 || ((trackID != 1) && (!fIonLetStore[l].IsPrimary())))
183 break;
184 }
185
186 // Just another type of ion/particle for our store...
187 if (l == fIonLetStore.size()) {
188 // Mass number
189 G4int A = particleDef->GetAtomicMass();
190
191 // Particle name
192 G4String fullName = particleDef->GetParticleName();
193 G4String name = fullName.substr(0, fullName.find("[")); // Cut excitation energy [x.y]
194
195 IonLet ion(trackID, PDGencoding, fullName, name, Z, A,
196 VoxelizedSensitiveDetector::GetInstance()->GetTotalVoxelNumber());
197 fIonLetStore.push_back(ion);
198 }
199
200 // Calculate ions LET and store them
201 fIonLetStore[l].Update(voxel, DE, DEElectrons, Lsn, DX);
202}
static VoxelizedSensitiveDetector * GetInstance()
Static method to retrieve a pointer to the only object existing in the simulation.
G4int GetThisVoxelNumber(G4int x, G4int j, G4int k) const
Method to get the absolute voxel index given its indexes in the three dimensions.

◆ GetTotalLETT()

const array_type RadioBio::LETAccumulable::GetTotalLETT ( ) const
inline

Definition at line 79 of file LETAccumulable.hh.

79{ return fTotalLETT; }

◆ GetTotalLETD()

const array_type RadioBio::LETAccumulable::GetTotalLETD ( ) const
inline

Definition at line 80 of file LETAccumulable.hh.

80{ return fTotalLETD; }

◆ GetDTotalLETT()

const array_type RadioBio::LETAccumulable::GetDTotalLETT ( ) const
inline

Definition at line 81 of file LETAccumulable.hh.

81{ return fDTotalLETT; }

◆ GetDTotalLETD()

const array_type RadioBio::LETAccumulable::GetDTotalLETD ( ) const
inline

Definition at line 82 of file LETAccumulable.hh.

82{ return fDTotalLETD; }

◆ GetIonLetStore()

const std::vector< IonLet > RadioBio::LETAccumulable::GetIonLetStore ( ) const
inline

Definition at line 84 of file LETAccumulable.hh.

84{ return fIonLetStore; }

◆ GetVerboseLevel()

G4int RadioBio::LETAccumulable::GetVerboseLevel ( ) const

Definition at line 206 of file LETAccumulable.cc.

207{
208 // return same level of LET class
210}
VRadiobiologicalQuantity * GetQuantity(G4String)
Definition Manager.cc:113
static Manager * GetInstance()
Definition Manager.cc:63

◆ Initialize()

void RadioBio::LETAccumulable::Initialize ( )
private

Definition at line 214 of file LETAccumulable.cc.

215{
216 if (GetVerboseLevel() > 0) {
217 G4cout << "LETAccumulable::Initialize(): " << G4endl;
218 }
219
221
222 fTotalLETT = array_type(0.0, voxNumber);
223 fTotalLETD = array_type(0.0, voxNumber);
224 fDTotalLETT = array_type(0.0, voxNumber);
225 fDTotalLETD = array_type(0.0, voxNumber);
226 fInitialized = true;
227}
std::valarray< G4double > array_type
G4int GetTotalVoxelNumber() const
Method to get the total voxel number.

Member Data Documentation

◆ fInitialized

G4bool RadioBio::LETAccumulable::fInitialized {false}
private

Definition at line 92 of file LETAccumulable.hh.

92{false};

◆ fTotalLETT

array_type RadioBio::LETAccumulable::fTotalLETT = {}
private

Definition at line 94 of file LETAccumulable.hh.

94{};

◆ fTotalLETD

array_type RadioBio::LETAccumulable::fTotalLETD = {}
private

Definition at line 95 of file LETAccumulable.hh.

95{};

◆ fDTotalLETT

array_type RadioBio::LETAccumulable::fDTotalLETT = {}
private

Definition at line 96 of file LETAccumulable.hh.

96{};

◆ fDTotalLETD

array_type RadioBio::LETAccumulable::fDTotalLETD = {}
private

Definition at line 97 of file LETAccumulable.hh.

97{};

◆ fIonLetStore

std::vector<IonLet> RadioBio::LETAccumulable::fIonLetStore = {}
private

Definition at line 99 of file LETAccumulable.hh.

99{};

The documentation for this class was generated from the following files:

Applications | User Support | Publications | Collaboration