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