Geant4-11
Public Member Functions | Protected Member Functions | Private Attributes
G4RKG3_Stepper Class Reference

#include <G4RKG3_Stepper.hh>

Inheritance diagram for G4RKG3_Stepper:
G4MagIntegratorStepper

Public Member Functions

G4double DistChord () const
 
 G4RKG3_Stepper (G4Mag_EqRhs *EqRhs)
 
G4EquationOfMotionGetEquationOfMotion ()
 
const G4EquationOfMotionGetEquationOfMotion () const
 
unsigned long GetfNoRHSCalls ()
 
G4int GetNumberOfStateVariables () const
 
G4int GetNumberOfVariables () const
 
G4int IntegrationOrder ()
 
G4int IntegratorOrder () const
 
G4bool IsFSAL () const
 
void NormalisePolarizationVector (G4double vec[12])
 
void NormaliseTangentVector (G4double vec[6])
 
void ResetfNORHSCalls ()
 
void RightHandSide (const G4double y[], G4double dydx[]) const
 
void RightHandSide (const G4double y[], G4double dydx[], G4double field[]) const
 
void SetEquationOfMotion (G4EquationOfMotion *newEquation)
 
void StepNoErr (const G4double tIn[8], const G4double dydx[6], G4double Step, G4double tOut[8], G4double B[3])
 
void Stepper (const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[], G4double yErr[])
 
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])
 
 ~G4RKG3_Stepper ()
 

Protected Member Functions

void SetFSAL (G4bool flag=true)
 
void SetIntegrationOrder (G4int order)
 

Private Attributes

G4ThreeVector BfldIn
 
G4EquationOfMotionfEquation_Rhs = nullptr
 
G4int fIntegrationOrder = -1
 
G4bool fIsFSAL = false
 
const G4int fNoIntegrationVariables = 0
 
unsigned long fNoRHSCalls = 0UL
 
const G4int fNoStateVariables = 0
 
G4ThreeVector fpInitial
 
G4ThreeVector fyFinal
 
G4ThreeVector fyInitial
 
G4ThreeVector fyMidPoint
 
G4double hStep = 0.0
 

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

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 196 of file G4RKG3_Stepper.cc.

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

References G4LineSection::Distline(), fyFinal, fyInitial, and fyMidPoint.

◆ GetEquationOfMotion() [1/2]

G4EquationOfMotion * G4MagIntegratorStepper::GetEquationOfMotion ( )
inlineinherited

◆ GetEquationOfMotion() [2/2]

const G4EquationOfMotion * G4MagIntegratorStepper::GetEquationOfMotion ( ) const
inlineinherited

◆ GetfNoRHSCalls()

unsigned long G4MagIntegratorStepper::GetfNoRHSCalls ( )
inlineinherited

◆ GetNumberOfStateVariables()

G4int G4MagIntegratorStepper::GetNumberOfStateVariables ( ) const
inlineinherited

◆ GetNumberOfVariables()

G4int G4MagIntegratorStepper::GetNumberOfVariables ( ) const
inlineinherited

◆ IntegrationOrder()

G4int G4MagIntegratorStepper::IntegrationOrder ( )
inlineinherited

◆ IntegratorOrder()

G4int G4RKG3_Stepper::IntegratorOrder ( ) const
inlinevirtual

Implements G4MagIntegratorStepper.

Definition at line 88 of file G4RKG3_Stepper.hh.

88{ return 4; }

◆ IsFSAL()

G4bool G4MagIntegratorStepper::IsFSAL ( ) const
inlineinherited

◆ NormalisePolarizationVector()

void G4MagIntegratorStepper::NormalisePolarizationVector ( G4double  vec[12])
inlineinherited

◆ NormaliseTangentVector()

void G4MagIntegratorStepper::NormaliseTangentVector ( G4double  vec[6])
inlineinherited

◆ ResetfNORHSCalls()

void G4MagIntegratorStepper::ResetfNORHSCalls ( )
inlineinherited

◆ RightHandSide() [1/2]

void G4MagIntegratorStepper::RightHandSide ( const G4double  y[],
G4double  dydx[] 
) const
inlineinherited

◆ RightHandSide() [2/2]

void G4MagIntegratorStepper::RightHandSide ( const G4double  y[],
G4double  dydx[],
G4double  field[] 
) const
inlineinherited

◆ SetEquationOfMotion()

void G4MagIntegratorStepper::SetEquationOfMotion ( G4EquationOfMotion newEquation)
inlineinherited

◆ SetFSAL()

void G4MagIntegratorStepper::SetFSAL ( G4bool  flag = true)
inlineprotectedinherited

◆ SetIntegrationOrder()

void G4MagIntegratorStepper::SetIntegrationOrder ( G4int  order)
inlineprotectedinherited

◆ StepNoErr()

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

Definition at line 130 of file G4RKG3_Stepper.cc.

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

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

References B(), BfldIn, G4EquationOfMotion::EvaluateRhsGivenB(), fpInitial, fyFinal, fyInitial, fyMidPoint, G4MagIntegratorStepper::GetEquationOfMotion(), hStep, and StepNoErr().

◆ 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 109 of file G4RKG3_Stepper.cc.

118{
119 G4Exception("G4RKG3_Stepper::StepWithEst()", "GeomField0001",
120 FatalException, "Method no longer used.");
121}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35

References FatalException, and G4Exception().

Field Documentation

◆ BfldIn

G4ThreeVector G4RKG3_Stepper::BfldIn
private

Definition at line 96 of file G4RKG3_Stepper.hh.

Referenced by Stepper().

◆ fEquation_Rhs

G4EquationOfMotion* G4MagIntegratorStepper::fEquation_Rhs = nullptr
privateinherited

Definition at line 124 of file G4MagIntegratorStepper.hh.

◆ fIntegrationOrder

G4int G4MagIntegratorStepper::fIntegrationOrder = -1
privateinherited

Definition at line 134 of file G4MagIntegratorStepper.hh.

◆ fIsFSAL

G4bool G4MagIntegratorStepper::fIsFSAL = false
privateinherited

Definition at line 136 of file G4MagIntegratorStepper.hh.

◆ fNoIntegrationVariables

const G4int G4MagIntegratorStepper::fNoIntegrationVariables = 0
privateinherited

Definition at line 125 of file G4MagIntegratorStepper.hh.

◆ fNoRHSCalls

unsigned long G4MagIntegratorStepper::fNoRHSCalls = 0UL
mutableprivateinherited

Definition at line 128 of file G4MagIntegratorStepper.hh.

◆ fNoStateVariables

const G4int G4MagIntegratorStepper::fNoStateVariables = 0
privateinherited

Definition at line 126 of file G4MagIntegratorStepper.hh.

◆ fpInitial

G4ThreeVector G4RKG3_Stepper::fpInitial
private

Definition at line 95 of file G4RKG3_Stepper.hh.

Referenced by Stepper().

◆ fyFinal

G4ThreeVector G4RKG3_Stepper::fyFinal
private

Definition at line 94 of file G4RKG3_Stepper.hh.

Referenced by DistChord(), and Stepper().

◆ fyInitial

G4ThreeVector G4RKG3_Stepper::fyInitial
private

Definition at line 92 of file G4RKG3_Stepper.hh.

Referenced by DistChord(), and Stepper().

◆ fyMidPoint

G4ThreeVector G4RKG3_Stepper::fyMidPoint
private

Definition at line 93 of file G4RKG3_Stepper.hh.

Referenced by DistChord(), and Stepper().

◆ hStep

G4double G4RKG3_Stepper::hStep = 0.0
private

Definition at line 97 of file G4RKG3_Stepper.hh.

Referenced by Stepper().


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