Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
32//
33// Created: Josh Xie (supported by Google Summer of Code 2014 )
34// Supervisors: Sandro Wenzel, John Apostolakis (CERN)
35// Adapted from G4G4TClassicalRK4 class
36// --------------------------------------------------------------------
37#include "G4ThreeVector.hh"
39#include "G4TMagErrorStepper.hh"
40
41template <class T_Equation, unsigned int N>
43 : public G4TMagErrorStepper<G4TClassicalRK4<T_Equation, N>, T_Equation, N>
44{
45 public: // with description
46 static constexpr G4double IntegratorCorrection = 1. / ((1 << 4) - 1);
47
48 G4TClassicalRK4(T_Equation* EqRhs, G4int numberOfVariables = 8);
49
50 virtual ~G4TClassicalRK4() { ; }
51
53 G4double dydx[])
54 {
55 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
56 }
57
58 // A stepper that does not know about errors.
59 // It is used by the MagErrorStepper stepper.
60
61 inline // __attribute__((always_inline))
62 void DumbStepper(const G4double yIn[],
63 const G4double dydx[],
64 G4double h,
65 G4double yOut[]);
66
67 public: // without description
68 G4int IntegratorOrder() const { return 4; }
69
72 // No copy constructor and assignment operator.
73
74 private:
75 // G4int fNumberOfVariables ; // is set default to 6 in constructor
76 G4double dydxm[N < 8 ? 8 : N];
77 G4double dydxt[N < 8 ? 8 : N];
78 G4double yt[N < 8 ? 8 : N];
79 // scratch space - not state
80
81 T_Equation* fEquation_Rhs;
82};
83
84template <class T_Equation, unsigned int N >
86G4TClassicalRK4(T_Equation* EqRhs, G4int numberOfVariables)
87 : G4TMagErrorStepper<G4TClassicalRK4<T_Equation, N>, T_Equation, N>(
88 EqRhs, numberOfVariables > 8 ? numberOfVariables : 8 )
89 , fEquation_Rhs(EqRhs)
90{
91 // unsigned int noVariables = std::max(numberOfVariables, 8); // For Time .. 7+1
92 if( dynamic_cast<G4EquationOfMotion*>(EqRhs) == nullptr )
93 {
94 G4Exception("G4TClassicalRK4: constructor", "GeomField0001",
95 FatalException, "Equation is not an G4EquationOfMotion.");
96 }
97}
98
99template <class T_Equation, unsigned int N >
100void
102 const G4double dydx[],
103 G4double h,
104 G4double yOut[])
105// Given values for the variables y[0,..,n-1] and their derivatives
106// dydx[0,...,n-1] known at x, use the classical 4th Runge-Kutta
107// method to advance the solution over an interval h and return the
108// incremented variables as yout[0,...,n-1], which not be a distinct
109// array from y. The user supplies the routine RightHandSide(x,y,dydx),
110// which returns derivatives dydx at x. The source is routine rk4 from
111// NRC p. 712-713 .
112{
113 G4double hh = h * 0.5, h6 = h / 6.0;
114
115 // Initialise time to t0, needed when it is not updated by the integration.
116 // [ Note: Only for time dependent fields (usually electric)
117 // is it neccessary to integrate the time.]
118 yt[7] = yIn[7];
119 yOut[7] = yIn[7];
120
121 for(unsigned int i = 0; i < N; ++i)
122 {
123 yt[i] = yIn[i] + hh * dydx[i]; // 1st Step K1=h*dydx
124 }
125 this->RightHandSideInl(yt, dydxt); // 2nd Step K2=h*dydxt
126
127 for(unsigned int i = 0; i < N; ++i)
128 {
129 yt[i] = yIn[i] + hh * dydxt[i];
130 }
131 this->RightHandSideInl(yt, dydxm); // 3rd Step K3=h*dydxm
132
133 for(unsigned int i = 0; i < N; ++i)
134 {
135 yt[i] = yIn[i] + h * dydxm[i];
136 dydxm[i] += dydxt[i]; // now dydxm=(K2+K3)/h
137 }
138 this->RightHandSideInl(yt, dydxt); // 4th Step K4=h*dydxt
139
140 for(unsigned int i = 0; i < N; ++i) // Final RK4 output
141 {
142 yOut[i] = yIn[i] + h6 * (dydx[i] + dydxt[i] +
143 2.0 * dydxm[i]); //+K1/6+K4/6+(K2+K3)/3
144 }
145 if(N == 12)
146 {
147 this->NormalisePolarizationVector(yOut);
148 }
149
150} // end of DumbStepper ....................................................
151
152// template <class T_Equation, unsigned int N >
153// G4TClassicalRK4<T_Equation,N>::
154
@ 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
void NormalisePolarizationVector(G4double vec[12])
G4TClassicalRK4(T_Equation *EqRhs, G4int numberOfVariables=8)
virtual ~G4TClassicalRK4()
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