![]() |
Institut für Astronomie und AstrophysikAbteilung AstronomieSand 1, D-72076 Tübingen, Germany |
![]() |
BARYCEN
Convert geocentric modified Julian date (MJD) to helio bary centric Julian
date.
This procedure correct for the extra light travel time between the Earth
and the bary center of the Sun.
bary_time = BARYCEN( date, ra, dec, /B1950, /TIME_DIFF,
ephemfile=ephemfile, orbitfile=orbitfile)
date - Modified Julian Date (= JD - 2400000.5), scalar or vector, MUST
be double precision.
The date describes an event at the geocentric time.
The format of date can be a vector which is recommended for
computing several dates because ephemeris file will be
accessed only once.
ra,dec - scalars giving right ascension and declination of the
observed object in DEGREES
Equinox is J2000 unless the /B1950 keyword is set.
The object direction is needed to take into account
the time delay due light travelling from center to
earth to the common solar bary center.
bary_time - solar bary centric Modified Julian Date.
If /TIME_DIFF is set, then BARYCEN() instead
returns the time difference in seconds of
barycentric-geocentric Julian date.
ephemfile - (string) The file which contains JPL ephemeris
information. It must have the JPLEPHREAD specific
format. Public available are
- JPLEPH.200 - JPL-DE200 which is the older but well
known ephemeris data file.
- JPLEPH.405 - JPL-DE405 is more recent and precise.
It is possible to build one's own ephemeris
file. Refer JPLEPHREAD procedure.
This file is not needed if the pinfo and pdata input
variables are set.
If this file is not given BARYCEN looks for the
environment variable "ASTRO_DATA" describing a
directory which must contain the file "JPLEPH.405".
If the ephemeris file is not found an error is raised.
orbit - 2-dim array of double containing the orbit position
information in respect of the geocenter. The position
must be given for all dates given with the DATE
input. Position must be given in FK5 coordintates (or
ICRF if using JPLEPH.405) in units of km.
Dimensions must be: array(3,n_elements(date)).
If orbit is not given the correction will be performed
at the center of the earth.
Note: as mentioned by the different authors the use of
the input date to define the orbit (and earth) vector
uses an circular argumentation; but the error due to
the small contribution of this correction is small.
Main steps of the correction are taken from
C. Marquardts description in
http://lheawww.gsfc.nasa.gov/users/craigm/bary/
which gives a good overview of the relevant aspects.
BARYCEN does not perform a dispersion time
correction.
The procedure first reads the jpl ephemeris file
taking the part which is needed for all dates to
convert. Then the coordinates of the earth for all
dates are computed within the FK5 (or ICRF when using
DE405 ephemerides which is mostly the same like FK5)
coordinate system, origin located at the solar system
bary center.
Additionally a unit vector is derived from the RA/DEC
position which is converted also to the FK5
coordinate system (J2000 equinox).
Finally the projection of this vector at the
coordinates of the earth yield the distance the light
must travel; ivided by the speed of light this gives the
time delay we are looking for.
We should mention that JPL DE200 ephemeris are
defined for the FK5 (J2000) coordinate system while
the DE405 is defined for the ICRS coordinate
system. They are defined close together but not
exact.
For the algorithm we always use FK5 coordinates
(i.e. computing all in J2000 and not e.g. in the
Hipparcos IRCS J1991.25 system). The difference lies
in the region of ~50-80mas.
Care must be taken when using the orbit input because
errors in the orbit position directly influence the
correction in a non-neglectable manner.
Internally all computations are performed in units of
km and seconds.
/B1950 - if set, then input coordinates are assumed to be in equinox
B1950 coordinates.
Default is the J2000 equinox.
/TIME_DIFF - if set, then HELIO_JD() returns the time difference
(helio bary centric JD - geocentric JD ) in seconds.
What is the barycentric Julian date of an observation of V402 Cygni
(J2000: RA = 20 9 7.8, Dec = 37 09 07) taken June 15, 1973 at 11:40 UT?
IDL> juldate, [1973,6,15,11,40], jd ; Get geocentric
; reduced Julian date
IDL> bjd = barycen( jd-0.5D0, $ ; convert to bary
ten(20,9,7.8)*15., $ ; center MJD
ten(37,9,7) )
==> bjd = 41848.988072569 (MJD)
Reading of orbit files could be performed with the fits
accessing function READORBIT().
!!! Care must be taken for the units of the result of the
!!! readorbit procedure. Satellite orbit files usually contain
!!! the spacecraft position in meters while here we are
!!! using kilometers!.
jplephread, jplephinterp,
jprecess, tdb2tdt
Processed RXTE lightcurves of GX301 (thank you, Ingo) with
barycen and checked result against fxbary product.
Maximal difference was proven to be in the range of
1.2572855e-06 sec while the averaged difference was
9.5243109e-08 sec. (17-09-2002, E. Goehler)
$Log: barycen.pro,v $
Revision 1.5 2003/07/09 08:34:24 goehler
Major change: date input now MJD instead RJD
Revision 1.4 2003/04/29 09:34:48 goehler
added warning concerning units in km when using the readorbit() function
Revision 1.3 2002/09/17 12:49:17 goehler
Fix of wrong used JPLEPHINTERP function in respect of velocities.
Tests yielded difference to fxbary in the range of less than 1usec.
Revision 1.2 2002/09/13 16:14:43 goehler
preliminary tests give good results (~10msec) for geocentric case. Use of
orbit files still bad.
Revision 1.1 2002/09/12 14:13:58 goehler
first version,
!!! still not tested !!!
[Home Page] [Software, Documentation] [IDL Documentation] [Quick Reference] [Feedback]