BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Public Member Functions | Data Fields | Static Public Attributes | Private Member Functions | Private Attributes | Static Private Attributes | Friends
BDSTrajectoryPoint Class Reference

A Point in a trajectory with extra information. More...

#include <BDSTrajectoryPoint.hh>

Inheritance diagram for BDSTrajectoryPoint:
Inheritance graph
Collaboration diagram for BDSTrajectoryPoint:
Collaboration graph

Public Member Functions

 BDSTrajectoryPoint ()
 Default constructor. More...
 
 BDSTrajectoryPoint (const G4Step *step, G4bool storeExtrasLocalIn, G4bool storeExtrasLinkIn, G4bool storeExtrasIonIn)
 
 BDSTrajectoryPoint (const G4Track *track, G4bool storeExtrasLocalIn, G4bool storeExtrasLinkIn, G4bool storeExtrasIonIn)
 
 BDSTrajectoryPoint (const BDSTrajectoryPoint &other)
 Implement copy constructor as we have pointers we own. More...
 
void * operator new (size_t)
 
void operator delete (void *aTrajectoryPoint)
 
int operator== (const BDSTrajectoryPoint &right) const
 
void DeleteExtraLocal ()
 
void DeleteExtraLinks ()
 
void DeleteExtraIon ()
 
G4bool IsScatteringPoint () const
 
G4bool NotTransportationLimitedStep () const
 Return true if step isn't defined by transportation processes. More...
 
G4int GetPreProcessType () const
 Accessor. More...
 
G4int GetPreProcessSubType () const
 Accessor. More...
 
G4int GetPostProcessType () const
 Accessor. More...
 
G4int GetPostProcessSubType () const
 Accessor. More...
 
G4double GetPreWeight () const
 Accessor. More...
 
G4double GetPostWeight () const
 Accessor. More...
 
G4double GetPreEnergy () const
 Accessor. More...
 
G4double GetPostEnergy () const
 Accessor. More...
 
G4double GetEnergyDeposit () const
 Accessor. More...
 
G4ThreeVector GetPreMomentum () const
 Accessor. More...
 
G4ThreeVector GetPostMomentum () const
 Accessor. More...
 
G4double GetPreS () const
 Accessor. More...
 
G4double GetPostS () const
 Accessor. More...
 
G4double GetPreGlobalTime () const
 Accessor. More...
 
G4double GetPostGlobalTime () const
 Accessor. More...
 
G4int GetBeamLineIndex () const
 Accessor. More...
 
BDSBeamlineGetBeamLine () const
 Accessor. More...
 
G4ThreeVector GetPrePosLocal () const
 Accessor. More...
 
G4ThreeVector GetPostPosLocal () const
 Accessor. More...
 
G4Material * GetMaterial () const
 Accessor. More...
 
void SetMaterial (G4Material *materialIn)
 
G4ThreeVector GetPositionLocal () const
 Accessor for the extra information local. More...
 
G4ThreeVector GetMomentumLocal () const
 Accessor. More...
 
G4int GetCharge () const
 Accessor for the extra information links. More...
 
G4double GetKineticEnergy () const
 Accessor for the extra information links. More...
 
G4int GetTurnsTaken () const
 Accessor for the extra information links. More...
 
G4double GetMass () const
 Accessor for the extra information links. More...
 
G4double GetRigidity () const
 Accessor for the extra information links. More...
 
G4bool GetIsIon () const
 Accessor for the extra information ions. More...
 
G4int GetIonA () const
 Accessor for the extra information ions. More...
 
G4int GetIonZ () const
 Accessor for the extra information ions. More...
 
G4int GetNElectrons () const
 Accessor for the extra information ions. More...
 
G4double PrePosR () const
 Return the transverse local radius in x,y. More...
 
G4double PostPosR () const
 Return the transverse local radius in x,y. More...
 
G4bool operator< (const BDSTrajectoryPoint &other) const
 Comparison operator. More...
 
G4bool operator> (const BDSTrajectoryPoint &other) const
 Comparison operator. More...
 
G4bool operator<= (const BDSTrajectoryPoint &other) const
 Comparison operator. More...
 
G4bool operator>= (const BDSTrajectoryPoint &other) const
 Comparison operator. More...
 

Static Public Member Functions

static G4bool IsScatteringPoint (const G4Step *step)
 Static function to determine whether a step corresponds to scattering. More...
 
static G4bool IsScatteringPoint (G4int postProcessType, G4int postProcessSubType, G4double totalEnergyDeposit)
 Static function to determine whether a step corresponds to scattering. More...
 

Data Fields

BDSTrajectoryPointLocalextraLocal
 
BDSTrajectoryPointLinkextraLink
 
BDSTrajectoryPointIonextraIon
 

Static Public Attributes

static G4double dEThresholdForScattering = 1e-8
 

Private Member Functions

void InitialiseVariables ()
 
void StoreExtrasLink (const G4Track *track)
 Utility function to prepare and fill extra link variables. More...
 
void StoreExtrasIon (const G4Track *track)
 Utility function to prepare and fill extra ion variables. More...
 

Private Attributes

G4int preProcessType
 Process type of pre-step point. More...
 
G4int preProcessSubType
 Process sub type of pre-step point. More...
 
G4int postProcessType
 Process type of post step point. More...
 
G4int postProcessSubType
 Process sub type of post step point. More...
 
G4double preWeight
 Weight associated with pre-step point. More...
 
G4double postWeight
 Weight associated with post step point. More...
 
G4double preEnergy
 Kinetic energy of pre-step point. More...
 
G4double postEnergy
 Kinetic energy of post step point. More...
 
G4ThreeVector preMomentum
 Momentum of pre-step point. More...
 
G4ThreeVector postMomentum
 Momentum of post-step point. More...
 
G4double energyDeposit
 Total energy deposited during step. More...
 
G4double preS
 Global curvilinear S coordinate of pre-step point. More...
 
G4double postS
 Global curvilinear S coordinate of post step point. More...
 
G4double preGlobalTime
 Time since event started of pre-step point. More...
 
G4double postGlobalTime
 Time since event started of post-step point. More...
 
G4int beamlineIndex
 Index to beam line element in the mass world beam line. More...
 
BDSBeamlinebeamline
 Beam line (if any) point belongs to (always mass world). More...
 
G4ThreeVector prePosLocal
 Local coordinates of pre-step point. More...
 
G4ThreeVector postPosLocal
 Local coordinates of post-step point. More...
 
G4Material * material
 Material point for pre-step point. More...
 

Static Private Attributes

static BDSAuxiliaryNavigatorauxNavigator = new BDSAuxiliaryNavigator()
 

Friends

std::ostream & operator<< (std::ostream &out, BDSTrajectoryPoint const &p)
 Output stream. More...
 

Detailed Description

A Point in a trajectory with extra information.

Author
S. Boogert

Definition at line 44 of file BDSTrajectoryPoint.hh.

Constructor & Destructor Documentation

◆ BDSTrajectoryPoint() [1/4]

BDSTrajectoryPoint::BDSTrajectoryPoint ( )

Default constructor.

Definition at line 56 of file BDSTrajectoryPoint.cc.

References InitialiseVariables().

Here is the call graph for this function:

◆ BDSTrajectoryPoint() [2/4]

BDSTrajectoryPoint::BDSTrajectoryPoint ( const G4Step *  step,
G4bool  storeExtrasLocalIn,
G4bool  storeExtrasLinkIn,
G4bool  storeExtrasIonIn 
)

◆ BDSTrajectoryPoint() [3/4]

BDSTrajectoryPoint::BDSTrajectoryPoint ( const G4Track *  track,
G4bool  storeExtrasLocalIn,
G4bool  storeExtrasLinkIn,
G4bool  storeExtrasIonIn 
)

◆ BDSTrajectoryPoint() [4/4]

BDSTrajectoryPoint::BDSTrajectoryPoint ( const BDSTrajectoryPoint other)

◆ ~BDSTrajectoryPoint()

BDSTrajectoryPoint::~BDSTrajectoryPoint ( )
virtual

Definition at line 218 of file BDSTrajectoryPoint.cc.

Member Function Documentation

◆ DeleteExtraIon()

void BDSTrajectoryPoint::DeleteExtraIon ( )
inline

Definition at line 75 of file BDSTrajectoryPoint.hh.

◆ DeleteExtraLinks()

void BDSTrajectoryPoint::DeleteExtraLinks ( )
inline

Definition at line 74 of file BDSTrajectoryPoint.hh.

◆ DeleteExtraLocal()

void BDSTrajectoryPoint::DeleteExtraLocal ( )
inline

Definition at line 73 of file BDSTrajectoryPoint.hh.

◆ GetBeamLine()

BDSBeamline * BDSTrajectoryPoint::GetBeamLine ( ) const
inline

Accessor.

Definition at line 111 of file BDSTrajectoryPoint.hh.

References beamline.

◆ GetBeamLineIndex()

G4int BDSTrajectoryPoint::GetBeamLineIndex ( ) const
inline

Accessor.

Definition at line 110 of file BDSTrajectoryPoint.hh.

References beamlineIndex.

◆ GetCharge()

G4int BDSTrajectoryPoint::GetCharge ( ) const
inline

Accessor for the extra information links.

Definition at line 126 of file BDSTrajectoryPoint.hh.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetEnergyDeposit()

G4double BDSTrajectoryPoint::GetEnergyDeposit ( ) const
inline

Accessor.

Definition at line 103 of file BDSTrajectoryPoint.hh.

References energyDeposit.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetIonA()

G4int BDSTrajectoryPoint::GetIonA ( ) const
inline

Accessor for the extra information ions.

Definition at line 135 of file BDSTrajectoryPoint.hh.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetIonZ()

G4int BDSTrajectoryPoint::GetIonZ ( ) const
inline

Accessor for the extra information ions.

Definition at line 136 of file BDSTrajectoryPoint.hh.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetIsIon()

G4bool BDSTrajectoryPoint::GetIsIon ( ) const
inline

Accessor for the extra information ions.

Definition at line 134 of file BDSTrajectoryPoint.hh.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetKineticEnergy()

G4double BDSTrajectoryPoint::GetKineticEnergy ( ) const
inline

Accessor for the extra information links.

Definition at line 127 of file BDSTrajectoryPoint.hh.

References preEnergy.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetMass()

G4double BDSTrajectoryPoint::GetMass ( ) const
inline

Accessor for the extra information links.

Definition at line 129 of file BDSTrajectoryPoint.hh.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetMaterial()

G4Material * BDSTrajectoryPoint::GetMaterial ( ) const
inline

Accessor.

Definition at line 114 of file BDSTrajectoryPoint.hh.

References material.

Referenced by BDSTrajectory::AppendStep(), and BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetMomentumLocal()

G4ThreeVector BDSTrajectoryPoint::GetMomentumLocal ( ) const
inline

Accessor.

Definition at line 122 of file BDSTrajectoryPoint.hh.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetNElectrons()

G4int BDSTrajectoryPoint::GetNElectrons ( ) const
inline

Accessor for the extra information ions.

Definition at line 137 of file BDSTrajectoryPoint.hh.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetPositionLocal()

G4ThreeVector BDSTrajectoryPoint::GetPositionLocal ( ) const
inline

Accessor for the extra information local.

Definition at line 121 of file BDSTrajectoryPoint.hh.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetPostEnergy()

G4double BDSTrajectoryPoint::GetPostEnergy ( ) const
inline

Accessor.

Definition at line 102 of file BDSTrajectoryPoint.hh.

References postEnergy.

◆ GetPostGlobalTime()

G4double BDSTrajectoryPoint::GetPostGlobalTime ( ) const
inline

Accessor.

Definition at line 109 of file BDSTrajectoryPoint.hh.

References postGlobalTime.

◆ GetPostMomentum()

G4ThreeVector BDSTrajectoryPoint::GetPostMomentum ( ) const
inline

Accessor.

Definition at line 105 of file BDSTrajectoryPoint.hh.

References postMomentum.

◆ GetPostPosLocal()

G4ThreeVector BDSTrajectoryPoint::GetPostPosLocal ( ) const
inline

Accessor.

Definition at line 113 of file BDSTrajectoryPoint.hh.

References postPosLocal.

◆ GetPostProcessSubType()

G4int BDSTrajectoryPoint::GetPostProcessSubType ( ) const
inline

Accessor.

Definition at line 98 of file BDSTrajectoryPoint.hh.

References postProcessSubType.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory(), and IsScatteringPoint().

Here is the caller graph for this function:

◆ GetPostProcessType()

G4int BDSTrajectoryPoint::GetPostProcessType ( ) const
inline

Accessor.

Definition at line 97 of file BDSTrajectoryPoint.hh.

References postProcessType.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory(), and IsScatteringPoint().

Here is the caller graph for this function:

◆ GetPostS()

G4double BDSTrajectoryPoint::GetPostS ( ) const
inline

Accessor.

Definition at line 107 of file BDSTrajectoryPoint.hh.

References postS.

◆ GetPostWeight()

G4double BDSTrajectoryPoint::GetPostWeight ( ) const
inline

Accessor.

Definition at line 100 of file BDSTrajectoryPoint.hh.

References postWeight.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetPreEnergy()

G4double BDSTrajectoryPoint::GetPreEnergy ( ) const
inline

Accessor.

Definition at line 101 of file BDSTrajectoryPoint.hh.

References preEnergy.

◆ GetPreGlobalTime()

G4double BDSTrajectoryPoint::GetPreGlobalTime ( ) const
inline

Accessor.

Definition at line 108 of file BDSTrajectoryPoint.hh.

References preGlobalTime.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetPreMomentum()

G4ThreeVector BDSTrajectoryPoint::GetPreMomentum ( ) const
inline

Accessor.

Definition at line 104 of file BDSTrajectoryPoint.hh.

References preMomentum.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetPrePosLocal()

G4ThreeVector BDSTrajectoryPoint::GetPrePosLocal ( ) const
inline

Accessor.

Definition at line 112 of file BDSTrajectoryPoint.hh.

References prePosLocal.

◆ GetPreProcessSubType()

G4int BDSTrajectoryPoint::GetPreProcessSubType ( ) const
inline

Accessor.

Definition at line 96 of file BDSTrajectoryPoint.hh.

References preProcessSubType.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetPreProcessType()

G4int BDSTrajectoryPoint::GetPreProcessType ( ) const
inline

Accessor.

Definition at line 95 of file BDSTrajectoryPoint.hh.

References preProcessType.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetPreS()

G4double BDSTrajectoryPoint::GetPreS ( ) const
inline

Accessor.

Definition at line 106 of file BDSTrajectoryPoint.hh.

References preS.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory(), and IsScatteringPoint().

Here is the caller graph for this function:

◆ GetPreWeight()

G4double BDSTrajectoryPoint::GetPreWeight ( ) const
inline

Accessor.

Definition at line 99 of file BDSTrajectoryPoint.hh.

References preWeight.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetRigidity()

G4double BDSTrajectoryPoint::GetRigidity ( ) const
inline

Accessor for the extra information links.

Definition at line 130 of file BDSTrajectoryPoint.hh.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ GetTurnsTaken()

G4int BDSTrajectoryPoint::GetTurnsTaken ( ) const
inline

Accessor for the extra information links.

Definition at line 128 of file BDSTrajectoryPoint.hh.

Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory().

Here is the caller graph for this function:

◆ InitialiseVariables()

void BDSTrajectoryPoint::InitialiseVariables ( )
private

Initialisation of variables in separate function to reduce duplication in multiple constructors.

Definition at line 225 of file BDSTrajectoryPoint.cc.

References beamline, beamlineIndex, energyDeposit, material, postEnergy, postGlobalTime, postMomentum, postPosLocal, postProcessSubType, postProcessType, postS, postWeight, preEnergy, preGlobalTime, preMomentum, prePosLocal, preProcessSubType, preProcessType, preS, and preWeight.

Referenced by BDSTrajectoryPoint().

Here is the caller graph for this function:

◆ IsScatteringPoint() [1/3]

G4bool BDSTrajectoryPoint::IsScatteringPoint ( ) const

Check to see if point is a scattering point (from a physics point of view). Uses static functions defined below. This is defined by whether the processes that defined the step length is non-transportation or the energy loss along step is greater than 1e-9 MeV, which could be due to continuous energy loss processes along step that might not limit the step length.

Definition at line 278 of file BDSTrajectoryPoint.cc.

References energyDeposit, GetPostProcessSubType(), GetPostProcessType(), GetPreS(), BDSProcessMap::GetProcessName(), BDSProcessMap::Instance(), IsScatteringPoint(), postProcessSubType, and postProcessType.

Referenced by BDSTrajectoryPrimary::AppendStep(), BDSTrajectory::FirstInteraction(), IsScatteringPoint(), BDSSDCollimator::ProcessHitsOrdered(), and BDSSDThinThing::ProcessHitsOrdered().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsScatteringPoint() [2/3]

G4bool BDSTrajectoryPoint::IsScatteringPoint ( const G4Step *  step)
static

Static function to determine whether a step corresponds to scattering.

Definition at line 322 of file BDSTrajectoryPoint.cc.

References IsScatteringPoint(), postProcessSubType, and postProcessType.

Here is the call graph for this function:

◆ IsScatteringPoint() [3/3]

G4bool BDSTrajectoryPoint::IsScatteringPoint ( G4int  postProcessType,
G4int  postProcessSubType,
G4double  totalEnergyDeposit 
)
static

Static function to determine whether a step corresponds to scattering.

Definition at line 344 of file BDSTrajectoryPoint.cc.

References dEThresholdForScattering, postProcessSubType, and postProcessType.

◆ NotTransportationLimitedStep()

G4bool BDSTrajectoryPoint::NotTransportationLimitedStep ( ) const

Return true if step isn't defined by transportation processes.

Definition at line 295 of file BDSTrajectoryPoint.cc.

References postProcessType, and preProcessType.

Referenced by BDSTrajectory::AppendStep().

Here is the caller graph for this function:

◆ operator delete()

void BDSTrajectoryPoint::operator delete ( void *  aTrajectoryPoint)
inline

Definition at line 213 of file BDSTrajectoryPoint.hh.

◆ operator new()

void * BDSTrajectoryPoint::operator new ( size_t  )
inline

Definition at line 206 of file BDSTrajectoryPoint.hh.

◆ operator<()

G4bool BDSTrajectoryPoint::operator< ( const BDSTrajectoryPoint other) const
inline

Comparison operator.

Definition at line 218 of file BDSTrajectoryPoint.hh.

References preGlobalTime, and preS.

◆ operator<=()

G4bool BDSTrajectoryPoint::operator<= ( const BDSTrajectoryPoint other) const
inline

Comparison operator.

Definition at line 152 of file BDSTrajectoryPoint.hh.

◆ operator==()

int BDSTrajectoryPoint::operator== ( const BDSTrajectoryPoint right) const
inline

Definition at line 71 of file BDSTrajectoryPoint.hh.

◆ operator>()

G4bool BDSTrajectoryPoint::operator> ( const BDSTrajectoryPoint other) const
inline

Comparison operator.

Definition at line 151 of file BDSTrajectoryPoint.hh.

◆ operator>=()

G4bool BDSTrajectoryPoint::operator>= ( const BDSTrajectoryPoint other) const
inline

Comparison operator.

Definition at line 153 of file BDSTrajectoryPoint.hh.

◆ PostPosR()

G4double BDSTrajectoryPoint::PostPosR ( ) const

Return the transverse local radius in x,y.

Definition at line 317 of file BDSTrajectoryPoint.cc.

References postPosLocal.

Referenced by BDSEventAction::IdentifyTrajectoriesForStorage().

Here is the caller graph for this function:

◆ PrePosR()

G4double BDSTrajectoryPoint::PrePosR ( ) const

Return the transverse local radius in x,y.

Definition at line 312 of file BDSTrajectoryPoint.cc.

References prePosLocal.

◆ SetMaterial()

void BDSTrajectoryPoint::SetMaterial ( G4Material *  materialIn)
inline

For initial points in a trajectory it is not possible to yet know the material so it has to be updated after the first step.

Definition at line 118 of file BDSTrajectoryPoint.hh.

References material.

◆ StoreExtrasIon()

void BDSTrajectoryPoint::StoreExtrasIon ( const G4Track *  track)
private

Utility function to prepare and fill extra ion variables.

Definition at line 265 of file BDSTrajectoryPoint.cc.

References BDS::IsIon().

Referenced by BDSTrajectoryPoint().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ StoreExtrasLink()

void BDSTrajectoryPoint::StoreExtrasLink ( const G4Track *  track)
private

Utility function to prepare and fill extra link variables.

Definition at line 252 of file BDSTrajectoryPoint.cc.

References BDSGlobalConstants::Instance(), BDS::IsFinite(), and BDS::Rigidity().

Referenced by BDSTrajectoryPoint().

Here is the call graph for this function:
Here is the caller graph for this function:

Friends And Related Function Documentation

◆ operator<<

std::ostream & operator<< ( std::ostream &  out,
BDSTrajectoryPoint const &  p 
)
friend

Output stream.

Definition at line 306 of file BDSTrajectoryPoint.cc.

Field Documentation

◆ auxNavigator

BDSAuxiliaryNavigator * BDSTrajectoryPoint::auxNavigator = new BDSAuxiliaryNavigator()
staticprivate

An auxiliary navigator to get curvilinear coordinates. Lots of points, but only need one navigator so make it static.

Definition at line 201 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint().

◆ beamline

BDSBeamline* BDSTrajectoryPoint::beamline
private

Beam line (if any) point belongs to (always mass world).

Definition at line 194 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetBeamLine(), and InitialiseVariables().

◆ beamlineIndex

G4int BDSTrajectoryPoint::beamlineIndex
private

Index to beam line element in the mass world beam line.

Definition at line 193 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetBeamLineIndex(), and InitialiseVariables().

◆ dEThresholdForScattering

G4double BDSTrajectoryPoint::dEThresholdForScattering = 1e-8
static

Threshold energy (in geant4 units) for considering a point a scattering point. In some cases the along step process such as multiple scattering won't show as the process that defined the step but may significantly degrade the energy. Use this value as a threshold over which we consider the step a 'scattering' one.

Definition at line 164 of file BDSTrajectoryPoint.hh.

Referenced by IsScatteringPoint().

◆ energyDeposit

G4double BDSTrajectoryPoint::energyDeposit
private

Total energy deposited during step.

Definition at line 188 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetEnergyDeposit(), InitialiseVariables(), and IsScatteringPoint().

◆ extraIon

BDSTrajectoryPointIon* BDSTrajectoryPoint::extraIon

Definition at line 158 of file BDSTrajectoryPoint.hh.

◆ extraLink

BDSTrajectoryPointLink* BDSTrajectoryPoint::extraLink

Definition at line 157 of file BDSTrajectoryPoint.hh.

◆ extraLocal

BDSTrajectoryPointLocal* BDSTrajectoryPoint::extraLocal

Definition at line 156 of file BDSTrajectoryPoint.hh.

◆ material

G4Material* BDSTrajectoryPoint::material
private

Material point for pre-step point.

Definition at line 197 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetMaterial(), InitialiseVariables(), and SetMaterial().

◆ postEnergy

G4double BDSTrajectoryPoint::postEnergy
private

Kinetic energy of post step point.

Definition at line 185 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPostEnergy(), and InitialiseVariables().

◆ postGlobalTime

G4double BDSTrajectoryPoint::postGlobalTime
private

Time since event started of post-step point.

Definition at line 192 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPostGlobalTime(), and InitialiseVariables().

◆ postMomentum

G4ThreeVector BDSTrajectoryPoint::postMomentum
private

Momentum of post-step point.

Definition at line 187 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPostMomentum(), and InitialiseVariables().

◆ postPosLocal

G4ThreeVector BDSTrajectoryPoint::postPosLocal
private

Local coordinates of post-step point.

Definition at line 196 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPostPosLocal(), InitialiseVariables(), and PostPosR().

◆ postProcessSubType

G4int BDSTrajectoryPoint::postProcessSubType
private

Process sub type of post step point.

Definition at line 180 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPostProcessSubType(), InitialiseVariables(), and IsScatteringPoint().

◆ postProcessType

G4int BDSTrajectoryPoint::postProcessType
private

Process type of post step point.

Definition at line 179 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPostProcessType(), InitialiseVariables(), IsScatteringPoint(), and NotTransportationLimitedStep().

◆ postS

G4double BDSTrajectoryPoint::postS
private

Global curvilinear S coordinate of post step point.

Definition at line 190 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPostS(), and InitialiseVariables().

◆ postWeight

G4double BDSTrajectoryPoint::postWeight
private

Weight associated with post step point.

Definition at line 183 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPostWeight(), and InitialiseVariables().

◆ preEnergy

G4double BDSTrajectoryPoint::preEnergy
private

Kinetic energy of pre-step point.

Definition at line 184 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetKineticEnergy(), GetPreEnergy(), and InitialiseVariables().

◆ preGlobalTime

G4double BDSTrajectoryPoint::preGlobalTime
private

Time since event started of pre-step point.

Definition at line 191 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPreGlobalTime(), InitialiseVariables(), and operator<().

◆ preMomentum

G4ThreeVector BDSTrajectoryPoint::preMomentum
private

Momentum of pre-step point.

Definition at line 186 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPreMomentum(), and InitialiseVariables().

◆ prePosLocal

G4ThreeVector BDSTrajectoryPoint::prePosLocal
private

Local coordinates of pre-step point.

Definition at line 195 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPrePosLocal(), InitialiseVariables(), and PrePosR().

◆ preProcessSubType

G4int BDSTrajectoryPoint::preProcessSubType
private

Process sub type of pre-step point.

Definition at line 178 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPreProcessSubType(), and InitialiseVariables().

◆ preProcessType

G4int BDSTrajectoryPoint::preProcessType
private

Process type of pre-step point.

Definition at line 177 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPreProcessType(), InitialiseVariables(), and NotTransportationLimitedStep().

◆ preS

G4double BDSTrajectoryPoint::preS
private

Global curvilinear S coordinate of pre-step point.

Definition at line 189 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPreS(), InitialiseVariables(), and operator<().

◆ preWeight

G4double BDSTrajectoryPoint::preWeight
private

Weight associated with pre-step point.

Definition at line 182 of file BDSTrajectoryPoint.hh.

Referenced by BDSTrajectoryPoint(), GetPreWeight(), and InitialiseVariables().


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