source: trunk/source/geometry/magneticfield/src/G4ConstRK4.cc@ 1199

Last change on this file since 1199 was 1058, checked in by garnier, 17 years ago

file release beta

File size: 7.7 KB
RevLine 
[1058]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: G4ConstRK4.cc,v 1.2 2008/10/29 14:17:42 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-02-cand-01 $
29//
30//
31// - 18.09.2008 - J.Apostolakis, T.Nikitina - Created
32// -------------------------------------------------------------------
33
34#include "G4ConstRK4.hh"
35#include "G4ThreeVector.hh"
36#include "G4LineSection.hh"
37
38//////////////////////////////////////////////////////////////////
39//
40// Constructor sets the number of variables (default = 8)
41
42G4ConstRK4::G4ConstRK4(G4Mag_EqRhs* EqRhs, G4int numberOfVariables)
43 : G4MagErrorStepper(EqRhs, numberOfVariables)
44{
45 if(numberOfVariables !=8 )
46 {
47 G4Exception("G4ConstRK4::G4ConstRK4()", "InvalidSetup", FatalException,
48 "Valid only for number of variables=8. Use another Stepper!");
49 }
50 else
51 {
52 fEq=EqRhs;
53 yMiddle= new G4double[8];
54 dydxMid= new G4double[8];
55 yInitial= new G4double[8];
56 yOneStep= new G4double[8];
57
58 dydxm = new G4double[8];
59 dydxt = new G4double[8];
60 yt = new G4double[8];
61 Field[0]=0.;Field[1]=0.;Field[2]=0.;
62 }
63}
64
65////////////////////////////////////////////////////////////////
66//
67// Destructor
68
69G4ConstRK4::~G4ConstRK4()
70{
71 delete [] yMiddle;
72 delete [] dydxMid;
73 delete [] yInitial;
74 delete [] yOneStep;
75 delete [] dydxm;
76 delete [] dydxt;
77 delete [] yt;
78}
79
80//////////////////////////////////////////////////////////////////////
81//
82// Given values for the variables y[0,..,n-1] and their derivatives
83// dydx[0,...,n-1] known at x, use the classical 4th Runge-Kutta
84// method to advance the solution over an interval h and return the
85// incremented variables as yout[0,...,n-1], which is not a distinct
86// array from y. The user supplies the routine RightHandSide(x,y,dydx),
87// which returns derivatives dydx at x. The source is routine rk4 from
88// NRC p. 712-713 .
89
90void G4ConstRK4::DumbStepper( const G4double yIn[],
91 const G4double dydx[],
92 G4double h,
93 G4double yOut[])
94{
95 G4double hh = h*0.5 , h6 = h/6.0 ;
96
97 // 1st Step K1=h*dydx
98 yt[5] = yIn[5] + hh*dydx[5] ;
99 yt[4] = yIn[4] + hh*dydx[4] ;
100 yt[3] = yIn[3] + hh*dydx[3] ;
101 yt[2] = yIn[2] + hh*dydx[2] ;
102 yt[1] = yIn[1] + hh*dydx[1] ;
103 yt[0] = yIn[0] + hh*dydx[0] ;
104 RightHandSideConst(yt,dydxt) ;
105
106 // 2nd Step K2=h*dydxt
107 yt[5] = yIn[5] + hh*dydxt[5] ;
108 yt[4] = yIn[4] + hh*dydxt[4] ;
109 yt[3] = yIn[3] + hh*dydxt[3] ;
110 yt[2] = yIn[2] + hh*dydxt[2] ;
111 yt[1] = yIn[1] + hh*dydxt[1] ;
112 yt[0] = yIn[0] + hh*dydxt[0] ;
113 RightHandSideConst(yt,dydxm) ;
114
115 // 3rd Step K3=h*dydxm
116 // now dydxm=(K2+K3)/h
117 yt[5] = yIn[5] + h*dydxm[5] ;
118 dydxm[5] += dydxt[5] ;
119 yt[4] = yIn[4] + h*dydxm[4] ;
120 dydxm[4] += dydxt[4] ;
121 yt[3] = yIn[3] + h*dydxm[3] ;
122 dydxm[3] += dydxt[3] ;
123 yt[2] = yIn[2] + h*dydxm[2] ;
124 dydxm[2] += dydxt[2] ;
125 yt[1] = yIn[1] + h*dydxm[1] ;
126 dydxm[1] += dydxt[1] ;
127 yt[0] = yIn[0] + h*dydxm[0] ;
128 dydxm[0] += dydxt[0] ;
129 RightHandSideConst(yt,dydxt) ;
130
131 // 4th Step K4=h*dydxt
132 yOut[5] = yIn[5]+h6*(dydx[5]+dydxt[5]+2.0*dydxm[5]);
133 yOut[4] = yIn[4]+h6*(dydx[4]+dydxt[4]+2.0*dydxm[4]);
134 yOut[3] = yIn[3]+h6*(dydx[3]+dydxt[3]+2.0*dydxm[3]);
135 yOut[2] = yIn[2]+h6*(dydx[2]+dydxt[2]+2.0*dydxm[2]);
136 yOut[1] = yIn[1]+h6*(dydx[1]+dydxt[1]+2.0*dydxm[1]);
137 yOut[0] = yIn[0]+h6*(dydx[0]+dydxt[0]+2.0*dydxm[0]);
138
139} // end of DumbStepper ....................................................
140
141////////////////////////////////////////////////////////////////
142//
143// Stepper
144
145void
146G4ConstRK4::Stepper( const G4double yInput[],
147 const G4double dydx[],
148 G4double hstep,
149 G4double yOutput[],
150 G4double yError [] )
151{
152 const G4int nvar = 8 ;
153 const G4int maxvar= 8;
154
155 G4int i;
156
157 // Correction for Richardson extrapolation
158 //
159 G4double correction = 1. / ( (1 << IntegratorOrder()) -1 );
160
161 // Saving yInput because yInput and yOutput can be aliases for same array
162
163 for (i=0;i<nvar;i++) { yInitial[i]=yInput[i]; }
164
165 yInitial[7]= yInput[7]; // Copy the time in case...even if not really needed
166 yMiddle[7] = yInput[7]; // Copy the time from initial value
167 yOneStep[7] = yInput[7]; // As it contributes to final value of yOutput ?
168 yOutput[7] = yInput[7]; // -> dumb stepper does it too for RK4
169 for (i=nvar;i<maxvar;i++) { yOutput[i]=yInput[i]; }
170 yError[7] = 0.0;
171
172 G4double halfStep = hstep * 0.5;
173
174 // Do two half steps
175 //
176 GetConstField(yInitial,Field);
177 DumbStepper (yInitial, dydx, halfStep, yMiddle);
178 RightHandSideConst(yMiddle, dydxMid);
179 DumbStepper (yMiddle, dydxMid, halfStep, yOutput);
180
181 // Store midpoint, chord calculation
182 //
183 fMidPoint = G4ThreeVector( yMiddle[0], yMiddle[1], yMiddle[2]);
184
185 // Do a full Step
186 //
187 DumbStepper(yInitial, dydx, hstep, yOneStep);
188 for(i=0;i<nvar;i++)
189 {
190 yError [i] = yOutput[i] - yOneStep[i] ;
191 yOutput[i] += yError[i]*correction ;
192 // Provides accuracy increased by 1 order via the
193 // Richardson extrapolation
194 }
195
196 fInitialPoint = G4ThreeVector( yInitial[0], yInitial[1], yInitial[2]);
197 fFinalPoint = G4ThreeVector( yOutput[0], yOutput[1], yOutput[2]);
198
199 return;
200}
201
202////////////////////////////////////////////////////////////////
203//
204// Estimate the maximum distance from the curve to the chord
205//
206// We estimate this using the distance of the midpoint to chord.
207// The method below is good only for angle deviations < 2 pi;
208// this restriction should not be a problem for the Runge Kutta methods,
209// which generally cannot integrate accurately for large angle deviations
210
211G4double G4ConstRK4::DistChord() const
212{
213 G4double distLine, distChord;
214
215 if (fInitialPoint != fFinalPoint)
216 {
217 distLine= G4LineSection::Distline( fMidPoint, fInitialPoint, fFinalPoint );
218 // This is a class method that gives distance of Mid
219 // from the Chord between the Initial and Final points
220 distChord = distLine;
221 }
222 else
223 {
224 distChord = (fMidPoint-fInitialPoint).mag();
225 }
226 return distChord;
227}
Note: See TracBrowser for help on using the repository browser.