Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BorisScheme.cc
Go to the documentation of this file.
1// ********************************************************************
2// * License and Disclaimer *
3// * *
4// * The Geant4 software is copyright of the Copyright Holders of *
5// * the Geant4 Collaboration. It is provided under the terms and *
6// * conditions of the Geant4 Software License, included in the file *
7// * LICENSE and available at http://cern.ch/geant4/license . These *
8// * include a list of copyright holders. *
9// * *
10// * Neither the authors of this software system, nor their employing *
11// * institutes,nor the agencies providing financial support for this *
12// * work make any representation or warranty, express or implied, *
13// * regarding this software system or assume any liability for its *
14// * use. Please see the license in the file LICENSE and URL above *
15// * for the full disclaimer and the limitation of liability. *
16// * *
17// * This code implementation is the result of the scientific and *
18// * technical work of the GEANT4 collaboration. *
19// * By using, copying, modifying or distributing the software (or *
20// * any work based on the software) you agree to acknowledge its *
21// * use in resulting scientific publications, and indicate your *
22// * acceptance of all terms of the Geant4 Software license. *
23// ********************************************************************
24//
25// G4BorisScheme implementation
26//
27// Author: Divyansh Tiwari, Google Summer of Code 2022
28// Supervision: John Apostolakis,Renee Fatemi, Soon Yung Jun
29// --------------------------------------------------------------------
30
31#include "G4BorisScheme.hh"
32#include "G4FieldUtils.hh"
33#include"G4SystemOfUnits.hh"
34#include "globals.hh"
36
37#include "G4EquationOfMotion.hh"
38//#include "G4EqMagElectricField.hh"
39
40using namespace field_utils;
41
43 G4int nvar )
44 : fEquation(equation), fnvar(nvar)
45{
46 if (nvar <= 0)
47 {
48 G4Exception("G4BorisScheme::G4BorisScheme()",
49 "GeomField0002", FatalException,
50 "Invalid number of variables; must be greater than zero!");
51 }
52}
53
54void G4BorisScheme::DoStep(const G4double restMass,const G4double charge, const G4double yIn[],
55 G4double yOut[], G4double hstep) const
56{
59
60 // Used the scheme described in the following paper:https://www.research-collection.ethz.ch/bitstream/handle/20.500.11850/153167/eth-5175-01.pdf?sequence=1
61 UpdatePosition(restMass, charge, yIn, yOut1Temp, hstep/2);
62 UpdateVelocity(restMass, charge, yOut1Temp, yOut2Temp, hstep);
63 UpdatePosition(restMass, charge, yOut2Temp, yOut, hstep/2);
64}
65
66void G4BorisScheme::UpdatePosition(const G4double restMass, const G4double /*charge*/, const G4double yIn[],
67 G4double yOut[], G4double hstep) const
68{
69 // Particle information
70 copy(yOut, yIn);
71
72 // Obtaining velocity
73 G4ThreeVector momentum_vec =G4ThreeVector(yIn[3],yIn[4],yIn[5]);
74 G4double momentum_mag = momentum_vec.mag();
75 G4ThreeVector momentum_dir =(1.0/momentum_mag)*momentum_vec;
76
77 G4double velocity_mag = momentum_mag*(c_l)/(std::sqrt(sqr(momentum_mag) +sqr(restMass)));
78 G4ThreeVector velocity = momentum_dir*velocity_mag;
79
80 //Obtaining the time step from the length step
81
82 hstep /= velocity_mag*CLHEP::m;
83
84 // Updating the Position
85 for(G4int i = 0; i <3; i++ )
86 {
87 G4double pos = yIn[i]/CLHEP::m;
88 pos += hstep*velocity[i];
89 yOut[i] = pos*CLHEP::m;
90 }
91}
92
93void G4BorisScheme::UpdateVelocity(const G4double restMass, const G4double charge, const G4double yIn[],
94 G4double yOut[], G4double hstep) const
95{
96 //Particle information
97 G4ThreeVector momentum_vec =G4ThreeVector(yIn[3],yIn[4],yIn[5]);
98 G4double momentum_mag = momentum_vec.mag();
99 G4ThreeVector momentum_dir =(1.0/momentum_mag)*momentum_vec;
100
101 G4double gamma = std::sqrt(sqr(momentum_mag) + sqr(restMass))/restMass;
102
103 G4double mass = (restMass/c_squared)/CLHEP::kg;
104
105 //Obtaining velocity
106
107 G4double velocity_mag = momentum_mag*(c_l)/(std::sqrt(sqr(momentum_mag) +sqr(restMass)));
108 G4ThreeVector velocity = momentum_dir*velocity_mag;
109
110 ////Obtaining the time step from the length step
111
112 hstep /= velocity_mag*CLHEP::m;
113
114 // Obtaining the field values
116 G4double fieldValue[6] ={0,0,0,0,0,0};
117 fEquation->EvaluateRhsReturnB(yIn, dydx, fieldValue);
118
119 //Initializing Vectors
122 copy(yOut, yIn);
123 for( G4int i = 0; i < 3; i++)
124 {
125 E[i] = fieldValue[i+3]/CLHEP::volt*CLHEP::meter;// FIXME - Check Units
126 B[i] = fieldValue[i]/CLHEP::tesla;
127 }
128
129 //Boris Algorithm
130 G4double qd = hstep*(charge/(2*mass*gamma));
131 G4ThreeVector h = qd*B;
132 G4ThreeVector u = velocity + qd*E;
133 G4double h_l = h[0]*h[0] + h[1]*h[1] + h[2]*h[2];
134 G4ThreeVector s_1 = (2*h)/(1 + h_l);
135 G4ThreeVector ud = u + (u + u.cross(h)).cross(s_1);
136 G4ThreeVector v_fi = ud +qd*E;
137 G4double v_mag = std::sqrt(v_fi.mag2());
138 G4ThreeVector v_dir = v_fi/v_mag;
139 G4double momen_mag = (restMass*v_mag)/(std::sqrt(c_l*c_l - v_mag*v_mag));
140 G4ThreeVector momen = momen_mag*v_dir;
141
142 // Storing the updated momentum
143 for(int i = 3; i < 6; i++)
144 {
145 yOut[i] = momen[i-3];
146 }
147}
148
149// ----------------------------------------------------------------------------------
150
151void G4BorisScheme::copy(G4double dst[], const G4double src[]) const
152{
153 std::memcpy(dst, src, sizeof(G4double) * fnvar);
154}
155
156// ----------------------------------------------------------------------------------
157// - Methods using the Boris Scheme Stepping to estimate integration error
158// ----------------------------------------------------------------------------------
160StepWithErrorEstimate(const G4double yIn[], G4double restMass, G4double charge, G4double hstep,
161 G4double yOut[], G4double yErr[]) const
162{
163 // Use two half-steps (comparing to a full step) to obtain output and error estimate
165 StepWithMidAndErrorEstimate(yIn, restMass, charge, hstep, yMid, yOut, yErr);
166}
167
168// ----------------------------------------------------------------------------------
169
171StepWithMidAndErrorEstimate(const G4double yIn[], G4double restMass, G4double charge, G4double hstep,
172 G4double yMid[], G4double yOut[], G4double yErr[]
173 ) const
174{
175 G4double halfStep= 0.5*hstep;
177
178 // In a single step
179 DoStep(restMass, charge, yIn, yOutAlt, hstep );
180
181 // Same, and also return mid-point evaluation
182 DoStep(restMass, charge, yIn, yMid, halfStep );
183 DoStep(restMass, charge, yMid, yOut, halfStep );
184
185 for( G4int i= 0; i<fnvar; i++ )
186 {
187 yErr[i] = yOutAlt[i] - yOut[i];
188 }
189}
G4double B(G4double temperature)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
double mag2() const
Hep3Vector cross(const Hep3Vector &) const
double mag() const
void UpdateVelocity(const G4double restMass, const G4double charge, const G4double yIn[], G4double yOut[], G4double hstep) const
void DoStep(G4double restMass, G4double charge, const G4double yIn[], G4double yOut[], G4double hstep) const
G4BorisScheme()=default
void UpdatePosition(const G4double restMass, const G4double charge, const G4double yIn[], G4double yOut[], G4double hstep) const
void StepWithErrorEstimate(const G4double yIn[], G4double restMass, G4double charge, G4double hstep, G4double yOut[], G4double yErr[]) const
void StepWithMidAndErrorEstimate(const G4double yIn[], G4double restMass, G4double charge, G4double hstep, G4double yMid[], G4double yOut[], G4double yErr[]) const
void EvaluateRhsReturnB(const G4double y[], G4double dydx[], G4double Field[]) const
T sqr(const T &x)
Definition: templates.hh:128