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

#include <G4RK547FEq3.hh>

+ Inheritance diagram for G4RK547FEq3:

Public Member Functions

 G4RK547FEq3 (G4EquationOfMotion *EqRhs, G4int integrationVariables=6)
 
 G4RK547FEq3 (const G4RK547FEq3 &)=delete
 
G4RK547FEq3operator= (const G4RK547FEq3 &)=delete
 
virtual void Stepper (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) override
 
void Stepper (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[], G4double dydxOutput[])
 
virtual G4double DistChord () const override
 
virtual G4int IntegratorOrder () 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
 
virtual void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])=0
 
virtual G4double DistChord () const =0
 
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
 
virtual G4int IntegratorOrder () const =0
 
G4int IntegrationOrder ()
 
G4EquationOfMotionGetEquationOfMotion ()
 
const G4EquationOfMotionGetEquationOfMotion () const
 
void SetEquationOfMotion (G4EquationOfMotion *newEquation)
 
unsigned long GetfNoRHSCalls ()
 
void ResetfNORHSCalls ()
 
G4bool IsFSAL () const
 

Additional Inherited Members

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

Detailed Description

Definition at line 44 of file G4RK547FEq3.hh.

Constructor & Destructor Documentation

◆ G4RK547FEq3() [1/2]

G4RK547FEq3::G4RK547FEq3 ( G4EquationOfMotion EqRhs,
G4int  integrationVariables = 6 
)

Definition at line 51 of file G4RK547FEq3.cc.

52 : G4MagIntegratorStepper(EqRhs, integrationVariables)
53{
54}

◆ G4RK547FEq3() [2/2]

G4RK547FEq3::G4RK547FEq3 ( const G4RK547FEq3 )
delete

Member Function Documentation

◆ DistChord()

G4double G4RK547FEq3::DistChord ( ) const
overridevirtual

Implements G4MagIntegratorStepper.

Definition at line 167 of file G4RK547FEq3.cc.

168{
170 makeStep(fyIn, fdydx, fhstep / 2., yMid);
171
172 const G4ThreeVector begin = makeVector(fyIn, Value3D::Position);
173 const G4ThreeVector mid = makeVector(yMid, Value3D::Position);
174 const G4ThreeVector end = makeVector(fyOut, Value3D::Position);
175
176 return G4LineSection::Distline(mid, begin, end);
177}
double G4double
Definition: G4Types.hh:83
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
G4ThreeVector makeVector(const ArrayType &array, Value3D value)

◆ IntegratorOrder()

virtual G4int G4RK547FEq3::IntegratorOrder ( ) const
inlineoverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 67 of file G4RK547FEq3.hh.

67{ return 4; }

◆ operator=()

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

◆ Stepper() [1/2]

void G4RK547FEq3::Stepper ( const G4double  yInput[],
const G4double  dydx[],
G4double  hstep,
G4double  yOutput[],
G4double  yError[] 
)
overridevirtual

Implements G4MagIntegratorStepper.

Definition at line 135 of file G4RK547FEq3.cc.

140{
141 copy(fyIn, yInput);
142 copy(fdydx, dydx);
143 fhstep = hstep;
144
145 makeStep(fyIn, fdydx, fhstep, fyOut, fdydxOut, yError);
146
147 copy(yOutput, fyOut);
148}
void copy(G4double dst[], const G4double src[], std::size_t size=G4FieldTrack::ncompSVEC)
Definition: G4FieldUtils.cc:98

◆ Stepper() [2/2]

void G4RK547FEq3::Stepper ( const G4double  yInput[],
const G4double  dydx[],
G4double  hstep,
G4double  yOutput[],
G4double  yError[],
G4double  dydxOutput[] 
)

Definition at line 150 of file G4RK547FEq3.cc.

156{
157 copy(fyIn, yInput);
158 copy(fdydx, dydx);
159 fhstep = hstep;
160
161 makeStep(fyIn, fdydx, fhstep, fyOut, fdydxOut, yError);
162
163 copy(yOutput, fyOut);
164 copy(dydxOutput, fdydxOut);
165}

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