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

#include <G4DoLoMcPriRK34.hh>

+ Inheritance diagram for G4DoLoMcPriRK34:

Public Member Functions

 G4DoLoMcPriRK34 (G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
 
 ~G4DoLoMcPriRK34 ()
 
 G4DoLoMcPriRK34 (const G4DoLoMcPriRK34 &)=delete
 
G4DoLoMcPriRK34operator= (const G4DoLoMcPriRK34 &)=delete
 
void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
 
void SetupInterpolation ()
 
void SetupInterpolate (const G4double yInput[], const G4double dydx[], const G4double Step)
 
void Interpolate (const G4double yInput[], const G4double dydx[], const G4double Step, G4double yOut[], G4double tau)
 
void Interpolate (G4double tau, G4double yOut[])
 
void interpolate (const G4double yInput[], const G4double dydx[], G4double yOut[], G4double Step, G4double tau)
 
G4double DistChord () const
 
G4int IntegratorOrder () const
 
- 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
 
virtual void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])=0
 
virtual G4double DistChord () const =0
 
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
 
virtual G4int IntegratorOrder () const =0
 
G4int IntegrationOrder ()
 
G4EquationOfMotionGetEquationOfMotion ()
 
const G4EquationOfMotionGetEquationOfMotion () const
 
void SetEquationOfMotion (G4EquationOfMotion *newEquation)
 
unsigned long GetfNoRHSCalls ()
 
void ResetfNORHSCalls ()
 
G4bool IsFSAL () const
 

Additional Inherited Members

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

Detailed Description

Definition at line 41 of file G4DoLoMcPriRK34.hh.

Constructor & Destructor Documentation

◆ G4DoLoMcPriRK34() [1/2]

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

Definition at line 37 of file G4DoLoMcPriRK34.cc.

40 : G4MagIntegratorStepper(EqRhs, noIntegrationVariables)
41{
42 const G4int numberOfVariables = noIntegrationVariables;
43
44 // New Chunk of memory being created for use by the Stepper
45
46 // aki - for storing intermediate RHS
47 //
48 ak2 = new G4double[numberOfVariables];
49 ak3 = new G4double[numberOfVariables];
50 ak4 = new G4double[numberOfVariables];
51 ak5 = new G4double[numberOfVariables];
52 ak6 = new G4double[numberOfVariables];
53
54 yTemp = new G4double[numberOfVariables] ;
55 yIn = new G4double[numberOfVariables] ;
56
57 fLastInitialVector = new G4double[numberOfVariables] ;
58 fLastFinalVector = new G4double[numberOfVariables] ;
59 fLastDyDx = new G4double[numberOfVariables];
60
61 fMidVector = new G4double[numberOfVariables];
62 fMidError = new G4double[numberOfVariables];
63 if( primary )
64 {
65 fAuxStepper = new G4DoLoMcPriRK34(EqRhs, numberOfVariables, !primary);
66 }
67}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85

◆ ~G4DoLoMcPriRK34()

G4DoLoMcPriRK34::~G4DoLoMcPriRK34 ( )

Definition at line 71 of file G4DoLoMcPriRK34.cc.

72{
73 // clear all previously allocated memory for Stepper and DistChord
74
75 delete [] ak2;
76 delete [] ak3;
77 delete [] ak4;
78 delete [] ak5;
79 delete [] ak6;
80
81 delete [] yTemp;
82 delete [] yIn;
83
84 delete [] fLastInitialVector;
85 delete [] fLastFinalVector;
86 delete [] fLastDyDx;
87 delete [] fMidVector;
88 delete [] fMidError;
89
90 delete fAuxStepper;
91}

◆ G4DoLoMcPriRK34() [2/2]

G4DoLoMcPriRK34::G4DoLoMcPriRK34 ( const G4DoLoMcPriRK34 )
delete

Member Function Documentation

◆ DistChord()

G4double G4DoLoMcPriRK34::DistChord ( ) const
virtual

Implements G4MagIntegratorStepper.

Definition at line 201 of file G4DoLoMcPriRK34.cc.

202{
203 G4double distLine, distChord;
204 G4ThreeVector initialPoint, finalPoint, midPoint;
205
206 // Store last initial and final points
207 // (they will be overwritten in self-Stepper call!)
208 //
209 initialPoint = G4ThreeVector( fLastInitialVector[0],
210 fLastInitialVector[1], fLastInitialVector[2] );
211 finalPoint = G4ThreeVector( fLastFinalVector[0],
212 fLastFinalVector[1], fLastFinalVector[2] );
213
214 // Do half a Step using StepNoErr
215
216 fAuxStepper->Stepper( fLastInitialVector, fLastDyDx, 0.5 * fLastStepLength,
217 fMidVector, fMidError );
218
219 midPoint = G4ThreeVector( fMidVector[0], fMidVector[1], fMidVector[2]);
220
221 // Use stored values of Initial and Endpoint + new Midpoint to evaluate
222 // distance of Chord
223 //
224 if (initialPoint != finalPoint)
225 {
226 distLine = G4LineSection::Distline( midPoint, initialPoint, finalPoint );
227 distChord = distLine;
228 }
229 else
230 {
231 distChord = (midPoint-initialPoint).mag();
232 }
233 return distChord;
234}
CLHEP::Hep3Vector G4ThreeVector
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)

◆ IntegratorOrder()

G4int G4DoLoMcPriRK34::IntegratorOrder ( ) const
inlinevirtual

Implements G4MagIntegratorStepper.

Definition at line 85 of file G4DoLoMcPriRK34.hh.

85{ return 3; }

◆ Interpolate() [1/2]

void G4DoLoMcPriRK34::Interpolate ( const G4double  yInput[],
const G4double  dydx[],
const G4double  Step,
G4double  yOut[],
G4double  tau 
)

Definition at line 255 of file G4DoLoMcPriRK34.cc.

260{
261 G4double bf1, bf2, bf3, bf4, bf5, bf6;
262
263 const G4int numberOfVariables = GetNumberOfVariables();
264
265 for(G4int i=0; i<numberOfVariables; ++i)
266 {
267 yIn[i]=yInput[i];
268 }
269
270 G4double tau_2 = tau*tau, tau_3 = tau*tau_2;
271
272 // Calculating the polynomials (coefficients for the respective stages)
273 //
274 bf1 = -(162.0*tau_3 - 504.0*tau_2 + 551.0*tau - 238.0)/238.0 ,
275 bf2 = 0.0 ,
276 bf3 = 27.0*tau*(27.0*tau_2 - 70.0*tau + 51.0 )/385.0 ,
277 bf4 = -27*tau*(27.0*tau_2 - 50.0*tau + 21.0)/85.0 ,
278 bf5 = 7.0*tau*(2232.0*tau_2 - 4166.0*tau + 1785.0 )/3278.0 ,
279 bf6 = tau*(tau - 1.0)*(387.0*tau - 238.0)/149.0 ;
280
281 for( G4int i=0; i<numberOfVariables; ++i)
282 {
283 yOut[i] = yIn[i] + Step*tau*(bf1*dydx[i] + bf2*ak2[i] + bf3*ak3[i]
284 + bf4*ak4[i] + bf5*ak5[i] + bf6*ak6[i] ) ;
285 }
286}
G4int GetNumberOfVariables() const

Referenced by Interpolate().

◆ interpolate()

void G4DoLoMcPriRK34::interpolate ( const G4double  yInput[],
const G4double  dydx[],
G4double  yOut[],
G4double  Step,
G4double  tau 
)

◆ Interpolate() [2/2]

void G4DoLoMcPriRK34::Interpolate ( G4double  tau,
G4double  yOut[] 
)

Definition at line 247 of file G4DoLoMcPriRK34.cc.

249{
250 Interpolate( fLastInitialVector, fLastDyDx, fLastStepLength, yOut, tau );
251}
void Interpolate(const G4double yInput[], const G4double dydx[], const G4double Step, G4double yOut[], G4double tau)

◆ operator=()

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

◆ SetupInterpolate()

void G4DoLoMcPriRK34::SetupInterpolate ( const G4double  yInput[],
const G4double  dydx[],
const G4double  Step 
)

Definition at line 240 of file G4DoLoMcPriRK34.cc.

243{
244 // Do Nothing
245}

◆ SetupInterpolation()

void G4DoLoMcPriRK34::SetupInterpolation ( )

Definition at line 236 of file G4DoLoMcPriRK34.cc.

237{
238}

◆ Stepper()

void G4DoLoMcPriRK34::Stepper ( const G4double  y[],
const G4double  dydx[],
G4double  h,
G4double  yout[],
G4double  yerr[] 
)
virtual

Implements G4MagIntegratorStepper.

Definition at line 98 of file G4DoLoMcPriRK34.cc.

103{
104 G4int i;
105
106 // The various constants defined on the basis of butcher tableu
107 //
108 const G4double b21 = 7.0/27.0 ,
109 b31 = 7.0/72.0 ,
110 b32 = 7.0/24.0 ,
111
112 b41 = 3043.0/3528.0 ,
113 b42 = -3757.0/1176.0 ,
114 b43 = 1445.0/441.0,
115
116 b51 = 17617.0/11662.0 ,
117 b52 = -4023.0/686.0 ,
118 b53 = 9372.0/1715.0 ,
119 b54 = -66.0/595.0 ,
120
121 b61 = 29.0/238.0 ,
122 b62 = 0.0 ,
123 b63 = 216.0/385.0 ,
124 b64 = 54.0/85.0 ,
125 b65 = -7.0/22.0 ,
126
127 dc1 = 363.0/2975.0 - b61 ,
128 dc2 = 0.0 - b62 ,
129 dc3 = 981.0/1750.0 - b63,
130 dc4 = 2709.0/4250.0 - b64 ,
131 dc5 = -3.0/10.0 - b65 ,
132 dc6 = -1.0/50.0 ; // end of declaration
133
134 const G4int numberOfVariables = GetNumberOfVariables();
135
136 // The number of variables to be integrated over
137 //
138 yOut[7] = yTemp[7] = yIn[7];
139
140 // Saving yInput because yInput and yOut can be aliases for same array
141 //
142 for(i=0; i<numberOfVariables; ++i)
143 {
144 yIn[i]=yInput[i];
145 }
146
147 // RightHandSide(yIn, DyDx) ; // 1st stage - Not doing, getting passed
148
149 for(i=0; i<numberOfVariables; ++i)
150 {
151 yTemp[i] = yIn[i] + b21*Step*DyDx[i] ;
152 }
153 RightHandSide(yTemp, ak2) ; // 2nd stage
154
155 for(i=0; i<numberOfVariables; ++i)
156 {
157 yTemp[i] = yIn[i] + Step*(b31*DyDx[i] + b32*ak2[i]) ;
158 }
159 RightHandSide(yTemp, ak3) ; // 3rd stage
160
161 for(i=0; i<numberOfVariables; ++i)
162 {
163 yTemp[i] = yIn[i] + Step*(b41*DyDx[i] + b42*ak2[i] + b43*ak3[i]) ;
164 }
165 RightHandSide(yTemp, ak4) ; // 4th stage
166
167 for(i=0; i<numberOfVariables; ++i)
168 {
169 yTemp[i] = yIn[i] + Step*(b51*DyDx[i] + b52*ak2[i]
170 + b53*ak3[i] + b54*ak4[i]) ;
171 }
172 RightHandSide(yTemp, ak5) ; // 5th stage
173
174 for(i=0; i<numberOfVariables; ++i)
175 {
176 yOut[i] = yIn[i] + Step*(b61*DyDx[i] + b62*ak2[i] + b63*ak3[i]
177 + b64*ak4[i] + b65*ak5[i]) ;
178 }
179 RightHandSide(yOut, ak6) ; // 6th and Final stage
180
181 for(i=0; i<numberOfVariables; ++i)
182 {
183
184 yErr[i] = Step*(dc1*DyDx[i] + dc2*ak2[i] + dc3*ak3[i] + dc4*ak4[i]
185 + dc5*ak5[i] + dc6*ak6[i] ) ;
186
187 // Store Input and Final values, for possible use in calculating chord
188 //
189 fLastInitialVector[i] = yIn[i] ;
190 fLastFinalVector[i] = yOut[i];
191 fLastDyDx[i] = DyDx[i];
192 }
193
194 fLastStepLength = Step;
195
196 return ;
197}
void RightHandSide(const G4double y[], G4double dydx[]) const

Referenced by DistChord().


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