Geant4-11
Static Public Member Functions | Private Member Functions
G4FieldTrackUpdator Class Reference

#include <G4FieldTrackUpdator.hh>

Static Public Member Functions

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

Private Member Functions

 G4FieldTrackUpdator ()
 
 ~G4FieldTrackUpdator ()
 

Detailed Description

Definition at line 40 of file G4FieldTrackUpdator.hh.

Constructor & Destructor Documentation

◆ G4FieldTrackUpdator()

G4FieldTrackUpdator::G4FieldTrackUpdator ( )
inlineprivate

Definition at line 49 of file G4FieldTrackUpdator.hh.

49{}

◆ ~G4FieldTrackUpdator()

G4FieldTrackUpdator::~G4FieldTrackUpdator ( )
inlineprivate

Definition at line 50 of file G4FieldTrackUpdator.hh.

50{}

Member Function Documentation

◆ CreateFieldTrack()

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

Definition at line 40 of file G4FieldTrackUpdator.cc.

41{
42 G4FieldTrack* ftrk = new G4FieldTrack(
43 trk->GetPosition(), trk->GetGlobalTime(), trk->GetMomentumDirection(),
47 0.0 // magnetic dipole moment to be implemented
48 );
49 return ftrk;
50}
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

References G4DynamicParticle::GetCharge(), G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4Track::GetKineticEnergy(), G4DynamicParticle::GetMass(), G4Track::GetMomentumDirection(), G4DynamicParticle::GetPolarization(), and G4Track::GetPosition().

◆ Update()

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

Definition at line 53 of file G4FieldTrackUpdator.cc.

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

References G4cerr, G4cout, G4endl, G4DynamicParticle::GetCharge(), G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4Track::GetKineticEnergy(), G4DynamicParticle::GetMagneticMoment(), G4DynamicParticle::GetMass(), G4FieldTrack::GetMomentum(), G4Track::GetMomentum(), G4Track::GetMomentumDirection(), G4DynamicParticle::GetParticleDefinition(), G4ParticleDefinition::GetPDGSpin(), G4DynamicParticle::GetPolarization(), G4Track::GetPosition(), G4Track::GetProperTime(), CLHEP::Hep3Vector::mag2(), G4FieldTrack::SetChargeAndMoments(), G4FieldTrack::SetPDGSpin(), G4FieldTrack::SetProperTimeOfFlight(), G4FieldTrack::SetRestMass(), G4FieldTrack::SetSpin(), and G4FieldTrack::UpdateState().

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


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