File indexing completed on 2025-01-31 09:22:17
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 #include "HadrontherapyDetectorSD.hh"
0030 #include "HadrontherapyDetectorHit.hh"
0031 #include "G4Step.hh"
0032 #include "G4VTouchable.hh"
0033 #include "G4TouchableHistory.hh"
0034 #include "G4SDManager.hh"
0035 #include "HadrontherapyMatrix.hh"
0036 #include "HadrontherapyLet.hh"
0037 #include "G4Track.hh"
0038 #include "G4SystemOfUnits.hh"
0039 #include "HadrontherapyMatrix.hh"
0040
0041
0042 #include "G4SteppingManager.hh"
0043 #include "G4TrackVector.hh"
0044 #include "HadrontherapySteppingAction.hh"
0045 #include "G4ios.hh"
0046 #include "G4SteppingManager.hh"
0047 #include "G4Track.hh"
0048 #include "G4Step.hh"
0049 #include "G4StepPoint.hh"
0050 #include "G4TrackStatus.hh"
0051 #include "G4TrackVector.hh"
0052 #include "G4ParticleDefinition.hh"
0053 #include "G4ParticleTypes.hh"
0054 #include "G4UserEventAction.hh"
0055 #include "G4TransportationManager.hh"
0056 #include "G4VSensitiveDetector.hh"
0057 #include "HadrontherapyRunAction.hh"
0058 #include "G4SystemOfUnits.hh"
0059 #include "HadrontherapyRBE.hh"
0060 #include <G4AccumulableManager.hh>
0061
0062
0063
0064 HadrontherapyDetectorSD::HadrontherapyDetectorSD(G4String name):
0065 G4VSensitiveDetector(name)
0066 {
0067 G4String HCname;
0068 collectionName.insert(HCname="HadrontherapyDetectorHitsCollection");
0069 HitsCollection = NULL;
0070 sensitiveDetectorName = name;
0071
0072 }
0073
0074
0075 HadrontherapyDetectorSD::~HadrontherapyDetectorSD()
0076 {}
0077
0078
0079 void HadrontherapyDetectorSD::Initialize(G4HCofThisEvent*)
0080 {
0081
0082 HitsCollection = new HadrontherapyDetectorHitsCollection(sensitiveDetectorName,
0083 collectionName[0]);
0084 }
0085
0086
0087 G4bool HadrontherapyDetectorSD::ProcessHits(G4Step* aStep, G4TouchableHistory* )
0088 {
0089
0090
0091 if (aStep -> GetPreStepPoint() -> GetPhysicalVolume() -> GetName() != "RODetectorZDivisionPhys") return false;
0092
0093
0094
0095 G4Track * theTrack = aStep -> GetTrack();
0096 G4double kineticEnergy = theTrack->GetKineticEnergy();
0097 G4ParticleDefinition *particleDef = theTrack -> GetDefinition();
0098
0099 G4String particleName = particleDef -> GetParticleName();
0100
0101
0102 G4int pdg = particleDef ->GetPDGEncoding();
0103
0104
0105 G4int trackID = theTrack -> GetTrackID();
0106
0107 G4double energyDeposit = aStep -> GetTotalEnergyDeposit();
0108
0109 G4double DX = aStep -> GetStepLength();
0110 G4int Z = particleDef-> GetAtomicNumber();
0111 G4int A = particleDef-> GetAtomicMass();
0112 G4StepPoint* PreStep = aStep->GetPreStepPoint();
0113
0114
0115 const G4VTouchable* touchable = aStep->GetPreStepPoint()->GetTouchable();
0116 G4int k = touchable->GetReplicaNumber(0);
0117 G4int i = touchable->GetReplicaNumber(2);
0118 G4int j = touchable->GetReplicaNumber(1);
0119
0120 G4TouchableHandle touchPreStep = PreStep->GetTouchableHandle();
0121 G4VPhysicalVolume* volumePre = touchPreStep->GetVolume();
0122 G4String namePre = volumePre->GetName();
0123
0124
0125
0126 HadrontherapyMatrix* matrix = HadrontherapyMatrix::GetInstance();
0127 HadrontherapyLet* let = HadrontherapyLet::GetInstance();
0128
0129 G4int* hitTrack = matrix -> GetHitTrack(i,j,k);
0130
0131
0132
0133 if (let)
0134 {
0135 if ( !(Z==0 && A==1) )
0136 {
0137 if( energyDeposit>0. && DX >0. )
0138 {
0139 if (pdg !=22 && pdg !=11)
0140 {
0141
0142
0143 G4double eKinPre = aStep -> GetPreStepPoint() -> GetKineticEnergy();
0144
0145 G4double eKinPost = aStep -> GetPostStepPoint() -> GetKineticEnergy();
0146
0147 G4double eKinMean = (eKinPre + eKinPost) * 0.5;
0148
0149
0150 const G4Material * materialStep = aStep -> GetPreStepPoint() -> GetMaterial();
0151
0152
0153 G4Step fstep = *theTrack -> GetStep();
0154
0155 const std::vector<const G4Track*> * secondary = fstep.GetSecondaryInCurrentStep();
0156
0157 size_t SecondarySize = (*secondary).size();
0158 G4double EnergySecondary = 0.;
0159
0160
0161 if (SecondarySize)
0162 {
0163 for (size_t numsec = 0; numsec< SecondarySize ; numsec ++)
0164 {
0165
0166 G4int PDGSecondary=(*secondary)[numsec]->GetDefinition()->GetPDGEncoding();
0167
0168 if(PDGSecondary == 11)
0169 {
0170
0171 EnergySecondary += (*secondary)[numsec]->GetKineticEnergy();
0172 }
0173 }
0174
0175 }
0176
0177
0178 let -> FillEnergySpectrum(trackID, particleDef, eKinMean, materialStep,
0179 energyDeposit,EnergySecondary,DX, i, j, k);
0180 }
0181 }
0182 }
0183 }
0184
0185
0186 if (matrix)
0187 {
0188
0189
0190
0191
0192 if ( *hitTrack != trackID )
0193 {
0194 *hitTrack = trackID;
0195
0196
0197
0198
0199 if ( Z >= 1) matrix -> Fill(trackID, particleDef, i, j, k, 0, true);
0200
0201 }
0202
0203 if(energyDeposit != 0)
0204 {
0205
0206
0207
0208
0209
0210
0211
0212 if ( Z>=1 )
0213 matrix -> Fill(trackID, particleDef, i, j, k, energyDeposit);
0214
0215
0216
0217 HadrontherapyDetectorHit* detectorHit = new HadrontherapyDetectorHit();
0218 detectorHit -> SetEdepAndPosition(i, j, k, energyDeposit);
0219 HitsCollection -> insert(detectorHit);
0220 }
0221 }
0222
0223 auto rbe = HadrontherapyRBE::GetInstance();
0224 if (rbe->IsCalculationEnabled())
0225 {
0226 if (!fRBEAccumulable)
0227 {
0228 fRBEAccumulable = dynamic_cast<HadrontherapyRBEAccumulable*>(G4AccumulableManager::Instance()->GetAccumulable("RBE"));
0229 if (!fRBEAccumulable)
0230 {
0231 G4Exception("HadrontherapyDetectorSD::ProcessHits", "NoAccumulable", FatalException, "Accumulable RBE not found.");
0232 }
0233 }
0234 if (A>0)
0235 {
0236 fRBEAccumulable->Accumulate(kineticEnergy / A, energyDeposit, DX, Z, i, j, k);
0237 }
0238 }
0239
0240
0241 return true;
0242 }
0243
0244
0245 void HadrontherapyDetectorSD::EndOfEvent(G4HCofThisEvent* HCE)
0246 {
0247
0248 static G4int HCID = -1;
0249 if(HCID < 0)
0250 {
0251 HCID = GetCollectionID(0);
0252 }
0253
0254 HCE -> AddHitsCollection(HCID,HitsCollection);
0255 }
0256