Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EqEMFieldWithSpin.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// G4EqEMFieldWithSpin implementation
27//
28// Created: Chris Gong & Peter Gumplinger, 30.08.2007
29// -------------------------------------------------------------------
30
33#include "G4ThreeVector.hh"
34#include "globals.hh"
36#include "G4SystemOfUnits.hh"
37
42
44
45void
47 G4double MomentumXc,
48 G4double particleMass)
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}
78
79void
81 const G4double Field[],
82 G4double dydx[] ) const
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}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
double z() const
double x() const
double y() const
G4double GetCharge() const
G4double GetMagneticDipoleMoment() const
G4double GetSpin() const
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
G4EquationOfMotion(G4Field *Field)
T sqr(const T &x)
Definition templates.hh:128