File indexing completed on 2025-04-10 08:06:18
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::fkNumberKinematicRegions>
0052 SteppingAction::fkArrayKinematicRegionNames = {"", "below 20 MeV", "above 20 MeV"};
0053
0054 const std::array<G4String, SteppingAction::fkNumberScoringPositions>
0055 SteppingAction::fkArrayScoringPositionNames = {"forward", "backward"};
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 iKinematicRegion, const G4int iScoringPosition,
0065 const G4int iParticleType)
0066 {
0067 G4int index = -1;
0068 if (iKinematicRegion >= 0 && iKinematicRegion < fkNumberKinematicRegions && iScoringPosition >= 0
0069 && iScoringPosition < fkNumberScoringPositions && iParticleType >= 0
0070 && iParticleType < fkNumberParticleTypes)
0071 {
0072 index = iKinematicRegion * fkNumberScoringPositions * fkNumberParticleTypes
0073 + iScoringPosition * 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 fIsFirstStepInScoringShell = true;
0101 fCubicVolumeScoringShell = 1.0;
0102 for (G4int i = 0; i < fkNumberCombinations; ++i) {
0103 fArraySumStepLengths[i] = 0.0;
0104 }
0105
0106
0107
0108
0109
0110
0111
0112
0113
0114
0115
0116
0117
0118
0119
0120
0121
0122
0123 }
0124
0125
0126
0127 void SteppingAction::UserSteppingAction(const G4Step* theStep)
0128 {
0129
0130 if (fIsFirstStepOfTheEvent) {
0131 if (theStep->GetTrack()->GetParentID() == 0) {
0132 fPrimaryParticleId = theStep->GetTrack()->GetDefinition()->GetPDGEncoding();
0133 fPrimaryParticleEnergy = theStep->GetPreStepPoint()->GetKineticEnergy();
0134 fPrimaryParticleDirection = theStep->GetPreStepPoint()->GetMomentumDirection();
0135 if (fRunPtr) {
0136 fRunPtr->SetPrimaryParticleId(fPrimaryParticleId);
0137 fRunPtr->SetPrimaryParticleEnergy(fPrimaryParticleEnergy);
0138 fRunPtr->SetPrimaryParticleDirection(fPrimaryParticleDirection);
0139 }
0140 fIsFirstStepOfTheEvent = false;
0141 }
0142 }
0143
0144 if (fIsFirstStepInTarget
0145 && theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiSphere")
0146 {
0147 fTargetMaterialName = theStep->GetPreStepPoint()->GetMaterial()->GetName();
0148 if (fRunPtr) fRunPtr->SetTargetMaterialName(fTargetMaterialName);
0149 fIsFirstStepInTarget = false;
0150 }
0151
0152 if (theStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() == "physiScoringShell") {
0153 if (fIsFirstStepInScoringShell) {
0154 fCubicVolumeScoringShell =
0155 theStep->GetTrack()->GetVolume()->GetLogicalVolume()->GetSolid()->GetCubicVolume();
0156 if (fRunPtr) fRunPtr->SetCubicVolumeScoringShell(fCubicVolumeScoringShell);
0157 fIsFirstStepInScoringShell = false;
0158 }
0159 G4double stepLength = theStep->GetTrack()->GetStepLength() * theStep->GetTrack()->GetWeight();
0160 G4int absPdg = theStep->GetTrack()->GetDefinition() == nullptr
0161 ? 0
0162 : std::abs(theStep->GetTrack()->GetDefinition()->GetPDGEncoding());
0163
0164
0165
0166
0167
0168
0169
0170
0171
0172
0173
0174
0175 G4int iKinematicRegion = theStep->GetPreStepPoint()->GetKineticEnergy() < 20.0 ? 1 : 2;
0176
0177
0178 G4int iScoringPosition =
0179 fPrimaryParticleDirection.dot(theStep->GetTrack()->GetPosition().unit()) > 0.0 ? 0 : 1;
0180 G4int iParticleType = -1;
0181 if (absPdg == 11)
0182 iParticleType = 1;
0183 else if (absPdg == 22)
0184 iParticleType = 2;
0185 else if (absPdg == 13)
0186 iParticleType = 3;
0187 else if (absPdg == 12 || absPdg == 14 || absPdg == 16)
0188 iParticleType = 4;
0189
0190 else if (absPdg == 111 || absPdg == 211)
0191 iParticleType = 5;
0192 else if (absPdg == 2112)
0193 iParticleType = 6;
0194 else if (absPdg == 2212)
0195 iParticleType = 7;
0196 else if (G4IonTable::IsIon(theStep->GetTrack()->GetDefinition()) ||
0197 G4IonTable::IsAntiIon(theStep->GetTrack()->GetDefinition()))
0198 iParticleType = 8;
0199 else if (absPdg < 1000)
0200 iParticleType = 9;
0201
0202 else if (absPdg > 1000)
0203 iParticleType = 10;
0204
0205
0206 G4int index = GetIndex(iKinematicRegion, iScoringPosition, iParticleType);
0207 fArraySumStepLengths[index] += stepLength;
0208
0209 index = GetIndex(iKinematicRegion, iScoringPosition, 0);
0210 fArraySumStepLengths[index] += stepLength;
0211
0212 index = GetIndex(0, iScoringPosition, iParticleType);
0213 fArraySumStepLengths[index] += stepLength;
0214
0215 index = GetIndex(0, iScoringPosition, 0);
0216 fArraySumStepLengths[index] += stepLength;
0217 if (fRunPtr) fRunPtr->SetSteppingArray(fArraySumStepLengths);
0218 }
0219 }
0220
0221