G4TDormandPrince45 is a templated version of G4DormandPrince745 5th order Runge-Kutta stepper.
More...
#include <G4TDormandPrince45.hh>
|
| | G4TDormandPrince45 (T_Equation *equation) |
| | G4TDormandPrince45 (T_Equation *equation, G4int numVar) |
| void | StepWithError (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) |
| void | Stepper (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) final |
| void | StepWithFinalDerivate (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[], G4double dydxOutput[]) |
| void | SetupInterpolation () |
| void | Interpolate (G4double tau, G4double yOut[]) const |
| G4double | DistChord () const final |
| G4int | IntegratorOrder () const override |
| G4StepperType | StepperType () const override |
| const field_utils::ShortState< N > & | GetYOut () const |
| void | Interpolate4thOrder (G4double yOut[], G4double tau) const |
| void | SetupInterpolation5thOrder () |
| void | Interpolate5thOrder (G4double yOut[], G4double tau) const |
| void | RightHandSideInl (const G4double y[], G4double dydx[]) |
| void | Stepper (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[], G4double dydxOutput[]) |
| T_Equation * | GetSpecificEquation () |
| | G4MagIntegratorStepper (G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false) |
| virtual | ~G4MagIntegratorStepper ()=default |
| | G4MagIntegratorStepper (const G4MagIntegratorStepper &)=delete |
| G4MagIntegratorStepper & | operator= (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 () |
| G4EquationOfMotion * | GetEquationOfMotion () |
| const G4EquationOfMotion * | GetEquationOfMotion () const |
| void | SetEquationOfMotion (G4EquationOfMotion *newEquation) |
| unsigned long | GetfNoRHSCalls () |
| void | ResetfNORHSCalls () |
| G4bool | IsFSAL () const |
| G4bool | isQSS () const |
| void | SetIsQSS (G4bool val) |
template<class T_Equation, unsigned int N = 6>
class G4TDormandPrince45< T_Equation, N >
G4TDormandPrince45 is a templated version of G4DormandPrince745 5th order Runge-Kutta stepper.
Definition at line 55 of file G4TDormandPrince45.hh.
◆ G4TDormandPrince45() [1/2]
template<class T_Equation, unsigned int N>
Definition at line 164 of file G4TDormandPrince45.hh.
167{
168
170 {
171 G4Exception(
"G4TDormandPrince745: constructor",
"GeomField0001",
173 }
174
175
176
177
178
179
180
181
182
183
184
185}
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)
G4TDormandPrince45 is a templated version of G4DormandPrince745 5th order Runge-Kutta stepper.
Referenced by G4TDormandPrince45().
◆ G4TDormandPrince45() [2/2]
template<class T_Equation, unsigned int N>
Definition at line 188 of file G4TDormandPrince45.hh.
191{
193 {
195 msg <<
"Equation has an incompatible number of variables." ;
196 msg <<
" template N = " <<
N
197 <<
" argument numVar = " <<
numVar ;
198
199 G4Exception(
"G4TCashKarpRKF45: constructor",
"GeomField0001",
201 }
203}
G4TDormandPrince45(T_Equation *equation)
◆ DistChord()
template<class T_Equation, unsigned int N>
Estimates the maximum distance of a chord from the true path over the segment last integrated.
Implements G4MagIntegratorStepper.
Definition at line 352 of file G4TDormandPrince45.hh.
353{
354
355
356
358 hf3 = 51252292925.0 / 65400821598.0,
359 hf4 = - 2691868925.0 / 45128329728.0,
360 hf5 = 187940372067.0 / 1594534317056.0,
361 hf6 = - 1776094331.0 / 19743644256.0,
362 hf7 = 11237099.0 / 235043384.0;
363
365
366 for(
unsigned int i = 0;
i < 3; ++
i)
367 {
368 mid[
i] = fyIn[
i] + 0.5 * fLastStepLength * (
371 }
372
375
377}
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
◆ GetSpecificEquation()
template<class T_Equation, unsigned int N = 6>
◆ GetYOut()
template<class T_Equation, unsigned int N = 6>
◆ IntegratorOrder()
template<class T_Equation, unsigned int N = 6>
◆ Interpolate()
template<class T_Equation, unsigned int N = 6>
Definition at line 83 of file G4TDormandPrince45.hh.
84 {
86 }
void Interpolate4thOrder(G4double yOut[], G4double tau) const
◆ Interpolate4thOrder()
template<class T_Equation, unsigned int N>
Definition at line 386 of file G4TDormandPrince45.hh.
388{
392
394 157015080.0 *
tau4 - 13107642775.0 *
tau3 + 34969693132.0 *
tau2 -
395 32272833064.0 *
tau + 11282082432.0);
396
398 15701508.0 *
tau3 - 914128567.0 *
tau2 + 2074956840.0 *
tau -
399 1323431896.0);
400
402 94209048.0 *
tau3 - 1518414297.0 *
tau2 + 2460397220.0 *
tau -
403 889289856.0);
404
406 52338360.0 *
tau3 - 451824525.0 *
tau2 + 687873124.0 *
tau -
407 259006536.0);
408
410 106151040.0 *
tau3 - 661884105.0 *
tau2 +
411 946554244.0 *
tau - 361440756.0);
412
414 8293050.0 *
tau2 - 82437520.0 *
tau + 44764047.0);
415
416 for(
unsigned int i = 0;
i <
N; ++
i)
417 {
418 yOut[
i] = fyIn[
i] + fLastStepLength *
tau * (
421 }
422}
Referenced by Interpolate().
◆ Interpolate5thOrder()
template<class T_Equation, unsigned int N>
Definition at line 481 of file G4TDormandPrince45.hh.
483{
484
485
487
488
490 bi[1][1] = -38039.0 / 7040.0,
491 bi[1][2] = 125923.0 / 10560.0,
492 bi[1][3] = -19683.0 / 1760.0,
493 bi[1][4] = 3303.0 / 880.0,
494
495
496
502
503
504
506 bi[3][1] = -12500.0 / 4081.0,
507 bi[3][2] = 205000.0 / 12243.0,
508 bi[3][3] = -90000.0 / 4081.0,
509 bi[3][4] = 36000.0 / 4081.0,
510
511
512
514 bi[4][1] = -3125.0 / 704.0,
515 bi[4][2] = 25625.0 / 1056.0,
516 bi[4][3] = -5625.0 / 176.0,
517 bi[4][4] = 1125.0 / 88.0,
518
519
520
522 bi[5][1] = 164025.0 / 74624.0,
523 bi[5][2] = -448335.0 / 37312.0,
524 bi[5][3] = 295245.0 / 18656.0,
525 bi[5][4] = -59049.0 / 9328.0,
526
527
528
530 bi[6][1] = -25.0 / 28.0,
531 bi[6][2] = 205.0 / 42.0,
532 bi[6][3] = -45.0 / 7.0,
533 bi[6][4] = 18.0 / 7.0,
534
535
536
538 bi[7][1] = -2.0 / 11.0,
539 bi[7][2] = 73.0 / 55.0,
540 bi[7][3] = -171.0 / 55.0,
541 bi[7][4] = 108.0 / 55.0,
542
543
544
546 bi[8][1] = 189.0 / 22.0,
547 bi[8][2] = -1593.0 / 55.0,
548 bi[8][3] = 3537.0 / 110.0,
549 bi[8][4] = -648.0 / 55.0,
550
551
552
554 bi[9][1] = 351.0 / 110.0,
555 bi[9][2] = -999.0 / 55.0,
556 bi[9][3] = 2943.0 / 110.0,
557 bi[9][4] = -648.0 / 55.0;
558
559
560
561
564
567 {
569 {
571 }
573 }
574
577 {
579 b[1] * fdydxIn[
i] +
b[2] * ak2[
i] +
b[3] * ak3[
i] +
580 b[4] * ak4[
i] +
b[5] * ak5[
i] +
b[6] * ak6[
i] +
581 b[7] * ak7[
i] +
b[8] * ak8[
i] +
b[9] * ak9[
i]
582 );
583 }
584}
◆ RightHandSideInl()
template<class T_Equation, unsigned int N = 6>
◆ SetupInterpolation()
template<class T_Equation, unsigned int N = 6>
◆ SetupInterpolation5thOrder()
template<class T_Equation, unsigned int N>
Definition at line 432 of file G4TDormandPrince45.hh.
433{
434
435
438 b83 = 8875.0 / 103032.0,
439 b84 = -125.0 / 1728.0,
440 b85 = 801.0 / 13568.0,
441 b86 = -13519.0 / 368064.0,
442 b87 = 11105.0 / 368064.0,
443
444 b91 = 632855.0 / 4478976.0,
446 b93 = 4146875.0 / 6491016.0,
447 b94 = 5490625.0 /14183424.0,
448 b95 = -15975.0 / 108544.0,
449 b96 = 8295925.0 / 220286304.0,
450 b97 = -1779595.0 / 62938944.0,
451 b98 = -805.0 / 4104.0;
452
454
455
456
457 for(
unsigned int i = 0;
i <
N; ++
i)
458 {
459 yTemp[
i] = fyIn[
i] + fLastStepLength * (
463 );
464 }
466
467 for(
unsigned int i = 0;
i <
N; ++
i)
468 {
469 yTemp[
i] = fyIn[
i] + fLastStepLength * (
473 );
474 }
476}
void RightHandSideInl(const G4double y[], G4double dydx[])
◆ Stepper() [1/2]
template<class T_Equation, unsigned int N>
The stepper for the Runge Kutta integration. The stepsize is fixed, with the step size given by 'h'. Integrates ODE starting values y[0 to 6]. Outputs yout[] and its estimated error yerr[].
- Parameters
-
| [in] | y | Starting values array of integration variables. |
| [in] | dydx | Derivatives array. |
| [in] | h | The given step size. |
| [out] | yout | Integration output. |
| [out] | yerr | The estimated error. |
Implements G4MagIntegratorStepper.
Definition at line 339 of file G4TDormandPrince45.hh.
344{
347
349}
void StepWithError(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[])
◆ Stepper() [2/2]
template<class T_Equation, unsigned int N = 6>
Definition at line 109 of file G4TDormandPrince45.hh.
113 {
115 }
void StepWithFinalDerivate(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[], G4double dydxOutput[])
◆ StepperType()
template<class T_Equation, unsigned int N = 6>
◆ StepWithError()
template<class T_Equation, unsigned int N>
Definition at line 226 of file G4TDormandPrince45.hh.
231{
232
233
235 b31 = 3.0 / 40.0,
b32 = 9.0 / 40.0,
236 b41 = 44.0 / 45.0,
b42 = -56.0 / 15.0,
b43 = 32.0/9.0,
237
238 b51 = 19372.0 / 6561.0,
b52 = -25360.0 / 2187.0,
b53 = 64448.0 / 6561.0,
239 b54 = -212.0 / 729.0,
240
241 b61 = 9017.0 / 3168.0 ,
b62 = -355.0 / 33.0,
242 b63 = 46732.0 / 5247.0,
b64 = 49.0 / 176.0,
243 b65 = -5103.0 / 18656.0,
244
245 b71 = 35.0 / 384.0,
b72 = 0.,
246 b73 = 500.0 / 1113.0,
b74 = 125.0 / 192.0,
247 b75 = -2187.0 / 6784.0,
b76 = 11.0 / 84.0,
248
249
250
251
252
253
254
255
256
257
258
259
260
261 dc1 = -(
b71 - 5179.0 / 57600.0),
263 dc3 = -(
b73 - 7571.0 / 16695.0),
264 dc4 = -(
b74 - 393.0 / 640.0),
265 dc5 = -(
b75 + 92097.0 / 339200.0),
266 dc6 = -(
b76 - 187.0 / 2100.0),
267 dc7 = -(- 1.0 / 40.0);
268
269
270
272
274
275
276
277 for(
unsigned int i = 0;
i <
N; ++
i)
278 {
281 }
283
284 for(
unsigned int i = 0;
i <
N; ++
i)
285 {
287 }
289
290 for(
unsigned int i = 0;
i <
N; ++
i)
291 {
294 }
296
297 for(
unsigned int i = 0;
i <
N; ++
i)
298 {
301 }
303
304 for(
unsigned int i = 0;
i <
N; ++
i)
305 {
309 }
311
312 for(
unsigned int i = 0;
i <
N; ++
i)
313 {
317 }
319
320 for(
unsigned int i = 0;
i <
N; ++
i)
321 {
326 ) + 1.5e-18;
327
328
329
332 }
333
334 fLastStepLength =
hstep;
335}
Referenced by Stepper(), and StepWithFinalDerivate().
◆ StepWithFinalDerivate()
template<class T_Equation, unsigned int N>
Definition at line 207 of file G4TDormandPrince45.hh.
213{
216}
void copy(G4double dst[], const G4double src[], std::size_t size=G4FieldTrack::ncompSVEC)
Referenced by Stepper().
◆ N8
template<class T_Equation, unsigned int N = 6>
The documentation for this class was generated from the following file: