.. include:: references.txt .. _astropy-coordinates-solarsystem: Solar System Ephemerides ************************ `astropy.coordinates` can calculate the |SkyCoord| of some of the major solar system objects. By default, it uses approximate orbital elements calculated using built-in `ERFA `_ routines, but it can also use more precise ones using the JPL ephemerides (which are derived from dynamical models). The default JPL ephemerides (DE430) provide predictions valid roughly for years between 1550 and 2650. The file is 115 MB and will need to be downloaded the first time you use this functionality, but will be cached after that. .. note:: Using JPL ephemerides requires that the `jplephem `_ package be installed. This is most easily achieved via ``pip install jplephem``, although whatever package management system you use might have it as well. Three functions are provided; :meth:`~astropy.coordinates.get_body`, :meth:`~astropy.coordinates.get_moon` and :meth:`~astropy.coordinates.get_body_barycentric`. The first two functions return |SkyCoord| objects in the `~astropy.coordinates.GCRS` frame, whilst the latter returns a `~astropy.coordinates.CartesianRepresentation` of the barycentric position of a body (i.e in the `~astropy.coordinates.ICRS` frame). Here is an example of using these functions with built-in ephemerides, i.e., without the need to download a large ephemerides file:: >>> from astropy.time import Time >>> from astropy.coordinates import solar_system_ephemeris, EarthLocation >>> from astropy.coordinates import get_body_barycentric, get_body, get_moon >>> t = Time("2014-09-22 23:22") >>> loc = EarthLocation.of_site('greenwich') # doctest: +REMOTE_DATA >>> with solar_system_ephemeris.set('builtin'): ... jup = get_body('jupiter', t, loc) # doctest: +REMOTE_DATA +IGNORE_OUTPUT >>> jup # doctest: +FLOAT_CMP +REMOTE_DATA Above, we used ``solar_system_ephemeris`` as a context, which sets the default ephemeris while in the ``with`` clause, and resets it at the end. To get more precise positions, one could use the ``de430`` ephemeris mentioned above, but between 1950 and 2050 one could also opt for the ``de432s`` ephemeris, which is stored in a smaller, ~10 MB, file (which will be downloaded and cached when the ephemeris is set): .. doctest-requires:: jplephem >>> solar_system_ephemeris.set('de432s') # doctest: +REMOTE_DATA, +IGNORE_OUTPUT >>> get_body('jupiter', t, loc) # doctest: +REMOTE_DATA, +FLOAT_CMP >>> get_moon(t, loc) # doctest: +REMOTE_DATA, +FLOAT_CMP >>> get_body_barycentric('moon', t) # doctest: +REMOTE_DATA, +FLOAT_CMP For one-off calculations with a given ephemeris, one can also pass it directly to the various functions: .. doctest-requires:: jplephem >>> get_body_barycentric('moon', t, ephemeris='de432s') ... # doctest: +REMOTE_DATA, +FLOAT_CMP >>> get_body_barycentric('moon', t, ephemeris='builtin') ... # doctest: +FLOAT_CMP For a list of the bodies for which positions can be calculated, do: .. note that we skip the next test if jplephem is not installed because if .. jplephem was not installed, we didn't change the science state higher up .. doctest-requires:: jplephem >>> solar_system_ephemeris.bodies # doctest: +REMOTE_DATA ('sun', 'mercury', 'venus', 'earth-moon-barycenter', 'earth', 'moon', 'mars', 'jupiter', 'saturn', 'uranus', 'neptune', 'pluto') >>> solar_system_ephemeris.set('builtin') >>> solar_system_ephemeris.bodies ('earth', 'sun', 'moon', 'mercury', 'venus', 'earth-moon-barycenter', 'mars', 'jupiter', 'saturn', 'uranus', 'neptune') .. note :: While the sun is included in the these ephemerides, it is important to recognize that `~astropy.coordinates.get_sun` always uses the built-in, polynomial model (as this requires no special download). So it is not safe to assume that ``get_body(time, 'sun')`` and ``get_sun(time)`` will give the same result.