Convert star position+velocity vector to catalog coordinates.
Status: support routine.
The specified pv-vector is the coordinate direction (and its rate of change) for the epoch at which the light leaving the star reached the solar-system barycenter.
The star data returned by this routine are "observables" for an imaginary observer at the solar-system barycenter. Proper motion and radial velocity are, strictly, in terms of barycentric coordinate time, TCB. For most practical applications, it is permissible to neglect the distinction between TCB and ordinary "proper" time on Earth (TT/TAI). The result will, as a rule, be limited by the intrinsic accuracy of the proper-motion and radial- velocity data; moreover, the supplied pv-vector is likely to be merely an intermediate result (for example generated by the routine STARPV), so that a change of time unit will cancel out overall.
In accordance with normal star-catalog conventions, the object's right ascension and declination are freed from the effects of secular aberration. The frame, which is aligned to the catalog equator and equinox, is Lorentzian and centered on the SSB.
Summarizing, the specified pv-vector is for most stars almost identical to the result of applying the standard geometrical "space motion" transformation to the catalog data. The differences, which are the subject of the Stumpff paper cited below, are:
(i) In stars with significant radial velocity and proper motion, the constantly changing light-time distorts the apparent proper motion. Note that this is a classical, not a relativistic, effect.
(ii) The transformation complies with special relativity.
Care is needed with units. The star coordinates are in radians and the proper motions in radians per Julian year, but the parallax is in arcseconds; the radial velocity is in km/s, but the pv-vector result is in au and au/day.
The proper motions are the rate of change of the right ascension and declination at the catalog epoch and are in radians per Julian year. The RA proper motion is in terms of coordinate angle, not true angle, and will thus be numerically larger at high declinations.
Straight-line motion at constant speed in the inertial frame is assumed. If the speed is greater than or equal to the speed of light, the routine aborts with an error status.
The inverse transformation is performed by the routine STARPV.
The pv
argument is documented as an input in the IAU routine,
But the velocity components are changed by this routine. In this
version, it is declared as intent(inout)
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(inout), | dimension(3,2) | :: | pv | pv-vector (au, au/day) [see Note 1] |
|
real(kind=wp), | intent(out) | :: | ra | right ascension (radians) [see Note 2] |
||
real(kind=wp), | intent(out) | :: | dec | declination (radians) [see Note 2] |
||
real(kind=wp), | intent(out) | :: | pmr | RA proper motion (radians/year) [see Note 2] |
||
real(kind=wp), | intent(out) | :: | pmd | Dec proper motion (radians/year) [see Note 2] |
||
real(kind=wp), | intent(out) | :: | px | parallax (arcsec) [see Note 2] |
||
real(kind=wp), | intent(out) | :: | rv | radial velocity (km/s, positive = receding) [see Note 2] |
||
integer, | intent(out) | :: | j | status [see Note 2]: * 0 = OK * -1 = superluminal speed (Note 5) * -2 = null position vector |
subroutine PVSTAR ( pv, ra, dec, pmr, pmd, px, rv, j )
implicit none
real(wp),dimension(3,2),intent(inout) :: pv !! pv-vector (au, au/day) [see Note 1]
real(wp),intent(out) :: ra !! right ascension (radians) [see Note 2]
real(wp),intent(out) :: dec !! declination (radians) [see Note 2]
real(wp),intent(out) :: pmr !! RA proper motion (radians/year) [see Note 2]
real(wp),intent(out) :: pmd !! Dec proper motion (radians/year) [see Note 2]
real(wp),intent(out) :: px !! parallax (arcsec) [see Note 2]
real(wp),intent(out) :: rv !! radial velocity (km/s, positive = receding) [see Note 2]
integer,intent(out) :: j !! status [see Note 2]:
!! * 0 = OK
!! * -1 = superluminal speed (Note 5)
!! * -2 = null position vector
! Julian years to days
real(wp),parameter :: y2d = 365.25_wp
! Radians to arcseconds
real(wp),parameter :: dr2as = 206264.8062470963551564734_wp
! Astronomical unit (m, IAU 2012)
real(wp),parameter :: aum = 149597870.7e3_wp
! Speed of light (au per day)
real(wp),parameter :: c = d2s*cmps/aum
real(wp) :: r, x(3), vr, ur(3), vt, ut(3), bett, betr, d, w, &
del, usr(3), ust(3), a, rad, decd, rd
! Isolate the radial component of the velocity (au/day, inertial).
call PN ( pv(1,1), r, x )
call PDP ( x, pv(1,2), vr )
call SXP ( vr, x, ur )
! Isolate the transverse component of the velocity (au/day, inertial).
call PMP ( pv(1,2), ur, ut )
call PM ( ut, vt )
! Special-relativity dimensionless parameters.
bett = vt / c
betr = vr / c
! The inertial-to-observed correction terms.
d = 1.0_wp + betr
w = betr*betr + bett*bett
if ( d==0.0_wp .or. w>=1.0_wp ) then
j = -1
return
end if
del = - w / ( sqrt(1.0_wp-w) + 1.0_wp )
! Apply relativistic correction factor to radial velocity component.
if ( betr/=0.0_wp ) then
w = ( betr-del ) / ( betr*d )
else
w = 1.0_wp
end if
call SXP ( w, ur, usr )
! Apply relativistic correction factor to tangential velocity component.
call SXP ( 1.0_wp/d, ut, ust )
! Combine the two to obtain the observed velocity vector (au/day).
call PPP ( usr, ust, pv(1,2) )
! Cartesian to spherical.
call PV2S ( pv, a, dec, r, rad, decd, rd )
if ( r == 0.0_wp ) then
j = -2
return
end if
! Return RA in range 0 to 2pi.
ra = ANP ( a )
! Return proper motions in radians per year.
pmr = rad * y2d
pmd = decd * y2d
! Return parallax in arcsec.
px = dr2as / r
! Return radial velocity in km/s.
rv = 1.0e-3_wp * rd * aum / d2s
! OK status.
j = 0
end subroutine PVSTAR