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

#include <G4TClassicalRK4.hh>

+ Inheritance diagram for G4TClassicalRK4< T_Equation, N >:

Public Member Functions

 G4TClassicalRK4 (T_Equation *EqRhs, G4int numberOfVariables=8)
 
virtual ~G4TClassicalRK4 ()
 
void RightHandSideInl (G4double y[], G4double dydx[])
 
void DumbStepper (const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[])
 
G4int IntegratorOrder () const
 
 G4TClassicalRK4 (const G4TClassicalRK4 &)=delete
 
G4TClassicalRK4operator= (const G4TClassicalRK4 &)=delete
 
- Public Member Functions inherited from G4TMagErrorStepper< G4TClassicalRK4< T_Equation, N >, T_Equation, N >
 G4TMagErrorStepper (T_Equation *EqRhs, G4int numberOfVariables, G4int numStateVariables=12)
 
virtual ~G4TMagErrorStepper ()
 
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
 
- 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 >

Definition at line 42 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 85 of file G4TClassicalRK4.hh.

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}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
#define N
Definition crc32.c:57

◆ ~G4TClassicalRK4()

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

Definition at line 50 of file G4TClassicalRK4.hh.

50{ ; }

◆ 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 101 of file G4TClassicalRK4.hh.

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 ....................................................
double G4double
Definition G4Types.hh:83
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

Implements G4MagIntegratorStepper.

Definition at line 68 of file G4TClassicalRK4.hh.

68{ 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 52 of file G4TClassicalRK4.hh.

54 {
55 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
56 }

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 46 of file G4TClassicalRK4.hh.


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