File indexing completed on 2025-11-06 10:05:07
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
0035 inline G4StepPoint* G4Step::GetPreStepPoint() const
0036 {
0037 return fpPreStepPoint;
0038 }
0039
0040 inline void G4Step::SetPreStepPoint(G4StepPoint* value)
0041 {
0042 delete fpPreStepPoint;
0043 fpPreStepPoint = value;
0044 }
0045
0046 inline G4StepPoint* G4Step::ResetPreStepPoint(G4StepPoint* value)
0047 {
0048 G4StepPoint* previousSP = fpPreStepPoint;
0049 fpPreStepPoint = value;
0050 return previousSP;
0051 }
0052
0053 inline G4StepPoint* G4Step::GetPostStepPoint() const
0054 {
0055 return fpPostStepPoint;
0056 }
0057
0058 inline void G4Step::SetPostStepPoint(G4StepPoint* value)
0059 {
0060 delete fpPostStepPoint;
0061 fpPostStepPoint = value;
0062 }
0063
0064 inline G4StepPoint* G4Step::ResetPostStepPoint(G4StepPoint* value)
0065 {
0066 G4StepPoint* previousSP = fpPostStepPoint;
0067 fpPostStepPoint = value;
0068 return previousSP;
0069 }
0070
0071 inline G4double G4Step::GetStepLength() const
0072 {
0073 return fStepLength;
0074 }
0075
0076 inline void G4Step::SetStepLength(G4double value)
0077 {
0078 fStepLength = value;
0079 }
0080
0081 inline G4ThreeVector G4Step::GetDeltaPosition() const
0082 {
0083 return fpPostStepPoint->GetPosition() - fpPreStepPoint->GetPosition();
0084 }
0085
0086 inline G4double G4Step::GetDeltaTime() const
0087 {
0088 return fpPostStepPoint->GetLocalTime() - fpPreStepPoint->GetLocalTime();
0089 }
0090
0091 inline G4double G4Step::GetTotalEnergyDeposit() const
0092 {
0093 return fTotalEnergyDeposit;
0094 }
0095
0096 inline void G4Step::SetTotalEnergyDeposit(G4double value)
0097 {
0098 fTotalEnergyDeposit = value;
0099 }
0100
0101 inline G4double G4Step::GetNonIonizingEnergyDeposit() const
0102 {
0103 return fNonIonizingEnergyDeposit;
0104 }
0105
0106 inline void G4Step::SetNonIonizingEnergyDeposit(G4double value)
0107 {
0108 fNonIonizingEnergyDeposit = value;
0109 }
0110
0111 inline void G4Step::AddTotalEnergyDeposit(G4double value)
0112 {
0113 fTotalEnergyDeposit += value;
0114 }
0115
0116 inline void G4Step::ResetTotalEnergyDeposit()
0117 {
0118 fTotalEnergyDeposit = 0.;
0119 fNonIonizingEnergyDeposit = 0.;
0120 }
0121
0122 inline void G4Step::AddNonIonizingEnergyDeposit(G4double value)
0123 {
0124 fNonIonizingEnergyDeposit += value;
0125 }
0126
0127 inline void G4Step::ResetNonIonizingEnergyDeposit()
0128 {
0129 fNonIonizingEnergyDeposit = 0.;
0130 }
0131
0132 inline G4double G4Step::GetDeltaEnergy() const
0133 {
0134 return fpPostStepPoint->GetKineticEnergy() - fpPreStepPoint->GetKineticEnergy();
0135 }
0136
0137 inline G4ThreeVector G4Step::GetDeltaMomentum() const
0138 {
0139 return fpPostStepPoint->GetMomentum() - fpPreStepPoint->GetMomentum();
0140 }
0141
0142 inline void G4Step::SetControlFlag(G4SteppingControl value)
0143 {
0144 fpSteppingControlFlag = value;
0145 }
0146
0147 inline G4SteppingControl G4Step::GetControlFlag() const
0148 {
0149 return fpSteppingControlFlag;
0150 }
0151
0152 inline void G4Step::CopyPostToPreStepPoint()
0153 {
0154
0155 *(fpPreStepPoint) = *(fpPostStepPoint);
0156 fpPostStepPoint->SetStepStatus(fUndefined);
0157
0158
0159 nSecondaryByLastStep = fSecondary->size();
0160 }
0161
0162
0163
0164
0165
0166
0167 #include "G4Track.hh"
0168
0169 inline G4Track* G4Step::GetTrack() const
0170 {
0171 return fpTrack;
0172 }
0173
0174 inline void G4Step::SetTrack(G4Track* value)
0175 {
0176 fpTrack = value;
0177 }
0178
0179 inline void G4Step::InitializeStep(G4Track* aValue)
0180 {
0181 fpTrack = aValue;
0182 fpTrack->SetStepLength(0.);
0183
0184
0185 fStepLength = 0.;
0186 fTotalEnergyDeposit = 0.;
0187 fNonIonizingEnergyDeposit = 0.;
0188
0189 nSecondaryByLastStep = 0;
0190
0191
0192
0193
0194
0195 fpPreStepPoint->SetSafety(0.0);
0196 fpPreStepPoint->SetProcessDefinedStep(nullptr);
0197 fpPreStepPoint->SetStepStatus(fUndefined);
0198
0199
0200
0201 fpPreStepPoint->SetWeight(fpTrack->GetWeight());
0202 fpPreStepPoint->SetPosition(fpTrack->GetPosition());
0203 fpPreStepPoint->SetLocalTime(fpTrack->GetLocalTime());
0204 fpPreStepPoint->SetGlobalTime(fpTrack->GetGlobalTime());
0205
0206
0207
0208 auto pParticle = fpTrack->GetDynamicParticle();
0209
0210 fpPreStepPoint->SetMass(pParticle->GetMass());
0211 fpPreStepPoint->SetCharge(pParticle->GetCharge());
0212 fpPreStepPoint->SetProperTime(pParticle->GetProperTime());
0213 fpPreStepPoint->SetKineticEnergy(pParticle->GetKineticEnergy());
0214 fpPreStepPoint->SetMomentumDirection(pParticle->GetMomentumDirection());
0215 fpPreStepPoint->SetPolarization(pParticle->GetPolarization());
0216
0217
0218
0219 auto lvol = fpTrack->GetTouchable()->GetVolume()->GetLogicalVolume();
0220
0221 fpPreStepPoint->SetTouchableHandle(fpTrack->GetTouchableHandle());
0222
0223 fpPreStepPoint->SetMaterial(lvol->GetMaterial());
0224 fpPreStepPoint->SetMaterialCutsCouple(lvol->GetMaterialCutsCouple());
0225 fpPreStepPoint->SetSensitiveDetector(lvol->GetSensitiveDetector());
0226
0227
0228 fpPreStepPoint->SetVelocity(fpTrack->CalculateVelocity());
0229
0230 (*fpPostStepPoint) = (*fpPreStepPoint);
0231 }
0232
0233 inline void G4Step::UpdateTrack()
0234 {
0235
0236
0237
0238 fpTrack->SetPosition(fpPostStepPoint->GetPosition());
0239 fpTrack->SetGlobalTime(fpPostStepPoint->GetGlobalTime());
0240 fpTrack->SetLocalTime(fpPostStepPoint->GetLocalTime());
0241
0242 auto* pParticle =
0243 (G4DynamicParticle*) (fpTrack->GetDynamicParticle());
0244 pParticle->SetMass(fpPostStepPoint->GetMass());
0245 pParticle->SetCharge(fpPostStepPoint->GetCharge());
0246 pParticle->SetProperTime(fpPostStepPoint->GetProperTime());
0247
0248 pParticle->SetMomentumDirection(fpPostStepPoint->GetMomentumDirection());
0249 pParticle->SetKineticEnergy(fpPostStepPoint->GetKineticEnergy());
0250 pParticle->SetPolarization(fpPostStepPoint->GetPolarization());
0251
0252 fpTrack->SetStepLength(fStepLength);
0253
0254
0255 fpTrack->SetNextTouchableHandle(fpPostStepPoint->GetTouchableHandle());
0256 fpTrack->SetWeight(fpPostStepPoint->GetWeight());
0257
0258
0259 fpTrack->SetVelocity(fpPostStepPoint->GetVelocity());
0260 }
0261
0262 inline std::size_t G4Step::GetNumberOfSecondariesInCurrentStep() const
0263 {
0264 return fSecondary->size() - nSecondaryByLastStep;
0265 }
0266
0267 inline const G4TrackVector* G4Step::GetSecondary() const
0268 {
0269 return fSecondary;
0270 }
0271
0272 inline G4TrackVector* G4Step::GetfSecondary()
0273 {
0274 return fSecondary;
0275 }
0276
0277 inline void G4Step::SetSecondary(G4TrackVector* value)
0278 {
0279 fSecondary = value;
0280 }
0281
0282 inline G4TrackVector* G4Step::NewSecondaryVector()
0283 {
0284 fSecondary = new G4TrackVector();
0285 return fSecondary;
0286 }
0287
0288 inline void G4Step::DeleteSecondaryVector()
0289 {
0290 if(fSecondary != nullptr)
0291 {
0292 fSecondary->clear();
0293 delete fSecondary;
0294 fSecondary = nullptr;
0295 }
0296 }
0297
0298 inline G4bool G4Step::IsFirstStepInVolume() const
0299 {
0300 return fFirstStepInVolume;
0301 }
0302
0303 inline G4bool G4Step::IsLastStepInVolume() const
0304 {
0305 return fLastStepInVolume;
0306 }
0307
0308 inline void G4Step::SetFirstStepFlag()
0309 {
0310 fFirstStepInVolume = true;
0311 }
0312
0313 inline void G4Step::ClearFirstStepFlag()
0314 {
0315 fFirstStepInVolume = false;
0316 }
0317
0318 inline void G4Step::SetLastStepFlag()
0319 {
0320 fLastStepInVolume = true;
0321 }
0322
0323 inline void G4Step::ClearLastStepFlag()
0324 {
0325 fLastStepInVolume = false;
0326 }
0327
0328 inline void
0329 G4Step::SetPointerToVectorOfAuxiliaryPoints(std::vector<G4ThreeVector>* vec)
0330 {
0331 fpVectorOfAuxiliaryPointsPointer = vec;
0332 }
0333
0334 inline std::vector<G4ThreeVector>*
0335 G4Step::GetPointerToVectorOfAuxiliaryPoints() const
0336 {
0337 return fpVectorOfAuxiliaryPointsPointer;
0338 }