callhorizons module¶
CALLHORIZONS - a Python interface to access JPL HORIZONS ephemerides and orbital elements.
This module provides a convenient python interface to the JPL HORIZONS system by directly accessing and parsing the HORIZONS website. Ephemerides can be obtained through get_ephemerides, orbital elements through get_elements. Function export2pyephem provides an interface to the PyEphem module.
michael.mommert (at) nau.edu, latest version: v1.0.1, 2016-07-19. This code is inspired by code created by Alex Hagen.
v1.0.3: ObsEclLon and ObsEclLat added to get_ephemerides v1.0.2: Python 3.5 compatibility implemented v1.0.1: get_ephemerides fixed v1.0: bugfixes completed, planets/satellites accessible, too v0.9: first release
-
class
callhorizons.
query
(targetname, smallbody=True, cap=True)[source]¶ -
__getitem__
(key)[source]¶ provides access to query data
Parameters: key : str/int
epoch index or property key
Returns: query data according to key
-
__init__
(targetname, smallbody=True, cap=True)[source]¶ Initialize query to Horizons
Parameters: targetname : str
HORIZONS-readable target number, name, or designation
smallbody : boolean
use
smallbody=False
if targetname is a planet or spacecraft (optional, default: True)cal : boolean
set to True to return the current apparition for comet targets.
-
__module__
= 'callhorizons'¶
-
dates
¶ returns list of epochs that have been queried (format ‘YYYY-MM-DD HH-MM-SS’)
-
dates_jd
¶ returns list of epochs that have been queried (Julian Dates)
-
export2pyephem
(center=u'500@10', equinox=2000.0)[source]¶ Call JPL HORIZONS website to obtain orbital elements based on the provided targetname, epochs, and center code and create a PyEphem (http://rhodesmill.org/pyephem/) object. This function requires PyEphem to be installed.
Parameters: center : str
center body (default: 500@10 = Sun)
equinox : float
equinox (default: 2000.0)
Examples
>>> import callhorizons >>> import numpy >>> import ephem
>>> ceres = callhorizons.query('Ceres') >>> ceres.set_epochrange('2016-02-23 00:00', '2016-02-24 00:00', '1h') >>> ceres_pyephem = ceres.export2pyephem()
>>> nau = ephem.Observer() # setup observer site >>> nau.lon = -111.653152/180.*numpy.pi >>> nau.lat = 35.184108/180.*numpy.pi >>> nau.elevation = 2100 # m >>> nau.date = '2015/10/5 01:23' # UT >>> print ('next rising: %s' % nau.next_rising(ceres_pyephem[0])) >>> print ('next transit: %s' % nau.next_transit(ceres_pyephem[0])) >>> print ('next setting: %s' % nau.next_setting(ceres_pyephem[0]))
-
fields
¶ returns list of available properties for all epochs
-
get_elements
(center=u'500@10')[source]¶ Call JPL HORIZONS website to obtain orbital elements based on the provided targetname, epochs, and center code. For valid center codes, please refer to http://ssd.jpl.nasa.gov/horizons.cgi
Parameters: center : str
center body (default: 500@10 = Sun)
Examples
>>> ceres = callhorizons.query('Ceres') >>> ceres.set_epochrange('2016-02-23 00:00', '2016-02-24 00:00', '1h') >>> print (ceres.get_elements(), 'epochs queried')
- The queried properties and their definitions are:
Property Definition targetname official number, name, designation [string] H absolute magnitude in V band (float, mag) G photometric slope parameter (float) datetime_jd epoch Julian Date (float) e eccentricity (float) p periapsis distance (float, au) a semi-major axis (float, au) incl inclination (float, deg) node longitude of Asc. Node (float, deg) argper argument of the perifocus (float, deg) Tp time of periapsis (float, Julian Date) meananomaly mean anomaly (float, deg) trueanomaly true anomaly (float, deg) period orbital period (float, Earth yr) Q apoapsis distance (float, au)
-
get_ephemerides
(observatory_code, airmass_lessthan=99, solar_elongation=(0, 180), skip_daylight=False)[source]¶ Call JPL HORIZONS website to obtain ephemerides based on the provided targetname, epochs, and observatory_code. For a list of valid observatory codes, refer to http://minorplanetcenter.net/iau/lists/ObsCodesF.html
Parameters: observatory_code : str/int
observer’s location code according to Minor Planet Center
airmass_lessthan : float
maximum airmass (optional, default: 99)
solar_elongation : tuple
permissible solar elongation range (optional, deg)
skip_daylight : boolean
crop daylight epoch during query (optional)
Examples
>>> ceres = callhorizons.query('Ceres') >>> ceres.set_epochrange('2016-02-23 00:00', '2016-02-24 00:00', '1h') >>> print (ceres.get_ephemerides(568), 'epochs queried')
- The queried properties and their definitions are:
Property Definition targetname official number, name, designation [string] H absolute magnitude in V band (float, mag) G photometric slope parameter (float) datetime epoch date and time (str, YYYY-MM-DD HH:MM:SS) datetime_jd epoch Julian Date (float) solar_presence information on Sun’s presence (str) lunar_presence information on Moon’s presence (str) RA target RA (float, J2000.0) DEC target DEC (float, J2000.0) RA_rate target rate RA (float, arcsec/s) DEC_rate target RA (float, arcsec/s, includes cos(DEC)) AZ Azimuth meas East(90) of North(0) (float, deg) EL Elevation (float, deg) airmass target optical airmass (float) magextinct V-mag extinction due airmass (float, mag) V V magnitude (comets: total mag) (float, mag) illumination fraction of illuminated disk (float) EclLon heliocentr. ecl. long. (float, deg, J2000.0) EclLat heliocentr. ecl. lat. (float, deg, J2000.0) ObsEclLon obscentr. ecl. long. (float, deg, J2000.0) ObsEclLat obscentr. ecl. lat. (float, deg, J2000.0) r heliocentric distance (float, au) r_rate heliocentric radial rate (float, km/s) delta distance from the observer (float, au) delta_rate obs-centric radial rate (float, km/s) lighttime one-way light time (float, s) elong solar elongation (float, deg) elongFlag app. position relative to Sun (str) alpha solar phase angle (float, deg) sunTargetPA PA of Sun->target vector (float, deg, EoN) velocityPA PA of velocity vector (float, deg, EoN) GlxLon galactic longitude (float, deg) GlxLat galactic latitude (float, deg) RA_3sigma 3sigma pos. unc. in RA (float, arcsec) DEC_3sigma 3sigma pos. unc. in DEC (float, arcsec)
-
isorbit_record
()[source]¶ True if targetname appears to be a comet orbit record number.
NAIF record numbers are 6 digits, begin with a ‘9’ and can change at any time.
-
parse_asteroid
()[source]¶ Parse targetname as if it were a asteroid.
Returns: des : string or None
The designation of the asteroid or None if targetname does not appear to be an asteroid name.
Examples
targetname des 1 1 (2) Pallas 2 (2001) Einstein 2001 2001 AT1 2001 AT1 (1714) Sy 1714 1714 SY 1714 SY # Note the near-confusion with (1714) 2014 MU69 2014 MU69 2017 AA 2017 AA
-
parse_comet
()[source]¶ Parse targetname as if it were a comet.
Returns: des : string or None
The designation of the comet or None if targetname does not appear to be a comet name. Note that comets starting with ‘X/’ are allowed, but this designation indicates a comet without an orbit, so query() should fail.
Examples
targetname des 1P/Halley 1P 3D/Biela 3D 9P/Tempel 1 9P 73P/Schwassmann Wachmann 3 C 73P # Note the missing “C”! 73P-C/Schwassmann Wachmann 3 C 73P-C 73P-BB 73P-BB 322P 322P X/1106 C1 X/1106 C1 P/1994 N2 (McNaught-Hartley) P/1994 N2 P/2001 YX127 (LINEAR) P/2001 YX127 C/-146 P1 C/-146 P1 C/2001 A2-A (LINEAR) C/2001 A2-A C/2013 US10 C/2013 US10 C/2015 V2 (Johnson) C/2015 V2
-
query
¶ returns URL that has been used in calling HORIZONS
-
set_discreteepochs
(discreteepochs)[source]¶ Set a list of discrete epochs, epochs have to be given as Julian Dates
Parameters: discreteepochs : list
list of floats or strings, maximum length: 15
Returns: None
Examples
>>> import callhorizons >>> ceres = callhorizons.query('Ceres') >>> ceres.set_discreteepochs([2457446.177083, 2457446.182343])
If more than 15 epochs are provided, the list will be cropped to 15 epochs.
-
set_epochrange
(start_epoch, stop_epoch, step_size)[source]¶ Set a range of epochs, all times are UT
Parameters: start_epoch : str
start epoch of the format ‘YYYY-MM-DD [HH-MM-SS]’
stop_epoch : str
final epoch of the format ‘YYYY-MM-DD [HH-MM-SS]’
step_size : str
epoch step size, e.g., ‘1d’ for 1 day, ‘10m’ for 10 minutes...
Returns: None
Examples
>>> import callhorizons >>> ceres = callhorizons.query('Ceres') >>> ceres.set_epochrange('2016-02-26', '2016-10-25', '1d')
Note that dates are mandatory; if no time is given, midnight is assumed.
-