File indexing completed on 2025-04-10 08:06:16
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 fAbsorberMaterialName = "";
0098 fActiveMaterialName = "";
0099 fIsFirstStepOfTheEvent = true;
0100 fIsFirstStepInAbsorberLayer = true;
0101 fIsFirstStepInActiveLayer = true;
0102 fIsFirstStepInScoringUpDown = true;
0103 fIsFirstStepInScoringSide = true;
0104 fCubicVolumeScoringUpDown = 1.0;
0105 fCubicVolumeScoringSide = 1.0;
0106 for (G4int i = 0; i < fkNumberCombinations; ++i) {
0107 fArraySumStepLengths[i] = 0.0;
0108 }
0109
0110
0111
0112
0113
0114
0115
0116
0117
0118
0119
0120
0121
0122
0123
0124
0125
0126
0127 }
0128
0129
0130
0131 void SteppingAction::UserSteppingAction(const G4Step* theStep)
0132 {
0133
0134 if (fIsFirstStepOfTheEvent) {
0135 if (theStep->GetTrack()->GetParentID() == 0) {
0136 fPrimaryParticleId = theStep->GetTrack()->GetDefinition()->GetPDGEncoding();
0137 fPrimaryParticleEnergy = theStep->GetPreStepPoint()->GetKineticEnergy();
0138 fPrimaryParticleDirection = theStep->GetPreStepPoint()->GetMomentumDirection();
0139 if (fRunPtr) {
0140 fRunPtr->SetPrimaryParticleId(fPrimaryParticleId);
0141 fRunPtr->SetPrimaryParticleEnergy(fPrimaryParticleEnergy);
0142 fRunPtr->SetPrimaryParticleDirection(fPrimaryParticleDirection);
0143 }
0144 fIsFirstStepOfTheEvent = false;
0145 }
0146 }
0147
0148 if (fIsFirstStepInAbsorberLayer
0149 && theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiAbsorber")
0150 {
0151 fAbsorberMaterialName = theStep->GetPreStepPoint()->GetMaterial()->GetName();
0152 if (fRunPtr) fRunPtr->SetAbsorberMaterialName(fAbsorberMaterialName);
0153 fIsFirstStepInAbsorberLayer = false;
0154 }
0155 if (fIsFirstStepInActiveLayer
0156 && theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiActive")
0157 {
0158 fActiveMaterialName = theStep->GetPreStepPoint()->GetMaterial()->GetName();
0159 if (fRunPtr) fRunPtr->SetActiveMaterialName(fActiveMaterialName);
0160 fIsFirstStepInActiveLayer = false;
0161 }
0162
0163 G4int iScoringVolume = -1;
0164 if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringDownstream") {
0165 iScoringVolume = 0;
0166 if (fIsFirstStepInScoringUpDown) {
0167 fCubicVolumeScoringUpDown =
0168 theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
0169 if (fRunPtr) fRunPtr->SetCubicVolumeScoringUpDown(fCubicVolumeScoringUpDown);
0170 fIsFirstStepInScoringUpDown = false;
0171 }
0172 }
0173 else if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringSide") {
0174 iScoringVolume = 1;
0175 if (fIsFirstStepInScoringSide) {
0176 fCubicVolumeScoringSide =
0177 theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
0178 if (fRunPtr) fRunPtr->SetCubicVolumeScoringSide(fCubicVolumeScoringSide);
0179 fIsFirstStepInScoringSide = false;
0180 }
0181 }
0182 else if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringUpstream") {
0183 iScoringVolume = 2;
0184 if (fIsFirstStepInScoringUpDown) {
0185 fCubicVolumeScoringUpDown =
0186 theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
0187 if (fRunPtr) fRunPtr->SetCubicVolumeScoringUpDown(fCubicVolumeScoringUpDown);
0188 fIsFirstStepInScoringUpDown = false;
0189 }
0190 }
0191 if (iScoringVolume >= 0) {
0192
0193
0194
0195 if (iScoringVolume == 2
0196 && fPrimaryParticleDirection.dot(theStep->GetPreStepPoint()->GetMomentumDirection()) > 0.0)
0197 return;
0198 G4double stepLength = theStep->GetTrack()->GetStepLength() * theStep->GetTrack()->GetWeight();
0199 G4int absPdg = theStep->GetTrack()->GetDefinition() == nullptr
0200 ? 0
0201 : std::abs(theStep->GetTrack()->GetDefinition()->GetPDGEncoding());
0202
0203
0204
0205
0206
0207
0208
0209
0210
0211
0212
0213
0214
0215
0216
0217 G4int iKinematicRegion = theStep->GetPreStepPoint()->GetKineticEnergy() < 20.0 ? 1 : 2;
0218 G4int iParticleType = -1;
0219 if (absPdg == 11)
0220 iParticleType = 1;
0221 else if (absPdg == 22)
0222 iParticleType = 2;
0223 else if (absPdg == 13)
0224 iParticleType = 3;
0225 else if (absPdg == 12 || absPdg == 14 || absPdg == 16)
0226 iParticleType = 4;
0227
0228 else if (absPdg == 111 || absPdg == 211)
0229 iParticleType = 5;
0230 else if (absPdg == 2112)
0231 iParticleType = 6;
0232 else if (absPdg == 2212)
0233 iParticleType = 7;
0234 else if (G4IonTable::IsIon(theStep->GetTrack()->GetDefinition()) ||
0235 G4IonTable::IsAntiIon(theStep->GetTrack()->GetDefinition()))
0236 iParticleType = 8;
0237 else if (absPdg < 1000)
0238 iParticleType = 9;
0239
0240 else if (absPdg > 1000)
0241 iParticleType = 10;
0242
0243
0244 G4int index = GetIndex(iScoringVolume, iKinematicRegion, iParticleType);
0245 fArraySumStepLengths[index] += stepLength;
0246
0247 index = GetIndex(iScoringVolume, iKinematicRegion, 0);
0248 fArraySumStepLengths[index] += stepLength;
0249
0250 index = GetIndex(iScoringVolume, 0, iParticleType);
0251 fArraySumStepLengths[index] += stepLength;
0252
0253 index = GetIndex(iScoringVolume, 0, 0);
0254 fArraySumStepLengths[index] += stepLength;
0255 if (fRunPtr) fRunPtr->SetSteppingArray(fArraySumStepLengths);
0256 }
0257 }
0258
0259