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

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

#include <G4MagErrorStepper.hh>

Inheritance diagram for G4MagErrorStepper:

Public Member Functions

 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
virtual void DumbStepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[])=0
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
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)

Additional Inherited Members

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

Detailed Description

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

Definition at line 49 of file G4MagErrorStepper.hh.

Constructor & Destructor Documentation

◆ G4MagErrorStepper() [1/2]

G4MagErrorStepper::G4MagErrorStepper ( G4EquationOfMotion * EqRhs,
G4int numberOfVariables,
G4int numStateVariables = 12 )

Constructor for G4MagErrorStepper.

Parameters
[in]EqRhsPointer to the provided equation of motion.
[in]numberOfVariablesThe number of integration variables.
[in]numberOfVariablesThe number of state variables.

Referenced by G4ClassicalRK4::G4ClassicalRK4(), G4ConstRK4::G4ConstRK4(), G4ExplicitEuler::G4ExplicitEuler(), G4ImplicitEuler::G4ImplicitEuler(), G4MagErrorStepper(), G4SimpleHeum::G4SimpleHeum(), G4SimpleRunge::G4SimpleRunge(), and operator=().

◆ ~G4MagErrorStepper()

G4MagErrorStepper::~G4MagErrorStepper ( )
override

Destructor.

Definition at line 34 of file G4MagErrorStepper.cc.

35{
36 delete [] yMiddle;
37 delete [] dydxMid;
38 delete [] yInitial;
39 delete [] yOneStep;
40}

◆ G4MagErrorStepper() [2/2]

G4MagErrorStepper::G4MagErrorStepper ( const G4MagErrorStepper & )
delete

Copy constructor and assignment operator not allowed.

Member Function Documentation

◆ DistChord()

G4double G4MagErrorStepper::DistChord ( ) const
overridevirtual

Estimates the maximum distance of curved solution and chord.

Implements G4MagIntegratorStepper.

Definition at line 100 of file G4MagErrorStepper.cc.

101{
102 // Estimate the maximum distance from the curve to the chord
103 //
104 // We estimate this using the distance of the midpoint to
105 // chord (the line between
106 //
107 // Method below is good only for angle deviations < 2 pi,
108 // This restriction should not a problem for the Runge cutta methods,
109 // which generally cannot integrate accurately for large angle deviations.
110
111 G4double distLine, distChord;
112
113 if (fInitialPoint != fFinalPoint)
114 {
115 distLine = G4LineSection::Distline(fMidPoint, fInitialPoint, fFinalPoint);
116 // This is a class method that gives distance of Mid
117 // from the Chord between the Initial and Final points.
118
119 distChord = distLine;
120 }
121 else
122 {
123 distChord = (fMidPoint-fInitialPoint).mag();
124 }
125
126 return distChord;
127}
double G4double
Definition G4Types.hh:83
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)

◆ DumbStepper()

virtual void G4MagErrorStepper::DumbStepper ( const G4double y[],
const G4double dydx[],
G4double h,
G4double yout[] )
pure virtual

Same as Stepper() function above, but should perform a 'dump' step without error calculation. To be implemented in concrete derived classes.

Implemented in G4ClassicalRK4, G4ConstRK4, G4ExplicitEuler, G4ImplicitEuler, G4SimpleHeum, and G4SimpleRunge.

Referenced by Stepper().

◆ operator=()

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

◆ Stepper()

void G4MagErrorStepper::Stepper ( const G4double y[],
const G4double dydx[],
G4double h,
G4double yout[],
G4double yerr[] )
overridevirtual

The stepper for the Runge Kutta integration. The stepsize is fixed, with the step size given by 'h'. Integrates ODE starting values y[0 to 6]. Outputs yout[] and its estimated error yerr[].

Parameters
[in]yStarting values array of integration variables.
[in]dydxDerivatives array.
[in]hThe given step size.
[out]youtIntegration output.
[out]yerrThe estimated error.

Implements G4MagIntegratorStepper.

Definition at line 42 of file G4MagErrorStepper.cc.

47{
48 const G4int nvar = GetNumberOfVariables();
49 const G4int maxvar = GetNumberOfStateVariables();
50
51 // correction for Richardson Extrapolation.
52 //
53 G4double correction = 1. / ( (1 << IntegratorOrder()) -1 );
54
55 // Saving yInput because yInput and yOutput can be aliases for same array
56 //
57 for(G4int i=0; i<nvar; ++i)
58 {
59 yInitial[i]=yInput[i];
60 }
61 yInitial[7] = yInput[7]; // Copy the time in case...even if not really needed
62 yMiddle[7] = yInput[7]; // Copy the time from initial value
63 yOneStep[7] = yInput[7]; // As it contributes to final value of yOutput ?
64 // yOutput[7] = yInput[7]; // -> dumb stepper does it too for RK4
65
66 for(G4int i=nvar; i<maxvar; ++i)
67 {
68 yOutput[i]=yInput[i];
69 }
70 // yError[7] = 0.0;
71
72 G4double halfStep = hstep * 0.5;
73
74 // Do two half steps
75 //
76 DumbStepper (yInitial, dydx, halfStep, yMiddle);
77 RightHandSide(yMiddle, dydxMid);
78 DumbStepper (yMiddle, dydxMid, halfStep, yOutput);
79
80 // Store midpoint, chord calculation
81 //
82 fMidPoint = G4ThreeVector( yMiddle[0], yMiddle[1], yMiddle[2]);
83
84 // Do a full Step
85 //
86 DumbStepper(yInitial, dydx, hstep, yOneStep);
87 for(G4int i=0; i<nvar; ++i)
88 {
89 yError [i] = yOutput[i] - yOneStep[i] ;
90 yOutput[i] += yError[i]*correction ;
91 // Provides accuracy increased by 1 order via Richardson Extrapolation
92 }
93
94 fInitialPoint = G4ThreeVector( yInitial[0], yInitial[1], yInitial[2]);
95 fFinalPoint = G4ThreeVector( yOutput[0], yOutput[1], yOutput[2]);
96
97 return;
98}
CLHEP::Hep3Vector G4ThreeVector
int G4int
Definition G4Types.hh:85
virtual void DumbStepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[])=0
G4int GetNumberOfVariables() const
void RightHandSide(const G4double y[], G4double dydx[]) const
G4int GetNumberOfStateVariables() const
virtual G4int IntegratorOrder() const =0

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