Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CashKarpRKF45.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 //
27 // $Id: G4CashKarpRKF45.cc 97598 2016-06-06 07:19:46Z gcosmo $
28 //
29 // The Cash-Karp Runge-Kutta-Fehlberg 4/5 method is an embedded fourth
30 // order method (giving fifth-order accuracy) for the solution of an ODE.
31 // Two different fourth order estimates are calculated; their difference
32 // gives an error estimate. [ref. Numerical Recipes in C, 2nd Edition]
33 // It is used to integrate the equations of the motion of a particle
34 // in a magnetic field.
35 //
36 // [ref. Numerical Recipes in C, 2nd Edition]
37 //
38 // -------------------------------------------------------------------
39 
40 #include "G4CashKarpRKF45.hh"
41 #include "G4LineSection.hh"
42 
44 //
45 // Constructor
46 
48  G4int noIntegrationVariables,
49  G4bool primary)
50  : G4MagIntegratorStepper(EqRhs, noIntegrationVariables),
51  fLastStepLength(0.), fAuxStepper(0)
52 {
53  const G4int numberOfVariables =
54  std::max( noIntegrationVariables,
55  ( ( (noIntegrationVariables-1)/4 + 1 ) * 4 ) );
56  // For better alignment with cache-line
57 
58  ak2 = new G4double[numberOfVariables] ;
59  ak3 = new G4double[numberOfVariables] ;
60  ak4 = new G4double[numberOfVariables] ;
61  ak5 = new G4double[numberOfVariables] ;
62  ak6 = new G4double[numberOfVariables] ;
63  // ak7 = 0;
64 
65  // Must ensure space extra 'state' variables exists - i.e. yIn[7]
66  const G4int numStateMax = std::max(GetNumberOfStateVariables(), 8);
67  const G4int numStateVars = std::max(noIntegrationVariables,
68  numStateMax );
69  // GetNumberOfStateVariables() );
70 
71  yTemp = new G4double[numStateVars] ;
72  yIn = new G4double[numStateVars] ;
73 
74  fLastInitialVector = new G4double[numStateVars] ;
75  fLastFinalVector = new G4double[numStateVars] ;
76  fLastDyDx = new G4double[numberOfVariables];
77 
78  fMidVector = new G4double[numStateVars];
79  fMidError = new G4double[numStateVars];
80  if( primary )
81  {
82  fAuxStepper = new G4CashKarpRKF45(EqRhs, numberOfVariables, !primary);
83  }
84 }
85 
87 //
88 // Destructor
89 
91 {
92  delete[] ak2;
93  delete[] ak3;
94  delete[] ak4;
95  delete[] ak5;
96  delete[] ak6;
97  // delete[] ak7;
98  delete[] yTemp;
99  delete[] yIn;
100 
101  delete[] fLastInitialVector;
102  delete[] fLastFinalVector;
103  delete[] fLastDyDx;
104  delete[] fMidVector;
105  delete[] fMidError;
106 
107  delete fAuxStepper;
108 }
109 
111 //
112 // Given values for n = 6 variables yIn[0,...,n-1]
113 // known at x, use the fifth-order Cash-Karp Runge-
114 // Kutta-Fehlberg-4-5 method to advance the solution over an interval
115 // Step and return the incremented variables as yOut[0,...,n-1]. Also
116 // return an estimate of the local truncation error yErr[] using the
117 // embedded 4th-order method. The user supplies routine
118 // RightHandSide(y,dydx), which returns derivatives dydx for y .
119 
120 void
122  const G4double dydx[],
123  G4double Step,
124  G4double yOut[],
125  G4double yErr[])
126 {
127  // const G4int nvar = 6 ;
128  // const G4double a2 = 0.2 , a3 = 0.3 , a4 = 0.6 , a5 = 1.0 , a6 = 0.875;
129  G4int i;
130 
131  const G4double b21 = 0.2 ,
132  b31 = 3.0/40.0 , b32 = 9.0/40.0 ,
133  b41 = 0.3 , b42 = -0.9 , b43 = 1.2 ,
134 
135  b51 = -11.0/54.0 , b52 = 2.5 , b53 = -70.0/27.0 ,
136  b54 = 35.0/27.0 ,
137 
138  b61 = 1631.0/55296.0 , b62 = 175.0/512.0 ,
139  b63 = 575.0/13824.0 , b64 = 44275.0/110592.0 ,
140  b65 = 253.0/4096.0 ,
141 
142  c1 = 37.0/378.0 , c3 = 250.0/621.0 , c4 = 125.0/594.0 ,
143  c6 = 512.0/1771.0 ,
144  dc5 = -277.0/14336.0 ;
145 
146  const G4double dc1 = c1 - 2825.0/27648.0 , dc3 = c3 - 18575.0/48384.0 ,
147  dc4 = c4 - 13525.0/55296.0 , dc6 = c6 - 0.25 ;
148 
149  // Initialise time to t0, needed when it is not updated by the integration.
150  // [ Note: Only for time dependent fields (usually electric)
151  // is it neccessary to integrate the time.]
152  yOut[7] = yTemp[7] = yIn[7];
153 
154  const G4int numberOfVariables= this->GetNumberOfVariables();
155  // The number of variables to be integrated over
156 
157  // Saving yInput because yInput and yOut can be aliases for same array
158 
159  for(i=0;i<numberOfVariables;i++)
160  {
161  yIn[i]=yInput[i];
162  }
163  // RightHandSide(yIn, dydx) ; // 1st Step
164 
165  for(i=0;i<numberOfVariables;i++)
166  {
167  yTemp[i] = yIn[i] + b21*Step*dydx[i] ;
168  }
169  RightHandSide(yTemp, ak2) ; // 2nd Step
170 
171  for(i=0;i<numberOfVariables;i++)
172  {
173  yTemp[i] = yIn[i] + Step*(b31*dydx[i] + b32*ak2[i]) ;
174  }
175  RightHandSide(yTemp, ak3) ; // 3rd Step
176 
177  for(i=0;i<numberOfVariables;i++)
178  {
179  yTemp[i] = yIn[i] + Step*(b41*dydx[i] + b42*ak2[i] + b43*ak3[i]) ;
180  }
181  RightHandSide(yTemp, ak4) ; // 4th Step
182 
183  for(i=0;i<numberOfVariables;i++)
184  {
185  yTemp[i] = yIn[i] + Step*(b51*dydx[i] + b52*ak2[i] + b53*ak3[i] +
186  b54*ak4[i]) ;
187  }
188  RightHandSide(yTemp, ak5) ; // 5th Step
189 
190  for(i=0;i<numberOfVariables;i++)
191  {
192  yTemp[i] = yIn[i] + Step*(b61*dydx[i] + b62*ak2[i] + b63*ak3[i] +
193  b64*ak4[i] + b65*ak5[i]) ;
194  }
195  RightHandSide(yTemp, ak6) ; // 6th Step
196 
197  for(i=0;i<numberOfVariables;i++)
198  {
199  // Accumulate increments with proper weights
200 
201  yOut[i] = yIn[i] + Step*(c1*dydx[i] + c3*ak3[i] + c4*ak4[i] + c6*ak6[i]) ;
202 
203  // Estimate error as difference between 4th and
204  // 5th order methods
205 
206  yErr[i] = Step*(dc1*dydx[i] + dc3*ak3[i] + dc4*ak4[i] +
207  dc5*ak5[i] + dc6*ak6[i]) ;
208 
209  // Store Input and Final values, for possible use in calculating chord
210  fLastInitialVector[i] = yIn[i] ;
211  fLastFinalVector[i] = yOut[i];
212  fLastDyDx[i] = dydx[i];
213  }
214  // NormaliseTangentVector( yOut ); // Not wanted
215 
216  fLastStepLength =Step;
217 
218  return ;
219 }
220 
222 
223 void
224 G4CashKarpRKF45::StepWithEst( const G4double*,
225  const G4double*,
226  G4double,
227  G4double*,
228  G4double&,
229  G4double&,
230  const G4double*,
231  G4double* )
232 {
233  G4Exception("G4CashKarpRKF45::StepWithEst()", "GeomField0001",
234  FatalException, "Method no longer used.");
235  return ;
236 }
237 
239 
241 {
242  G4double distLine, distChord;
243  G4ThreeVector initialPoint, finalPoint, midPoint;
244 
245  // Store last initial and final points (they will be overwritten in self-Stepper call!)
246  initialPoint = G4ThreeVector( fLastInitialVector[0],
247  fLastInitialVector[1], fLastInitialVector[2]);
248  finalPoint = G4ThreeVector( fLastFinalVector[0],
249  fLastFinalVector[1], fLastFinalVector[2]);
250 
251  // Do half a step using StepNoErr
252 
253  fAuxStepper->Stepper( fLastInitialVector, fLastDyDx, 0.5 * fLastStepLength,
254  fMidVector, fMidError );
255 
256  midPoint = G4ThreeVector( fMidVector[0], fMidVector[1], fMidVector[2]);
257 
258  // Use stored values of Initial and Endpoint + new Midpoint to evaluate
259  // distance of Chord
260 
261 
262  if (initialPoint != finalPoint)
263  {
264  distLine = G4LineSection::Distline( midPoint, initialPoint, finalPoint );
265  distChord = distLine;
266  }
267  else
268  {
269  distChord = (midPoint-initialPoint).mag();
270  }
271  return distChord;
272 }
273 
274 
CLHEP::Hep3Vector G4ThreeVector
G4CashKarpRKF45(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
int G4int
Definition: G4Types.hh:78
G4int GetNumberOfVariables() const
bool G4bool
Definition: G4Types.hh:79
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double DistChord() const
G4int GetNumberOfStateVariables() const
void RightHandSide(const double y[], double dydx[])
double G4double
Definition: G4Types.hh:76