4. API Documentation (C)

This documents the C API for REBOUNDx. The main reference point in the documentation is Implemented Effects ,which has descriptions for each effect and its parameters, citations, and links to examples.

REBOUNDx API definition.

Author

Dan Tamayo tamayo.daniel@gmail.com, Hanno Rein

4.1. LICENSE

Copyright (c) 2019 Dan Tamayo, Hanno Rein

This file is part of reboundx.

reboundx is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

reboundx is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with rebound. If not, see http://www.gnu.org/licenses/.

Main REBOUNDx Functions

struct rebx_extras *rebx_attach(struct reb_simulation *sim)

Adds REBOUNDx functionality to a passed REBOUND simulation.

Parameters:

sim – Pointer to the reb_simulation on which to add REBOUNDx functionality.

Returns:

Pointer to a rebx_extras structure.

void rebx_detach(struct reb_simulation *sim, struct rebx_extras *rebx)

Detaches REBOUNDx from simulation, resetting all the simulation’s function pointers that REBOUNDx has set.

This does not free the memory allocated by REBOUNDx (call rebx_free).

Parameters:

sim – Pointer to the simulation from which to remove REBOUNDx

void rebx_extras_cleanup(struct reb_simulation *sim)
void rebx_free(struct rebx_extras *rebx)

Frees all memory allocated by REBOUNDx instance.

Should be called after simulation is done if memory is a concern.

Parameters:

rebx – The rebx_extras pointer returned from the initial call to rebx_attach.

int rebx_remove_force(struct rebx_extras *rebx, struct rebx_force *force)
int rebx_remove_operator(struct rebx_extras *rebx, struct rebx_operator *operator)
void rebx_output_binary(struct rebx_extras *rebx, char *filename)

Save a binary file with all the effects in the simulation, as well as all particle and effect parameters.

Parameters:
  • rebx – Pointer to the rebx_extras instance

  • filename – Filename to which to save the binary file.

struct rebx_extras *rebx_create_extras_from_binary(struct reb_simulation *sim, const char *const filename)

Reads a REBOUNDx binary file, loads all effects and parameters.

Parameters:
  • sim – Pointer to the simulation to which the effects and parameters should be added.

  • filename – Filename of the saved binary file.

void rebx_init_extras_from_binary(struct rebx_extras *rebx, const char *const filename, enum rebx_input_binary_messages *warnings)

Similar to rebx_create_extras_from_binary(), but takes an extras instance (must be attached to a simulation) and allows for manual message handling.

Parameters:
  • rebx – Pointer to a rebx_extras instance to be updated.

  • filename – Filename of the saved binary file.

  • warnings – Pointer to an array of warnings to be populated during loading.

Functions for manipulating effects in REBOUNDx

int rebx_add_operator(struct rebx_extras *rebx, struct rebx_operator *operator)

Main function for adding effects in REBOUNDx.

Parameters:
  • rebx – Pointer to the rebx_extras instance

  • name – Name of the effect we want to add

Returns:

Returns a pointer to a rebx_effect structure for the effect.

int rebx_add_operator_step(struct rebx_extras *rebx, struct rebx_operator *operator, const double dt_fraction, enum rebx_timing timing)
int rebx_add_force(struct rebx_extras *rebx, struct rebx_force *force)
struct rebx_operator *rebx_load_operator(struct rebx_extras *const rebx, const char *name)
struct rebx_force *rebx_load_force(struct rebx_extras *const rebx, const char *name)
struct rebx_operator *rebx_create_operator(struct rebx_extras *const rebx, const char *name)
struct rebx_force *rebx_create_force(struct rebx_extras *const rebx, const char *name)
struct rebx_effect *rebx_add_custom_force(struct rebx_extras *rebx, const char *name, void (*custom_force)(struct reb_simulation *const sim, struct rebx_effect *const effect, struct reb_particle *const particles, const int N), const int force_is_velocity_dependent)

Function for adding a custom force in REBOUNDx.

Parameters:
  • rebx – Pointer to the rebx_extras instance

  • name – String with the name of the custom effect

  • custom_force – User-implemented function that updates the accelerations of particles.

  • force_is_velocity_dependent – Should be set to 1 if the custom force uses particle velocities, 0 otherwise

Returns:

Returns a pointer to a rebx_effect structure for the effect

struct rebx_effect *rebx_add_custom_operator(struct rebx_extras *rebx, const char *name, void (*custom_operator)(struct reb_simulation *const sim, struct rebx_effect *const effect, const double dt, enum rebx_timing timing))

Function for adding a custom post_timestep_modification in REBOUNDx.

Parameters:
  • rebx – Pointer to the rebx_extras instance

  • name – String with the name of the custom effect

  • custom_ptm – User-implemented function that updates particles.

Returns:

Returns a pointer to a rebx_effect structure for the effect.

struct rebx_force *rebx_get_force(struct rebx_extras *const rebx, const char *const name)

Get a pointer to a force by name.

Parameters:
  • rebx – Pointer to the rebx_extras instance

  • effect_name – Name of the force (string)

Returns:

Pointer to the corresponding rebx_force structure, or NULL if not found.

struct rebx_operator *rebx_get_operator(struct rebx_extras *const rebx, const char *const name)

Get a pointer to an operator by name.

Parameters:
  • rebx – Pointer to the rebx_extras instance

  • effect_name – Name of the operator (string)

Returns:

Pointer to the corresponding rebx_operator structure, or NULL if not found.

Functions for accessing and modifying particle and effect parameters.

int rebx_remove_param(struct rebx_node **apptr, const char *const param_name)

Removes a parameter from a particle or effect.

Parameters:
  • object – Pointer to the particle or effect we want to remove a parameter from.

  • param_name – Name of the parameter we want to remove.

Returns:

1 if parameter found and successfully removed, 0 otherwise.

void *rebx_get_param(struct rebx_extras *const rebx, struct rebx_node *ap, const char *const param_name)

Gets a parameter from a particle or effect.

Parameters:
Returns:

A void pointer to the parameter. NULL if not found.

struct rebx_param *rebx_get_param_struct(struct rebx_extras *const rebx, struct rebx_node *ap, const char *const param_name)
void rebx_set_param_pointer(struct rebx_extras *const rebx, struct rebx_node **apptr, const char *const param_name, void *val)
void rebx_set_param_double(struct rebx_extras *const rebx, struct rebx_node **apptr, const char *const param_name, double val)
void rebx_set_param_int(struct rebx_extras *const rebx, struct rebx_node **apptr, const char *const param_name, int val)
void rebx_set_param_uint32(struct rebx_extras *const rebx, struct rebx_node **apptr, const char *const param_name, uint32_t val)
void rebx_set_param_vec3d(struct rebx_extras *const rebx, struct rebx_node **apptr, const char *const param_name, struct reb_vec3d val)
void rebx_register_param(struct rebx_extras *const rebx, const char *name, enum rebx_param_type type)

General utility functions

struct reb_vec3d rebx_tools_spin_angular_momentum(struct rebx_extras *const rebx)

Calculate spin angular momentum in the simulation of any bodies with spin parameters set (moment of inertia I and angular rotation frequency vector Omega).

Parameters:

rebx – Pointer to the rebx_extras instance

void rebx_simulation_irotate(struct rebx_extras *const rebx, const struct reb_rotation q)

Convenience Functions for Various Effects

double rebx_rad_calc_beta(const double G, const double c, const double source_mass, const double source_luminosity, const double radius, const double density, const double Q_pr)

Calculates beta, the ratio between the radiation pressure force and the gravitational force from the star.

Parameters:
  • G – Gravitational constant.

  • c – Speed of light.

  • source_mass – Mass of the source body.

  • source_luminosity – Luminosity of radiation source.

  • radius – Particle physical radius.

  • density – density of particle.

  • Q_pr – Radiation pressure coefficient (Burns et al. 1979).

Returns:

Beta parameter (double).

double rebx_rad_calc_particle_radius(const double G, const double c, const double source_mass, const double source_luminosity, const double beta, const double density, const double Q_pr)

Calculates the particle radius from physical parameters and beta, the ratio of radiation to gravitational forces from the star.

Parameters:
  • G – Gravitational constant.

  • c – Speed of light.

  • source_mass – Mass of the source body.

  • source_luminosity – Luminosity of radiation source.

  • beta – ratio of radiation force to gravitational force from the radiation source body.

  • density – density of particle.

  • Q_pr – Radiation pressure coefficient (Burns et al. 1979).

Returns:

Particle radius (double).

void rebx_spin_initialize_ode(struct rebx_extras *const rebx, struct rebx_force *const effect)

Count how many particles have their moment of inertia and spin set, and initialize corresponding spin ODEs.

Must be called after setting the moment of inertia and spin of all particles you want to evolve. Attaches spin_ode param to passed effect struct.

Parameters:
  • rebx – Pointer to the rebx_extras instance.

  • effect – (rebx_force) Force structure to which to attach spin_ode object.

double rebx_central_force_Acentral(const struct reb_particle p, const struct reb_particle primary, const double pomegadot, const double gamma)

Calculates the Aradial parameter for central_force effect required for a particle to have a particular pericenter precession rate.

Parameters:
  • p – Particle whose pericenter precession rate we want to match.

  • primary – Central particle for the central force (to which we add the Acentral and gammacentral parameters).

  • pomegadot – Pericenter precession rate we want to obtain.

  • gamma – Index of the central force law.

Returns:

Acentral Normalization to add to the central particle.

double rebx_gr_hamiltonian(struct rebx_extras *const rebx, const struct rebx_force *const gr)

Calculates the hamiltonian for gr, including the classical Hamiltonian.

Assumes there is only one source particle (with gr_source set to 1)

Parameters:
  • rebx – pointer to the REBOUNDx extras instance.

  • gr – Force structure returned by rebx_load_force

Returns:

Total Hamiltonian, including classical Hamiltonian (double).

double rebx_gr_full_hamiltonian(struct rebx_extras *const rebx, const struct rebx_force *const gr_full)

Calculates the hamiltonian for gr_full, including the classical Hamiltonian.

Parameters:
  • rebx – pointer to the REBOUNDx extras instance.

  • gr_full – Force structure returned by rebx_load_force

Returns:

Total Hamiltonian, including classical Hamiltonian (double).

double rebx_gr_potential_potential(struct rebx_extras *const rebx, const struct rebx_force *const gr_potential)

Calculates the potential for gr_potential, including the classical Hamiltonian.

Parameters:
  • rebx – pointer to the REBOUNDx extras instance.

  • gr_potential – Force structure returned by rebx_load_force

Returns:

Total Hamiltonian, including classical Hamiltonian (double).

double rebx_tides_constant_time_lag_potential(struct rebx_extras *const rebx)

Calculates the potential for the conservative piece of the tides_constant_time_lag effect.

Will be conserved if tctl_tau = 0

Parameters:

rebx – pointer to the REBOUNDx extras instance.

Returns:

Potential corresponding to tides_constant_time_lag effect.

double rebx_tides_spin_energy(struct rebx_extras *const rebx)

Calculates the total energy associated with bodies’ spin and their tidal interactions in tides_spin effect. This includes the spin kinetic energy plus gravitational potential between the tidally and rotationally induced quadrupoles and other bodies (treated as point masses). You should add sim.energy() to include bodies’ overall kinetic energy and point-mass gravitational potential energy.

Parameters:

rebx – pointer to the REBOUNDx extras instance.

Returns:

Energy associated with tides_spin effect.

double rebx_central_force_potential(struct rebx_extras *const rebx)

Calculates the potential for central_force effect.

Parameters:

rebx – pointer to the REBOUNDx extras instance.

Returns:

Potential corresponding to central_force effect.

double rebx_gravitational_harmonics_potential(struct rebx_extras *const rebx)

Calculates the potential for all particles with additional gravity field harmonics beyond the monopole (i.e., J2, J4).

Parameters:

rebx – pointer to the REBOUNDx extras instance.

Returns:

Potential corresponding to the effect from all particles of their additional gravity field harmonics

Specialized functions for accessing and modifying particle and effect parameters.

void *rebx_get_param_check(struct reb_simulation *sim, struct rebx_node *ap, const char *const param_name, enum rebx_param_type param_type)

Gets the full rebx_param structure for a particular parameter, rather than just the pointer to the contents.

This can be useful to check properties of the parameter, like the param_type or shape.

Returns a void pointer to the parameter just like rebx_get_param, but additionally checks that the param_type matches what is expected.

Effects should use this function rather than rebx_get_param to ensure that the user appropriately set parameters if working from Python.

Parameters:
  • object – Pointer to the particle or effect that holds the parameter.

  • param_name – Name of the parameter we want to get (see Effects page at http://reboundx.readthedocs.org)

  • object – Pointer to the particle or effect that holds the parameter.

  • param_name – Name of the parameter we want to get (see Effects page at http://reboundx.readthedocs.org)

Returns:

Pointer to the rebx_param structure that holds the parameter. NULL if not found.

Returns:

Void pointer to the parameter. NULL if not found or type does not match (will write error to stderr).

struct rebx_param *rebx_get_or_add_param(struct rebx_extras *const rebx, struct rebx_node **apptr, const char *const param_name)

REBOUND Stepper Functions

void rebx_ias15_step(struct reb_simulation *const sim, struct rebx_operator *const operator, const double dt)

Executes a step for passed time dt using the IAS15 integrator in REBOUND.

Parameters:
  • sim – Pointer to the simulation to step.

  • operator – Unused pointer (kept for consistency with other operators). Can pass NULL.

  • dt – timestep for which to step in simulation time units.

void rebx_kepler_step(struct reb_simulation *const sim, struct rebx_operator *const operator, const double dt)

Executes a Kepler step for passed time dt using the WHFast integrator in REBOUND.

Will use the coordinates and other options set in sim.ri_whfast

Parameters:
  • sim – Pointer to the simulation to step.

  • operator – Unused pointer (kept for consistency with other operators). Can pass NULL.

  • dt – timestep for which to step in simulation time units.

void rebx_jump_step(struct reb_simulation *const sim, struct rebx_operator *const operator, const double dt)

Executes a jump step for passed time dt using the WHFast integrator in REBOUND.

Will use the coordinates and other options set in sim.ri_whfast

Parameters:
  • sim – Pointer to the simulation to step.

  • operator – Unused pointer (kept for consistency with other operators). Can pass NULL.

  • dt – timestep for which to step in simulation time units.

void rebx_interaction_step(struct reb_simulation *const sim, struct rebx_operator *const operator, const double dt)

Executes an interaction step for passed time dt using the WHFast integrator in REBOUND.

Will use the coordinates and other options set in sim.ri_whfast

Parameters:
  • sim – Pointer to the simulation to step.

  • operator – Unused pointer (kept for consistency with other operators). Can pass NULL.

  • dt – timestep for which to step in simulation time units.

void rebx_drift_step(struct reb_simulation *const sim, struct rebx_operator *const operator, const double dt)

Executes a drift step for passed time dt using the leapfrog integrator in REBOUND.

Parameters:
  • sim – Pointer to the simulation to step.

  • operator – Unused pointer (kept for consistency with other operators). Can pass NULL.

  • dt – timestep for which to step in simulation time units.

void rebx_kick_step(struct reb_simulation *const sim, struct rebx_operator *const operator, const double dt)

Executes a kick step for passed time dt using the leapfrog integrator in REBOUND.

Parameters:
  • sim – Pointer to the simulation to step.

  • operator – Unused pointer (kept for consistency with other operators). Can pass NULL.

  • dt – timestep for which to step in simulation time units.

Interpolation Routines

struct rebx_interpolator *rebx_create_interpolator(struct rebx_extras *const rebx, const int Nvalues, const double *times, const double *values, enum rebx_interpolation_type interpolation)

Takes an array of times and corresponding array of values and returns a structure that allows interpolation of values at arbitrary times.

See the parameter interpolation examples in C and Python.

Parameters:
  • rebx – pointer to the REBOUNDx extras instance.

  • Nvalues – Length of times and values arrays (must be equal for both).

  • times – Array of times at which the corresponding values are supplied.

  • values – Array of values at each corresponding time.

  • interpolation – Enum specifying the interpolation method. Defaults to spline.

Returns:

Pointer to a rebx_interpolator structure. Call rebx_interpolate to get values.

void rebx_free_interpolator(struct rebx_interpolator *const interpolator)

Frees the memory for a rebx_interpolator structure.

double rebx_interpolate(struct rebx_extras *const rebx, struct rebx_interpolator *const interpolator, const double time)

Interpolate value at arbitrary times.

Need to first rebx_create_interpolator with an array of times and corresponding values to interpolate between. See parameter interpolation examples.

Parameters:
  • rebx – Pointer to the REBOUNDx extras instance.

  • interpolator – Pointer to the rebx_interpolator structure to interpolate from.

  • time – Time at which to interpolate value.

Returns:

Interpolated value at passed time.

Testing Functions

FILE *rebx_input_inspect_binary(const char *const filename, enum rebx_input_binary_messages *warnings)

Loads a binary file, reads the header, and gives back the file pointer for manual reading.

Parameters:
  • filename – File to open

  • warnings – Pointer to warnings enum to store warnings that come up

Returns:

Returns a pointer to the binary file at the position following the header

struct rebx_binary_field rebx_input_read_binary_field(FILE *inf)

Read the next field in a binary file.

Parameters:

inf – Pointer to the input file

Returns:

Returns rebx_binary_field struct, initialized to 0 if read fails

void rebx_input_skip_binary_field(FILE *inf, long field_size)

Skip forward in binary file.

Parameters:
  • inf – Pointer to the input file

  • field_size – Length by which to skip from current file position

Defines

M_PI
REBXGITHASH

Enums

enum rebx_param_type

Data types available in REBOUNDx for parameters.

Values:

enumerator REBX_TYPE_NONE
enumerator REBX_TYPE_DOUBLE
enumerator REBX_TYPE_INT
enumerator REBX_TYPE_POINTER
enumerator REBX_TYPE_FORCE
enumerator REBX_TYPE_UINT32
enumerator REBX_TYPE_ORBIT
enumerator REBX_TYPE_ODE
enumerator REBX_TYPE_VEC3D
enum REBX_COORDINATES

Different coordinate systems.

Values:

enumerator REBX_COORDINATES_JACOBI

Jacobi coordinates (default)

enumerator REBX_COORDINATES_BARYCENTRIC

Coordinates referenced to pos/vel of system’s center of mass.

enumerator REBX_COORDINATES_PARTICLE

Coordinates relative to pos/vel of a particular particle.

enum rebx_timing

Flag for whether steps should happen before or after the timestep.

Values:

enumerator REBX_TIMING_PRE

Pre timestep.

enumerator REBX_TIMING_POST

Post timestep.

enum rebx_force_type

Force types.

Values:

enumerator REBX_FORCE_NONE

Uninitialized default.

enumerator REBX_FORCE_POS

Force derivable from a position-dependent potential.

enumerator REBX_FORCE_VEL

velocity (or pos and vel) dependent force

enum rebx_operator_type

Operator types.

Values:

enumerator REBX_OPERATOR_NONE

Uninitialized default.

enumerator REBX_OPERATOR_UPDATER

operator that modifies x,v or m,

enumerator REBX_OPERATOR_RECORDER

operator that leaves state unchanged. Just records

enum rebx_binary_field_type

Different fields for binary files.

Values:

enumerator REBX_BINARY_FIELD_TYPE_NONE
enumerator REBX_BINARY_FIELD_TYPE_OPERATOR
enumerator REBX_BINARY_FIELD_TYPE_PARTICLE
enumerator REBX_BINARY_FIELD_TYPE_REBX_STRUCTURE
enumerator REBX_BINARY_FIELD_TYPE_PARAM
enumerator REBX_BINARY_FIELD_TYPE_NAME
enumerator REBX_BINARY_FIELD_TYPE_PARAM_TYPE
enumerator REBX_BINARY_FIELD_TYPE_PARAM_VALUE
enumerator REBX_BINARY_FIELD_TYPE_END
enumerator REBX_BINARY_FIELD_TYPE_PARTICLE_INDEX
enumerator REBX_BINARY_FIELD_TYPE_REBX_INTEGRATOR
enumerator REBX_BINARY_FIELD_TYPE_FORCE_TYPE
enumerator REBX_BINARY_FIELD_TYPE_OPERATOR_TYPE
enumerator REBX_BINARY_FIELD_TYPE_STEP
enumerator REBX_BINARY_FIELD_TYPE_STEP_DT_FRACTION
enumerator REBX_BINARY_FIELD_TYPE_REGISTERED_PARAM
enumerator REBX_BINARY_FIELD_TYPE_ADDITIONAL_FORCE
enumerator REBX_BINARY_FIELD_TYPE_PARAM_LIST
enumerator REBX_BINARY_FIELD_TYPE_REGISTERED_PARAMETERS
enumerator REBX_BINARY_FIELD_TYPE_ALLOCATED_FORCES
enumerator REBX_BINARY_FIELD_TYPE_ALLOCATED_OPERATORS
enumerator REBX_BINARY_FIELD_TYPE_ADDITIONAL_FORCES
enumerator REBX_BINARY_FIELD_TYPE_PRE_TIMESTEP_MODIFICATIONS
enumerator REBX_BINARY_FIELD_TYPE_POST_TIMESTEP_MODIFICATIONS
enumerator REBX_BINARY_FIELD_TYPE_PARTICLES
enumerator REBX_BINARY_FIELD_TYPE_FORCE
enumerator REBX_BINARY_FIELD_TYPE_SNAPSHOT
enum rebx_input_binary_messages

Possible errors that might occur during binary file reading.

Values:

enumerator REBX_INPUT_BINARY_WARNING_NONE
enumerator REBX_INPUT_BINARY_ERROR_NOFILE
enumerator REBX_INPUT_BINARY_ERROR_CORRUPT
enumerator REBX_INPUT_BINARY_ERROR_NO_MEMORY
enumerator REBX_INPUT_BINARY_ERROR_REBX_NOT_LOADED
enumerator REBX_INPUT_BINARY_ERROR_REGISTERED_PARAM_NOT_LOADED
enumerator REBX_INPUT_BINARY_WARNING_PARAM_NOT_LOADED
enumerator REBX_INPUT_BINARY_WARNING_PARTICLE_PARAMS_NOT_LOADED
enumerator REBX_INPUT_BINARY_WARNING_FORCE_NOT_LOADED
enumerator REBX_INPUT_BINARY_WARNING_OPERATOR_NOT_LOADED
enumerator REBX_INPUT_BINARY_WARNING_STEP_NOT_LOADED
enumerator REBX_INPUT_BINARY_WARNING_ADDITIONAL_FORCE_NOT_LOADED
enumerator REBX_INPUT_BINARY_WARNING_FIELD_UNKNOWN
enumerator REBX_INPUT_BINARY_WARNING_LIST_UNKNOWN
enumerator REBX_INPUT_BINARY_WARNING_PARAM_VALUE_NULL
enumerator REBX_INPUT_BINARY_WARNING_VERSION
enumerator REBX_INPUT_BINARY_WARNING_FORCE_PARAM_NOT_LOADED
enum rebx_integrator

Different schemes for integrating across the interaction step.

Values:

enumerator REBX_INTEGRATOR_NONE
enumerator REBX_INTEGRATOR_IMPLICIT_MIDPOINT
enumerator REBX_INTEGRATOR_RK4
enumerator REBX_INTEGRATOR_EULER
enumerator REBX_INTEGRATOR_RK2
enum rebx_interpolation_type

Different interpolation options.

Values:

enumerator REBX_INTERPOLATION_NONE
enumerator REBX_INTERPOLATION_SPLINE

Functions

void rebx_error(struct rebx_extras *rebx, const char *const msg)

Variables

const char *rebx_build_str

Date and time build string.

const char *rebx_version_str

Version string.

const char *rebx_githash_str

Current git hash.

struct rebx_node
#include <reboundx.h>

Node structure for all REBOUNDx linked lists.

Public Members

void *object

Pointer to object (param, force, step, etc)

struct rebx_node *next

Pointer to next node in list.

struct rebx_param
#include <reboundx.h>

Main structure used for all parameters added to objects.

Public Members

char *name

For searching linked lists and informative errors.

enum rebx_param_type type

Needed to cast value.

void *value

Pointer to parameter value.

struct rebx_force
#include <reboundx.h>

Structure for REBOUNDx forces.

Public Members

char *name

For searching linked lists and informative errors.

struct rebx_node *ap

Additional parameters linked list.

struct reb_simulation *sim

Pointer to attached sim. Needed for error checks.

enum rebx_force_type force_type

Force type for internal logic.

void (*update_accelerations)(struct reb_simulation *const sim, struct rebx_force *const force, struct reb_particle *const particles, const int N)

Function pointer to add additional accelerations.

struct rebx_operator
#include <reboundx.h>

Structure for REBOUNDx operators.

Public Members

char *name

For searching linked lists and informative errors.

struct rebx_node *ap

Additional parameters linked list.

struct reb_simulation *sim

Pointer to attached sim. Needed for error checks.

enum rebx_operator_type operator_type

Operator type for internal logic.

void (*step_function)(struct reb_simulation *sim, struct rebx_operator *operator, const double dt)

Function pointer to execute step.

struct rebx_step
#include <reboundx.h>

Structure for a REBOUNDx step.

A step is just a combination of an operator with a fraction of a timestep (see Sec. 6 of REBOUNdx paper). Can use same operator for different steps of different lengths to build higher order splitting schemes.

Public Members

struct rebx_operator * operator

Pointer to operator to use.

double dt_fraction

Fraction of sim.dt to use each time it’s called.

struct rebx_binary_field
#include <reboundx.h>

Structure used as building block to save and load binary files.

Public Members

enum rebx_binary_field_type type

Type of object.

long size

Size in bytes of the object data (not including this structure). So you can skip ahead.

struct rebx_interpolator

Public Members

enum rebx_interpolation_type interpolation
double *times
double *values
int Nvalues
double *y2
int klo
struct rebx_extras
#include <reboundx.h>

Main REBOUNDx structure.

These fields are used internally by REBOUNDx and generally should not be changed manually by the user. Use the API instead.

Public Members

struct reb_simulation *sim

Pointer to the simulation REBOUNDx is linked to.

struct rebx_node *additional_forces

Linked list of extra forces.

struct rebx_node *pre_timestep_modifications

Linked list of rebx_steps to apply before each timestep.

struct rebx_node *post_timestep_modifications

Linked list of rebx_steps to apply after each timestep.

struct rebx_node *registered_params

Linked list of rebx_params with all the parameter names registered with their type (for type safety)

struct rebx_node *allocated_forces

For memory management.

struct rebx_node *allocated_operators

For memory management.