Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4EqEMFieldWithSpin.cc
이 파일의 문서화 페이지로 가기
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: G4EqEMFieldWithSpin.cc 95822 2016-02-26 08:04:51Z gcosmo $
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 
38 #include "G4EqEMFieldWithSpin.hh"
40 #include "G4ThreeVector.hh"
41 #include "globals.hh"
42 #include "G4PhysicalConstants.hh"
43 #include "G4SystemOfUnits.hh"
44 
46  : G4EquationOfMotion( emField ), charge(0.), mass(0.), magMoment(0.),
47  spin(0.), fElectroMagCof(0.), fMassCof(0.), omegac(0.),
48  anomaly(0.0011659208), beta(0.), gamma(0.)
49 {
50 }
51 
53 {
54 }
55 
56 void
58  G4double MomentumXc,
59  G4double particleMass)
60 {
61  charge = particleCharge.GetCharge();
62  mass = particleMass;
63  magMoment = particleCharge.GetMagneticDipoleMoment();
64  spin = particleCharge.GetSpin();
65 
67  fMassCof = mass*mass;
68 
69  omegac = (eplus/mass)*c_light;
70 
71  G4double muB = 0.5*eplus*hbar_Planck/(mass/c_squared);
72 
73  G4double g_BMT;
74  if ( spin != 0. ) g_BMT = (std::abs(magMoment)/muB)/spin;
75  else g_BMT = 2.;
76 
77  anomaly = (g_BMT - 2.)/2.;
78 
79  G4double E = std::sqrt(sqr(MomentumXc)+sqr(mass));
80  beta = MomentumXc/E;
81  gamma = E/mass;
82 }
83 
84 void
86  const G4double Field[],
87  G4double dydx[] ) const
88 {
89 
90  // Components of y:
91  // 0-2 dr/ds,
92  // 3-5 dp/ds - momentum derivatives
93  // 9-11 dSpin/ds = (1/beta) dSpin/dt - spin derivatives
94 
95  // The BMT equation, following J.D.Jackson, Classical
96  // Electrodynamics, Second Edition,
97  // dS/dt = (e/mc) S \cross
98  // [ (g/2-1 +1/\gamma) B
99  // -(g/2-1)\gamma/(\gamma+1) (\beta \cdot B)\beta
100  // -(g/2-\gamma/(\gamma+1) \beta \cross E ]
101  // where
102  // S = \vec{s}, where S^2 = 1
103  // B = \vec{B}
104  // \beta = \vec{\beta} = \beta \vec{u} with u^2 = 1
105  // E = \vec{E}
106 
107  G4double pSquared = y[3]*y[3] + y[4]*y[4] + y[5]*y[5] ;
108 
109  G4double Energy = std::sqrt( pSquared + fMassCof );
110  G4double cof2 = Energy/c_light ;
111 
112  G4double pModuleInverse = 1.0/std::sqrt(pSquared) ;
113 
114  G4double inverse_velocity = Energy * pModuleInverse / c_light;
115 
116  G4double cof1 = fElectroMagCof*pModuleInverse ;
117 
118  dydx[0] = y[3]*pModuleInverse ;
119  dydx[1] = y[4]*pModuleInverse ;
120  dydx[2] = y[5]*pModuleInverse ;
121 
122  dydx[3] = cof1*(cof2*Field[3] + (y[4]*Field[2] - y[5]*Field[1])) ;
123 
124  dydx[4] = cof1*(cof2*Field[4] + (y[5]*Field[0] - y[3]*Field[2])) ;
125 
126  dydx[5] = cof1*(cof2*Field[5] + (y[3]*Field[1] - y[4]*Field[0])) ;
127 
128  dydx[6] = dydx[8] = 0.;//not used
129 
130  // Lab Time of flight
131  dydx[7] = inverse_velocity;
132 
133  G4ThreeVector BField(Field[0],Field[1],Field[2]);
134  G4ThreeVector EField(Field[3],Field[4],Field[5]);
135 
136  EField /= c_light;
137 
138  G4ThreeVector u(y[3], y[4], y[5]);
139  u *= pModuleInverse;
140 
141  G4double udb = anomaly*beta*gamma/(1.+gamma) * (BField * u);
142  G4double ucb = (anomaly+1./gamma)/beta;
143  G4double uce = anomaly + 1./(gamma+1.);
144 
145  G4ThreeVector Spin(y[9],y[10],y[11]);
146 
147  G4double pcharge;
148  if (charge == 0.) pcharge = 1.;
149  else pcharge = charge;
150 
151  G4ThreeVector dSpin(0.,0.,0.);
152  if (Spin.mag2() != 0.) {
153  dSpin =
154  pcharge*omegac*( ucb*(Spin.cross(BField))-udb*(Spin.cross(u))
155  // from Jackson
156  // -uce*Spin.cross(u.cross(EField)) );
157  // but this form has one less operation
158  - uce*(u*(Spin*EField) - EField*(Spin*u)) );
159  }
160 
161  dydx[ 9] = dSpin.x();
162  dydx[10] = dSpin.y();
163  dydx[11] = dSpin.z();
164 
165  return ;
166 }
static constexpr double hbar_Planck
G4double GetSpin() const
Float_t y
Definition: compare.C:6
double z() const
Double_t beta
static constexpr double c_squared
double G4double
Definition: G4Types.hh:76
double mag2() const
Hep3Vector cross(const Hep3Vector &) const
static constexpr double eplus
Definition: G4SIunits.hh:199
static constexpr double c_light
G4EqEMFieldWithSpin(G4ElectroMagneticField *emField)
G4double GetMagneticDipoleMoment() const
double x() const
G4double GetCharge() const
T sqr(const T &x)
Definition: templates.hh:145
void SetChargeMomentumMass(G4ChargeState particleCharge, G4double MomentumXc, G4double mass)
double y() const
void EvaluateRhsGivenB(const G4double y[], const G4double Field[], G4double dydx[]) const