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

#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 ()
 
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
 
- 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

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

Definition at line 46 of file G4TMagErrorStepper.hh.

Constructor & Destructor Documentation

◆ G4TMagErrorStepper()

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

51 : G4MagIntegratorStepper(EqRhs, numberOfVariables, numStateVariables)
52 , fEquation_Rhs(EqRhs)
53 {
54 // G4int nvar = std::max(this->GetNumberOfVariables(), 8);
55 }

◆ ~G4TMagErrorStepper()

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

Definition at line 57 of file G4TMagErrorStepper.hh.

57{ ; }

Member Function Documentation

◆ DistChord()

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

Implements G4MagIntegratorStepper.

Definition at line 153 of file G4TMagErrorStepper.hh.

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

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

60 {
61 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
62 }

Referenced by G4TMagErrorStepper< T_Stepper, T_Equation, N >::Stepper().

◆ Stepper()

template<class T_Stepper , class T_Equation , unsigned int N>
void G4TMagErrorStepper< T_Stepper, T_Equation, N >::Stepper ( const G4double  yInput[],
const G4double  dydx[],
G4double  hstep,
G4double  yOutput[],
G4double  yError[] 
)
inlinefinaloverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 94 of file G4TMagErrorStepper.hh.

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

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