File indexing completed on 2025-01-18 09:59: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 void G4Step::SetControlFlag(G4SteppingControl value)
0133 {
0134 fpSteppingControlFlag = value;
0135 }
0136
0137 inline G4SteppingControl G4Step::GetControlFlag() const
0138 {
0139 return fpSteppingControlFlag;
0140 }
0141
0142 inline void G4Step::CopyPostToPreStepPoint()
0143 {
0144
0145 *(fpPreStepPoint) = *(fpPostStepPoint);
0146 fpPostStepPoint->SetStepStatus(fUndefined);
0147
0148
0149 nSecondaryByLastStep = fSecondary->size();
0150 }
0151
0152
0153
0154
0155
0156
0157 #include "G4Track.hh"
0158
0159 inline G4Track* G4Step::GetTrack() const
0160 {
0161 return fpTrack;
0162 }
0163
0164 inline void G4Step::SetTrack(G4Track* value)
0165 {
0166 fpTrack = value;
0167 }
0168
0169 inline void G4Step::InitializeStep(G4Track* aValue)
0170 {
0171 fpTrack = aValue;
0172 fpTrack->SetStepLength(0.);
0173
0174
0175 fStepLength = 0.;
0176 fTotalEnergyDeposit = 0.;
0177 fNonIonizingEnergyDeposit = 0.;
0178
0179 nSecondaryByLastStep = 0;
0180
0181
0182
0183
0184
0185 fpPreStepPoint->SetSafety(0.0);
0186 fpPreStepPoint->SetProcessDefinedStep(nullptr);
0187 fpPreStepPoint->SetStepStatus(fUndefined);
0188
0189
0190
0191 fpPreStepPoint->SetWeight(fpTrack->GetWeight());
0192 fpPreStepPoint->SetPosition(fpTrack->GetPosition());
0193 fpPreStepPoint->SetLocalTime(fpTrack->GetLocalTime());
0194 fpPreStepPoint->SetGlobalTime(fpTrack->GetGlobalTime());
0195
0196
0197
0198 auto pParticle = fpTrack->GetDynamicParticle();
0199
0200 fpPreStepPoint->SetMass(pParticle->GetMass());
0201 fpPreStepPoint->SetCharge(pParticle->GetCharge());
0202 fpPreStepPoint->SetProperTime(pParticle->GetProperTime());
0203 fpPreStepPoint->SetKineticEnergy(pParticle->GetKineticEnergy());
0204 fpPreStepPoint->SetMomentumDirection(pParticle->GetMomentumDirection());
0205 fpPreStepPoint->SetPolarization(pParticle->GetPolarization());
0206
0207
0208
0209 auto lvol = fpTrack->GetTouchable()->GetVolume()->GetLogicalVolume();
0210
0211 fpPreStepPoint->SetTouchableHandle(fpTrack->GetTouchableHandle());
0212
0213 fpPreStepPoint->SetMaterial(lvol->GetMaterial());
0214 fpPreStepPoint->SetMaterialCutsCouple(lvol->GetMaterialCutsCouple());
0215 fpPreStepPoint->SetSensitiveDetector(lvol->GetSensitiveDetector());
0216
0217
0218 fpPreStepPoint->SetVelocity(fpTrack->CalculateVelocity());
0219
0220 (*fpPostStepPoint) = (*fpPreStepPoint);
0221 }
0222
0223 inline void G4Step::UpdateTrack()
0224 {
0225
0226
0227
0228 fpTrack->SetPosition(fpPostStepPoint->GetPosition());
0229 fpTrack->SetGlobalTime(fpPostStepPoint->GetGlobalTime());
0230 fpTrack->SetLocalTime(fpPostStepPoint->GetLocalTime());
0231
0232 auto* pParticle =
0233 (G4DynamicParticle*) (fpTrack->GetDynamicParticle());
0234 pParticle->SetMass(fpPostStepPoint->GetMass());
0235 pParticle->SetCharge(fpPostStepPoint->GetCharge());
0236 pParticle->SetProperTime(fpPostStepPoint->GetProperTime());
0237
0238 pParticle->SetMomentumDirection(fpPostStepPoint->GetMomentumDirection());
0239 pParticle->SetKineticEnergy(fpPostStepPoint->GetKineticEnergy());
0240 pParticle->SetPolarization(fpPostStepPoint->GetPolarization());
0241
0242 fpTrack->SetStepLength(fStepLength);
0243
0244
0245 fpTrack->SetNextTouchableHandle(fpPostStepPoint->GetTouchableHandle());
0246 fpTrack->SetWeight(fpPostStepPoint->GetWeight());
0247
0248
0249 fpTrack->SetVelocity(fpPostStepPoint->GetVelocity());
0250 }
0251
0252 inline std::size_t G4Step::GetNumberOfSecondariesInCurrentStep() const
0253 {
0254 return fSecondary->size() - nSecondaryByLastStep;
0255 }
0256
0257 inline const G4TrackVector* G4Step::GetSecondary() const
0258 {
0259 return fSecondary;
0260 }
0261
0262 inline G4TrackVector* G4Step::GetfSecondary()
0263 {
0264 return fSecondary;
0265 }
0266
0267 inline void G4Step::SetSecondary(G4TrackVector* value)
0268 {
0269 fSecondary = value;
0270 }
0271
0272 inline G4TrackVector* G4Step::NewSecondaryVector()
0273 {
0274 fSecondary = new G4TrackVector();
0275 return fSecondary;
0276 }
0277
0278 inline void G4Step::DeleteSecondaryVector()
0279 {
0280 if(fSecondary != nullptr)
0281 {
0282 fSecondary->clear();
0283 delete fSecondary;
0284 fSecondary = nullptr;
0285 }
0286 }
0287
0288 inline G4bool G4Step::IsFirstStepInVolume() const
0289 {
0290 return fFirstStepInVolume;
0291 }
0292
0293 inline G4bool G4Step::IsLastStepInVolume() const
0294 {
0295 return fLastStepInVolume;
0296 }
0297
0298 inline void G4Step::SetFirstStepFlag()
0299 {
0300 fFirstStepInVolume = true;
0301 }
0302
0303 inline void G4Step::ClearFirstStepFlag()
0304 {
0305 fFirstStepInVolume = false;
0306 }
0307
0308 inline void G4Step::SetLastStepFlag()
0309 {
0310 fLastStepInVolume = true;
0311 }
0312
0313 inline void G4Step::ClearLastStepFlag()
0314 {
0315 fLastStepInVolume = false;
0316 }
0317
0318 inline void
0319 G4Step::SetPointerToVectorOfAuxiliaryPoints(std::vector<G4ThreeVector>* vec)
0320 {
0321 fpVectorOfAuxiliaryPointsPointer = vec;
0322 }
0323
0324 inline std::vector<G4ThreeVector>*
0325 G4Step::GetPointerToVectorOfAuxiliaryPoints() const
0326 {
0327 return fpVectorOfAuxiliaryPointsPointer;
0328 }