source: trunk/source/geometry/magneticfield/src/G4RKG3_Stepper.cc@ 1036

Last change on this file since 1036 was 921, checked in by garnier, 17 years ago

en test de gl2ps. Problemes de libraries

File size: 7.2 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: G4RKG3_Stepper.cc,v 1.15 2007/08/21 10:17:41 tnikitin Exp $
28// GEANT4 tag $Name: geant4-09-02-cand-01 $
29//
30// -------------------------------------------------------------------
31
32#include "G4RKG3_Stepper.hh"
33#include "G4LineSection.hh"
34#include "G4Mag_EqRhs.hh"
35
36G4RKG3_Stepper::G4RKG3_Stepper(G4Mag_EqRhs *EqRhs)
37 : G4MagIntegratorStepper(EqRhs,6)
38{
39
40 fPtrMagEqOfMot=EqRhs;
41
42}
43
44G4RKG3_Stepper::~G4RKG3_Stepper()
45{
46}
47
48void G4RKG3_Stepper::Stepper( const G4double yInput[7],
49 const G4double dydx[7],
50 G4double Step,
51 G4double yOut[7],
52 G4double yErr[])
53{
54 G4double B[3];
55 G4int nvar = 6 ;
56 G4int i;
57 G4double by15 = 1. / 15. ; // was 0.066666666 ;
58
59 G4double yTemp[7], dydxTemp[6], yIn[7] ;
60 // Saving yInput because yInput and yOut can be aliases for same array
61 for(i=0;i<nvar;i++) yIn[i]=yInput[i];
62
63 G4double h = Step * 0.5;
64 hStep=Step;
65 // Do two half steps
66
67 StepNoErr(yIn, dydx,h, yTemp,B) ;
68
69 //Store Bfld for DistChord Calculation
70 for(i=0;i<3;i++)BfldIn[i]=B[i];
71
72 // RightHandSide(yTemp,dydxTemp) ;
73
74 GetEquationOfMotion()->EvaluateRhsGivenB(yTemp,B,dydxTemp) ;
75 StepNoErr(yTemp,dydxTemp,h,yOut,B);
76
77 // Store midpoint, chord calculation
78
79 fyMidPoint = G4ThreeVector( yTemp[0], yTemp[1], yTemp[2]);
80
81 // Do a full Step
82
83 h *= 2 ;
84 StepNoErr(yIn,dydx,h,yTemp,B); // ,beTemp2) ;
85 for(i=0;i<nvar;i++)
86 {
87 yErr[i] = yOut[i] - yTemp[i] ;
88 yOut[i] += yErr[i]*by15 ; // Provides 5th order of accuracy
89 }
90
91 //Store values for DistChord method
92
93 fyInitial = G4ThreeVector( yIn[0], yIn[1], yIn[2]);
94 fpInitial = G4ThreeVector( yIn[3], yIn[4], yIn[5]);
95 fyFinal = G4ThreeVector( yOut[0], yOut[1], yOut[2]);
96
97 // NormaliseTangentVector( yOut ); // Deleted
98}
99
100// ---------------------------------------------------------------------------
101
102// Integrator for RK from G3 with evaluation of error in solution and delta
103// geometry based on naive similarity with the case of uniform magnetic field.
104// B1[3] is input and is the first magnetic field values
105// B2[3] is output and is the final magnetic field values.
106
107void G4RKG3_Stepper::StepWithEst( const G4double*,
108 const G4double*,
109 G4double,
110 G4double*,
111 G4double&,
112 G4double&,
113 const G4double*,
114 G4double* )
115
116{
117 G4Exception("G4RKG3_Stepper::StepWithEst()", "ObsoleteMethod",
118 FatalException, "Method no longer used.");
119}
120
121// -----------------------------------------------------------------
122
123
124// Integrator RK Stepper from G3 with only two field evaluation per Step.
125// It is used in propagation initial Step by small substeps after solution
126// error and delta geometry considerations. B[3] is magnetic field which
127// is passed from substep to substep.
128
129void G4RKG3_Stepper::StepNoErr(const G4double tIn[7],
130 const G4double dydx[7],
131 G4double Step,
132 G4double tOut[7],
133 G4double B[3] ) // const
134
135{
136
137 // Copy and edit the routine above, to delete alpha2, beta2, ...
138 G4double K1[7],K2[7],K3[7],K4[7] ;
139 G4double tTemp[7], yderiv[6] ;
140
141 // Need Momentum value to give correct values to the coefficients in equation
142 // Integration on unit velocity, but tIn[3,4,5] is momentum
143 G4double mom,inverse_mom;
144 G4int i ;
145 const G4double c1=0.5,c2=0.125,c3=1./6.;
146
147 // GetEquationOfMotion()->EvaluateRhsReturnB(tIn,dydx,B1) ;
148 // Correction for momentum not a velocity
149 // Need the protection !!! must be not zero
150 mom=std::sqrt(tIn[3]*tIn[3]+tIn[4]*tIn[4]+tIn[5]*tIn[5]);
151 inverse_mom=1./mom;
152 for(i=0;i<3;i++)
153 {
154 K1[i] = Step * dydx[i+3]*inverse_mom;
155 tTemp[i] = tIn[i] + Step*(c1*tIn[i+3]*inverse_mom + c2*K1[i]) ;
156 tTemp[i+3] = tIn[i+3] + c1*K1[i]*mom ;
157
158 }
159
160 GetEquationOfMotion()->EvaluateRhsReturnB(tTemp,yderiv,B) ;
161
162
163 for(i=0;i<3;i++)
164 {
165 K2[i] = Step * yderiv[i+3]*inverse_mom;
166 tTemp[i+3] = tIn[i+3] + c1*K2[i]*mom ;
167 }
168
169 // Given B, calculate yderiv !
170 GetEquationOfMotion()->EvaluateRhsGivenB(tTemp,B,yderiv) ;
171
172 for(i=0;i<3;i++)
173 {
174 K3[i] = Step * yderiv[i+3]*inverse_mom;
175 tTemp[i] = tIn[i] + Step*(tIn[i+3]*inverse_mom + c1*K3[i]) ;
176 tTemp[i+3] = tIn[i+3] + K3[i]*mom ;
177 }
178
179
180 // Calculates y-deriv(atives) & returns B too!
181 GetEquationOfMotion()->EvaluateRhsReturnB(tTemp,yderiv,B) ;
182
183
184 for(i=0;i<3;i++) // Output trajectory vector
185 {
186 K4[i] = Step * yderiv[i+3]*inverse_mom;
187 tOut[i] = tIn[i] + Step*(tIn[i+3]*inverse_mom+ (K1[i] + K2[i] + K3[i])*c3) ;
188 tOut[i+3] = tIn[i+3] + mom*(K1[i] + 2*K2[i] + 2*K3[i] +K4[i])*c3 ;
189 }
190
191 // NormaliseTangentVector( tOut );
192
193
194}
195
196
197// ---------------------------------------------------------------------------
198
199 G4double G4RKG3_Stepper::DistChord() const
200 {
201 // Soon: must check whether h/R > 2 pi !!
202 // Method below is good only for < 2 pi
203 G4double distChord,distLine;
204
205 if (fyInitial != fyFinal) {
206 distLine= G4LineSection::Distline(fyMidPoint,fyInitial,fyFinal );
207
208 distChord = distLine;
209 }else{
210 distChord = (fyMidPoint-fyInitial).mag();
211 }
212
213
214 return distChord;
215
216 }
217
Note: See TracBrowser for help on using the repository browser.