File indexing completed on 2025-04-15 08:09:25
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
0034 #include "SteppingAction.hh"
0035
0036 #include "Run.hh"
0037
0038 #include "G4IonTable.hh"
0039 #include "G4LossTableManager.hh"
0040 #include "G4ParticleDefinition.hh"
0041 #include "G4ParticleTypes.hh"
0042 #include "G4Step.hh"
0043 #include "G4StepPoint.hh"
0044 #include "G4SystemOfUnits.hh"
0045 #include "G4TouchableHistory.hh"
0046 #include "G4Track.hh"
0047 #include "G4VPhysicalVolume.hh"
0048 #include "G4VSolid.hh"
0049 #include "G4VTouchable.hh"
0050
0051 const std::array<G4String, SteppingAction::fkNumberScoringVolumes>
0052 SteppingAction::fkArrayScoringVolumeNames = {"downstream", "side", "upstream"};
0053
0054 const std::array<G4String, SteppingAction::fkNumberKinematicRegions>
0055 SteppingAction::fkArrayKinematicRegionNames = {"", "below 20 MeV", "above 20 MeV"};
0056
0057 const std::array<G4String, SteppingAction::fkNumberParticleTypes>
0058 SteppingAction::fkArrayParticleTypeNames = {"all", "electron", "gamma", "muon",
0059 "neutrino", "pion", "neutron", "proton",
0060 "ion", "otherMeson", "otherBaryon"};
0061
0062
0063
0064 G4int SteppingAction::GetIndex(const G4int iScoringVolume, const G4int iKinematicRegion,
0065 const G4int iParticleType)
0066 {
0067 G4int index = -1;
0068 if (iScoringVolume >= 0 && iScoringVolume < fkNumberScoringVolumes && iKinematicRegion >= 0
0069 && iKinematicRegion < fkNumberKinematicRegions && iParticleType >= 0
0070 && iParticleType < fkNumberParticleTypes)
0071 {
0072 index = iScoringVolume * fkNumberKinematicRegions * fkNumberParticleTypes
0073 + iKinematicRegion * fkNumberParticleTypes + iParticleType;
0074 }
0075 if (index < 0 || index >= fkNumberCombinations) {
0076 G4cerr << "SteppingAction::GetIndex : WRONG index=" << index << " set it to 0 !" << G4endl;
0077 index = 0;
0078 }
0079 return index;
0080 }
0081
0082
0083
0084 SteppingAction::SteppingAction() : G4UserSteppingAction()
0085 {
0086 Initialize();
0087 }
0088
0089
0090
0091 void SteppingAction::Initialize()
0092 {
0093
0094 fPrimaryParticleId = 0;
0095 fPrimaryParticleEnergy = 0.0;
0096 fPrimaryParticleDirection = G4ThreeVector(0.0, 0.0, 1.0);
0097 fTargetMaterialName = "";
0098 fIsFirstStepOfTheEvent = true;
0099 fIsFirstStepInTarget = true;
0100 fIsFirstStepInScoringUpDown = true;
0101 fIsFirstStepInScoringSide = true;
0102 fCubicVolumeScoringUpDown = 1.0;
0103 fCubicVolumeScoringSide = 1.0;
0104 for (G4int i = 0; i < fkNumberCombinations; ++i) {
0105 fArraySumStepLengths[i] = 0.0;
0106 }
0107
0108
0109
0110
0111
0112
0113
0114
0115
0116
0117
0118
0119
0120
0121
0122
0123
0124
0125 }
0126
0127
0128
0129 void SteppingAction::UserSteppingAction(const G4Step* theStep)
0130 {
0131
0132 if (fIsFirstStepOfTheEvent) {
0133 if (theStep->GetTrack()->GetParentID() == 0) {
0134 fPrimaryParticleId = theStep->GetTrack()->GetDefinition()->GetPDGEncoding();
0135 fPrimaryParticleEnergy = theStep->GetPreStepPoint()->GetKineticEnergy();
0136 fPrimaryParticleDirection = theStep->GetPreStepPoint()->GetMomentumDirection();
0137 if (fRunPtr) {
0138 fRunPtr->SetPrimaryParticleId(fPrimaryParticleId);
0139 fRunPtr->SetPrimaryParticleEnergy(fPrimaryParticleEnergy);
0140 fRunPtr->SetPrimaryParticleDirection(fPrimaryParticleDirection);
0141 }
0142 fIsFirstStepOfTheEvent = false;
0143 }
0144 }
0145
0146 if (fIsFirstStepInTarget
0147 && theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiLayer")
0148 {
0149 fTargetMaterialName = theStep->GetPreStepPoint()->GetMaterial()->GetName();
0150 if (fRunPtr) fRunPtr->SetTargetMaterialName(fTargetMaterialName);
0151 fIsFirstStepInTarget = false;
0152 }
0153
0154 G4int iScoringVolume = -1;
0155 if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringDownstream") {
0156 iScoringVolume = 0;
0157 if (fIsFirstStepInScoringUpDown) {
0158 fCubicVolumeScoringUpDown =
0159 theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
0160 if (fRunPtr) fRunPtr->SetCubicVolumeScoringUpDown(fCubicVolumeScoringUpDown);
0161 fIsFirstStepInScoringUpDown = false;
0162 }
0163 }
0164 else if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringSide") {
0165 iScoringVolume = 1;
0166 if (fIsFirstStepInScoringSide) {
0167 fCubicVolumeScoringSide =
0168 theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
0169 if (fRunPtr) fRunPtr->SetCubicVolumeScoringSide(fCubicVolumeScoringSide);
0170 fIsFirstStepInScoringSide = false;
0171 }
0172 }
0173 else if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringUpstream") {
0174 iScoringVolume = 2;
0175 if (fIsFirstStepInScoringUpDown) {
0176 fCubicVolumeScoringUpDown =
0177 theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
0178 if (fRunPtr) fRunPtr->SetCubicVolumeScoringUpDown(fCubicVolumeScoringUpDown);
0179 fIsFirstStepInScoringUpDown = false;
0180 }
0181 }
0182 if (iScoringVolume >= 0) {
0183
0184
0185
0186 if (iScoringVolume == 2
0187 && fPrimaryParticleDirection.dot(theStep->GetPreStepPoint()->GetMomentumDirection()) > 0.0)
0188 {
0189 return;
0190 }
0191 G4double stepLength = theStep->GetTrack()->GetStepLength() * theStep->GetTrack()->GetWeight();
0192 G4int absPdg = theStep->GetTrack()->GetDefinition() == nullptr
0193 ? 0
0194 : std::abs(theStep->GetTrack()->GetDefinition()->GetPDGEncoding());
0195
0196
0197
0198
0199
0200
0201
0202
0203
0204
0205
0206
0207
0208
0209
0210 G4int iKinematicRegion = theStep->GetPreStepPoint()->GetKineticEnergy() < 20.0 ? 1 : 2;
0211 G4int iParticleType = -1;
0212 if (absPdg == 11)
0213 iParticleType = 1;
0214 else if (absPdg == 22)
0215 iParticleType = 2;
0216 else if (absPdg == 13)
0217 iParticleType = 3;
0218 else if (absPdg == 12 || absPdg == 14 || absPdg == 16)
0219 iParticleType = 4;
0220
0221 else if (absPdg == 111 || absPdg == 211)
0222 iParticleType = 5;
0223 else if (absPdg == 2112)
0224 iParticleType = 6;
0225 else if (absPdg == 2212)
0226 iParticleType = 7;
0227 else if (G4IonTable::IsIon(theStep->GetTrack()->GetDefinition()) ||
0228 G4IonTable::IsAntiIon(theStep->GetTrack()->GetDefinition()))
0229 iParticleType = 8;
0230 else if (absPdg < 1000)
0231 iParticleType = 9;
0232
0233 else if (absPdg > 1000)
0234 iParticleType = 10;
0235
0236
0237 G4int index = GetIndex(iScoringVolume, iKinematicRegion, iParticleType);
0238 fArraySumStepLengths[index] += stepLength;
0239
0240 index = GetIndex(iScoringVolume, iKinematicRegion, 0);
0241 fArraySumStepLengths[index] += stepLength;
0242
0243 index = GetIndex(iScoringVolume, 0, iParticleType);
0244 fArraySumStepLengths[index] += stepLength;
0245
0246 index = GetIndex(iScoringVolume, 0, 0);
0247 fArraySumStepLengths[index] += stepLength;
0248 if (fRunPtr) fRunPtr->SetSteppingArray(fArraySumStepLengths);
0249 }
0250 }
0251
0252