Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EqEMFieldWithEDM.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// 19.02.2009 Kevin Lynch, based on G4EqEMFieldWithSpin
33// 06.11.2009 Hiromi Iinuma see:
34// http://hypernews.slac.stanford.edu/HyperNews/geant4/get/emfields/161.html
35//
36// -------------------------------------------------------------------
37
38#include "G4EqEMFieldWithEDM.hh"
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), eta(0.), pcharge(0.), E(0.),
48 gamma(0.), beta(0.)
49{
50}
51
53{
54}
55
56void
58 G4double MomentumXc,
59 G4double particleMass)
60{
61 fElectroMagCof = eplus*particleCharge*c_light ;
62 fMassCof = particleMass*particleMass ;
63
64 omegac = (eplus/particleMass)*c_light;
65
66 pcharge = particleCharge;
67
68 E = std::sqrt(sqr(MomentumXc)+sqr(particleMass));
69 beta = MomentumXc/E;
70 gamma = E/particleMass;
71
72}
73
74void
76 const G4double Field[],
77 G4double dydx[] ) const
78{
79
80 // Components of y:
81 // 0-2 dr/ds,
82 // 3-5 dp/ds - momentum derivatives
83 // 9-11 dSpin/ds = (1/beta) dSpin/dt - spin derivatives
84
85 // The BMT equation, following J.D.Jackson, Classical
86 // Electrodynamics, Second Edition, with additions for EDM
87 // evolution from
88 // M.Nowakowski, et.al. Eur.J.Phys.26, pp 545-560, (2005)
89 // or
90 // Silenko, Phys.Rev.ST Accel.Beams 9:034003, (2006)
91
92 // dS/dt = (e/m) S \cross
93 // MDM: [ (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 //
97 // EDM: eta/2( E - gamma/(gamma+1) \beta (\beta \cdot E)
98 // + \beta \cross B ) ]
99 //
100 // where
101 // S = \vec{s}, where S^2 = 1
102 // B = \vec{B}
103 // \beta = \vec{\beta} = \beta \vec{u} with u^2 = 1
104 // E = \vec{E}
105
106 G4double pSquared = y[3]*y[3] + y[4]*y[4] + y[5]*y[5] ;
107
108 G4double Energy = std::sqrt( pSquared + fMassCof );
109 G4double cof2 = Energy/c_light ;
110
111 G4double pModuleInverse = 1.0/std::sqrt(pSquared) ;
112
113 G4double inverse_velocity = Energy * pModuleInverse / c_light;
114
115 G4double cof1 = fElectroMagCof*pModuleInverse ;
116
117 dydx[0] = y[3]*pModuleInverse ;
118 dydx[1] = y[4]*pModuleInverse ;
119 dydx[2] = y[5]*pModuleInverse ;
120
121 dydx[3] = cof1*(cof2*Field[3] + (y[4]*Field[2] - y[5]*Field[1])) ;
122
123 dydx[4] = cof1*(cof2*Field[4] + (y[5]*Field[0] - y[3]*Field[2])) ;
124
125 dydx[5] = cof1*(cof2*Field[5] + (y[3]*Field[1] - y[4]*Field[0])) ;
126
127 dydx[6] = dydx[8] = 0.;//not used
128
129 // Lab Time of flight
130 dydx[7] = inverse_velocity;
131
132 G4ThreeVector BField(Field[0],Field[1],Field[2]);
133 G4ThreeVector EField(Field[3],Field[4],Field[5]);
134
135 EField /= c_light;
136
137 G4ThreeVector u(y[3], y[4], y[5]);
138 u *= pModuleInverse;
139
140 G4double udb = anomaly*beta*gamma/(1.+gamma) * (BField * u);
141 G4double ucb = (anomaly+1./gamma)/beta;
142 G4double uce = anomaly + 1./(gamma+1.);
143 G4double ude = beta*gamma/(1.+gamma)*(EField*u);
144
145 G4ThreeVector Spin(y[9],y[10],y[11]);
146
147 G4ThreeVector dSpin
148 = pcharge*omegac*( ucb*(Spin.cross(BField))-udb*(Spin.cross(u))
149 // from Jackson
150 // -uce*Spin.cross(u.cross(EField)) )
151 // but this form has one less operation
152 - uce*(u*(Spin*EField) - EField*(Spin*u))
153 + eta/2.*(Spin.cross(EField) - ude*(Spin.cross(u))
154 // +Spin.cross(u.cross(Bfield))
155 + (u*(Spin*BField) - BField*(Spin*u)) ) );
156
157 dydx[ 9] = dSpin.x();
158 dydx[10] = dSpin.y();
159 dydx[11] = dSpin.z();
160
161 return ;
162}
double G4double
Definition: G4Types.hh:64
double z() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
void SetChargeMomentumMass(G4double particleCharge, G4double MomentumXc, G4double mass)
void EvaluateRhsGivenB(const G4double y[], const G4double Field[], G4double dydx[]) const
G4EqEMFieldWithEDM(G4ElectroMagneticField *emField)
T sqr(const T &x)
Definition: templates.hh:145