Back to home page

EIC code displayed by LXR

 
 

    


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

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/DoseAccumulable.cc
0028 /// \brief Implementation of the RadioBio::DoseAccumulable class
0029 
0030 #include "DoseAccumulable.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 DoseAccumulable::DoseAccumulable() : VRadiobiologicalAccumulable("Dose") {}
0049 
0050 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0051 
0052 void DoseAccumulable::Merge(const G4VAccumulable& rhs)
0053 {
0054   if (GetVerboseLevel() > 1) {
0055     G4cout << "DoseAccumulable::Merge()" << G4endl;
0056   }
0057   const DoseAccumulable& other = dynamic_cast<const DoseAccumulable&>(rhs);
0058 
0059   // Merges the counter
0060   fEnDep += other.GetEnDeposit();
0061 }
0062 
0063 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0064 
0065 void DoseAccumulable::Reset()
0066 {
0067   if (GetVerboseLevel() > 0) {
0068     G4cout << "DoseAccumulable::Reset()" << G4endl;
0069   }
0070   if (fInitialized) {
0071     fEnDep = 0.0;
0072   }
0073   else {
0074     Initialize();
0075   }
0076 }
0077 
0078 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0079 
0080 // To accumulate given the hit
0081 void DoseAccumulable::Accumulate(Hit* hit)
0082 {
0083   // Calculation done only if energy is deposited
0084   if (hit->GetDeltaE() <= 0.) return;
0085 
0086   if (GetVerboseLevel() > 1) {
0087     G4cout << "DoseAccumulable::Accumulate()" << G4endl;
0088   }
0089 
0090   // Get voxel number
0091   G4int xIndex = hit->GetXindex();
0092   G4int yIndex = hit->GetYindex();
0093   G4int zIndex = hit->GetZindex();
0094 
0095   G4int voxel =
0096     VoxelizedSensitiveDetector::GetInstance()->GetThisVoxelNumber(xIndex, yIndex, zIndex);
0097 
0098   // Get deposited energy
0099   G4double DE = hit->GetDeltaE();
0100 
0101   // Update total dose
0102   fEnDep[voxel] += DE;
0103 }
0104 
0105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0106 
0107 G4int DoseAccumulable::GetVerboseLevel() const
0108 {
0109   // Returns same verbosity of Dose class
0110   return Manager::GetInstance()->GetQuantity("Dose")->GetVerboseLevel();
0111 }
0112 
0113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0114 
0115 void DoseAccumulable::Initialize()
0116 {
0117   if (GetVerboseLevel() > 0) {
0118     G4cout << "DoseAccumulable::Initialize(): " << G4endl;
0119   }
0120 
0121   fEnDep = array_type(0.0, VoxelizedSensitiveDetector::GetInstance()->GetTotalVoxelNumber());
0122 
0123   fInitialized = true;
0124 }
0125 
0126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0127 
0128 }  // namespace RadioBio