Skyfield: HomeTable of ContentsAPI Reference

API Reference

Quick links to the sections below:

Opening files

# File you already have.

from skyfield.api import load_file
planets = load_file('~/Downloads/de405.bsp')
load_file(path) Open a file on your local drive, using its extension to guess its type.
# File you want Skyfield to download automatically.

from skyfield.api import load
ts = load.timescale()
planets = load('de405.bsp')
Loader(directory[, verbose, expire]) A tool for downloading and opening astronomical data files.
Loader.build_url(filename) Return the URL Skyfield will try downloading for a given filename.
Loader.days_old(filename) Return how recently filename was modified, measured in days.
Loader.download(url[, filename, backup]) Download a file, even if it’s already on disk; return its path.
Loader.path_to(filename) Return the path to filename in this loader’s directory.
Loader.timescale([delta_t, builtin]) Return a Timescale built using official Earth rotation data.
Loader.tle_file(url[, reload, filename, ts, …]) Load and parse a TLE file, returning a list of Earth satellites.

Time scales

A script will typically start by building a single Skyfield Timescale to use for all date and time conversions:

from skyfield import api
ts = api.load.timescale()

Its methods are:

Timescale.now() Return the current date and time as a Time object.
Timescale.from_datetime(datetime) Return a Time for a Python datetime.
Timescale.from_datetimes(datetime_list) Return a Time for a list of Python datetimes.
Timescale.utc(year[, month, day, hour, …]) Build a Time from a UTC Calendar date.
Timescale.tai([year, month, day, hour, …]) Build a Time from an International Atomic Time Calendar date.
Timescale.tai_jd(jd[, fraction]) Build a Time from an International Atomic Time Julian date.
Timescale.tt([year, month, day, hour, …]) Build a Time from a Terrestrial Time Calendar date.
Timescale.tt_jd(jd[, fraction]) Build a Time from a Terrestrial Time Julian date.
Timescale.J(year) Build a Time from a Terrestrial Time Julian year or array.
Timescale.tdb([year, month, day, hour, …]) Build a Time from a Barycentric Dynamical Time Calendar date.
Timescale.tdb_jd(jd[, fraction]) Build Time from a Barycentric Dynamical Time Julian date.
Timescale.ut1([year, month, day, hour, …]) Build a Time from a UT1 Universal Time Calendar date.
Timescale.ut1_jd(jd) Build a Time from a UT1 Universal Time Julian date.
Timescale.from_astropy(t) Build a Skyfield Time from an AstroPy time object.

Time objects

The Time class is Skyfield’s way of representing either a single time, or a whole array of times. The same time can be represented in several different time scales.

t.tai International Atomic Time (TAI) as a Julian date.
t.tt Terrestrial Time (TT) as a Julian date.
t.J Terrestrial Time (TT) as floating point Julian years.
t.tdb Barycentric Dynamical Time (TDB) as a Julian date.
t.ut1 Universal Time (UT1) as a Julian date.

A couple of offsets between time scales are also available.

t.delta_t Difference TT − UT1 in seconds.
t.dut1 Difference UT1 − UTC in seconds.

Other time scales and conversions are available through its methods.

Time.utc_jpl() Convert to a string like A.D.2014-Jan-18 01:35:37.5000 UT..
Time.utc_iso([delimiter, places]) Convert to an ISO 8601 string like 2014-01-18T01:35:38Z in UTC.
Time.utc_strftime([format]) Format the UTC time using a Python datetime formatting string.
Time.utc_datetime() Convert to a Python datetime in UTC.
Time.utc_datetime_and_leap_second() Convert to a Python datetime in UTC, plus a leap second value.
Time.astimezone(tz) Convert to a Python datetime in a particular timezone tz.
Time.astimezone_and_leap_second(tz) Convert to a Python datetime and leap second in a timezone.
Time.toordinal() Return the proleptic Gregorian ordinal of the UTC date.
Time.tai_calendar() TAI as a (year, month, day, hour, minute, second) Calendar date.
Time.tt_calendar() TT as a (year, month, day, hour, minute, second) Calendar date.
Time.tdb_calendar() TDB as a (year, month, day, hour, minute, second) Calendar date.
Time.ut1_calendar() UT1 as a (year, month, day, hour, minute, second) Calendar date.
Time.tai_strftime([format]) Format TAI with a datetime strftime() format string.
Time.tt_strftime([format]) Format TT with a datetime strftime() format string.
Time.tdb_strftime([format]) Format TDB with a datetime strftime() format string.
Time.ut1_strftime([format]) Format UT1 with a datetime strftime() format string.
Time.M 3×3 rotation matrix: ICRS → equinox of this date.
Time.MT 3×3 rotation matrix: equinox of this date → ICRS.
Time.gmst Greenwich Mean Sidereal Time as decimal hours.
Time.gast Greenwich Apparent Sidereal Time as decimal hours.
Time.nutation_matrix() Compute the 3×3 nutation matrix N for this date.
Time.precession_matrix() Compute the 3×3 precession matrix P for this date.

Time utilities

compute_calendar_date(jd_integer[, …]) Convert Julian day jd_integer into a calendar (year, month, day).

Vector functions

The common API shared by planets, Earth locations, and Earth satellites.

VectorFunction Given a time, computes a corresponding position.
VectorFunction.at(t) At time t, compute the target’s position relative to the center.

Either adding two vector functions v1 + v2 or subtracting them v1 - v2 produces a new function of time that, when invoked with .at(t), returns the sum or difference of the vectors returned by the two functions.

Planetary ephemerides

By downloading a SpiceKernel file, Skyfield users can build vector functions predicting the positions of the Moon, Sun, and planets. See Planets, and Choosing an Ephemeris.

SpiceKernel(path) Ephemeris file in NASA .bsp format.
SpiceKernel.close() Close this ephemeris file.
SpiceKernel.comments() Return the comments string of this kernel.
SpiceKernel.names() Return all target names that are valid with this kernel.
SpiceKernel.decode(name) Translate a target name into its integer code.

Kernels also support lookup using the Python kernel['Mars'] syntax, in which case they return a function of time that returns vectors from the Solar System barycenter to the named body.

Planetary magnitudes

skyfield.magnitudelib.planetary_magnitude(position)

Given the position of a planet, return its visual magnitude.

This prototype function — which so far only supports Mercury, Venus, Earth, Jupiter, and Uranus — computes the visual magnitude of a planet, given its position relative to an observer.

>>> from skyfield.api import load
>>> from skyfield.magnitudelib import planetary_magnitude
>>> ts = load.timescale()
>>> t = ts.utc(2020, 7, 31)
>>> eph = load('de421.bsp')
>>> astrometric = eph['earth'].at(t).observe(eph['jupiter barycenter'])
>>> print('%.2f' % planetary_magnitude(astrometric))
-2.73

The routine does not yet take into account whether the observer is facing the equator or poles of Uranus, so will only be accurate to within about 0.1 magnitudes.

Planetary reference frames

PlanetaryConstants Planetary constants manager.
Frame Planetary constants frame, for building rotation matrices.

Almanac

Routines to search for events like sunrise, sunset, and Moon phase.

phase_angle(ephemeris, body, t) Compute the phase angle of a body viewed from Earth.
fraction_illuminated(ephemeris, body, t) Compute the illuminated fraction of a body viewed from Earth.
seasons(ephemeris) Build a function of time that returns the quarter of the year.
sunrise_sunset(ephemeris, topos) Build a function of time that returns whether the sun is up.
dark_twilight_day(ephemeris, topos) Build a function of time returning whether it is dark, twilight, or day.
moon_phases(ephemeris) Build a function of time that returns the moon phase 0 through 3.

Topocentric locations

You can create a vector function that computes the location of any position on the Earth’s surface.

Topos([latitude, longitude, …]) A vector function that knows the position of a place on Earth.

Kepler orbits

See Kepler Orbits for computing the positions of comets, asteroids, and other minor planets.

Earth satellites

By downloading TLE satellite element sets, Skyfield users can build vector functions that predict their positions. See Earth Satellites.

EarthSatellite(line1, line2[, name, ts]) An Earth satellite loaded from a TLE file and propagated with SGP4.
EarthSatellite.from_satrec(satrec, ts) Build an EarthSatellite from a raw sgp4 Satrec object.

Stars and other distant objects

Star The position in the sky of a star or other fixed object.

Astronomical positions

The ICRF three-dimensional position vector serves as the base class for all of the following position classes. Each class represents an (x,y,z) .position and .velocity in the International Terrestrial Reference Frame (ITRF), an inertial system that is an update to J2000 and that does not rotate with the Earth itself.

ICRF An (x,y,z) position and velocity oriented to the ICRF axes.
Barycentric An (x,y,z) position measured from the Solar System barycenter.
Astrometric An astrometric (x,y,z) position relative to a particular observer.
Apparent An apparent [x y z] position relative to a particular observer.
Geocentric An (x,y,z) position measured from the center of the Earth.

Positions are usually generated by the at(t) method of a vector function, rather than being constructed manually. But you can also build a position directly from a raw vector, or from right ascension and declination coordinates with skyfield.positionlib.position_of_radec().

position_of_radec(ra_hours, dec_degrees[, …]) Build a position object from a right ascension and declination.

Position methods and attributes

All position objects have three attributes which provide access to their raw data:

ICRF.t The Time of the position.
ICRF.position A Distance array giving x, y, z.
ICRF.velocity A Velocity array giving ẋ, ẏ, ż.

If a position lacks a velocity, then the attribute is simply the value None.

All positions support a basic set of methods:

ICRF.distance() Compute the distance from the origin to this position.
ICRF.speed() Compute the magnitude of the velocity vector.
ICRF.radec([epoch]) Compute equatorial (RA, declination, distance)
ICRF.separation_from(another_icrf) Return the angle between this position and another.
ICRF.ecliptic_xyz([epoch]) Compute J2000 ecliptic position vector (x,y,z).
ICRF.ecliptic_velocity() Compute J2000 ecliptic velocity vector (x_dot, y_dot, z_dot)
ICRF.ecliptic_latlon([epoch]) Compute J2000 ecliptic coordinates (lat, lon, distance)
ICRF.galactic_xyz() Compute galactic coordinates (x,y,z)
ICRF.galactic_latlon() Compute galactic coordinates (lat, lon, distance)
ICRF.is_sunlit(ephemeris) Return whether a position in Earth orbit is in sunlight.
ICRF.from_altaz([alt, az, alt_degrees, …]) Generate an Apparent position from an altitude and azimuth.

Position methods specific to one class

Barycentric.observe(body) Compute the Astrometric position of a body from this location.
Astrometric.apparent() Compute an Apparent position for this body.
Apparent.altaz([temperature_C, pressure_mbar]) Compute (alt, az, distance) relative to the observer’s horizon
Geocentric.subpoint() Return the latitude and longitude directly beneath this position.

Constellations

skyfield.api.load_constellation_map()

Load Skyfield’s constellation boundaries and return a lookup function.

Skyfield carries an internal map of constellation boundaries that is optimized for quick position lookup. Call this function to load the map and return a function mapping position to constellation name.

>>> from skyfield.api import position_of_radec, load_constellation_map
>>> constellation_at = load_constellation_map()
>>> north_pole = position_of_radec(0, 90)
>>> constellation_at(north_pole)
'UMi'

If you pass an array of positions, you’ll receive an array of names.

skyfield.data.stellarium.parse_constellations(lines)

Return a list of constellation outlines.

Each constellation outline is a list of edges, each of which is drawn between a pair of specific stars:

[
    (name, [(star1, star2), (star3, star4), ...]),
    (name, [(star1, star2), (star3, star4), ...]),
    ...
]

Each name is a 3-letter constellation abbreviation; each star is an integer Hipparcos catalog number. See Drawing a finder chart for comet NEOWISE for an example of how to combine this data with the Hipparcos star catalog to draw constellation lines on a chart.

skyfield.data.stellarium.parse_star_names(lines)

Return the names in a Stellarium star_names.fab file.

Returns a list of named tuples, each of which offers a .hip attribute with a Hipparcos catalog number and a .name attribute with the star name. Do not depend on the tuple having only length two; additional fields may be added in the future.

Searching

skyfield.searchlib.find_discrete()

Find the times at which a discrete function of time changes value.

This routine is used to find instantaneous events like sunrise, transits, and the seasons. See Searching for the dates of astronomical events for how to use it yourself.

skyfield.searchlib.find_maxima()

Find the local maxima in the values returned by a function of time.

This routine is used to find events like highest altitude and maximum elongation. See Searching for the dates of astronomical events for how to use it yourself.

skyfield.searchlib.find_minima()

Find the local minima in the values returned by a function of time.

This routine is used to find events like minimum elongation. See Searching for the dates of astronomical events for how to use it yourself.

Osculating orbital elements

This routine returns osculating orbital elements for an object’s instantaneous position and velocity.

osculating_elements_of(position[, …]) Produce the osculating orbital elements for a position.
OsculatingElements.apoapsis_distance Distance object
OsculatingElements.argument_of_latitude Angle object
OsculatingElements.argument_of_periapsis Angle object
OsculatingElements.eccentric_anomaly Angle object
OsculatingElements.eccentricity numpy.ndarray
OsculatingElements.inclination Angle object
OsculatingElements.longitude_of_ascending_node Angle object
OsculatingElements.longitude_of_periapsis Angle object
OsculatingElements.mean_anomaly Angle object
OsculatingElements.mean_longitude Angle object
OsculatingElements.mean_motion_per_day Angle object
OsculatingElements.periapsis_distance Distance object
OsculatingElements.periapsis_time Time object
OsculatingElements.period_in_days numpy.ndarray
OsculatingElements.semi_latus_rectum Distance object
OsculatingElements.semi_major_axis Distance object
OsculatingElements.semi_minor_axis Distance object
OsculatingElements.time Time object
OsculatingElements.true_anomaly Angle object
OsculatingElements.true_longitude Angle object

Units

Distance Distance measure.
Distance.au Astronomical Units.
Distance.km Kilometers.
Distance.m Meters.
Velocity Velocity measure.
Velocity.au_per_d Astronomical Units.
Velocity.km_per_s Kilometers.
Angle Angle measure.
Angle.degrees Degrees of arc (360 in a complete circle).
Angle.hours Hours of arc (24 in a complete circle).
Angle.radians Radians (τ = 2π in a complete circle).

All three kinds of quantity support one or more methods.

Distance.length() Compute the length when this is an x,y,z vector.
Distance.light_seconds() Return the length of this vector in light seconds.
Distance.to(unit) Convert this distance to the given AstroPy unit.
Velocity.to(unit) Convert this velocity to the given AstroPy unit.
Angle.arcminutes() Return the angle in arcminutes.
Angle.arcseconds() Return the angle in arcseconds.
Angle.mas() Return the angle in milliarcseconds.
Angle.to(unit) Convert this angle to the given AstroPy unit.
Angle.hms([warn]) Convert to a tuple (hours, minutes, seconds).
Angle.signed_hms([warn]) Convert to a tuple (sign, hours, minutes, seconds).
Angle.hstr([places, warn]) Convert to a string like 12h 07m 30.00s.
Angle.dms([warn]) Convert to a tuple (degrees, minutes, seconds).
Angle.signed_dms([warn]) Convert to a tuple (sign, degrees, minutes, seconds).
Angle.dstr([places, warn]) Convert to a string like 181deg 52' 30.0".

Trigonometry

position_angle_of(anglepair1, anglepair2) Return the position angle of one position with respect to another.