Geant4-11
G4FieldTrackUpdator.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4FieldTrackUpdator class implementation
27//
28// Author: M.Asai, 28 April 2006
29// --------------------------------------------------------------------
30
31#include "globals.hh"
33#include "G4ThreeVector.hh"
34#include "G4Track.hh"
35#include "G4DynamicParticle.hh"
36#include "G4TrackStatus.hh"
37#include "G4FieldTrack.hh"
38
39//---------------------------------------------------------------------
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}
51
52//---------------------------------------------------------------------
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
G4double GetMass() const
G4double GetCharge() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetMagneticMoment() const
const G4ThreeVector & GetPolarization() const
static void Update(G4FieldTrack *, const G4Track *)
static G4FieldTrack * CreateFieldTrack(const G4Track *)
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
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4double GetProperTime() const
G4ThreeVector GetMomentum() const
const G4DynamicParticle * GetDynamicParticle() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const