Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TClassicalRK4< T_Equation, N > Class Template Reference

G4TClassicalRK4 is a templated version of G4ClassicalRK4 4th order Runge-Kutta stepper. More...

#include <G4TClassicalRK4.hh>

Inheritance diagram for G4TClassicalRK4< T_Equation, N >:

Public Member Functions

 G4TClassicalRK4 (T_Equation *EqRhs, G4int numberOfVariables=8)
 ~G4TClassicalRK4 () override=default
 G4TClassicalRK4 (const G4TClassicalRK4 &)=delete
G4TClassicalRK4operator= (const G4TClassicalRK4 &)=delete
void RightHandSideInl (G4double y[], G4double dydx[])
void DumbStepper (const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[])
G4int IntegratorOrder () const
Public Member Functions inherited from G4TMagErrorStepper< G4TClassicalRK4< T_Equation, N >, T_Equation, N >
 G4TMagErrorStepper (T_Equation *EqRhs, G4int numberOfVariables, G4int numStateVariables=12)
virtual ~G4TMagErrorStepper ()=default
G4TMagErrorStepperoperator= (const G4TMagErrorStepper &)=delete
void RightHandSide (G4double y[], G4double dydx[])
void Stepper (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) override final
G4double DistChord () const override final
G4StepperType StepperType () const override
Public Member Functions inherited from G4MagIntegratorStepper
 G4MagIntegratorStepper (G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)
virtual ~G4MagIntegratorStepper ()=default
 G4MagIntegratorStepper (const G4MagIntegratorStepper &)=delete
G4MagIntegratorStepperoperator= (const G4MagIntegratorStepper &)=delete
void NormaliseTangentVector (G4double vec[6])
void NormalisePolarizationVector (G4double vec[12])
void RightHandSide (const G4double y[], G4double dydx[]) const
void RightHandSide (const G4double y[], G4double dydx[], G4double field[]) const
G4int GetNumberOfVariables () const
G4int GetNumberOfStateVariables () const
G4int IntegrationOrder ()
G4EquationOfMotionGetEquationOfMotion ()
const G4EquationOfMotionGetEquationOfMotion () const
void SetEquationOfMotion (G4EquationOfMotion *newEquation)
unsigned long GetfNoRHSCalls ()
void ResetfNORHSCalls ()
G4bool IsFSAL () const
G4bool isQSS () const
void SetIsQSS (G4bool val)

Static Public Attributes

static constexpr G4double IntegratorCorrection = 1. / ((1 << 4) - 1)

Additional Inherited Members

Protected Member Functions inherited from G4MagIntegratorStepper
void SetIntegrationOrder (G4int order)
void SetFSAL (G4bool flag=true)

Detailed Description

template<class T_Equation, unsigned int N>
class G4TClassicalRK4< T_Equation, N >

G4TClassicalRK4 is a templated version of G4ClassicalRK4 4th order Runge-Kutta stepper.

Definition at line 49 of file G4TClassicalRK4.hh.

Constructor & Destructor Documentation

◆ G4TClassicalRK4() [1/2]

template<class T_Equation, unsigned int N>
G4TClassicalRK4< T_Equation, N >::G4TClassicalRK4 ( T_Equation * EqRhs,
G4int numberOfVariables = 8 )

Definition at line 88 of file G4TClassicalRK4.hh.

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}
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4TClassicalRK4 is a templated version of G4ClassicalRK4 4th order Runge-Kutta stepper.
G4TMagErrorStepper(T_Equation *EqRhs, G4int numberOfVariables, G4int numStateVariables=12)

Referenced by G4TClassicalRK4(), G4TClassicalRK4(), and operator=().

◆ ~G4TClassicalRK4()

template<class T_Equation, unsigned int N>
G4TClassicalRK4< T_Equation, N >::~G4TClassicalRK4 ( )
overridedefault

◆ G4TClassicalRK4() [2/2]

template<class T_Equation, unsigned int N>
G4TClassicalRK4< T_Equation, N >::G4TClassicalRK4 ( const G4TClassicalRK4< T_Equation, N > & )
delete

Member Function Documentation

◆ DumbStepper()

template<class T_Equation, unsigned int N>
void G4TClassicalRK4< T_Equation, N >::DumbStepper ( const G4double yIn[],
const G4double dydx[],
G4double h,
G4double yOut[] )
inline

Definition at line 104 of file G4TClassicalRK4.hh.

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 {
151 }
152}
void NormalisePolarizationVector(G4double vec[12])
void RightHandSideInl(G4double y[], G4double dydx[])

◆ IntegratorOrder()

template<class T_Equation, unsigned int N>
G4int G4TClassicalRK4< T_Equation, N >::IntegratorOrder ( ) const
inlinevirtual

Returns the order of the integrator, i.e. its error behaviour is of the order O(h^order).

Implements G4MagIntegratorStepper.

Definition at line 75 of file G4TClassicalRK4.hh.

75{ return 4; }

◆ operator=()

template<class T_Equation, unsigned int N>
G4TClassicalRK4 & G4TClassicalRK4< T_Equation, N >::operator= ( const G4TClassicalRK4< T_Equation, N > & )
delete

◆ RightHandSideInl()

template<class T_Equation, unsigned int N>
void G4TClassicalRK4< T_Equation, N >::RightHandSideInl ( G4double y[],
G4double dydx[] )
inline

Definition at line 62 of file G4TClassicalRK4.hh.

63 {
64 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
65 }

Referenced by DumbStepper().

Member Data Documentation

◆ IntegratorCorrection

template<class T_Equation, unsigned int N>
G4double G4TClassicalRK4< T_Equation, N >::IntegratorCorrection = 1. / ((1 << 4) - 1)
staticconstexpr

Definition at line 53 of file G4TClassicalRK4.hh.


The documentation for this class was generated from the following file: