Geant4 9.6.0
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 ()
 
void Stepper (const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[], G4double yErr[])
 
G4double DistChord () const
 
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
 
- Public Member Functions inherited from G4MagIntegratorStepper
 G4MagIntegratorStepper (G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12)
 
virtual ~G4MagIntegratorStepper ()
 
virtual void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])=0
 
virtual G4double DistChord () const =0
 
virtual void ComputeRightHandSide (const G4double y[], G4double dydx[])
 
void NormaliseTangentVector (G4double vec[6])
 
void NormalisePolarizationVector (G4double vec[12])
 
void RightHandSide (const double y[], double dydx[])
 
G4int GetNumberOfVariables () const
 
G4int GetNumberOfStateVariables () const
 
virtual G4int IntegratorOrder () const =0
 
G4EquationOfMotionGetEquationOfMotion ()
 
void SetEquationOfMotion (G4EquationOfMotion *newEquation)
 

Detailed Description

Definition at line 51 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), hStep(0.), fPtrMagEqOfMot(EqRhs)
37{
38}

◆ ~G4RKG3_Stepper()

G4RKG3_Stepper::~G4RKG3_Stepper ( )

Definition at line 40 of file G4RKG3_Stepper.cc.

41{
42}

Member Function Documentation

◆ DistChord()

G4double G4RKG3_Stepper::DistChord ( ) const
virtual

Implements G4MagIntegratorStepper.

Definition at line 197 of file G4RKG3_Stepper.cc.

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

◆ IntegratorOrder()

G4int G4RKG3_Stepper::IntegratorOrder ( ) const
inlinevirtual

Implements G4MagIntegratorStepper.

Definition at line 95 of file G4RKG3_Stepper.hh.

95{ 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 126 of file G4RKG3_Stepper.cc.

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

Implements G4MagIntegratorStepper.

Definition at line 44 of file G4RKG3_Stepper.cc.

49{
50 G4double B[3];
51 G4int nvar = 6 ;
52 G4int i;
53 G4double by15 = 1. / 15. ; // was 0.066666666 ;
54
55 G4double yTemp[8], dydxTemp[6], yIn[8] ;
56 // Saving yInput because yInput and yOut can be aliases for same array
57 for(i=0;i<nvar;i++) yIn[i]=yInput[i];
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 for(i=0;i<3;i++)BfldIn[i]=B[i];
68
69 // RightHandSide(yTemp,dydxTemp) ;
70
71 GetEquationOfMotion()->EvaluateRhsGivenB(yTemp,B,dydxTemp) ;
72 StepNoErr(yTemp,dydxTemp,h,yOut,B);
73
74 // Store midpoint, chord calculation
75
76 fyMidPoint = G4ThreeVector( yTemp[0], yTemp[1], yTemp[2]);
77
78 // Do a full Step
79
80 h *= 2 ;
81 StepNoErr(yIn,dydx,h,yTemp,B);
82 for(i=0;i<nvar;i++)
83 {
84 yErr[i] = yOut[i] - yTemp[i] ;
85 yOut[i] += yErr[i]*by15 ; // Provides 5th order of accuracy
86 }
87
88 //Store values for DistChord method
89
90 fyInitial = G4ThreeVector( yIn[0], yIn[1], yIn[2]);
91 fpInitial = G4ThreeVector( yIn[3], yIn[4], yIn[5]);
92 fyFinal = G4ThreeVector( yOut[0], yOut[1], yOut[2]);
93
94 // NormaliseTangentVector( yOut ); // Deleted
95}
CLHEP::Hep3Vector G4ThreeVector
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 104 of file G4RKG3_Stepper.cc.

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

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