API documentation (Python)

An N-body integrator package for python.

class rebound.SimulationArchive(filename, setup=None, setup_args=(), rebxfilename=None)[source]

SimulationArchive Class.

The SimulationArchive is a binary file that includes all settings, constants as well as initial conditions to reproduce a simulation exactly (down to the last bit). It is also possible to add snapshots of the current state as the simulation runs. This allows for fast access to long running simulations by making use of the SimulationArchive binary file. For a full discussion of the functionality see the paper by Rein & Tamayo 2017.

Examples

Here is a simple example:

>>> sa = rebound.SimulationArchive("archive.bin")
>>> sim = sa.getSimulation(t=1e6)
>>> print(sim.particles[1])
>>> for sim in sa:
>>>     print(sim.t, sim.particles[1].e)

Methods

estimateTime(t, tbefore=None)[source]

This function estimates the time needed to integrate the simulation exactly to time the t starting from the nearest snapshot or the current status of the simulation (whichever is smaller).

If an array is passed as an argument, the function will estimate the time it will take to integrate to all the times in the array. The array is assumed to be sorted. The function assumes a simulation can be reused to get to the next requested time if this is faster than reloading the simulation from the nearest snapshot.

Note that the estimates are based on the runtime of the original simulation. If the original simulation was run on a different machine, the estimated runtime may differ.

Returns:An approximation of the runtime (in seconds) required to integrate to time t.
getSimulation(t, mode='snapshot', keep_unsynchronized=1)[source]

This function returns a simulation object at (or close to) the requested time t.

Returns:

A rebound.Simulation object. This object will be modified

the next time getSimulation is called. Making any manual

changes to this object could have unforseen consequences.

Examples

Here is a simple example on how to load a simulation from a Simulation Archive file with the getSimulation method. As the mode argument is set to close, the simulation will be integrated from the nearest snapshot to the request time.

>>> sa = rebound.SimulationArchive("archive.bin")
>>> sim = sa.getSimulation(t=1e6, mode="close")
>>> print(sim.particles[1])
>>> for sim in sa:
>>>     print(sim.t, sim.particles[1].e)
getSimulations(times, **kwargs)[source]

A generator to quickly access many simulations. The arguments are the same as for getSimulation.

class rebound.Simulation[source]

REBOUND Simulation Object.

This object encapsulated an entire REBOUND simulation. It is an abstraction of the C struct reb_simulation. You can create mutiple REBOUND simulations (the c library is thread safe).

Examples

Most simulation parameters can be directly changed with the property syntax:

>>> sim = rebound.Simulation()
>>> sim.G = 1.                  # Sets the graviational constant (default 1)
>>> sim.softening = 1.          # Sets the graviational softening parameter (default 0)
>>> sim.testparticle_type = 1   # Allows massive particles to feel influence from testparticles (default 0)
>>> sim.dt = 0.1                # Sets the timestep (will change for adaptive integrators such as IAS15).
>>> sim.t = 0.                  # Sets the current simulation time (default 0)
>>> print(sim.N)                # Gets the current number of particles
>>> print(sim.N_active)         # Gets the current number of active particles

Attributes

Methods

N_real

Get the current number of real (i.e. no variational/shadow) particles in the simulation.

add(particle=None, **kwargs)[source]

Adds a particle to REBOUND. Accepts one of the following:

  1. A single Particle structure.
  2. The particle’s mass and a set of cartesian coordinates: m,x,y,z,vx,vy,vz.
  3. The primary as a Particle structure, the particle’s mass and a set of orbital elements: primary,m,a,anom,e,omega,inv,Omega,MEAN (see Orbit for the definition of orbital elements).
  4. A name of an object (uses NASA Horizons to look up coordinates)
  5. A list of particles or names.
add_particles_ascii(s)[source]

Adds particles from an ASCII string.

Parameters:

s : string

One particle per line. Each line should include particle’s mass, radius, position and velocity.

add_variation(order=1, first_order=None, first_order_2=None, testparticle=-1)[source]

This function adds a set of variational particles to the simulation.

If there are N real particles in the simulation, this functions adds N additional variational particles. To see how many particles (real and variational) are in a simulation, use 'sim.N'. To see how many variational particles are in a simulation use 'sim.N_var'.

Currently Leapfrog, WHFast and IAS15 support first order variational equations. IAS15 also supports second order variational equations.

Parameters:

order : integer, optional

By default the function adds a set of first order variational particles to the simulation. Set this flag to 2 for second order.

first_order : Variation, optional

Second order variational equations depend on their corresponding first order variational equations. This parameter expects the Variation object corresponding to the first order variational equations.

first_order_2 : Variation, optional

Same as first_order. But allows to set two different indicies to calculate off-diagonal elements. If omitted, then first_order will be used for both first order equations.

testparticle : int, optional

If set to a value >= 0, then only one variational particle will be added and be treated as a test particle.

Returns:

Returns Variation object (a copy–you can only modify it through its particles property or vary method).

additional_forces

Get or set a function pointer for calculating additional forces in the simulation.

The argument can be a python function or something that can be cast to a C function of type CFUNCTYPE(None,POINTER(Simulaton)). If the forces are velocity dependent, the flag force_is_velocity_dependent needs to be set to 1. Otherwise the particle structures might contain incorrect velocity values.

boundary

Get or set the boundary module.

Available boundary modules are:

  • 'none' (default)
  • 'open'
  • 'periodic'
  • 'shear'

Check the online documentation for a full description of each of the modules.

calculate_angular_momentum()[source]

Returns a list of the three (x,y,z) components of the total angular momentum of all particles in the simulation.

calculate_com(first=0, last=None)[source]

Returns the center of momentum for all particles in the simulation.

Parameters:

first: int, optional

If first is specified, only calculate the center of momentum starting from index=``first``.

last : int or None, optional

If last is specified only calculate the center of momentum up to (but excluding) index=``last``. Same behavior as Python’s range function.

Examples

>>> sim = rebound.Simulation()
>>> sim.add(m=1, x=-20)
>>> sim.add(m=1, x=-10)
>>> sim.add(m=1, x=0)
>>> sim.add(m=1, x=10)
>>> sim.add(m=1, x=20)
>>> com = sim.calculate_com()
>>> com.x
0.0 
>>> com = sim.calculate_com(first=2,last=4) # Considers indices 2,3
>>> com.x
5.0
calculate_energy()[source]

Returns the sum of potential and kinetic energy of all particles in the simulation.

calculate_lyapunov()[source]

Return the current Lyapunov Characteristic Number (LCN). Note that you need to call init_megno() before the start of the simulation. To get a timescale (the Lyapunov timescale), take the inverse of this quantity.

calculate_megno()[source]

Return the current MEGNO value. Note that you need to call init_megno() before the start of the simulation.

calculate_orbits(heliocentric=False, barycentric=False)[source]

Calculate orbital parameters for all partices in the simulation. By default this functions returns the orbits in Jacobi coordinates.

If MEGNO is enabled, variational particles will be ignored.

Parameters:

heliocentric : bool, optional

Set the parameter heliocentric to True to return orbits referenced to sim.particles[0].

barycentric : bool, optional

Set the parameter barycentric to True to return orbits referenced to the system’s barycenter.

Returns:

Returns an array of Orbits of length N-1.

coefficient_of_restitution

Get or set a function pointer that defined the coefficient of restitution.

collision

Get or set the collision module.

Available collision modules are:

  • 'none' (default)
  • 'direct'
  • 'tree'

Check the online documentation for a full description of each of the modules.

collision_resolve

Get or set a function pointer for collision resolving routine.

Possible options for setting:
  1. Function pointer
  2. “merge”: two colliding particles will merge)
  3. “harsphere”: two colliding particles will bounce of using a set coefficient of restitution
configure_box(boxsize, root_nx=1, root_ny=1, root_nz=1)[source]

Initialize the simulation box.

This function only needs to be called it boundary conditions other than “none” are used. In such a case the boxsize must be known and is set with this function.

Parameters:

boxsize : float, optional

The size of one root box.

root_nx, root_ny, root_nz : int, optional

The number of root boxes in each direction. The total size of the simulation box will be root_nx * boxsize, root_ny * boxsize and root_nz * boxsize. By default there will be exactly one root box in each direction.

configure_ghostboxes(nghostx=0, nghosty=0, nghostz=0)[source]

Initialize the ghost boxes.

This function only needs to be called it boundary conditions other than “none” or “open” are used. In such a case the number of ghostboxes must be known and is set with this function.

Parameters:

nghostx, nghosty, nghostz : int

The number of ghost boxes in each direction. All values default to 0 (no ghost boxes).

convert_particle_units(*args)[source]

Will convert the units for the simulation (i.e. convert G) as well as the particles’ cartesian elements. Must have set sim.units ahead of calling this function so REBOUND knows what units to convert from.

Parameters:3 strings corresponding to units of time, length and mass. Can be in any order and aren’t case sensitive. You can add new units to rebound/rebound/units.py
estimateSimulationArchiveSize(tmax)[source]

This function estimates the SimulationArchive file size (in bytes) before a simulation is run. This is useful to check if the interval results in a resonable filesize.

Note that the simulation setup needs to be complete, that is with all the particles present, before this function can return meaningful results.

free_particle_ap

Get or set a function pointer for freeing the ap pointer whenever sim.remove is called.

classmethod from_archive(filename, snapshot=-1)[source]

Loads a REBOUND simulation from a SimulationArchive binary file. It uses the last snapshot in that file unless otherwise specified. This function is only really useful for restarting a simulation. For full access to the Simulation Archive, use the SimulationArchive class.

After loading the REBOUND simulation from file, you need to reset any function pointers manually. Please read all documentation before using this function as it has several requirements to work correctly.

Returns:A rebound.Simulation object.
classmethod from_file(filename)[source]

Loads a REBOUND simulation from a file.

After loading the REBOUND simulation from file, you need to reset any function pointers manually.

Returns:A rebound.Simulation object.

Examples

The following example creates a simulation, saves it to a file and then creates a copy of the simulation store in the binary file.

>>> sim = rebound.Simulation()
>>> sim.add(m=1.)
>>> sim.add(m=1.e-3,x=1.,vy=1.)
>>> sim.save("simulation.bin")
>>> sim_copy = rebound.Simulation.from_file("simulation.bin")
getWidget(**kwargs)[source]

Wrapper function that returns a new widget attached to this simulation.

Widgets provide real-time 3D visualizations from within an Jupyter notebook. See the Widget class for more details on the possible arguments.

Returns:A rebound.Widget object.

Examples

>>> sim = rebound.Simulation()
>>> sim.add(m=1.)
>>> sim.add(m=1.e-3,x=1.,vy=1.)
>>> sim.getWidget()
gravity

Get or set the gravity module.

Available gravity modules are:

  • 'none'
  • 'basic' (default)
  • 'compensated'
  • 'tree'

Check the online documentation for a full description of each of the modules.

heartbeat

Set a function pointer for a heartbeat function. The heartbeat function is called every timestep and can be used to monitor long simulations, check for stalled simulations and output debugging information.

The argument can be a python function or something that can be cast to a C function or a python function.

The function called will receive a pointer to the simulation object as its argument.

Examples

>>> import rebound
>>> sim = rebound.Simulation()
>>> sim.add(m=1.)
>>> sim.add(m=1e-3, a=1)
>>> def heartbeat(sim):
>>>     # sim is a pointer to the simulation object,
>>>     # thus use contents to access object data.
>>>     # See ctypes documentation for details.
>>>     print(sim.contents.t)
>>> sim.heartbeat = heartbeat
>>> sim.integrate(1.)
initSimulationArchive(filename, interval=None, interval_walltime=None)[source]

This function initializes the Simulation Archive so that binary data can be outputted to the SimulationArchive file during the simulation.

Examples

The following example creates a simulation, then initializes the Simulation Archive and integrates it forward in time.

>>> sim = rebound.Simulation()
>>> sim.add(m=1.)
>>> sim.add(m=1.e-3,x=1.,vy=1.)
>>> sim.initSimulationArchive("sa.bin",interval=1000.)
>>> sim.integrate(1e8)
init_megno()[source]

This function initialises the chaos indicator MEGNO particles and enables their integration.

MEGNO is short for Mean Exponential Growth of Nearby orbits. It can be used to test if a system is chaotic or not. In the backend, the integrator is integrating an additional set of particles using the variational equation. Note that variational equations are better suited for this than shadow particles. MEGNO is currently only supported in the IAS15 and WHFast integrators.

This function also needs to be called if you are interested in the Lyapunov exponent as it is calculate with the help of MEGNO. See Rein and Tamayo 2015 for details on the implementation.

For more information on MENGO see e.g. http://dx.doi.org/10.1051/0004-6361:20011189

integrate(tmax, exact_finish_time=1)[source]

Main integration function. Call this function when you have setup your simulation and want to integrate it forward (or backward) in time. The function might be called many times to integrate the simulation in steps and create outputs in-between steps.

Parameters:

tmax : float

The final time of your simulation. If the current time is 100, and tmax=200, then after the calling the integrate routine, the time has advanced to t=200. If tmax is larger than or equal to the current time, no integration will be performed.

exact_finish_time: int, optional

This argument determines whether REBOUND should try to finish at the exact time (tmax) you give it or if it is allowed to overshoot. Overshooting could happen if one starts at t=0, has a timestep of dt=10 and wants to integrate to tmax=25. With exact_finish_time=1, the integrator will choose the last timestep such that t is exactly 25 after the integration, otherwise t=30. Note that changing the timestep does affect the accuracy of symplectic integrators negatively.

Examples

The typical usage is as follows. Note the use of np.linspace to create equally spaced outputs. Using np.logspace can be used to easily produce logarithmically spaced outputs.

>>> import numpy as np
>>> for time in np.linspace(0,100.,10): 
>>>     sim.integrate(time)
>>>     perform_output(sim)
integrator

Get or set the intergrator module.

Available integrators are:

  • 'ias15' (default)
  • 'whfast'
  • 'whfasthelio'
  • 'sei'
  • 'leapfrog'
  • 'hermes'
  • 'none'

Check the online documentation for a full description of each of the integrators.

integrator_synchronize()[source]

Call this function if safe-mode is disabled and you need synchronize particle positions and velocities between timesteps.

move_to_com()[source]

This function moves all particles in the simulation to a center of momentum frame. In that frame, the center of mass is at the origin and does not move. It makes sense to call this function at the beginning of the integration, especially for the high accuray integrators IAS15 and WHFast.

particles

Returns a Particles object that allows users to access particles like a dictionary using indices, hashes, or strings.

The Particles object uses pointers and thus the contents update as the simulation progresses. Note that the pointers could change, for example when a particle is added or removed from the simulation.

particles_ascii(prec=8)[source]

Returns an ASCII string with all particles’ masses, radii, positions and velocities.

Parameters:

prec : int, optional

Number of digits after decimal point. Default 8.

post_timestep_modifications

Get or set a function pointer for post-timestep modifications.

The argument can be a python function or something that can be cast to a C function or a python function.

pre_timestep_modifications

Get or set a function pointer for pre-timestep modifications.

The argument can be a python function or something that can be cast to a C function or a python function.

refreshWidgets()[source]

This function manually refreshed all widgets attached to this simulation.

You want to call this function if any particle data has been manually changed.

remove(index=None, hash=None, keepSorted=True)[source]

Removes a particle from the simulation.

Parameters:

index : int, optional

Specify particle to remove by index.

hash : c_uint32 or string, optional

Specifiy particle to remove by hash (if a string is passed, the corresponding hash is calculated).

keepSorted : bool, optional

By default, remove preserves the order of particles in the particles array. Might set it to zero in cases with many particles and many removals to speed things up.

save(filename)[source]

Save the entire REBOUND simulation to a binary file.

serialize_particle_data(**kwargs)[source]

Fast way to access serialized particle data via numpy arrays.

This function can directly set the values of numpy arrays to current particle data. This is significantly faster than accessing particle data via sim.particles as all the copying is done on the C side. No memory is allocated by this function. It expects correctly sized numpy arrays as arguments. The argument name indicates what kind of particle data is written to the array.

Possible argument names are “hash”, “m”, “r”, “xyz”, “vxvyvz”. The datatype for the “hash” array needs to be uint32. The other arrays expect a datatype of float64. The lengths of “hash”, “m”, “r” arrays need to be at least sim.N. The lengths of xyz and vxvyvz need to be at least 3*sim.N. Exceptions are raised otherwise.

Note that this routine is only intended for special use cases where speed is an issue. For normal use, it is recommended to access particle data via the sim.particles array. Be aware of potential issues that arrise by directly accesing the memory of numpy arrays (see numpy documentation for more details).

Examples

This sets an array to the xyz positions of all particles:

>>> import numpy as np
>>> a = np.zeros((sim.N,3),dtype="float64")
>>> sim.serialize_particle_data(xyz=a)
>>> print(a)

To get all current radii of particles:

>>> a = np.zeros(sim.N,dtype="float64")
>>> sim.serialize_particle_data(r=a)
>>> print(a)

To get all current radii and hashes of particles:

>>> a = np.zeros(sim.N,dtype="float64")
>>> b = np.zeros(sim.N,dtype="uint32")
>>> sim.serialize_particle_data(r=a,hash=b)
>>> print(a,b)
status()[source]

Prints a summary of the current status of the simulation.

step()[source]

Perform exactly one integration step with REBOUND. This function is rarely needed. Instead, use integrate().

tree_update()[source]

Call this function to update the tree structure manually after removing particles.

units

Tuple of the units for length, time and mass. Can be set in any order, and strings are not case-sensitive. See ipython_examples/Units.ipynb for more information. You can check the units’ exact values and add Additional units in rebound/rebound/units.py. Units should be set before adding particles to the simulation (will give error otherwise).

Examples

>>> sim = rebound.Simulation()
>>> sim.units = ('yr', 'AU', 'Msun')
class rebound.Orbit[source]

A class containing orbital parameters for a particle. This is an abstraction of the reb_orbit data structure in C.

When using the various REBOUND functions using Orbits, all angles are in radians. The following image illustrated the most important angles used. In REBOUND the reference direction is the positive x direction, the reference plane is the xy plane.

_images/orbit.png

Image from wikipedia. CC-BY-SA-3.

Attributes

d (float) radial distance from reference
v (float) velocity relative to central object’s velocity
h (float) specific angular momentum
P (float) orbital period (negative if hyperbolic)
n (float) mean motion (negative if hyperbolic)
a (float) semimajor axis
e (float) eccentricity
inc (float) inclination
Omega (float) longitude of ascending node
omega (float) argument of pericenter
pomega (float) longitude of pericenter
f (float) true anomaly
M (float) mean anomaly
l (float) mean longitude = Omega + omega + M
theta (float) true longitude = Omega + omega + f
T (float) time of pericenter passage
rebound.OrbitPlot(sim, figsize=None, lim=None, limz=None, Narc=100, unitlabel=None, color=False, periastron=False, trails=True, lw=1.0, slices=False)[source]

Convenience function for plotting instantaneous orbits.

Parameters:

slices : bool, optional

Plot all three slices if set to True. Default is False and plots orbits only in the xy plane.

figsize : tuple of float, optional

Tuple defining the figure size (default: (5,5))

lim : float, optional

Limit for axes (default: None = automatically determined)

limz : float, optional

Limit for z axis, only used if slices=True (default: None = automatically determined)

unitlabel : str, optional

String describing the units, shown on axis labels (default: None)

color : bool, optional

Enable color (default: False)

periastron : bool, optional

Draw a marker at periastron (default: False)

trails : bool, optional

Draw trails instead of solid lines (default: False)

lw : float, optional

Linewidth (default: 1.)

Returns:

fig

A matplotlib figure

Examples

The following example illustrates a typical use case.

>>> sim = rebound.Simulation()
>>> sim.add(m=1)
>>> sim.add(a=1)
>>> fig = rebound.OrbitPlot(sim)
>>> fig.savefig("image.png") # save figure to file
>>> fig.show() # show figure on screen
class rebound.Particle(simulation=None, particle=None, m=None, x=None, y=None, z=None, vx=None, vy=None, vz=None, primary=None, a=None, P=None, e=None, inc=None, Omega=None, omega=None, pomega=None, f=None, M=None, l=None, theta=None, T=None, r=None, date=None, variation=None, variation2=None, h=None, k=None, ix=None, iy=None, hash=0)[source]

The main REBOUND particle data structure. This is an abstraction of the reb_particle structure in C. The Particle fields are set at the end of simulation.py to avoid circular references.

Attributes

x, y, z (float) Particle positions
vx, vy, vz (float) Particle velocities
ax, ay, az (float) Particle accelerations
m (float) Particle mass
r (float) Particle radius
lastcollision (float) Last time the particle had a physical collision (if checking for collisions)
c (c_void_p (C void pointer)) Pointer to the cell the particle is currently in (if using tree code)
ap (c_void_p (C void pointer)) Pointer to additional parameters one might want to add to particles
_sim (POINTER(rebound.Simulation)) Internal pointer to the parent simulation (used in C version of REBOUND)

Methods

__init__(simulation=None, particle=None, m=None, x=None, y=None, z=None, vx=None, vy=None, vz=None, primary=None, a=None, P=None, e=None, inc=None, Omega=None, omega=None, pomega=None, f=None, M=None, l=None, theta=None, T=None, r=None, date=None, variation=None, variation2=None, h=None, k=None, ix=None, iy=None, hash=0)[source]

Initializes a Particle structure. Rather than explicitly creating a Particle structure, users may use the add() member function of a Simulation instance, which will both create a Particle and then add it to the simulation with one function call.

This function accepts either cartesian positions and velocities, classical orbital elements together with the reference Particle (the primary), as well as orbital parameters defined by Pal (2009).

For convenience, optional keywords that are not passed default to zero (mass, cartesian and orbital elements).

Whenever initializing a particle from orbital elements, one must specify either the semimajor axis or the period of the orbit.

For classical orbital paramerers, one can specify the longitude of the ascending node by passing Omega, to specify the pericenter one can pass either omega or pomega (not both), and for the longitude/anomaly one can pass one of f, M, l or theta. See ipython_examples/OrbitalElements.ipynb for examples. See also Murray & Dermott Solar System Dynamics for formal definitions of angles in orbital mechanics.

All angles should be specified in radians.

Parameters:

simulation : Simulation

Simulation instance associated with this particle (Required if passing orbital elements or setting up a variation).

particle : Particle, optional

If a particle is passed, a copy of that particle is returned. If a variational particle is initialized, then particle is original particle that will be varied.

m : float

Mass (Default: 0)

x, y, z : float

Positions in Cartesian coordinates (Default: 0)

vx, vy, vz : float

Velocities in Cartesian coordinates (Default: 0)

primary : Particle

Primary body for converting orbital elements to cartesian (Default: center of mass of the particles in the passed simulation, i.e., this will yield Jacobi coordinates as one progressively adds particles)

a : float

Semimajor axis (a or P required if passing orbital elements)

P : float

Orbital period (a or P required if passing orbital elements)

e : float

Eccentricity (Default: 0)

inc : float

Inclination (Default: 0)

Omega : float

Longitude of ascending node (Default: 0)

omega : float

Argument of pericenter (Default: 0)

pomega : float

Longitude of pericenter (Default: 0)

f : float

True anomaly (Default: 0)

M : float

Mean anomaly (Default: 0)

l : float

Mean longitude (Default: 0)

theta : float

True longitude (Default: 0)

T : float

Time of pericenter passage

h : float

h variable, see Pal (2009) for a definition (Default: 0)

k : float

k variable, see Pal (2009) for a definition (Default: 0)

ix : float

ix variable, see Pal (2009) for a definition (Default: 0)

iy : float

iy variable, see Pal (2009) for a definition (Default: 0)

r : float

Particle radius (only used for collisional simulations)

date : string

For consistency with adding particles through horizons. Not used here.

variation : string (Default: None)

Set this string to the name of an orbital parameter to initialize the particle as a variational particle. Can be one of the following: m, a, e, inc, omega, Omega, f, k, h, lambda, ix, iy.

variation2 : string (Default: None)

Set this string to the name of a second orbital parameter to initialize the particle as a second order variational particle. Only used for second order variational equations. Can be one of the following: m, a, e, inc, omega, Omega, f, k, h, lambda, ix, iy.

hash : c_uint32

Unsigned integer identifier for particle. Can pass an integer directly, or a string that will be converted to a hash. User is responsible for assigning unique hashes.

Examples

>>> sim = rebound.Simulation()
>>> sim.add(m=1.)
>>> p1 = rebound.Particle(simulation=sim, m=0.001, a=0.5, e=0.01)
>>> p2 = rebound.Particle(simulation=sim, m=0.0, x=1., vy=1.)
>>> p3 = rebound.Particle(simulation=sim, m=0.001, a=1.5, h=0.1, k=0.2, l=0.1)
calculate_orbit(primary=None, G=None)[source]

Returns a rebound.Orbit object with the keplerian orbital elements corresponding to the particle around the passed primary (rebound.Particle) If no primary is passed, defaults to Jacobi coordinates.

Parameters:

primary : rebound.Particle

Central body (Optional. Default uses Jacobi coordinates)

G : float

Gravitational constant (Optional. Default takes G from simulation in which particle is in)

Returns:

A rebound.Orbit object

Examples

>>> sim = rebound.Simulation()
>>> sim.add(m=1.)
>>> sim.add(x=1.,vy=1.)
>>> orbit = sim.particles[1].calculate_orbit(sim.particles[0])
>>> print(orbit.e) # gives the eccentricity
copy()[source]

Returns a deep copy of the particle. The particle is not added to any simulation by default.

hash

Get or set the particle’s hash. If set to a string, the corresponding integer hash is calculated.

vxyz

Get or set the xyz velocity coordinates of the particle.

xyz

Get or set the xyz position coordinates of the particle.

exception rebound.SimulationError[source]

The simulation exited with a generic error.

exception rebound.Encounter[source]

The simulation exited because a close encounter has been detected. You may want to search for the pair of bodies which have the smallest distance.

exception rebound.Escape[source]

The simulation exited because a particle has been se encounter has been detected. You may want to search for the particle with the largest distance from the origin and remove it from the simulation.

exception rebound.NoParticles[source]

The simulation exited because no particles are left in the simulation.

exception rebound.ParticleNotFound[source]

Particle was not found in the simulation.

class rebound.InterruptiblePool(processes=None, initializer=None, initargs=(), **kwargs)[source]

A modified version of multiprocessing.pool.Pool that has better behavior with regard to KeyboardInterrupts in the map() method.

Parameters:
  • processes – (optional) The number of worker processes to use; defaults to the number of CPUs.
  • initializer – (optional) Either None, or a callable that will be invoked by each worker process when it starts.
  • initargs – (optional) Arguments for initializer; it will be called as initializer(*initargs).
  • kwargs – (optional) Extra arguments. Python 2.7 supports a maxtasksperchild parameter.

Methods

map(func, iterable, chunksize=None)[source]

Equivalent of map() built-in, without swallowing KeyboardInterrupt.

Parameters:
  • func – The function to apply to the items.
  • iterable – An iterable of items that will have func applied to them.
class rebound.Variation[source]

REBOUND Variational Configuration Object.

This object encapsulated the configuration of one set of variational equations in a REBOUND simulation. It is an abstraction of the C struct reb_variational_configuration.

None of the fields in this struct should be changed after it has been initialized.

One rebound simulation can include any number of first and second order variational equations.

Note that variations are only encoded as particles for convenience. A variational particle’s position and velocity should be interpreted as a derivative, i.e. how much that position orvelocity varies with respect to the first or second-order variation. See ipython_examples/VariationalEquations.ipynb and Rein and Tamayo (2016) for details.

Examples

>>> sim = rebound.Simulation()          # Create a simulation
>>> sim.add(m=1.)                       # Add a star
>>> sim.add(m=1.e-3, a=1.)              #     a planet
>>> var_config = sim.add_variation()    # Add a set of variational equations. 
>>> var_config.particles[1].x = 1.      # Initialize the variational particle corresponding to the planet
>>> sim.integrate(100.)                 # Integrate the simulation
>>> print(var_config.particles[0].vx)   # Print the velocity of the variational particle corresponding to the star

Attributes

Methods

particles

Access the variational particles corresponding to this set of variational equations.

The function returns a list of particles which are sorted in the same way as those in sim.particles

The particles are pointers and thus can be modified.

If there are N real particles, this function will also return a list of N particles (all of which are variational particles).

vary(particle_index, variation, variation2=None, primary=None)[source]

This function can be used to initialize the variational particles that are part of a Variation.

Note that rather than using this convenience function, one can also directly manipulate the particles’ coordinate using the following syntax:

>>> var = sim.add_variation()
>>> var.particles[0].x = 1.

The vary() function is useful for initializing variations corresponding to changes in one of the orbital parameters for a particle on a bound Keplerian orbit.

The function supports both first and second order variations in the following classical orbital parameters:

a, e, inc, omega, Omega, f
as well as the Pal (2009) coordinates:
a, h, k, ix, iy, lambda

and in both cases the mass m of the particle. The advantage of the Pal coordinate system is that all derivatives are well behaved (infinitely differentiable). Classical orbital parameters on the other hand exhibit coordinate singularities, for example when e=0.

The following example initializes the variational particles corresponding to a change in the semi-major axis of the particle with index 1:

>>> var = sim.add_variation()
>>> var.vary(1,"a")
Parameters:

particle_index : int

The index of the particle that should be varied. The index starts at 0 and runs through N-1. The first particle added to a simulation receives the index 0, the second 1, and the on.

variation : string

This parameter determines which orbital parameter is varied.

variation2: string, optional

This is only used for second order variations which can depend on two varying parameters. If omitted, then it is assumed that the parameter variation is variation2.

primary: Particle, optional

By default variational particles are created in the Heliocentric frame. Set this parameter to use any other particles as a primary (e.g. the center of mass).

class rebound.reb_simulation_integrator_whfast[source]

This class is an abstraction of the C-struct reb_simulation_integrator_whfast. It controls the behaviour of the symplectic WHFast integrator described in Rein and Tamayo (2015).

This struct should be accessed via the simulation class only. Here is an example:

>>> sim = rebound.Simulation()
>>> sim.ri_whfast.corrector =  11
Variables:
  • corrector (int) – The order of the symplectic corrector in the WHFast integrator. By default the symplectic correctors are turned off (=0). For high accuracy simulation set this value to 11. For more details read Rein and Tamayo (2015).
  • recalculate_jacobi_this_timestep (int) – Sets a flag that tells WHFast that the particles have changed. Setting this flag to 1 (default 0) triggers the WHFast integrator to recalculate Jacobi coordinates. This is needed if the user changes the particle position, velocity or mass inbetween timesteps. After every timestep the flag is set back to 0, so if you continuously update the particles manually, you need to set this flag to 1 after every timestep.
  • safe_mode (int) – If safe_mode is 1 (default) particles can be modified between timesteps and particle velocities and positions are always synchronised. If you set safe_mode to 0, the speed and accuracy of WHFast improve. However, make sure you are aware of the consequences. Read the iPython tutorial on advanced WHFast usage to learn more.
class rebound.reb_simulation_integrator_sei[source]

This class is an abstraction of the C-struct reb_simulation_integrator_sei. It controls the behaviour of the symplectic SEI integrator for shearing sheet calculations. It is described in Rein and Tremaine (2011).

This struct should be accessed via the simulation class only. Here is an example:

>>> sim = rebound.Simulation()
>>> sim.ri_sei.OMEGA =  1.58
Variables:
  • OMEGA (float) – The epicyclic frequency OMEGA. For simulations making use of shearing sheet boundary conditions, REBOUND needs to know the epicyclic frequency. By default OMEGA is 1. For more details read Rein and Tremaine 2011.
  • OMEGAZ (float) – The z component of the epicyclic frequency OMEGA. By default it is assuming OMEGAZ is the same as OMEGA.