API documentation (Python)

An N-body integrator package for python.

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

SimulationArchive Class.

The SimulationArchive is a binary file format which includes all settings, constants as well as particle positions and velocities. This makes it possible to reproduce a simulation exactly (down to the last bit). The SimulationArchive allows you to add an arbitrary number of snapshots. Simulations can be reconstructed from these snapshots. Since version 2 of the SimulationArchive (Spring 2018), you can change anything inbetween snapshots, icluding settings like the integrator, the timestep, the number of particles. The file format is efficient in that only data that changed is stored in the SimulationArchive file. This is all done automatically. All the user has to do is call the function to create a snapshot. The SimulationArchive thus allows for fast access to any long running simulations. 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)
Attributes
auto_interval

Structure/Union member

auto_step

Structure/Union member

auto_walltime

Structure/Union member

nblobs

Structure/Union member

offset

Structure/Union member

size_first

Structure/Union member

size_snapshot

Structure/Union member

t

Structure/Union member

version

Structure/Union member

Methods

getBezierPaths(self[, origin])

This function returns array that can be used as a Cubic Bezier Path in matplotlib.

getSimulation(self, t[, mode, …])

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

getSimulations(self, times, \*\*kwargs)

A generator to quickly access many simulations.

__init__(self, filename, setup=None, setup_args=(), process_warnings=True)[source]
getBezierPaths(self, origin=None)[source]

This function returns array that can be used as a Cubic Bezier Path in matplotlib. The function returns two arrays, the first one contains the verticies for each particles and has the shape (Nvert, Nparticles, 2) where Nvert is the number of verticies. The second array returned describes the type of verticies to be used with matplotlib’s Patch class.

Examples

The following example reads in a SimulationArchive and plots the trajectories as Cubic Bezier Curves. It also plots the actual datapoints stored in the SimulationArchive. Note that the SimulationArchive needs to have enough datapoints to allow for smooth and reasonable orbits.

>>> from matplotlib.path import Path
>>> import matplotlib.patches as patches
>>> sa = rebound.SimulationArchive("test.bin")
>>> verts, codes = sa.getBezierPaths(origin=0)
>>> fig, ax = plt.subplots()
>>> for j in range(sa[0].N):
>>>     path = Path(verts[:,j,:], codes)
>>>     patch = patches.PathPatch(path, facecolor='none')
>>>     ax.add_patch(patch)
>>>     ax.scatter(verts[::3,j,0],verts[::3,j,1])
>>> ax.set_aspect('equal')
>>> ax.autoscale_view()
getSimulation(self, t, mode='snapshot', keep_unsynchronized=1)[source]

This function returns a simulation object at (or close to) the requested time t. Everytime this function is called a new simulation object is created.

Returns
A rebound.Simulation object. Everytime the function gets called
a new object gets created.

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(self, times, **kwargs)[source]

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

class rebound.Simulation(filename=None, snapshot=None)[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

By calling rebound.Simulation() as shown above, you create a new simulation object 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("simulation.bin")

Similarly, you can create a simulation, from a simulation archive by specifying the snapshot you want to load.

>>> sim = rebound.Simulation("archive.bin", 34)

or

>>> sim = rebound.Simulation(filename="archive.bin", snapshot=34)
Attributes
G

Structure/Union member

N

Structure/Union member

N_active

Structure/Union member

N_lookup

Structure/Union member

N_real

Get the current number of real (i.e.

N_var

Structure/Union member

additional_forces

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

allocatedN_lookup

Structure/Union member

allocated_N

Structure/Union member

boundary

Get or set the boundary module.

boxsize

Structure/Union member

boxsize_max

Structure/Union member

coefficient_of_restitution

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

collision

Get or set the collision module.

collision_resolve

Get or set a function pointer for collision resolving routine.

collision_resolve_keep_sorted

Structure/Union member

collisions

Structure/Union member

collisions_Nlog

Structure/Union member

collisions_allocatedN

Structure/Union member

collisions_plog

Structure/Union member

display_data

Structure/Union member

dt

Structure/Union member

dt_last_done

Structure/Union member

energy_offset

Structure/Union member

exact_finish_time

Structure/Union member

exit_max_distance

Structure/Union member

exit_min_distance

Structure/Union member

extras

Structure/Union member

force_is_velocity_dependent

Structure/Union member

free_particle_ap

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

gravity

Get or set the gravity module.

gravity_cs

Structure/Union member

gravity_cs_allocatedN

Structure/Union member

gravity_ignore

Structure/Union member

hash_ctr

Structure/Union member

heartbeat

Set a function pointer for a heartbeat function.

integrator

Get or set the intergrator module.

max_radius

Structure/Union member

messages

Structure/Union member

minimum_collision_velocity

Structure/Union member

nghostx

Structure/Union member

nghosty

Structure/Union member

nghostz

Structure/Union member

opening_angle2

Structure/Union member

particles

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

post_timestep_modifications

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

pre_timestep_modifications

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

python_unit_l

Structure/Union member

python_unit_m

Structure/Union member

python_unit_t

Structure/Union member

ri_eos

Structure/Union member

ri_ias15

Structure/Union member

ri_janus

Structure/Union member

ri_mercurius

Structure/Union member

ri_saba

Structure/Union member

ri_sei

Structure/Union member

ri_whfast

Structure/Union member

root_n

Structure/Union member

root_nx

Structure/Union member

root_ny

Structure/Union member

root_nz

Structure/Union member

root_size

Structure/Union member

save_messages

Structure/Union member

simulationarchive_auto_interval

Structure/Union member

simulationarchive_auto_step

Structure/Union member

simulationarchive_auto_walltime

Structure/Union member

simulationarchive_filename

Returns the current SimulationArchive filename in use.

simulationarchive_next

Structure/Union member

simulationarchive_next_step

Structure/Union member

simulationarchive_size_first

Structure/Union member

simulationarchive_size_snapshot

Structure/Union member

simulationarchive_version

Structure/Union member

softening

Structure/Union member

steps_done

Structure/Union member

t

Structure/Union member

testparticle_type

Structure/Union member

track_energy_offset

Structure/Union member

units

Tuple of the units for length, time and mass.

usleep

Structure/Union member

var_config

Structure/Union member

var_config_N

Structure/Union member

walltime

Structure/Union member

Methods

add(self[, particle])

Adds a particle to REBOUND.

add_particles_ascii(self, s)

Adds particles from an ASCII string.

add_variation(self[, order, first_order, …])

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

automateSimulationArchive(self, filename[, …])

This function automates taking snapshots during a simulationusing the Simulation Archive.

calculate_angular_momentum(self)

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

calculate_com(self[, first, last])

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

calculate_energy(self)

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

calculate_lyapunov(self)

Return the current Lyapunov Characteristic Number (LCN).

calculate_megno(self)

Return the current MEGNO value.

calculate_orbits(self[, primary, …])

Calculate orbital parameters for all partices in the simulation.

configure_box(self, boxsize[, root_nx, …])

Initialize the simulation box.

configure_ghostboxes(self[, nghostx, …])

Initialize the ghost boxes.

convert_particle_units(self, \*args)

Will convert the units for the simulation (i.e.

copy(self)

Returns a deep copy of a REBOUND simulation.

from_archive(filename[, snapshot])

rebound.Simulation.from_archive(filename,snapshot) is deprecated and will be removed in the futute.

from_file(filename)

rebound.Simulation.from_file(filename) is deprecated and will be removed in the futute.

getWidget(self, \*\*kwargs)

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

init_megno(self[, seed])

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

integrate(self, tmax[, exact_finish_time])

Main integration function.

integrator_reset(self)

Call this function to reset temporary integrator variables

integrator_synchronize(self)

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

move_to_com(self)

This function moves all particles in the simulation to a center of momentum frame.

move_to_hel(self)

This function moves all particles in the simulation to the heliocentric frame.

particles_ascii(self[, prec])

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

refreshWidgets(self)

This function manually refreshed all widgets attached to this simulation.

remove(self[, index, hash, keepSorted])

Removes a particle from the simulation.

save(self, filename)

Save the entire REBOUND simulation to a binary file.

serialize_particle_data(self, \*\*kwargs)

Fast way to access serialized particle data via numpy arrays.

set_serialized_particle_data(self, \*\*kwargs)

Fast way to set serialized particle data via numpy arrays.

simulationarchive_snapshot(self, filename[, …])

Take a snapshot and save it to a SimulationArchive file.

status(self)

Prints a summary of the current status of the simulation.

step(self)

Perform exactly one integration step with REBOUND.

tree_update(self)

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

multiply

process_messages

update_units

property N_real

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

__init__(self, filename=None, snapshot=None)[source]

Initialize self. See help(type(self)) for accurate signature.

add(self, 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(self, s)[source]

Adds particles from an ASCII string.

Parameters
sstring

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

add_variation(self, 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
orderinteger, 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_orderVariation, 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_2Variation, 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.

testparticleint, 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).
property 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.

automateSimulationArchive(self, filename, interval=None, walltime=None, step=None, deletefile=False)[source]

This function automates taking snapshots during a simulationusing the Simulation Archive. Instead of using this function, one can also call simulationarchive_snapshot() manually to create snapshots.

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.automateSimulationArchive("sa.bin",interval=1000.)
>>> sim.integrate(1e8)

The SimulationArchive can later be read in using the following syntax:

>>> sa = rebound.SimulationArchive("sa.bin")
>>> sim = sa[0]   # get the first snapshot in the SA file (initial conditions)
>>> sim = sa[-1]  # get the last snapshot in the SA file
property 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(self)[source]

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

calculate_com(self, 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``.

lastint 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(self)[source]

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

calculate_lyapunov(self)[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(self)[source]

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

calculate_orbits(self, primary=None, jacobi_masses=False, heliocentric=None, barycentric=None)[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
primaryrebound.Particle, optional

Set the primary against which to reference the osculating orbit. Default(use Jacobi center of mass)

jacobi_masses: bool

Whether to use jacobi primary mass in orbit calculation. (Default: False)

heliocentric: bool, DEPRECATED

To calculate heliocentric elements, pass primary=sim.particles[0]

barycentricbool, DEPRECATED

To calculate barycentric elements, pass primary=sim.calculate_com()

Returns
Returns an array of Orbits of length N-1.
property coefficient_of_restitution

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

property collision

Get or set the collision module.

Available collision modules are:

  • 'none' (default)

  • 'direct'

  • 'tree'

  • 'mercurius'

  • 'direct'

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

property 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. “hardsphere”: two colliding particles will bounce of using a set coefficient of restitution

configure_box(self, 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
boxsizefloat, optional

The size of one root box.

root_nx, root_ny, root_nzint, 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(self, 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, nghostzint

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

convert_particle_units(self, *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
copy(self)[source]

Returns a deep copy of a REBOUND simulation. You need to reset any function pointers on the copy.

Returns
A rebound.Simulation object.
property 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]

rebound.Simulation.from_archive(filename,snapshot) is deprecated and will be removed in the futute. Use rebound.Simulation(filename,snapshot) instead

classmethod from_file(filename)[source]

rebound.Simulation.from_file(filename) is deprecated and will be removed in the futute. Use rebound.Simulation(filename) instead

getWidget(self, **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()
property 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.

property 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.)
init_megno(self, seed=None)[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 MEGNO see e.g. http://dx.doi.org/10.1051/0004-6361:20011189

integrate(self, 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
tmaxfloat

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)
property integrator

Get or set the intergrator module.

Available integrators include:

  • 'IAS15' (default)

  • 'WHFast'

  • 'SEI'

  • 'LEAPFROG'

  • 'JANUS'

  • 'MERCURIUS'

  • 'WHCKL'

  • 'WHCKM'

  • 'WHCKC'

  • 'SABA4'

  • 'SABACL4'

  • 'SABACM4'

  • 'SABA(10,6,4)'

  • 'EOS'

  • 'none'

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

integrator_reset(self)[source]

Call this function to reset temporary integrator variables

integrator_synchronize(self)[source]

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

move_to_com(self)[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.

move_to_hel(self)[source]

This function moves all particles in the simulation to the heliocentric frame. Note that the simulation will not stay in the heliocentric frame during integrations as the heliocentric frame is not an innertial frame.

property 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(self, prec=8)[source]

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

Parameters
precint, optional

Number of digits after decimal point. Default 8.

property 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.

property 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(self)[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(self, index=None, hash=None, keepSorted=True)[source]

Removes a particle from the simulation.

Parameters
indexint, optional

Specify particle to remove by index.

hashc_uint32 or string, optional

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

keepSortedbool, 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(self, filename)[source]

Save the entire REBOUND simulation to a binary file.

serialize_particle_data(self, **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”, and “xyzvxvyvz”. 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. The length of “xyzvxvyvz” arrays need to be 6*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)
set_serialized_particle_data(self, **kwargs)[source]

Fast way to set serialized particle data via numpy arrays. This is the inverse of Simulation.serialize_particle_data() and uses the same syntax

property simulationarchive_filename

Returns the current SimulationArchive filename in use. Do not set manually. Use sim.automateSimulationArchive() instead

simulationarchive_snapshot(self, filename, deletefile=False)[source]

Take a snapshot and save it to a SimulationArchive file. If the file does not exist yet, a new one will be created. If the file does exist, a snapshot will be appended.

status(self)[source]

Prints a summary of the current status of the simulation.

step(self)[source]

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

tree_update(self)[source]

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

property 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
dfloat

radial distance from reference

vfloat

velocity relative to central object’s velocity

hfloat

specific angular momentum

Pfloat

orbital period (negative if hyperbolic)

nfloat

mean motion (negative if hyperbolic)

afloat

semimajor axis

efloat

eccentricity

incfloat

inclination

Omegafloat

longitude of ascending node

omegafloat

argument of pericenter

pomegafloat

longitude of pericenter

ffloat

true anomaly

Mfloat

mean anomaly

lfloat

mean longitude = Omega + omega + M

thetafloat

true longitude = Omega + omega + f

Tfloat

time of pericenter passage

rhillfloat

Hill radius ( =a*pow(m/(3M),1./3.) )

rebound.OrbitPlot(sim, figsize=None, fancy=False, slices=0, xlim=None, ylim=None, unitlabel=None, color=False, periastron=False, orbit_type='trail', lw=1.0, plotparticles=[], primary=None, Narc=128)[source]

Convenience function for plotting instantaneous orbits.

Parameters
figsizetuple of float, optional

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

fancybool (default: False)

Changes various settings to create a fancy looking plot

slicesfloat, optional

Default is 0, showing the orbits in the xy plane only. Set to a value between 0 and 1 to create three plots, showing the orbits from different directions. The value corresponds to the size of the additional plots relative to the main plot.

xlimtuple of float, optional

Limits for x axes (default: None = automatically determined)

ylimtuple of float, optional

Limits for y axes (default: None = automatically determined)

unitlabelstr, optional

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

colorbool, str or list, optional

By default plots are black and white. If set to True, plots use a color cycle. If a string or list of strings, e.g. [‘red’, ‘cyan’], will cycle between passed colors.

periastronbool, optional

Draw a marker at periastron (default: False)

orbit_typestr, optional

This argument determines the type of orbit show. By default it shows the orbit as a trailing and fading line (“trail”). Other object are: “solid”, None.

lwfloat, optional

Linewidth used in plots (default: 1.)

plotparticleslist, optional

List of particles to plot. Can be a list of any valid keys for accessing sim.particles, i.e., integer indices or hashes (default: plot all particles)

primaryrebound.Particle, optional

Primary to use for the osculating orbit (default: Jacobi center of mass)

Narcint, optional

Number of points used in an orbit. Increase this number for highly eccentric orbits. (default: 128)

Returns
fig, ax_main, (ax_top, ax_right)

The function return the matplotlib figure as well as the axes (three axes if slices>0.)

Examples

The following example illustrates a typical use case.

>>> sim = rebound.Simulation()
>>> sim.add(m=1)
>>> sim.add(a=1)
>>> fig, ax_main = 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, jacobi_masses=False)[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, zfloat

Particle positions

vx, vy, vzfloat

Particle velocities

ax, ay, azfloat

Particle accelerations

mfloat

Particle mass

rfloat

Particle radius

lastcollisionfloat

Last time the particle had a physical collision (if checking for collisions)

cc_void_p (C void pointer)

Pointer to the cell the particle is currently in (if using tree code)

hashc_uint32

Get or set the particle’s hash.

apc_void_p (C void pointer)

Pointer to additional parameters one might want to add to particles

_simPOINTER(rebound.Simulation)

Internal pointer to the parent simulation (used in C version of REBOUND)

a, e, inc, Omega, omega, ffloat

(Kepler Elements) Semi-major axis, eccentricity, inclination, longitude of the ascending node, argument of periapsis, and true anomaly respectively. The Keplerian Elements are in Jacobi coordinates (with mu = G*Minc, where Minc is the total mass from index 0 to the particle’s index, inclusive).

Methods

calculate_orbit(self[, primary, G])

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 (with mu = G*Minc, where Minc is the total mass from index 0 to the particle’s index, inclusive).

copy(self)

Returns a deep copy of the particle.

sample_orbit(self[, Npts, primary, …])

Returns a nested list of xyz positions along the osculating orbit of the particle.

__init__(self, 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, jacobi_masses=False)[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
simulationSimulation

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

particleParticle, 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.

mfloat

Mass (Default: 0)

x, y, zfloat

Positions in Cartesian coordinates (Default: 0)

vx, vy, vzfloat

Velocities in Cartesian coordinates (Default: 0)

primaryParticle

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)

afloat

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

Pfloat

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

efloat

Eccentricity (Default: 0)

incfloat

Inclination (Default: 0)

Omegafloat

Longitude of ascending node (Default: 0)

omegafloat

Argument of pericenter (Default: 0)

pomegafloat

Longitude of pericenter (Default: 0)

ffloat

True anomaly (Default: 0)

Mfloat

Mean anomaly (Default: 0)

lfloat

Mean longitude (Default: 0)

thetafloat

True longitude (Default: 0)

Tfloat

Time of pericenter passage

hfloat

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

kfloat

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

ixfloat

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

iyfloat

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

rfloat

Particle radius (only used for collisional simulations)

datestring

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

variationstring (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.

variation2string (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.

hashc_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.

jacobi_masses: bool

Whether to use jacobi primary mass in orbit initialization. Particle mass will still be set to physical value (Default: False)

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)
>>> p4 = rebound.Particle(simulation=sim, m=0.001, a=1.5, omega=”uniform”) # omega will be a random number between 0 and 2pi
calculate_orbit(self, 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 (with mu = G*Minc, where Minc is the total mass from index 0 to the particle’s index, inclusive).

Parameters
primaryrebound.Particle

Central body (Optional. Default uses Jacobi coordinates)

Gfloat

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(self)[source]

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

property hash

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

sample_orbit(self, Npts=100, primary=None, trailing=True, timespan=None, useTrueAnomaly=None)[source]

Returns a nested list of xyz positions along the osculating orbit of the particle. If primary is not passed, returns xyz positions along the Jacobi osculating orbit (with mu = G*Minc, where Minc is the total mass from index 0 to the particle’s index, inclusive).

Parameters
Nptsint, optional

Number of points along the orbit to return (default: 100)

primaryrebound.Particle, optional

Primary to use for the osculating orbit (default: Jacobi center of mass)

trailing: bool, optional

Whether to return points stepping backwards in time (True) or forwards (False). (default: True)

timespan: float, optional

Return points (for the osculating orbit) from the current position to timespan (forwards or backwards in time depending on trailing keyword). Defaults to the orbital period for bound orbits, and to the rough time it takes the orbit to move by the current distance from the primary for a hyperbolic orbit. Implementation currently only supports this option if useTrueAnomaly=False.

useTrueAnomaly: bool, optional

Will sample equally spaced points in true anomaly if True, otherwise in mean anomaly. Latter might be better for hyperbolic orbits, where true anomaly can stay near the limiting value for a long time, and then switch abruptly at pericenter. (Default: True for bound orbits, False for unbound orbits)

property vxyz

Get or set the xyz velocity coordinates of the particle.

property 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.Collision[source]

The simulation exited because a collision has been detected. You may want to search for which particles have a lastcollision time equal to the simulation time.

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

apply(self, func[, args, kwds])

Equivalent of func(*args, **kwds).

apply_async(self, func[, args, kwds, …])

Asynchronous version of apply() method.

imap(self, func, iterable[, chunksize])

Equivalent of map() – can be MUCH slower than Pool.map().

imap_unordered(self, func, iterable[, chunksize])

Like imap() method but ordering of results is arbitrary.

map(self, func, iterable[, chunksize])

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

map_async(self, func, iterable[, chunksize, …])

Asynchronous version of map() method.

starmap(self, func, iterable[, chunksize])

Like map() method but the elements of the iterable are expected to be iterables as well and will be unpacked as arguments.

starmap_async(self, func, iterable[, …])

Asynchronous version of starmap() method.

Process

close

join

terminate

__init__(self, processes=None, initializer=None, initargs=(), **kwargs)[source]

Initialize self. See help(type(self)) for accurate signature.

map(self, 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
index

Structure/Union member

index_1st_order_a

Structure/Union member

index_1st_order_b

Structure/Union member

order

Structure/Union member

particles

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

testparticle

Structure/Union member

Methods

vary(self, particle_index, variation[, …])

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

property 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(self, 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_indexint

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.

variationstring

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) and in Rein, Tamayo, and Brown (2019).

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

>>> sim = rebound.Simulation()
>>> sim.integrator = "whfast"
>>> sim.ri_whfast.corrector =  11
>>> sim.ri_whfast.kernel = "lazy"
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 or 17. For more details read Rein and Tamayo (2015).

  • corrector2 (int) – Second correctors (C2 of Wisdom et al 1996). By default the second symplectic correctors are turned off (=0). Set to 1 to turn them on.

  • kernel (int/string) – Kernel option. Set to 0 for the default WH kernel (standard kick step). Other options are “modifiedkick” (1), “composition” (2), “lazy” (3).

  • recalculate_coordinates_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/heliocenctric 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.

  • coordinates (string) – Sets the internal coordinate system that WHFast is using. By default it uses 'jacobi' (=0) coordinates. Other options are 'democraticheliocentric' (=1) and 'whds' (=2). See Hernandez and Dehnen (2017) for more information.

  • 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.

Attributes
coordinates

Get or set the internal coordinate system.

corrector

Structure/Union member

corrector2

Structure/Union member

is_synchronized

Structure/Union member

keep_unsynchronized

Structure/Union member

kernel

Get or set the WHFast Kernel.

recalculate_coordinates_this_timestep

Structure/Union member

safe_mode

Structure/Union member

property coordinates

Get or set the internal coordinate system.

Available coordinate systems are:

  • 'jacobi' (default)

  • 'democraticheliocentric'

  • 'whds'

property kernel

Get or set the WHFast Kernel.

Available kernels are:

  • 'default' (standard WH kernel, kick)

  • 'modifiedkick' (modified kick for newtonian gravity)

  • 'composition' (Wisdom’s composition method)

  • 'lazy' (Lazy implementer’s method)

class rebound.reb_simulation_integrator_saba[source]

This class is an abstraction of the C-struct reb_simulation_integrator_saba. It controls the behaviour of the SABA integrator family. See Rein, Tamayo, and Brown (2019) for more details.

Variables
  • type (str) – Set the type of SABA integrator manually. The type can also be set by setting the integrator field in the REBOUND simulation.

  • safe_mode (int) – This variable acts the same as for WHFast. 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 the integrator will improve. However, make sure you are aware of the consequences. Read the iPython tutorial on advanced WHFast usage to learn more.

Example usage:

>>> sim = rebound.Simulation()
>>> sim.integrator = "SABA(10,6,4)"
>>> sim = rebound.Simulation()
>>> sim.integrator = "SABA"
>>> sim.ri_saba.type = "(10,6,4)"
>>> sim = rebound.Simulation()
>>> sim.integrator = "SABACL4"
>>> sim.ri_saba.safe_mode = 0
Attributes
is_synchronized

Structure/Union member

keep_unsynchronized

Structure/Union member

safe_mode

Structure/Union member

type

Get or set the type of SABA integrator.

property type

Get or set the type of SABA integrator.

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.

Attributes
OMEGA

Structure/Union member

OMEGAZ

Structure/Union member

class rebound.reb_simulation_integrator_mercurius[source]

This class is an abstraction of the C-struct reb_simulation_integrator_mercurius. It controls the behaviour of the MERCURIUS integrator. See Rein et al. (2019) for more details.

Variables

hillfac (float) – Switching radius in units of the hill radius.

Example usage:

>>> sim = rebound.Simulation()
>>> sim.integrator = "mercurius"
>>> sim.ri_mercurius.hillfac = 3.
Attributes
L

Structure/Union member

hillfac

Structure/Union member

mode

Structure/Union member

recalculate_coordinates_this_timestep

Structure/Union member

recalculate_dcrit_this_timestep

Structure/Union member

safe_mode

Structure/Union member