REBOUND can add particles to simulations by obtaining ephemerides from NASA's powerful HORIZONS database. HORIZONS supports many different options, and we will certainly not try to cover everything here. This is meant to serve as an introduction to the basics, beyond what's in Churyumov-Gerasimenko.ipynb. If you catch any errors, or would either like to expand on this documentation or improve REBOUND's HORIZONS interface (rebound/horizons.py
), please do fork the repository and send us a pull request.
Adding particles
When we add particles by passing a string, REBOUND queries the HORIZONS database and takes the first dataset HORIZONS offers. For the Sun, moons, and small bodies, this will typically return the body itself. For planets, it will return the barycenter of the system (for moonless planets like Venus it will say barycenter but there is no distinction). If you want the planet specifically, you have to use, e.g., "NAME=Pluto" rather than "Pluto". In all cases, REBOUND will print out the name of the HORIZONS entry it's using.
You can also add bodies using their integer NAIF IDs: NAIF IDs. Note that because of the number of small bodies (asteroids etc.) we have discovered, this convention only works for large objects. For small bodies, instead use "NAME=name" (see the SMALL BODIES section in the HORIZONS Documentation).
import rebound
sim = rebound.Simulation()
sim.add("Sun")
## Other examples:
# sim.add("Venus")
# sim.add("399")
# sim.add("Europa")
# sim.add("NAME=Ida")
# sim.add("Pluto")
# sim.add("NAME=Pluto")
sim.status()
Searching NASA Horizons for 'Sun'... Found: Sun (10) --------------------------------- REBOUND version: 3.17.4 REBOUND built on: Sep 28 2021 11:18:40 Number of particles: 1 Selected integrator: ias15 Simulation time: 0.0000000000000000e+00 Current timestep: 0.001000 --------------------------------- <rebound.particle.Particle object at 0x7f8222fc2240, m=0.9999999999950272 x=-0.008209560798576862 y=0.004129411861725172 z=0.0001581586597591073 vx=-0.0002624462976527945 vy=-0.0004663115919721187 vz=9.942935240423656e-06> ---------------------------------
Currently, HORIZONS does not have any mass information for solar system bodies. rebound/horizons.py
has a hard-coded list provided by Jon Giorgini (10 May 2015) that includes the planets, their barycenters (total mass of planet plus moons), and the largest moons. If REBOUND doesn't find the corresponding mass for an object from this list (like for the asteroid Ida below), it will print a warning message. If you need the body's mass for your simulation, you can set it manually, e.g. (see Units.ipynb for an overview of using different units):
sim.add("NAME=Ida")
print(sim.particles[-1]) # Ida before setting the mass
sim.particles[-1].m = 2.1e-14 # Setting mass of Ida in Solar masses
print(sim.particles[-1]) # Ida after setting the mass
Searching NASA Horizons for 'NAME=Ida'... Found: 243 Ida (A884 SB) <rebound.particle.Particle object at 0x7f8222fc2240, m=0.0 x=-2.2120992198130702 y=1.7367783911044365 z=0.001902768696505959 vx=-0.3911361032281711 vy=-0.45991860843972715 vz=-0.011860531651599532> <rebound.particle.Particle object at 0x7f8222fc2240, m=2.1e-14 x=-2.2120992198130702 y=1.7367783911044365 z=0.001902768696505959 vx=-0.3911361032281711 vy=-0.45991860843972715 vz=-0.011860531651599532>
/home/user/rebound/rebound/horizons.py:146: RuntimeWarning: Warning: Mass cannot be retrieved from NASA HORIZONS. Set to 0. warnings.warn("Warning: Mass cannot be retrieved from NASA HORIZONS. Set to 0.", RuntimeWarning)
Time
By default, REBOUND queries HORIZONS for objects' current positions. Specifically, it caches the current time the first time you call rebound.add
, and gets the corresponding ephemeris. All subsequent calls to rebound.add
will then use that initial cached time to make sure you get a synchronized set of ephemerides.
You can also explicitly pass REBOUND the time at which you would like the particles ephemerides:
sim = rebound.Simulation()
date = "2005-06-30 15:24" # You can also use Julian Days. For example: date = "JD2458327.500000"
sim.add("Venus")
sim.add("Venus", date=date)
sim.status()
Searching NASA Horizons for 'Venus'... Found: Venus Barycenter (299) (chosen from query 'Venus') Searching NASA Horizons for 'Venus'... Found: Venus Barycenter (299) (chosen from query 'Venus') --------------------------------- REBOUND version: 3.17.4 REBOUND built on: Sep 28 2021 11:18:40 Number of particles: 2 Selected integrator: ias15 Simulation time: 0.0000000000000000e+00 Current timestep: 0.001000 --------------------------------- <rebound.particle.Particle object at 0x7f8222fc2240, m=2.447838287784771e-06 x=0.4025719758008419 y=-0.5962868903087453 z=-0.03178627392421859 vx=0.962311037745445 vy=0.6593297118253567 vz=-0.04648020412614457> <rebound.particle.Particle object at 0x7f8222fc2c40, m=2.447838287784771e-06 x=-0.6612400840226685 y=0.26970377066078294 z=0.04197305542931168 vx=-0.445669425641083 vy=-1.0953912971910784 vz=0.010725839790529726> ---------------------------------
We see that the two Venus positions are different. The first call cached the current time, but since the second call specified a date, it overrode the default. Any time you pass a date, it will overwrite the default cached time, so:
sim = rebound.Simulation()
date = "2005-06-30 15:24" # You can also use Julian Days. For example: date = "JD2458327.500000"
sim.add("Venus", date=date)
sim.add("Earth")
sim.status()
Searching NASA Horizons for 'Venus'... Found: Venus Barycenter (299) (chosen from query 'Venus') Searching NASA Horizons for 'Earth'... Found: Earth-Moon Barycenter (3) (chosen from query 'Earth') --------------------------------- REBOUND version: 3.17.4 REBOUND built on: Sep 28 2021 11:18:40 Number of particles: 2 Selected integrator: ias15 Simulation time: 0.0000000000000000e+00 Current timestep: 0.001000 --------------------------------- <rebound.particle.Particle object at 0x7f8222fc2cc0, m=2.447838287784771e-06 x=-0.6612400840226685 y=0.26970377066078294 z=0.04197305542931168 vx=-0.445669425641083 vy=-1.0953912971910784 vz=0.010725839790529726> <rebound.particle.Particle object at 0x7f8222fc2d40, m=3.0404326489511185e-06 x=0.9897025414054398 y=0.09460299451532682 z=0.00015007541573282246 vx=-0.10681915648565657 vy=0.9918343709972148 vz=-3.7394344962530276e-05> ---------------------------------
would set up a simulation with Venus and Earth, all synchronized to 2005-06-30 15:24. All dates should either be passed in the format Year-Month-Day Hour:Minute or in JDxxxxxxx.xxxx for a date in Julian Days.
REBOUND takes these absolute times to the nearest minute, since at the level of seconds you have to worry about exactly what time system you're using, and small additional perturbations probably start to matter. For reference HORIZONS interprets all times for ephemerides as Coordinate (or Barycentric Dynamical) Time.
Reference Frame
REBOUND queries for particles' positions and velocities relative to the Sun:
sim = rebound.Simulation()
sim.add("Sun")
sim.status()
Searching NASA Horizons for 'Sun'... Found: Sun (10) --------------------------------- REBOUND version: 3.17.4 REBOUND built on: Sep 28 2021 11:18:40 Number of particles: 1 Selected integrator: ias15 Simulation time: 0.0000000000000000e+00 Current timestep: 0.001000 --------------------------------- <rebound.particle.Particle object at 0x7f8222fc2240, m=0.9999999999950272 x=-0.008209560798576862 y=0.004129411861725172 z=0.0001581586597591073 vx=-0.0002624462976527945 vy=-0.0004663115919721187 vz=9.942935240423656e-06> ---------------------------------
The reference plane is the ecliptic (Earth's orbital plane) of J2000 (Jan. 1st 2000 12:00 GMT), with the x-axis along the ascending node of the ecliptic and the Earth's mean equator (also at J2000).