Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NystromRK4 Class Reference

#include <G4NystromRK4.hh>

Inheritance diagram for G4NystromRK4:
Collaboration 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 G4double P[], G4double 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, bool isFSAL=false)
 
virtual ~G4MagIntegratorStepper ()
 
void NormaliseTangentVector (G4double vec[6])
 
void NormalisePolarizationVector (G4double vec[12])
 
void RightHandSide (const double y[], double dydx[])
 
G4int GetNumberOfVariables () const
 
G4int GetNumberOfStateVariables () const
 
G4int IntegrationOrder ()
 
G4EquationOfMotionGetEquationOfMotion ()
 
void SetEquationOfMotion (G4EquationOfMotion *newEquation)
 
unsigned long GetfNoRHSCalls ()
 
void ResetfNORHSCalls ()
 
bool IsFSAL ()
 

Additional Inherited Members

- Protected Member Functions inherited from G4MagIntegratorStepper
void SetIntegrationOrder (int order)
 
void SetFSAL (bool flag=true)
 

Detailed Description

Definition at line 51 of file G4NystromRK4.hh.

Constructor & Destructor Documentation

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 }
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, bool isFSAL=false)
G4NystromRK4::~G4NystromRK4 ( )

Definition at line 63 of file G4NystromRK4.cc.

64 {
65 }

Member Function Documentation

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

Reimplemented from G4MagIntegratorStepper.

Definition at line 210 of file G4NystromRK4.cc.

211 {
212  G4double P4vec[4]= { P[0], P[1], P[2], P[7] }; // Time is P[7]
213  getField(P4vec);
214  m_mom = std::sqrt(P[3]*P[3]+P[4]*P[4]+P[5]*P[5]) ;
215  m_imom = 1./m_mom ;
216  m_cof = m_fEq->FCof()*m_imom ;
217  m_cachedMom = true ; // Caching the value
218  dPdS[0] = P[3]*m_imom ; // dx /ds
219  dPdS[1] = P[4]*m_imom ; // dy /ds
220  dPdS[2] = P[5]*m_imom ; // dz /ds
221  dPdS[3] = m_cof*(P[4]*m_lastField[2]-P[5]*m_lastField[1]) ; // dPx/ds
222  dPdS[4] = m_cof*(P[5]*m_lastField[0]-P[3]*m_lastField[2]) ; // dPy/ds
223  dPdS[5] = m_cof*(P[3]*m_lastField[1]-P[4]*m_lastField[0]) ; // dPz/ds
224 }
static double P[]
G4double FCof() const
Definition: G4Mag_EqRhs.hh:84
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4NystromRK4::DistChord ( ) const
virtual

Implements G4MagIntegratorStepper.

Definition at line 186 of file G4NystromRK4.cc.

187 {
188  G4double ax = m_fPoint[0]-m_iPoint[0];
189  G4double ay = m_fPoint[1]-m_iPoint[1];
190  G4double az = m_fPoint[2]-m_iPoint[2];
191  G4double dx = m_mPoint[0]-m_iPoint[0];
192  G4double dy = m_mPoint[1]-m_iPoint[1];
193  G4double dz = m_mPoint[2]-m_iPoint[2];
194  G4double d2 = (ax*ax+ay*ay+az*az) ;
195 
196  if(d2!=0.) {
197  G4double ds = (ax*dx+ay*dy+az*dz)/d2;
198  dx -= (ds*ax) ;
199  dy -= (ds*ay) ;
200  dz -= (ds*az) ;
201  }
202  return std::sqrt(dx*dx+dy*dy+dz*dz);
203 }
static const G4double d2
double G4double
Definition: G4Types.hh:76
G4double G4NystromRK4::GetDistanceForConstantField ( ) const
inline

Definition at line 112 of file G4NystromRK4.hh.

113 {
114  return m_magdistance;
115 }
G4int G4NystromRK4::IntegratorOrder ( ) const
inlinevirtual

Implements G4MagIntegratorStepper.

Definition at line 73 of file G4NystromRK4.hh.

73 {return 4;}
void G4NystromRK4::SetDistanceForConstantField ( G4double  length)
inline

Definition at line 106 of file G4NystromRK4.hh.

107 {
108  m_magdistance= length;
109  m_magdistance2 = length*length;
110 }
void G4NystromRK4::Stepper ( const G4double  P[],
const G4double  dPdS[],
G4double  step,
G4double  Po[],
G4double  Err[] 
)
virtual

Implements G4MagIntegratorStepper.

Definition at line 73 of file G4NystromRK4.cc.

74 {
75  const G4double perMillion = 1.0e-6;
76  G4double R[4] = { P[0], P[1] , P[2], P[7] }; // x, y, z, t
77  G4double A[3] = {dPdS[0], dPdS[1], dPdS[2]};
78 
79  m_iPoint[0]=R[0]; m_iPoint[1]=R[1]; m_iPoint[2]=R[2];
80 
81  constexpr G4double one_sixth= 1./6.;
82  const G4double S = Step ;
83  const G4double S5 = .5*Step ;
84  const G4double S4 = .25*Step ;
85  const G4double S6 = Step * one_sixth; // Step / 6.;
86 
87  // Ensure that the location and cached field value are correct
88  getField( R );
89 
90  // Ensure that the momentum is set correctly.
91 
92  // - Quick check momentum magnitude (squared) against previous value
93  G4double newmom2 = (P[3]*P[3]+P[4]*P[4]+P[5]*P[5]);
94  G4double oldmom2 = m_mom * m_mom;
95  if( std::fabs(newmom2 - oldmom2) > perMillion * oldmom2 ) {
96  m_mom = sqrt(newmom2) ;
97  m_imom = 1./m_mom;
98  m_cof = m_fEq->FCof()*m_imom;
99  }
100 
101 #ifdef G4DEBUG_FIELD
102  CheckCachedMomemtum( P, m_mom );
103  CheckFieldPosition( P, m_fldPosition );
104 #endif
105 
106  // Point 1
107  //
108  G4double K1[3] = { m_imom*dPdS[3], m_imom*dPdS[4], m_imom*dPdS[5] };
109 
110  // Point2
111  //
112  G4double p[4] = {R[0]+S5*(A[0]+S4*K1[0]),
113  R[1]+S5*(A[1]+S4*K1[1]),
114  R[2]+S5*(A[2]+S4*K1[2]),
115  P[7] };
116  getField(p);
117 
118  G4double A2[3] = {A[0]+S5*K1[0],A[1]+S5*K1[1],A[2]+S5*K1[2]};
119  G4double K2[3] = {(A2[1]*m_lastField[2]-A2[2]*m_lastField[1])*m_cof,
120  (A2[2]*m_lastField[0]-A2[0]*m_lastField[2])*m_cof,
121  (A2[0]*m_lastField[1]-A2[1]*m_lastField[0])*m_cof};
122 
123  m_mPoint[0]=p[0]; m_mPoint[1]=p[1]; m_mPoint[2]=p[2];
124 
125  // Point 3 with the same magnetic field
126  //
127  G4double A3[3] = {A[0]+S5*K2[0],A[1]+S5*K2[1],A[2]+S5*K2[2]};
128  G4double K3[3] = {(A3[1]*m_lastField[2]-A3[2]*m_lastField[1])*m_cof,
129  (A3[2]*m_lastField[0]-A3[0]*m_lastField[2])*m_cof,
130  (A3[0]*m_lastField[1]-A3[1]*m_lastField[0])*m_cof};
131 
132  // Point 4
133  //
134  p[0] = R[0]+S*(A[0]+S5*K3[0]);
135  p[1] = R[1]+S*(A[1]+S5*K3[1]);
136  p[2] = R[2]+S*(A[2]+S5*K3[2]);
137 
138  getField(p);
139 
140  G4double A4[3] = {A[0]+S*K3[0],A[1]+S*K3[1],A[2]+S*K3[2]};
141  G4double K4[3] = {(A4[1]*m_lastField[2]-A4[2]*m_lastField[1])*m_cof,
142  (A4[2]*m_lastField[0]-A4[0]*m_lastField[2])*m_cof,
143  (A4[0]*m_lastField[1]-A4[1]*m_lastField[0])*m_cof};
144 
145  // New position
146  //
147  Po[0] = P[0]+S*(A[0]+S6*(K1[0]+K2[0]+K3[0]));
148  Po[1] = P[1]+S*(A[1]+S6*(K1[1]+K2[1]+K3[1]));
149  Po[2] = P[2]+S*(A[2]+S6*(K1[2]+K2[2]+K3[2]));
150 
151  m_fPoint[0]=Po[0]; m_fPoint[1]=Po[1]; m_fPoint[2]=Po[2];
152 
153  // New direction
154  //
155  Po[3] = A[0]+S6*(K1[0]+K4[0]+2.*(K2[0]+K3[0]));
156  Po[4] = A[1]+S6*(K1[1]+K4[1]+2.*(K2[1]+K3[1]));
157  Po[5] = A[2]+S6*(K1[2]+K4[2]+2.*(K2[2]+K3[2]));
158 
159  // Errors
160  //
161  Err[3] = S*std::fabs(K1[0]-K2[0]-K3[0]+K4[0]);
162  Err[4] = S*std::fabs(K1[1]-K2[1]-K3[1]+K4[1]);
163  Err[5] = S*std::fabs(K1[2]-K2[2]-K3[2]+K4[2]);
164  Err[0] = S*Err[3] ;
165  Err[1] = S*Err[4] ;
166  Err[2] = S*Err[5] ;
167  Err[3]*= m_mom ;
168  Err[4]*= m_mom ;
169  Err[5]*= m_mom ;
170 
171  // Normalize momentum
172  //
173  G4double normF = m_mom/std::sqrt(Po[3]*Po[3]+Po[4]*Po[4]+Po[5]*Po[5]);
174  Po [3]*=normF; Po[4]*=normF; Po[5]*=normF;
175 
176  // Pass Energy, time unchanged -- time is not integrated !!
177  Po[6]=P[6]; Po[7]=P[7];
178 }
static constexpr double perMillion
Definition: G4SIunits.hh:334
double S(double temp)
const char * p
Definition: xmltok.h:285
static double P[]
double A(double temperature)
G4double FCof() const
Definition: G4Mag_EqRhs.hh:84
Definition: Step.hh:41
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:


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