Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TDormandPrince45< T_Equation, N > Class Template Reference

#include <G4TDormandPrince45.hh>

+ Inheritance diagram for G4TDormandPrince45< T_Equation, N >:

Public Member Functions

 G4TDormandPrince45 (T_Equation *equation)
 
 G4TDormandPrince45 (T_Equation *equation, G4int numVar)
 
void StepWithError (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[])
 
virtual void Stepper (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) override 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
 
virtual G4double DistChord () const override final
 
virtual G4int IntegratorOrder () 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 ()
 
- 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
 

Static Public Attributes

static constexpr int N8 = N > 8 ? N : 8
 

Additional Inherited Members

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

Detailed Description

template<class T_Equation, unsigned int N = 6>
class G4TDormandPrince45< T_Equation, N >

Definition at line 48 of file G4TDormandPrince45.hh.

Constructor & Destructor Documentation

◆ G4TDormandPrince45() [1/2]

template<class T_Equation , unsigned int N>
G4TDormandPrince45< T_Equation, N >::G4TDormandPrince45 ( T_Equation *  equation)

Definition at line 166 of file G4TDormandPrince45.hh.

167 : G4MagIntegratorStepper(dynamic_cast<G4EquationOfMotion*>(equation), N )
168 , fEquation_Rhs(equation)
169{
170 // assert( dynamic_cast<G4EquationOfMotion*>(equation) != nullptr );
171 if( dynamic_cast<G4EquationOfMotion*>(equation) == nullptr )
172 {
173 G4Exception("G4TDormandPrince745: constructor", "GeomField0001",
174 FatalException, "T_Equation is not an G4EquationOfMotion.");
175 }
176
177 /***
178 assert( equation->GetNumberOfVariables == N );
179 if( equation->GetNumberOfVariables != N ){
180 G4ExceptionDescription msg;
181 msg << "Equation has an incompatible number of variables." ;
182 msg << " template N = " << N << " equation-Nvar= "
183 << equation->GetNumberOfVariables;
184 G4Exception("G4TCashKarpRKF45: constructor", "GeomField0001",
185 FatalException, msg );
186 } ****/
187}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35

◆ G4TDormandPrince45() [2/2]

template<class T_Equation , unsigned int N>
G4TDormandPrince45< T_Equation, N >::G4TDormandPrince45 ( T_Equation *  equation,
G4int  numVar 
)

Definition at line 190 of file G4TDormandPrince45.hh.

192{
193 if( numVar != G4int(N)){
195 msg << "Equation has an incompatible number of variables." ;
196 msg << " template N = " << N
197 << " argument numVar = " << numVar ;
198 // << " equation-Nvar= " << equation->GetNumberOfVariables(); // --> Expected later
199 G4Exception("G4TCashKarpRKF45: constructor", "GeomField0001",
201 }
202 assert( numVar == N );
203}
@ FatalErrorInArgument
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
int G4int
Definition: G4Types.hh:85

Member Function Documentation

◆ DistChord()

template<class T_Equation , unsigned int N>
G4double G4TDormandPrince45< T_Equation, N >::DistChord
finaloverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 354 of file G4TDormandPrince45.hh.

355{
356 // Coefficients were taken from Some Practical Runge-Kutta Formulas
357 // by Lawrence F. Shampine, page 149, c*
358 //
359 const G4double hf1 = 6025192743.0 / 30085553152.0,
360 hf3 = 51252292925.0 / 65400821598.0,
361 hf4 = - 2691868925.0 / 45128329728.0,
362 hf5 = 187940372067.0 / 1594534317056.0,
363 hf6 = - 1776094331.0 / 19743644256.0,
364 hf7 = 11237099.0 / 235043384.0;
365
366 G4ThreeVector mid;
367
368 for(unsigned int i = 0; i < 3; ++i)
369 {
370 mid[i] = fyIn[i] + 0.5 * fLastStepLength * (
371 hf1 * fdydxIn[i] + hf3 * ak3[i] +
372 hf4 * ak4[i] + hf5 * ak5[i] + hf6 * ak6[i] + hf7 * ak7[i]);
373 }
374
377
378 return G4LineSection::Distline(mid, begin, end);
379}
double G4double
Definition: G4Types.hh:83
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
G4ThreeVector makeVector(const ArrayType &array, Value3D value)

◆ GetSpecificEquation()

template<class T_Equation , unsigned int N = 6>
T_Equation * G4TDormandPrince45< T_Equation, N >::GetSpecificEquation ( )
inline

Definition at line 111 of file G4TDormandPrince45.hh.

111{ return fEquation_Rhs; }

◆ GetYOut()

template<class T_Equation , unsigned int N = 6>
const field_utils::ShortState< N > & G4TDormandPrince45< T_Equation, N >::GetYOut ( ) const
inline

Definition at line 88 of file G4TDormandPrince45.hh.

88{ return fyOut; }

◆ IntegratorOrder()

template<class T_Equation , unsigned int N = 6>
virtual G4int G4TDormandPrince45< T_Equation, N >::IntegratorOrder ( ) const
inlineoverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 86 of file G4TDormandPrince45.hh.

86{ return 4; }

◆ Interpolate()

template<class T_Equation , unsigned int N = 6>
void G4TDormandPrince45< T_Equation, N >::Interpolate ( G4double  tau,
G4double  yOut[] 
) const
inline

Definition at line 78 of file G4TDormandPrince45.hh.

79 {
80 Interpolate4thOrder(yOut, tau);
81 }
void Interpolate4thOrder(G4double yOut[], G4double tau) const

◆ Interpolate4thOrder()

template<class T_Equation , unsigned int N>
void G4TDormandPrince45< T_Equation, N >::Interpolate4thOrder ( G4double  yOut[],
G4double  tau 
) const

Definition at line 388 of file G4TDormandPrince45.hh.

390{
391 // const G4int numberOfVariables = GetNumberOfVariables();
392
393 const G4double tau2 = tau * tau,
394 tau3 = tau * tau2,
395 tau4 = tau2 * tau2;
396
397 const G4double bf1 = 1.0 / 11282082432.0 * (
398 157015080.0 * tau4 - 13107642775.0 * tau3 + 34969693132.0 * tau2 -
399 32272833064.0 * tau + 11282082432.0);
400
401 const G4double bf3 = - 100.0 / 32700410799.0 * tau * (
402 15701508.0 * tau3 - 914128567.0 * tau2 + 2074956840.0 * tau -
403 1323431896.0);
404
405 const G4double bf4 = 25.0 / 5641041216.0 * tau * (
406 94209048.0 * tau3 - 1518414297.0 * tau2 + 2460397220.0 * tau -
407 889289856.0);
408
409 const G4double bf5 = - 2187.0 / 199316789632.0 * tau * (
410 52338360.0 * tau3 - 451824525.0 * tau2 + 687873124.0 * tau -
411 259006536.0);
412
413 const G4double bf6 = 11.0 / 2467955532.0 * tau * (
414 106151040.0 * tau3 - 661884105.0 * tau2 +
415 946554244.0 * tau - 361440756.0);
416
417 const G4double bf7 = 1.0 / 29380423.0 * tau * (1.0 - tau) * (
418 8293050.0 * tau2 - 82437520.0 * tau + 44764047.0);
419
420 for(unsigned int i = 0; i < N; ++i)
421 {
422 yOut[i] = fyIn[i] + fLastStepLength * tau * (
423 bf1 * fdydxIn[i] + bf3 * ak3[i] + bf4 * ak4[i] +
424 bf5 * ak5[i] + bf6 * ak6[i] + bf7 * ak7[i]);
425 }
426}

Referenced by G4TDormandPrince45< T_Equation, N >::Interpolate().

◆ Interpolate5thOrder()

template<class T_Equation , unsigned int N>
void G4TDormandPrince45< T_Equation, N >::Interpolate5thOrder ( G4double  yOut[],
G4double  tau 
) const

Definition at line 486 of file G4TDormandPrince45.hh.

488{
489 // Define the coefficients for the polynomials
490 //
491 G4double bi[10][5];
492
493 // COEFFICIENTS OF bi[1]
494 bi[1][0] = 1.0,
495 bi[1][1] = -38039.0 / 7040.0,
496 bi[1][2] = 125923.0 / 10560.0,
497 bi[1][3] = -19683.0 / 1760.0,
498 bi[1][4] = 3303.0 / 880.0,
499 // --------------------------------------------------------
500 //
501 // COEFFICIENTS OF bi[2]
502 bi[2][0] = 0.0,
503 bi[2][1] = 0.0,
504 bi[2][2] = 0.0,
505 bi[2][3] = 0.0,
506 bi[2][4] = 0.0,
507 // --------------------------------------------------------
508 //
509 // COEFFICIENTS OF bi[3]
510 bi[3][0] = 0.0,
511 bi[3][1] = -12500.0 / 4081.0,
512 bi[3][2] = 205000.0 / 12243.0,
513 bi[3][3] = -90000.0 / 4081.0,
514 bi[3][4] = 36000.0 / 4081.0,
515 // --------------------------------------------------------
516 //
517 // COEFFICIENTS OF bi[4]
518 bi[4][0] = 0.0,
519 bi[4][1] = -3125.0 / 704.0,
520 bi[4][2] = 25625.0 / 1056.0,
521 bi[4][3] = -5625.0 / 176.0,
522 bi[4][4] = 1125.0 / 88.0,
523 // --------------------------------------------------------
524 //
525 // COEFFICIENTS OF bi[5]
526 bi[5][0] = 0.0,
527 bi[5][1] = 164025.0 / 74624.0,
528 bi[5][2] = -448335.0 / 37312.0,
529 bi[5][3] = 295245.0 / 18656.0,
530 bi[5][4] = -59049.0 / 9328.0,
531 // --------------------------------------------------------
532 //
533 // COEFFICIENTS OF bi[6]
534 bi[6][0] = 0.0,
535 bi[6][1] = -25.0 / 28.0,
536 bi[6][2] = 205.0 / 42.0,
537 bi[6][3] = -45.0 / 7.0,
538 bi[6][4] = 18.0 / 7.0,
539 // --------------------------------------------------------
540 //
541 // COEFFICIENTS OF bi[7]
542 bi[7][0] = 0.0,
543 bi[7][1] = -2.0 / 11.0,
544 bi[7][2] = 73.0 / 55.0,
545 bi[7][3] = -171.0 / 55.0,
546 bi[7][4] = 108.0 / 55.0,
547 // --------------------------------------------------------
548 //
549 // COEFFICIENTS OF bi[8]
550 bi[8][0] = 0.0,
551 bi[8][1] = 189.0 / 22.0,
552 bi[8][2] = -1593.0 / 55.0,
553 bi[8][3] = 3537.0 / 110.0,
554 bi[8][4] = -648.0 / 55.0,
555 // --------------------------------------------------------
556 //
557 // COEFFICIENTS OF bi[9]
558 bi[9][0] = 0.0,
559 bi[9][1] = 351.0 / 110.0,
560 bi[9][2] = -999.0 / 55.0,
561 bi[9][3] = 2943.0 / 110.0,
562 bi[9][4] = -648.0 / 55.0;
563 // --------------------------------------------------------
564
565 // Calculating the polynomials
566
567 G4double b[10];
568 std::memset(b, 0.0, sizeof(b));
569
570 G4double tauPower = 1.0;
571 for(G4int j = 0; j <= 4; ++j)
572 {
573 for(G4int iStage = 1; iStage <= 9; ++iStage)
574 {
575 b[iStage] += bi[iStage][j] * tauPower;
576 }
577 tauPower *= tau;
578 }
579
580 // const G4int numberOfVariables = GetNumberOfVariables();
581 const G4double stepLen = fLastStepLength * tau;
582 for(G4int i = 0; i < N; ++i)
583 {
584 yOut[i] = fyIn[i] + stepLen * (
585 b[1] * fdydxIn[i] + b[2] * ak2[i] + b[3] * ak3[i] +
586 b[4] * ak4[i] + b[5] * ak5[i] + b[6] * ak6[i] +
587 b[7] * ak7[i] + b[8] * ak8[i] + b[9] * ak9[i]
588 );
589 }
590}

◆ RightHandSideInl()

template<class T_Equation , unsigned int N = 6>
void G4TDormandPrince45< T_Equation, N >::RightHandSideInl ( const G4double  y[],
G4double  dydx[] 
)
inline

Definition at line 96 of file G4TDormandPrince45.hh.

98 {
99 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
100 }

◆ SetupInterpolation()

template<class T_Equation , unsigned int N = 6>
void G4TDormandPrince45< T_Equation, N >::SetupInterpolation ( )
inline

Definition at line 76 of file G4TDormandPrince45.hh.

76{}

◆ SetupInterpolation5thOrder()

template<class T_Equation , unsigned int N>
void G4TDormandPrince45< T_Equation, N >::SetupInterpolation5thOrder

Definition at line 436 of file G4TDormandPrince45.hh.

437{
438 // Coefficients for the additional stages
439 //
440 const G4double b81 = 6245.0 / 62208.0,
441 b82 = 0.0,
442 b83 = 8875.0 / 103032.0,
443 b84 = -125.0 / 1728.0,
444 b85 = 801.0 / 13568.0,
445 b86 = -13519.0 / 368064.0,
446 b87 = 11105.0 / 368064.0,
447
448 b91 = 632855.0 / 4478976.0,
449 b92 = 0.0,
450 b93 = 4146875.0 / 6491016.0,
451 b94 = 5490625.0 /14183424.0,
452 b95 = -15975.0 / 108544.0,
453 b96 = 8295925.0 / 220286304.0,
454 b97 = -1779595.0 / 62938944.0,
455 b98 = -805.0 / 4104.0;
456
457 // const G4int numberOfVariables = GetNumberOfVariables();
459
460 // Evaluate the extra stages
461 //
462 for(unsigned int i = 0; i < N; ++i)
463 {
464 yTemp[i] = fyIn[i] + fLastStepLength * (
465 b81 * fdydxIn[i] + b82 * ak2[i] + b83 * ak3[i] +
466 b84 * ak4[i] + b85 * ak5[i] + b86 * ak6[i] +
467 b87 * ak7[i]
468 );
469 }
470 RightHandSideInl(yTemp, ak8); // 8th Stage
471
472 for(unsigned int i = 0; i < N; ++i)
473 {
474 yTemp[i] = fyIn[i] + fLastStepLength * (
475 b91 * fdydxIn[i] + b92 * ak2[i] + b93 * ak3[i] +
476 b94 * ak4[i] + b95 * ak5[i] + b96 * ak6[i] +
477 b97 * ak7[i] + b98 * ak8[i]
478 );
479 }
480 RightHandSideInl(yTemp, ak9); // 9th Stage
481}
void RightHandSideInl(const G4double y[], G4double dydx[])
G4double[N] ShortState
Definition: G4FieldUtils.hh:48

◆ Stepper() [1/2]

template<class T_Equation , unsigned int N>
void G4TDormandPrince45< T_Equation, N >::Stepper ( const G4double  yInput[],
const G4double  dydx[],
G4double  hstep,
G4double  yOutput[],
G4double  yError[] 
)
inlinefinaloverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 340 of file G4TDormandPrince45.hh.

345{
346 assert( yOutput != yInput );
347 assert( yError != yInput );
348
349 StepWithError( yInput, dydx, Step, yOutput, yError);
350}
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>
void G4TDormandPrince45< T_Equation, N >::Stepper ( const G4double  yInput[],
const G4double  dydx[],
G4double  hstep,
G4double  yOutput[],
G4double  yError[],
G4double  dydxOutput[] 
)
inline

Definition at line 103 of file G4TDormandPrince45.hh.

106 {
107 StepWithFinalDerivate(yInput, dydx, hstep,
108 yOutput, yError, dydxOutput);
109 }
void StepWithFinalDerivate(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[], G4double dydxOutput[])

◆ StepWithError()

template<class T_Equation , unsigned int N>
void G4TDormandPrince45< T_Equation, N >::StepWithError ( const G4double  yInput[],
const G4double  dydx[],
G4double  hstep,
G4double  yOutput[],
G4double  yError[] 
)
inline

Definition at line 226 of file G4TDormandPrince45.hh.

231{
232 // The parameters of the Butcher tableu
233 //
234 constexpr G4double b21 = 0.2,
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 //Sum of columns, sum(bij) = ei
250 // e1 = 0. ,
251 // e2 = 1.0/5.0 ,
252 // e3 = 3.0/10.0 ,
253 // e4 = 4.0/5.0 ,
254 // e5 = 8.0/9.0 ,
255 // e6 = 1.0 ,
256 // e7 = 1.0 ,
257
258 // Difference between the higher and the lower order method coeff. :
259 // b7j are the coefficients of higher order
260
261 dc1 = -(b71 - 5179.0 / 57600.0),
262 dc2 = -(b72 - .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 // const G4int numberOfVariables = GetNumberOfVariables();
270 // The number of variables to be integrated over
272
273 yOut[7] = yTemp[7] = fyIn[7] = yInput[7]; // Pass along the time - used in RightHandSide
274
275 // Saving yInput because yInput and yOut can be aliases for same array
276 //
277 for(unsigned int i = 0; i < N; ++i)
278 {
279 fyIn[i] = yInput[i];
280 yTemp[i] = yInput[i] + b21 * hstep * dydx[i];
281 }
282 RightHandSideInl(yTemp, ak2); // 2nd stage
283
284 for(unsigned int i = 0; i < N; ++i)
285 {
286 yTemp[i] = fyIn[i] + hstep * (b31 * dydx[i] + b32 * ak2[i]);
287 }
288 RightHandSideInl(yTemp, ak3); // 3rd stage
289
290 for(unsigned int i = 0; i < N; ++i)
291 {
292 yTemp[i] = fyIn[i] + hstep * (
293 b41 * dydx[i] + b42 * ak2[i] + b43 * ak3[i]);
294 }
295 RightHandSideInl(yTemp, ak4); // 4th stage
296
297 for(unsigned int i = 0; i < N; ++i)
298 {
299 yTemp[i] = fyIn[i] + hstep * (
300 b51 * dydx[i] + b52 * ak2[i] + b53 * ak3[i] + b54 * ak4[i]);
301 }
302 RightHandSideInl(yTemp, ak5); // 5th stage
303
304 for(unsigned int i = 0; i < N; ++i)
305 {
306 yTemp[i] = fyIn[i] + hstep * (
307 b61 * dydx[i] + b62 * ak2[i] +
308 b63 * ak3[i] + b64 * ak4[i] + b65 * ak5[i]);
309 }
310 RightHandSideInl(yTemp, ak6); // 6th stage
311
312 for(unsigned int i = 0; i < N; ++i)
313 {
314 yOut[i] = fyIn[i] + hstep * (
315 b71 * dydx[i] + b72 * ak2[i] + b73 * ak3[i] +
316 b74 * ak4[i] + b75 * ak5[i] + b76 * ak6[i]);
317 }
318 RightHandSideInl(yOut, ak7); // 7th and Final stage
319
320 for(unsigned int i = 0; i < N; ++i)
321 {
322 yErr[i] = hstep * (
323 dc1 * dydx[i] + dc2 * ak2[i] +
324 dc3 * ak3[i] + dc4 * ak4[i] +
325 dc5 * ak5[i] + dc6 * ak6[i] + dc7 * ak7[i]
326 ) + 1.5e-18;
327
328 // Store Input and Final values, for possible use in calculating chord
329 //
330 fyOut[i] = yOut[i];
331 fdydxIn[i] = dydx[i];
332 }
333
334 fLastStepLength = hstep;
335}

◆ StepWithFinalDerivate()

template<class T_Equation , unsigned int N>
void G4TDormandPrince45< T_Equation, N >::StepWithFinalDerivate ( const G4double  yInput[],
const G4double  dydx[],
G4double  hstep,
G4double  yOutput[],
G4double  yError[],
G4double  dydxOutput[] 
)
inline

Definition at line 206 of file G4TDormandPrince45.hh.

213{
214 StepWithError(yInput, dydx, hstep, yOutput, yError);
215 field_utils::copy(dydxOutput, ak7, N);
216}
void copy(G4double dst[], const G4double src[], size_t size=G4FieldTrack::ncompSVEC)
Definition: G4FieldUtils.cc:98

Referenced by G4TDormandPrince45< T_Equation, N >::Stepper().

Member Data Documentation

◆ N8

template<class T_Equation , unsigned int N = 6>
constexpr int G4TDormandPrince45< T_Equation, N >::N8 = N > 8 ? N : 8
staticconstexpr

Definition at line 113 of file G4TDormandPrince45.hh.


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