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

#include <G4FieldTrack.hh>

Public Types

enum  { ncompSVEC = 12 }
 

Public Member Functions

 G4FieldTrack (const G4ThreeVector &pPosition, G4double LaboratoryTimeOfFlight, const G4ThreeVector &pMomentumDirection, G4double kineticEnergy, G4double restMass_c2, G4double charge, const G4ThreeVector &polarization, G4double magnetic_dipole_moment=0.0, G4double curve_length=0.0, G4double PDGspin=-1.0)
 
 G4FieldTrack (char)
 
 G4FieldTrack (const G4ThreeVector &pPosition, const G4ThreeVector &pMomentumDirection, G4double curve_length, G4double kineticEnergy, const G4double restMass_c2, G4double velocity, G4double LaboratoryTimeOfFlight=0.0, G4double ProperTimeOfFlight=0.0, const G4ThreeVector *pPolarization=nullptr, G4double PDGspin=-1.0)
 
 ~G4FieldTrack ()
 
 G4FieldTrack (const G4FieldTrack &pFieldTrack)
 
G4FieldTrackoperator= (const G4FieldTrack &rStVec)
 
 G4FieldTrack (G4FieldTrack &&from)
 
G4FieldTrackoperator= (G4FieldTrack &&from)
 
void UpdateState (const G4ThreeVector &pPosition, G4double LaboratoryTimeOfFlight, const G4ThreeVector &pMomentumDirection, G4double kineticEnergy)
 
void UpdateFourMomentum (G4double kineticEnergy, const G4ThreeVector &momentumDirection)
 
void SetChargeAndMoments (G4double charge, G4double magnetic_dipole_moment=DBL_MAX, G4double electric_dipole_moment=DBL_MAX, G4double magnetic_charge=DBL_MAX)
 
void SetPDGSpin (G4double pdgSpin)
 
G4double GetPDGSpin ()
 
G4ThreeVector GetMomentum () const
 
G4ThreeVector GetPosition () const
 
const G4ThreeVectorGetMomentumDir () const
 
G4ThreeVector GetMomentumDirection () const
 
G4double GetCurveLength () const
 
G4ThreeVector GetPolarization () const
 
void SetPolarization (const G4ThreeVector &vecPol)
 
const G4ChargeStateGetChargeState () const
 
G4double GetLabTimeOfFlight () const
 
G4double GetProperTimeOfFlight () const
 
G4double GetKineticEnergy () const
 
G4double GetCharge () const
 
G4double GetRestMass () const
 
void SetPosition (G4ThreeVector nPos)
 
void SetMomentum (G4ThreeVector nMomDir)
 
void SetMomentumDir (G4ThreeVector nMomDir)
 
void SetRestMass (G4double Mass_c2)
 
void SetCurveLength (G4double nCurve_s)
 
void SetKineticEnergy (G4double nEnergy)
 
void SetLabTimeOfFlight (G4double tofLab)
 
void SetProperTimeOfFlight (G4double tofProper)
 
void DumpToArray (G4double valArr[ncompSVEC]) const
 
void LoadFromArray (const G4double valArr[ncompSVEC], G4int noVarsIntegrated)
 
void InitialiseSpin (const G4ThreeVector &vecPolarization)
 
G4ThreeVector GetSpin () const
 
void SetSpin (G4ThreeVector vSpin)
 

Friends

std::ostream & operator<< (std::ostream &os, const G4FieldTrack &SixVec)
 

Detailed Description

Definition at line 44 of file G4FieldTrack.hh.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
ncompSVEC 

Definition at line 144 of file G4FieldTrack.hh.

144{ ncompSVEC = 12 };

Constructor & Destructor Documentation

◆ G4FieldTrack() [1/5]

G4FieldTrack::G4FieldTrack ( const G4ThreeVector pPosition,
G4double  LaboratoryTimeOfFlight,
const G4ThreeVector pMomentumDirection,
G4double  kineticEnergy,
G4double  restMass_c2,
G4double  charge,
const G4ThreeVector polarization,
G4double  magnetic_dipole_moment = 0.0,
G4double  curve_length = 0.0,
G4double  PDGspin = -1.0 
)

Definition at line 81 of file G4FieldTrack.cc.

91: fDistanceAlongCurve(curve_length),
92 fKineticEnergy(kineticEnergy),
93 fRestMass_c2(restMass_c2),
94 fLabTimeOfFlight(LaboratoryTimeOfFlight),
95 fProperTimeOfFlight(0.),
96 // fMomentumDir(pMomentumDirection),
97 fChargeState( charge, magnetic_dipole_moment, pdgSpin )
98 // fChargeState( charge, magnetic_dipole_moment ) ,
99 // fPDGSpin( pdgSpin )
100{
101 UpdateFourMomentum( kineticEnergy, pMomentumDirection );
102 // Sets momentum direction as well.
103
104 SetPosition( pPosition );
105 SetPolarization( vecPolarization );
106}
void UpdateFourMomentum(G4double kineticEnergy, const G4ThreeVector &momentumDirection)
void SetPolarization(const G4ThreeVector &vecPol)
void SetPosition(G4ThreeVector nPos)

◆ G4FieldTrack() [2/5]

G4FieldTrack::G4FieldTrack ( char  )

Definition at line 136 of file G4FieldTrack.cc.

137 : fKineticEnergy(0.), fRestMass_c2(0.), fLabTimeOfFlight(0.),
138 fProperTimeOfFlight(0.), fChargeState( DBL_MAX , DBL_MAX, -1 )
139{
140 G4ThreeVector Zero(0.0, 0.0, 0.0);
141 SetCurvePnt( Zero, Zero, 0.0 );
142 SetPolarization( Zero );
143 // fInitialMomentumMag = 0.00; // Invalid
144 // fLastMomentumMag = 0.0;
145}
#define DBL_MAX
Definition: templates.hh:62

◆ G4FieldTrack() [3/5]

G4FieldTrack::G4FieldTrack ( const G4ThreeVector pPosition,
const G4ThreeVector pMomentumDirection,
G4double  curve_length,
G4double  kineticEnergy,
const G4double  restMass_c2,
G4double  velocity,
G4double  LaboratoryTimeOfFlight = 0.0,
G4double  ProperTimeOfFlight = 0.0,
const G4ThreeVector pPolarization = nullptr,
G4double  PDGspin = -1.0 
)

Definition at line 108 of file G4FieldTrack.cc.

118 : fDistanceAlongCurve(curve_length),
119 fKineticEnergy(kineticEnergy),
120 fRestMass_c2(restMass_c2),
121 fLabTimeOfFlight(pLaboratoryTimeOfFlight),
122 fProperTimeOfFlight(pProperTimeOfFlight),
123 fChargeState( DBL_MAX, DBL_MAX, -1.0 ) // charge not set
124{
125 UpdateFourMomentum( kineticEnergy, pMomentumDirection );
126 // Sets momentum direction as well.
127
128 SetPosition( pPosition );
129 fChargeState.SetPDGSpin( pdgSpin );
130
131 G4ThreeVector PolarVec(0.0, 0.0, 0.0);
132 if( pPolarization ) { PolarVec= *pPolarization; }
133 SetPolarization( PolarVec );
134}
void SetPDGSpin(G4double spin)

◆ ~G4FieldTrack()

G4FieldTrack::~G4FieldTrack ( )

◆ G4FieldTrack() [4/5]

G4FieldTrack::G4FieldTrack ( const G4FieldTrack pFieldTrack)
inline

◆ G4FieldTrack() [5/5]

G4FieldTrack::G4FieldTrack ( G4FieldTrack &&  from)
inline

Member Function Documentation

◆ DumpToArray()

◆ GetCharge()

G4double G4FieldTrack::GetCharge ( ) const
inline

◆ GetChargeState()

const G4ChargeState * G4FieldTrack::GetChargeState ( ) const
inline

◆ GetCurveLength()

◆ GetKineticEnergy()

◆ GetLabTimeOfFlight()

◆ GetMomentum()

◆ GetMomentumDir()

◆ GetMomentumDirection()

◆ GetPDGSpin()

G4double G4FieldTrack::GetPDGSpin ( )
inline

◆ GetPolarization()

G4ThreeVector G4FieldTrack::GetPolarization ( ) const
inline

◆ GetPosition()

G4ThreeVector G4FieldTrack::GetPosition ( ) const
inline

Referenced by G4ImportanceProcess::AlongStepGetPhysicalInteractionLength(), G4WeightCutOffProcess::AlongStepGetPhysicalInteractionLength(), G4WeightWindowProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldScoringProcess::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4ParallelGeometriesLimiterProcess::AlongStepGetPhysicalInteractionLength(), G4FastSimulationManagerProcess::AlongStepGetPhysicalInteractionLength(), G4ChordFinder::ApproxCurvePointS(), G4ChordFinder::ApproxCurvePointV(), G4VIntersectionLocator::CheckAndReEstimateEndpoint(), G4PathFinder::ComputeStep(), G4ITPathFinder::ComputeStep(), G4PropagatorInField::ComputeStep(), G4PathFinder::DoNextCurvedStep(), G4ITPathFinder::DoNextCurvedStep(), G4PathFinder::DoNextLinearStep(), G4ITPathFinder::DoNextLinearStep(), G4BrentLocator::EstimateIntersectionPoint(), G4MultiLevelLocator::EstimateIntersectionPoint(), G4SimpleLocator::EstimateIntersectionPoint(), G4PathFinder::Locate(), G4DriverReporter::PrintStat_Aux(), G4MagInt_Driver::PrintStat_Aux(), G4OldMagIntDriver::PrintStat_Aux(), G4PropagatorInField::printStatus(), G4VIntersectionLocator::printStatus(), G4MagInt_Driver::PrintStatus(), G4OldMagIntDriver::PrintStatus(), G4DriverReporter::PrintStatus(), G4MagInt_Driver::QuickAdvance(), G4OldMagIntDriver::QuickAdvance(), G4VIntersectionLocator::ReEstimateEndpoint(), and G4PathFinder::ReLocate().

◆ GetProperTimeOfFlight()

G4double G4FieldTrack::GetProperTimeOfFlight ( ) const
inline

◆ GetRestMass()

G4double G4FieldTrack::GetRestMass ( ) const
inline

◆ GetSpin()

◆ InitialiseSpin()

void G4FieldTrack::InitialiseSpin ( const G4ThreeVector vecPolarization)
inline

◆ LoadFromArray()

void G4FieldTrack::LoadFromArray ( const G4double  valArr[ncompSVEC],
G4int  noVarsIntegrated 
)

Definition at line 175 of file G4FieldTrack.cc.

177{
178 // Fill the variables not integrated with zero -- so it's clear !!
179 //
180 G4double valArr[ncompSVEC];
181 for(G4int i=0; i<noVarsIntegrated; ++i)
182 {
183 valArr[i] = valArrIn[i];
184 }
185 for(G4int i=noVarsIntegrated; i<ncompSVEC; ++i)
186 {
187 valArr[i] = 0.0;
188 }
189
190 SixVector[0] = valArr[0];
191 SixVector[1] = valArr[1];
192 SixVector[2] = valArr[2];
193 SixVector[3] = valArr[3];
194 SixVector[4] = valArr[4];
195 SixVector[5] = valArr[5];
196
197 G4ThreeVector Momentum(valArr[3],valArr[4],valArr[5]);
198
199 G4double momentum_square= Momentum.mag2();
200 fMomentumDir= Momentum.unit();
201
202 fKineticEnergy = momentum_square
203 / (std::sqrt(momentum_square+fRestMass_c2*fRestMass_c2)
204 + fRestMass_c2 );
205 // The above equation is stable for small and large momenta
206
207 // The following components may or may not be
208 // integrated over -- integration is optional
209 // fKineticEnergy = valArr[6];
210
211 fLabTimeOfFlight = valArr[7];
212 fProperTimeOfFlight = valArr[8];
213 G4ThreeVector vecPolarization= G4ThreeVector(valArr[9],valArr[10],valArr[11]);
214 SetPolarization( vecPolarization );
215
216 // fMomentumDir=G4ThreeVector(valArr[13],valArr[14],valArr[15]);
217 // fDistanceAlongCurve= valArr[];
218}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85

Referenced by G4MagInt_Driver::AccurateAdvance(), G4OldMagIntDriver::AccurateAdvance(), G4MagInt_Driver::PrintStatus(), G4OldMagIntDriver::PrintStatus(), G4DriverReporter::PrintStatus(), G4MagInt_Driver::QuickAdvance(), and G4OldMagIntDriver::QuickAdvance().

◆ operator=() [1/2]

G4FieldTrack & G4FieldTrack::operator= ( const G4FieldTrack rStVec)
inline

◆ operator=() [2/2]

G4FieldTrack & G4FieldTrack::operator= ( G4FieldTrack &&  from)
inline

◆ SetChargeAndMoments()

void G4FieldTrack::SetChargeAndMoments ( G4double  charge,
G4double  magnetic_dipole_moment = DBL_MAX,
G4double  electric_dipole_moment = DBL_MAX,
G4double  magnetic_charge = DBL_MAX 
)

Definition at line 147 of file G4FieldTrack.cc.

152{
153 fChargeState.SetChargesAndMoments( charge,
154 magnetic_dipole_moment,
155 electric_dipole_moment,
156 magnetic_charge );
157
158 // NOTE: Leaves Spin unchanged !
159 //
160 // G4double pdgSpin= fChargeState.GetSpin();
161 // New Property of ChargeState (not well documented! )
162
163 // IDEA: Improve the implementation using handles
164 // -- and handle to the old one (which can be shared by other copies) and
165 // must not be left to hang loose
166 //
167 // fpChargeState= new G4ChargeState( charge, magnetic_dipole_moment,
168 // electric_dipole_moment, magnetic_charge );
169}
void SetChargesAndMoments(G4double charge, G4double magnetic_dipole_moment, G4double electric_dipole_moment, G4double magnetic_charge)

Referenced by G4FieldTrackUpdator::Update().

◆ SetCurveLength()

◆ SetKineticEnergy()

void G4FieldTrack::SetKineticEnergy ( G4double  nEnergy)
inline

◆ SetLabTimeOfFlight()

void G4FieldTrack::SetLabTimeOfFlight ( G4double  tofLab)
inline

◆ SetMomentum()

void G4FieldTrack::SetMomentum ( G4ThreeVector  nMomDir)
inline

◆ SetMomentumDir()

void G4FieldTrack::SetMomentumDir ( G4ThreeVector  nMomDir)
inline

◆ SetPDGSpin()

void G4FieldTrack::SetPDGSpin ( G4double  pdgSpin)
inline

◆ SetPolarization()

void G4FieldTrack::SetPolarization ( const G4ThreeVector vecPol)
inline

Referenced by G4FieldTrack(), and LoadFromArray().

◆ SetPosition()

◆ SetProperTimeOfFlight()

void G4FieldTrack::SetProperTimeOfFlight ( G4double  tofProper)
inline

◆ SetRestMass()

void G4FieldTrack::SetRestMass ( G4double  Mass_c2)
inline

◆ SetSpin()

void G4FieldTrack::SetSpin ( G4ThreeVector  vSpin)
inline

◆ UpdateFourMomentum()

void G4FieldTrack::UpdateFourMomentum ( G4double  kineticEnergy,
const G4ThreeVector momentumDirection 
)
inline

Referenced by G4FieldTrack().

◆ UpdateState()

void G4FieldTrack::UpdateState ( const G4ThreeVector pPosition,
G4double  LaboratoryTimeOfFlight,
const G4ThreeVector pMomentumDirection,
G4double  kineticEnergy 
)
inline

Friends And Related Function Documentation

◆ operator<<

std::ostream & operator<< ( std::ostream &  os,
const G4FieldTrack SixVec 
)
friend

Definition at line 33 of file G4FieldTrack.cc.

34{
35 const G4double* SixV = SixVec.SixVector;
36 const G4int precPos= 9; // For position
37 const G4int precEp= 9; // For Energy / momentum
38 const G4int precLen= 12; // For Length along track
39 const G4int precSpin= 9; // For polarisation
40 const G4int precTime= 6; // For time of flight
41 const G4int oldpr= os.precision(precPos);
42 os << " ( ";
43 os << " X= " << SixV[0] << " " << SixV[1] << " "
44 << SixV[2] << " "; // Position
45 os.precision(precEp);
46 os << " P= " << SixV[3] << " " << SixV[4] << " "
47 << SixV[5] << " "; // Momentum
48 os << " Pmag= "
49 << G4ThreeVector(SixV[3], SixV[4], SixV[5]).mag(); // mom magnitude
50 os << " Ekin= " << SixVec.fKineticEnergy ;
51 os.precision(precLen);
52 os << " l= " << SixVec.GetCurveLength();
53 os.precision(6);
54 os << " m0= " << SixVec.fRestMass_c2;
55 os << " (Pdir-1)= " << SixVec.fMomentumDir.mag()-1.0;
56 if( SixVec.fLabTimeOfFlight > 0.0 )
57 {
58 os.precision(precTime);
59 }
60 else
61 {
62 os.precision(3);
63 }
64 os << " t_lab= " << SixVec.fLabTimeOfFlight;
65 os << " t_proper= " << SixVec.fProperTimeOfFlight ;
66 G4ThreeVector pol= SixVec.GetPolarization();
67 if( pol.mag2() > 0.0 )
68 {
69 os.precision(precSpin);
70 os << " PolV= " << pol; // SixVec.GetPolarization();
71 }
72 else
73 {
74 os << " PolV= (0,0,0) ";
75 }
76 os << " ) ";
77 os.precision(oldpr);
78 return os;
79}
double mag2() const
double mag() const
G4double GetCurveLength() const
G4ThreeVector GetPolarization() const

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