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

#include <G4EqEMFieldWithSpin.hh>

+ Inheritance diagram for G4EqEMFieldWithSpin:

Public Member Functions

 G4EqEMFieldWithSpin (G4ElectroMagneticField *emField)
 
 ~G4EqEMFieldWithSpin () override
 
void SetChargeMomentumMass (G4ChargeState particleCharge, G4double MomentumXc, G4double mass) override
 
void EvaluateRhsGivenB (const G4double y[], const G4double Field[], G4double dydx[]) const override
 
void SetAnomaly (G4double a)
 
G4double GetAnomaly () const
 
- Public Member Functions inherited from G4EquationOfMotion
 G4EquationOfMotion (G4Field *Field)
 
virtual ~G4EquationOfMotion ()
 
virtual void EvaluateRhsGivenB (const G4double y[], const G4double B[3], G4double dydx[]) const =0
 
void RightHandSide (const G4double y[], G4double dydx[]) const
 
void EvaluateRhsReturnB (const G4double y[], G4double dydx[], G4double Field[]) const
 
void GetFieldValue (const G4double Point[4], G4double Field[]) const
 
const G4FieldGetFieldObj () const
 
G4FieldGetFieldObj ()
 
void SetFieldObj (G4Field *pField)
 

Detailed Description

Definition at line 43 of file G4EqEMFieldWithSpin.hh.

Constructor & Destructor Documentation

◆ G4EqEMFieldWithSpin()

G4EqEMFieldWithSpin::G4EqEMFieldWithSpin ( G4ElectroMagneticField * emField)

Definition at line 38 of file G4EqEMFieldWithSpin.cc.

39 : G4EquationOfMotion( emField )
40{
41}
G4EquationOfMotion(G4Field *Field)

◆ ~G4EqEMFieldWithSpin()

G4EqEMFieldWithSpin::~G4EqEMFieldWithSpin ( )
overridedefault

Member Function Documentation

◆ EvaluateRhsGivenB()

void G4EqEMFieldWithSpin::EvaluateRhsGivenB ( const G4double y[],
const G4double Field[],
G4double dydx[] ) const
override

Definition at line 80 of file G4EqEMFieldWithSpin.cc.

83{
84
85 // Components of y:
86 // 0-2 dr/ds,
87 // 3-5 dp/ds - momentum derivatives
88 // 9-11 dSpin/ds = (1/beta) dSpin/dt - spin derivatives
89
90 // The BMT equation, following J.D.Jackson, Classical
91 // Electrodynamics, Second Edition,
92 // dS/dt = (e/mc) S \cross
93 // [ (g/2-1 +1/\gamma) B
94 // -(g/2-1)\gamma/(\gamma+1) (\beta \cdot B)\beta
95 // -(g/2-\gamma/(\gamma+1) \beta \cross E ]
96 // where
97 // S = \vec{s}, where S^2 = 1
98 // B = \vec{B}
99 // \beta = \vec{\beta} = \beta \vec{u} with u^2 = 1
100 // E = \vec{E}
101
102 G4double pSquared = y[3]*y[3] + y[4]*y[4] + y[5]*y[5] ;
103
104 G4double Energy = std::sqrt( pSquared + fMassCof );
105 G4double cof2 = Energy/c_light ;
106
107 G4double pModuleInverse = 1.0/std::sqrt(pSquared) ;
108
109 G4double inverse_velocity = Energy * pModuleInverse / c_light;
110
111 G4double cof1 = fElectroMagCof*pModuleInverse ;
112
113 dydx[0] = y[3]*pModuleInverse ;
114 dydx[1] = y[4]*pModuleInverse ;
115 dydx[2] = y[5]*pModuleInverse ;
116
117 dydx[3] = cof1*(cof2*Field[3] + (y[4]*Field[2] - y[5]*Field[1])) ;
118
119 dydx[4] = cof1*(cof2*Field[4] + (y[5]*Field[0] - y[3]*Field[2])) ;
120
121 dydx[5] = cof1*(cof2*Field[5] + (y[3]*Field[1] - y[4]*Field[0])) ;
122
123 dydx[6] = dydx[8] = 0.;//not used
124
125 // Lab Time of flight
126 dydx[7] = inverse_velocity;
127
128 G4ThreeVector BField(Field[0],Field[1],Field[2]);
129 G4ThreeVector EField(Field[3],Field[4],Field[5]);
130
131 EField /= c_light;
132
133 G4ThreeVector u(y[3], y[4], y[5]);
134 u *= pModuleInverse;
135
136 G4double udb = anomaly*beta*gamma/(1.+gamma) * (BField * u);
137 G4double ucb = (anomaly+1./gamma)/beta;
138 G4double uce = anomaly + 1./(gamma+1.);
139
140 G4ThreeVector Spin(y[9],y[10],y[11]);
141
142 G4double pcharge;
143 if (charge == 0.)
144 {
145 pcharge = 1.;
146 }
147 else
148 {
149 pcharge = charge;
150 }
151
152 G4ThreeVector dSpin(0.,0.,0.);
153 if (Spin.mag2() != 0.)
154 {
155 dSpin = pcharge*omegac*( ucb*(Spin.cross(BField))-udb*(Spin.cross(u))
156 // from Jackson
157 // -uce*Spin.cross(u.cross(EField)) );
158 // but this form has one less operation
159 - uce*(u*(Spin*EField) - EField*(Spin*u)) );
160 }
161
162 dydx[ 9] = dSpin.x();
163 dydx[10] = dSpin.y();
164 dydx[11] = dSpin.z();
165
166 return;
167}
double G4double
Definition G4Types.hh:83

◆ GetAnomaly()

G4double G4EqEMFieldWithSpin::GetAnomaly ( ) const
inline

Definition at line 61 of file G4EqEMFieldWithSpin.hh.

61{ return anomaly; }

◆ SetAnomaly()

void G4EqEMFieldWithSpin::SetAnomaly ( G4double a)
inline

Definition at line 60 of file G4EqEMFieldWithSpin.hh.

60{ anomaly = a; }

◆ SetChargeMomentumMass()

void G4EqEMFieldWithSpin::SetChargeMomentumMass ( G4ChargeState particleCharge,
G4double MomentumXc,
G4double mass )
overridevirtual

Implements G4EquationOfMotion.

Definition at line 46 of file G4EqEMFieldWithSpin.cc.

49{
50 charge = particleCharge.GetCharge();
51 mass = particleMass;
52 magMoment = particleCharge.GetMagneticDipoleMoment();
53 spin = particleCharge.GetSpin();
54
55 fElectroMagCof = eplus*charge*c_light ;
56 fMassCof = mass*mass;
57
58 omegac = (eplus/mass)*c_light;
59
60 G4double muB = 0.5*eplus*hbar_Planck/(mass/c_squared);
61
62 G4double g_BMT;
63 if ( spin != 0. )
64 {
65 g_BMT = (std::abs(magMoment)/muB)/spin;
66 }
67 else
68 {
69 g_BMT = 2.;
70 }
71
72 anomaly = (g_BMT - 2.)/2.;
73
74 G4double E = std::sqrt(sqr(MomentumXc)+sqr(mass));
75 beta = MomentumXc/E;
76 gamma = E/mass;
77}
G4double GetCharge() const
G4double GetMagneticDipoleMoment() const
G4double GetSpin() const
T sqr(const T &x)
Definition templates.hh:128

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