Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ClassicalRK4 Class Reference

#include <G4ClassicalRK4.hh>

+ Inheritance diagram for G4ClassicalRK4:

Public Member Functions

 G4ClassicalRK4 (G4EquationOfMotion *EquationMotion, G4int numberOfVariables=6)
 
 ~G4ClassicalRK4 () override
 
 G4ClassicalRK4 (const G4ClassicalRK4 &)=delete
 
G4ClassicalRK4operator= (const G4ClassicalRK4 &)=delete
 
void DumbStepper (const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[]) override
 
G4int IntegratorOrder () const override
 
- Public Member Functions inherited from G4MagErrorStepper
 G4MagErrorStepper (G4EquationOfMotion *EqRhs, G4int numberOfVariables, G4int numStateVariables=12)
 
 ~G4MagErrorStepper () override
 
 G4MagErrorStepper (const G4MagErrorStepper &)=delete
 
G4MagErrorStepperoperator= (const G4MagErrorStepper &)=delete
 
void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
 
G4double DistChord () 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)
 

Additional Inherited Members

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

Detailed Description

Definition at line 40 of file G4ClassicalRK4.hh.

Constructor & Destructor Documentation

◆ G4ClassicalRK4() [1/2]

G4ClassicalRK4::G4ClassicalRK4 ( G4EquationOfMotion * EquationMotion,
G4int numberOfVariables = 6 )

Definition at line 38 of file G4ClassicalRK4.cc.

40 : G4MagErrorStepper(EqRhs, numberOfVariables)
41{
42 unsigned int noVariables= std::max(numberOfVariables,8); // For Time .. 7+1
43
44 dydxm = new G4double[noVariables];
45 dydxt = new G4double[noVariables];
46 yt = new G4double[noVariables];
47}
double G4double
Definition G4Types.hh:83
G4MagErrorStepper(G4EquationOfMotion *EqRhs, G4int numberOfVariables, G4int numStateVariables=12)

◆ ~G4ClassicalRK4()

G4ClassicalRK4::~G4ClassicalRK4 ( )
override

Definition at line 53 of file G4ClassicalRK4.cc.

54{
55 delete [] dydxm;
56 delete [] dydxt;
57 delete [] yt;
58}

◆ G4ClassicalRK4() [2/2]

G4ClassicalRK4::G4ClassicalRK4 ( const G4ClassicalRK4 & )
delete

Member Function Documentation

◆ DumbStepper()

void G4ClassicalRK4::DumbStepper ( const G4double yIn[],
const G4double dydx[],
G4double h,
G4double yOut[] )
overridevirtual

Implements G4MagErrorStepper.

Definition at line 71 of file G4ClassicalRK4.cc.

75{
76 const G4int nvar = GetNumberOfVariables(); // fNumberOfVariables();
77 G4int i;
78 G4double hh = h*0.5, h6 = h/6.0;
79
80 // Initialise time to t0, needed when it is not updated by the integration.
81 // [ Note: Only for time dependent fields (usually electric)
82 // is it neccessary to integrate the time.]
83 yt[7] = yIn[7];
84 yOut[7] = yIn[7];
85
86 for(i=0; i<nvar; ++i)
87 {
88 yt[i] = yIn[i] + hh*dydx[i] ; // 1st Step K1=h*dydx
89 }
90 RightHandSide(yt,dydxt) ; // 2nd Step K2=h*dydxt
91
92 for(i=0; i<nvar; ++i)
93 {
94 yt[i] = yIn[i] + hh*dydxt[i] ;
95 }
96 RightHandSide(yt,dydxm) ; // 3rd Step K3=h*dydxm
97
98 for(i=0; i<nvar; ++i)
99 {
100 yt[i] = yIn[i] + h*dydxm[i] ;
101 dydxm[i] += dydxt[i] ; // now dydxm=(K2+K3)/h
102 }
103 RightHandSide(yt,dydxt) ; // 4th Step K4=h*dydxt
104
105 for(i=0; i<nvar; ++i) // Final RK4 output
106 {
107 yOut[i] = yIn[i]+h6*(dydx[i]+dydxt[i]+2.0*dydxm[i]); //+K1/6+K4/6+(K2+K3)/3
108 }
109 if ( nvar == 12 ) { NormalisePolarizationVector ( yOut ); }
110
111} // end of DumbStepper ....................................................
int G4int
Definition G4Types.hh:85
void NormalisePolarizationVector(G4double vec[12])
G4int GetNumberOfVariables() const
void RightHandSide(const G4double y[], G4double dydx[]) const

◆ IntegratorOrder()

G4int G4ClassicalRK4::IntegratorOrder ( ) const
inlineoverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 68 of file G4ClassicalRK4.hh.

68{ return 4; }

◆ operator=()

G4ClassicalRK4 & G4ClassicalRK4::operator= ( const G4ClassicalRK4 & )
delete

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