Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TMagErrorStepper< T_Stepper, T_Equation, N > Class Template Reference

G4TMagErrorStepper is a templated version of G4MagErrorStepper. More...

#include <G4TMagErrorStepper.hh>

Inheritance diagram for G4TMagErrorStepper< T_Stepper, T_Equation, N >:

Public Member Functions

 G4TMagErrorStepper (T_Equation *EqRhs, G4int numberOfVariables, G4int numStateVariables=12)
virtual ~G4TMagErrorStepper ()=default
 G4TMagErrorStepper (const G4TMagErrorStepper &)=delete
G4TMagErrorStepperoperator= (const G4TMagErrorStepper &)=delete
void RightHandSide (G4double y[], G4double dydx[])
void Stepper (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) override final
G4double DistChord () const override final
G4StepperType StepperType () 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
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

template<class T_Stepper, class T_Equation, unsigned int N>
class G4TMagErrorStepper< T_Stepper, T_Equation, N >

G4TMagErrorStepper is a templated version of G4MagErrorStepper.

Definition at line 49 of file G4TMagErrorStepper.hh.

Constructor & Destructor Documentation

◆ G4TMagErrorStepper() [1/2]

template<class T_Stepper, class T_Equation, unsigned int N>
G4TMagErrorStepper< T_Stepper, T_Equation, N >::G4TMagErrorStepper ( T_Equation * EqRhs,
G4int numberOfVariables,
G4int numStateVariables = 12 )
inline

Definition at line 53 of file G4TMagErrorStepper.hh.

56 , fEquation_Rhs(EqRhs) { ; }
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)
G4TMagErrorStepper is a templated version of G4MagErrorStepper.

Referenced by G4TMagErrorStepper(), and operator=().

◆ ~G4TMagErrorStepper()

template<class T_Stepper, class T_Equation, unsigned int N>
virtual G4TMagErrorStepper< T_Stepper, T_Equation, N >::~G4TMagErrorStepper ( )
virtualdefault

◆ G4TMagErrorStepper() [2/2]

template<class T_Stepper, class T_Equation, unsigned int N>
G4TMagErrorStepper< T_Stepper, T_Equation, N >::G4TMagErrorStepper ( const G4TMagErrorStepper< T_Stepper, T_Equation, N > & )
delete

Member Function Documentation

◆ DistChord()

template<class T_Stepper, class T_Equation, unsigned int N>
G4double G4TMagErrorStepper< T_Stepper, T_Equation, N >::DistChord ( ) const
inlinefinaloverridevirtual

Estimates the maximum distance of a chord from the true path over the segment last integrated.

Implements G4MagIntegratorStepper.

Definition at line 154 of file G4TMagErrorStepper.hh.

155{
156 // Estimate the maximum distance from the curve to the chord
157 //
158 // We estimate this using the distance of the midpoint to
159 // chord (the line between
160 //
161 // Method below is good only for angle deviations < 2 pi,
162 // This restriction should not a problem for the Runge cutta methods,
163 // which generally cannot integrate accurately for large angle deviations.
165
166 if(fInitialPoint != fFinalPoint)
167 {
168 distLine = G4LineSection::Distline(fMidPoint, fInitialPoint, fFinalPoint);
169 // This is a class method that gives distance of Mid
170 // from the Chord between the Initial and Final points.
171
173 }
174 else
175 {
176 distChord = (fMidPoint - fInitialPoint).mag();
177 }
178
179 return distChord;
180}
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)

◆ operator=()

template<class T_Stepper, class T_Equation, unsigned int N>
G4TMagErrorStepper & G4TMagErrorStepper< T_Stepper, T_Equation, N >::operator= ( const G4TMagErrorStepper< T_Stepper, T_Equation, N > & )
delete

◆ RightHandSide()

template<class T_Stepper, class T_Equation, unsigned int N>
void G4TMagErrorStepper< T_Stepper, T_Equation, N >::RightHandSide ( G4double y[],
G4double dydx[] )
inline

Definition at line 63 of file G4TMagErrorStepper.hh.

64 {
65 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
66 }

Referenced by Stepper().

◆ Stepper()

template<class T_Stepper, class T_Equation, unsigned int N>
void G4TMagErrorStepper< T_Stepper, T_Equation, N >::Stepper ( const G4double y[],
const G4double dydx[],
G4double h,
G4double yout[],
G4double yerr[] )
inlinefinaloverridevirtual

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 95 of file G4TMagErrorStepper.hh.

105{
106 const unsigned int maxvar = GetNumberOfStateVariables();
107
108 // Saving yInput because yInput and yOutput can be aliases for same array
109 for(unsigned int i = 0; i < N; ++i)
110 yInitial[i] = yInput[i];
111 yInitial[7] =
112 yInput[7]; // Copy the time in case ... even if not really needed
113 yMiddle[7] = yInput[7]; // Copy the time from initial value
114 yOneStep[7] = yInput[7]; // As it contributes to final value of yOutput ?
115 // yOutput[7] = yInput[7]; // -> dumb stepper does it too for RK4
116 for(unsigned int i = N; i < maxvar; ++i)
117 yOutput[i] = yInput[i];
118
119 G4double halfStep = hstep * 0.5;
120
121 // Do two half steps
122
123 static_cast<T_Stepper*>(this)->DumbStepper(yInitial, dydx, halfStep,
124 yMiddle);
125 this->RightHandSide(yMiddle, dydxMid);
126 static_cast<T_Stepper*>(this)->DumbStepper(yMiddle, dydxMid, halfStep,
127 yOutput);
128
129 // Store midpoint, chord calculation
130
131 fMidPoint = G4ThreeVector(yMiddle[0], yMiddle[1], yMiddle[2]);
132
133 // Do a full Step
134 static_cast<T_Stepper*>(this)->DumbStepper(yInitial, dydx, hstep, yOneStep);
135 for(unsigned int i = 0; i < N; ++i)
136 {
137 yError[i] = yOutput[i] - yOneStep[i];
138 yOutput[i] +=
139 yError[i] *
140 T_Stepper::IntegratorCorrection; // Provides accuracy increased
141 // by 1 order via the
142 // Richardson Extrapolation
143 }
144
145 fInitialPoint = G4ThreeVector(yInitial[0], yInitial[1], yInitial[2]);
146 fFinalPoint = G4ThreeVector(yOutput[0], yOutput[1], yOutput[2]);
147
148 return;
149}
CLHEP::Hep3Vector G4ThreeVector
G4int GetNumberOfStateVariables() const
void RightHandSide(G4double y[], G4double dydx[])

◆ StepperType()

template<class T_Stepper, class T_Equation, unsigned int N>
G4StepperType G4TMagErrorStepper< T_Stepper, T_Equation, N >::StepperType ( ) const
inlineoverridevirtual

Returns the stepper type ID ('kUserStepper'). This function should be overriden in derived classes.

Reimplemented from G4MagIntegratorStepper.

Definition at line 72 of file G4TMagErrorStepper.hh.

72{ return kTMagErrorStepper; }

Referenced by StepperType().


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