Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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//
27// $Id$
28//
29//
30// This is the standard right-hand side for equation of motion.
31//
32// 30.08.2007 Chris Gong, Peter Gumplinger
33// 14.02.2009 Kevin Lynch
34// 06.11.2009 Hiromi Iinuma
35//
36// -------------------------------------------------------------------
37
40#include "G4ThreeVector.hh"
41#include "globals.hh"
43#include "G4SystemOfUnits.hh"
44
46 : G4EquationOfMotion( emField ), fElectroMagCof(0.), fMassCof(0.),
47 omegac(0.), anomaly(0.0011659208), pcharge(0.), E(0.), gamma(0.), beta(0.)
48{
49}
50
52{
53}
54
55void
57 G4double MomentumXc,
58 G4double particleMass)
59{
60 fElectroMagCof = eplus*particleCharge*c_light ;
61 fMassCof = particleMass*particleMass ;
62
63 omegac = (eplus/particleMass)*c_light;
64
65 pcharge = particleCharge;
66
67 E = std::sqrt(sqr(MomentumXc)+sqr(particleMass));
68 beta = MomentumXc/E;
69 gamma = E/particleMass;
70
71}
72
73void
75 const G4double Field[],
76 G4double dydx[] ) const
77{
78
79 // Components of y:
80 // 0-2 dr/ds,
81 // 3-5 dp/ds - momentum derivatives
82 // 9-11 dSpin/ds = (1/beta) dSpin/dt - spin derivatives
83
84 // The BMT equation, following J.D.Jackson, Classical
85 // Electrodynamics, Second Edition,
86 // dS/dt = (e/mc) S \cross
87 // [ (g/2-1 +1/\gamma) B
88 // -(g/2-1)\gamma/(\gamma+1) (\beta \cdot B)\beta
89 // -(g/2-\gamma/(\gamma+1) \beta \cross E ]
90 // where
91 // S = \vec{s}, where S^2 = 1
92 // B = \vec{B}
93 // \beta = \vec{\beta} = \beta \vec{u} with u^2 = 1
94 // E = \vec{E}
95
96 G4double pSquared = y[3]*y[3] + y[4]*y[4] + y[5]*y[5] ;
97
98 G4double Energy = std::sqrt( pSquared + fMassCof );
99 G4double cof2 = Energy/c_light ;
100
101 G4double pModuleInverse = 1.0/std::sqrt(pSquared) ;
102
103 G4double inverse_velocity = Energy * pModuleInverse / c_light;
104
105 G4double cof1 = fElectroMagCof*pModuleInverse ;
106
107 dydx[0] = y[3]*pModuleInverse ;
108 dydx[1] = y[4]*pModuleInverse ;
109 dydx[2] = y[5]*pModuleInverse ;
110
111 dydx[3] = cof1*(cof2*Field[3] + (y[4]*Field[2] - y[5]*Field[1])) ;
112
113 dydx[4] = cof1*(cof2*Field[4] + (y[5]*Field[0] - y[3]*Field[2])) ;
114
115 dydx[5] = cof1*(cof2*Field[5] + (y[3]*Field[1] - y[4]*Field[0])) ;
116
117 dydx[6] = dydx[8] = 0.;//not used
118
119 // Lab Time of flight
120 dydx[7] = inverse_velocity;
121
122 G4ThreeVector BField(Field[0],Field[1],Field[2]);
123 G4ThreeVector EField(Field[3],Field[4],Field[5]);
124
125 EField /= c_light;
126
127 G4ThreeVector u(y[3], y[4], y[5]);
128 u *= pModuleInverse;
129
130 G4double udb = anomaly*beta*gamma/(1.+gamma) * (BField * u);
131 G4double ucb = (anomaly+1./gamma)/beta;
132 G4double uce = anomaly + 1./(gamma+1.);
133
134 G4ThreeVector Spin(y[9],y[10],y[11]);
135
136 G4ThreeVector dSpin
137 = pcharge*omegac*( ucb*(Spin.cross(BField))-udb*(Spin.cross(u))
138 // from Jackson
139 // -uce*Spin.cross(u.cross(EField)) );
140 // but this form has one less operation
141 - uce*(u*(Spin*EField) - EField*(Spin*u)) );
142
143 dydx[ 9] = dSpin.x();
144 dydx[10] = dSpin.y();
145 dydx[11] = dSpin.z();
146
147 return ;
148}
double G4double
Definition: G4Types.hh:64
double z() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
G4EqEMFieldWithSpin(G4ElectroMagneticField *emField)
void EvaluateRhsGivenB(const G4double y[], const G4double Field[], G4double dydx[]) const
void SetChargeMomentumMass(G4double particleCharge, G4double MomentumXc, G4double mass)
T sqr(const T &x)
Definition: templates.hh:145