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

#include <G4NystromRK4.hh>

+ Inheritance diagram for G4NystromRK4:

Public Member Functions

 G4NystromRK4 (G4Mag_EqRhs *EquationMotion, G4double distanceConstField=0.0)
 
 ~G4NystromRK4 ()
 
void Stepper (const G4double P[], const G4double dPdS[], G4double step, G4double Po[], G4double Err[])
 
virtual void ComputeRightHandSide (const double P[], double dPdS[])
 
void SetDistanceForConstantField (G4double length)
 
G4double GetDistanceForConstantField () const
 
G4int IntegratorOrder () const
 
G4double DistChord () 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 G4NystromRK4.hh.

Constructor & Destructor Documentation

◆ G4NystromRK4()

G4NystromRK4::G4NystromRK4 ( G4Mag_EqRhs EquationMotion,
G4double  distanceConstField = 0.0 
)

Definition at line 41 of file G4NystromRK4.cc.

42 : G4MagIntegratorStepper(magEqRhs, 6), // number of variables
43 m_fEq( magEqRhs ),
44 m_magdistance( distanceConstField ),
45 m_cof( 0.0 ),
46 m_mom( 0.0 ),
47 m_imom( 0.0 ),
48 m_cachedMom( false )
49{
50 m_fldPosition[0] = m_iPoint[0] = m_fPoint[0] = m_mPoint[0] = 9.9999999e+99 ;
51 m_fldPosition[1] = m_iPoint[1] = m_fPoint[1] = m_mPoint[1] = 9.9999999e+99 ;
52 m_fldPosition[2] = m_iPoint[2] = m_fPoint[2] = m_mPoint[2] = 9.9999999e+99 ;
53 m_fldPosition[3] = -9.9999999e+99;
54 m_lastField[0] = m_lastField[1] = m_lastField[2] = 0.0;
55
56 m_magdistance2 = distanceConstField*distanceConstField;
57}

◆ ~G4NystromRK4()

G4NystromRK4::~G4NystromRK4 ( )

Definition at line 63 of file G4NystromRK4.cc.

64{
65}

Member Function Documentation

◆ ComputeRightHandSide()

void G4NystromRK4::ComputeRightHandSide ( const double  P[],
double  dPdS[] 
)
virtual

Reimplemented from G4MagIntegratorStepper.

Definition at line 196 of file G4NystromRK4.cc.

197{
198 G4double P4vec[4]= { P[0], P[1], P[2], P[7] }; // Time is P[7]
199 getField(P4vec);
200 m_mom = std::sqrt(P[3]*P[3]+P[4]*P[4]+P[5]*P[5]) ;
201 m_imom = 1./m_mom ;
202 m_cof = m_fEq->FCof()*m_imom ;
203 m_cachedMom = true ; // Caching the value
204 dPdS[0] = P[3]*m_imom ; // dx /ds
205 dPdS[1] = P[4]*m_imom ; // dy /ds
206 dPdS[2] = P[5]*m_imom ; // dz /ds
207 dPdS[3] = m_cof*(P[4]*m_lastField[2]-P[5]*m_lastField[1]) ; // dPx/ds
208 dPdS[4] = m_cof*(P[5]*m_lastField[0]-P[3]*m_lastField[2]) ; // dPy/ds
209 dPdS[5] = m_cof*(P[3]*m_lastField[1]-P[4]*m_lastField[0]) ; // dPz/ds
210}
double G4double
Definition: G4Types.hh:64
G4double FCof() const
Definition: G4Mag_EqRhs.hh:83

◆ DistChord()

G4double G4NystromRK4::DistChord ( ) const
virtual

Implements G4MagIntegratorStepper.

Definition at line 172 of file G4NystromRK4.cc.

173{
174 G4double ax = m_fPoint[0]-m_iPoint[0];
175 G4double ay = m_fPoint[1]-m_iPoint[1];
176 G4double az = m_fPoint[2]-m_iPoint[2];
177 G4double dx = m_mPoint[0]-m_iPoint[0];
178 G4double dy = m_mPoint[1]-m_iPoint[1];
179 G4double dz = m_mPoint[2]-m_iPoint[2];
180 G4double d2 = (ax*ax+ay*ay+az*az) ;
181
182 if(d2!=0.) {
183 G4double ds = (ax*dx+ay*dy+az*dz)/d2;
184 dx -= (ds*ax) ;
185 dy -= (ds*ay) ;
186 dz -= (ds*az) ;
187 }
188 return std::sqrt(dx*dx+dy*dy+dz*dz);
189}

◆ GetDistanceForConstantField()

G4double G4NystromRK4::GetDistanceForConstantField ( ) const
inline

Definition at line 108 of file G4NystromRK4.hh.

109{
110 return m_magdistance;
111}

◆ IntegratorOrder()

G4int G4NystromRK4::IntegratorOrder ( ) const
inlinevirtual

Implements G4MagIntegratorStepper.

Definition at line 73 of file G4NystromRK4.hh.

73{return 4;}

◆ SetDistanceForConstantField()

void G4NystromRK4::SetDistanceForConstantField ( G4double  length)
inline

Definition at line 102 of file G4NystromRK4.hh.

103{
104 m_magdistance= length;
105 m_magdistance2 = length*length;
106}

◆ Stepper()

void G4NystromRK4::Stepper ( const G4double  P[],
const G4double  dPdS[],
G4double  step,
G4double  Po[],
G4double  Err[] 
)
virtual

Implements G4MagIntegratorStepper.

Definition at line 72 of file G4NystromRK4.cc.

74{
75 G4double R[3] = { P[0], P[1] , P[2]};
76 G4double A[3] = {dPdS[0], dPdS[1], dPdS[2]};
77
78 m_iPoint[0]=R[0]; m_iPoint[1]=R[1]; m_iPoint[2]=R[2];
79
80 const G4double one_sixth= 1./6.;
81 G4double S = Step ;
82 G4double S5 = .5*Step ;
83 G4double S4 = .25*Step ;
84 G4double S6 = Step * one_sixth; // Step / 6.;
85
86
87 // John A added, in order to emulate effect of call to changed/derived RHS
88 // m_mom = sqrt(P[3]*P[3]+P[4]*P[4]+P[5]*P[5]);
89 // m_imom = 1./m_mom;
90 // m_cof = m_fEq->FCof()*m_imom;
91
92 // Point 1
93 //
94 G4double K1[3] = { m_imom*dPdS[3], m_imom*dPdS[4], m_imom*dPdS[5] };
95
96 // Point2
97 //
98 G4double p[4] = {R[0]+S5*(A[0]+S4*K1[0]),
99 R[1]+S5*(A[1]+S4*K1[1]),
100 R[2]+S5*(A[2]+S4*K1[2]),
101 P[7] };
102 getField(p);
103
104 G4double A2[3] = {A[0]+S5*K1[0],A[1]+S5*K1[1],A[2]+S5*K1[2]};
105 G4double K2[3] = {(A2[1]*m_lastField[2]-A2[2]*m_lastField[1])*m_cof,
106 (A2[2]*m_lastField[0]-A2[0]*m_lastField[2])*m_cof,
107 (A2[0]*m_lastField[1]-A2[1]*m_lastField[0])*m_cof};
108
109 m_mPoint[0]=p[0]; m_mPoint[1]=p[1]; m_mPoint[2]=p[2];
110
111 // Point 3 with the same magnetic field
112 //
113 G4double A3[3] = {A[0]+S5*K2[0],A[1]+S5*K2[1],A[2]+S5*K2[2]};
114 G4double K3[3] = {(A3[1]*m_lastField[2]-A3[2]*m_lastField[1])*m_cof,
115 (A3[2]*m_lastField[0]-A3[0]*m_lastField[2])*m_cof,
116 (A3[0]*m_lastField[1]-A3[1]*m_lastField[0])*m_cof};
117
118 // Point 4
119 //
120 p[0] = R[0]+S*(A[0]+S5*K3[0]);
121 p[1] = R[1]+S*(A[1]+S5*K3[1]);
122 p[2] = R[2]+S*(A[2]+S5*K3[2]);
123
124 getField(p);
125
126 G4double A4[3] = {A[0]+S*K3[0],A[1]+S*K3[1],A[2]+S*K3[2]};
127 G4double K4[3] = {(A4[1]*m_lastField[2]-A4[2]*m_lastField[1])*m_cof,
128 (A4[2]*m_lastField[0]-A4[0]*m_lastField[2])*m_cof,
129 (A4[0]*m_lastField[1]-A4[1]*m_lastField[0])*m_cof};
130
131 // New position
132 //
133 Po[0] = P[0]+S*(A[0]+S6*(K1[0]+K2[0]+K3[0]));
134 Po[1] = P[1]+S*(A[1]+S6*(K1[1]+K2[1]+K3[1]));
135 Po[2] = P[2]+S*(A[2]+S6*(K1[2]+K2[2]+K3[2]));
136
137 m_fPoint[0]=Po[0]; m_fPoint[1]=Po[1]; m_fPoint[2]=Po[2];
138
139 // New direction
140 //
141 Po[3] = A[0]+S6*(K1[0]+K4[0]+2.*(K2[0]+K3[0]));
142 Po[4] = A[1]+S6*(K1[1]+K4[1]+2.*(K2[1]+K3[1]));
143 Po[5] = A[2]+S6*(K1[2]+K4[2]+2.*(K2[2]+K3[2]));
144
145 // Errors
146 //
147 Err[3] = S*std::fabs(K1[0]-K2[0]-K3[0]+K4[0]);
148 Err[4] = S*std::fabs(K1[1]-K2[1]-K3[1]+K4[1]);
149 Err[5] = S*std::fabs(K1[2]-K2[2]-K3[2]+K4[2]);
150 Err[0] = S*Err[3] ;
151 Err[1] = S*Err[4] ;
152 Err[2] = S*Err[5] ;
153 Err[3]*= m_mom ;
154 Err[4]*= m_mom ;
155 Err[5]*= m_mom ;
156
157 // Normalize momentum
158 //
159 G4double normF = m_mom/std::sqrt(Po[3]*Po[3]+Po[4]*Po[4]+Po[5]*Po[5]);
160 Po [3]*=normF; Po[4]*=normF; Po[5]*=normF;
161
162 // Pass Energy, time unchanged -- time is not integrated !!
163 Po[6]=P[6]; Po[7]=P[7];
164}

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