Skyfield: HomeTable of ContentsChangelogAPI Reference

API Reference — Astronomical Positions

See Positions for a detailed guide to these various kinds of position that Skyfield can compute, and to the selection of coordinate systems that can be used to express them.

Generic ICRF position

class skyfield.positionlib.ICRF(position_au, velocity_au_per_d=None, t=None, center=None, target=None)

An (x,y,z) position and velocity oriented to the ICRF axes.

The International Celestial Reference Frame (ICRF) is a permanent reference frame that is the replacement for J2000. Their axes agree to within 0.02 arcseconds. It also supersedes older equinox-based systems like B1900 and B1950.

Each instance of this class provides a .position vector and a .velocity vector that specify (x,y,z) coordinates along the axes of the ICRF. A specific time .t might be specified or might be None.

t

The Time coordinate of this position.

position

A Distance object offering the position’s (x,y,z) coordinates.

velocity

A Velocity object offering the velocity’s (dx/dt,dy/dt,dz/dt) coordinates.

This attribute will have the value None if no velocity was specified for this position.

classmethod from_time_and_frame_vectors(t, frame, distance, velocity)

Constructor: build a position from two vectors in a reference frame.

distance()

Compute the distance from the origin to this position.

The return value is a Distance that prints itself out in astronomical units (au) but that also offers attributes au, km, and m if you want to access its magnitude as a number.

>>> v = ICRF([1, 1, 0])
>>> print(v.distance())
1.41421 au
speed()

Compute the magnitude of the velocity vector.

>>> v = ICRF([0, 0, 0], [1, 2, 3])
>>> print(v.speed())
3.74166 au/day
light_time

Length of this vector in days of light travel time.

radec(epoch=None)

Compute equatorial RA, declination, and distance.

When called without a parameter, this returns standard ICRF right ascension and declination:

>>> from skyfield.api import load
>>> ts = load.timescale()
>>> t = ts.utc(2020, 5, 13, 10, 32)
>>> eph = load('de421.bsp')
>>> astrometric = eph['earth'].at(t).observe(eph['sun'])
>>> ra, dec, distance = astrometric.radec()
>>> print(ra, dec, sep='\n')
03h 21m 47.67s
+18deg 28' 55.3"

If you instead want the coordinates referenced to the dynamical system defined by the Earth’s true equator and equinox, provide a specific epoch time.

>>> ra, dec, distance = astrometric.apparent().radec(epoch='date')
>>> print(ra, dec, sep='\n')
03h 22m 54.73s
+18deg 33' 04.5"

To get J2000.0 coordinates, simply pass ts.J2000.

hadec()

Compute hour angle, declination, and distance.

Returns a tuple of two Angle objects plus the Distance to the target. The angles are the hour angle (±12 hours) east or west of your meridian along the ITRS celestial equator, and the declination (±90 degees) above or below it. This only works for positions whose center is a geographic location; otherwise, there is no local meridian from which to measure the hour angle.

Because this declination is measured from the plane of the Earth’s physical geographic equator, it will be slightly different than the declination returned by radec() if you have loaded a Polar Motion file.

The coordinates are not adjusted for atmospheric refraction near the horizon.

altaz(temperature_C=None, pressure_mbar='standard')

Compute (alt, az, distance) relative to the observer’s horizon

The altitude returned is an Angle measured in degrees above the horizon, while the azimuth Angle measures east along the horizon from geographic north (so 0 degrees means north, 90 is east, 180 is south, and 270 is west).

By default, Skyfield does not adjust the altitude for atmospheric refraction. If you want Skyfield to estimate how high the atmosphere might lift the body’s image, give the argument temperature_C either the temperature in degrees centigrade, or the string 'standard' (in which case 10°C is used).

When calculating refraction, Skyfield uses the observer’s elevation above sea level to estimate the atmospheric pressure. If you want to override that value, simply provide a number through the pressure_mbar parameter.

separation_from(another_icrf)

Return the angle between this position and another.

>>> from skyfield.api import load
>>> ts = load.timescale()
>>> t = ts.utc(2020, 4, 18)
>>> eph = load('de421.bsp')
>>> sun, venus, earth = eph['sun'], eph['venus'], eph['earth']
>>> e = earth.at(t)
>>> s = e.observe(sun)
>>> v = e.observe(venus)
>>> print(s.separation_from(v))
43deg 23' 23.1"

You can also compute separations across an array of positions.

>>> t = ts.utc(2020, 4, [18, 19, 20])
>>> e = earth.at(t)
>>> print(e.observe(sun).separation_from(e.observe(venus)))
3 values from 43deg 23' 23.1" to 42deg 49' 46.6"
cirs_xyz(epoch)

Compute cartesian CIRS coordinates at a given epoch (x,y,z).

Calculate coordinates in the Celestial Intermediate Reference System (CIRS), a dynamical coordinate system referenced to the Celestial Intermediate Origin (CIO). As this is a dynamical system it must be calculated at a specific epoch.

cirs_radec(epoch)

Get spherical CIRS coordinates at a given epoch (ra, dec, distance).

Calculate coordinates in the Celestial Intermediate Reference System (CIRS), a dynamical coordinate system referenced to the Celestial Intermediate Origin (CIO). As this is a dynamical system it must be calculated at a specific epoch.

frame_xyz(frame)

Return this position as an (x,y,z) vector in a reference frame.

Returns a Distance object giving the (x,y,z) of this position in the given frame. See Coordinates in other reference frames.

frame_xyz_and_velocity(frame)

Return (x,y,z) position and velocity vectors in a reference frame.

Returns two vectors in the given coordinate frame: a Distance providing an (x,y,z) position and a Velocity giving (xdot,ydot,zdot) velocity. See Coordinates in other reference frames.

frame_latlon(frame)

Return longitude, latitude, and distance in the given frame.

Returns a 3-element tuple giving the latitude and longitude as Angle objects and the range to the target as a Distance. See Coordinates in other reference frames.

frame_latlon_and_rates(frame)

Return a reference frame longitude, latitude, range, and rates.

Return a 6-element tuple of 3 coordinates and 3 rates-of-change for this position in the given reference frame:

If the reference frame is the ICRS, or is J2000, or otherwise involves the celestial equator and pole, then the latitude and longitude returned will measure what are more commonly called “declination” and “right ascension”. Note that right ascension is usually expressed as hours (24 in a circle), rather than in the degrees that this routine will return.

to_skycoord(unit=None)

Convert this distance to an AstroPy SkyCoord object.

Currently, this will only work with Skyfield positions whose center is the Solar System barycenter or else the geocenter.

phase_angle(sun)

Return this position’s phase angle: the angle Sun-target-observer.

Given a Sun object (which you can build by loading an ephemeris and looking up eph['Sun']), return the Angle from the body’s point of view between light arriving from the Sun and the light departing toward the observer. This angle is 0° if the observer is in the same direction as the Sun and sees the body as fully illuminated, and 180° if the observer is behind the body and sees only its dark side.

New in version 1.42.

fraction_illuminated(sun)

Return the fraction of the target’s disc that is illuminated.

Given a Sun object (which you can build by loading an ephemeris and looking up eph['Sun']), compute what fraction from 0.0 to 1.0 of this target’s disc is illuminated, under the assumption that the target is a sphere.

New in version 1.42.

is_sunlit(ephemeris)

Return whether a position in Earth orbit is in sunlight.

Returns True or False, or an array of such values, to indicate whether this position is in sunlight or is blocked by the Earth’s shadow. It should work with positions produced either by calling at() on a satellite object, or by calling at() on the relative position sat - topos of a satellite with respect to an Earth observer’s position. See Find when a satellite is in sunlight.

is_behind_earth()

Return whether the Earth blocks the view of this object.

For a position centered on an Earth-orbiting satellite, return whether the target is in eclipse behind the disc of the Earth. See Find whether the Earth blocks a satellite’s view.

from_altaz(alt=None, az=None, alt_degrees=None, az_degrees=None, distance=<Distance 0.1 au>)

Generate an Apparent position from an altitude and azimuth.

The altitude and azimuth can each be provided as an Angle object, or else as a number of degrees provided as either a float or a tuple of degrees, arcminutes, and arcseconds:

alt=Angle(...), az=Angle(...)
alt_degrees=23.2289, az_degrees=142.1161
alt_degrees=(23, 13, 44.1), az_degrees=(142, 6, 58.1)

The distance should be a Distance object, if provided; otherwise a default of 0.1 au is used.

Position measured from the Solar System barycenter

class skyfield.positionlib.Barycentric(position_au, velocity_au_per_d=None, t=None, center=None, target=None)

An (x,y,z) position measured from the Solar System barycenter.

Skyfield generates a Barycentric position measured from the gravitational center of the Solar System whenever you ask a body for its location at a particular time:

>>> t = ts.utc(2003, 8, 29)
>>> mars.at(t)
<Barycentric BCRS position and velocity at date t center=0 target=499>

This class’s .position and .velocity are (x,y,z) vectors in the Barycentric Celestial Reference System (BCRS), the modern replacement for J2000 coordinates measured from the Solar System Barycenter.

This class inherits the methods of is parent class ICRF as well as the orientation of its axes in space.

observe(body)

Compute the Astrometric position of a body from this location.

To compute the body’s astrometric position, it is first asked for its position at the time t of this position itself. The distance to the body is then divided by the speed of light to find how long it takes its light to arrive. Finally, the light travel time is subtracted from t and the body is asked for a series of increasingly exact positions to learn where it was when it emitted the light that is now reaching this position.

>>> earth.at(t).observe(mars)
<Astrometric ICRS position and velocity at date t center=399 target=499>

Astrometric position relative to an observer

class skyfield.positionlib.Astrometric(position_au, velocity_au_per_d=None, t=None, center=None, target=None)

An astrometric (x,y,z) position relative to a particular observer.

The astrometric position of a body is its position relative to an observer, adjusted for light-time delay. It is the position of the body back when it emitted (or reflected) the light that is now reaching the observer’s eye or telescope. Astrometric positions are usually generated in Skyfield by calling the Barycentric method observe(), which performs the light-time correction.

Both the .position and .velocity are (x,y,z) vectors oriented along the axes of the ICRF, the modern replacement for the J2000 reference frame.

It is common to either call .radec() (with no argument) on an astrometric position to generate an astrometric place right ascension and declination with respect to the ICRF axes, or else to call .apparent() to generate an Apparent position.

This class inherits the methods of is parent class ICRF as well as the orientation of its axes in space.

altaz()

Compute (alt, az, distance) relative to the observer’s horizon

The altitude returned is an Angle measured in degrees above the horizon, while the azimuth Angle measures east along the horizon from geographic north (so 0 degrees means north, 90 is east, 180 is south, and 270 is west).

By default, Skyfield does not adjust the altitude for atmospheric refraction. If you want Skyfield to estimate how high the atmosphere might lift the body’s image, give the argument temperature_C either the temperature in degrees centigrade, or the string 'standard' (in which case 10°C is used).

When calculating refraction, Skyfield uses the observer’s elevation above sea level to estimate the atmospheric pressure. If you want to override that value, simply provide a number through the pressure_mbar parameter.

apparent()

Compute an Apparent position for this body.

This applies two effects to the position that arise from relativity and shift slightly where the other body will appear in the sky: the deflection that the image will experience if its light passes close to large masses in the Solar System, and the aberration of light caused by the observer’s own velocity.

>>> earth.at(t).observe(mars).apparent()
<Apparent GCRS position and velocity at date t center=399 target=499>

These transforms convert the position from the BCRS reference frame of the Solar System barycenter and to the reference frame of the observer. In the specific case of an Earth observer, the output reference frame is the GCRS.

Apparent position relative to an observer

class skyfield.positionlib.Apparent(position_au, velocity_au_per_d=None, t=None, center=None, target=None)

An apparent (x,y,z) position relative to a particular observer.

This class’s vectors provide the position and velocity of a body relative to an observer, adjusted to predict where the body’s image will really appear (hence “apparent”) in the sky:

  • Light-time delay, as already present in an Astrometric position.
  • Deflection: gravity bends light, and thus the image of a distant object, as the light passes massive objects like Jupiter, Saturn, and the Sun. For an observer on the Earth’s surface or in Earth orbit, the slight deflection by the gravity of the Earth itself is also included.
  • Aberration: incoming light arrives slanted because of the observer’s motion through space.

These positions are usually produced in Skyfield by calling the apparent() method of an Astrometric object.

Both the .position and .velocity are (x,y,z) vectors oriented along the axes of the ICRF, the modern replacement for the J2000 reference frame. If the observer is at the geocenter, they are more specifically GCRS coordinates. Two common coordinates that this vector can generate are:

  • Proper place: call .radec() without arguments to compute right ascension and declination with respect to the fixed axes of the ICRF.
  • Apparent place, the most popular option: call .radec('date') to generate right ascension and declination with respect to the equator and equinox of date.

This class inherits the methods of is parent class ICRF as well as the orientation of its axes in space.

Geocentric position relative to the Earth

class skyfield.positionlib.Geocentric(position_au, velocity_au_per_d=None, t=None, center=None, target=None)

An (x,y,z) position measured from the center of the Earth.

A geocentric position is the difference between the position of the Earth at a given instant and the position of a target body at the same instant, without accounting for light-travel time or the effect of relativity on the light itself.

Its .position and .velocity vectors have (x,y,z) axes that are those of the Geocentric Celestial Reference System (GCRS), an inertial system that is an update to J2000 and that does not rotate with the Earth itself.

This class inherits the methods of is parent class ICRF as well as the orientation of its axes in space.

itrf_xyz()

Deprecated; instead, call .frame_xyz(itrs). See Coordinates in other reference frames.

subpoint()

Deprecated; instead, call either iers2010.subpoint(pos) or wgs84.subpoint(pos).

Building a position from right ascension and declination

skyfield.positionlib.position_of_radec(ra_hours, dec_degrees, distance_au=206264806247096.38, epoch=None, t=None, center=None, target=None)

Build a position object from a right ascension and declination.

If a specific distance_au is not provided, Skyfield returns a position vector a gigaparsec in length. This puts the position at a great enough distance that it will stand at the same right ascension and declination from any viewing position in the Solar System, to very high precision (within a few hundredths of a microarcsecond).

If an epoch is specified, the input coordinates are understood to be in the dynamical system of that particular date. Otherwise, they will be assumed to be ICRS (the modern replacement for J2000).

New in version 1.21: This replaces a deprecated function position_from_radec() whose distance argument was not as well designed.