File indexing completed on 2025-02-23 09:22:20
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030 #include "LETAccumulable.hh"
0031
0032 #include "G4EmCalculator.hh"
0033 #include "G4LogicalVolume.hh"
0034 #include "G4ParticleDefinition.hh"
0035
0036 #include "Hit.hh"
0037 #include "Manager.hh"
0038 #include "VoxelizedSensitiveDetector.hh"
0039 #include <G4SystemOfUnits.hh>
0040
0041 #include <tuple>
0042
0043 namespace RadioBio
0044 {
0045
0046
0047
0048 LETAccumulable::LETAccumulable() : VRadiobiologicalAccumulable("LET") {}
0049
0050
0051
0052 void LETAccumulable::Merge(const G4VAccumulable& rhs)
0053 {
0054 if (GetVerboseLevel() > 0) {
0055 G4cout << "LETAccumulable::Merge()" << G4endl;
0056 }
0057 const LETAccumulable& other = dynamic_cast<const LETAccumulable&>(rhs);
0058
0059
0060 fTotalLETT += other.GetTotalLETT();
0061 fDTotalLETT += other.GetDTotalLETT();
0062 fDTotalLETD += other.GetDTotalLETD();
0063 fTotalLETD += other.GetTotalLETD();
0064
0065
0066 auto rhsIonLetStore = other.GetIonLetStore();
0067
0068
0069 for (unsigned int l = 0; l < rhsIonLetStore.size(); l++) {
0070 G4int PDGencoding = rhsIonLetStore[l].GetPDGencoding();
0071 size_t q;
0072
0073 for (q = 0; q < fIonLetStore.size(); q++) {
0074
0075 if (fIonLetStore[q].GetPDGencoding() == PDGencoding) {
0076 if (rhsIonLetStore[l].IsPrimary() == fIonLetStore[q].IsPrimary()) break;
0077 }
0078 }
0079
0080 if (q == fIonLetStore.size())
0081 fIonLetStore.push_back(rhsIonLetStore[l]);
0082 else
0083 fIonLetStore[q].Merge(&(rhsIonLetStore[l]));
0084 }
0085 }
0086
0087
0088
0089 void LETAccumulable::Reset()
0090 {
0091 if (GetVerboseLevel() > 0) {
0092 G4cout << "LETAccumulable::Reset()" << G4endl;
0093 }
0094 if (fInitialized) {
0095 fTotalLETT = 0.0;
0096 fDTotalLETT = 0.0;
0097 fDTotalLETD = 0.0;
0098 fTotalLETD = 0.0;
0099 fIonLetStore.clear();
0100 }
0101 else {
0102 Initialize();
0103 }
0104 }
0105
0106
0107
0108
0109 void LETAccumulable::Accumulate(Hit* hit)
0110 {
0111 if (GetVerboseLevel() > 1) {
0112 G4cout << "LETAccumulable::Accumulate()" << G4endl;
0113 }
0114
0115
0116 if (hit->GetDeltaE() == 0. || hit->GetTrackLength() == 0.) return;
0117
0118
0119 G4int Z = hit->GetPartType()->GetAtomicNumber();
0120 if (Z < 1) return;
0121
0122 G4int PDGencoding = hit->GetPartType()->GetPDGEncoding();
0123
0124 if (PDGencoding == 22 || PDGencoding == 11) return;
0125
0126
0127 PDGencoding -= PDGencoding % 10;
0128
0129
0130 G4int xIndex = hit->GetXindex();
0131 G4int yIndex = hit->GetYindex();
0132 ;
0133 G4int zIndex = hit->GetZindex();
0134 ;
0135
0136 G4int voxel =
0137 VoxelizedSensitiveDetector::GetInstance()->GetThisVoxelNumber(xIndex, yIndex, zIndex);
0138
0139
0140 G4double ekinMean = hit->GetEkinMean();
0141
0142
0143 const G4ParticleDefinition* particleDef = hit->GetPartType();
0144
0145
0146 G4Material* mat = hit->GetPhysicalVolume()->GetLogicalVolume()->GetMaterial();
0147
0148
0149 G4EmCalculator emCal;
0150
0151 G4double Lsn = emCal.ComputeElectronicDEDX(ekinMean, particleDef, mat);
0152
0153
0154
0155
0156 G4double DE = hit->GetDeltaE();
0157
0158
0159 G4double DEElectrons = hit->GetElectronEnergy();
0160
0161
0162 G4double DX = hit->GetTrackLength();
0163
0164
0165 G4int trackID = hit->GetTrackID();
0166
0167
0168
0169 fTotalLETD[voxel] += (DE + DEElectrons) * Lsn;
0170
0171 fDTotalLETD[voxel] += DE + DEElectrons;
0172
0173 fTotalLETT[voxel] += DX * Lsn;
0174
0175 fDTotalLETT[voxel] += DX;
0176
0177 size_t l;
0178 for (l = 0; l < fIonLetStore.size(); l++) {
0179
0180 if (fIonLetStore[l].GetPDGencoding() == PDGencoding)
0181 if (((trackID == 1) && (fIonLetStore[l].IsPrimary()))
0182 || ((trackID != 1) && (!fIonLetStore[l].IsPrimary())))
0183 break;
0184 }
0185
0186
0187 if (l == fIonLetStore.size()) {
0188
0189 G4int A = particleDef->GetAtomicMass();
0190
0191
0192 G4String fullName = particleDef->GetParticleName();
0193 G4String name = fullName.substr(0, fullName.find("["));
0194
0195 IonLet ion(trackID, PDGencoding, fullName, name, Z, A,
0196 VoxelizedSensitiveDetector::GetInstance()->GetTotalVoxelNumber());
0197 fIonLetStore.push_back(ion);
0198 }
0199
0200
0201 fIonLetStore[l].Update(voxel, DE, DEElectrons, Lsn, DX);
0202 }
0203
0204
0205
0206 G4int LETAccumulable::GetVerboseLevel() const
0207 {
0208
0209 return Manager::GetInstance()->GetQuantity("LET")->GetVerboseLevel();
0210 }
0211
0212
0213
0214 void LETAccumulable::Initialize()
0215 {
0216 if (GetVerboseLevel() > 0) {
0217 G4cout << "LETAccumulable::Initialize(): " << G4endl;
0218 }
0219
0220 G4int voxNumber = VoxelizedSensitiveDetector::GetInstance()->GetTotalVoxelNumber();
0221
0222 fTotalLETT = array_type(0.0, voxNumber);
0223 fTotalLETD = array_type(0.0, voxNumber);
0224 fDTotalLETT = array_type(0.0, voxNumber);
0225 fDTotalLETD = array_type(0.0, voxNumber);
0226 fInitialized = true;
0227 }
0228
0229
0230
0231 }