Geant4-11
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
39 : G4EquationOfMotion( emField ), charge(0.), mass(0.), magMoment(0.),
40 spin(0.), fElectroMagCof(0.), fMassCof(0.), omegac(0.),
41 anomaly(0.0011659208), beta(0.), gamma(0.)
42{
43}
44
46{
47}
48
49void
51 G4double MomentumXc,
52 G4double particleMass)
53{
54 charge = particleCharge.GetCharge();
55 mass = particleMass;
56 magMoment = particleCharge.GetMagneticDipoleMoment();
57 spin = particleCharge.GetSpin();
58
61
63
65
66 G4double g_BMT;
67 if ( spin != 0. ) g_BMT = (std::abs(magMoment)/muB)/spin;
68 else g_BMT = 2.;
69
70 anomaly = (g_BMT - 2.)/2.;
71
72 G4double E = std::sqrt(sqr(MomentumXc)+sqr(mass));
73 beta = MomentumXc/E;
74 gamma = E/mass;
75}
76
77void
79 const G4double Field[],
80 G4double dydx[] ) const
81{
82
83 // Components of y:
84 // 0-2 dr/ds,
85 // 3-5 dp/ds - momentum derivatives
86 // 9-11 dSpin/ds = (1/beta) dSpin/dt - spin derivatives
87
88 // The BMT equation, following J.D.Jackson, Classical
89 // Electrodynamics, Second Edition,
90 // dS/dt = (e/mc) S \cross
91 // [ (g/2-1 +1/\gamma) B
92 // -(g/2-1)\gamma/(\gamma+1) (\beta \cdot B)\beta
93 // -(g/2-\gamma/(\gamma+1) \beta \cross E ]
94 // where
95 // S = \vec{s}, where S^2 = 1
96 // B = \vec{B}
97 // \beta = \vec{\beta} = \beta \vec{u} with u^2 = 1
98 // E = \vec{E}
99
100 G4double pSquared = y[3]*y[3] + y[4]*y[4] + y[5]*y[5] ;
101
102 G4double Energy = std::sqrt( pSquared + fMassCof );
103 G4double cof2 = Energy/c_light ;
104
105 G4double pModuleInverse = 1.0/std::sqrt(pSquared) ;
106
107 G4double inverse_velocity = Energy * pModuleInverse / c_light;
108
109 G4double cof1 = fElectroMagCof*pModuleInverse ;
110
111 dydx[0] = y[3]*pModuleInverse ;
112 dydx[1] = y[4]*pModuleInverse ;
113 dydx[2] = y[5]*pModuleInverse ;
114
115 dydx[3] = cof1*(cof2*Field[3] + (y[4]*Field[2] - y[5]*Field[1])) ;
116
117 dydx[4] = cof1*(cof2*Field[4] + (y[5]*Field[0] - y[3]*Field[2])) ;
118
119 dydx[5] = cof1*(cof2*Field[5] + (y[3]*Field[1] - y[4]*Field[0])) ;
120
121 dydx[6] = dydx[8] = 0.;//not used
122
123 // Lab Time of flight
124 dydx[7] = inverse_velocity;
125
126 G4ThreeVector BField(Field[0],Field[1],Field[2]);
127 G4ThreeVector EField(Field[3],Field[4],Field[5]);
128
129 EField /= c_light;
130
131 G4ThreeVector u(y[3], y[4], y[5]);
132 u *= pModuleInverse;
133
134 G4double udb = anomaly*beta*gamma/(1.+gamma) * (BField * u);
135 G4double ucb = (anomaly+1./gamma)/beta;
136 G4double uce = anomaly + 1./(gamma+1.);
137
138 G4ThreeVector Spin(y[9],y[10],y[11]);
139
140 G4double pcharge;
141 if (charge == 0.)
142 {
143 pcharge = 1.;
144 }
145 else
146 {
147 pcharge = charge;
148 }
149
150 G4ThreeVector dSpin(0.,0.,0.);
151 if (Spin.mag2() != 0.)
152 {
153 dSpin = pcharge*omegac*( ucb*(Spin.cross(BField))-udb*(Spin.cross(u))
154 // from Jackson
155 // -uce*Spin.cross(u.cross(EField)) );
156 // but this form has one less operation
157 - uce*(u*(Spin*EField) - EField*(Spin*u)) );
158 }
159
160 dydx[ 9] = dSpin.x();
161 dydx[10] = dSpin.y();
162 dydx[11] = dSpin.z();
163
164 return;
165}
static constexpr double eplus
Definition: G4SIunits.hh:184
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)
void EvaluateRhsGivenB(const G4double y[], const G4double Field[], G4double dydx[]) const
void SetChargeMomentumMass(G4ChargeState particleCharge, G4double MomentumXc, G4double mass)
float c_light
Definition: hepunit.py:256
float hbar_Planck
Definition: hepunit.py:263
float c_squared
Definition: hepunit.py:257
T sqr(const T &x)
Definition: templates.hh:128