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

G4RKG3_Stepper implements a Runga-Kutta integrator stepper used in Geant-3. More...

#include <G4RKG3_Stepper.hh>

Inheritance diagram for G4RKG3_Stepper:

Public Member Functions

 G4RKG3_Stepper (G4Mag_EqRhs *EqRhs)
 ~G4RKG3_Stepper () override=default
 G4RKG3_Stepper (const G4RKG3_Stepper &)=delete
G4RKG3_Stepperoperator= (const G4RKG3_Stepper &)=delete
void Stepper (const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[], G4double yErr[]) override
G4double DistChord () const override
void StepNoErr (const G4double tIn[8], const G4double dydx[6], G4double Step, G4double tOut[8], G4double B[3])
void StepWithEst (const G4double tIn[8], const G4double dydx[6], G4double Step, G4double tOut[8], G4double &alpha2, G4double &beta2, const G4double B1[3], G4double B2[3])
G4int IntegratorOrder () const override
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
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

G4RKG3_Stepper implements a Runga-Kutta integrator stepper used in Geant-3.

Definition at line 48 of file G4RKG3_Stepper.hh.

Constructor & Destructor Documentation

◆ G4RKG3_Stepper() [1/2]

G4RKG3_Stepper::G4RKG3_Stepper ( G4Mag_EqRhs * EqRhs)

Constructor for G4RKG3_Stepper.

Parameters
[in]EqRhsPointer to the provided equation of motion.

Definition at line 35 of file G4RKG3_Stepper.cc.

36 : G4MagIntegratorStepper(EqRhs,6)
37{
38}
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)

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

◆ ~G4RKG3_Stepper()

G4RKG3_Stepper::~G4RKG3_Stepper ( )
overridedefault

Default Destructor.

◆ G4RKG3_Stepper() [2/2]

G4RKG3_Stepper::G4RKG3_Stepper ( const G4RKG3_Stepper & )
delete

Copy constructor and assignment operator not allowed.

Member Function Documentation

◆ DistChord()

G4double G4RKG3_Stepper::DistChord ( ) const
overridevirtual

Returns the distance from chord line.

Implements G4MagIntegratorStepper.

Definition at line 192 of file G4RKG3_Stepper.cc.

193{
194 // Soon: must check whether h/R > 2 pi !!
195 // Method below is good only for < 2 pi
196
197 G4double distChord,distLine;
198
199 if (fyInitial != fyFinal)
200 {
201 distLine = G4LineSection::Distline(fyMidPoint,fyInitial,fyFinal);
202 distChord = distLine;
203 }
204 else
205 {
206 distChord = (fyMidPoint-fyInitial).mag();
207 }
208
209 return distChord;
210}
double G4double
Definition G4Types.hh:83
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)

◆ IntegratorOrder()

G4int G4RKG3_Stepper::IntegratorOrder ( ) const
inlineoverridevirtual

Returns the order, 4, of integration.

Implements G4MagIntegratorStepper.

Definition at line 123 of file G4RKG3_Stepper.hh.

123{ return 4; }

◆ operator=()

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

◆ StepNoErr()

void G4RKG3_Stepper::StepNoErr ( const G4double tIn[8],
const G4double dydx[6],
G4double Step,
G4double tOut[8],
G4double B[3] )

Integrator of Runge-Kutta Stepper from Geant-3 with only two field evaluation per Step. It is used in propagating the initial Step by small substeps after solution error and delta geometry considerations. B[3] is magnetic field which is passed from substep to substep.

Definition at line 126 of file G4RKG3_Stepper.cc.

132{
133
134 // Copy and edit the routine above, to delete alpha2, beta2, ...
135 //
136 G4double K1[7], K2[7], K3[7], K4[7];
137 G4double tTemp[8]={0.0}, yderiv[6]={0.0};
138
139 // Need Momentum value to give correct values to the coefficients in
140 // equation. Integration on unit velocity, but tIn[3,4,5] is momentum
141
142 G4double mom, inverse_mom;
143 const G4double c1=0.5, c2=0.125, c3=1./6.;
144
145 // Correction for momentum not a velocity
146 // Need the protection !!! must be not zero
147 //
148 mom = std::sqrt(tIn[3]*tIn[3]+tIn[4]*tIn[4]+tIn[5]*tIn[5]);
149 inverse_mom = 1./mom;
150 for(auto i=0; i<3; ++i)
151 {
152 K1[i] = Step * dydx[i+3]*inverse_mom;
153 tTemp[i] = tIn[i] + Step*(c1*tIn[i+3]*inverse_mom + c2*K1[i]) ;
154 tTemp[i+3] = tIn[i+3] + c1*K1[i]*mom ;
155 }
156
157 GetEquationOfMotion()->EvaluateRhsReturnB(tTemp,yderiv,B) ;
158
159 for(auto i=0; i<3; ++i)
160 {
161 K2[i] = Step * yderiv[i+3]*inverse_mom;
162 tTemp[i+3] = tIn[i+3] + c1*K2[i]*mom ;
163 }
164
165 // Given B, calculate yderiv !
166 //
167 GetEquationOfMotion()->EvaluateRhsGivenB(tTemp,B,yderiv) ;
168
169 for(auto i=0; i<3; ++i)
170 {
171 K3[i] = Step * yderiv[i+3]*inverse_mom;
172 tTemp[i] = tIn[i] + Step*(tIn[i+3]*inverse_mom + c1*K3[i]) ;
173 tTemp[i+3] = tIn[i+3] + K3[i]*mom ;
174 }
175
176 // Calculates y-deriv(atives) & returns B too!
177 //
178 GetEquationOfMotion()->EvaluateRhsReturnB(tTemp,yderiv,B) ;
179
180 for(auto i=0; i<3; ++i) // Output trajectory vector
181 {
182 K4[i] = Step * yderiv[i+3]*inverse_mom;
183 tOut[i] = tIn[i] + Step*(tIn[i+3]*inverse_mom+ (K1[i]+K2[i]+K3[i])*c3) ;
184 tOut[i+3] = tIn[i+3] + mom*(K1[i] + 2*K2[i] + 2*K3[i] +K4[i])*c3 ;
185 }
186 tOut[6] = tIn[6];
187 tOut[7] = tIn[7];
188}
G4double B(G4double temperature)
virtual void EvaluateRhsGivenB(const G4double y[], const G4double B[3], G4double dydx[]) const =0
void EvaluateRhsReturnB(const G4double y[], G4double dydx[], G4double Field[]) const
G4EquationOfMotion * GetEquationOfMotion()

Referenced by Stepper().

◆ Stepper()

void G4RKG3_Stepper::Stepper ( const G4double yIn[],
const G4double dydx[],
G4double h,
G4double yOut[],
G4double yErr[] )
overridevirtual

The stepper for the Runge Kutta integration. Method provided, even if less efficient. The stepsize is fixed, with the step size given by 'h'. Integrates ODE starting values yInput[0 to 6]. Outputs yOut[] and its estimated error yErr[].

Parameters
[in]yInStarting 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 40 of file G4RKG3_Stepper.cc.

45{
46 G4double B[3];
47 G4int nvar = 6 ;
48 G4double by15 = 1. / 15. ; // was 0.066666666 ;
49
50 G4double yTemp[8], dydxTemp[6], yIn[8];
51
52 // Saving yInput because yInput and yOut can be aliases for same array
53 //
54 for(G4int i=0; i<nvar; ++i)
55 {
56 yIn[i]=yInput[i];
57 }
58 yIn[6] = yInput[6];
59 yIn[7] = yInput[7];
60 G4double h = Step * 0.5;
61 hStep = Step;
62 // Do two half steps
63
64 StepNoErr(yIn, dydx,h, yTemp,B) ;
65
66 // Store Bfld for DistChord Calculation
67 //
68 for(auto i=0; i<3; ++i)
69 {
70 BfldIn[i] = B[i];
71 }
72 // RightHandSide(yTemp,dydxTemp) ;
73
74 GetEquationOfMotion()->EvaluateRhsGivenB(yTemp,B,dydxTemp) ;
75 StepNoErr(yTemp,dydxTemp,h,yOut,B);
76
77 // Store midpoint, chord calculation
78
79 fyMidPoint = G4ThreeVector(yTemp[0], yTemp[1], yTemp[2]);
80
81 // Do a full Step
82 //
83 h *= 2 ;
84 StepNoErr(yIn,dydx,h,yTemp,B);
85 for(G4int i=0; i<nvar; ++i)
86 {
87 yErr[i] = yOut[i] - yTemp[i] ;
88 yOut[i] += yErr[i]*by15 ; // Provides 5th order of accuracy
89 }
90
91 // Store values for DistChord method
92 //
93 fyInitial = G4ThreeVector( yIn[0], yIn[1], yIn[2]);
94 fpInitial = G4ThreeVector( yIn[3], yIn[4], yIn[5]);
95 fyFinal = G4ThreeVector( yOut[0], yOut[1], yOut[2]);
96}
CLHEP::Hep3Vector G4ThreeVector
int G4int
Definition G4Types.hh:85
void StepNoErr(const G4double tIn[8], const G4double dydx[6], G4double Step, G4double tOut[8], G4double B[3])

◆ StepperType()

G4StepperType G4RKG3_Stepper::StepperType ( ) const
inlineoverridevirtual

Returns the stepper type-ID, "kRKG3Stepper".

Reimplemented from G4MagIntegratorStepper.

Definition at line 128 of file G4RKG3_Stepper.hh.

128{ return kRKG3Stepper; }
@ kRKG3Stepper
G4RKG3_Stepper.

◆ StepWithEst()

void G4RKG3_Stepper::StepWithEst ( const G4double tIn[8],
const G4double dydx[6],
G4double Step,
G4double tOut[8],
G4double & alpha2,
G4double & beta2,
const G4double B1[3],
G4double B2[3] )

Integrator for Runge-Kutta from Geant-3 with evaluation of error in solution and delta geometry based on naive similarity with the case of uniform magnetic field. B1[3] is in input and is the first magnetic field values B2[3] is the output and is the final magnetic field values.

Definition at line 105 of file G4RKG3_Stepper.cc.

114{
115 G4Exception("G4RKG3_Stepper::StepWithEst()", "GeomField0001",
116 FatalException, "Method no longer used.");
117}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)

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