File indexing completed on 2025-02-23 09:20:59
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
0031
0032
0033 #include "SteppingAction.hh"
0034
0035 #include "DetectorConstruction.hh"
0036 #include "EventAction.hh"
0037 #include "HistoManager.hh"
0038 #include "Run.hh"
0039
0040 #include "G4PhysicalConstants.hh"
0041 #include "G4Positron.hh"
0042 #include "G4RunManager.hh"
0043
0044
0045
0046 SteppingAction::SteppingAction(DetectorConstruction* det, EventAction* evt)
0047 : fDetector(det), fEventAct(evt)
0048 {}
0049
0050
0051
0052 void SteppingAction::UserSteppingAction(const G4Step* aStep)
0053 {
0054
0055 const G4StepPoint* prePoint = aStep->GetPreStepPoint();
0056
0057
0058
0059 G4VPhysicalVolume* volume = prePoint->GetTouchableHandle()->GetVolume();
0060
0061 const G4Material* mat = volume->GetLogicalVolume()->GetMaterial();
0062 if (mat == fDetector->GetWorldMaterial()) return;
0063
0064 const G4StepPoint* endPoint = aStep->GetPostStepPoint();
0065 const G4ParticleDefinition* particle = aStep->GetTrack()->GetDefinition();
0066
0067
0068
0069 G4int absorNum = prePoint->GetTouchableHandle()->GetCopyNumber(0);
0070 G4int layerNum = prePoint->GetTouchableHandle()->GetCopyNumber(1);
0071
0072
0073 Run* run = static_cast<Run*>(G4RunManager::GetRunManager()->GetNonConstCurrentRun());
0074
0075
0076 G4double edep = aStep->GetTotalEnergyDeposit() * aStep->GetTrack()->GetWeight();
0077
0078
0079 G4double stepl = 0.;
0080 if (particle->GetPDGCharge() != 0.) {
0081 stepl = aStep->GetStepLength();
0082 run->AddChargedStep();
0083 }
0084 else {
0085 run->AddNeutralStep();
0086 }
0087
0088
0089
0090
0091 fEventAct->SumEnergy(absorNum, edep, stepl);
0092
0093
0094 if (edep > 0.) {
0095 G4AnalysisManager::Instance()->FillH1(kMaxAbsor + absorNum, G4double(layerNum + 1), edep);
0096 }
0097
0098
0099
0100 G4int Idnow = (fDetector->GetNbOfAbsor()) * layerNum + absorNum;
0101 G4int plane;
0102
0103
0104 if (endPoint->GetStepStatus() == fGeomBoundary) {
0105 G4ThreeVector position = endPoint->GetPosition();
0106 G4ThreeVector direction = endPoint->GetMomentumDirection();
0107 G4double sizeYZ = 0.5 * fDetector->GetCalorSizeYZ();
0108 G4double Eflow = endPoint->GetKineticEnergy();
0109 if (particle == G4Positron::Positron()) Eflow += 2 * electron_mass_c2;
0110 if ((std::abs(position.y()) >= sizeYZ) || (std::abs(position.z()) >= sizeYZ))
0111 run->SumLateralEleak(Idnow, Eflow);
0112 else if (direction.x() >= 0.)
0113 run->SumEnergyFlow(plane = Idnow + 1, Eflow);
0114 else
0115 run->SumEnergyFlow(plane = Idnow, -Eflow);
0116 }
0117
0118
0119
0120
0121
0122
0123 }
0124
0125
0126
0127 G4double SteppingAction::BirksAttenuation(const G4Step* aStep)
0128 {
0129
0130
0131
0132 const G4Material* material = aStep->GetTrack()->GetMaterial();
0133 G4double birk1 = material->GetIonisation()->GetBirksConstant();
0134 G4double destep = aStep->GetTotalEnergyDeposit();
0135 G4double stepl = aStep->GetStepLength();
0136 G4double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
0137
0138 G4double response = destep;
0139 if (birk1 * destep * stepl * charge != 0.) {
0140 response = destep / (1. + birk1 * destep / stepl);
0141 }
0142 return response;
0143 }
0144
0145