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/SD.cc
0028 /// \brief Implementation of the RadioBio::SD class
0029 
0030 #include "SD.hh"
0031 
0032 #include "G4AccumulableManager.hh"
0033 #include "G4HCofThisEvent.hh"
0034 #include "G4SDManager.hh"
0035 #include "G4Step.hh"
0036 #include "G4ThreeVector.hh"
0037 
0038 #include "Hit.hh"
0039 #include "Manager.hh"
0040 #include "VRadiobiologicalAccumulable.hh"
0041 
0042 namespace RadioBio
0043 {
0044 
0045 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0046 
0047 SD::SD(const G4String& name, const G4String& hitsCollectionName) : G4VSensitiveDetector(name)
0048 {
0049   collectionName.insert(hitsCollectionName);
0050 }
0051 
0052 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0053 
0054 void SD::Initialize(G4HCofThisEvent* hce)
0055 {
0056   // Create hits collection
0057   fHitsCollection = new RadioBioHitsCollection(SensitiveDetectorName, collectionName[0]);
0058 
0059   // Add this collection in hce
0060   G4int hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
0061   hce->AddHitsCollection(hcID, fHitsCollection);
0062 }
0063 
0064 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0065 
0066 G4bool SD::ProcessHits(G4Step* aStep, G4TouchableHistory*)
0067 {
0068   // FIXME why the namespace specification? Should not be necessary...
0069   RadioBio::Hit* newHit = new RadioBio::Hit();
0070 
0071   // Get the pre-step kinetic energy
0072   G4double eKinPre = aStep->GetPreStepPoint()->GetKineticEnergy();
0073   // Get the post-step kinetic energy
0074   G4double eKinPost = aStep->GetPostStepPoint()->GetKineticEnergy();
0075   // Get the step average kinetic energy
0076   G4double eKinMean = (eKinPre + eKinPost) * 0.5;
0077 
0078   const std::vector<const G4Track*>* secondary = aStep->GetSecondaryInCurrentStep();
0079 
0080   size_t SecondarySize = (*secondary).size();
0081   G4double EnergySecondary = 0.;
0082 
0083   // Get secondary electrons energy deposited
0084   if (SecondarySize)  // Calculate only secondary particles
0085   {
0086     for (size_t numsec = 0; numsec < SecondarySize; numsec++) {
0087       // Get the PDG code of every secondaty particles in current step
0088       G4int PDGSecondary = (*secondary)[numsec]->GetDefinition()->GetPDGEncoding();
0089 
0090       if (PDGSecondary == 11)  // Calculate only secondary electrons
0091       {
0092         // Calculate the energy deposit of secondary electrons in current step
0093         EnergySecondary += (*secondary)[numsec]->GetKineticEnergy();
0094       }
0095     }
0096   }
0097 
0098   // Update the hit with the necessary quantities
0099   newHit->SetTrackID(aStep->GetTrack()->GetTrackID());
0100   newHit->SetPartType(aStep->GetTrack()->GetParticleDefinition());
0101   newHit->SetEkinMean(eKinMean);
0102   newHit->SetDeltaE(aStep->GetTotalEnergyDeposit());
0103   newHit->SetEinit(aStep->GetPreStepPoint()->GetKineticEnergy());
0104   newHit->SetTrackLength(aStep->GetStepLength());
0105   newHit->SetElectronEnergy(EnergySecondary);
0106   newHit->SetPhysicalVolume(aStep->GetPreStepPoint()->GetPhysicalVolume());
0107   newHit->SetVoxelIndexes(aStep->GetPreStepPoint()->GetTouchableHandle());
0108 
0109   // Insert the hit in the hitcollection
0110   fHitsCollection->insert(newHit);
0111 
0112   // Accumulables are accumulated only if calculation is enabled
0113   for (G4int i = 0; i < G4AccumulableManager::Instance()->GetNofAccumulables(); ++i) {
0114     // Get the accumulable from the proper manager
0115     G4VAccumulable* GenAcc = G4AccumulableManager::Instance()->GetAccumulable(i);
0116 
0117     // Get the quantity from the proper manager using the name
0118     auto q = Manager::GetInstance()->GetQuantity(GenAcc->GetName());
0119 
0120     VRadiobiologicalAccumulable* radioAcc = dynamic_cast<VRadiobiologicalAccumulable*>(GenAcc);
0121 
0122     // If the dynamic_cast did not work, this means that accumulable
0123     // was not a VRadiobiologicalAccumulable
0124     if (radioAcc == nullptr) continue;
0125 
0126     // If calculation is enabled, accumulate
0127     if (q->IsCalculationEnabled()) radioAcc->Accumulate(newHit);
0128   }
0129 
0130   return true;
0131 }
0132 
0133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0134 
0135 void SD::EndOfEvent(G4HCofThisEvent*)
0136 {
0137   if (verboseLevel > 1) {
0138     G4int nofHits = fHitsCollection->entries();
0139     G4cout << G4endl << "-------->Hits Collection: in this event they are " << nofHits
0140            << " hits in the detector slices: " << G4endl;
0141     for (G4int i = 0; i < nofHits; i++)
0142       (*fHitsCollection)[i]->Print();
0143   }
0144 }
0145 
0146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0147 
0148 }  // namespace RadioBio