source: trunk/source/geometry/magneticfield/src/G4MagHelicalStepper.cc@ 1326

Last change on this file since 1326 was 1231, checked in by garnier, 16 years ago

update from CVS

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