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