[1274] | 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: G4EqEMFieldWithEDM.cc,v 1.2 2009/11/06 22:31:54 gum Exp $ |
---|
| 28 | // GEANT4 tag $Name: $ |
---|
| 29 | // |
---|
| 30 | // |
---|
| 31 | // This is the standard right-hand side for equation of motion. |
---|
| 32 | // |
---|
| 33 | // 19.02.2009 Kevin Lynch, based on G4EqEMFieldWithSpin |
---|
| 34 | // 06.11.2009 Hiromi Iinuma see: |
---|
| 35 | // http://hypernews.slac.stanford.edu/HyperNews/geant4/get/emfields/161.html |
---|
| 36 | // |
---|
| 37 | // ------------------------------------------------------------------- |
---|
| 38 | |
---|
| 39 | #include "G4EqEMFieldWithEDM.hh" |
---|
| 40 | #include "G4ElectroMagneticField.hh" |
---|
| 41 | #include "G4ThreeVector.hh" |
---|
| 42 | #include "globals.hh" |
---|
| 43 | |
---|
| 44 | G4EqEMFieldWithEDM::G4EqEMFieldWithEDM(G4ElectroMagneticField *emField ) |
---|
| 45 | : G4EquationOfMotion( emField ) |
---|
| 46 | { |
---|
| 47 | anomaly = 0.0011659208; |
---|
| 48 | eta = 0.; |
---|
| 49 | } |
---|
| 50 | |
---|
| 51 | G4EqEMFieldWithEDM::~G4EqEMFieldWithEDM() |
---|
| 52 | { |
---|
| 53 | } |
---|
| 54 | |
---|
| 55 | void |
---|
| 56 | G4EqEMFieldWithEDM::SetChargeMomentumMass(G4double particleCharge, // e+ units |
---|
| 57 | G4double MomentumXc, |
---|
| 58 | G4double particleMass) |
---|
| 59 | { |
---|
| 60 | fElectroMagCof = eplus*particleCharge*c_light ; |
---|
| 61 | fMassCof = particleMass*particleMass ; |
---|
| 62 | |
---|
| 63 | omegac = 0.105658387*GeV/particleMass * 2.837374841e-3*(rad/cm/kilogauss); |
---|
| 64 | |
---|
| 65 | ParticleCharge = particleCharge; |
---|
| 66 | |
---|
| 67 | E = std::sqrt(sqr(MomentumXc)+sqr(particleMass)); |
---|
| 68 | beta = MomentumXc/E; |
---|
| 69 | gamma = E/particleMass; |
---|
| 70 | |
---|
| 71 | } |
---|
| 72 | |
---|
| 73 | void |
---|
| 74 | G4EqEMFieldWithEDM::EvaluateRhsGivenB(const G4double y[], |
---|
| 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, with additions for EDM |
---|
| 86 | // evolution from |
---|
| 87 | // M.Nowakowski, et.al. Eur.J.Phys.26, pp 545-560, (2005) |
---|
| 88 | // or |
---|
| 89 | // Silenko, Phys.Rev.ST Accel.Beams 9:034003, (2006) |
---|
| 90 | |
---|
| 91 | // dS/dt = (e/m) S \cross |
---|
| 92 | // MDM: [ (g/2-1 +1/\gamma) B |
---|
| 93 | // -(g/2-1)\gamma/(\gamma+1) (\beta \cdot B)\beta |
---|
| 94 | // -(g/2-\gamma/(\gamma+1) \beta \cross E |
---|
| 95 | // |
---|
| 96 | // EDM: eta/2( E - gamma/(gamma+1) \beta (\beta \cdot E) |
---|
| 97 | // + \beta \cross B ) ] |
---|
| 98 | // |
---|
| 99 | // where |
---|
| 100 | // S = \vec{s}, where S^2 = 1 |
---|
| 101 | // B = \vec{B} |
---|
| 102 | // \beta = \vec{\beta} = \beta \vec{u} with u^2 = 1 |
---|
| 103 | // E = \vec{E} |
---|
| 104 | |
---|
| 105 | G4double pSquared = y[3]*y[3] + y[4]*y[4] + y[5]*y[5] ; |
---|
| 106 | |
---|
| 107 | G4double Energy = std::sqrt( pSquared + fMassCof ); |
---|
| 108 | G4double cof2 = Energy/c_light ; |
---|
| 109 | |
---|
| 110 | G4double pModuleInverse = 1.0/std::sqrt(pSquared) ; |
---|
| 111 | |
---|
| 112 | G4double inverse_velocity = Energy * pModuleInverse / c_light; |
---|
| 113 | |
---|
| 114 | G4double cof1 = fElectroMagCof*pModuleInverse ; |
---|
| 115 | |
---|
| 116 | dydx[0] = y[3]*pModuleInverse ; |
---|
| 117 | dydx[1] = y[4]*pModuleInverse ; |
---|
| 118 | dydx[2] = y[5]*pModuleInverse ; |
---|
| 119 | |
---|
| 120 | dydx[3] = cof1*(cof2*Field[3] + (y[4]*Field[2] - y[5]*Field[1])) ; |
---|
| 121 | |
---|
| 122 | dydx[4] = cof1*(cof2*Field[4] + (y[5]*Field[0] - y[3]*Field[2])) ; |
---|
| 123 | |
---|
| 124 | dydx[5] = cof1*(cof2*Field[5] + (y[3]*Field[1] - y[4]*Field[0])) ; |
---|
| 125 | |
---|
| 126 | dydx[6] = dydx[8] = 0.;//not used |
---|
| 127 | |
---|
| 128 | // Lab Time of flight |
---|
| 129 | dydx[7] = inverse_velocity; |
---|
| 130 | |
---|
| 131 | G4ThreeVector BField(Field[0],Field[1],Field[2]); |
---|
| 132 | G4ThreeVector EField(Field[3],Field[4],Field[5]); |
---|
| 133 | |
---|
| 134 | EField /= c_light; |
---|
| 135 | |
---|
| 136 | G4ThreeVector u(y[3], y[4], y[5]); |
---|
| 137 | u *= pModuleInverse; |
---|
| 138 | |
---|
| 139 | G4double udb = anomaly*beta*gamma/(1.+gamma) * (BField * u); |
---|
| 140 | G4double ucb = (anomaly+1./gamma)/beta; |
---|
| 141 | G4double uce = anomaly + 1./(gamma+1.); |
---|
| 142 | G4double ude = beta*gamma/(1.+gamma)*(EField*u); |
---|
| 143 | |
---|
| 144 | G4ThreeVector Spin(y[9],y[10],y[11]); |
---|
| 145 | |
---|
| 146 | G4ThreeVector dSpin |
---|
| 147 | = ParticleCharge*omegac*( ucb*(Spin.cross(BField))-udb*(Spin.cross(u)) |
---|
| 148 | // from Jackson |
---|
| 149 | // -uce*Spin.cross(u.cross(EField)) ) |
---|
| 150 | // but this form has one less operation |
---|
| 151 | - uce*(u*(Spin*EField) - EField*(Spin*u)) |
---|
| 152 | +eta/2.*(Spin.cross(EField) - ude*(Spin.cross(u)) |
---|
| 153 | // +Spin.cross(u.cross(Bfield)) |
---|
| 154 | + (u*(Spin*BField) - BField*(Spin*u)) |
---|
| 155 | ) |
---|
| 156 | ); |
---|
| 157 | |
---|
| 158 | dydx[ 9] = dSpin.x(); |
---|
| 159 | dydx[10] = dSpin.y(); |
---|
| 160 | dydx[11] = dSpin.z(); |
---|
| 161 | |
---|
| 162 | return ; |
---|
| 163 | } |
---|