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

#include <G4Track.hh>

Public Member Functions

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

Detailed Description

Definition at line 65 of file G4Track.hh.

Constructor & Destructor Documentation

◆ G4Track() [1/3]

G4Track::G4Track ( )

Definition at line 62 of file G4Track.cc.

63 : fVelocity(c_light)
64 , fpDynamicParticle(new G4DynamicParticle())
65{}

Referenced by CopyTrackInfo(), G4Track(), operator delete(), operator!=(), operator=(), and operator==().

◆ 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)
52 , fVelocity(c_light)
53{
54 fpDynamicParticle = (apValueDynamicParticle) != nullptr
55 ? apValueDynamicParticle : new G4DynamicParticle();
56 // check if the particle type is Optical Photon
57 is_OpticalPhoton =
58 (fpDynamicParticle->GetDefinition()->GetPDGEncoding() == -22);
59}

◆ G4Track() [3/3]

G4Track::G4Track ( const G4Track & right,
G4bool copyTouchables = true )

Definition at line 68 of file G4Track.cc.

69 : fVelocity(c_light)
70{
71 fCopyTouchables = copyTouchables;
72 *this = right;
73 fCopyTouchables = true;
74}

◆ ~G4Track()

G4Track::~G4Track ( )

Definition at line 77 of file G4Track.cc.

78{
79 delete fpDynamicParticle;
80 delete fpUserInformation;
81 ClearAuxiliaryTrackInformation();
82}

Member Function Documentation

◆ AddTrackLength()

void G4Track::AddTrackLength ( const G4double aValue)
inline

◆ CalculateVelocity()

◆ CalculateVelocityForOpticalPhoton()

G4double G4Track::CalculateVelocityForOpticalPhoton ( ) const

Definition at line 166 of file G4Track.cc.

167{
168 G4double velocity = c_light;
169
170 G4Material* mat = nullptr;
171 G4bool update_groupvel = false;
172 if(fpStep != nullptr)
173 {
174 mat = this->GetMaterial(); // Fix for repeated volumes
175 }
176 else
177 {
178 if(fpTouchable)
179 {
180 mat = fpTouchable->GetVolume()->GetLogicalVolume()->GetMaterial();
181 }
182 }
183 // check if previous step is in the same volume
184 // and get new GROUPVELOCITY table if necessary
185 if((mat != nullptr) && ((mat != prev_mat) || (groupvel == nullptr)))
186 {
187 groupvel = nullptr;
188 if(mat->GetMaterialPropertiesTable() != nullptr)
190 update_groupvel = true;
191 }
192 prev_mat = mat;
193
194 if(groupvel != nullptr)
195 {
196 // light velocity = c/(rindex+d(rindex)/d(log(E_phot)))
197 // values stored in GROUPVEL material properties vector
198 velocity = prev_velocity;
199
200 // check if momentum is same as in the previous step
201 // and calculate group velocity if necessary
202 G4double current_momentum = fpDynamicParticle->GetTotalMomentum();
203 if(update_groupvel || (current_momentum != prev_momentum))
204 {
205 velocity = groupvel->Value(current_momentum);
206 prev_velocity = velocity;
207 prev_momentum = current_momentum;
208 }
209 }
210
211 return velocity;
212}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
G4Material * GetMaterial() const

Referenced by G4ITTransportation::AlongStepDoIt(), and G4ParticleChange::UpdateStepForAlongStep().

◆ CopyTrackInfo()

void G4Track::CopyTrackInfo ( const G4Track & right,
G4bool copyTouchables = true )

Definition at line 158 of file G4Track.cc.

159{
160 fCopyTouchables = copyTouchables;
161 *this = right;
162 fCopyTouchables = true;
163}

Referenced by G4WorkerSubEvtRunManager::DoWork().

◆ GetAuxiliaryTrackInformation()

G4VAuxiliaryTrackInformation * G4Track::GetAuxiliaryTrackInformation ( G4int id) const

Definition at line 235 of file G4Track.cc.

236{
237 if(fpAuxiliaryTrackInformationMap == nullptr)
238 return nullptr;
239
240 auto itr = fpAuxiliaryTrackInformationMap->find(id);
241 if(itr == fpAuxiliaryTrackInformationMap->cend())
242 return nullptr;
243 return (*itr).second;
244}

◆ 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
inline

◆ GetCurrentStepNumber()

◆ GetDefinition()

◆ GetDynamicParticle()

const G4DynamicParticle * G4Track::GetDynamicParticle ( ) const
inline

Referenced by G4AdjointAlongStepWeightCorrection::AlongStepDoIt(), G4AdjointProcessEquivalentToDirectProcess::AlongStepDoIt(), G4ContinuousGainOfEnergy::AlongStepDoIt(), G4DynamicParticleIonisation::AlongStepDoIt(), G4ErrorEnergyLoss::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4Transportation::AlongStepDoIt(), G4VEnergyLossProcess::AlongStepDoIt(), G4AdjointProcessEquivalentToDirectProcess::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4ITTransportation::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4TransportationWithMsc::AlongStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4ITStepProcessor::ApplyProductionCut(), G4AdjointProcessEquivalentToDirectProcess::AtRestDoIt(), G4DecayWithSpin::AtRestDoIt(), G4AdjointProcessEquivalentToDirectProcess::AtRestGetPhysicalInteractionLength(), G4Decay::AtRestGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4BetheBlochIonGasModel::ChargeSquareRatio(), G4BraggIonGasModel::ChargeSquareRatio(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4ITTransportation::ComputeStep(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4LowEWentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanAdjointMscModel::ComputeTruePathLengthLimit(), G4UrbanMscModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4FieldTrackUpdator::CreateFieldTrack(), G4PionDecayMakeSpin::DaughterPolarization(), G4VRadioactiveDecay::DecayAnalog(), G4Decay::DecayIt(), G4RadioactiveDecay::DecayIt(), G4UnknownDecay::DecayIt(), G4VRadioactiveDecay::DecayIt(), G4ErrorEnergyLoss::GetContinuousStepLimit(), G4hImpactIonisation::GetContinuousStepLimit(), G4AnnihiToMuPair::GetMeanFreePath(), G4Decay::GetMeanFreePath(), G4ElNeutrinoNucleusProcess::GetMeanFreePath(), G4GammaConversionToMuons::GetMeanFreePath(), G4HadronicProcess::GetMeanFreePath(), G4hImpactIonisation::GetMeanFreePath(), G4MuNeutrinoNucleusProcess::GetMeanFreePath(), G4MuonicAtomDecay::GetMeanFreePath(), G4NeutrinoElectronProcess::GetMeanFreePath(), G4OpAbsorption::GetMeanFreePath(), G4OpMieHG::GetMeanFreePath(), G4OpRayleigh::GetMeanFreePath(), G4OpWLS2::GetMeanFreePath(), G4OpWLS::GetMeanFreePath(), G4SynchrotronRadiation::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetMeanFreePath(), G4TauNeutrinoNucleusProcess::GetMeanFreePath(), G4VRadioactiveDecay::GetMeanFreePath(), G4VXTRenergyLoss::GetMeanFreePath(), G4XrayReflection::GetMeanFreePath(), G4Decay::GetMeanLifeTime(), G4MuonicAtomDecay::GetMeanLifeTime(), G4SynchrotronRadiationInMat::GetPhotonEnergy(), G4Scintillation::GetScintillationYieldByParticleType(), G4HadProjectile::Initialise(), G4FastStep::Initialize(), G4ParticleChange::Initialize(), G4ParticleChangeForDecay::Initialize(), G4ParticleChangeForLoss::InitializeForAlongStep(), G4VEmProcess::MeanFreePath(), G4VEnergyLossProcess::MeanFreePath(), G4SmartTrackStack::PopFromStack(), G4AdjointProcessEquivalentToDirectProcess::PostStepDoIt(), G4AnnihiToMuPair::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4DecayWithSpin::PostStepDoIt(), G4DynamicParticleIonisation::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4NeutronGeneralProcess::PostStepDoIt(), G4NuVacOscProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpMieHG::PostStepDoIt(), G4OpRayleigh::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4XrayReflection::PostStepDoIt(), G4AdjointProcessEquivalentToDirectProcess::PostStepGetPhysicalInteractionLength(), G4Cerenkov::PostStepGetPhysicalInteractionLength(), G4Decay::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4GammaGeneralProcess::PostStepGetPhysicalInteractionLength(), G4MaxTimeCuts::PostStepGetPhysicalInteractionLength(), G4MinEkineCuts::PostStepGetPhysicalInteractionLength(), G4UnknownDecay::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VITDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4ErrorFreeTrajState::PropagateError(), 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(), G4DNAGeneralIonIonisationModel::StartTracking(), G4DNAIonChargeDecreaseModel::StartTracking(), G4DNAIonChargeIncreaseModel::StartTracking(), G4GoudsmitSaundersonMscModel::StartTracking(), G4HadronicProcess::StartTracking(), G4UrbanAdjointMscModel::StartTracking(), G4UrbanMscModel::StartTracking(), and G4FieldTrackUpdator::Update().

◆ GetGlobalTime()

G4double G4Track::GetGlobalTime ( ) const
inline

Referenced by G4ITTrackHolder::_PushTrack(), G4ITTransportation::AlongStepDoIt(), G4Transportation::AlongStepDoIt(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4ITTransportation::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4DecayWithSpin::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4FastStep::CheckIt(), G4VParticleChange::CheckSecondary(), G4StackChecker::ClassifyNewTrack(), G4ITStepProcessor::ComputeInteractionLength(), G4ITTransportation::ComputeStep(), G4ITModelProcessor::ComputeTrackReaction(), G4FieldTrackUpdator::CreateFieldTrack(), G4VRadioactiveDecay::DecayAnalog(), G4Decay::DecayIt(), G4DNAMolecularDissociation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4UnknownDecay::DecayIt(), G4DNABrownianTransportation::Diffusion(), G4ChannelingFastSimModel::DoIt(), G4HadronicProcess::FillResult(), G4SynchrotronRadiation::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetMeanFreePath(), G4FastStep::Initialize(), G4ParticleChange::Initialize(), G4ParticleChangeForDecay::Initialize(), G4FastSimHitMaker::make(), GFlashHitMaker::make(), G4DNAIRT::MakeReaction(), G4DNAMakeReaction::MakeReaction(), G4DNAMolecularReaction::MakeReaction(), G4CoherentPairProduction::PostStepDoIt(), G4DNABrownianTransportation::PostStepDoIt(), G4DNAScavengerProcess::PostStepDoIt(), G4DNASecondOrderReaction::PostStepDoIt(), G4DynamicParticleIonisation::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4UCNAbsorption::PostStepDoIt(), G4UCNMultiScattering::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4DNAPolyNucleotideReactionProcess::PostStepGetPhysicalInteractionLength(), G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4MaxTimeCuts::PostStepGetPhysicalInteractionLength(), G4NeutronGeneralProcess::PostStepGetPhysicalInteractionLength(), G4NeutronKiller::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4ITSteppingVerbose::PreStepVerbose(), G4ITTrackHolder::PushDelayed(), G4TrackingInformation::RecordCurrentPositionNTime(), G4DNAIRT::Sampling(), G4FieldTrackUpdator::Update(), G4FastStep::UpdateStepForAtRest(), and G4FastStep::UpdateStepForPostStep().

◆ GetKineticEnergy()

G4double G4Track::GetKineticEnergy ( ) const
inline

Referenced by G4AdjointForcedInteractionForGamma::AlongStepDoIt(), G4ErrorEnergyLoss::AlongStepDoIt(), G4ITTransportation::AlongStepDoIt(), G4VMultipleScattering::AlongStepDoIt(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4ITTransportation::AlongStepGetPhysicalInteractionLength(), G4TransportationWithMsc::AlongStepGetPhysicalInteractionLength(), G4VMultipleScattering::AlongStepGetPhysicalInteractionLength(), G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4ITStepProcessor::ApplyProductionCut(), G4VEmModel::ChargeSquareRatio(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4FastStep::CheckIt(), G4VParticleChange::CheckSecondary(), G4StackChecker::ClassifyNewTrack(), G4FieldTrackUpdator::CreateFieldTrack(), G4PhysChemIO::FormattedText::CreateSolvatedElectron(), G4PhysChemIO::G4Analysis::CreateSolvatedElectron(), G4ITStepProcessor::DealWithSecondaries(), G4ChannelingFastSimModel::DoIt(), G4HadronicProcess::DumpState(), G4HadronicProcess::FillResult(), G4RichTrajectory::G4RichTrajectory(), G4SmoothTrajectory::G4SmoothTrajectory(), G4Trajectory::G4Trajectory(), G4AdjointAlongStepWeightCorrection::GetContinuousStepLimit(), G4ContinuousGainOfEnergy::GetContinuousStepLimit(), G4ErrorEnergyLoss::GetContinuousStepLimit(), G4MicroElecCapture::GetMeanFreePath(), G4NuVacOscProcess::GetMeanFreePath(), G4PhononDownconversion::GetMeanFreePath(), G4PhononScattering::GetMeanFreePath(), G4VAdjointReverseReaction::GetMeanFreePath(), G4VRadioactiveDecay::GetMeanFreePath(), G4VTransitionRadiation::GetMeanFreePath(), G4VRadioactiveDecay::GetMeanLifeTime(), G4Scintillation::GetScintillationYieldByParticleType(), G4ParticleChangeForLoss::InitializeForAlongStep(), G4ParticleChangeForGamma::InitializeForPostStep(), G4VEmProcess::MeanFreePath(), G4VEnergyLossProcess::MeanFreePath(), G4ChannelingFastSimModel::ModelTrigger(), GFlashShowerModel::ModelTrigger(), G4AdjointForcedInteractionForGamma::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4LowECapture::PostStepDoIt(), G4MicroElecCapture::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4NuVacOscProcess::PostStepDoIt(), G4PhononReflection::PostStepDoIt(), G4PhononScattering::PostStepDoIt(), G4SpecialCuts::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4WeightWindowProcess::PostStepDoIt(), G4AdjointForcedInteractionForGamma::PostStepGetPhysicalInteractionLength(), G4GammaGeneralProcess::PostStepGetPhysicalInteractionLength(), G4HadronicProcess::PostStepGetPhysicalInteractionLength(), G4LowECapture::PostStepGetPhysicalInteractionLength(), G4NeutronKiller::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4TransportationLogger::ReportLoopingTrack(), G4FieldTrackUpdator::Update(), G4ParticleChange::UpdateStepForAlongStep(), and G4AdjointSteppingAction::UserSteppingAction().

◆ GetLocalTime()

◆ GetLogicalVolumeAtVertex()

const G4LogicalVolume * G4Track::GetLogicalVolumeAtVertex ( ) const
inline

◆ GetMaterial()

G4Material * G4Track::GetMaterial ( ) const
inline

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

◆ GetMaterialCutsCouple()

const G4MaterialCutsCouple * G4Track::GetMaterialCutsCouple ( ) const
inline

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

◆ GetMomentum()

◆ GetMomentumDirection()

const G4ThreeVector & G4Track::GetMomentumDirection ( ) const
inline

Referenced by G4DNABrownianTransportation::AlongStepDoIt(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4VParticleChange::CheckSecondary(), G4StackChecker::ClassifyNewTrack(), G4DNABrownianTransportation::ComputeGeomLimit(), G4VMscModel::ComputeGeomLimit(), G4DNABrownianTransportation::ComputeStep(), G4FieldTrackUpdator::CreateFieldTrack(), G4HadronicProcess::DumpState(), G4HadronicProcess::FillResult(), G4CoherentPairProduction::GetMeanFreePath(), G4VFastSimSensitiveDetector::Hit(), G4VGFlashSensitiveDetector::Hit(), G4BOptnForceCommonTruncatedExp::Initialize(), G4ParticleChangeForGamma::InitializeForPostStep(), G4ParticleChangeForLoss::InitializeForPostStep(), G4CoupledTransportation::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4PhononReflection::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4Transportation::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4VTransitionRadiation::PostStepDoIt(), G4CoupledTransportation::StartTracking(), G4FastSimulationManagerProcess::StartTracking(), G4ImportanceProcess::StartTracking(), G4ParallelGeometriesLimiterProcess::StartTracking(), G4ParallelWorldProcess::StartTracking(), G4ParallelWorldScoringProcess::StartTracking(), G4VPhononProcess::StartTracking(), G4WeightCutOffProcess::StartTracking(), G4WeightWindowProcess::StartTracking(), and G4FieldTrackUpdator::Update().

◆ GetNextMaterial()

G4Material * G4Track::GetNextMaterial ( ) const
inline

◆ GetNextMaterialCutsCouple()

const G4MaterialCutsCouple * G4Track::GetNextMaterialCutsCouple ( ) const
inline

◆ GetNextTouchable()

const G4VTouchable * G4Track::GetNextTouchable ( ) const
inline

◆ GetNextTouchableHandle()

◆ GetNextVolume()

◆ GetOriginTouchable()

const G4VTouchable * G4Track::GetOriginTouchable ( ) const
inline

◆ GetOriginTouchableHandle()

const G4TouchableHandle & G4Track::GetOriginTouchableHandle ( ) const
inline

◆ GetParentID()

◆ GetParentResonanceDef()

const G4ParticleDefinition * G4Track::GetParentResonanceDef ( ) const
inline

◆ GetParentResonanceID()

G4int G4Track::GetParentResonanceID ( ) const
inline

◆ GetParentResonanceMass()

G4double G4Track::GetParentResonanceMass ( ) const
inline

◆ GetParentResonanceName()

G4String G4Track::GetParentResonanceName ( ) const
inline

◆ GetParentResonancePDGEncoding()

G4int G4Track::GetParentResonancePDGEncoding ( ) const
inline

◆ GetParticleDefinition()

◆ GetPolarization()

◆ GetPosition()

const G4ThreeVector & G4Track::GetPosition ( ) const
inline

Referenced by G4AdjointForcedInteractionForGamma::AlongStepDoIt(), G4DNABrownianTransportation::AlongStepDoIt(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4ITTransportation::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4eplusAnnihilation::AtRestDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4DNAIndependentReactionTimeStepper::CalculateStep(), G4StackChecker::ClassifyNewTrack(), G4DNABrownianTransportation::ComputeStep(), G4ITTransportation::ComputeStep(), G4FieldTrackUpdator::CreateFieldTrack(), G4DNAMultipleIonisationManager::CreateMultipleIonisedWaterMolecule(), G4DNAChemistryManager::CreateSolvatedElectron(), G4PhysChemIO::FormattedText::CreateSolvatedElectron(), G4PhysChemIO::G4Analysis::CreateSolvatedElectron(), G4DNAChemistryManager::CreateWaterMolecule(), G4PhysChemIO::FormattedText::CreateWaterMolecule(), G4PhysChemIO::G4Analysis::CreateWaterMolecule(), G4VRadioactiveDecay::DecayAnalog(), G4Decay::DecayIt(), G4DNAMolecularDissociation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4UnknownDecay::DecayIt(), G4HadronicProcess::DumpState(), G4HadronicProcess::FillResult(), G4DNASmoluchowskiReactionModel::FindReaction(), G4SmoothTrajectory::G4SmoothTrajectory(), G4Trajectory::G4Trajectory(), G4CoherentPairProduction::GetMeanFreePath(), G4SynchrotronRadiation::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetMeanFreePath(), G4XrayReflection::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetPhotonEnergy(), G4GFlashSpot::GetPosition(), G4IT::GetPosition(), G4TDNAOneStepThermalizationModel< DNA::Penetration::Meesungnoen2002 >::GetRmean(), G4DiffusionControlledReactionModel::GetTimeToEncounter(), G4BOptnForceCommonTruncatedExp::Initialize(), G4FastStep::Initialize(), G4ParticleChange::Initialize(), G4FastSimHitMaker::make(), G4DNAIRT::MakeReaction(), G4DNAMolecularReaction::MakeReaction(), G4CoherentPairProduction::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4DNAScavengerProcess::PostStepDoIt(), G4DNASecondOrderReaction::PostStepDoIt(), G4DynamicParticleIonisation::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4Transportation::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength(), G4ErrorMagFieldLimitProcess::PostStepGetPhysicalInteractionLength(), G4ITSteppingVerbose::PostStepVerbose(), G4ITSteppingVerbose::PreStepVerbose(), G4ErrorFreeTrajState::PropagateError(), G4TrackingInformation::RecordCurrentPositionNTime(), G4TransportationLogger::ReportLoopingTrack(), G4DNAIRT::Sampling(), G4CoupledTransportation::StartTracking(), G4FastSimulationManagerProcess::StartTracking(), G4ImportanceProcess::StartTracking(), G4ParallelGeometriesLimiterProcess::StartTracking(), G4ParallelWorldProcess::StartTracking(), G4ParallelWorldScoringProcess::StartTracking(), G4WeightCutOffProcess::StartTracking(), G4WeightWindowProcess::StartTracking(), G4ErrorFreeTrajParam::Update(), G4ErrorFreeTrajState::Update(), G4FieldTrackUpdator::Update(), and G4DNAMakeReaction::UpdatePositionForReaction().

◆ GetProperTime()

◆ GetStep()

◆ GetStepLength()

◆ GetTotalEnergy()

◆ GetTouchable()

◆ GetTouchableHandle()

◆ GetTrackID()

G4int G4Track::GetTrackID ( ) const
inline

Referenced by G4ITTrackHolder::_PushTrack(), G4DNABrownianTransportation::AlongStepDoIt(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4DNAIndependentReactionTimeStepper::CalculateStep(), G4DNAIRTMoleculeEncounterStepper::CalculateStep(), G4DNAMoleculeEncounterStepper::CalculateStep(), G4DNAPolyNucleotideReactionProcess::CalculateTimeStep(), G4ITReactionSet::CanAddThisReaction(), G4FastList< typename ObjectW::type >::CheckFlag(), G4StackChecker::ClassifyNewTrack(), G4ITStepProcessor::ComputeInteractionLength(), G4DNABrownianTransportation::ComputeStep(), G4ITModelProcessor::ComputeTrackReaction(), G4DNAMultipleIonisationManager::CreateMultipleIonisedWaterMolecule(), G4DNAChemistryManager::CreateSolvatedElectron(), G4PhysChemIO::FormattedText::CreateSolvatedElectron(), G4PhysChemIO::G4Analysis::CreateSolvatedElectron(), G4DNAChemistryManager::CreateWaterMolecule(), G4PhysChemIO::FormattedText::CreateWaterMolecule(), G4PhysChemIO::G4Analysis::CreateWaterMolecule(), G4DNAMolecularDissociation::DecayIt(), G4DNABrownianTransportation::Diffusion(), G4ITStepProcessor::DoIt(), G4HadronicProcess::DumpState(), G4HadronicProcess::FillResult(), G4DNAIndependentReactionTimeStepper::FindReaction(), G4RichTrajectory::G4RichTrajectory(), G4SmoothTrajectory::G4SmoothTrajectory(), G4Trajectory::G4Trajectory(), G4VAdjointReverseReaction::GetMeanFreePath(), G4Scintillation::GetScintillationYieldByParticleType(), G4PSPassageCellCurrent::IsPassed(), G4PSPassageCellFlux::IsPassed(), G4PSPassageTrackLength::IsPassed(), G4ITTrackHolder::KillTracks(), compTrackPerID::operator()(), G4StackManager::PopNextTrack(), G4Cerenkov::PostStepDoIt(), G4CoherentPairProduction::PostStepDoIt(), G4DNABrownianTransportation::PostStepDoIt(), G4DNAScavengerProcess::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4AdjointForcedInteractionForGamma::PostStepGetPhysicalInteractionLength(), G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength(), G4ITSteppingVerbose::PostStepVerbose(), G4ITSteppingVerbose::PreStepVerbose(), G4PSPopulation::ProcessHits(), G4StackManager::PushOneTrack(), G4DNAIRT::Sampling(), G4ITStepProcessor::Stepping(), G4ITTrackingInteractivity::TrackBanner(), G4VITSteppingVerbose::TrackBanner(), and G4ITSteppingVerbose::TrackingEnded().

◆ GetTrackLength()

◆ GetTrackStatus()

◆ GetUserInformation()

◆ GetVelocity()

◆ GetVertexKineticEnergy()

G4double G4Track::GetVertexKineticEnergy ( ) const
inline

◆ GetVertexMomentumDirection()

const G4ThreeVector & G4Track::GetVertexMomentumDirection ( ) const
inline

◆ GetVertexPosition()

const G4ThreeVector & G4Track::GetVertexPosition ( ) const
inline

◆ GetVolume()

G4VPhysicalVolume * G4Track::GetVolume ( ) const
inline

Referenced by G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4FastSimulationManagerProcess::AlongStepGetPhysicalInteractionLength(), G4ImportanceProcess::AlongStepGetPhysicalInteractionLength(), G4ITTransportation::AlongStepGetPhysicalInteractionLength(), G4ParallelGeometriesLimiterProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldScoringProcess::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4WeightCutOffProcess::AlongStepGetPhysicalInteractionLength(), G4WeightWindowProcess::AlongStepGetPhysicalInteractionLength(), G4DecayWithSpin::AtRestDoIt(), G4FastSimulationManagerProcess::AtRestGetPhysicalInteractionLength(), G4DNABrownianTransportation::ComputeGeomLimit(), G4RadioactiveDecay::DecayIt(), G4VRadioactiveDecay::DecayIt(), G4HadronicProcess::DumpState(), G4Channeling::GetMeanFreePath(), G4CoherentPairProduction::GetMeanFreePath(), G4SynchrotronRadiation::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetMeanFreePath(), G4VTransitionRadiation::GetMeanFreePath(), G4VXTRenergyLoss::GetMeanFreePath(), G4XrayReflection::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetPhotonEnergy(), G4BOptnForceCommonTruncatedExp::Initialize(), G4Channeling::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4VTransitionRadiation::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength(), G4FastSimulationManagerProcess::PostStepGetPhysicalInteractionLength(), G4LowECapture::PostStepGetPhysicalInteractionLength(), G4MaxTimeCuts::PostStepGetPhysicalInteractionLength(), G4MinEkineCuts::PostStepGetPhysicalInteractionLength(), G4StepLimiter::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4TransportationLogger::ReportLoopingTrack(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4PolarizedIonisationModel::SampleSecondaries(), and G4VPhononProcess::StartTracking().

◆ GetWeight()

◆ HasParentResonance()

G4bool G4Track::HasParentResonance ( ) const
inline

◆ IncrementCurrentStepNumber()

void G4Track::IncrementCurrentStepNumber ( )
inline

◆ IsBelowThreshold()

G4bool G4Track::IsBelowThreshold ( ) const
inline

◆ IsGoodForTracking()

G4bool G4Track::IsGoodForTracking ( ) const
inline

◆ 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 85 of file G4Track.cc.

86{
87 if(this != &right)
88 {
89 fPosition = right.fPosition;
90 fGlobalTime = right.fGlobalTime;
91 fLocalTime = right.fLocalTime;
92 fTrackLength = right.fTrackLength;
93 fWeight = right.fWeight;
94 fStepLength = right.fStepLength;
95
96 // additional fields required for geometrical splitting
97 if(fCopyTouchables) {
98 fpTouchable = right.fpTouchable;
99 fpNextTouchable = right.fpNextTouchable;
100 fpOriginTouchable = right.fpOriginTouchable;
101 }
102
103 // Track ID (and Parent ID) is not copied and set to zero for new track
104 fTrackID = 0;
105 fParentID = 0;
106
107 // CurrentStepNumber is set to be 0
108 fCurrentStepNumber = 0;
109
110 // Creator model ID
111 fCreatorModelID = right.fCreatorModelID;
112
113 // Parent resonance
114 fParentResonanceDef = right.fParentResonanceDef;
115 fParentResonanceID = right.fParentResonanceID;
116
117 // velocity information
118 fVelocity = right.fVelocity;
119
120 // dynamic particle information
121 delete fpDynamicParticle;
122 fpDynamicParticle = new G4DynamicParticle(*(right.fpDynamicParticle));
123
124 // track status and flags for tracking
125 fTrackStatus = right.fTrackStatus;
126 fBelowThreshold = right.fBelowThreshold;
127 fGoodForTracking = right.fGoodForTracking;
128
129 // Step information (Step Length, Step Number, pointer to the Step,)
130 // are not copied
131 fpStep = nullptr;
132
133 // vertex information
134 fVtxPosition = right.fVtxPosition;
135 fpLVAtVertex = right.fpLVAtVertex;
136 fVtxKineticEnergy = right.fVtxKineticEnergy;
137 fVtxMomentumDirection = right.fVtxMomentumDirection;
138
139 // CreatorProcess and UserInformation are not copied
140 fpCreatorProcess = nullptr;
141 delete fpUserInformation;
142 fpUserInformation = nullptr;
143
144 prev_mat = right.prev_mat;
145 groupvel = right.groupvel;
146 prev_velocity = right.prev_velocity;
147 prev_momentum = right.prev_momentum;
148
149 is_OpticalPhoton = right.is_OpticalPhoton;
150 useGivenVelocity = right.useGivenVelocity;
151
152 ClearAuxiliaryTrackInformation();
153 }
154 return *this;
155}

◆ operator==()

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

◆ RemoveAuxiliaryTrackInformation() [1/2]

void G4Track::RemoveAuxiliaryTrackInformation ( G4int id)

Definition at line 247 of file G4Track.cc.

248{
249 if(fpAuxiliaryTrackInformationMap != nullptr &&
250 fpAuxiliaryTrackInformationMap->find(id) != fpAuxiliaryTrackInformationMap->cend())
251 {
252 fpAuxiliaryTrackInformationMap->erase(id);
253 }
254}

Referenced by RemoveAuxiliaryTrackInformation().

◆ RemoveAuxiliaryTrackInformation() [2/2]

void G4Track::RemoveAuxiliaryTrackInformation ( G4String & name)

Definition at line 257 of file G4Track.cc.

258{
259 if(fpAuxiliaryTrackInformationMap != nullptr)
260 {
263 }
264}
int G4int
Definition G4Types.hh:85
static G4int GetModelID(const G4int modelIndex)
void RemoveAuxiliaryTrackInformation(G4int id)
Definition G4Track.cc:247

◆ SetAuxiliaryTrackInformation()

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

Definition at line 215 of file G4Track.cc.

217{
218 if(fpAuxiliaryTrackInformationMap == nullptr)
219 {
220 fpAuxiliaryTrackInformationMap =
221 new std::map<G4int, G4VAuxiliaryTrackInformation*>;
222 }
224 {
226 ED << "Process/model ID <" << id << "> is invalid.";
227 G4Exception("G4VAuxiliaryTrackInformation::G4VAuxiliaryTrackInformation()",
228 "TRACK0982", FatalException, ED);
229 }
230 (*fpAuxiliaryTrackInformationMap)[id] = info;
231}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
static G4int GetModelIndex(const G4int modelID)

Referenced by G4eplusAnnihilation::AtRestDoIt().

◆ SetBelowThresholdFlag()

void G4Track::SetBelowThresholdFlag ( G4bool value = true)
inline

◆ SetCreatorModelID()

◆ SetCreatorProcess()

void G4Track::SetCreatorProcess ( const G4VProcess * aValue)
inline

◆ SetGlobalTime()

void G4Track::SetGlobalTime ( const G4double aValue)
inline

◆ SetGoodForTrackingFlag()

◆ SetKineticEnergy()

◆ SetLocalTime()

void G4Track::SetLocalTime ( const G4double aValue)
inline

◆ SetLogicalVolumeAtVertex()

void G4Track::SetLogicalVolumeAtVertex ( const G4LogicalVolume * )
inline

◆ SetMomentumDirection()

void G4Track::SetMomentumDirection ( const G4ThreeVector & aValue)
inline

◆ SetNextTouchableHandle()

void G4Track::SetNextTouchableHandle ( const G4TouchableHandle & apValue)
inline

◆ SetOriginTouchableHandle()

void G4Track::SetOriginTouchableHandle ( const G4TouchableHandle & apValue)
inline

◆ SetParentID()

◆ SetParentResonanceDef()

void G4Track::SetParentResonanceDef ( const G4ParticleDefinition * parent)
inline

◆ SetParentResonanceID()

void G4Track::SetParentResonanceID ( const G4int parentID)
inline

◆ SetPolarization()

void G4Track::SetPolarization ( const G4ThreeVector & aValue)
inline

◆ SetPosition()

◆ SetProperTime()

void G4Track::SetProperTime ( const G4double aValue)
inline

◆ SetStep()

void G4Track::SetStep ( const G4Step * aValue)
inline

◆ SetStepLength()

void G4Track::SetStepLength ( G4double value)
inline

◆ SetTouchableHandle()

◆ SetTrackID()

void G4Track::SetTrackID ( const G4int aValue)
inline

◆ SetTrackStatus()

◆ SetUserInformation()

void G4Track::SetUserInformation ( G4VUserTrackInformation * aValue) const
inline

◆ SetVelocity()

void G4Track::SetVelocity ( G4double val)
inline

◆ SetVertexKineticEnergy()

void G4Track::SetVertexKineticEnergy ( const G4double aValue)
inline

◆ SetVertexMomentumDirection()

void G4Track::SetVertexMomentumDirection ( const G4ThreeVector & aValue)
inline

◆ SetVertexPosition()

void G4Track::SetVertexPosition ( const G4ThreeVector & aValue)
inline

◆ SetWeight()

◆ UseGivenVelocity() [1/2]

G4bool G4Track::UseGivenVelocity ( ) const
inline

◆ UseGivenVelocity() [2/2]

void G4Track::UseGivenVelocity ( G4bool val)
inline

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