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 "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
0046
0047 SD::SD(const G4String& name, const G4String& hitsCollectionName) : G4VSensitiveDetector(name)
0048 {
0049 collectionName.insert(hitsCollectionName);
0050 }
0051
0052
0053
0054 void SD::Initialize(G4HCofThisEvent* hce)
0055 {
0056
0057 fHitsCollection = new RadioBioHitsCollection(SensitiveDetectorName, collectionName[0]);
0058
0059
0060 G4int hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
0061 hce->AddHitsCollection(hcID, fHitsCollection);
0062 }
0063
0064
0065
0066 G4bool SD::ProcessHits(G4Step* aStep, G4TouchableHistory*)
0067 {
0068
0069 RadioBio::Hit* newHit = new RadioBio::Hit();
0070
0071
0072 G4double eKinPre = aStep->GetPreStepPoint()->GetKineticEnergy();
0073
0074 G4double eKinPost = aStep->GetPostStepPoint()->GetKineticEnergy();
0075
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
0084 if (SecondarySize)
0085 {
0086 for (size_t numsec = 0; numsec < SecondarySize; numsec++) {
0087
0088 G4int PDGSecondary = (*secondary)[numsec]->GetDefinition()->GetPDGEncoding();
0089
0090 if (PDGSecondary == 11)
0091 {
0092
0093 EnergySecondary += (*secondary)[numsec]->GetKineticEnergy();
0094 }
0095 }
0096 }
0097
0098
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
0110 fHitsCollection->insert(newHit);
0111
0112
0113 for (G4int i = 0; i < G4AccumulableManager::Instance()->GetNofAccumulables(); ++i) {
0114
0115 G4VAccumulable* GenAcc = G4AccumulableManager::Instance()->GetAccumulable(i);
0116
0117
0118 auto q = Manager::GetInstance()->GetQuantity(GenAcc->GetName());
0119
0120 VRadiobiologicalAccumulable* radioAcc = dynamic_cast<VRadiobiologicalAccumulable*>(GenAcc);
0121
0122
0123
0124 if (radioAcc == nullptr) continue;
0125
0126
0127 if (q->IsCalculationEnabled()) radioAcc->Accumulate(newHit);
0128 }
0129
0130 return true;
0131 }
0132
0133
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
0147
0148 }