Geant4 11.2.2
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 40 of file G4FieldTrackUpdator.hh.

Member Function Documentation

◆ CreateFieldTrack()

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

Definition at line 40 of file G4FieldTrackUpdator.cc.

41{
42 return new G4FieldTrack(
43 trk->GetPosition(), trk->GetGlobalTime(), trk->GetMomentumDirection(),
47 0.0 // magnetic dipole moment to be implemented
48 );
49}
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 52 of file G4FieldTrackUpdator.cc.

53{
54 const G4DynamicParticle* ptDynamicParticle = trk->GetDynamicParticle();
55
56 // The following properties must be updated 1) for each new track, and
57 ftrk->SetRestMass(ptDynamicParticle->GetMass());
58 // 2) Since ion can lose/gain electrons, this must be done at every step
59
60 ftrk->UpdateState(trk->GetPosition(), trk->GetGlobalTime(),
62
63#ifdef G4CHECK
64 if((trk->GetMomentum() - ftrk->GetMomentum()).mag2() >
65 1.e-16 * trk->GetMomentum().mag2())
66 {
67 G4cerr << "ERROR> G4FieldTrackUpdator sees *Disagreement* in momentum "
68 << G4endl;
69 G4cout << " FTupdator: Tracking Momentum= " << trk->GetMomentum()
70 << G4endl;
71 G4cout << " FTupdator: FldTrack Momentum= " << ftrk->GetMomentum()
72 << G4endl;
73 G4cout << " FTupdator: FldTrack-Tracking= "
74 << ftrk->GetMomentum() - trk->GetMomentum() << G4endl;
75 }
76#endif
77
79
80 ftrk->SetChargeAndMoments(ptDynamicParticle->GetCharge(),
81 ptDynamicParticle->GetMagneticMoment());
82 ftrk->SetPDGSpin(ptDynamicParticle->GetParticleDefinition()->GetPDGSpin());
83 // The charge can change during tracking
84 ftrk->SetSpin(ptDynamicParticle->GetPolarization());
85}
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double mag2() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetMagneticMoment() const
void SetProperTimeOfFlight(G4double tofProper)
void SetPDGSpin(G4double pdgSpin)
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 SetSpin(const G4ThreeVector &vSpin)
void SetRestMass(G4double Mass_c2)
G4ThreeVector GetMomentum() const
G4double GetProperTime() const
G4ThreeVector GetMomentum() const

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


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