Skyfield: HomeTable of ContentsAPI Reference

API Reference

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.path_to(filename) Return the path to filename in this loader’s directory.
Loader.timescale([delta_t]) Open or download three time scale files, returning a Timescale.
Loader.tle(url[, reload]) Load and parse a satellite TLE file.

Time scales

A Skyfield Timescale object is typically built at the beginning of each program:

from skyfield import api
ts = api.load_timescale()

It downloads and parses the data tables necessary to correctly convert between Universal Time and the more stable time scales used by astronomers.

Timescale.now() Return the current date and time as a Time object.
Timescale.utc(year[, month, day, hour, …]) Return the Time corresponding to a specific moment in UTC.
Timescale.tai([year, month, day, hour, …]) Return the Time corresponding to a specific moment in TAI.
Timescale.tt([year, month, day, hour, …]) Return the Time corresponding to a specific moment in TT.
Timescale.tdb([year, month, day, hour, …]) Return the Time corresponding to a specific moment in TDB.
Timescale.ut1([year, month, day, hour, …]) Return the Time corresponding to a specific moment in UT1.
Timescale.from_astropy(t) Return a Skyfield time corresponding to the AstroPy time t.

Time objects

The Time class is Skyfield’s way of representing either a single time, or a whole array of times. Each object has four floating point attributes that present the time in several basic 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 decimal 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 an A.D.2014-Jan-18 01:35:37.5000 UT string..
Time.utc_iso([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 date 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 pytz provided timezone.
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() Return TAI as a tuple (year, month, day, hour, minute, second).
Time.tt_calendar() Return TT as a tuple (year, month, day, hour, minute, second).

Vector Functions

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

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

Planetary Ephemerides

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

SpiceKernel(path) Ephemeris file in NASA .bsp format.
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.
SpiceKernel.__getitem__(target) Return a vector function for computing the location of target.

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.

Earth Satellites

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

EarthSatellite(line1, line2[, name, ts]) An Earth satellite loaded from a TLE file and propagated with SGP4.

Stars and other distant objects

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

Astronomical positions

The ICRF class serves as the base for all other positions classes, which each share its axes but have more specific meanings. Positions are usually returned by calling at(t) on a vector function rather than being constructed manually.

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

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_position() Compute J2000 ecliptic coordinates (x, y, z)
ICRF.ecliptic_latlon() Compute J2000 ecliptic coordinates (lat, lon, distance)
ICRF.galactic_position() Compute galactic coordinates (x, y, z)
ICRF.galactic_latlon() Compute galactic coordinates (lat, lon, distance)
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.

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.to(unit) Convert this distance to the given AstroPy unit.
Velocity.to(unit) Convert this velocity to the given AstroPy unit.
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 (degrees, hours, minutes, seconds).
Angle.dstr([places, warn]) Convert to a string like 181deg 52' 30.0".