Geant4-11
Public Types | Public Member Functions | Private Member Functions | Private Attributes
G4Track Class Reference

#include <G4Track.hh>

Public Types

using ProfilerConfig = G4ProfilerConfig< G4ProfileType::Track >
 

Public Member Functions

void AddTrackLength (const G4double aValue)
 
G4double CalculateVelocity () const
 
G4double CalculateVelocityForOpticalPhoton () const
 
void CopyTrackInfo (const G4Track &)
 
 G4Track ()
 
 G4Track (const G4Track &)
 
 G4Track (G4DynamicParticle *apValueDynamicParticle, G4double aValueTime, const G4ThreeVector &aValuePosition)
 
G4VAuxiliaryTrackInformationGetAuxiliaryTrackInformation (G4int id) const
 
std::map< G4int, G4VAuxiliaryTrackInformation * > * GetAuxiliaryTrackInformationMap () const
 
G4int GetCreatorModelID () const
 
G4int GetCreatorModelIndex () const
 
const G4String GetCreatorModelName () const
 
const G4VProcessGetCreatorProcess () const
 
G4int GetCurrentStepNumber () const
 
G4ParticleDefinitionGetDefinition () const
 
const G4DynamicParticleGetDynamicParticle () const
 
G4double GetGlobalTime () const
 
G4double GetKineticEnergy () const
 
G4double GetLocalTime () const
 
const G4LogicalVolumeGetLogicalVolumeAtVertex () const
 
G4MaterialGetMaterial () const
 
const G4MaterialCutsCoupleGetMaterialCutsCouple () const
 
G4ThreeVector GetMomentum () const
 
const G4ThreeVectorGetMomentumDirection () const
 
G4MaterialGetNextMaterial () const
 
const G4MaterialCutsCoupleGetNextMaterialCutsCouple () const
 
const G4VTouchableGetNextTouchable () const
 
const G4TouchableHandleGetNextTouchableHandle () const
 
G4VPhysicalVolumeGetNextVolume () const
 
const G4VTouchableGetOriginTouchable () const
 
const G4TouchableHandleGetOriginTouchableHandle () const
 
G4int GetParentID () const
 
const G4ParticleDefinitionGetParticleDefinition () const
 
const G4ThreeVectorGetPolarization () const
 
const G4ThreeVectorGetPosition () const
 
G4double GetProperTime () const
 
const G4StepGetStep () const
 
G4double GetStepLength () const
 
G4double GetTotalEnergy () const
 
const G4VTouchableGetTouchable () const
 
const G4TouchableHandleGetTouchableHandle () const
 
G4int GetTrackID () const
 
G4double GetTrackLength () const
 
G4TrackStatus GetTrackStatus () const
 
G4VUserTrackInformationGetUserInformation () const
 
G4double GetVelocity () const
 
G4double GetVertexKineticEnergy () const
 
const G4ThreeVectorGetVertexMomentumDirection () const
 
const G4ThreeVectorGetVertexPosition () const
 
G4VPhysicalVolumeGetVolume () const
 
G4double GetWeight () const
 
void IncrementCurrentStepNumber ()
 
G4bool IsBelowThreshold () const
 
G4bool IsGoodForTracking () const
 
void operator delete (void *aTrack)
 
void * operator new (std::size_t)
 
G4bool operator!= (const G4Track &)
 
G4Trackoperator= (const G4Track &)
 
G4bool operator== (const G4Track &)
 
void RemoveAuxiliaryTrackInformation (G4int id)
 
void RemoveAuxiliaryTrackInformation (G4String &name)
 
void SetAuxiliaryTrackInformation (G4int id, G4VAuxiliaryTrackInformation *info) const
 
void SetBelowThresholdFlag (G4bool value=true)
 
void SetCreatorModelID (const G4int id)
 
void SetCreatorProcess (const G4VProcess *aValue)
 
void SetGlobalTime (const G4double aValue)
 
void SetGoodForTrackingFlag (G4bool value=true)
 
void SetKineticEnergy (const G4double aValue)
 
void SetLocalTime (const G4double aValue)
 
void SetLogicalVolumeAtVertex (const G4LogicalVolume *)
 
void SetMomentumDirection (const G4ThreeVector &aValue)
 
void SetNextTouchableHandle (const G4TouchableHandle &apValue)
 
void SetOriginTouchableHandle (const G4TouchableHandle &apValue)
 
void SetParentID (const G4int aValue)
 
void SetPolarization (const G4ThreeVector &aValue)
 
void SetPosition (const G4ThreeVector &aValue)
 
void SetProperTime (const G4double aValue)
 
void SetStep (const G4Step *aValue)
 
void SetStepLength (G4double value)
 
void SetTouchableHandle (const G4TouchableHandle &apValue)
 
void SetTrackID (const G4int aValue)
 
void SetTrackStatus (const G4TrackStatus aTrackStatus)
 
void SetUserInformation (G4VUserTrackInformation *aValue) const
 
void SetVelocity (G4double val)
 
void SetVertexKineticEnergy (const G4double aValue)
 
void SetVertexMomentumDirection (const G4ThreeVector &aValue)
 
void SetVertexPosition (const G4ThreeVector &aValue)
 
void SetWeight (G4double aValue)
 
G4bool UseGivenVelocity () const
 
void UseGivenVelocity (G4bool val)
 
 ~G4Track ()
 

Private Member Functions

void ClearAuxiliaryTrackInformation ()
 

Private Attributes

G4bool fBelowThreshold = false
 
G4int fCreatorModelID = -1
 
G4int fCurrentStepNumber = 0
 
G4double fGlobalTime = 0.0
 
G4bool fGoodForTracking = false
 
G4double fLocalTime = 0.0
 
G4int fParentID = 0
 
std::map< G4int, G4VAuxiliaryTrackInformation * > * fpAuxiliaryTrackInformationMap = nullptr
 
const G4VProcessfpCreatorProcess = nullptr
 
G4DynamicParticlefpDynamicParticle = nullptr
 
const G4LogicalVolumefpLVAtVertex = nullptr
 
G4TouchableHandle fpNextTouchable
 
G4TouchableHandle fpOriginTouchable
 
G4ThreeVector fPosition
 
const G4StepfpStep = nullptr
 
G4TouchableHandle fpTouchable
 
G4VUserTrackInformationfpUserInformation = nullptr
 
G4double fStepLength = 0.0
 
G4int fTrackID = 0
 
G4double fTrackLength = 0.0
 
G4TrackStatus fTrackStatus = fAlive
 
G4double fVelocity = 0.0
 
G4double fVtxKineticEnergy = 0.0
 
G4ThreeVector fVtxMomentumDirection
 
G4ThreeVector fVtxPosition
 
G4double fWeight = 1.0
 
G4MaterialPropertyVectorgroupvel = nullptr
 
G4bool is_OpticalPhoton = false
 
G4Materialprev_mat = nullptr
 
G4double prev_momentum = 0.0
 
G4double prev_velocity = 0.0
 
G4bool useGivenVelocity = false
 

Detailed Description

Definition at line 66 of file G4Track.hh.

Member Typedef Documentation

◆ ProfilerConfig

Definition at line 70 of file G4Track.hh.

Constructor & Destructor Documentation

◆ G4Track() [1/3]

G4Track::G4Track ( )

Definition at line 62 of file G4Track.cc.

65{}
G4double fVelocity
Definition: G4Track.hh:271
G4DynamicParticle * fpDynamicParticle
Definition: G4Track.hh:278
float c_light
Definition: hepunit.py:256

◆ G4Track() [2/3]

G4Track::G4Track ( G4DynamicParticle apValueDynamicParticle,
G4double  aValueTime,
const G4ThreeVector aValuePosition 
)

Definition at line 47 of file G4Track.cc.

50 : fPosition(aValuePosition)
51 , fGlobalTime(aValueTime)
53{
54 fpDynamicParticle = (apValueDynamicParticle)
55 ? apValueDynamicParticle : new G4DynamicParticle();
56 // check if the particle type is Optical Photon
59}
G4ParticleDefinition * GetDefinition() const
G4ThreeVector fPosition
Definition: G4Track.hh:262
G4double fGlobalTime
Definition: G4Track.hh:264
G4bool is_OpticalPhoton
Definition: G4Track.hh:330

References fpDynamicParticle, G4DynamicParticle::GetDefinition(), G4ParticleDefinition::GetPDGEncoding(), and is_OpticalPhoton.

◆ G4Track() [3/3]

G4Track::G4Track ( const G4Track right)

Definition at line 68 of file G4Track.cc.

70{
71 *this = right;
72}

◆ ~G4Track()

G4Track::~G4Track ( )

Definition at line 75 of file G4Track.cc.

76{
77 delete fpDynamicParticle;
78 delete fpUserInformation;
80}
G4VUserTrackInformation * fpUserInformation
Definition: G4Track.hh:303
void ClearAuxiliaryTrackInformation()
Definition: G4Track.cc:257

References ClearAuxiliaryTrackInformation(), fpDynamicParticle, and fpUserInformation.

Member Function Documentation

◆ AddTrackLength()

void G4Track::AddTrackLength ( const G4double  aValue)

◆ CalculateVelocity()

G4double G4Track::CalculateVelocity ( ) const

◆ CalculateVelocityForOpticalPhoton()

G4double G4Track::CalculateVelocityForOpticalPhoton ( ) const

Definition at line 156 of file G4Track.cc.

157{
158 G4double velocity = c_light;
159
160 G4Material* mat = nullptr;
161 G4bool update_groupvel = false;
162 if(fpStep != nullptr)
163 {
164 mat = this->GetMaterial(); // Fix for repeated volumes
165 }
166 else
167 {
168 if(fpTouchable != 0)
169 {
171 }
172 }
173 // check if previous step is in the same volume
174 // and get new GROUPVELOCITY table if necessary
175 if((mat != nullptr) && ((mat != prev_mat) || (groupvel == nullptr)))
176 {
177 groupvel = nullptr;
178 if(mat->GetMaterialPropertiesTable() != nullptr)
180 update_groupvel = true;
181 }
182 prev_mat = mat;
183
184 if(groupvel != nullptr)
185 {
186 // light velocity = c/(rindex+d(rindex)/d(log(E_phot)))
187 // values stored in GROUPVEL material properties vector
188 velocity = prev_velocity;
189
190 // check if momentum is same as in the previous step
191 // and calculate group velocity if necessary
192 G4double current_momentum = fpDynamicParticle->GetTotalMomentum();
193 if(update_groupvel || (current_momentum != prev_momentum))
194 {
195 velocity = groupvel->Value(current_momentum);
196 prev_velocity = velocity;
197 prev_momentum = current_momentum;
198 }
199 }
200
201 return velocity;
202}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4double GetTotalMomentum() const
G4Material * GetMaterial() const
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
G4double Value(const G4double energy, std::size_t &lastidx) const
G4double prev_velocity
Definition: G4Track.hh:307
const G4Step * fpStep
Definition: G4Track.hh:290
G4double prev_momentum
Definition: G4Track.hh:308
G4Material * GetMaterial() const
G4Material * prev_mat
Definition: G4Track.hh:305
G4TouchableHandle fpTouchable
Definition: G4Track.hh:273
G4MaterialPropertyVector * groupvel
Definition: G4Track.hh:306
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34

References source.hepunit::c_light, fpDynamicParticle, fpStep, fpTouchable, G4VPhysicalVolume::GetLogicalVolume(), G4LogicalVolume::GetMaterial(), GetMaterial(), G4Material::GetMaterialPropertiesTable(), G4MaterialPropertiesTable::GetProperty(), G4DynamicParticle::GetTotalMomentum(), G4VTouchable::GetVolume(), groupvel, kGROUPVEL, prev_mat, prev_momentum, prev_velocity, and G4PhysicsVector::Value().

Referenced by G4ITTransportation::AlongStepDoIt().

◆ ClearAuxiliaryTrackInformation()

void G4Track::ClearAuxiliaryTrackInformation ( )
private

Definition at line 257 of file G4Track.cc.

258{
259 if(fpAuxiliaryTrackInformationMap == nullptr)
260 return;
261 for(auto itr = fpAuxiliaryTrackInformationMap->cbegin();
262 itr != fpAuxiliaryTrackInformationMap->cend(); ++itr)
263 {
264 delete(*itr).second;
265 }
268}
std::map< G4int, G4VAuxiliaryTrackInformation * > * fpAuxiliaryTrackInformationMap
Definition: G4Track.hh:312

References fpAuxiliaryTrackInformationMap.

Referenced by operator=(), and ~G4Track().

◆ CopyTrackInfo()

void G4Track::CopyTrackInfo ( const G4Track right)

Definition at line 150 of file G4Track.cc.

151{
152 *this = right;
153}

◆ GetAuxiliaryTrackInformation()

G4VAuxiliaryTrackInformation * G4Track::GetAuxiliaryTrackInformation ( G4int  id) const

◆ GetAuxiliaryTrackInformationMap()

std::map< G4int, G4VAuxiliaryTrackInformation * > * G4Track::GetAuxiliaryTrackInformationMap ( ) const
inline

◆ GetCreatorModelID()

G4int G4Track::GetCreatorModelID ( ) const
inline

◆ GetCreatorModelIndex()

G4int G4Track::GetCreatorModelIndex ( ) const
inline

◆ GetCreatorModelName()

const G4String G4Track::GetCreatorModelName ( ) const
inline

◆ GetCreatorProcess()

const G4VProcess * G4Track::GetCreatorProcess ( ) const

◆ GetCurrentStepNumber()

G4int G4Track::GetCurrentStepNumber ( ) const

Referenced by G4CoupledTransportation::AlongStepDoIt(), G4ImportanceProcess::AlongStepGetPhysicalInteractionLength(), G4WeightCutOffProcess::AlongStepGetPhysicalInteractionLength(), G4WeightWindowProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldScoringProcess::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4ParallelGeometriesLimiterProcess::AlongStepGetPhysicalInteractionLength(), G4FastSimulationManagerProcess::AlongStepGetPhysicalInteractionLength(), G4RichTrajectory::AppendStep(), G4DNABrownianTransportation::Diffusion(), export_G4Track(), G4RichTrajectoryPoint::G4RichTrajectoryPoint(), G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength(), G4AdjointForcedInteractionForGamma::PostStepGetPhysicalInteractionLength(), G4PSCellCharge::ProcessHits(), G4PSMinKinEAtGeneration::ProcessHits(), G4PSNofSecondary::ProcessHits(), G4BOptrForceCollision::ProposeOccurenceBiasingOperation(), G4TransportationLogger::ReportLoopingTrack(), G4ITStepProcessor::SetInitialStep(), G4SteppingManager::SetInitialStep(), G4ITSteppingVerbose::StepInfo(), G4SteppingVerbose::StepInfo(), G4SteppingVerboseWithUnits::StepInfo(), G4ITSteppingVerbose::StepInfoForLeadingTrack(), G4SteppingVerbose::TrackingStarted(), G4SteppingVerboseWithUnits::TrackingStarted(), G4ITSteppingVerbose::VerboseTrack(), G4SteppingVerbose::VerboseTrack(), and G4SteppingVerboseWithUnits::VerboseTrack().

◆ GetDefinition()

G4ParticleDefinition * G4Track::GetDefinition ( ) const

Referenced by G4SDParticleFilter::Accept(), G4VAtomDeexcitation::AlongStepDeexcitation(), G4AdjointForcedInteractionForGamma::AlongStepDoIt(), G4ITStepProcessor::ApplyProductionCut(), G4SteppingManager::ApplyProductionCut(), G4eplusAnnihilation::AtRestDoIt(), G4FastSimulationManager::AtRestGetFastSimulationManagerTrigger(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4ParticleChange::CheckIt(), G4ParticleChangeForDecay::CheckIt(), G4ParticleChangeForGamma::CheckIt(), G4ParticleChangeForLoss::CheckIt(), G4ParticleChangeForMSC::CheckIt(), GFlashShowerModel::CheckParticleDefAndContainment(), G4VParticleChange::CheckSecondary(), G4StackChecker::ClassifyNewTrack(), G4ITStepProcessor::DealWithSecondaries(), GFlashShowerModel::DoIt(), G4ITStepProcessor::DoStepping(), export_G4Track(), G4MuonicAtomDecay::FillResult(), G4SmoothTrajectory::G4SmoothTrajectory(), G4Trajectory::G4Trajectory(), G4VAdjointReverseReaction::GetMeanFreePath(), G4VTransitionRadiation::GetMeanFreePath(), G4ITStepProcessor::GetProcessInfo(), G4SteppingManager::GetProcessNumber(), G4ErrorPropagator::InitG4Track(), G4SteppingManager::InvokeAlongStepDoItProcs(), G4SteppingManager::InvokeAtRestDoItProcs(), G4SteppingManager::InvokePSDIP(), GFlashShowerModel::ModelTrigger(), G4hImpactIonisation::PostStepDoIt(), G4FastSimulationManager::PostStepGetFastSimulationManagerTrigger(), G4MinEkineCuts::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4LowECapture::PostStepGetPhysicalInteractionLength(), G4PSNofSecondary::ProcessHits(), G4TrackingManager::ProcessOneTrack(), G4ErrorPropagator::Propagate(), G4BOptrForceCollision::ProposeNonPhysicsBiasingOperation(), G4ChannelingOptrChangeCrossSection::ProposeOccurenceBiasingOperation(), G4BOptrForceCollision::ProposeOccurenceBiasingOperation(), G4DecayWithSpin::Spin_Precession(), G4EnergySplitter::SplitEnergyInVolumes(), G4ITTrackingManager::StartTracking(), G4VEnergyLossProcess::StartTracking(), G4TrackingManager::TrackBanner(), G4ITTrackingInteractivity::TrackBanner(), G4VITSteppingVerbose::TrackBanner(), G4ITSteppingVerbose::TrackingEnded(), and G4AdjointSteppingAction::UserSteppingAction().

◆ GetDynamicParticle()

const G4DynamicParticle * G4Track::GetDynamicParticle ( ) const

Referenced by G4AdjointAlongStepWeightCorrection::AlongStepDoIt(), G4ContinuousGainOfEnergy::AlongStepDoIt(), G4VEnergyLossProcess::AlongStepDoIt(), G4ErrorEnergyLoss::AlongStepDoIt(), G4CoupledTransportation::AlongStepDoIt(), G4Transportation::AlongStepDoIt(), G4AdjointProcessEquivalentToDirectProcess::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4VContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4AdjointProcessEquivalentToDirectProcess::AlongStepGetPhysicalInteractionLength(), G4ITTransportation::AlongStepGetPhysicalInteractionLength(), G4EmBiasingManager::ApplyDirectionalSplitting(), G4ITStepProcessor::ApplyProductionCut(), G4SteppingManager::ApplyProductionCut(), G4EmBiasingManager::ApplySplitting(), G4DecayWithSpin::AtRestDoIt(), G4AdjointProcessEquivalentToDirectProcess::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4Decay::AtRestGetPhysicalInteractionLength(), G4AdjointProcessEquivalentToDirectProcess::AtRestGetPhysicalInteractionLength(), G4BraggIonGasModel::ChargeSquareRatio(), G4BetheBlochIonGasModel::ChargeSquareRatio(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4PolarizedAnnihilation::ComputeSaturationFactor(), G4PolarizedCompton::ComputeSaturationFactor(), G4PolarizedIonisation::ComputeSaturationFactor(), G4ITTransportation::ComputeStep(), G4UrbanAdjointMscModel::ComputeTruePathLengthLimit(), G4LowEWentzelVIModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4UrbanMscModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4FieldTrackUpdator::CreateFieldTrack(), G4PionDecayMakeSpin::DaughterPolarization(), G4RadioactiveDecay::DecayAnalog(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4MuonicAtomDecay::DecayIt(), export_G4Track(), G4ErrorEnergyLoss::GetContinuousStepLimit(), G4hImpactIonisation::GetContinuousStepLimit(), G4AnnihiToMuPair::GetMeanFreePath(), G4MuonicAtomDecay::GetMeanFreePath(), G4Decay::GetMeanFreePath(), G4GammaConversionToMuons::GetMeanFreePath(), G4VXTRenergyLoss::GetMeanFreePath(), G4NeutrinoElectronProcess::GetMeanFreePath(), G4HadronicProcess::GetMeanFreePath(), G4ElNeutrinoNucleusProcess::GetMeanFreePath(), G4MuNeutrinoNucleusProcess::GetMeanFreePath(), G4OpAbsorption::GetMeanFreePath(), G4OpMieHG::GetMeanFreePath(), G4OpRayleigh::GetMeanFreePath(), G4OpWLS::GetMeanFreePath(), G4OpWLS2::GetMeanFreePath(), G4RadioactiveDecay::GetMeanFreePath(), G4hImpactIonisation::GetMeanFreePath(), G4SynchrotronRadiation::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetMeanFreePath(), G4MuonicAtomDecay::GetMeanLifeTime(), G4Decay::GetMeanLifeTime(), G4RadioactiveDecay::GetMeanLifeTime(), G4SynchrotronRadiationInMat::GetPhotonEnergy(), G4Scintillation::GetScintillationYieldByParticleType(), G4HadProjectile::Initialise(), G4FastStep::Initialize(), G4ParticleChange::Initialize(), G4ParticleChangeForDecay::Initialize(), G4ParticleChangeForLoss::InitializeForAlongStep(), G4ParticleChangeForLoss::InitializeForPostStep(), G4VEmProcess::MeanFreePath(), G4VEnergyLossProcess::MeanFreePath(), G4DNARuddIonisationExtendedModel::PartialCrossSection(), G4DNARuddIonisationModel::PartialCrossSection(), G4SmartTrackStack::PopFromStack(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4DecayWithSpin::PostStepDoIt(), G4AnnihiToMuPair::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpMieHG::PostStepDoIt(), G4OpRayleigh::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4AdjointProcessEquivalentToDirectProcess::PostStepDoIt(), G4Cerenkov::PostStepGetPhysicalInteractionLength(), G4MaxTimeCuts::PostStepGetPhysicalInteractionLength(), G4MinEkineCuts::PostStepGetPhysicalInteractionLength(), G4VITDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4GammaGeneralProcess::PostStepGetPhysicalInteractionLength(), G4Decay::PostStepGetPhysicalInteractionLength(), G4UnknownDecay::PostStepGetPhysicalInteractionLength(), G4AdjointProcessEquivalentToDirectProcess::PostStepGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4ErrorFreeTrajState::PropagateError(), G4ErrorFreeTrajState::PropagateErrorIoni(), G4ErrorFreeTrajState::PropagateErrorMSC(), G4SmartTrackStack::PushToStack(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointhIonisationModel::RapidSampleSecondaries(), G4AdjointBremsstrahlungModel::SampleSecondaries(), G4AdjointComptonModel::SampleSecondaries(), G4AdjointeIonisationModel::SampleSecondaries(), G4AdjointhIonisationModel::SampleSecondaries(), G4AdjointIonIonisationModel::SampleSecondaries(), G4AdjointPhotoElectricModel::SampleSecondaries(), G4GammaGeneralProcess::SelectHadProcess(), G4SteppingManager::SetInitialStep(), G4AdjointProcessEquivalentToDirectProcess::StartTracking(), G4UrbanAdjointMscModel::StartTracking(), G4GoudsmitSaundersonMscModel::StartTracking(), G4UrbanMscModel::StartTracking(), G4FieldTrackUpdator::Update(), and G4ParticleChangeForTransport::UpdateStepForAlongStep().

◆ GetGlobalTime()

G4double G4Track::GetGlobalTime ( ) const

Referenced by G4ITTrackHolder::_PushTrack(), G4ParticleChangeForGamma::AddSecondary(), G4ITTransportation::AlongStepDoIt(), G4CoupledTransportation::AlongStepDoIt(), G4Transportation::AlongStepDoIt(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4ITTransportation::AlongStepGetPhysicalInteractionLength(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4DecayWithSpin::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4DNAIRTMoleculeEncounterStepper::CheckAndRecordResults(), G4DNAMoleculeEncounterStepper::CheckAndRecordResults(), G4DNAIndependentReactionTimeStepper::CheckAndRecordResults(), G4FastStep::CheckIt(), G4ParticleChange::CheckIt(), G4VParticleChange::CheckSecondary(), G4StackChecker::ClassifyNewTrack(), G4ITStepProcessor::ComputeInteractionLength(), G4ITTransportation::ComputeStep(), G4ITModelProcessor::ComputeTrackReaction(), G4FieldTrackUpdator::CreateFieldTrack(), G4VPhononProcess::CreateSecondary(), G4RadioactiveDecay::DecayAnalog(), G4DNAMolecularDissociation::DecayIt(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4Radioactivation::DecayIt(), G4MuonicAtomDecay::DecayIt(), G4DNABrownianTransportation::Diffusion(), export_G4Track(), G4HadronicProcess::FillResult(), G4MuonicAtomDecay::FillResult(), G4SynchrotronRadiation::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetMeanFreePath(), G4FastStep::Initialize(), G4ParticleChange::Initialize(), G4ParticleChangeForDecay::Initialize(), G4DNAIRT::MakeReaction(), G4DNAIRT_geometries::MakeReaction(), G4DNAMakeReaction::MakeReaction(), G4DNAMolecularReaction::MakeReaction(), G4DNAElectronHoleRecombination::MakeReaction(), G4DNASecondOrderReaction::PostStepDoIt(), G4DNAScavengerProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNAbsorption::PostStepDoIt(), G4UCNMultiScattering::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4DNABrownianTransportation::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4MaxTimeCuts::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4NeutronKiller::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength(), G4ITSteppingVerbose::PreStepVerbose(), G4ITTrackHolder::PushDelayed(), G4TrackingInformation::RecordCurrentPositionNTime(), G4DNAIRT::Sampling(), G4DNAIRT_geometries::Sampling(), G4FieldTrackUpdator::Update(), G4FastStep::UpdateStepForAtRest(), G4FastStep::UpdateStepForPostStep(), G4ITSteppingVerbose::VerboseTrack(), G4SteppingVerbose::VerboseTrack(), G4SteppingVerboseWithUnits::VerboseTrack(), and G4Molecule::~G4Molecule().

◆ GetKineticEnergy()

G4double G4Track::GetKineticEnergy ( ) const

Referenced by G4VMultipleScattering::AlongStepDoIt(), G4ErrorEnergyLoss::AlongStepDoIt(), G4AdjointForcedInteractionForGamma::AlongStepDoIt(), G4ITTransportation::AlongStepDoIt(), G4VMultipleScattering::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4ITTransportation::AlongStepGetPhysicalInteractionLength(), G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4ITStepProcessor::ApplyProductionCut(), G4SteppingManager::ApplyProductionCut(), G4VEmModel::ChargeSquareRatio(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4FastStep::CheckIt(), G4ParticleChange::CheckIt(), G4ParticleChangeForDecay::CheckIt(), G4ParticleChangeForGamma::CheckIt(), G4ParticleChangeForLoss::CheckIt(), G4ParticleChangeForMSC::CheckIt(), G4VParticleChange::CheckSecondary(), G4StackChecker::ClassifyNewTrack(), G4FieldTrackUpdator::CreateFieldTrack(), G4PhysChemIO::FormattedText::CreateSolvatedElectron(), G4PhysChemIO::G4Analysis::CreateSolvatedElectron(), G4ITStepProcessor::DealWithSecondaries(), G4HadronicProcess::DumpState(), G4MuonicAtomDecay::DumpState(), G4ExceptionHandler::DumpTrackInfo(), GFlashShowerModel::ElectronDoIt(), export_G4Track(), G4HadronicProcess::FillResult(), G4MuonicAtomDecay::FillResult(), G4RichTrajectory::G4RichTrajectory(), G4RichTrajectoryPoint::G4RichTrajectoryPoint(), G4SmoothTrajectory::G4SmoothTrajectory(), G4Trajectory::G4Trajectory(), G4ErrorEnergyLoss::GetContinuousStepLimit(), G4AdjointAlongStepWeightCorrection::GetContinuousStepLimit(), G4ContinuousGainOfEnergy::GetContinuousStepLimit(), G4PhononDownconversion::GetMeanFreePath(), G4PhononScattering::GetMeanFreePath(), G4VAdjointReverseReaction::GetMeanFreePath(), G4VTransitionRadiation::GetMeanFreePath(), G4Scintillation::GetScintillationYieldByParticleType(), G4ErrorPropagator::InitG4Track(), G4ParticleChangeForLoss::InitializeForAlongStep(), G4ParticleChangeForGamma::InitializeForPostStep(), G4ParticleChangeForLoss::InitializeForPostStep(), G4ITStepProcessor::InvokeAlongStepDoItProcs(), G4SteppingManager::InvokeAlongStepDoItProcs(), G4SteppingManager::InvokeAtRestDoItProcs(), G4SteppingManager::InvokePSDIP(), G4UCNBoundaryProcess::InvokeSD(), G4PhononDownconversion::MakeLTSecondaries(), G4PhononDownconversion::MakeTTSecondaries(), G4VEmProcess::MeanFreePath(), G4VEnergyLossProcess::MeanFreePath(), GFlashShowerModel::ModelTrigger(), G4SpecialCuts::PostStepDoIt(), G4WeightWindowProcess::PostStepDoIt(), G4PhononReflection::PostStepDoIt(), G4PhononScattering::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4AdjointForcedInteractionForGamma::PostStepDoIt(), G4LowECapture::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutronKiller::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4GammaGeneralProcess::PostStepGetPhysicalInteractionLength(), G4AdjointForcedInteractionForGamma::PostStepGetPhysicalInteractionLength(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4LowECapture::PostStepGetPhysicalInteractionLength(), G4ErrorPropagator::Propagate(), G4TransportationLogger::ReportLoopingTrack(), G4ITStepProcessor::SetInitialStep(), G4SteppingManager::SetInitialStep(), G4ITSteppingVerbose::StepInfo(), G4SteppingVerbose::StepInfo(), G4SteppingVerboseWithUnits::StepInfo(), G4ITSteppingVerbose::StepInfoForLeadingTrack(), G4SteppingVerbose::TrackingStarted(), G4SteppingVerboseWithUnits::TrackingStarted(), G4FieldTrackUpdator::Update(), G4AdjointSteppingAction::UserSteppingAction(), G4ITSteppingVerbose::VerboseTrack(), G4SteppingVerbose::VerboseTrack(), and G4SteppingVerboseWithUnits::VerboseTrack().

◆ GetLocalTime()

G4double G4Track::GetLocalTime ( ) const

◆ GetLogicalVolumeAtVertex()

const G4LogicalVolume * G4Track::GetLogicalVolumeAtVertex ( ) const

Referenced by export_G4Track().

◆ GetMaterial()

G4Material * G4Track::GetMaterial ( ) const

Referenced by G4ErrorEnergyLoss::AlongStepDoIt(), G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4VContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4HadronStoppingProcess::AtRestDoIt(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), CalculateVelocityForOpticalPhoton(), G4VEmModel::ChargeSquareRatio(), G4PolarizedAnnihilation::ComputeSaturationFactor(), G4PolarizedCompton::ComputeSaturationFactor(), G4PolarizedIonisation::ComputeSaturationFactor(), G4DNABrownianTransportation::ComputeStep(), G4DNABrownianTransportation::Diffusion(), G4HadronicProcess::DumpState(), G4MuonicAtomDecay::DumpState(), export_G4Track(), G4DNAElectronHoleRecombination::FindReactant(), G4ErrorEnergyLoss::GetContinuousStepLimit(), G4AnnihiToMuPair::GetMeanFreePath(), G4GammaConversionToMuons::GetMeanFreePath(), G4NeutrinoElectronProcess::GetMeanFreePath(), G4HadronicProcess::GetMeanFreePath(), G4ElNeutrinoNucleusProcess::GetMeanFreePath(), G4MuNeutrinoNucleusProcess::GetMeanFreePath(), G4OpAbsorption::GetMeanFreePath(), G4OpMieHG::GetMeanFreePath(), G4OpRayleigh::GetMeanFreePath(), G4OpWLS::GetMeanFreePath(), G4OpWLS2::GetMeanFreePath(), G4UCNAbsorption::GetMeanFreePath(), G4UCNLoss::GetMeanFreePath(), G4UCNMultiScattering::GetMeanFreePath(), G4Scintillation::GetScintillationYieldByParticleType(), G4HadProjectile::Initialise(), G4NeutrinoElectronProcess::PostStepDoIt(), G4AnnihiToMuPair::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4OpMieHG::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4VTransitionRadiation::PostStepDoIt(), G4Cerenkov::PostStepGetPhysicalInteractionLength(), G4VITDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4Decay::PostStepGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), and G4ElementSelector::SelectZandA().

◆ GetMaterialCutsCouple()

const G4MaterialCutsCouple * G4Track::GetMaterialCutsCouple ( ) const

Referenced by G4VMultipleScattering::AlongStepDoIt(), G4NuclearStopping::AlongStepDoIt(), G4AdjointForcedInteractionForGamma::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4VMultipleScattering::AlongStepGetPhysicalInteractionLength(), G4EmBiasingManager::ApplyDirectionalSplitting(), G4EmBiasingManager::ApplyRangeCut(), G4EmBiasingManager::ApplySplitting(), G4eplusAnnihilation::AtRestDoIt(), G4UrbanAdjointMscModel::ComputeTruePathLengthLimit(), G4LowEWentzelVIModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4UrbanMscModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4hImpactIonisation::GetContinuousStepLimit(), G4AdjointAlongStepWeightCorrection::GetContinuousStepLimit(), G4ContinuousGainOfEnergy::GetContinuousStepLimit(), G4hImpactIonisation::GetMeanFreePath(), G4VAdjointReverseReaction::GetMeanFreePath(), G4VEmProcess::MeanFreePath(), G4VEnergyLossProcess::MeanFreePath(), G4AdjointForcedInteractionForGamma::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4Cerenkov::PostStepGetPhysicalInteractionLength(), G4MinEkineCuts::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4GammaGeneralProcess::PostStepGetPhysicalInteractionLength(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointhIonisationModel::RapidSampleSecondaries(), G4AdjointBremsstrahlungModel::SampleSecondaries(), G4AdjointPhotoElectricModel::SampleSecondaries(), and G4EmSaturation::VisibleEnergyDepositionAtAStep().

◆ GetMomentum()

G4ThreeVector G4Track::GetMomentum ( ) const

◆ GetMomentumDirection()

const G4ThreeVector & G4Track::GetMomentumDirection ( ) const

Referenced by G4DNABrownianTransportation::AlongStepDoIt(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4VParticleChange::CheckSecondary(), G4StackChecker::ClassifyNewTrack(), G4VMscModel::ComputeGeomLimit(), G4DNABrownianTransportation::ComputeGeomLimit(), G4DNABrownianTransportation::ComputeStep(), G4FieldTrackUpdator::CreateFieldTrack(), G4HadronicProcess::DumpState(), G4MuonicAtomDecay::DumpState(), G4ExceptionHandler::DumpTrackInfo(), GFlashShowerModel::ElectronDoIt(), export_G4Track(), G4HadronicProcess::FillResult(), G4VFastSimSensitiveDetector::Hit(), G4VGFlashSensitiveDetector::Hit(), G4BOptnForceCommonTruncatedExp::Initialize(), G4ParticleChangeForGamma::InitializeForPostStep(), G4ParticleChangeForLoss::InitializeForPostStep(), G4PhononReflection::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4VTransitionRadiation::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4Transportation::PostStepDoIt(), G4ITStepProcessor::SetInitialStep(), G4SteppingManager::SetInitialStep(), G4SamplingPostStepAction::Split(), G4ParallelGeometriesLimiterProcess::StartTracking(), G4ImportanceProcess::StartTracking(), G4WeightCutOffProcess::StartTracking(), G4WeightWindowProcess::StartTracking(), G4FastSimulationManagerProcess::StartTracking(), G4ParallelWorldProcess::StartTracking(), G4ParallelWorldScoringProcess::StartTracking(), G4CoupledTransportation::StartTracking(), G4VPhononProcess::StartTracking(), G4FieldTrackUpdator::Update(), G4ParallelWorldScoringProcess::Verbose(), G4ScoreSplittingProcess::Verbose(), G4ITSteppingVerbose::VerboseTrack(), G4SteppingVerbose::VerboseTrack(), and G4SteppingVerboseWithUnits::VerboseTrack().

◆ GetNextMaterial()

G4Material * G4Track::GetNextMaterial ( ) const

◆ GetNextMaterialCutsCouple()

const G4MaterialCutsCouple * G4Track::GetNextMaterialCutsCouple ( ) const

◆ GetNextTouchable()

const G4VTouchable * G4Track::GetNextTouchable ( ) const

◆ GetNextTouchableHandle()

const G4TouchableHandle & G4Track::GetNextTouchableHandle ( ) const

◆ GetNextVolume()

G4VPhysicalVolume * G4Track::GetNextVolume ( ) const

◆ GetOriginTouchable()

const G4VTouchable * G4Track::GetOriginTouchable ( ) const

◆ GetOriginTouchableHandle()

const G4TouchableHandle & G4Track::GetOriginTouchableHandle ( ) const

◆ GetParentID()

G4int G4Track::GetParentID ( ) const

◆ GetParticleDefinition()

const G4ParticleDefinition * G4Track::GetParticleDefinition ( ) const

◆ GetPolarization()

const G4ThreeVector & G4Track::GetPolarization ( ) const

◆ GetPosition()

const G4ThreeVector & G4Track::GetPosition ( ) const

Referenced by G4ParticleChangeForGamma::AddSecondary(), G4DNABrownianTransportation::AlongStepDoIt(), G4AdjointForcedInteractionForGamma::AlongStepDoIt(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4ITTransportation::AlongStepGetPhysicalInteractionLength(), G4EmBiasingManager::ApplyDirectionalSplitting(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4Molecule::BuildTrack(), G4DNAIndependentReactionTimeStepper::CalculateStep(), G4DNAIRTMoleculeEncounterStepper::CheckAndRecordResults(), G4DNAMoleculeEncounterStepper::CheckAndRecordResults(), G4ParticleChange::CheckIt(), G4ParticleChangeForDecay::CheckIt(), G4ParticleChangeForGamma::CheckIt(), G4ParticleChangeForLoss::CheckIt(), G4ParticleChangeForMSC::CheckIt(), G4VParticleChange::CheckSecondary(), G4StackChecker::ClassifyNewTrack(), G4ITTransportation::ComputeStep(), G4DNABrownianTransportation::ComputeStep(), G4FieldTrackUpdator::CreateFieldTrack(), G4VPhononProcess::CreateSecondary(), G4PhysChemIO::FormattedText::CreateSolvatedElectron(), G4PhysChemIO::G4Analysis::CreateSolvatedElectron(), G4DNAChemistryManager::CreateSolvatedElectron(), G4DNAChemistryManager::CreateWaterMolecule(), G4PhysChemIO::FormattedText::CreateWaterMolecule(), G4PhysChemIO::G4Analysis::CreateWaterMolecule(), G4RadioactiveDecay::DecayAnalog(), G4DNAMolecularDissociation::DecayIt(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4Radioactivation::DecayIt(), G4MuonicAtomDecay::DecayIt(), G4HadronicProcess::DumpState(), G4MuonicAtomDecay::DumpState(), GFlashShowerModel::ElectronDoIt(), export_G4Track(), G4HadronicProcess::FillResult(), G4MuonicAtomDecay::FillResult(), G4DNAElectronHoleRecombination::FindReactant(), G4DNASmoluchowskiReactionModel::FindReaction(), G4SmoothTrajectory::G4SmoothTrajectory(), G4Trajectory::G4Trajectory(), G4SynchrotronRadiation::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetPhotonEnergy(), G4GFlashSpot::GetPosition(), G4IT::GetPosition(), G4DNAPartiallyDiffusionControlled::GetTimeToEncounter(), G4DNATotallyDiffusionControlled::GetTimeToEncounter(), G4FastStep::Initialize(), G4ParticleChange::Initialize(), G4BOptnForceCommonTruncatedExp::Initialize(), G4FastSimHitMaker::make(), G4DNAIRT::MakeReaction(), G4DNAIRT_geometries::MakeReaction(), G4DNAMakeReaction::MakeReaction(), G4DNAMolecularReaction::MakeReaction(), G4DNAElectronHoleRecombination::MakeReaction(), G4IT::operator[](), G4DNASecondOrderReaction::PostStepDoIt(), G4DNAScavengerProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4Transportation::PostStepDoIt(), G4ErrorMagFieldLimitProcess::PostStepGetPhysicalInteractionLength(), G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength(), G4ITSteppingVerbose::PostStepVerbose(), G4ITSteppingVerbose::PreStepVerbose(), G4ErrorFreeTrajState::PropagateError(), G4TrackingInformation::RecordCurrentPositionNTime(), G4TransportationLogger::ReportLoopingTrack(), G4TDNAOneStepThermalizationModel< MODEL >::SampleSecondaries(), G4DNAIRT::Sampling(), G4DNAIRT_geometries::Sampling(), G4FastTrack::SetCurrentTrack(), G4ITStepProcessor::SetInitialStep(), G4SteppingManager::SetInitialStep(), G4ParallelGeometriesLimiterProcess::StartTracking(), G4ImportanceProcess::StartTracking(), G4WeightCutOffProcess::StartTracking(), G4WeightWindowProcess::StartTracking(), G4FastSimulationManagerProcess::StartTracking(), G4ParallelWorldProcess::StartTracking(), G4ParallelWorldScoringProcess::StartTracking(), G4CoupledTransportation::StartTracking(), G4ITSteppingVerbose::StepInfo(), G4SteppingVerbose::StepInfo(), G4SteppingVerboseWithUnits::StepInfo(), G4ITSteppingVerbose::StepInfoForLeadingTrack(), G4SteppingVerbose::TrackingStarted(), G4SteppingVerboseWithUnits::TrackingStarted(), G4ITSteppingVerbose::TrackingStarted(), G4ErrorFreeTrajParam::Update(), G4ErrorFreeTrajState::Update(), G4FieldTrackUpdator::Update(), G4DNAMakeReaction::UpdatePositionForReaction(), G4ITSteppingVerbose::VerboseTrack(), G4SteppingVerbose::VerboseTrack(), G4SteppingVerboseWithUnits::VerboseTrack(), and G4Molecule::~G4Molecule().

◆ GetProperTime()

G4double G4Track::GetProperTime ( ) const

◆ GetStep()

const G4Step * G4Track::GetStep ( ) const

Referenced by G4AdjointForcedInteractionForGamma::AlongStepDoIt(), G4VMscModel::ComputeGeomLimit(), G4DNABrownianTransportation::ComputeGeomLimit(), G4ITModelProcessor::ComputeTrackReaction(), G4UrbanAdjointMscModel::ComputeTruePathLengthLimit(), G4LowEWentzelVIModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4UrbanMscModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), export_G4Track(), G4DNASmoluchowskiReactionModel::FindReaction(), G4NeutrinoElectronProcess::GetMeanFreePath(), G4ElNeutrinoNucleusProcess::GetMeanFreePath(), G4MuNeutrinoNucleusProcess::GetMeanFreePath(), G4Channeling::GetPost(), G4Channeling::GetPre(), G4ErrorPropagator::MakeOneStep(), G4ParallelWorldProcess::PostStepDoIt(), G4Channeling::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength(), G4AdjointForcedInteractionForGamma::PostStepGetPhysicalInteractionLength(), G4ITSteppingVerbose::PostStepVerbose(), G4PSDoseDeposit::ProcessHits(), G4ErrorFreeTrajState::PropagateError(), G4ErrorFreeTrajState::PropagateErrorIoni(), G4ErrorFreeTrajState::PropagateErrorMSC(), G4BOptrForceCollision::ProposeNonPhysicsBiasingOperation(), G4BOptrForceCollision::ProposeOccurenceBiasingOperation(), G4ITStepProcessor::SetTrack(), G4ParallelWorldProcess::StartTracking(), G4ScoreSplittingProcess::StartTracking(), G4ITSteppingVerbose::TrackingStarted(), and G4Channeling::UpdateParameters().

◆ GetStepLength()

G4double G4Track::GetStepLength ( ) const

◆ GetTotalEnergy()

G4double G4Track::GetTotalEnergy ( ) const

◆ GetTouchable()

const G4VTouchable * G4Track::GetTouchable ( ) const

◆ GetTouchableHandle()

const G4TouchableHandle & G4Track::GetTouchableHandle ( ) const

◆ GetTrackID()

G4int G4Track::GetTrackID ( ) const

Referenced by G4ITTrackHolder::_PushTrack(), G4DNABrownianTransportation::AlongStepDoIt(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4DNAIRTMoleculeEncounterStepper::CalculateStep(), G4DNAMoleculeEncounterStepper::CalculateStep(), G4DNAIndependentReactionTimeStepper::CalculateStep(), G4ITReactionSet::CanAddThisReaction(), G4DNAIRTMoleculeEncounterStepper::CheckAndRecordResults(), G4DNAMoleculeEncounterStepper::CheckAndRecordResults(), G4DNAIndependentReactionTimeStepper::CheckAndRecordResults(), G4StackChecker::ClassifyNewTrack(), G4ITStepProcessor::ComputeInteractionLength(), G4DNABrownianTransportation::ComputeStep(), G4ITModelProcessor::ComputeTrackReaction(), G4PhysChemIO::FormattedText::CreateSolvatedElectron(), G4PhysChemIO::G4Analysis::CreateSolvatedElectron(), G4DNAChemistryManager::CreateSolvatedElectron(), G4DNAChemistryManager::CreateWaterMolecule(), G4PhysChemIO::FormattedText::CreateWaterMolecule(), G4PhysChemIO::G4Analysis::CreateWaterMolecule(), G4ITStepProcessor::DealWithSecondaries(), G4DNAMolecularDissociation::DecayIt(), G4DNABrownianTransportation::Diffusion(), G4ITStepProcessor::DoIt(), G4EventManager::DoProcessing(), G4ITStepProcessor::DoStepping(), G4HadronicProcess::DumpState(), G4MuonicAtomDecay::DumpState(), G4ExceptionHandler::DumpTrackInfo(), export_G4Track(), G4HadronicProcess::FillResult(), G4DNAIndependentReactionTimeStepper::FindReaction(), G4SmoothTrajectory::G4SmoothTrajectory(), G4Trajectory::G4Trajectory(), G4VAdjointReverseReaction::GetMeanFreePath(), G4Scintillation::GetScintillationYieldByParticleType(), G4DNAIndependentReactionTimeStepper::GetTimeToEncounter(), G4ITStepProcessor::InvokeAlongStepDoItProcs(), G4SteppingManager::InvokeAlongStepDoItProcs(), G4SteppingManager::InvokeAtRestDoItProcs(), G4SteppingManager::InvokePSDIP(), G4PSPassageCellCurrent::IsPassed(), G4PSPassageCellFlux::IsPassed(), G4PSPassageTrackLength::IsPassed(), G4ITTrackHolder::KillTracks(), compTrackPerID::operator()(), G4StackManager::PopNextTrack(), G4DNAScavengerProcess::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4DNABrownianTransportation::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4AdjointForcedInteractionForGamma::PostStepGetPhysicalInteractionLength(), G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength(), G4ITSteppingVerbose::PostStepVerbose(), G4ITSteppingVerbose::PreStepVerbose(), G4PSPopulation::ProcessHits(), G4StackManager::PushOneTrack(), G4DNAIRT::Sampling(), G4DNAIRT_geometries::Sampling(), G4ITStepProcessor::SetTrack(), G4ITSteppingVerbose::StepInfo(), G4ITSteppingVerbose::StepInfoForLeadingTrack(), G4ITStepProcessor::Stepping(), G4TrackingManager::TrackBanner(), G4ITTrackingInteractivity::TrackBanner(), G4VITSteppingVerbose::TrackBanner(), G4ITSteppingVerbose::TrackingEnded(), G4ITSteppingVerbose::TrackingStarted(), G4ITSteppingVerbose::VerboseTrack(), G4SteppingVerbose::VerboseTrack(), and G4SteppingVerboseWithUnits::VerboseTrack().

◆ GetTrackLength()

G4double G4Track::GetTrackLength ( ) const

◆ GetTrackStatus()

G4TrackStatus G4Track::GetTrackStatus ( ) const

Referenced by G4BiasingProcessInterface::AlongStepDoIt(), G4DNAIRTMoleculeEncounterStepper::CheckAndRecordResults(), G4DNAMoleculeEncounterStepper::CheckAndRecordResults(), G4DNAIndependentReactionTimeStepper::CheckAndRecordResults(), G4ErrorPropagator::CheckIfLastStep(), G4RadioactiveDecay::DecayAnalog(), G4Decay::DecayIt(), G4MuonicAtomDecay::DecayIt(), G4StackManager::DefaultClassification(), G4ITStepProcessor::DoDefinePhysicalStepLength(), G4ITStepProcessor::DoIt(), G4EventManager::DoProcessing(), G4ITStepProcessor::DoStepping(), G4BOptrForceCollision::EndTracking(), export_G4Track(), G4ITStepProcessor::ExtractDoItData(), G4ITStepProcessor::ExtractILData(), G4DNAMolecularReaction::FindReaction(), G4DNAIndependentReactionTimeStepper::FindReaction(), G4Scintillation::GetScintillationYieldByParticleType(), G4ParticleChangeForNothing::Initialize(), G4ParticleChangeForLoss::InitializeForAlongStep(), G4ParticleChangeForGamma::InitializeForPostStep(), G4ParticleChangeForLoss::InitializeForPostStep(), G4ITStepProcessor::InvokeAlongStepDoItProcs(), G4SteppingManager::InvokeAlongStepDoItProcs(), G4ITStepProcessor::InvokePostStepDoItProcs(), G4SteppingManager::InvokePostStepDoItProcs(), G4ITStepProcessor::InvokeTransportationProc(), G4ErrorPropagator::MakeSteps(), G4ImportanceProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4Decay::PostStepDoIt(), G4DecayWithSpin::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4VTransitionRadiation::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4Transportation::PostStepDoIt(), G4PSTermination::ProcessHits(), G4TrackingManager::ProcessOneTrack(), G4ITTrackHolder::PushToKill(), G4ITStepProcessor::SetInitialStep(), G4SteppingManager::SetInitialStep(), G4SteppingManager::Stepping(), G4ParticleChangeForMSC::UpdateStepForAlongStep(), G4ITSteppingVerbose::VerboseTrack(), G4SteppingVerbose::VerboseTrack(), and G4SteppingVerboseWithUnits::VerboseTrack().

◆ GetUserInformation()

G4VUserTrackInformation * G4Track::GetUserInformation ( ) const

Referenced by GetIT().

◆ GetVelocity()

G4double G4Track::GetVelocity ( ) const

◆ GetVertexKineticEnergy()

G4double G4Track::GetVertexKineticEnergy ( ) const

◆ GetVertexMomentumDirection()

const G4ThreeVector & G4Track::GetVertexMomentumDirection ( ) const

◆ GetVertexPosition()

const G4ThreeVector & G4Track::GetVertexPosition ( ) const

◆ GetVolume()

G4VPhysicalVolume * G4Track::GetVolume ( ) const

Referenced by G4ImportanceProcess::AlongStepGetPhysicalInteractionLength(), G4WeightCutOffProcess::AlongStepGetPhysicalInteractionLength(), G4WeightWindowProcess::AlongStepGetPhysicalInteractionLength(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldScoringProcess::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4ParallelGeometriesLimiterProcess::AlongStepGetPhysicalInteractionLength(), G4FastSimulationManagerProcess::AlongStepGetPhysicalInteractionLength(), G4ITTransportation::AlongStepGetPhysicalInteractionLength(), G4DecayWithSpin::AtRestDoIt(), G4FastSimulationManagerProcess::AtRestGetPhysicalInteractionLength(), G4DNABrownianTransportation::ComputeGeomLimit(), G4PolarizedAnnihilation::ComputeSaturationFactor(), G4PolarizedCompton::ComputeSaturationFactor(), G4PolarizedIonisation::ComputeSaturationFactor(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4HadronicProcess::DumpState(), G4MuonicAtomDecay::DumpState(), export_G4Track(), G4Channeling::GetMatData(), G4Channeling::GetMeanFreePath(), G4VXTRenergyLoss::GetMeanFreePath(), G4SynchrotronRadiation::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetMeanFreePath(), G4VTransitionRadiation::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetPhotonEnergy(), G4BOptnForceCommonTruncatedExp::Initialize(), G4VEnergyLossProcess::IsRegionForCubcutProcessor(), G4ScoreSplittingProcess::PostStepDoIt(), G4Channeling::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4VTransitionRadiation::PostStepDoIt(), G4MaxTimeCuts::PostStepGetPhysicalInteractionLength(), G4MinEkineCuts::PostStepGetPhysicalInteractionLength(), G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength(), G4FastSimulationManagerProcess::PostStepGetPhysicalInteractionLength(), G4StepLimiter::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4LowECapture::PostStepGetPhysicalInteractionLength(), G4ErrorFreeTrajState::PropagateErrorIoni(), G4ErrorFreeTrajState::PropagateErrorMSC(), G4TransportationLogger::ReportLoopingTrack(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4PolarizedIonisationModel::SampleSecondaries(), G4ITStepProcessor::SetInitialStep(), G4SteppingManager::SetInitialStep(), G4VPhononProcess::StartTracking(), G4SteppingVerboseWithUnits::StepInfo(), G4SteppingVerboseWithUnits::TrackingStarted(), and G4Channeling::UpdateParameters().

◆ GetWeight()

G4double G4Track::GetWeight ( ) const

Referenced by G4EmBiasingManager::ApplyDirectionalSplitting(), G4BOptnForceFreeFlight::ApplyFinalStateBiasing(), G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4EmBiasingManager::ApplySplitting(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4RadioactiveDecay::DecayAnalog(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4MuonicAtomDecay::DecayIt(), export_G4Track(), G4MuonicAtomDecay::FillResult(), G4FastStep::Initialize(), G4ParticleChangeForLoss::InitializeForAlongStep(), G4ParticleChangeForGamma::InitializeForPostStep(), G4ParticleChangeForLoss::InitializeForPostStep(), G4ImportanceProcess::PostStepDoIt(), G4WeightCutOffProcess::PostStepDoIt(), G4WeightWindowProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4AdjointTrackingAction::PreUserTrackingAction(), G4BOptrForceCollision::ProposeNonPhysicsBiasingOperation(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointhIonisationModel::RapidSampleSecondaries(), G4AdjointBremsstrahlungModel::SampleSecondaries(), G4AdjointComptonModel::SampleSecondaries(), G4AdjointeIonisationModel::SampleSecondaries(), G4AdjointhIonisationModel::SampleSecondaries(), G4AdjointIonIonisationModel::SampleSecondaries(), G4AdjointPhotoElectricModel::SampleSecondaries(), G4ParticleChangeForOccurenceBiasing::StealSecondaries(), and G4AdjointSteppingAction::UserSteppingAction().

◆ IncrementCurrentStepNumber()

void G4Track::IncrementCurrentStepNumber ( )

◆ IsBelowThreshold()

G4bool G4Track::IsBelowThreshold ( ) const

◆ IsGoodForTracking()

G4bool G4Track::IsGoodForTracking ( ) const

◆ operator delete()

void G4Track::operator delete ( void *  aTrack)
inline

◆ operator new()

void * G4Track::operator new ( std::size_t  )
inline

◆ operator!=()

G4bool G4Track::operator!= ( const G4Track )
inline

◆ operator=()

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

Definition at line 83 of file G4Track.cc.

84{
85 if(this != &right)
86 {
87 fPosition = right.fPosition;
89 fLocalTime = right.fLocalTime;
91 fWeight = right.fWeight;
93
94 // additional fields required for geometrical splitting
98
99 // Track ID (and Parent ID) is not copied and set to zero for new track
100 fTrackID = 0;
101 fParentID = 0;
102
103 // CurrentStepNumber is set to be 0
105
106 // Creator model ID
108
109 // velocity information
110 fVelocity = right.fVelocity;
111
112 // dynamic particle information
113 delete fpDynamicParticle;
115
116 // track status and flags for tracking
120
121 // Step information (Step Length, Step Number, pointer to the Step,)
122 // are not copied
123 fpStep = nullptr;
124
125 // vertex information
130
131 // CreatorProcess and UserInformation are not copied
132 fpCreatorProcess = nullptr;
133 delete fpUserInformation;
134 fpUserInformation = nullptr;
135
136 prev_mat = right.prev_mat;
137 groupvel = right.groupvel;
140
143
145 }
146 return *this;
147}
G4double fLocalTime
Definition: G4Track.hh:266
G4TrackStatus fTrackStatus
Definition: G4Track.hh:279
G4double fWeight
Definition: G4Track.hh:287
G4double fVtxKineticEnergy
Definition: G4Track.hh:296
G4int fParentID
Definition: G4Track.hh:320
G4double fStepLength
Definition: G4Track.hh:281
G4int fCreatorModelID
Definition: G4Track.hh:317
G4TouchableHandle fpOriginTouchable
Definition: G4Track.hh:275
const G4VProcess * fpCreatorProcess
Definition: G4Track.hh:300
G4ThreeVector fVtxMomentumDirection
Definition: G4Track.hh:294
G4bool useGivenVelocity
Definition: G4Track.hh:332
G4int fTrackID
Definition: G4Track.hh:321
const G4LogicalVolume * fpLVAtVertex
Definition: G4Track.hh:298
G4double fTrackLength
Definition: G4Track.hh:268
G4int fCurrentStepNumber
Definition: G4Track.hh:314
G4bool fBelowThreshold
Definition: G4Track.hh:323
G4TouchableHandle fpNextTouchable
Definition: G4Track.hh:274
G4ThreeVector fVtxPosition
Definition: G4Track.hh:292
G4bool fGoodForTracking
Definition: G4Track.hh:326

References ClearAuxiliaryTrackInformation(), fBelowThreshold, fCreatorModelID, fCurrentStepNumber, fGlobalTime, fGoodForTracking, fLocalTime, fParentID, fpCreatorProcess, fpDynamicParticle, fpLVAtVertex, fpNextTouchable, fpOriginTouchable, fPosition, fpStep, fpTouchable, fpUserInformation, fStepLength, fTrackID, fTrackLength, fTrackStatus, fVelocity, fVtxKineticEnergy, fVtxMomentumDirection, fVtxPosition, fWeight, groupvel, is_OpticalPhoton, prev_mat, prev_momentum, prev_velocity, and useGivenVelocity.

◆ operator==()

G4bool G4Track::operator== ( const G4Track )
inline

◆ RemoveAuxiliaryTrackInformation() [1/2]

void G4Track::RemoveAuxiliaryTrackInformation ( G4int  id)

Definition at line 237 of file G4Track.cc.

238{
239 if(fpAuxiliaryTrackInformationMap != nullptr &&
241 {
243 }
244}

References fpAuxiliaryTrackInformationMap.

Referenced by RemoveAuxiliaryTrackInformation().

◆ RemoveAuxiliaryTrackInformation() [2/2]

void G4Track::RemoveAuxiliaryTrackInformation ( G4String name)

Definition at line 247 of file G4Track.cc.

248{
249 if(fpAuxiliaryTrackInformationMap != nullptr)
250 {
253 }
254}
int G4int
Definition: G4Types.hh:85
static G4int GetModelID(const G4int modelIndex)
void RemoveAuxiliaryTrackInformation(G4int id)
Definition: G4Track.cc:237
const char * name(G4int ptype)

References fpAuxiliaryTrackInformationMap, G4PhysicsModelCatalog::GetModelID(), G4InuclParticleNames::name(), and RemoveAuxiliaryTrackInformation().

◆ SetAuxiliaryTrackInformation()

void G4Track::SetAuxiliaryTrackInformation ( G4int  id,
G4VAuxiliaryTrackInformation info 
) const

Definition at line 205 of file G4Track.cc.

207{
208 if(fpAuxiliaryTrackInformationMap == nullptr)
209 {
211 new std::map<G4int, G4VAuxiliaryTrackInformation*>;
212 }
214 {
216 ED << "Process/model ID <" << id << "> is invalid.";
217 G4Exception("G4VAuxiliaryTrackInformation::G4VAuxiliaryTrackInformation()",
218 "TRACK0982", FatalException, ED);
219 }
220 (*fpAuxiliaryTrackInformationMap)[id] = info;
221}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static G4int GetModelIndex(const G4int modelID)

References FatalException, fpAuxiliaryTrackInformationMap, G4Exception(), and G4PhysicsModelCatalog::GetModelIndex().

Referenced by G4eplusAnnihilation::AtRestDoIt(), G4Channeling::GetTrackData(), G4BOptrForceCollision::OperationApplied(), and G4BOptrForceCollision::ProposeNonPhysicsBiasingOperation().

◆ SetBelowThresholdFlag()

void G4Track::SetBelowThresholdFlag ( G4bool  value = true)

◆ SetCreatorModelID()

void G4Track::SetCreatorModelID ( const G4int  id)
inline

◆ SetCreatorProcess()

void G4Track::SetCreatorProcess ( const G4VProcess aValue)

◆ SetGlobalTime()

void G4Track::SetGlobalTime ( const G4double  aValue)

◆ SetGoodForTrackingFlag()

void G4Track::SetGoodForTrackingFlag ( G4bool  value = true)

◆ SetKineticEnergy()

void G4Track::SetKineticEnergy ( const G4double  aValue)

◆ SetLocalTime()

void G4Track::SetLocalTime ( const G4double  aValue)

◆ SetLogicalVolumeAtVertex()

void G4Track::SetLogicalVolumeAtVertex ( const G4LogicalVolume )

◆ SetMomentumDirection()

void G4Track::SetMomentumDirection ( const G4ThreeVector aValue)

◆ SetNextTouchableHandle()

void G4Track::SetNextTouchableHandle ( const G4TouchableHandle apValue)

◆ SetOriginTouchableHandle()

void G4Track::SetOriginTouchableHandle ( const G4TouchableHandle apValue)

◆ SetParentID()

void G4Track::SetParentID ( const G4int  aValue)

◆ SetPolarization()

void G4Track::SetPolarization ( const G4ThreeVector aValue)

◆ SetPosition()

void G4Track::SetPosition ( const G4ThreeVector aValue)

◆ SetProperTime()

void G4Track::SetProperTime ( const G4double  aValue)

◆ SetStep()

void G4Track::SetStep ( const G4Step aValue)

◆ SetStepLength()

void G4Track::SetStepLength ( G4double  value)

◆ SetTouchableHandle()

void G4Track::SetTouchableHandle ( const G4TouchableHandle apValue)

◆ SetTrackID()

void G4Track::SetTrackID ( const G4int  aValue)

◆ SetTrackStatus()

void G4Track::SetTrackStatus ( const G4TrackStatus  aTrackStatus)

◆ SetUserInformation()

void G4Track::SetUserInformation ( G4VUserTrackInformation aValue) const

◆ SetVelocity()

void G4Track::SetVelocity ( G4double  val)

◆ SetVertexKineticEnergy()

void G4Track::SetVertexKineticEnergy ( const G4double  aValue)

◆ SetVertexMomentumDirection()

void G4Track::SetVertexMomentumDirection ( const G4ThreeVector aValue)

◆ SetVertexPosition()

void G4Track::SetVertexPosition ( const G4ThreeVector aValue)

◆ SetWeight()

void G4Track::SetWeight ( G4double  aValue)

◆ UseGivenVelocity() [1/2]

G4bool G4Track::UseGivenVelocity ( ) const

◆ UseGivenVelocity() [2/2]

void G4Track::UseGivenVelocity ( G4bool  val)

Field Documentation

◆ fBelowThreshold

G4bool G4Track::fBelowThreshold = false
private

Definition at line 323 of file G4Track.hh.

Referenced by operator=().

◆ fCreatorModelID

G4int G4Track::fCreatorModelID = -1
private

Definition at line 317 of file G4Track.hh.

Referenced by operator=().

◆ fCurrentStepNumber

G4int G4Track::fCurrentStepNumber = 0
private

Definition at line 314 of file G4Track.hh.

Referenced by operator=().

◆ fGlobalTime

G4double G4Track::fGlobalTime = 0.0
private

Definition at line 264 of file G4Track.hh.

Referenced by operator=().

◆ fGoodForTracking

G4bool G4Track::fGoodForTracking = false
private

Definition at line 326 of file G4Track.hh.

Referenced by operator=().

◆ fLocalTime

G4double G4Track::fLocalTime = 0.0
private

Definition at line 266 of file G4Track.hh.

Referenced by operator=().

◆ fParentID

G4int G4Track::fParentID = 0
private

Definition at line 320 of file G4Track.hh.

Referenced by operator=().

◆ fpAuxiliaryTrackInformationMap

std::map<G4int, G4VAuxiliaryTrackInformation*>* G4Track::fpAuxiliaryTrackInformationMap = nullptr
mutableprivate

◆ fpCreatorProcess

const G4VProcess* G4Track::fpCreatorProcess = nullptr
private

Definition at line 300 of file G4Track.hh.

Referenced by operator=().

◆ fpDynamicParticle

G4DynamicParticle* G4Track::fpDynamicParticle = nullptr
private

Definition at line 278 of file G4Track.hh.

Referenced by CalculateVelocityForOpticalPhoton(), G4Track(), operator=(), and ~G4Track().

◆ fpLVAtVertex

const G4LogicalVolume* G4Track::fpLVAtVertex = nullptr
private

Definition at line 298 of file G4Track.hh.

Referenced by operator=().

◆ fpNextTouchable

G4TouchableHandle G4Track::fpNextTouchable
private

Definition at line 274 of file G4Track.hh.

Referenced by operator=().

◆ fpOriginTouchable

G4TouchableHandle G4Track::fpOriginTouchable
private

Definition at line 275 of file G4Track.hh.

Referenced by operator=().

◆ fPosition

G4ThreeVector G4Track::fPosition
private

Definition at line 262 of file G4Track.hh.

Referenced by operator=().

◆ fpStep

const G4Step* G4Track::fpStep = nullptr
private

Definition at line 290 of file G4Track.hh.

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

◆ fpTouchable

G4TouchableHandle G4Track::fpTouchable
private

Definition at line 273 of file G4Track.hh.

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

◆ fpUserInformation

G4VUserTrackInformation* G4Track::fpUserInformation = nullptr
mutableprivate

Definition at line 303 of file G4Track.hh.

Referenced by operator=(), and ~G4Track().

◆ fStepLength

G4double G4Track::fStepLength = 0.0
private

Definition at line 281 of file G4Track.hh.

Referenced by operator=().

◆ fTrackID

G4int G4Track::fTrackID = 0
private

Definition at line 321 of file G4Track.hh.

Referenced by operator=().

◆ fTrackLength

G4double G4Track::fTrackLength = 0.0
private

Definition at line 268 of file G4Track.hh.

Referenced by operator=().

◆ fTrackStatus

G4TrackStatus G4Track::fTrackStatus = fAlive
private

Definition at line 279 of file G4Track.hh.

Referenced by operator=().

◆ fVelocity

G4double G4Track::fVelocity = 0.0
private

Definition at line 271 of file G4Track.hh.

Referenced by operator=().

◆ fVtxKineticEnergy

G4double G4Track::fVtxKineticEnergy = 0.0
private

Definition at line 296 of file G4Track.hh.

Referenced by operator=().

◆ fVtxMomentumDirection

G4ThreeVector G4Track::fVtxMomentumDirection
private

Definition at line 294 of file G4Track.hh.

Referenced by operator=().

◆ fVtxPosition

G4ThreeVector G4Track::fVtxPosition
private

Definition at line 292 of file G4Track.hh.

Referenced by operator=().

◆ fWeight

G4double G4Track::fWeight = 1.0
private

Definition at line 287 of file G4Track.hh.

Referenced by operator=().

◆ groupvel

G4MaterialPropertyVector* G4Track::groupvel = nullptr
mutableprivate

Definition at line 306 of file G4Track.hh.

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

◆ is_OpticalPhoton

G4bool G4Track::is_OpticalPhoton = false
private

Definition at line 330 of file G4Track.hh.

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

◆ prev_mat

G4Material* G4Track::prev_mat = nullptr
mutableprivate

Definition at line 305 of file G4Track.hh.

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

◆ prev_momentum

G4double G4Track::prev_momentum = 0.0
mutableprivate

Definition at line 308 of file G4Track.hh.

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

◆ prev_velocity

G4double G4Track::prev_velocity = 0.0
mutableprivate

Definition at line 307 of file G4Track.hh.

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

◆ useGivenVelocity

G4bool G4Track::useGivenVelocity = false
private

Definition at line 332 of file G4Track.hh.

Referenced by operator=().


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