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

#include <G4Step.hh>

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 ()
 
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 60 of file G4Step.hh.

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}

Referenced by G4Step(), and operator=().

◆ ~G4Step()

G4Step::~G4Step ( )

Definition at line 47 of file G4Step.cc.

48{
49 delete fpPreStepPoint;
50 delete fpPostStepPoint;
51
52 delete secondaryInCurrentStep;
53 delete fSecondary;
54}

◆ G4Step() [2/2]

G4Step::G4Step ( const G4Step & right)

Definition at line 57 of file G4Step.cc.

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

Member Function Documentation

◆ AddNonIonizingEnergyDeposit()

◆ AddTotalEnergyDeposit()

◆ ClearFirstStepFlag()

void G4Step::ClearFirstStepFlag ( )
inline

◆ ClearLastStepFlag()

void G4Step::ClearLastStepFlag ( )
inline

◆ CopyPostToPreStepPoint()

void G4Step::CopyPostToPreStepPoint ( )
inline

◆ DeleteSecondaryVector()

void G4Step::DeleteSecondaryVector ( )
inline

◆ GetControlFlag()

G4SteppingControl G4Step::GetControlFlag ( ) const
inline

◆ GetDeltaEnergy()

G4double G4Step::GetDeltaEnergy ( ) const
inline

◆ GetDeltaMomentum()

G4ThreeVector G4Step::GetDeltaMomentum ( ) const
inline

◆ GetDeltaPosition()

G4ThreeVector G4Step::GetDeltaPosition ( ) const
inline

◆ GetDeltaTime()

◆ GetfSecondary()

G4TrackVector * G4Step::GetfSecondary ( )
inline

◆ GetNonIonizingEnergyDeposit()

◆ GetNumberOfSecondariesInCurrentStep()

std::size_t G4Step::GetNumberOfSecondariesInCurrentStep ( ) const
inline

◆ GetPointerToVectorOfAuxiliaryPoints()

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

◆ GetPostStepPoint()

G4StepPoint * G4Step::GetPostStepPoint ( ) const
inline

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

◆ GetPreStepPoint()

G4StepPoint * G4Step::GetPreStepPoint ( ) const
inline

Referenced by G4SDChargedFilter::Accept(), G4SDKineticEnergyFilter::Accept(), G4SDNeutralFilter::Accept(), G4VAtomDeexcitation::AlongStepDeexcitation(), G4ITTransportation::AlongStepDoIt(), G4NuclearStopping::AlongStepDoIt(), G4Transportation::AlongStepDoIt(), G4ClonedRichTrajectory::AppendStep(), G4RichTrajectory::AppendStep(), G4eplusAnnihilation::AtRestDoIt(), G4VReadOutGeometry::CheckROVolume(), G4VPrimitiveScorer::ComputeCurrentSolid(), G4DNABrownianTransportation::ComputeGeomLimit(), G4VMscModel::ComputeGeomLimit(), G4NIELCalculator::ComputeNIEL(), G4VPrimitiveScorer::ComputeSolid(), G4DNABrownianTransportation::ComputeStep(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4LowEWentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanAdjointMscModel::ComputeTruePathLengthLimit(), G4UrbanMscModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4ParallelWorldProcess::CopyStep(), G4AdjointCrossSurfChecker::CrossingAnInterfaceBetweenTwoVolumes(), G4AdjointCrossSurfChecker::CrossingASphere(), G4CellScoreComposer::EstimatorCalculation(), G4DNASmoluchowskiReactionModel::FindReaction(), G4VReadOutGeometry::FindROTouchable(), G4ClonedRichTrajectoryPoint::G4ClonedRichTrajectoryPoint(), G4RichTrajectoryPoint::G4RichTrajectoryPoint(), 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(), G4VPrimitiveScorer::GetIndex(), G4ElNeutrinoNucleusProcess::GetMeanFreePath(), G4MuNeutrinoNucleusProcess::GetMeanFreePath(), G4NeutrinoElectronProcess::GetMeanFreePath(), G4NuVacOscProcess::GetMeanFreePath(), G4TauNeutrinoNucleusProcess::GetMeanFreePath(), G4Scintillation::GetScintillationYieldByParticleType(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolume(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolumeByExtSurface(), G4VFastSimSensitiveDetector::Hit(), G4VGFlashSensitiveDetector::Hit(), G4PSPassageCellCurrent::IsPassed(), G4PSPassageCellFlux::IsPassed(), G4PSPassageTrackLength::IsPassed(), G4PSCylinderSurfaceCurrent::IsSelectedSurface(), G4PSCylinderSurfaceFlux::IsSelectedSurface(), G4PSFlatSurfaceCurrent::IsSelectedSurface(), G4PSFlatSurfaceFlux::IsSelectedSurface(), G4PSSphereSurfaceCurrent::IsSelectedSurface(), G4PSSphereSurfaceFlux::IsSelectedSurface(), G4ElectronIonPair::MeanNumberOfIonsAlongStep(), G4Cerenkov::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4ImportanceProcess::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4NuVacOscProcess::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4VTransitionRadiation::PostStepDoIt(), 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(), G4EnergySplitter::SplitEnergyInVolumes(), G4ParallelWorldProcess::StartTracking(), G4ScoreSplittingProcess::StartTracking(), G4ParticleChange::UpdateStepForAlongStep(), G4ParticleChangeForLoss::UpdateStepForAlongStep(), G4ParticleChangeForTransport::UpdateStepForAlongStep(), G4VParticleChange::UpdateStepForAlongStep(), G4MSSteppingAction::UserSteppingAction(), G4ParallelWorldScoringProcess::Verbose(), and G4ScoreSplittingProcess::Verbose().

◆ GetSecondary()

const G4TrackVector * G4Step::GetSecondary ( ) const
inline

◆ GetSecondaryInCurrentStep()

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

Definition at line 163 of file G4Step.cc.

164{
165 secondaryInCurrentStep->clear();
166 std::size_t nSecondary = fSecondary->size();
167 for(std::size_t i = nSecondaryByLastStep; i < nSecondary; ++i)
168 {
169 secondaryInCurrentStep->push_back((*fSecondary)[i]);
170 }
171 return secondaryInCurrentStep;
172}

Referenced by G4NIELCalculator::RecoilEnergy().

◆ GetStepLength()

G4double G4Step::GetStepLength ( ) const
inline

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

◆ GetTotalEnergyDeposit()

◆ GetTrack()

◆ InitializeStep()

void G4Step::InitializeStep ( G4Track * aValue)
inline

◆ IsFirstStepInVolume()

G4bool G4Step::IsFirstStepInVolume ( ) const
inline

◆ IsLastStepInVolume()

G4bool G4Step::IsLastStepInVolume ( ) const
inline

◆ NewSecondaryVector()

G4TrackVector * G4Step::NewSecondaryVector ( )
inline

◆ operator=()

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

Definition at line 100 of file G4Step.cc.

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

◆ ResetNonIonizingEnergyDeposit()

void G4Step::ResetNonIonizingEnergyDeposit ( )
inline

◆ ResetPostStepPoint()

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

◆ ResetPreStepPoint()

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

◆ ResetTotalEnergyDeposit()

void G4Step::ResetTotalEnergyDeposit ( )
inline

◆ SetControlFlag()

void G4Step::SetControlFlag ( G4SteppingControl StepControlFlag)
inline

◆ SetFirstStepFlag()

void G4Step::SetFirstStepFlag ( )
inline

◆ SetLastStepFlag()

void G4Step::SetLastStepFlag ( )
inline

◆ SetNonIonizingEnergyDeposit()

void G4Step::SetNonIonizingEnergyDeposit ( G4double value)
inline

◆ SetPointerToVectorOfAuxiliaryPoints()

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

◆ SetPostStepPoint()

void G4Step::SetPostStepPoint ( G4StepPoint * value)
inline

◆ SetPreStepPoint()

void G4Step::SetPreStepPoint ( G4StepPoint * value)
inline

◆ SetSecondary()

void G4Step::SetSecondary ( G4TrackVector * value)
inline

◆ SetStepLength()

◆ SetTotalEnergyDeposit()

void G4Step::SetTotalEnergyDeposit ( G4double value)
inline

◆ SetTrack()

void G4Step::SetTrack ( G4Track * value)
inline

◆ UpdateTrack()

void G4Step::UpdateTrack ( )
inline

Member Data Documentation

◆ fNonIonizingEnergyDeposit

G4double G4Step::fNonIonizingEnergyDeposit = 0.0
protected

Definition at line 182 of file G4Step.hh.

Referenced by G4Step(), and operator=().

◆ fTotalEnergyDeposit

G4double G4Step::fTotalEnergyDeposit = 0.0
protected

Definition at line 179 of file G4Step.hh.

Referenced by G4Step(), and operator=().


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