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

#include <G4RKG3_Stepper.hh>

+ Inheritance diagram for G4RKG3_Stepper:

Public Member Functions

 G4RKG3_Stepper (G4Mag_EqRhs *EqRhs)
 
 ~G4RKG3_Stepper () override
 
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
 
- 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

Definition at line 43 of file G4RKG3_Stepper.hh.

Constructor & Destructor Documentation

◆ G4RKG3_Stepper()

G4RKG3_Stepper::G4RKG3_Stepper ( G4Mag_EqRhs * EqRhs)

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)

◆ ~G4RKG3_Stepper()

G4RKG3_Stepper::~G4RKG3_Stepper ( )
overridedefault

Member Function Documentation

◆ DistChord()

G4double G4RKG3_Stepper::DistChord ( ) const
overridevirtual

Implements G4MagIntegratorStepper.

Definition at line 194 of file G4RKG3_Stepper.cc.

195{
196 // Soon: must check whether h/R > 2 pi !!
197 // Method below is good only for < 2 pi
198
199 G4double distChord,distLine;
200
201 if (fyInitial != fyFinal)
202 {
203 distLine = G4LineSection::Distline(fyMidPoint,fyInitial,fyFinal);
204 distChord = distLine;
205 }
206 else
207 {
208 distChord = (fyMidPoint-fyInitial).mag();
209 }
210
211 return distChord;
212}
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

Implements G4MagIntegratorStepper.

Definition at line 86 of file G4RKG3_Stepper.hh.

86{ return 4; }

◆ StepNoErr()

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

Definition at line 128 of file G4RKG3_Stepper.cc.

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

Implements G4MagIntegratorStepper.

Definition at line 42 of file G4RKG3_Stepper.cc.

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

◆ 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] )

Definition at line 107 of file G4RKG3_Stepper.cc.

116{
117 G4Exception("G4RKG3_Stepper::StepWithEst()", "GeomField0001",
118 FatalException, "Method no longer used.");
119}
@ 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: