![]() |
SMILE
v2.5
Schwarzschild Modelling Interactive expLoratory Environment
|
15-th order implicit Runge-Kutta scheme from Rein & Spiegel, 2014, MNRAS (adapted from Rebound). More...
#include <odeint.h>


Public Member Functions | |
| COdeIntegratorIAS15 (COdeSystem *_odeSystem, double _acc) | |
| virtual void | integrateToTime (double timeEnd) |
| perform integration | |
| virtual double | getInterpolatedSolution (unsigned int c, double t) const |
| return interpolated value of c-th variable at time t, where t must lie within current timestep interval; called from various runtime functions that need to obtain x(t) at arbitrary times | |
Public Member Functions inherited from smile::CBasicOdeIntegrator | |
| CBasicOdeIntegrator (COdeSystem *_odeSystem) | |
| double | getTime () const |
| return the time to which the integration has proceeded so far | |
Static Public Member Functions | |
| static const char * | myName () |
Private Member Functions | |
| void | predict_next_step (double ratio, int N3, vectord _e[7], vectord _b[7]) |
| Helper function for resetting the b and e coefficients. | |
| int | integrator_ias15_step () |
| Actual integration for one timestep. | |
Private Attributes | |
| double | integrator_epsilon |
| accuracy for the integrator | |
| unsigned int | nfcn |
| statistics: number of force evaluations | |
| unsigned int | nstep |
| unsigned int | naccpt |
| unsigned int | nrejct |
| statistics: number of timesteps (total, accepted, rejected) | |
| OdeStateType | stateCurr |
| ODE system state at the end of timestep. | |
| OdeStateType | statePrev |
| ODE system state at the beginning of timestep. | |
| OdeStateType | deriv |
| temp.buffer for computing rhs of equations of motion | |
| vectord | at |
| Temp.buffer for accelerations. | |
| vectord | a0 |
| Accelerations at the beginning of timestep. | |
| vectord | compsum |
| Temp.buffer for compensated summation. | |
| vectord | g [7] |
| vectord | b [7] |
| vectord | e [7] |
| coefficients of approximation in the integrator | |
| vectord | br [7] |
| vectord | er [7] |
| coefficients of approximation at the beginning of timestep | |
| double | dt |
| current timestep length | |
| double | dt_last_success |
| Last accepted timestep (corresponding to br and er) | |
Additional Inherited Members | |
Public Types inherited from smile::CBasicOdeIntegrator | |
| enum | STEPPERKIND { SK_DEFAULT, SK_LEAPFROG_NB, SK_DOP853, SK_IAS15, SK_HERMITE, SK_ODEINT_CK5, SK_ODEINT_DP5, SK_ODEINT_BS3, SK_ODEINT_BS, SK_ODEINT_RK4, SK_ODEINT_SYMPL4 } |
| list of various ODE stepper types | |
Protected Attributes inherited from smile::CBasicOdeIntegrator | |
| COdeSystem * | odeSystem |
| double | timePrev |
| double | timeCurr |
15-th order implicit Runge-Kutta scheme from Rein & Spiegel, 2014, MNRAS (adapted from Rebound).
It has its own accuracy parameter, with typical values 10^-4..10^-3 providing relative accuracy in the range 10^-15(almost machine precision)..10^-10. It works only for "Standard Hamiltonian" systems (position+velocity variables).
1.8.8