Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TClassicalRK4.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// G4TClassicalRK4
27//
28// Class description:
29//
30// Templated version of G4ClassicalRK4.
31// Adapted from G4TClassicalRK4 class.
32
33// Author: Josh Xie (CERN, Google Summer of Code 2014), June 2014
34// Supervisors: Sandro Wenzel, John Apostolakis (CERN)
35// --------------------------------------------------------------------
36#ifndef G4TCLASSICALRK4_HH
37#define G4TCLASSICALRK4_HH
38
39#include "G4ThreeVector.hh"
41#include "G4TMagErrorStepper.hh"
42
43/**
44 * @brief G4TClassicalRK4 is a templated version of G4ClassicalRK4
45 * 4th order Runge-Kutta stepper.
46 */
47
48template <class T_Equation, unsigned int N>
49class G4TClassicalRK4 : public G4TMagErrorStepper<G4TClassicalRK4<T_Equation, N>, T_Equation, N>
50{
51 public:
52
53 static constexpr G4double IntegratorCorrection = 1. / ((1 << 4) - 1);
54
55 G4TClassicalRK4(T_Equation* EqRhs, G4int numberOfVariables = 8);
56
57 ~G4TClassicalRK4() override = default;
58
61
63 {
64 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
65 }
66
67 // A stepper that does not know about errors.
68 // It is used by the MagErrorStepper stepper.
69
70 inline void DumbStepper(const G4double yIn[],
71 const G4double dydx[],
72 G4double h,
73 G4double yOut[]);
74
75 G4int IntegratorOrder() const { return 4; }
76
77 private:
78
79 G4double dydxm[N < 8 ? 8 : N];
80 G4double dydxt[N < 8 ? 8 : N];
81 G4double yt[N < 8 ? 8 : N];
82 // scratch space - not state
83
84 T_Equation* fEquation_Rhs;
85};
86
87template <class T_Equation, unsigned int N >
89G4TClassicalRK4(T_Equation* EqRhs, G4int numberOfVariables)
90 : G4TMagErrorStepper<G4TClassicalRK4<T_Equation, N>, T_Equation, N>(
91 EqRhs, numberOfVariables > 8 ? numberOfVariables : 8 )
92 , fEquation_Rhs(EqRhs)
93{
94 // unsigned int noVariables = std::max(numberOfVariables, 8); // For Time .. 7+1
95 if( dynamic_cast<G4EquationOfMotion*>(EqRhs) == nullptr )
96 {
97 G4Exception("G4TClassicalRK4: constructor", "GeomField0001",
98 FatalException, "Equation is not an G4EquationOfMotion.");
99 }
100}
101
102template <class T_Equation, unsigned int N >
103void
105 const G4double dydx[],
106 G4double h,
107 G4double yOut[])
108// Given values for the variables y[0,..,n-1] and their derivatives
109// dydx[0,...,n-1] known at x, use the classical 4th Runge-Kutta
110// method to advance the solution over an interval h and return the
111// incremented variables as yout[0,...,n-1], which not be a distinct
112// array from y. The user supplies the routine RightHandSide(x,y,dydx),
113// which returns derivatives dydx at x. The source is routine rk4 from
114// NRC p. 712-713 .
115{
116 G4double hh = h * 0.5, h6 = h / 6.0;
117
118 // Initialise time to t0, needed when it is not updated by the integration.
119 // [ Note: Only for time dependent fields (usually electric)
120 // is it neccessary to integrate the time.]
121 yt[7] = yIn[7];
122 yOut[7] = yIn[7];
123
124 for(unsigned int i = 0; i < N; ++i)
125 {
126 yt[i] = yIn[i] + hh * dydx[i]; // 1st Step K1=h*dydx
127 }
128 this->RightHandSideInl(yt, dydxt); // 2nd Step K2=h*dydxt
129
130 for(unsigned int i = 0; i < N; ++i)
131 {
132 yt[i] = yIn[i] + hh * dydxt[i];
133 }
134 this->RightHandSideInl(yt, dydxm); // 3rd Step K3=h*dydxm
135
136 for(unsigned int i = 0; i < N; ++i)
137 {
138 yt[i] = yIn[i] + h * dydxm[i];
139 dydxm[i] += dydxt[i]; // now dydxm=(K2+K3)/h
140 }
141 this->RightHandSideInl(yt, dydxt); // 4th Step K4=h*dydxt
142
143 for(unsigned int i = 0; i < N; ++i) // Final RK4 output
144 {
145 yOut[i] = yIn[i] + h6 * (dydx[i] + dydxt[i] +
146 2.0 * dydxm[i]); //+K1/6+K4/6+(K2+K3)/3
147 }
148 if(N == 12)
149 {
150 this->NormalisePolarizationVector(yOut);
151 }
152}
153
154#endif
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4EquationOfMotion is the abstract base class for the right hand size of the equation of motion of a ...
void NormalisePolarizationVector(G4double vec[12])
G4TClassicalRK4(T_Equation *EqRhs, G4int numberOfVariables=8)
~G4TClassicalRK4() override=default
void DumbStepper(const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[])
static constexpr G4double IntegratorCorrection
void RightHandSideInl(G4double y[], G4double dydx[])
G4TClassicalRK4 & operator=(const G4TClassicalRK4 &)=delete
G4TClassicalRK4(const G4TClassicalRK4 &)=delete
G4int IntegratorOrder() const
G4TMagErrorStepper(T_Equation *EqRhs, G4int numberOfVariables, G4int numStateVariables=12)
#define N
Definition crc32.c:57