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

#include <G4FieldTrackUpdator.hh>

Static Public Member Functions

static G4FieldTrackCreateFieldTrack (const G4Track *)
 
static void Update (G4FieldTrack *, const G4Track *)
 

Detailed Description

Definition at line 43 of file G4FieldTrackUpdator.hh.

Member Function Documentation

◆ CreateFieldTrack()

G4FieldTrack * G4FieldTrackUpdator::CreateFieldTrack ( const G4Track trk)
static

Definition at line 45 of file G4FieldTrackUpdator.cc.

46{
47 G4FieldTrack* ftrk = new G4FieldTrack(
48 trk->GetPosition(),
49 trk->GetGlobalTime(),
51 trk->GetKineticEnergy(),
55 0.0 // magnetic dipole moment to be implemented
56 );
57 return ftrk;
58}
G4double GetMass() const
G4double GetCharge() const
const G4ThreeVector & GetPolarization() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const

◆ Update()

void G4FieldTrackUpdator::Update ( G4FieldTrack ftrk,
const G4Track trk 
)
static

Definition at line 60 of file G4FieldTrackUpdator.cc.

61{
62 ftrk->UpdateState(
63 trk->GetPosition(),
64 trk->GetGlobalTime(),
66 trk->GetKineticEnergy()
67 );
68 const G4DynamicParticle* ptDynamicParticle= trk->GetDynamicParticle();
69
70 ftrk->SetChargeAndMoments( ptDynamicParticle->GetCharge() );
71 // The charge can change during tracking
72 ftrk->SetSpin( ptDynamicParticle->GetPolarization() );
73
74 // The following properties must be updated ONCE for each new track (at least)
75 ftrk->SetRestMass(ptDynamicParticle->GetMass());
76}
void SetSpin(G4ThreeVector nSpin)
void SetChargeAndMoments(G4double charge, G4double magnetic_dipole_moment=DBL_MAX, G4double electric_dipole_moment=DBL_MAX, G4double magnetic_charge=DBL_MAX)
void UpdateState(const G4ThreeVector &pPosition, G4double LaboratoryTimeOfFlight, const G4ThreeVector &pMomentumDirection, G4double kineticEnergy)
void SetRestMass(G4double Mass_c2)

Referenced by G4ImportanceProcess::AlongStepGetPhysicalInteractionLength(), G4WeightCutOffProcess::AlongStepGetPhysicalInteractionLength(), G4WeightWindowProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldScoringProcess::AlongStepGetPhysicalInteractionLength(), and G4FastSimulationManagerProcess::AlongStepGetPhysicalInteractionLength().


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