Geant4 11.4.0
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 (CERN, Google Summer of Code 2022), 05.11.2022
28// Supervision: John Apostolakis (CERN), Renee Fatemi, Soon Yung Jun (FNAL)
29// --------------------------------------------------------------------
30
31#include "G4BorisScheme.hh"
32#include "G4FieldUtils.hh"
33#include"G4SystemOfUnits.hh"
34#include "globals.hh"
36
37#include "G4EquationOfMotion.hh"
38
39using namespace field_utils;
40
42 : fEquation(equation), fnvar(nvar)
43{
44 if (nvar <= 0)
45 {
46 G4Exception("G4BorisScheme::G4BorisScheme()",
47 "GeomField0002", FatalException,
48 "Invalid number of variables; must be greater than zero!");
49 }
50}
51
52void G4BorisScheme::DoStep(const G4double restMass,const G4double charge,
53 const G4double yIn[], G4double yOut[],
54 G4double hstep) const
55{
58
59 // Used the scheme described in the following paper:
60 // 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,
67 const G4double /*crg*/, const G4double yIn[],
68 G4double yOut[], G4double hstep) const
69{
70 // Particle information
71 copy(yOut, yIn);
72
73 // Obtaining velocity
74 G4ThreeVector momentum_vec =G4ThreeVector(yIn[3],yIn[4],yIn[5]);
75 G4double momentum_mag = momentum_vec.mag();
76 G4ThreeVector momentum_dir =(1.0/momentum_mag)*momentum_vec;
77
78 G4double velocity_mag = momentum_mag*(c_l)/(std::sqrt(sqr(momentum_mag) +sqr(restMass)));
79 G4ThreeVector velocity = momentum_dir*velocity_mag;
80
81 // Obtaining the time step from the length step
82
83 hstep /= velocity_mag*CLHEP::m;
84
85 // Updating the Position
86 for(G4int i = 0; i <3; i++ )
87 {
88 G4double pos = yIn[i]/CLHEP::m;
89 pos += hstep*velocity[i];
90 yOut[i] = pos*CLHEP::m;
91 }
92}
93
94void G4BorisScheme::UpdateVelocity(const G4double restMass,
95 const G4double charge, const G4double yIn[],
96 G4double yOut[], G4double hstep) const
97{
98 // Particle information
99
100 G4ThreeVector momentum_vec = G4ThreeVector(yIn[3],yIn[4],yIn[5]);
101 G4double momentum_mag = momentum_vec.mag();
102 G4ThreeVector momentum_dir =(1.0/momentum_mag)*momentum_vec;
103
104 G4double gamma = std::sqrt(sqr(momentum_mag) + sqr(restMass))/restMass;
105
106 G4double mass = (restMass/c_squared)/CLHEP::kg;
107
108 // Obtaining velocity
109
110 G4double velocity_mag = momentum_mag*(c_l)/(std::sqrt(sqr(momentum_mag) +sqr(restMass)));
111 G4ThreeVector velocity = momentum_dir*velocity_mag;
112
113 // Obtaining the time step from the length step
114
115 hstep /= velocity_mag*CLHEP::m;
116
117 // Obtaining the field values
118
120 G4double fieldValue[6] ={0,0,0,0,0,0};
121 fEquation->EvaluateRhsReturnB(yIn, dydx, fieldValue);
122
123 // Initializing Vectors
124
127 copy(yOut, yIn);
128 for( G4int i = 0; i < 3; i++)
129 {
130 E[i] = fieldValue[i+3]/CLHEP::volt*CLHEP::meter;// FIXME - Check Units
131 B[i] = fieldValue[i]/CLHEP::tesla;
132 }
133
134 // Boris Algorithm
135
136 G4double qd = hstep*(charge/(2*mass*gamma));
137 G4ThreeVector h = qd*B;
138 G4ThreeVector u = velocity + qd*E;
139 G4double h_l = h[0]*h[0] + h[1]*h[1] + h[2]*h[2];
140 G4ThreeVector s_1 = (2*h)/(1 + h_l);
141 G4ThreeVector ud = u + (u + u.cross(h)).cross(s_1);
142 G4ThreeVector v_fi = ud +qd*E;
143 G4double v_mag = std::sqrt(v_fi.mag2());
144 G4ThreeVector v_dir = v_fi/v_mag;
145 G4double momen_mag = (restMass*v_mag)/(std::sqrt(c_l*c_l - v_mag*v_mag));
146 G4ThreeVector momen = momen_mag*v_dir;
147
148 // Storing the updated momentum
149
150 for(G4int i = 3; i < 6; i++)
151 {
152 yOut[i] = momen[i-3];
153 }
154}
155
156// ----------------------------------------------------------------------------------
157
158void G4BorisScheme::copy(G4double dst[], const G4double src[]) const
159{
160 std::memcpy(dst, src, sizeof(G4double) * fnvar);
161}
162
163// ----------------------------------------------------------------------------------
164// - Methods using the Boris Scheme Stepping to estimate integration error
165// ----------------------------------------------------------------------------------
167StepWithErrorEstimate(const G4double yIn[], G4double restMass,
168 G4double charge, G4double hstep,
169 G4double yOut[], G4double yErr[]) const
170{
171 // Use two half-steps (comparing to a full step) to obtain output
172 // and error estimate
173
175 StepWithMidAndErrorEstimate(yIn, restMass, charge, hstep, yMid, yOut, yErr);
176}
177
178// ----------------------------------------------------------------------------------
179
181StepWithMidAndErrorEstimate(const G4double yIn[], G4double restMass,
182 G4double charge, G4double hstep, G4double yMid[],
183 G4double yOut[], G4double yErr[]) const
184{
185 G4double halfStep= 0.5*hstep;
187
188 // In a single step
189 DoStep(restMass, charge, yIn, yOutAlt, hstep );
190
191 // Same, and also return mid-point evaluation
192 DoStep(restMass, charge, yIn, yMid, halfStep );
193 DoStep(restMass, charge, yMid, yOut, halfStep );
194
195 for( G4int i = 0; i < fnvar; i++ )
196 {
197 yErr[i] = yOutAlt[i] - yOut[i];
198 }
199}
G4double B(G4double temperature)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
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 DoStep(G4double restMass, G4double charge, const G4double yIn[], G4double yOut[], G4double hstep) const
G4BorisScheme()=default
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
G4EquationOfMotion is the abstract base class for the right hand size of the equation of motion of a ...
field_utils is a helper namespace, including simple methods to extract vectors from arrays in convent...
void copy(G4double dst[], const G4double src[], std::size_t size=G4FieldTrack::ncompSVEC)
T sqr(const T &x)
Definition templates.hh:128