Geant4-11
G4TCashKarpRKF45.hh
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// G4TCashKarpRKF45
27//
28// Class description:
29//
30// Templated version of Cash-Karp 4th/5th order embedded stepper
31//
32// Knowing the type (class) of the equation of motion enables a non-
33// virtual call of its methods.
34// As an embedded 5th order method, it requires fewer field evaluations
35// (1 initial + 5 others per step = 6 per step) than ClassicalRK4 and
36// also non-embedded methods of the same order.
37//
38// Can be used to enable use of non-virtual calls for field, equation,
39// and stepper - potentially with inlined methods.
40//
41// Created: Josh Xie June 2014 (supported by Google Summer of Code 2014 )
42// Supervisors: Sandro Wenzel, John Apostolakis (CERN)
43//
44// Adapted from G4CashKarpRKF45 class
45// --------------------------------------------------------------------
46// Original description (G4CashKarpRKF45):
47// The Cash-Karp Runge-Kutta-Fehlberg 4/5 method is an embedded fourth
48// order method (giving fifth-order accuracy) for the solution of an ODE.
49// Two different fourth order estimates are calculated; their difference
50// gives an error estimate. [ref. Numerical Recipes in C, 2nd Edition]
51// Used to integrate the equations of motion of a particle in a field.
52// Original Authors: J.Apostolakis, V.Grichine - 30.01.1997
53
54#ifndef G4T_CASH_KARP_RKF45_HH
55#define G4T_CASH_KARP_RKF45_HH
56
57#include <cassert>
58
59#include "G4LineSection.hh"
61
62template <class T_Equation, unsigned int N = 6 >
64{
65 public:
66
67 G4TCashKarpRKF45(T_Equation* EqRhs, // G4int noIntegrationVariables = 6,
68 G4bool primary = true);
69
70 virtual ~G4TCashKarpRKF45();
71
72 inline void
73 StepWithError(const G4double yInput[], // * __restrict__ yInput,
74 const G4double dydx[], // * __restrict__ dydx,
75 G4double Step,
76 G4double yOut[], // * __restrict__ yOut,
77 G4double yErr[] ); // * __restrict__ yErr);
78
79 virtual void Stepper(const G4double yInput[],
80 const G4double dydx[],
81 G4double hstep,
82 G4double yOutput[],
83 G4double yError[]) override final;
84
85 // __attribute__((always_inline))
86 void RightHandSideInl( const G4double y[], // * __restrict__ y,
87 G4double dydx[] ) // * __restrict__ dydx )
88 {
89 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
90 }
91
92 inline G4double DistChord() const override;
93
94 inline G4int IntegratorOrder() const override { return 4; }
95
96 private:
99 // private copy constructor and assignment operator.
100
101 private:
102 G4double ak2[N], ak3[N], ak4[N], ak5[N], ak6[N], ak7[N], yTemp[N], yIn[N];
103 // scratch space
104
111 // for DistChord calculations
112
114 // ... or G4TCashKarpRKF45<T_Equation, N>* fAuxStepper;
115 T_Equation* fEquation_Rhs;
116};
117
119//
120// Constructor
121//
122template <class T_Equation, unsigned int N >
124 G4bool primary)
125 : G4MagIntegratorStepper(dynamic_cast<G4EquationOfMotion*>(EqRhs), N )
126 , fEquation_Rhs(EqRhs)
127{
128 if( dynamic_cast<G4EquationOfMotion*>(EqRhs) == nullptr )
129 {
130 G4Exception("G4TCashKarpRKF45: constructor", "GeomField0001",
131 FatalException, "Equation is not an G4EquationOfMotion.");
132 }
133
135 fLastFinalVector = new G4double[N];
136 fLastDyDx = new G4double[N];
137
138 fMidVector = new G4double[N];
139 fMidError = new G4double[N];
140
141 if(primary)
142 {
143 fAuxStepper = new G4TCashKarpRKF45<T_Equation, N> (EqRhs, !primary);
144 }
145}
146
147template <class T_Equation, unsigned int N >
149{
150 delete[] fLastInitialVector;
151 delete[] fLastFinalVector;
152 delete[] fLastDyDx;
153 delete[] fMidVector;
154 delete[] fMidError;
155
156 delete fAuxStepper;
157}
158
160//
161// Given values for n = 6 variables yIn[0,...,n-1]
162// known at x, use the fifth-order Cash-Karp Runge-
163// Kutta-Fehlberg-4-5 method to advance the solution over an interval
164// Step and return the incremented variables as yOut[0,...,n-1]. Also
165// return an estimate of the local truncation error yErr[] using the
166// embedded 4th-order method. The equation's method is called (inline)
167// via RightHandSideInl(y,dydx), which returns derivatives dydx for y .
168//
169template <class T_Equation, unsigned int N >
170inline void
172 const G4double* dydx,
173 G4double Step,
174 G4double * yOut,
175 G4double * yErr)
176{
177 // const G4double a2 = 0.2 , a3 = 0.3 , a4 = 0.6 , a5 = 1.0 , a6 = 0.875;
178
179 const G4double b21 = 0.2, b31 = 3.0 / 40.0, b32 = 9.0 / 40.0, b41 = 0.3,
180 b42 = -0.9, b43 = 1.2,
181
182 b51 = -11.0 / 54.0, b52 = 2.5, b53 = -70.0 / 27.0,
183 b54 = 35.0 / 27.0,
184
185 b61 = 1631.0 / 55296.0, b62 = 175.0 / 512.0,
186 b63 = 575.0 / 13824.0, b64 = 44275.0 / 110592.0,
187 b65 = 253.0 / 4096.0,
188
189 c1 = 37.0 / 378.0, c3 = 250.0 / 621.0, c4 = 125.0 / 594.0,
190 c6 = 512.0 / 1771.0, dc5 = -277.0 / 14336.0;
191
192 const G4double dc1 = c1 - 2825.0 / 27648.0, dc3 = c3 - 18575.0 / 48384.0,
193 dc4 = c4 - 13525.0 / 55296.0, dc6 = c6 - 0.25;
194
195 // Initialise time to t0, needed when it is not updated by the integration.
196 // [ Note: Only for time dependent fields (usually electric)
197 // is it neccessary to integrate the time.]
198 // yOut[7] = yTemp[7] = yIn[7];
199
200 // Saving yInput because yInput and yOut can be aliases for same array
201 for(unsigned int i = 0; i < N; ++i)
202 {
203 yIn[i] = yInput[i];
204 }
205 // RightHandSideInl(yIn, dydx) ; // 1st Step
206
207 for(unsigned int i = 0; i < N; ++i)
208 {
209 yTemp[i] = yIn[i] + b21 * Step * dydx[i];
210 }
211 this->RightHandSideInl(yTemp, ak2); // 2nd Step
212
213 for(unsigned int i = 0; i < N; ++i)
214 {
215 yTemp[i] = yIn[i] + Step * (b31 * dydx[i] + b32 * ak2[i]);
216 }
217 this->RightHandSideInl(yTemp, ak3); // 3rd Step
218
219 for(unsigned int i = 0; i < N; ++i)
220 {
221 yTemp[i] = yIn[i] + Step * (b41 * dydx[i] + b42 * ak2[i] + b43 * ak3[i]);
222 }
223 this->RightHandSideInl(yTemp, ak4); // 4th Step
224
225 for(unsigned int i = 0; i < N; ++i)
226 {
227 yTemp[i] = yIn[i] + Step * (b51 * dydx[i] + b52 * ak2[i] + b53 * ak3[i] +
228 b54 * ak4[i]);
229 }
230 this->RightHandSideInl(yTemp, ak5); // 5th Step
231
232 for(unsigned int i = 0; i < N; ++i)
233 {
234 yTemp[i] = yIn[i] + Step * (b61 * dydx[i] + b62 * ak2[i] + b63 * ak3[i] +
235 b64 * ak4[i] + b65 * ak5[i]);
236 }
237 this->RightHandSideInl(yTemp, ak6); // 6th Step
238
239 for(unsigned int i = 0; i < N; ++i)
240 {
241 // Accumulate increments with proper weights
242
243 yOut[i] = yIn[i] +
244 Step * (c1 * dydx[i] + c3 * ak3[i] + c4 * ak4[i] + c6 * ak6[i]);
245 }
246 for(unsigned int i = 0; i < N; ++i)
247 {
248 // Estimate error as difference between 4th and
249 // 5th order methods
250
251 yErr[i] = Step * (dc1 * dydx[i] + dc3 * ak3[i] + dc4 * ak4[i] +
252 dc5 * ak5[i] + dc6 * ak6[i]);
253 }
254 for(unsigned int i = 0; i < N; ++i)
255 {
256 // Store Input and Final values, for possible use in calculating chord
257 fLastInitialVector[i] = yIn[i];
258 fLastFinalVector[i] = yOut[i];
259 fLastDyDx[i] = dydx[i];
260 }
261 // NormaliseTangentVector( yOut ); // Not wanted
262
263 fLastStepLength = Step;
264
265 return;
266}
267
268template <class T_Equation, unsigned int N >
269inline G4double
271{
272 G4double distLine, distChord;
273 G4ThreeVector initialPoint, finalPoint, midPoint;
274
275 // Store last initial and final points (they will be overwritten in
276 // self-Stepper call!)
277 initialPoint = G4ThreeVector(fLastInitialVector[0], fLastInitialVector[1],
278 fLastInitialVector[2]);
279 finalPoint = G4ThreeVector(fLastFinalVector[0], fLastFinalVector[1],
280 fLastFinalVector[2]);
281
282 // Do half a step using StepNoErr
283
284 fAuxStepper->G4TCashKarpRKF45::Stepper(fLastInitialVector, fLastDyDx,
285 0.5 * fLastStepLength, fMidVector,
286 fMidError);
287
288 midPoint = G4ThreeVector(fMidVector[0], fMidVector[1], fMidVector[2]);
289
290 // Use stored values of Initial and Endpoint + new Midpoint to evaluate
291 // distance of Chord
292
293 if(initialPoint != finalPoint)
294 {
295 distLine = G4LineSection::Distline(midPoint, initialPoint, finalPoint);
296 distChord = distLine;
297 }
298 else
299 {
300 distChord = (midPoint - initialPoint).mag();
301 }
302 return distChord;
303}
304
305template <class T_Equation, unsigned int N >
306inline void
308 const G4double dydx[],
309 G4double Step,
310 G4double yOutput[],
311 G4double yError[])
312{
313 assert( yOutput != yInput );
314 assert( yError != yInput );
315
316 StepWithError( yInput, dydx, Step, yOutput, yError);
317}
318
319
320
321#endif /* G4TCashKARP_RKF45_hh */
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static const G4double ak2
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
G4double * fLastFinalVector
G4double * fLastInitialVector
G4TCashKarpRKF45(const G4TCashKarpRKF45 &)
G4int IntegratorOrder() const override
G4TCashKarpRKF45 * fAuxStepper
G4TCashKarpRKF45(T_Equation *EqRhs, G4bool primary=true)
G4double DistChord() const override
G4TCashKarpRKF45 & operator=(const G4TCashKarpRKF45 &)
void RightHandSideInl(const G4double y[], G4double dydx[])
void StepWithError(const G4double yInput[], const G4double dydx[], G4double Step, G4double yOut[], G4double yErr[])
T_Equation * fEquation_Rhs
virtual ~G4TCashKarpRKF45()
virtual void Stepper(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) override final