Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:22:20

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 //
0027 /// \file radiobiology/src/LETAccumulable.cc
0028 /// \brief Implementation of the RadioBio::LETAccumulable class
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0047 
0048 LETAccumulable::LETAccumulable() : VRadiobiologicalAccumulable("LET") {}
0049 
0050 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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   // Merges standard counters
0060   fTotalLETT += other.GetTotalLETT();
0061   fDTotalLETT += other.GetDTotalLETT();
0062   fDTotalLETD += other.GetDTotalLETD();
0063   fTotalLETD += other.GetTotalLETD();
0064 
0065   // Merges ion counters
0066   auto rhsIonLetStore = other.GetIonLetStore();
0067 
0068   // Loop over rhs ions
0069   for (unsigned int l = 0; l < rhsIonLetStore.size(); l++) {
0070     G4int PDGencoding = rhsIonLetStore[l].GetPDGencoding();
0071     size_t q;
0072     // Loop over lhs ions to find the one
0073     for (q = 0; q < fIonLetStore.size(); q++) {
0074       // Check if primary status is the same
0075       if (fIonLetStore[q].GetPDGencoding() == PDGencoding) {
0076         if (rhsIonLetStore[l].IsPrimary() == fIonLetStore[q].IsPrimary()) break;
0077       }
0078     }
0079 
0080     if (q == fIonLetStore.size())  // Just copy the IonLet element in lhs store
0081       fIonLetStore.push_back(rhsIonLetStore[l]);
0082     else  // Merges rhs data with lhs ones
0083       fIonLetStore[q].Merge(&(rhsIonLetStore[l]));
0084   }
0085 }
0086 
0087 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0107 
0108 // To accumulate given the hit
0109 void LETAccumulable::Accumulate(Hit* hit)
0110 {
0111   if (GetVerboseLevel() > 1) {
0112     G4cout << "LETAccumulable::Accumulate()" << G4endl;
0113   }
0114 
0115   // No need to accumulate if no energy deposit or track length is zero.
0116   if (hit->GetDeltaE() == 0. || hit->GetTrackLength() == 0.) return;
0117 
0118   // Atomic number
0119   G4int Z = hit->GetPartType()->GetAtomicNumber();
0120   if (Z < 1) return;  // Calculate only protons and ions
0121 
0122   G4int PDGencoding = hit->GetPartType()->GetPDGEncoding();
0123 
0124   if (PDGencoding == 22 || PDGencoding == 11) return;  // do not accumulate for electrons or gammas
0125 
0126   // Get simple Particle data Group ID without excitation level
0127   PDGencoding -= PDGencoding % 10;
0128 
0129   // Get voxel number
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   // Get mean kinetic energy
0140   G4double ekinMean = hit->GetEkinMean();
0141 
0142   // Get particle definition
0143   const G4ParticleDefinition* particleDef = hit->GetPartType();
0144 
0145   // Get material
0146   G4Material* mat = hit->GetPhysicalVolume()->GetLogicalVolume()->GetMaterial();
0147 
0148   // ICRU stopping power calculation
0149   G4EmCalculator emCal;
0150   // Use the mean kinetic energy of ions in a step to calculate ICRU stopping power
0151   G4double Lsn = emCal.ComputeElectronicDEDX(ekinMean, particleDef, mat);
0152   // G4cout << "ERASE ME. Lsn = " << Lsn << G4endl;
0153   // G4cout << "ERASE ME. ekinMean = " << ekinMean << G4endl;
0154 
0155   // Get deposited energy
0156   G4double DE = hit->GetDeltaE();
0157 
0158   // Get deposited energy considering secondary electrons
0159   G4double DEElectrons = hit->GetElectronEnergy();
0160 
0161   // Get track length
0162   G4double DX = hit->GetTrackLength();
0163 
0164   // Get track ID
0165   G4int trackID = hit->GetTrackID();
0166 
0167   // Total LET calculation...
0168   // Total dose LET Numerator, including secondary electrons energy deposit
0169   fTotalLETD[voxel] += (DE + DEElectrons) * Lsn;
0170   // Total dose LET Denominator, including secondary electrons energy deposit
0171   fDTotalLETD[voxel] += DE + DEElectrons;
0172   // Total track LET Numerator
0173   fTotalLETT[voxel] += DX * Lsn;
0174   // Total track LET Denominator
0175   fDTotalLETT[voxel] += DX;
0176 
0177   size_t l;
0178   for (l = 0; l < fIonLetStore.size(); l++) {
0179     // Judge species of the current particle and store it
0180     if (fIonLetStore[l].GetPDGencoding() == PDGencoding)
0181       if (((trackID == 1) && (fIonLetStore[l].IsPrimary()))
0182           || ((trackID != 1) && (!fIonLetStore[l].IsPrimary())))
0183         break;
0184   }
0185 
0186   // Just another type of ion/particle for our store...
0187   if (l == fIonLetStore.size()) {
0188     // Mass number
0189     G4int A = particleDef->GetAtomicMass();
0190 
0191     // Particle name
0192     G4String fullName = particleDef->GetParticleName();
0193     G4String name = fullName.substr(0, fullName.find("["));  // Cut excitation energy [x.y]
0194 
0195     IonLet ion(trackID, PDGencoding, fullName, name, Z, A,
0196                VoxelizedSensitiveDetector::GetInstance()->GetTotalVoxelNumber());
0197     fIonLetStore.push_back(ion);
0198   }
0199 
0200   // Calculate ions LET and store them
0201   fIonLetStore[l].Update(voxel, DE, DEElectrons, Lsn, DX);
0202 }
0203 
0204 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0205 
0206 G4int LETAccumulable::GetVerboseLevel() const
0207 {
0208   // return same level of LET class
0209   return Manager::GetInstance()->GetQuantity("LET")->GetVerboseLevel();
0210 }
0211 
0212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0230 
0231 }  // namespace RadioBio