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