![]()  | 
  
    SMILE
    v2.5
    
   Schwarzschild Modelling Interactive expLoratory Environment 
   | 
 
8th order Runge-Kutta integrator from Hairer,Norsett&Wanner. More...
#include <odeint.h>


Public Member Functions | |
| COdeIntegratorDOP853 (COdeSystem *_odeSystem, double _accAbs, double _accRel) | |
| virtual void | integrateToTime (double timeEnd) | 
| perform integration  | |
| virtual double | getInterpolatedSolution (unsigned int c, double t) const | 
| computes interpolated coordinate at time t within current timestep  | |
  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 | |
| int | dop853 (double tstart, OdeStateType &y, double tend, double atoler, double rtoler) | 
| driver routine of the DOP853 integrator  More... | |
| void | fcn (unsigned int, double t, const OdeStateType &vars, OdeStateType &derivs) const | 
| shortcut for calling the function that provides the derivative of the ODE system  | |
| double | hinit (unsigned int n, double x, OdeStateType &y, double posneg, OdeStateType &f0, OdeStateType &f1, OdeStateType &yy1, int iord, double hmax, double atoler, double rtoler) | 
| internal integration initialization routine  | |
Private Attributes | |
| double | accAbs | 
| double | accRel | 
| relative and absolute 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)  | |
| double | hout | 
| current timestep  | |
| OdeStateType | rcont1 | 
| temporary storage for dense output interpolation  | |
| OdeStateType | rcont2 | 
| OdeStateType | rcont3 | 
| OdeStateType | rcont4 | 
| OdeStateType | rcont5 | 
| OdeStateType | rcont6 | 
| OdeStateType | rcont7 | 
| OdeStateType | rcont8 | 
| OdeStateType | stateCurr | 
| current values of variables inside the integrator  | |
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 | 
8th order Runge-Kutta integrator from Hairer,Norsett&Wanner.
It is suitable for integrating orbits in a smooth potential (everything except the treecode potential).
      
  | 
  private | 
driver routine of the DOP853 integrator
| tstart | initial t-value | 
| y | initial values for variables | 
| tend | final t-value (tend-tstart may be positive or negative) | 
| atoler | absolute error tolerance | 
| rtoler | relative error tolerance | 
 1.8.8