source: trunk/source/geometry/magneticfield/src/G4EqEMFieldWithSpin.cc@ 852

Last change on this file since 852 was 850, checked in by garnier, 17 years ago

geant4.8.2 beta

File size: 4.5 KB
Line 
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,v 1.2 2008/04/24 12:33:08 tnikitin Exp $
28// GEANT4 tag $Name: HEAD $
29//
30//
31// This is the standard right-hand side for equation of motion.
32//
33// The only case another is required is when using a moving reference
34// frame ... or extending the class to include additional Forces,
35// eg an electric field
36//
37// 30.08.2007 Chris Gong, Peter Gumplinger
38//
39// -------------------------------------------------------------------
40
41#include "G4EqEMFieldWithSpin.hh"
42#include "G4ThreeVector.hh"
43#include "globals.hh"
44
45G4EqEMFieldWithSpin::G4EqEMFieldWithSpin(G4ElectroMagneticField *emField )
46 : G4EquationOfMotion( emField ) { anomaly = 1.165923e-3; }
47
48void
49G4EqEMFieldWithSpin::SetChargeMomentumMass(G4double particleCharge, // e+ units
50 G4double MomentumXc,
51 G4double particleMass)
52{
53 fElectroMagCof = eplus*particleCharge*c_light ;
54 fMassCof = particleMass*particleMass ;
55
56 omegac = 0.105658387*GeV/particleMass * 2.837374841e-3*(rad/cm/kilogauss);
57
58 ParticleCharge = particleCharge;
59
60 E = std::sqrt(sqr(MomentumXc)+sqr(particleMass));
61 beta = MomentumXc/E;
62 gamma = E/particleMass;
63}
64
65
66
67void
68G4EqEMFieldWithSpin::EvaluateRhsGivenB(const G4double y[],
69 const G4double Field[],
70 G4double dydx[] ) const
71{
72
73 // Components of y:
74 // 0-2 dr/ds,
75 // 3-5 dp/ds - momentum derivatives
76
77 G4double pSquared = y[3]*y[3] + y[4]*y[4] + y[5]*y[5] ;
78
79 G4double Energy = std::sqrt( pSquared + fMassCof );
80 G4double cof2 = Energy/c_light ;
81
82 G4double pModuleInverse = 1.0/std::sqrt(pSquared) ;
83
84 // G4double inverse_velocity = Energy * c_light * pModuleInverse;
85 G4double inverse_velocity = Energy * pModuleInverse / c_light;
86
87 G4double cof1 = fElectroMagCof*pModuleInverse ;
88
89 // G4double vDotE = y[3]*Field[3] + y[4]*Field[4] + y[5]*Field[5] ;
90
91
92 dydx[0] = y[3]*pModuleInverse ;
93 dydx[1] = y[4]*pModuleInverse ;
94 dydx[2] = y[5]*pModuleInverse ;
95
96 dydx[3] = cof1*(cof2*Field[3] + (y[4]*Field[2] - y[5]*Field[1])) ;
97
98 dydx[4] = cof1*(cof2*Field[4] + (y[5]*Field[0] - y[3]*Field[2])) ;
99
100 dydx[5] = cof1*(cof2*Field[5] + (y[3]*Field[1] - y[4]*Field[0])) ;
101
102 dydx[6] = dydx[8] = 0.;//not used
103
104 // Lab Time of flight
105 dydx[7] = inverse_velocity;
106
107 G4ThreeVector BField(Field[0],Field[1],Field[2]);
108
109 G4ThreeVector u(y[3], y[4], y[5]);
110 u *= pModuleInverse;
111
112 G4double udb = anomaly*beta*gamma/(1.+gamma) * (BField * u);
113 G4double ucb = (anomaly+1./gamma)/beta;
114
115 G4ThreeVector Spin(y[9],y[10],y[11]);
116 G4ThreeVector dSpin;
117
118 dSpin = ParticleCharge*omegac*(ucb*(Spin.cross(BField))-udb*(Spin.cross(u)));
119
120 dydx[ 9] = dSpin.x();
121 dydx[10] = dSpin.y();
122 dydx[11] = dSpin.z();
123
124 return ;
125}
Note: See TracBrowser for help on using the repository browser.