Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MagHelicalStepper.cc
Go to the documentation of this file.
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// G4MagHelicalStepper implementation
27//
28// Given a purely magnetic field a better approach than adding a straight line
29// (as in the normal runge-kutta-methods) is to add helix segments to the
30// current position
31//
32// Author: John Apostolakis (CERN), 05.11.1998
33// --------------------------------------------------------------------
34
37#include "G4SystemOfUnits.hh"
38#include "G4LineSection.hh"
39#include "G4Mag_EqRhs.hh"
40
41// Constant for determining unit conversion when using normal as integrand.
42//
43const G4double G4MagHelicalStepper::fUnitConstant = 0.299792458*(GeV/(tesla*m));
44
46 : G4MagIntegratorStepper(EqRhs, 6), // integrate over 6 variables only !!
47 // position & velocity
48 fPtrMagEqOfMot(EqRhs)
49{
50}
51
52void
54 const G4ThreeVector& Bfld,
55 G4double h,
56 G4double yHelix[],
57 G4double yHelix2[] )
58{
59 // const G4int nvar = 6;
60
61 // OLD const G4double approc_limit = 0.05;
62 // OLD approc_limit = 0.05 gives max.error=x^5/5!=(0.05)^5/5!=2.6*e-9
63 // NEW approc_limit = 0.005 gives max.error=x^5/5!=2.6*e-14
64
65 const G4double approc_limit = 0.005;
66 G4ThreeVector Bnorm, B_x_P, vperp, vpar;
67
68 G4double B_d_P;
69 G4double B_v_P;
70 G4double Theta;
71 G4double R_1;
72 G4double R_Helix;
73 G4double CosT2, SinT2, CosT, SinT;
74 G4ThreeVector positionMove, endTangent;
75
76 G4double Bmag = Bfld.mag();
77 const G4double* pIn = yIn+3;
78 G4ThreeVector initVelocity = G4ThreeVector( pIn[0], pIn[1], pIn[2]);
79 G4double velocityVal = initVelocity.mag();
80 G4ThreeVector initTangent = (1.0/velocityVal) * initVelocity;
81
82 R_1 = GetInverseCurve(velocityVal,Bmag);
83
84 // for too small magnetic fields there is no curvature
85 // (include momentum here) FIXME
86
87 if( (std::fabs(R_1) < 1e-10)||(Bmag<1e-12) )
88 {
89 LinearStep( yIn, h, yHelix );
90
91 // Store and/or calculate parameters for chord distance
92
93 SetAngCurve(1.);
94 SetCurve(h);
95 SetRadHelix(0.);
96 }
97 else
98 {
99 Bnorm = (1.0/Bmag)*Bfld;
100
101 // calculate the direction of the force
102
103 B_x_P = Bnorm.cross(initTangent);
104
105 // parallel and perp vectors
106
107 B_d_P = Bnorm.dot(initTangent); // this is the fraction of P parallel to B
108
109 vpar = B_d_P * Bnorm; // the component parallel to B
110 vperp= initTangent - vpar; // the component perpendicular to B
111
112 B_v_P = std::sqrt( 1 - B_d_P * B_d_P); // Fraction of P perp to B
113
114 // calculate the stepping angle
115
116 Theta = R_1 * h; // * B_v_P;
117
118 // Trigonometrix
119
120 if( std::fabs(Theta) > approc_limit )
121 {
122 SinT = std::sin(Theta);
123 CosT = std::cos(Theta);
124 }
125 else
126 {
127 G4double Theta2 = Theta*Theta;
128 G4double Theta3 = Theta2 * Theta;
129 G4double Theta4 = Theta2 * Theta2;
130 SinT = Theta - 1.0/6.0 * Theta3;
131 CosT = 1 - 0.5 * Theta2 + 1.0/24.0 * Theta4;
132 }
133
134 // the actual "rotation"
135
136 G4double R = 1.0 / R_1;
137
138 positionMove = R * ( SinT * vperp + (1-CosT) * B_x_P) + h * vpar;
139 endTangent = CosT * vperp + SinT * B_x_P + vpar;
140
141 // Store the resulting position and tangent
142
143 yHelix[0] = yIn[0] + positionMove.x();
144 yHelix[1] = yIn[1] + positionMove.y();
145 yHelix[2] = yIn[2] + positionMove.z();
146 yHelix[3] = velocityVal * endTangent.x();
147 yHelix[4] = velocityVal * endTangent.y();
148 yHelix[5] = velocityVal * endTangent.z();
149
150 // Store 2*h step Helix if exist
151
152 if(yHelix2 != nullptr)
153 {
154 SinT2 = 2.0 * SinT * CosT;
155 CosT2 = 1.0 - 2.0 * SinT * SinT;
156 endTangent = (CosT2 * vperp + SinT2 * B_x_P + vpar);
157 positionMove = R * ( SinT2 * vperp + (1-CosT2) * B_x_P) + h*2 * vpar;
158
159 yHelix2[0] = yIn[0] + positionMove.x();
160 yHelix2[1] = yIn[1] + positionMove.y();
161 yHelix2[2] = yIn[2] + positionMove.z();
162 yHelix2[3] = velocityVal * endTangent.x();
163 yHelix2[4] = velocityVal * endTangent.y();
164 yHelix2[5] = velocityVal * endTangent.z();
165 }
166
167 // Store and/or calculate parameters for chord distance
168
169 G4double ptan=velocityVal*B_v_P;
170
171 G4double particleCharge = fPtrMagEqOfMot->FCof() / (eplus*c_light);
172 R_Helix =std::abs( ptan/(fUnitConstant * particleCharge*Bmag));
173
174 SetAngCurve(std::abs(Theta));
175 SetCurve(std::abs(R));
176 SetRadHelix(R_Helix);
177 }
178}
179
180// Use the midpoint method to get an error estimate and correction
181// modified from G4ClassicalRK4: W.Wander <wwc@mit.edu> 12/09/97
182//
183void
185 const G4double*,
186 G4double hstep,
187 G4double yOut[],
188 G4double yErr[] )
189{
190 const G4int nvar = 6;
191
192 // correction for Richardson Extrapolation.
193 // G4double correction = 1. / ( (1 << IntegratorOrder()) -1 );
194
195 G4double yTemp[8], yIn[8] ;
196 G4ThreeVector Bfld_initial, Bfld_midpoint;
197
198 // Saving yInput because yInput and yOut can be aliases for same array
199 //
200 for(G4int i=0; i<nvar; ++i)
201 {
202 yIn[i]=yInput[i];
203 }
204
205 G4double h = hstep * 0.5;
206
207 MagFieldEvaluate(yIn, Bfld_initial) ;
208
209 // Do two half steps
210 //
211 DumbStepper(yIn, Bfld_initial, h, yTemp);
212 MagFieldEvaluate(yTemp, Bfld_midpoint) ;
213 DumbStepper(yTemp, Bfld_midpoint, h, yOut);
214
215 // Do a full Step
216 //
217 h = hstep ;
218 DumbStepper(yIn, Bfld_initial, h, yTemp);
219
220 // Error estimation
221 //
222 for(G4int i=0; i<nvar; ++i)
223 {
224 yErr[i] = yOut[i] - yTemp[i] ;
225 }
226
227 return;
228}
229
232{
233 // Check whether h/R > pi !!
234 // Method DistLine is good only for < pi
235
236 G4double Ang=GetAngCurve();
237 if(Ang<=pi)
238 {
239 return GetRadHelix()*(1-std::cos(0.5*Ang));
240 }
241
242 if(Ang<twopi)
243 {
244 return GetRadHelix()*(1+std::cos(0.5*(twopi-Ang)));
245 }
246
247 // return Diameter of projected circle
248 return 2*GetRadHelix();
249}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
double z() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const
double mag() const
void SetCurve(const G4double Curve)
void SetRadHelix(const G4double Rad)
virtual void DumbStepper(const G4double y[], G4ThreeVector Bfld, G4double h, G4double yout[])=0
G4double GetRadHelix() const
void MagFieldEvaluate(const G4double y[], G4ThreeVector &Bfield)
G4double DistChord() const override
G4MagHelicalStepper(G4Mag_EqRhs *EqRhs)
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
G4double GetInverseCurve(const G4double Momentum, const G4double Bmag)
void LinearStep(const G4double yIn[], G4double h, G4double yHelix[]) const
void AdvanceHelix(const G4double yIn[], const G4ThreeVector &Bfld, G4double h, G4double yHelix[], G4double yHelix2[]=nullptr)
void SetAngCurve(const G4double Ang)
G4double GetAngCurve() const
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)
G4Mag_EqRhs is the "standard" equation of motion of a particle in a pure magnetic field.