Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MagIntegratorStepper Class Referenceabstract

G4MagIntegratorStepper is an abstract base class for integrator of particle's equation of motion, used in tracking in space dependent magnetic field. More...

#include <G4MagIntegratorStepper.hh>

Inheritance diagram for G4MagIntegratorStepper:

Public Member Functions

 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
virtual G4StepperType StepperType () 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)

Protected Member Functions

void SetIntegrationOrder (G4int order)
void SetFSAL (G4bool flag=true)

Detailed Description

G4MagIntegratorStepper is an abstract base class for integrator of particle's equation of motion, used in tracking in space dependent magnetic field.

Definition at line 58 of file G4MagIntegratorStepper.hh.

Constructor & Destructor Documentation

◆ G4MagIntegratorStepper() [1/2]

G4MagIntegratorStepper::G4MagIntegratorStepper ( G4EquationOfMotion * Equation,
G4int numIntegrationVariables,
G4int numStateVariables = 12,
G4bool isFSAL = false )

Constructor for G4MagIntegratorStepper.

Parameters
[in]EquationPointer to the provided equation of motion.
[in]numIntegrationVariablesThe number of integration variables.
[in]numStateVariablesThe number of state variables.
[in]isFSALFlag to indicate if it is an FSAL (First Same As Last) type driver.

Definition at line 36 of file G4MagIntegratorStepper.cc.

41 : fEquation_Rhs(Equation),
42 fNoIntegrationVariables(num_integration_vars),
43 fNoStateVariables(std::max(num_state_vars,8)),
44 fIsFSAL(isFSAL)
45{
46 if( Equation == nullptr )
47 {
48 G4Exception( "G4MagIntegratorStepper::G4MagIntegratorStepper", "GeomField0003",
49 FatalErrorInArgument, "Must have non-null equation." );
50 }
51}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)

Referenced by G4BogackiShampine23::G4BogackiShampine23(), G4BogackiShampine45::G4BogackiShampine45(), G4CashKarpRKF45::G4CashKarpRKF45(), G4DoLoMcPriRK34::G4DoLoMcPriRK34(), G4DormandPrince745::G4DormandPrince745(), G4DormandPrinceRK56::G4DormandPrinceRK56(), G4DormandPrinceRK78::G4DormandPrinceRK78(), G4MagHelicalStepper::G4MagHelicalStepper(), G4MagIntegratorStepper(), G4NystromRK4::G4NystromRK4(), G4QSStepper::G4QSStepper(), G4RK547FEq1::G4RK547FEq1(), G4RK547FEq2::G4RK547FEq2(), G4RK547FEq3::G4RK547FEq3(), G4RKG3_Stepper::G4RKG3_Stepper(), G4TCashKarpRKF45< T_Equation, N >::G4TCashKarpRKF45(), G4TDormandPrince45< T_Equation, N >::G4TDormandPrince45(), G4TMagErrorStepper< T_Stepper, T_Equation, N >::G4TMagErrorStepper(), G4TsitourasRK45::G4TsitourasRK45(), operator=(), and G4HelixMixedStepper::SetupStepper().

◆ ~G4MagIntegratorStepper()

virtual G4MagIntegratorStepper::~G4MagIntegratorStepper ( )
virtualdefault

Default virtual Destructor.

◆ G4MagIntegratorStepper() [2/2]

G4MagIntegratorStepper::G4MagIntegratorStepper ( const G4MagIntegratorStepper & )
delete

Copy constructor and assignment operator not allowed.

Member Function Documentation

◆ DistChord()

◆ GetEquationOfMotion() [1/2]

G4EquationOfMotion * G4MagIntegratorStepper::GetEquationOfMotion ( )
inline

Methods returning the pointer to the associated equation of motion. As some steppers (e.g. RKG3) require other methods of Eq_Rhs this function allows for access to them.

Referenced by G4DormandPrince745::GetSpecificEquation(), G4NystromRK4::SetDistanceForConstantField(), G4RKG3_Stepper::StepNoErr(), and G4RKG3_Stepper::Stepper().

◆ GetEquationOfMotion() [2/2]

const G4EquationOfMotion * G4MagIntegratorStepper::GetEquationOfMotion ( ) const
inline

◆ GetfNoRHSCalls()

unsigned long G4MagIntegratorStepper::GetfNoRHSCalls ( )
inline

Methods for counting/resetting the number of calls to RHS method(s).

◆ GetNumberOfStateVariables()

◆ GetNumberOfVariables()

◆ IntegrationOrder()

G4int G4MagIntegratorStepper::IntegrationOrder ( )
inline

Replacement method - using new data member.

◆ IntegratorOrder()

◆ IsFSAL()

G4bool G4MagIntegratorStepper::IsFSAL ( ) const
inline

Returns true if the stepper is of FSAL (First Same As Last) type.

◆ isQSS()

G4bool G4MagIntegratorStepper::isQSS ( ) const
inline

Returns true if the stepper is of QSS (Quantum State Simulation) type.

Referenced by G4ChordFinder::G4ChordFinder().

◆ NormalisePolarizationVector()

void G4MagIntegratorStepper::NormalisePolarizationVector ( G4double vec[12])
inline

◆ NormaliseTangentVector()

void G4MagIntegratorStepper::NormaliseTangentVector ( G4double vec[6])
inline

Simple utility function to (re)normalise 'unit velocity' vector.

◆ operator=()

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

◆ ResetfNORHSCalls()

void G4MagIntegratorStepper::ResetfNORHSCalls ( )
inline

◆ RightHandSide() [1/2]

◆ RightHandSide() [2/2]

void G4MagIntegratorStepper::RightHandSide ( const G4double y[],
G4double dydx[],
G4double field[] ) const
inline

Calculates 'dydx' and 'field' at point 'y'.

◆ SetEquationOfMotion()

void G4MagIntegratorStepper::SetEquationOfMotion ( G4EquationOfMotion * newEquation)
inline

Setter for the equation of motion.

◆ SetFSAL()

void G4MagIntegratorStepper::SetFSAL ( G4bool flag = true)
inlineprotected

◆ SetIntegrationOrder()

void G4MagIntegratorStepper::SetIntegrationOrder ( G4int order)
inlineprotected

Setters for the integration order and FSAL type.

Referenced by G4BogackiShampine23::G4BogackiShampine23().

◆ SetIsQSS()

void G4MagIntegratorStepper::SetIsQSS ( G4bool val)
inline

◆ Stepper()

virtual void G4MagIntegratorStepper::Stepper ( const G4double y[],
const G4double dydx[],
G4double h,
G4double yout[],
G4double yerr[] )
pure virtual

◆ StepperType()


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