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

#include <G4CashKarpRKF45.hh>

+ Inheritance diagram for G4CashKarpRKF45:

Public Member Functions

 G4CashKarpRKF45 (G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
 
 ~G4CashKarpRKF45 () override
 
 G4CashKarpRKF45 (const G4CashKarpRKF45 &)=delete
 
G4CashKarpRKF45operator= (const G4CashKarpRKF45 &)=delete
 
void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
 
G4double DistChord () const override
 
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 44 of file G4CashKarpRKF45.hh.

Constructor & Destructor Documentation

◆ G4CashKarpRKF45() [1/2]

G4CashKarpRKF45::G4CashKarpRKF45 ( G4EquationOfMotion * EqRhs,
G4int numberOfVariables = 6,
G4bool primary = true )

Definition at line 47 of file G4CashKarpRKF45.cc.

50 : G4MagIntegratorStepper(EqRhs, noIntegrationVariables)
51{
52 const G4int numberOfVariables =
53 std::max( noIntegrationVariables,
54 ( ( (noIntegrationVariables-1)/4 + 1 ) * 4 ) );
55 // For better alignment with cache-line
56
57 ak2 = new G4double[numberOfVariables] ;
58 ak3 = new G4double[numberOfVariables] ;
59 ak4 = new G4double[numberOfVariables] ;
60 ak5 = new G4double[numberOfVariables] ;
61 ak6 = new G4double[numberOfVariables] ;
62 // ak7 = 0;
63
64 // Must ensure space extra 'state' variables exists - i.e. yIn[7]
65 const G4int numStateMax = std::max(GetNumberOfStateVariables(), 8);
66 const G4int numStateVars = std::max(noIntegrationVariables,
67 numStateMax );
68 // GetNumberOfStateVariables() );
69
70 yTemp = new G4double[numStateVars] ;
71 yIn = new G4double[numStateVars] ;
72
73 fLastInitialVector = new G4double[numStateVars] ;
74 fLastFinalVector = new G4double[numStateVars] ;
75 fLastDyDx = new G4double[numberOfVariables];
76
77 fMidVector = new G4double[numStateVars];
78 fMidError = new G4double[numStateVars];
79 if( primary )
80 {
81 fAuxStepper = new G4CashKarpRKF45(EqRhs, numberOfVariables, !primary);
82 }
83}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4CashKarpRKF45(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)
G4int GetNumberOfStateVariables() const

Referenced by G4CashKarpRKF45().

◆ ~G4CashKarpRKF45()

G4CashKarpRKF45::~G4CashKarpRKF45 ( )
override

Definition at line 89 of file G4CashKarpRKF45.cc.

90{
91 delete [] ak2;
92 delete [] ak3;
93 delete [] ak4;
94 delete [] ak5;
95 delete [] ak6;
96 // delete [] ak7;
97 delete [] yTemp;
98 delete [] yIn;
99
100 delete [] fLastInitialVector;
101 delete [] fLastFinalVector;
102 delete [] fLastDyDx;
103 delete [] fMidVector;
104 delete [] fMidError;
105
106 delete fAuxStepper;
107}

◆ G4CashKarpRKF45() [2/2]

G4CashKarpRKF45::G4CashKarpRKF45 ( const G4CashKarpRKF45 & )
delete

Member Function Documentation

◆ DistChord()

G4double G4CashKarpRKF45::DistChord ( ) const
overridevirtual

Implements G4MagIntegratorStepper.

Definition at line 238 of file G4CashKarpRKF45.cc.

239{
240 G4double distLine, distChord;
241 G4ThreeVector initialPoint, finalPoint, midPoint;
242
243 // Store last initial and final points
244 // (they will be overwritten in self-Stepper call!)
245 //
246 initialPoint = G4ThreeVector( fLastInitialVector[0],
247 fLastInitialVector[1], fLastInitialVector[2]);
248 finalPoint = G4ThreeVector( fLastFinalVector[0],
249 fLastFinalVector[1], fLastFinalVector[2]);
250
251 // Do half a step using StepNoErr
252 //
253 fAuxStepper->Stepper( fLastInitialVector, fLastDyDx,
254 0.5 * fLastStepLength, fMidVector, fMidError );
255
256 midPoint = G4ThreeVector( fMidVector[0], fMidVector[1], fMidVector[2]);
257
258 // Use stored values of Initial and Endpoint + new Midpoint to evaluate
259 // distance of Chord
260 //
261 if (initialPoint != finalPoint)
262 {
263 distLine = G4LineSection::Distline( midPoint, initialPoint, finalPoint );
264 distChord = distLine;
265 }
266 else
267 {
268 distChord = (midPoint-initialPoint).mag();
269 }
270 return distChord;
271}
CLHEP::Hep3Vector G4ThreeVector
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)

◆ IntegratorOrder()

G4int G4CashKarpRKF45::IntegratorOrder ( ) const
inlineoverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 64 of file G4CashKarpRKF45.hh.

64{ return 4; }

◆ operator=()

G4CashKarpRKF45 & G4CashKarpRKF45::operator= ( const G4CashKarpRKF45 & )
delete

◆ Stepper()

void G4CashKarpRKF45::Stepper ( const G4double y[],
const G4double dydx[],
G4double h,
G4double yout[],
G4double yerr[] )
overridevirtual

Implements G4MagIntegratorStepper.

Definition at line 120 of file G4CashKarpRKF45.cc.

125{
126 // const G4int nvar = 6 ;
127 // const G4double a2 = 0.2 , a3 = 0.3 , a4 = 0.6 , a5 = 1.0 , a6 = 0.875;
128 G4int i;
129
130 const G4double b21 = 0.2 ,
131 b31 = 3.0/40.0 , b32 = 9.0/40.0 ,
132 b41 = 0.3 , b42 = -0.9 , b43 = 1.2 ,
133
134 b51 = -11.0/54.0 , b52 = 2.5 , b53 = -70.0/27.0 ,
135 b54 = 35.0/27.0 ,
136
137 b61 = 1631.0/55296.0 , b62 = 175.0/512.0 ,
138 b63 = 575.0/13824.0 , b64 = 44275.0/110592.0 ,
139 b65 = 253.0/4096.0 ,
140
141 c1 = 37.0/378.0 , c3 = 250.0/621.0 , c4 = 125.0/594.0 ,
142 c6 = 512.0/1771.0 , dc5 = -277.0/14336.0 ;
143
144 const G4double dc1 = c1 - 2825.0/27648.0 , dc3 = c3 - 18575.0/48384.0 ,
145 dc4 = c4 - 13525.0/55296.0 , dc6 = c6 - 0.25 ;
146
147 // Initialise time to t0, needed when it is not updated by the integration.
148 // [ Note: Only for time dependent fields (usually electric)
149 // is it neccessary to integrate the time.]
150 yOut[7] = yTemp[7] = yIn[7] = yInput[7];
151
152 const G4int numberOfVariables= this->GetNumberOfVariables();
153 // The number of variables to be integrated over
154
155 // Saving yInput because yInput and yOut can be aliases for same array
156
157 for(i=0; i<numberOfVariables; ++i)
158 {
159 yIn[i]=yInput[i];
160 }
161 // RightHandSide(yIn, dydx) ; // 1st Step
162
163 for(i=0; i<numberOfVariables; ++i)
164 {
165 yTemp[i] = yIn[i] + b21*Step*dydx[i] ;
166 }
167 RightHandSide(yTemp, ak2) ; // 2nd Step
168
169 for(i=0; i<numberOfVariables; ++i)
170 {
171 yTemp[i] = yIn[i] + Step*(b31*dydx[i] + b32*ak2[i]) ;
172 }
173 RightHandSide(yTemp, ak3) ; // 3rd Step
174
175 for(i=0; i<numberOfVariables; ++i)
176 {
177 yTemp[i] = yIn[i] + Step*(b41*dydx[i] + b42*ak2[i] + b43*ak3[i]) ;
178 }
179 RightHandSide(yTemp, ak4) ; // 4th Step
180
181 for(i=0; i<numberOfVariables; ++i)
182 {
183 yTemp[i] = yIn[i] + Step*(b51*dydx[i]
184 + b52*ak2[i] + b53*ak3[i] + b54*ak4[i]) ;
185 }
186 RightHandSide(yTemp, ak5) ; // 5th Step
187
188 for(i=0; i<numberOfVariables; ++i)
189 {
190 yTemp[i] = yIn[i] + Step*(b61*dydx[i]
191 + b62*ak2[i] + b63*ak3[i] + b64*ak4[i] + b65*ak5[i]) ;
192 }
193 RightHandSide(yTemp, ak6) ; // 6th Step
194
195 for(i=0; i<numberOfVariables; ++i)
196 {
197 // Accumulate increments with proper weights
198 //
199 yOut[i] = yIn[i] + Step*(c1*dydx[i] + c3*ak3[i] + c4*ak4[i] + c6*ak6[i]) ;
200
201 // Estimate error as difference between 4th and 5th order methods
202 //
203 yErr[i] = Step*(dc1*dydx[i]
204 + dc3*ak3[i] + dc4*ak4[i] + dc5*ak5[i] + dc6*ak6[i]) ;
205
206 // Store Input and Final values, for possible use in calculating chord
207 //
208 fLastInitialVector[i] = yIn[i] ;
209 fLastFinalVector[i] = yOut[i];
210 fLastDyDx[i] = dydx[i];
211 }
212 // NormaliseTangentVector( yOut ); // Not wanted
213
214 fLastStepLength = Step;
215
216 return;
217}
G4int GetNumberOfVariables() const
void RightHandSide(const G4double y[], G4double dydx[]) const

Referenced by DistChord().


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