Along with
integrated solar system ephemerides, the PixInsight 1.8.6 platform includes new functionalities for the reduction of positions of solar system bodies and stars. In this document I am going to describe the JavaScript implementation at a basic level, sufficient for the practical calculation of astrometric, proper and apparent positions, among others. These reduced positions can be used for a variety of powerful applications, such as high-precision astrometry, drawing annotations on plate solved images, solar system simulations, and generation of high-level ephemerides such as occultations of stars by asteroids and the Moon, just to name a few.
The JavaScript implementation, available through the PJSR runtime in the PixInsight core application, uses the PCL/C++ implementation as its back end. For information on the entire functionality available with detailed descriptions of the implemented algorithms, see the corresponding sections of the
official PCL documentation:
C++ documentation for the integrated solar system ephemerides system.C++ documentation for the Position class.
StarPositionThis object is a simple structure to define positional data of a star. Positions and proper motions must be referred to ICRS/J2000.0. If out-of-range coordinates are specified, their values will be constrained to their proper ranges: right ascension to [0°,360°) and declination to [-90°,+90°].
Number StarPosition.alphaICRS right ascension in degrees. Zero if not specified.
Number StarPosition.deltaICRS declination in degrees. Zero if not specified.
Date StarPosition.epochEpoch of catalog coordinates. The fundamental J2000.0 epoch will be assumed if not specified.
Number StarPosition.muAlphaProper motion in right ascension, in milliarcseconds per year, multiplied by the cosine of the declination. Zero if not specified.
Number StarPosition.muDeltaProper motion in declination, in mas/year. Zero if not specified.
Number StarPosition.parallaxParallax in arcseconds. Zero if not specified.
Number StarPosition.radialVelocityRadial velocity in km/s, positive away from Earth. Zero if not specified.
ObserverPositionGeodetic coordinates of a terrestrial observer.
This structure provides the data necessary to calculate topocentric places of solar system bodies and stars. Typically, the values used should be WGS84 coordinates (for example, as distributed by GPS) or ITRF coordinates—both systems close together at the level of a few centimeters.
Geographic longitudes grow eastward, so longitudes are positive to the east and negative to the west of the central meridian. For ease of calculation, the
ObserverPosition.longitude property converts assigned longitudes to the [0°,360°) range.
If out-of-range coordinates are specified, their values will be constrained to their proper ranges: longitude to [0°,360°), and latitude to [-90°,+90°].
Boolean ObserverPosition.cioBasedWhether to compute CIO-based positions (if true) or equinox-based positions (if false). False by default, that is, equinox-based topocentric coordinates are computed by default.
Number ObserverPosition.equatorialRadiusEquatorial radius of the reference spheroid in meters. By default this property is the IAU 2009/2012 equatorial radius of Earth: 6,378,136.6 m.
Number ObserverPosition.flatteningFlattening of the reference spheroid. The default value is the IERS 2010 flattening of Earth: 1/298.25642.
Number ObserverPosition.heightGeodetic height in meters. Zero by default.
Number ObserverPosition.latitudeGeodetic latitude in degrees. Zero if not specified.
Number ObserverPosition.longitudeGeodetic longitude in degrees. Zero if not specified.
Vector ObserverPosition.regionalCenterGeocentric rectangular coordinates of the center of the regional spheroid, in meters. A zero 3-D vector (that is, the geocenter) by default.
PositionThis object implements algorithms for reduction of positions of solar system bodies and stars. It allows for calculation of geometric, astrometric, proper, apparent and intermediate places, including geocentric and topocentric coordinates.
The implemented vector astrometry and ephemeris calculation algorithms are rigorous and compliant with current IAU and IERS resolutions. Both equinox-based and CIO-based paradigms have been implemented for calculation of positions that depend on Earth's rotation. The apparent and intermediate places include the following corrections:
- Light-travel time for solar system bodies.
- Space motion for stars, including parallax, radial velocity and proper motions, with corrections for the relativistic Doppler effect that takes into account the change in light-travel time for stars.
- Relativistic deflection of light due to solar gravitation.
- Aberration of light, relativistic model.
- Frame bias, precession and nutation. (equinox-based and CIO-based).
- Accurate topocentric places with polar motion corrections.
Vector components are expressed in astronomical units (au) for stars and all solar system bodies except the Moon, for which positions are given in kilometers.
As of writing this document (January 2019), the IAU 2006/2000A precession-nutation theory is implemented (adjusted model with corrections to nutation angles, IAU 2006/2000A
R). The standard fundamental ephemerides are JPL's DE438/LE438.
Constructor new Position( Date t[, String timescale="TT"] )
new Position( String isoTime[, String timescale="TT"] )
new Position( Number jd1[, Number jd2=0[, String timescale="TT"]] )Constructs a new object for reduction of positions at the specified time point, which can be provided as a
Date object, as an ISO 8601 date/time string, or as a Julian date with one or two components, JD=jd1+jd2. The timescale argument can be one of the following:
TT | Terrestrial Time. This is the default timescale. |
TDB | Barycentric Dynamical Time. |
Teph | Ephemeris time, as defined by JPL DE/LE numerical integrations. For all purposes, this is equivalent to TDB. |
UTC | Coordinated Universal Time. |
TAI | Atomic International Time. |
UT1 | Universal Time. |
UT | Universal Time (same as UT1). |
Properties ObserverPosition Position.observerThis property defines the geodetic positional and reference data necessary for calculation of topocentric positions of solar system objects and stars.
By default, an instance of the
Position object calculates geocentric coordinates. After assigning a value to this property, subsequently calculated geometric, proper, astrometric, apparent and intermediate places will be topocentric, that is, will be referred to the location of the observer with respect to the center of the Earth.
By assigning a value to this property, all previously computed positional data will be erased with the exception of fundamental ephemerides and existing bias-precession-nutation matrices, which can always be preserved.
The
ObserverPosition.cioBased property of the assigned object selects the celestial system and paradigm to be used in the calculation of positions that depend on Earth's rotation. If the
cioBased property is true, equinox-based places cannot be calculated, and any of the apparent() methods will throw an exception if called. Conversely, if
cioBased is false, CIO-based places cannot be calculated and no call to an intermediate() method will be allowed.
If polar motion corrections are enabled, the position of the Celestial Intermediate Pole (CIP) with respect to the ITRS is interpolated from the global CIP_ITRS database, if it provides data for the current time of calculation. In such case, polar motion is taken into account in the calculation of the observer's geocentric position and velocity. For the geocentric velocity a standard constant value for the rotation rate of the Earth is used; the velocity component due to precession and nutation is not taken into account since its effect is negligible. See the
polarMotionEnabled property for more information and additional references.
Boolean Position.polarMotionEnabledThis property is true if topocentric places take into account polar motion corrections to compute the geocentric position and velocity of the observer. This involves calculation of CIP (Celestial Intermediate Pole) coordinates with respect to the ITRS, as well as access to a database of CIP/ITRS positions, which is part of the official PixInsight distribution and is updated regularly as new observational data become available.
Polar motion introduces changes at the mas level for calculation of topocentric coordinates of the Moon. For the rest of objects, the effect of polar motion corrections is completely negligible. For topocentric positions of the Moon, polar motion can be necessary to achieve the highest accuracy, but in such case one may also have to take into account a regional geoid referred to the Earth's reference ellipsoid. See the
ObserverPosition object.
Number Position.epsAMean obliquity of the ecliptic, in radians.
Methods Vector Position.apparent( EphemerisHandle H )Computes the apparent place of a solar system body. The returned vector is the apparent place of the object defined by the specified ephemeris handle
H in geocentric or topocentric rectangular equatorial coordinates, calculated for the instant defined by this object in the TT timescale. The implemented reduction algorithm includes the following corrections:
- Light-travel time.
- Relativistic deflection of light due to solar gravitation (except for the Sun, the Moon, and any object closer from Earth than the Sun at the time of observation.
- Aberration of light, including relativistic terms.
- Frame bias, precession and nutation. The origin of right ascension is the true equinox of date.
Vector Position.apparent( StarPosition S )Computes the apparent place of a star. The returned vector is the apparent place calculated for the positional star data defined by the specified object
S. Its components are geocentric or topocentric rectangular equatorial coordinates, calculated for the instant defined by this object in the TT timescale. The implemented reduction algorithm includes the following corrections:
- Space motion, including parallax, radial velocity and proper motions, with corrections for the relativistic Doppler effect.
- Relativistic deflection of light due to solar gravitation.
- Aberration of light, including relativistic terms.
- Frame bias, precession and nutation. The origin of right ascension is the true equinox of date.
Number|null Position.apparentVisualMagnitude( EphemerisHandle H )Returns the observed visual magnitude of a solar system body.
For objects with known H and G values (absolute magnitude and slope parameters, respectively; see
EphemerisHandle.H and
EphemerisHandle.G), the apparent visual magnitude is calculated applying the algorithm for minor planets described in Bowell et al. (1989). See also The Explanatory Supplement, Section 10.4.3.
For Mercury, Venus, Mars, Jupiter, Saturn and Neptune, we apply the equations described in the following paper:
For Saturn, we compute the apparent visual magnitude taking into account the planet's rings.
For Uranus, Pluto and the Galilean satellites of Jupiter, data from various sources are taken from Table 10.6 of the Explanatory Supplement.
If the required data are not available, or if no algorithm is known for the calculation of the apparent visual magnitude of the specified object, this method returns
null.
A
null value is also returned when the phase angle of the object at the time of calculation is beyond the limits of the set of observations used to generate the underlying models. For Mercury, apparent magnitudes are only calculated for phase angles 2° ?
i ? 170°. For Venus, the magnitude is only calculated for 0° <
i ? 163°.7. The valid range for Mars is
i ? 50°.
See also
Position.canComputeApparentVisualMagnitude() Vector Position.astrometric( EphemerisHandle H )Computes the astrometric place of a solar system body. The returned vector is the astrometric place of the object defined by the specified ephemeris handle
H in geocentric or topocentric rectangular equatorial coordinates, calculated for the instant defined by this object in the TT timescale. The implemented reduction algorithm includes the following corrections:
- Light-travel time.
- Relativistic deflection of light due to solar gravitation (except for the Sun, the Moon, and any object closer from Earth than the Sun at the time of observation).
An astrometric place does not include annual aberration, nutation and precession corrections. Hence it is referred to an 'hybrid' reference system, but similar to GCRS J2000.0.
Vector Position.astrometric( StarPosition S )Computes the astrometric place of a star. The returned vector is the astrometric place calculated for the positional star data defined by the specified object
S. Its components are geocentric or topocentric rectangular equatorial coordinates, calculated for the instant defined by this object in the TT timescale. The implemented reduction algorithm includes the following corrections:
- Space motion, including parallax, radial velocity and proper motions, with corrections for the relativistic Doppler effect.
- Relativistic deflection of light due to solar gravitation.
An astrometric place does not include annual aberration, nutation and precession corrections. Hence it is referred to an 'hybrid' reference system, but similar to GCRS J2000.0.
Boolean Position.canComputeApparentVisualMagnitude( EphemerisHandle H )Returns true only if the apparent visual magnitude of the object represented by the specified handle
H can be calculated with the current implementation, at the calculation time defined by this
Position object.
Currently apparent visual magnitudes can be calculated for the following solar system bodies:
- Objects providing valid H and G parameters (absolute magnitude and slope coefficient). This is true for most asteroids included in standard XEPH files.
- Mercury, Venus, Mars, Jupiter, Saturn, Uranus, Neptune and Pluto.
- The four Galilean satellites of Jupiter: Io, Europa, Ganymede and Callisto.
Vector Position.geometric( EphemerisHandle H )Computes the geometric position of a solar system body. The components of the returned vector are the geocentric or topocentric rectangular equatorial coordinates for the instant of calculation defined by this object in the TT timescale, accounting for light-travel time, for the body defined by the specified ephemeris handle
H.
The implemented reduction algorithm includes just the correction for light-travel time, but no corrections for light deflection, annual aberration, nutation, or precession. The position so calculated allows to plot the specified body directly on an existing sky chart referred to GCRS/J2000.0. Note however, that for generation of new graphical representations for a given date using star catalog data, astrometric or proper places should be used instead.
Vector Position.geometric( StarPosition S )Computes the geometric position of a star. The components of the returned vector are the geocentric or topocentric rectangular equatorial coordinates for the instant of calculation defined by this object in the TT timescale, for the positional star data defined by the specified object
S.
The implemented reduction algorithm includes just the corrections for space motion: parallax, radial velocity and proper motions, when the corresponding data items are nonzero in the specified object S. The space motion vector includes terms to account for the relativistic Doppler effect.
void Position.initCIOBasedParameters()Calculates all parameters and data structures necessary for CIO-based reduction of positions. This method starts by calling
initEquinoxBasedParameters(), so it implicitly calculates all equinox-based parameters. Then it calculates the CIO-based combined bias-precession-nutation matrix.
Since all of these items depend exclusively on time, they are computed only once in the first call to this function, and subsequent calls will have no effect.
Normally, you don't have to call this function directly because it is invoked automatically when necessary by the different position reduction routines.
void Position.initEquinoxBasedParameters()Calculates all parameters and data structures necessary for equinox-based reduction of positions. This method calculates the following structures:
- Precession+bias angles, IAU 2006 precession model, Fukushima-Williams parameterization.
- Mean obliquity of the ecliptic, IAU 2006 precession model.
- Nutation angles, IAU 2006/2000AR nutation model.
- Combined bias-precession-nutation matrix, equinox-based.
- Position of the Celestial Intermediate Pole (CIP).
- Celestial Intermediate Origin (CIO) locator.
- Equation of the origins (EO).
- Earth rotation angle (ERA) for the UT1 time of calculation.
- Greenwich Apparent Sidereal Time (GAST), IAU 2006.
Since all of these items depend exclusively on time, they are computed only once in the first call to this function, and subsequent calls will have no effect.
Normally, you don't have to call this function directly because it is invoked automatically when necessary by the different position reduction routines.
Vector Position.intermediate( EphemerisHandle H )Computes the intermediate place of a solar system body. The returned vector is the intermediate place of the object defined by the specified ephemeris handle
H in geocentric or topocentric rectangular equatorial coordinates, calculated for the instant defined by this object in the TT timescale. The implemented reduction algorithm includes the following corrections:
- Light-travel time.
- Relativistic deflection of light due to solar gravitation (except for the Sun, the Moon, and any object closer from Earth than the Sun at the time of observation.
- Aberration of light, including relativistic terms.
- Frame bias, precession and nutation. The origin of right ascension is the Celestial Intermediate Origin (CIO), following the IAU recommendations since January 2003.
Vector Position.intermediate( StarPosition S )Computes the intermediate place of a star. The returned vector is the intermediate place calculated for the positional star data defined by the specified object
S. Its components are geocentric or topocentric rectangular equatorial coordinates, calculated for the instant defined by this object in the TT timescale. The implemented reduction algorithm includes the following corrections:
- Space motion, including parallax, radial velocity and proper motions, with corrections for the relativistic Doppler effect.
- Relativistic deflection of light due to solar gravitation.
- Aberration of light, including relativistic terms.
- Frame bias, precession and nutation. The origin of right ascension is the Celestial Intermediate Origin (CIO), following the IAU recommendations since January 2003.
Array Position.nutationAngles()Returns the nutation angles in radians as an
Array object
P, where
P[0] is the nutation in longitude and
P[1] is the nutation in obliquity. This method calls
initEquinoxBasedParameters() to ensure that all data required for equinox-based reduction of positions is available.
Vector Position.proper( EphemerisHandle H )Computes the proper place of a solar system body. The returned vector is the proper place of the object defined by the specified ephemeris handle
H in geocentric or topocentric rectangular equatorial coordinates, calculated for the instant defined by this object in the TT timescale. The implemented reduction algorithm includes the following corrections:
- Light-travel time.
- Relativistic deflection of light due to solar gravitation (except for the Sun, the Moon, and any object closer from Earth than the Sun at the time of observation.
- Aberration of light, including relativistic terms.
A proper place does not include nutation and precession corrections, hence it is referred to the reference coordinate system: GCRS J2000.0.
Vector Position.proper( StarPosition S )Computes the proper place of a star. The returned vector is the proper place calculated for the positional star data defined by the specified object
S. Its components are geocentric or topocentric rectangular equatorial coordinates, calculated for the instant defined by this object in the TT timescale. The implemented reduction algorithm includes the following corrections:
- Space motion, including parallax, radial velocity and proper motions, with corrections for the relativistic Doppler effect.
- Relativistic deflection of light due to solar gravitation.
- Aberration of light, including relativistic terms.
A proper place does not include nutation and precession corrections, hence it is referred to the reference coordinate system: GCRS J2000.0.
Vector Position.true( EphemerisHandle H )Computes the true position of a solar system body. The components of the returned vector are the geocentric or topocentric rectangular equatorial coordinates for the calculation instant defined by this object in the TT timescale,
without accounting for light-travel time, for the body defined by the specified ephemeris handle
H.
This function calls
Position.geometric( EphemerisHandle ) internally to compute, if necessary, the geometric position with correction for light time, that is, no separate calculation routine is implemented for true positions. The returned vector is only useful to compute the true geocentric or topocentric distance, and for verification purposes.
Vector Position.true( StarPosition S )Computes the geometric position of a star. This function has been implemented for completeness. It is a synonym for
Position.geometric( StarPosition ). There are no known 'true' positions of stars, since their light-travel time is implicitly included in the space motion equations.
Number Position.trueDistance( EphemerisHandle H )Computes the true distance of a solar system body. The true distance is the actual distance from the body to the observer, geocentric or topocentric, at the instant of calculation. This excludes the light-travel time correction.
Number Position.trueDistance( StarPosition S )Computes the true distance of a star. The returned value is just the norm of the geometric position vector, that is:
Position.geometric( S ).l2norm()This should be an actual distance in au only for stars with known parallaxes. For stars where the parallax is unknown or undefined, this value is meaningless because in such cases position vectors are unit vectors, whose components are also known as
direction cosines.Static Methods Vector Position.equatorialToEcliptic( Vector r, Number eps )Conversion from rectangular equatorial to rectangular ecliptic coordinates.
r is a vector of rectangular equatorial coordinates.
eps is the oliquity of the ecliptic in radians.
Returns a vector whose components are the rectangular ecliptic coordinates corresponding to the specified equatorial position
r at the epoch where the specified obliquity has been calculated.
Vector Position.equatorialToEcliptic( Vector r, Number se, Number ce )This method is identical to the previous one, namely
Position.equatorialToEcliptic( Vector, Number ), but the obliquity of the ecliptic is specified by its sine
se and cosine
ce.
Vector Position.icrsEquatorialToGalactic( Vector r )Conversion from ICRS rectangular equatorial to rectangular galactic coordinates.
r is a vector of rectangular equatorial coordinates in the ICRS. Returns a vector whose components are the calculated rectangular galactic coordinates.
In this routine we adopt the proposed ICRS coordinates of the galactic pole in:
- Jia-Cheng Liu, Zi Zhu, and Hong Zhang, Reconsidering the galactic coordinate system, Astronomy & Astrophysics manuscript no. AA2010, October 26, 2018.
The applied conventional definitions are as follows. The ICRS equatorial coordinates of the zero point of galactic coordinates are:
alpha = 17
h45
m40
s.0400
delta = –29°00'28".138
The equatorial coordinates of the galactic pole, coherent with the ICRS, are:
alpha
p = 12
h51
m36
s.7151981
delta
p = +27°06'11".193172
Note that these definitions are not consistent with the conventional values currently accepted by the IAU. The current (as of October 2018) galactic coordinate system was defined by the IAU in 1959 in the FK4 B1950.0 reference system.