Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Step Class Reference

#include <G4Step.hh>

Public Types

using ProfilerConfig = G4ProfilerConfig< G4ProfileType::Step >
 

Public Member Functions

 G4Step ()
 
 ~G4Step ()
 
 G4Step (const G4Step &)
 
G4Stepoperator= (const G4Step &)
 
G4TrackGetTrack () const
 
void SetTrack (G4Track *value)
 
G4StepPointGetPreStepPoint () const
 
void SetPreStepPoint (G4StepPoint *value)
 
G4StepPointResetPreStepPoint (G4StepPoint *value=nullptr)
 
G4StepPointGetPostStepPoint () const
 
void SetPostStepPoint (G4StepPoint *value)
 
G4StepPointResetPostStepPoint (G4StepPoint *value=nullptr)
 
G4double GetStepLength () const
 
void SetStepLength (G4double value)
 
G4double GetTotalEnergyDeposit () const
 
void SetTotalEnergyDeposit (G4double value)
 
G4double GetNonIonizingEnergyDeposit () const
 
void SetNonIonizingEnergyDeposit (G4double value)
 
G4SteppingControl GetControlFlag () const
 
void SetControlFlag (G4SteppingControl StepControlFlag)
 
void AddTotalEnergyDeposit (G4double value)
 
void ResetTotalEnergyDeposit ()
 
void AddNonIonizingEnergyDeposit (G4double value)
 
void ResetNonIonizingEnergyDeposit ()
 
G4bool IsFirstStepInVolume () const
 
G4bool IsLastStepInVolume () const
 
void SetFirstStepFlag ()
 
void ClearFirstStepFlag ()
 
void SetLastStepFlag ()
 
void ClearLastStepFlag ()
 
G4ThreeVector GetDeltaPosition () const
 
G4double GetDeltaTime () const
 
G4ThreeVector GetDeltaMomentum () const
 
G4double GetDeltaEnergy () const
 
void InitializeStep (G4Track *aValue)
 
void UpdateTrack ()
 
void CopyPostToPreStepPoint ()
 
G4PolylineCreatePolyline () const
 
void SetPointerToVectorOfAuxiliaryPoints (std::vector< G4ThreeVector > *vec)
 
std::vector< G4ThreeVector > * GetPointerToVectorOfAuxiliaryPoints () const
 
std::size_t GetNumberOfSecondariesInCurrentStep () const
 
const std::vector< const G4Track * > * GetSecondaryInCurrentStep () const
 
const G4TrackVectorGetSecondary () const
 
G4TrackVectorGetfSecondary ()
 
G4TrackVectorNewSecondaryVector ()
 
void DeleteSecondaryVector ()
 
void SetSecondary (G4TrackVector *value)
 

Protected Attributes

G4double fTotalEnergyDeposit = 0.0
 
G4double fNonIonizingEnergyDeposit = 0.0
 

Detailed Description

Definition at line 61 of file G4Step.hh.

Member Typedef Documentation

◆ ProfilerConfig

Constructor & Destructor Documentation

◆ G4Step() [1/2]

G4Step::G4Step ( )

Definition at line 38 of file G4Step.cc.

39{
40 fpPreStepPoint = new G4StepPoint();
41 fpPostStepPoint = new G4StepPoint();
42
43 secondaryInCurrentStep = new std::vector<const G4Track*>;
44}

◆ ~G4Step()

G4Step::~G4Step ( )

Definition at line 47 of file G4Step.cc.

48{
49 delete fpPreStepPoint;
50 delete fpPostStepPoint;
51
52 secondaryInCurrentStep->clear();
53 delete secondaryInCurrentStep;
54
55 if(fSecondary != nullptr)
56 {
57 fSecondary->clear();
58 delete fSecondary;
59 }
60}

◆ G4Step() [2/2]

G4Step::G4Step ( const G4Step right)

Definition at line 63 of file G4Step.cc.

66 , fStepLength(right.fStepLength)
67 , fpTrack(right.fpTrack)
68 , fpSteppingControlFlag(right.fpSteppingControlFlag)
69 , fFirstStepInVolume(right.fFirstStepInVolume)
70 , fLastStepInVolume(right.fLastStepInVolume)
71 , nSecondaryByLastStep(right.nSecondaryByLastStep)
72 , secondaryInCurrentStep(right.secondaryInCurrentStep)
73 , fpVectorOfAuxiliaryPointsPointer(right.fpVectorOfAuxiliaryPointsPointer)
74{
75 if(right.fpPreStepPoint != nullptr)
76 {
77 fpPreStepPoint = new G4StepPoint(*(right.fpPreStepPoint));
78 }
79 else
80 {
81 fpPreStepPoint = new G4StepPoint();
82 }
83 if(right.fpPostStepPoint != nullptr)
84 {
85 fpPostStepPoint = new G4StepPoint(*(right.fpPostStepPoint));
86 }
87 else
88 {
89 fpPostStepPoint = new G4StepPoint();
90 }
91
92 if(right.fSecondary != nullptr)
93 {
94 fSecondary = new G4TrackVector(*(right.fSecondary));
95 }
96 else
97 {
98 fSecondary = new G4TrackVector();
99 }
100
101 // secondaryInCurrentStep is cleared
102 secondaryInCurrentStep = new std::vector<const G4Track*>;
103}
std::vector< G4Track * > G4TrackVector
G4double fNonIonizingEnergyDeposit
Definition: G4Step.hh:187
G4double fTotalEnergyDeposit
Definition: G4Step.hh:184

Member Function Documentation

◆ AddNonIonizingEnergyDeposit()

◆ AddTotalEnergyDeposit()

◆ ClearFirstStepFlag()

void G4Step::ClearFirstStepFlag ( )

◆ ClearLastStepFlag()

◆ CopyPostToPreStepPoint()

void G4Step::CopyPostToPreStepPoint ( )

◆ CreatePolyline()

G4Polyline * G4Step::CreatePolyline ( ) const

◆ DeleteSecondaryVector()

void G4Step::DeleteSecondaryVector ( )

◆ GetControlFlag()

G4SteppingControl G4Step::GetControlFlag ( ) const

◆ GetDeltaEnergy()

G4double G4Step::GetDeltaEnergy ( ) const

Definition at line 185 of file G4Step.cc.

186{
187 static G4ThreadLocal G4bool isFirstTime = true;
188 if(isFirstTime)
189 {
190 isFirstTime = false;
191#ifdef G4VERBOSE
192 G4Exception("G4Step::GetDeltaEnergy()", "Warning", JustWarning,
193 "This method is obsolete and will be removed soon");
194#endif
195 }
196
197 return fpPostStepPoint->GetKineticEnergy() -
198 fpPreStepPoint->GetKineticEnergy();
199}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
bool G4bool
Definition: G4Types.hh:86
G4double GetKineticEnergy() const
#define G4ThreadLocal
Definition: tls.hh:77

◆ GetDeltaMomentum()

G4ThreeVector G4Step::GetDeltaMomentum ( ) const

Definition at line 169 of file G4Step.cc.

170{
171 static G4ThreadLocal G4bool isFirstTime = true;
172 if(isFirstTime)
173 {
174 isFirstTime = false;
175#ifdef G4VERBOSE
176 G4Exception("G4Step::GetDeltaMomentum()", "Warning", JustWarning,
177 "This method is obsolete and will be removed soon");
178#endif
179 }
180
181 return fpPostStepPoint->GetMomentum() - fpPreStepPoint->GetMomentum();
182}
G4ThreeVector GetMomentum() const

◆ GetDeltaPosition()

G4ThreeVector G4Step::GetDeltaPosition ( ) const

◆ GetDeltaTime()

◆ GetfSecondary()

◆ GetNonIonizingEnergyDeposit()

◆ GetNumberOfSecondariesInCurrentStep()

std::size_t G4Step::GetNumberOfSecondariesInCurrentStep ( ) const

◆ GetPointerToVectorOfAuxiliaryPoints()

std::vector< G4ThreeVector > * G4Step::GetPointerToVectorOfAuxiliaryPoints ( ) const
inline

◆ GetPostStepPoint()

G4StepPoint * G4Step::GetPostStepPoint ( ) const

Referenced by G4VAtomDeexcitation::AlongStepDeexcitation(), G4AdjointAlongStepWeightCorrection::AlongStepDoIt(), G4ContinuousGainOfEnergy::AlongStepDoIt(), G4VMultipleScattering::AlongStepDoIt(), G4NuclearStopping::AlongStepDoIt(), G4AdjointForcedInteractionForGamma::AlongStepDoIt(), G4Transportation::AlongStepDoIt(), G4RichTrajectory::AppendStep(), G4SmoothTrajectory::AppendStep(), G4Trajectory::AppendStep(), G4BOptnForceFreeFlight::ApplyFinalStateBiasing(), G4DecayWithSpin::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4NIELCalculator::ComputeNIEL(), G4DNABrownianTransportation::ComputeStep(), G4VITSteppingVerbose::CopyState(), G4ParallelWorldProcess::CopyStep(), G4AdjointCrossSurfChecker::CrossingAGivenRegisteredSurface(), G4AdjointCrossSurfChecker::CrossingAnInterfaceBetweenTwoVolumes(), G4AdjointCrossSurfChecker::CrossingASphere(), G4ITStepProcessor::DoDefinePhysicalStepLength(), G4ITStepProcessor::DoStepping(), G4ImportanceProcess::G4ImportanceProcess(), G4ParallelWorldProcess::G4ParallelWorldProcess(), G4ParallelWorldScoringProcess::G4ParallelWorldScoringProcess(), G4RichTrajectoryPoint::G4RichTrajectoryPoint(), G4ScoreSplittingProcess::G4ScoreSplittingProcess(), G4SteppingManager::G4SteppingManager(), G4WeightCutOffProcess::G4WeightCutOffProcess(), G4WeightWindowProcess::G4WeightWindowProcess(), G4Scintillation::GetScintillationYieldByParticleType(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolume(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolumeByExtSurface(), G4ParticleChangeForMSC::InitialiseMSC(), G4ITStepProcessor::InvokeAtRestDoItProcs(), G4ITStepProcessor::InvokePSDIP(), G4PSPassageCellCurrent::IsPassed(), G4PSPassageCellFlux::IsPassed(), G4PSPassageTrackLength::IsPassed(), G4PSFlatSurfaceCurrent::IsSelectedSurface(), G4PSFlatSurfaceFlux::IsSelectedSurface(), G4PSSphereSurfaceCurrent::IsSelectedSurface(), G4PSSphereSurfaceFlux::IsSelectedSurface(), G4PSCylinderSurfaceCurrent::IsSelectedSurface(), G4PSCylinderSurfaceFlux::IsSelectedSurface(), G4ErrorPropagator::MakeOneStep(), G4ImportanceProcess::PostStepDoIt(), G4WeightCutOffProcess::PostStepDoIt(), G4WeightWindowProcess::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4Channeling::PostStepDoIt(), G4PhononReflection::PostStepDoIt(), G4PhononScattering::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNAbsorption::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4UCNMultiScattering::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4ITSteppingVerbose::PostStepVerbose(), G4PSCellCharge::ProcessHits(), G4PSCylinderSurfaceFlux::ProcessHits(), G4PSFlatSurfaceFlux::ProcessHits(), G4PSNofCollision::ProcessHits(), G4PSSphereSurfaceFlux::ProcessHits(), G4PSTrackCounter::ProcessHits(), G4PSVolumeFlux::ProcessHits(), G4ElectronIonPair::ResidualeChargePostStep(), G4ElectronIonPair::SampleIonsAlongStep(), G4NeutronGeneralProcess::SelectedProcess(), G4GammaGeneralProcess::SelectedProcess(), G4ITStepProcessor::SetupMembers(), G4ITSteppingVerbose::ShowStep(), G4SteppingVerbose::ShowStep(), G4SteppingVerboseWithUnits::ShowStep(), G4ParallelWorldProcess::StartTracking(), G4ScoreSplittingProcess::StartTracking(), G4ITSteppingVerbose::StepInfo(), G4SteppingVerbose::StepInfo(), G4SteppingVerboseWithUnits::StepInfo(), G4ITSteppingVerbose::StepInfoForLeadingTrack(), G4SteppingManager::Stepping(), G4ParallelWorldProcess::SwitchMaterial(), G4ParticleChangeForOccurenceBiasing::UpdateStepForAlongStep(), G4VParticleChange::UpdateStepForAlongStep(), G4ParticleChangeForLoss::UpdateStepForAlongStep(), G4ParticleChangeForMSC::UpdateStepForAlongStep(), G4ParticleChangeForTransport::UpdateStepForAlongStep(), G4ParticleChange::UpdateStepForAlongStep(), G4ParticleChangeForGamma::UpdateStepForAtRest(), G4FastStep::UpdateStepForAtRest(), G4VParticleChange::UpdateStepForAtRest(), G4ParticleChangeForDecay::UpdateStepForAtRest(), G4ParticleChange::UpdateStepForAtRest(), G4ParticleChangeForOccurenceBiasing::UpdateStepForPostStep(), G4FastStep::UpdateStepForPostStep(), G4VParticleChange::UpdateStepForPostStep(), G4ParticleChangeForDecay::UpdateStepForPostStep(), G4ParticleChangeForGamma::UpdateStepForPostStep(), G4ParticleChangeForLoss::UpdateStepForPostStep(), G4ParticleChangeForTransport::UpdateStepForPostStep(), G4ParticleChange::UpdateStepForPostStep(), G4AdjointSteppingAction::UserSteppingAction(), G4ParallelWorldScoringProcess::Verbose(), and G4ScoreSplittingProcess::Verbose().

◆ GetPreStepPoint()

G4StepPoint * G4Step::GetPreStepPoint ( ) const

Referenced by G4SDChargedFilter::Accept(), G4SDKineticEnergyFilter::Accept(), G4SDNeutralFilter::Accept(), G4VAtomDeexcitation::AlongStepDeexcitation(), G4NuclearStopping::AlongStepDoIt(), G4ITTransportation::AlongStepDoIt(), G4Transportation::AlongStepDoIt(), G4TransportationWithMsc::AlongStepGetPhysicalInteractionLength(), G4RichTrajectory::AppendStep(), G4VReadOutGeometry::CheckROVolume(), G4VPrimitiveScorer::ComputeCurrentSolid(), G4VMscModel::ComputeGeomLimit(), G4DNABrownianTransportation::ComputeGeomLimit(), G4NIELCalculator::ComputeNIEL(), G4VPrimitiveScorer::ComputeSolid(), G4DNABrownianTransportation::ComputeStep(), G4UrbanAdjointMscModel::ComputeTruePathLengthLimit(), G4LowEWentzelVIModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4UrbanMscModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4VITSteppingVerbose::CopyState(), G4ParallelWorldProcess::CopyStep(), G4AdjointCrossSurfChecker::CrossingAnInterfaceBetweenTwoVolumes(), G4AdjointCrossSurfChecker::CrossingASphere(), G4CellScoreComposer::EstimatorCalculation(), G4DNASmoluchowskiReactionModel::FindReaction(), G4VReadOutGeometry::FindROTouchable(), G4ImportanceProcess::G4ImportanceProcess(), G4ParallelWorldProcess::G4ParallelWorldProcess(), G4ParallelWorldScoringProcess::G4ParallelWorldScoringProcess(), G4RichTrajectoryPoint::G4RichTrajectoryPoint(), G4ScoreSplittingProcess::G4ScoreSplittingProcess(), G4SteppingManager::G4SteppingManager(), G4WeightCutOffProcess::G4WeightCutOffProcess(), G4WeightWindowProcess::G4WeightWindowProcess(), G4VPrimitiveScorer::GetIndex(), G4PSCellCharge3D::GetIndex(), G4PSCellFlux3D::GetIndex(), G4PSCylinderSurfaceCurrent3D::GetIndex(), G4PSCylinderSurfaceFlux3D::GetIndex(), G4PSDoseDeposit3D::GetIndex(), G4PSEnergyDeposit3D::GetIndex(), G4PSFlatSurfaceCurrent3D::GetIndex(), G4PSFlatSurfaceFlux3D::GetIndex(), G4PSMinKinEAtGeneration3D::GetIndex(), G4PSNofCollision3D::GetIndex(), G4PSNofSecondary3D::GetIndex(), G4PSNofStep3D::GetIndex(), G4PSPassageCellCurrent3D::GetIndex(), G4PSPassageCellFlux3D::GetIndex(), G4PSPassageTrackLength3D::GetIndex(), G4PSPopulation3D::GetIndex(), G4PSSphereSurfaceCurrent3D::GetIndex(), G4PSSphereSurfaceFlux3D::GetIndex(), G4PSStepChecker3D::GetIndex(), G4PSTermination3D::GetIndex(), G4PSTrackCounter3D::GetIndex(), G4PSTrackLength3D::GetIndex(), G4PSVolumeFlux3D::GetIndex(), G4NeutrinoElectronProcess::GetMeanFreePath(), G4ElNeutrinoNucleusProcess::GetMeanFreePath(), G4MuNeutrinoNucleusProcess::GetMeanFreePath(), G4TauNeutrinoNucleusProcess::GetMeanFreePath(), G4Scintillation::GetScintillationYieldByParticleType(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolume(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolumeByExtSurface(), G4VFastSimSensitiveDetector::Hit(), G4VGFlashSensitiveDetector::Hit(), G4ITStepProcessor::InitDefineStep(), G4PSPassageCellCurrent::IsPassed(), G4PSPassageCellFlux::IsPassed(), G4PSPassageTrackLength::IsPassed(), G4PSFlatSurfaceCurrent::IsSelectedSurface(), G4PSFlatSurfaceFlux::IsSelectedSurface(), G4PSSphereSurfaceCurrent::IsSelectedSurface(), G4PSSphereSurfaceFlux::IsSelectedSurface(), G4PSCylinderSurfaceCurrent::IsSelectedSurface(), G4PSCylinderSurfaceFlux::IsSelectedSurface(), G4ElectronIonPair::MeanNumberOfIonsAlongStep(), G4ImportanceProcess::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4VTransitionRadiation::PostStepDoIt(), G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength(), G4AdjointForcedInteractionForGamma::PostStepGetPhysicalInteractionLength(), G4PSCellCharge::ProcessHits(), G4PSCellFlux::ProcessHits(), G4PSCylinderSurfaceCurrent::ProcessHits(), G4PSCylinderSurfaceFlux::ProcessHits(), G4PSDoseDeposit::ProcessHits(), G4PSEnergyDeposit::ProcessHits(), G4PSFlatSurfaceCurrent::ProcessHits(), G4PSFlatSurfaceFlux::ProcessHits(), G4PSMinKinEAtGeneration::ProcessHits(), G4PSNofCollision::ProcessHits(), G4PSNofSecondary::ProcessHits(), G4PSPassageCellCurrent::ProcessHits(), G4PSPassageCellFlux::ProcessHits(), G4PSPopulation::ProcessHits(), G4PSSphereSurfaceCurrent::ProcessHits(), G4PSSphereSurfaceFlux::ProcessHits(), G4PSTermination::ProcessHits(), G4PSTrackCounter::ProcessHits(), G4PSTrackLength::ProcessHits(), G4PSVolumeFlux::ProcessHits(), G4ErrorFreeTrajState::PropagateError(), G4TransportationLogger::ReportLoopingTrack(), G4ElectronIonPair::SampleIonsAlongStep(), G4ITStepProcessor::SetupMembers(), G4ITSteppingVerbose::ShowStep(), G4SteppingVerbose::ShowStep(), G4SteppingVerboseWithUnits::ShowStep(), G4EnergySplitter::SplitEnergyInVolumes(), G4ParallelWorldProcess::StartTracking(), G4ScoreSplittingProcess::StartTracking(), G4SteppingManager::Stepping(), G4VParticleChange::UpdateStepForAlongStep(), G4ParticleChangeForLoss::UpdateStepForAlongStep(), G4ParticleChangeForTransport::UpdateStepForAlongStep(), G4ParticleChange::UpdateStepForAlongStep(), G4MSSteppingAction::UserSteppingAction(), G4ParallelWorldScoringProcess::Verbose(), and G4ScoreSplittingProcess::Verbose().

◆ GetSecondary()

◆ GetSecondaryInCurrentStep()

const std::vector< const G4Track * > * G4Step::GetSecondaryInCurrentStep ( ) const

Definition at line 202 of file G4Step.cc.

203{
204 secondaryInCurrentStep->clear();
205 std::size_t nSecondary = fSecondary->size();
206 for(std::size_t i = nSecondaryByLastStep; i < nSecondary; ++i)
207 {
208 secondaryInCurrentStep->push_back((*fSecondary)[i]);
209 }
210 return secondaryInCurrentStep;
211}

Referenced by G4NIELCalculator::RecoilEnergy(), and G4SteppingVerboseWithUnits::StepInfo().

◆ GetStepLength()

G4double G4Step::GetStepLength ( ) const

Referenced by G4VAtomDeexcitation::AlongStepDeexcitation(), G4AdjointAlongStepWeightCorrection::AlongStepDoIt(), G4ContinuousGainOfEnergy::AlongStepDoIt(), G4VEnergyLossProcess::AlongStepDoIt(), G4VMultipleScattering::AlongStepDoIt(), G4ErrorEnergyLoss::AlongStepDoIt(), G4BiasingProcessInterface::AlongStepDoIt(), G4NuclearStopping::AlongStepDoIt(), G4AdjointForcedInteractionForGamma::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4RayTrajectory::AppendStep(), G4NIELCalculator::ComputeNIEL(), G4ParallelWorldProcess::CopyStep(), G4ITStepProcessor::DoStepping(), G4CellScoreComposer::EstimatorCalculation(), G4PSPassageCellFlux::IsPassed(), G4PSPassageTrackLength::IsPassed(), G4ImportanceProcess::PostStepDoIt(), G4WeightWindowProcess::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4DNABrownianTransportation::PostStepDoIt(), G4BiasingProcessInterface::PostStepDoIt(), G4VTransitionRadiation::PostStepDoIt(), G4MultiFunctionalDetector::ProcessHits(), G4PSCellFlux::ProcessHits(), G4PSNofStep::ProcessHits(), G4PSTrackLength::ProcessHits(), G4ErrorFreeTrajState::PropagateError(), G4TransportationLogger::ReportLoopingTrack(), G4EnergySplitter::SplitEnergyInVolumes(), G4ITSteppingVerbose::StepInfo(), G4SteppingVerbose::StepInfo(), G4SteppingVerboseWithUnits::StepInfo(), G4ITSteppingVerbose::StepInfoForLeadingTrack(), G4SteppingManager::Stepping(), G4SteppingVerbose::TrackingStarted(), G4SteppingVerboseWithUnits::TrackingStarted(), G4BOptnForceCommonTruncatedExp::UpdateForStep(), G4MSSteppingAction::UserSteppingAction(), G4ParallelWorldScoringProcess::Verbose(), G4ScoreSplittingProcess::Verbose(), and G4EmSaturation::VisibleEnergyDepositionAtAStep().

◆ GetTotalEnergyDeposit()

◆ GetTrack()

◆ InitializeStep()

void G4Step::InitializeStep ( G4Track aValue)

◆ IsFirstStepInVolume()

G4bool G4Step::IsFirstStepInVolume ( ) const

◆ IsLastStepInVolume()

G4bool G4Step::IsLastStepInVolume ( ) const

◆ NewSecondaryVector()

◆ operator=()

G4Step & G4Step::operator= ( const G4Step right)

Definition at line 106 of file G4Step.cc.

107{
108 if(this != &right)
109 {
112 fStepLength = right.fStepLength;
113 fpTrack = right.fpTrack;
114 fpSteppingControlFlag = right.fpSteppingControlFlag;
115 fFirstStepInVolume = right.fFirstStepInVolume;
116 fLastStepInVolume = right.fLastStepInVolume;
117 nSecondaryByLastStep = right.nSecondaryByLastStep;
118 secondaryInCurrentStep = right.secondaryInCurrentStep;
119 fpVectorOfAuxiliaryPointsPointer = right.fpVectorOfAuxiliaryPointsPointer;
120
121 delete fpPreStepPoint;
122
123 if(right.fpPreStepPoint != nullptr)
124 {
125 fpPreStepPoint = new G4StepPoint(*(right.fpPreStepPoint));
126 }
127 else
128 {
129 fpPreStepPoint = new G4StepPoint();
130 }
131
132 delete fpPostStepPoint;
133
134 if(right.fpPostStepPoint != nullptr)
135 {
136 fpPostStepPoint = new G4StepPoint(*(right.fpPostStepPoint));
137 }
138 else
139 {
140 fpPostStepPoint = new G4StepPoint();
141 }
142
143 if(fSecondary != nullptr)
144 {
145 fSecondary->clear();
146 delete fSecondary;
147 }
148 if(right.fSecondary != nullptr)
149 {
150 fSecondary = new G4TrackVector(*(right.fSecondary));
151 }
152 else
153 {
154 fSecondary = new G4TrackVector();
155 }
156
157 // secondaryInCurrentStep is not copied
158 if(secondaryInCurrentStep != nullptr)
159 {
160 secondaryInCurrentStep->clear();
161 delete secondaryInCurrentStep;
162 }
163 secondaryInCurrentStep = new std::vector<const G4Track*>;
164 }
165 return *this;
166}

◆ ResetNonIonizingEnergyDeposit()

void G4Step::ResetNonIonizingEnergyDeposit ( )

◆ ResetPostStepPoint()

G4StepPoint * G4Step::ResetPostStepPoint ( G4StepPoint value = nullptr)

◆ ResetPreStepPoint()

G4StepPoint * G4Step::ResetPreStepPoint ( G4StepPoint value = nullptr)

◆ ResetTotalEnergyDeposit()

void G4Step::ResetTotalEnergyDeposit ( )

◆ SetControlFlag()

◆ SetFirstStepFlag()

◆ SetLastStepFlag()

◆ SetNonIonizingEnergyDeposit()

void G4Step::SetNonIonizingEnergyDeposit ( G4double  value)

◆ SetPointerToVectorOfAuxiliaryPoints()

void G4Step::SetPointerToVectorOfAuxiliaryPoints ( std::vector< G4ThreeVector > *  vec)
inline

◆ SetPostStepPoint()

void G4Step::SetPostStepPoint ( G4StepPoint value)

◆ SetPreStepPoint()

void G4Step::SetPreStepPoint ( G4StepPoint value)

◆ SetSecondary()

void G4Step::SetSecondary ( G4TrackVector value)

◆ SetStepLength()

◆ SetTotalEnergyDeposit()

◆ SetTrack()

◆ UpdateTrack()

Member Data Documentation

◆ fNonIonizingEnergyDeposit

G4double G4Step::fNonIonizingEnergyDeposit = 0.0
protected

Definition at line 187 of file G4Step.hh.

Referenced by operator=().

◆ fTotalEnergyDeposit

G4double G4Step::fTotalEnergyDeposit = 0.0
protected

Definition at line 184 of file G4Step.hh.

Referenced by operator=().


The documentation for this class was generated from the following files: