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

#include <G4RK547FEq1.hh>

+ Inheritance diagram for G4RK547FEq1:

Public Member Functions

 G4RK547FEq1 (G4EquationOfMotion *EqRhs, G4int integrationVariables=6)
 
 G4RK547FEq1 (const G4RK547FEq1 &)=delete
 
G4RK547FEq1operator= (const G4RK547FEq1 &)=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 G4RK547FEq1.hh.

Constructor & Destructor Documentation

◆ G4RK547FEq1() [1/2]

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

Definition at line 51 of file G4RK547FEq1.cc.

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

◆ G4RK547FEq1() [2/2]

G4RK547FEq1::G4RK547FEq1 ( const G4RK547FEq1 )
delete

Member Function Documentation

◆ DistChord()

G4double G4RK547FEq1::DistChord ( ) const
overridevirtual

Implements G4MagIntegratorStepper.

Definition at line 164 of file G4RK547FEq1.cc.

165{
167 makeStep(fyIn, fdydx, fhstep / 2., yMid);
168
169 const G4ThreeVector begin = makeVector(fyIn, Value3D::Position);
170 const G4ThreeVector mid = makeVector(yMid, Value3D::Position);
171 const G4ThreeVector end = makeVector(fyOut, Value3D::Position);
172
173 return G4LineSection::Distline(mid, begin, end);
174}
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 G4RK547FEq1::IntegratorOrder ( ) const
inlineoverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 67 of file G4RK547FEq1.hh.

67{ return 4; }

◆ operator=()

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

◆ Stepper() [1/2]

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

Implements G4MagIntegratorStepper.

Definition at line 132 of file G4RK547FEq1.cc.

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

◆ Stepper() [2/2]

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

Definition at line 147 of file G4RK547FEq1.cc.

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

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