Integrator Structures

group IntegratorStructs

Enums

enum REB_EOS_TYPE

Available opperator splitting methods for phi0 and phi1 in EOS integrators.

Values:

REB_EOS_LF = 0x00
REB_EOS_LF4 = 0x01
REB_EOS_LF6 = 0x02
REB_EOS_LF8 = 0x03
REB_EOS_LF4_2 = 0x04
REB_EOS_LF8_6_4 = 0x05
REB_EOS_PLF7_6_4 = 0x06
REB_EOS_PMLF4 = 0x07
REB_EOS_PMLF6 = 0x08
struct reb_simulation_integrator_ias15
#include <rebound.h>

This structure contains variables and pointer used by the IAS15 integrator.

Public Members

double epsilon

This parameter controls the accuracy of the integrator.

Set to 0 to make IAS15 a non-adaptive integrator. The default value is: 1e-9.

double min_dt

The minimum allowed timestep.

The default value is 0 (no minimal timestep). Set a finite value to this variable if the IAS15 integrator has problems and the timestep becomes excessively small.

unsigned int epsilon_global

Flag that determines how relative acceleration error is estimated.

If set to 1, estimate the fractional error by max(acceleration_error)/max(acceleration), where max is take over all particles. If set to 0, estimate the fractional error by max(acceleration_error/acceleration).

struct reb_simulation_integrator_mercurius
#include <rebound.h>

This structure contains variables and pointer used by the MERCURIUS integrator.

Public Members

double (*L)(const struct reb_simulation *const r, double d, double dcrit)

This is a function pointer to the force switching function used.

If NULL (the default), the MERCURY switching function will be used. The argument d is the distance between two particles. The argument dcrit is the maximum of the critical distances of the two particles. The return value is a scalar between 0 and 1. If it always returns 1, then the integrator becomes the standard Wisdom-Holman integrator.

double hillfac

Switching distance in units of the Hill radius.

The switching distances for particles are calculated automastically based on multiple criteria. One criteria calculates the Hill radius of particles and then multiplies it with the hillfac variable.

unsigned int recalculate_coordinates_this_timestep

Setting this flag to one will recalculate heliocentric coordinates from the particle structure at the beginning of the next timestep.

After one timestep, the flag gets set back to 0. If you want to change particles after every timestep, you also need to set this flag to 1 before every timestep. Default is 0.

unsigned int recalculate_dcrit_this_timestep

Setting this flag to one will recalculate the critical switchover distances dcrit at the the beginning of the next timestep.

After one timestep, the flag gets set back to 0. If you want to recalculate dcrit at every every timestep, you also need to set this flag to 1 before every timestep. Default is 0.

unsigned int safe_mode

If this flag is set (the default), the integrator will recalculate heliocentric coordinates and synchronize after every timestep, to avoid problems with outputs or particle modifications between timesteps.

Setting it to 0 will result in a speedup, but care must be taken to synchronize and recalculate coordinates when needed.

unsigned int is_synchronized

Flag to determine if current particle structure is synchronized.

unsigned int mode

Internal. 0 if WH is operating, 1 if IAS15 is operating.

unsigned int encounterN

Number of particles currently having an encounter.

unsigned int encounterNactive

Number of particles currently having an encounter.

unsigned int allocatedN

Current size of allocated internal arrays.

unsigned int allocatedN_additionalforces

Current size of allocated internal particles_backup_additionalforces array.

unsigned int dcrit_allocatedN

Current size of dcrit arrays.

double *dcrit

Switching radii for particles.

struct reb_particle* REBOUND_RESTRICT reb_simulation_integrator_mercurius::particles_backup

Internal array, contains coordinates before Kepler step for encounter prediction.

struct reb_particle* REBOUND_RESTRICT reb_simulation_integrator_mercurius::particles_backup_additionalforces

Internal array, contains coordinates before Kepler step for encounter prediction.

int *encounter_map

Map to represent which particles are integrated with ias15.

struct reb_vec3d com_pos

Used internally to keep track of the centre of mass during the timestep.

struct reb_vec3d com_vel

Used internally to keep track of the centre of mass during the timestep.

struct reb_simulation_integrator_sei
#include <rebound.h>

This structure contains variables used by the SEI integrator.

This is where the user sets the orbital frequency OMEGA for shearing sheet simulations.

Public Members

double OMEGA

Epicyclic/orbital frequency.

double OMEGAZ

Epicyclic frequency in vertical direction.

struct reb_simulation_integrator_saba
#include <rebound.h>

This structure contains variables used by the SABA integrator.

Public Types

enum [anonymous]

SABA type.

Available types include: SABA1, SABA2, SABA3, SABA4, SABACM1, SABACM2, SABACM3, SABACM4, SABACL1, SABACL2, SABACL3, SABACL4, SABA(10,4), SABA(8,6,4), SABA(10,6,4), SABAH(8,4,4), SABAH(8,6,4), and SABAH(10,6,4).

Values:

REB_SABA_1 = 0x0
REB_SABA_2 = 0x1
REB_SABA_3 = 0x2
REB_SABA_4 = 0x3
REB_SABA_CM_1 = 0x100
REB_SABA_CM_2 = 0x101
REB_SABA_CM_3 = 0x102
REB_SABA_CM_4 = 0x103
REB_SABA_CL_1 = 0x200
REB_SABA_CL_2 = 0x201
REB_SABA_CL_3 = 0x202
REB_SABA_CL_4 = 0x203
REB_SABA_10_4 = 0x4
REB_SABA_8_6_4 = 0x5
REB_SABA_10_6_4 = 0x6
REB_SABA_H_8_4_4 = 0x7
REB_SABA_H_8_6_4 = 0x8
REB_SABA_H_10_6_4 = 0x9

Public Members

reb_simulation_integrator_saba::[anonymous] type

SABA type.

Available types include: SABA1, SABA2, SABA3, SABA4, SABACM1, SABACM2, SABACM3, SABACM4, SABACL1, SABACL2, SABACL3, SABACL4, SABA(10,4), SABA(8,6,4), SABA(10,6,4), SABAH(8,4,4), SABAH(8,6,4), and SABAH(10,6,4).

unsigned int safe_mode

Safe_mode has the same functionality as in WHFast.

unsigned int is_synchronized

Flag to determine if current particle structure is synchronized.

unsigned int keep_unsynchronized

Flaf that determines if the inertial coordinates generated are discared in subsequent timesteps (Jacobi coordinates are used instead).

Danger zone! Only use this flag if you are absolutely sure what you are doing. This is intended for simulation which have to be reproducible on a bit by bit basis.

struct reb_simulation_integrator_whfast
#include <rebound.h>

This structure contains variables used by the WHFast integrator.

Public Types

enum [anonymous]

This variable determines the kernel of the WHFast integrator.

  • 0 (default): Uses a standard WH kick step

  • 1: uses the exact modified kick (for Newtonian gravity)

  • 2: uses the composition kernel

  • 3: uses the lazy implementer’s modified kick

Values:

REB_WHFAST_KERNEL_DEFAULT = 0
REB_WHFAST_KERNEL_MODIFIEDKICK = 1
REB_WHFAST_KERNEL_COMPOSITION = 2
REB_WHFAST_KERNEL_LAZY = 3
enum [anonymous]

Chooses the coordinate system for the WHFast algorithm. Default is Jacobi Coordinates.

Values:

REB_WHFAST_COORDINATES_JACOBI = 0

Jacobi coordinates (default)

REB_WHFAST_COORDINATES_DEMOCRATICHELIOCENTRIC = 1

Democratic Heliocentric coordinates.

REB_WHFAST_COORDINATES_WHDS = 2

WHDS coordinates (Hernandez and Dehnen, 2017)

Public Members

unsigned int corrector

This variable turns on/off different first symplectic correctors for WHFast.

These correctors remove terms of order O(eps*dt^2)

  • 0 (default): turns off all first correctors

  • 3: uses third order (two-stage) first corrector

  • 5: uses fifth order (four-stage) first corrector

  • 7: uses seventh order (six-stage) first corrector

  • 11: uses eleventh order (ten-stage) first corrector

  • 17: uses 17th order (16-stage) first corrector

unsigned int corrector2

This variable turns on/off the second symplectic correctors for WHFast.

  • 0 (default): turns off second correctors

  • 1: uses second corrector

reb_simulation_integrator_whfast::[anonymous] kernel

This variable determines the kernel of the WHFast integrator.

  • 0 (default): Uses a standard WH kick step

  • 1: uses the exact modified kick (for Newtonian gravity)

  • 2: uses the composition kernel

  • 3: uses the lazy implementer’s modified kick

reb_simulation_integrator_whfast::[anonymous] coordinates

Chooses the coordinate system for the WHFast algorithm. Default is Jacobi Coordinates.

unsigned int recalculate_coordinates_this_timestep

Setting this flag to one will recalculate Jacobi/heliocentric coordinates from the particle structure in the next timestep.

After the timestep, the flag gets set back to 0. If you want to change particles after every timestep, you also need to set this flag to 1 before every timestep. Default is 0.

unsigned int safe_mode

If this flag is set (the default), whfast will recalculate jacobi/heliocentric coordinates and synchronize every timestep, to avoid problems with outputs or particle modifications between timesteps.

Setting it to 0 will result in a speedup, but care must be taken to synchronize and recalculate jacobi coordinates when needed. See AdvWHFast.ipynb in the python_tutorials folder (navigate to it on github if you don’t have ipython notebook installed). The explanation is general, and the python and C flags have the same names.

struct reb_particle* REBOUND_RESTRICT reb_simulation_integrator_whfast::p_jh

Jacobi/heliocentric coordinates.

This array contains the Jacobi/heliocentric coordinates of all particles. It is automatically filled and updated by WHfast. Access this array with caution.

struct reb_particle* REBOUND_RESTRICT reb_simulation_integrator_whfast::p_temp

Internal temporary array used for lazy implementer’s kernel method.

unsigned int keep_unsynchronized

Generate inertial coordinates at the end of the integration, but do not change the Jacobi/heliocentric coordinates.

Danger zone! Only use this flag if you are absolutely sure what you are doing. This is intended for simulation which have to be reproducible on a bit by bit basis.

struct reb_simulation_integrator_eos
#include <rebound.h>

This structure contains variables used by the EOS integrator family.

Public Members

REB_EOS_TYPE phi0

Outer opperator splitting scheme.

REB_EOS_TYPE phi1

Inner opperator splitting scheme.

unsigned int safe_mode

If set to 0, always combine drift steps at the beginning and end of phi0. If set to 1, n needs to be bigger than 1.

unsigned int is_synchronized

Flag to indicate if the drift step at the end of the last timestep has been taken.

struct reb_simulation_integrator_janus

Public Members

double scale_pos

Scale of the problem. Positions get divided by this number before the conversion to an integer.

double scale_vel

Scale of the problem. Velocities get divided by this number before the conversion to an integer.

unsigned int order

Order of the scheme. Default is 6.

unsigned int recalculate_integer_coordinates_this_timestep

If this flag is set, then janus will recalculate integer coordinates at the next timestep.