Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
RadioBio::SD Class Reference

RadioBio sensitive detector class. More...

#include <Doxymodules_medical.h>

Inheritance diagram for RadioBio::SD:
G4VSensitiveDetector

Public Member Functions

 SD (const G4String &name, const G4String &hitsCollectionName)
 
 ~SD () override=default
 
void Initialize (G4HCofThisEvent *hitCollection) override
 
G4bool ProcessHits (G4Step *step, G4TouchableHistory *history) override
 
void EndOfEvent (G4HCofThisEvent *hitCollection) override
 

Private Attributes

RadioBioHitsCollectionfHitsCollection = nullptr
 

Detailed Description

RadioBio sensitive detector class.

The hits are accounted in hits in ProcessHits() function which is called by Geant4 kernel at each step. A hit is created with each step with non zero energy deposit.

Definition at line 145 of file Doxymodules_medical.h.

Constructor & Destructor Documentation

◆ SD()

RadioBio::SD::SD ( const G4String name,
const G4String hitsCollectionName 
)

Definition at line 47 of file SD.cc.

48{
49 collectionName.insert(hitsCollectionName);
50}

◆ ~SD()

RadioBio::SD::~SD ( )
overridedefault

Member Function Documentation

◆ Initialize()

void RadioBio::SD::Initialize ( G4HCofThisEvent hitCollection)
override

Definition at line 54 of file SD.cc.

55{
56 // Create hits collection
57 fHitsCollection = new RadioBioHitsCollection(SensitiveDetectorName, collectionName[0]);
58
59 // Add this collection in hce
60 G4int hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
61 hce->AddHitsCollection(hcID, fHitsCollection);
62}
RadioBioHitsCollection * fHitsCollection
Definition SD.hh:62
G4THitsCollection< Hit > RadioBioHitsCollection
Definition Hit.hh:110

◆ ProcessHits()

G4bool RadioBio::SD::ProcessHits ( G4Step step,
G4TouchableHistory history 
)
override

Definition at line 66 of file SD.cc.

67{
68 // FIXME why the namespace specification? Should not be necessary...
69 RadioBio::Hit* newHit = new RadioBio::Hit();
70
71 // Get the pre-step kinetic energy
72 G4double eKinPre = aStep->GetPreStepPoint()->GetKineticEnergy();
73 // Get the post-step kinetic energy
74 G4double eKinPost = aStep->GetPostStepPoint()->GetKineticEnergy();
75 // Get the step average kinetic energy
76 G4double eKinMean = (eKinPre + eKinPost) * 0.5;
77
78 const std::vector<const G4Track*>* secondary = aStep->GetSecondaryInCurrentStep();
79
80 size_t SecondarySize = (*secondary).size();
81 G4double EnergySecondary = 0.;
82
83 // Get secondary electrons energy deposited
84 if (SecondarySize) // Calculate only secondary particles
85 {
86 for (size_t numsec = 0; numsec < SecondarySize; numsec++) {
87 // Get the PDG code of every secondaty particles in current step
88 G4int PDGSecondary = (*secondary)[numsec]->GetDefinition()->GetPDGEncoding();
89
90 if (PDGSecondary == 11) // Calculate only secondary electrons
91 {
92 // Calculate the energy deposit of secondary electrons in current step
93 EnergySecondary += (*secondary)[numsec]->GetKineticEnergy();
94 }
95 }
96 }
97
98 // Update the hit with the necessary quantities
99 newHit->SetTrackID(aStep->GetTrack()->GetTrackID());
100 newHit->SetPartType(aStep->GetTrack()->GetParticleDefinition());
101 newHit->SetEkinMean(eKinMean);
102 newHit->SetDeltaE(aStep->GetTotalEnergyDeposit());
103 newHit->SetEinit(aStep->GetPreStepPoint()->GetKineticEnergy());
104 newHit->SetTrackLength(aStep->GetStepLength());
105 newHit->SetElectronEnergy(EnergySecondary);
106 newHit->SetPhysicalVolume(aStep->GetPreStepPoint()->GetPhysicalVolume());
107 newHit->SetVoxelIndexes(aStep->GetPreStepPoint()->GetTouchableHandle());
108
109 // Insert the hit in the hitcollection
110 fHitsCollection->insert(newHit);
111
112 // Accumulables are accumulated only if calculation is enabled
113 for (G4int i = 0; i < G4AccumulableManager::Instance()->GetNofAccumulables(); ++i) {
114 // Get the accumulable from the proper manager
115 G4VAccumulable* GenAcc = G4AccumulableManager::Instance()->GetAccumulable(i);
116
117 // Get the quantity from the proper manager using the name
118 auto q = Manager::GetInstance()->GetQuantity(GenAcc->GetName());
119
120 VRadiobiologicalAccumulable* radioAcc = dynamic_cast<VRadiobiologicalAccumulable*>(GenAcc);
121
122 // If the dynamic_cast did not work, this means that accumulable
123 // was not a VRadiobiologicalAccumulable
124 if (radioAcc == nullptr) continue;
125
126 // If calculation is enabled, accumulate
127 if (q->IsCalculationEnabled()) radioAcc->Accumulate(newHit);
128 }
129
130 return true;
131}
Detector hit class.
void SetTrackID(G4int track)
Definition Hit.hh:71
void SetEinit(G4double e)
Definition Hit.hh:75
void SetElectronEnergy(G4double elEnergy)
Definition Hit.hh:77
void SetPartType(const G4ParticleDefinition *part)
Definition Hit.hh:72
void SetDeltaE(double DeltaE)
Definition Hit.hh:74
void SetEkinMean(double EkinMean)
Definition Hit.hh:73
void SetVoxelIndexes(const G4TouchableHandle &TH)
Definition Hit.cc:109
void SetTrackLength(G4double x)
Definition Hit.hh:76
void SetPhysicalVolume(G4VPhysicalVolume *PV)
Definition Hit.hh:78
VRadiobiologicalQuantity * GetQuantity(G4String)
Definition Manager.cc:113
static Manager * GetInstance()
Definition Manager.cc:63

◆ EndOfEvent()

void RadioBio::SD::EndOfEvent ( G4HCofThisEvent hitCollection)
override

Definition at line 135 of file SD.cc.

136{
137 if (verboseLevel > 1) {
138 G4int nofHits = fHitsCollection->entries();
139 G4cout << G4endl << "-------->Hits Collection: in this event they are " << nofHits
140 << " hits in the detector slices: " << G4endl;
141 for (G4int i = 0; i < nofHits; i++)
142 (*fHitsCollection)[i]->Print();
143 }
144}

Member Data Documentation

◆ fHitsCollection

RadioBioHitsCollection* RadioBio::SD::fHitsCollection = nullptr
private

Definition at line 62 of file SD.hh.


The documentation for this class was generated from the following files:

Applications | User Support | Publications | Collaboration