1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 3.2//EN">
<!--Converted with LaTeX2HTML 97.1 (release) (July 13th, 1997)
by Nikos Drakos (nikos@cbl.leeds.ac.uk), CBLU, University of Leeds
* revised and updated by: Marcus Hennecke, Ross Moore, Herb Swan
* with significant contributions from:
Jens Lippman, Marek Rouchal, Martin Wilck and others -->
<HTML>
<HEAD>
<TITLE>Ephemerides</TITLE>
<META NAME="description" CONTENT="Ephemerides">
<META NAME="keywords" CONTENT="sun67">
<META NAME="resource-type" CONTENT="document">
<META NAME="distribution" CONTENT="global">
<META HTTP-EQUIV="Content-Type" CONTENT="text/html; charset=iso_8859_1">
<LINK REL="STYLESHEET" HREF="sun67.css">
<LINK REL="next" HREF="node225.html">
<LINK REL="previous" HREF="node223.html">
<LINK REL="up" HREF="node197.html">
<LINK REL="next" HREF="node225.html">
</HEAD>
<BODY >
<BR> <HR>
<A NAME="tex2html2693" HREF="node225.html">
<IMG WIDTH="37" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="next" SRC="next_motif.gif"></A>
<A NAME="tex2html2691" HREF="node197.html">
<IMG WIDTH="26" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="up" SRC="up_motif.gif"></A>
<A NAME="tex2html2685" HREF="node223.html">
<IMG WIDTH="63" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="previous" SRC="previous_motif.gif"></A> <A HREF="sun67.html#stardoccontents"><IMG ALIGN="BOTTOM" BORDER="0"
SRC="contents_motif.gif"></A>
<BR>
<B> Next:</B> <A NAME="tex2html2694" HREF="node225.html">Radial Velocity and Light-Time Corrections</A>
<BR>
<B>Up:</B> <A NAME="tex2html2692" HREF="node197.html">EXPLANATION AND EXAMPLES</A>
<BR>
<B> Previous:</B> <A NAME="tex2html2686" HREF="node223.html">Geocentric Coordinates</A>
<BR> <HR> <P>
<P><!--End of Navigation Panel-->
<H2><A NAME="SECTION000518000000000000000">
Ephemerides</A>
</H2>
SLALIB includes routines for generating positions and
velocities of Solar-System bodies. The accuracy objectives are
modest, and the SLALIB facilities do not attempt
to compete with precomputed ephemerides such as
those provided by JPL, or with models containing
thousands of terms. It is also worth noting
that SLALIB's very accurate star coordinate conversion
routines are not strictly applicable to solar-system cases,
though they are adequate for most practical purposes.
<P>
Earth/Sun ephemerides can be generated using the routine
sla_EVP,
which predicts Earth position and velocity with respect to both the
solar-system barycentre and the
Sun. Maximum velocity error is 0.42 metres per second; maximum
heliocentric position error is 1600 km (about <IMG WIDTH="17" HEIGHT="17" ALIGN="BOTTOM" BORDER="0"
SRC="img316.gif"
ALT="$2\hspace{-0.05em}^{'\hspace{-0.1em}'}$">), with
barycentric position errors about 4 times worse.
(The Sun's position as
seen from the Earth can, of course, be obtained simply by
reversing the signs of the Cartesian components of the
Earth:Sun vector.)
<P>
Geocentric Moon ephemerides are available from
sla_DMOON,
which predicts the Moon's position and velocity with respect to
the Earth's centre. Direction accuracy is usually better than
10 km (<IMG WIDTH="17" HEIGHT="18" ALIGN="BOTTOM" BORDER="0"
SRC="img131.gif"
ALT="$5\hspace{-0.05em}^{'\hspace{-0.1em}'}$">) and distance accuracy a little worse.
<P>
Lower-precision but faster predictions for the Sun and Moon
can be made by calling
sla_EARTH
and
sla_MOON.
Both are single precision and accept dates in the form of
year, day-in-year and fraction of day
(starting from a calendar date you need to call
sla_CLYD
or
sla_CALYD
to get the required year and day).
The
sla_EARTH
routine returns the heliocentric position and velocity
of the Earth's centre for the mean equator and
equinox of date. The accuracy is better than 20,000 km in position
and 10 metres per second in speed.
The
position and velocity of the Moon with respect to the
Earth's centre for the mean equator and ecliptic of date
can be obtained by calling
sla_MOON.
The positional accuracy is better than <IMG WIDTH="25" HEIGHT="18" ALIGN="BOTTOM" BORDER="0"
SRC="img82.gif"
ALT="$30\hspace{-0.05em}^{'\hspace{-0.1em}'}$"> in direction
and 1000 km in distance.
<P>
Approximate ephemerides for all the major planets
can be generated by calling
sla_PLANET
or
sla_RDPLAN. These routines offer arcminute accuracy (much
better for the inner planets and for Pluto) over a span of several
millennia (but only <IMG WIDTH="39" HEIGHT="25" ALIGN="MIDDLE" BORDER="0"
SRC="img317.gif"
ALT="$\pm100$"> years for Pluto).
The routine
sla_PLANET produces heliocentric position and
velocity in the form of equatorial <IMG WIDTH="106" HEIGHT="29" ALIGN="MIDDLE" BORDER="0"
SRC="img51.gif"
ALT="$[\,x,y,z,\dot{x},\dot{y},\dot{z}\,]$"> for the
mean equator and equinox of J2000. The vectors
produced by
sla_PLANET
can be used in a variety of ways according to the
requirements of the application concerned. The routine
sla_RDPLAN
uses
sla_PLANET
and
sla_DMOON
to deal with the common case of predicting
a planet's apparent <IMG WIDTH="42" HEIGHT="29" ALIGN="MIDDLE" BORDER="0"
SRC="img3.gif"
ALT="$[\,\alpha,\delta\,]$"> and angular size as seen by a
terrestrial observer.
<P>
Note that in predicting the position in the sky of a solar-system body
it is necessary to allow for geocentric parallax. This correction
is <I>essential</I> in the case of the Moon, where the observer's
position on the Earth can affect the Moon's <IMG WIDTH="42" HEIGHT="29" ALIGN="MIDDLE" BORDER="0"
SRC="img3.gif"
ALT="$[\,\alpha,\delta\,]$"> by up to
<IMG WIDTH="18" HEIGHT="13" ALIGN="BOTTOM" BORDER="0"
SRC="img318.gif"
ALT="$1^\circ$">. The calculation can most conveniently be done by calling
sla_PVOBS and subtracting the resulting 6-vector from the
one produced by
sla_DMOON, as is demonstrated by the following example:
<P><PRE>
* Demonstrate the size of the geocentric parallax correction
* in the case of the Moon. The test example is for the AAT,
* before midnight, in summer, near first quarter.
IMPLICIT NONE
CHARACTER NAME*40,SH,SD
INTEGER J,I,IHMSF(4),IDMSF(4)
DOUBLE PRECISION SLONGW,SLAT,H,DJUTC,FDUTC,DJUT1,DJTT,STL,
: RMATN(3,3),PMM(6),PMT(6),RM,DM,PVO(6),TL
DOUBLE PRECISION sla_DTT,sla_GMST,sla_EQEQX,sla_DRANRM
* Get AAT longitude and latitude in radians and height in metres
CALL sla_OBS(0,'AAT',NAME,SLONGW,SLAT,H)
* UTC (1992 January 13, 11 13 59) to MJD
CALL sla_CLDJ(1992,1,13,DJUTC,J)
CALL sla_DTF2D(11,13,59.0D0,FDUTC,J)
DJUTC=DJUTC+FDUTC
* UT1 (UT1-UTC value of -0.152 sec is from IERS Bulletin B)
DJUT1=DJUTC+(-0.152D0)/86400D0
* TT
DJTT=DJUTC+sla_DTT(DJUTC)/86400D0
* Local apparent sidereal time
STL=sla_GMST(DJUT1)-SLONGW+sla_EQEQX(DJTT)
* Geocentric position/velocity of Moon (mean of date)
CALL sla_DMOON(DJTT,PMM)
* Nutation to true equinox of date
CALL sla_NUT(DJTT,RMATN)
CALL sla_DMXV(RMATN,PMM,PMT)
CALL sla_DMXV(RMATN,PMM(4),PMT(4))
* Report geocentric HA,Dec
CALL sla_DCC2S(PMT,RM,DM)
CALL sla_DR2TF(2,sla_DRANRM(STL-RM),SH,IHMSF)
CALL sla_DR2AF(1,DM,SD,IDMSF)
WRITE (*,'(1X,'' geocentric:'',2X,A,I2.2,2I3.2,''.'',I2.2,'//
: '1X,A,I2.2,2I3.2,''.'',I1)')
: SH,IHMSF,SD,IDMSF
* Geocentric position of observer (true equator and equinox of date)
CALL sla_PVOBS(SLAT,H,STL,PVO)
* Place origin at observer
DO I=1,6
PMT(I)=PMT(I)-PVO(I)
END DO
* Allow for planetary aberration
TL=499.004782D0*SQRT(PMT(1)**2+PMT(2)**2+PMT(3)**2)
DO I=1,3
PMT(I)=PMT(I)-TL*PMT(I+3)
END DO
* Report topocentric HA,Dec
CALL sla_DCC2S(PMT,RM,DM)
CALL sla_DR2TF(2,sla_DRANRM(STL-RM),SH,IHMSF)
CALL sla_DR2AF(1,DM,SD,IDMSF)
WRITE (*,'(1X,''topocentric:'',2X,A,I2.2,2I3.2,''.'',I2.2,'//
: '1X,A,I2.2,2I3.2,''.'',I1)')
: SH,IHMSF,SD,IDMSF
END
</PRE>
<P>
The output produced is as follows:
<P><PRE>
geocentric: +03 06 55.59 +15 03 39.0
topocentric: +03 09 23.79 +15 40 51.5
</PRE>
<P>(An easier but
less instructive method of estimating the topocentric apparent place of the
Moon is to call the routine
sla_RDPLAN.)
<P>
As an example of using
sla_PLANET,
the following program estimates the geocentric separation
between Venus and Jupiter during a close conjunction
in 2BC, which is a star-of-Bethlehem candidate:
<P><PRE>
* Compute time and minimum geocentric apparent separation
* between Venus and Jupiter during the close conjunction of 2 BC.
IMPLICIT NONE
DOUBLE PRECISION SEPMIN,DJD0,FD,DJD,DJDM,DF,PV(6),RMATP(3,3),
: PVM(6),PVE(6),TL,RV,DV,RJ,DJ,SEP
INTEGER IHOUR,IMIN,J,I,IHMIN,IMMIN
DOUBLE PRECISION sla_EPJ,sla_DSEP
* Search for closest approach on the given day
DJD0=1720859.5D0
SEPMIN=1D10
DO IHOUR=20,22
DO IMIN=0,59
CALL sla_DTF2D(IHOUR,IMIN,0D0,FD,J)
* Julian date and MJD
DJD=DJD0+FD
DJDM=DJD-2400000.5D0
* Earth to Moon (mean of date)
CALL sla_DMOON(DJDM,PV)
* Precess Moon position to J2000
CALL sla_PRECL(sla_EPJ(DJDM),2000D0,RMATP)
CALL sla_DMXV(RMATP,PV,PVM)
* Sun to Earth-Moon Barycentre (mean J2000)
CALL sla_PLANET(DJDM,3,PVE,J)
* Correct from EMB to Earth
DO I=1,3
PV(I)=PVE(I)-0.012150581D0*PVM(I)
END DO
* Sun to Venus
CALL sla_PLANET(DJDM,2,PV,J)
* Earth to Venus
DO I=1,6
PV(I)=PV(I)-PVE(I)
END DO
* Light time to Venus (sec)
TL=499.004782D0*SQRT((PV(1)-PVE(1))**2+
: (PV(2)-PVE(2))**2+
: (PV(3)-PVE(3))**2)
* Extrapolate backwards in time by that much
DO I=1,3
PV(I)=PV(I)-TL*PV(I+3)
END DO
* To RA,Dec
CALL sla_DCC2S(PV,RV,DV)
* Same for Jupiter
CALL sla_PLANET(DJDM,5,PV,J)
DO I=1,6
PV(I)=PV(I)-PVE(I)
END DO
TL=499.004782D0*SQRT((PV(1)-PVE(1))**2+
: (PV(2)-PVE(2))**2+
: (PV(3)-PVE(3))**2)
DO I=1,3
PV(I)=PV(I)-TL*PV(I+3)
END DO
CALL sla_DCC2S(PV,RJ,DJ)
* Separation (arcsec)
SEP=sla_DSEP(RV,DV,RJ,DJ)
* Keep if smallest so far
IF (SEP.LT.SEPMIN) THEN
IHMIN=IHOUR
IMMIN=IMIN
SEPMIN=SEP
END IF
END DO
END DO
* Report
WRITE (*,'(1X,I2.2,'':'',I2.2,F6.1)') IHMIN,IMMIN,
: 206264.8062D0*SEPMIN
END
</PRE>
<P>
The output produced (the Ephemeris Time on the day in question, and
the closest approach in arcseconds) is as follows:
<P><PRE>
21:19 33.7
</PRE>
<P>
For comparison, accurate predictions based on the JPL DE102 ephemeris
give a separation about <IMG WIDTH="17" HEIGHT="18" ALIGN="BOTTOM" BORDER="0"
SRC="img319.gif"
ALT="$8\hspace{-0.05em}^{'\hspace{-0.1em}'}$"> less than
the above estimate, occurring about half an hour earlier
(see <I>Sky and Telescope,</I> April 1987, p357).
<P>
The following program demonstrates
sla_RDPLAN.
<PRE>
* For a given date, time and geographical location, output
* a table of planetary positions and diameters.
IMPLICIT NONE
CHARACTER PNAMES(0:9)*7,B*80,S
INTEGER I,NP,IY,J,IM,ID,IHMSF(4),IDMSF(4)
DOUBLE PRECISION R2AS,FD,DJM,ELONG,PHI,RA,DEC,DIAM
PARAMETER (R2AS=206264.80625D0)
DATA PNAMES / 'Sun','Mercury','Venus','Moon','Mars','Jupiter',
: 'Saturn','Uranus','Neptune', 'Pluto' /
* Loop until 'end' typed
B=' '
DO WHILE (B.NE.'END'.AND.B.NE.'end')
* Get date, time and observer's location
PRINT *,'Date? (Y,M,D, Gregorian)'
READ (*,'(A)') B
IF (B.NE.'END'.AND.B.NE.'end') THEN
I=1
CALL sla_INTIN(B,I,IY,J)
CALL sla_INTIN(B,I,IM,J)
CALL sla_INTIN(B,I,ID,J)
PRINT *,'Time? (H,M,S, dynamical)'
READ (*,'(A)') B
I=1
CALL sla_DAFIN(B,I,FD,J)
FD=FD*2.3873241463784300365D0
CALL sla_CLDJ(IY,IM,ID,DJM,J)
DJM=DJM+FD
PRINT *,'Longitude? (D,M,S, east +ve)'
READ (*,'(A)') B
I=1
CALL sla_DAFIN(B,I,ELONG,J)
PRINT *,'Latitude? (D,M,S, (geodetic)'
READ (*,'(A)') B
I=1
CALL sla_DAFIN(B,I,PHI,J)
* Loop planet by planet
DO NP=0,8
* Get RA,Dec and diameter
CALL sla_RDPLAN(DJM,NP,ELONG,PHI,RA,DEC,DIAM)
* One line of report
CALL sla_DR2TF(2,RA,S,IHMSF)
CALL sla_DR2AF(1,DEC,S,IDMSF)
WRITE (*,
: '(1X,A,2X,3I3.2,''.'',I2.2,2X,A,I2.2,2I3.2,''.'',I1,F8.1)')
: PNAMES(NP),IHMSF,S,IDMSF,R2AS*DIAM
* Next planet
END DO
PRINT *,' '
END IF
* Next case
END DO
END
</PRE>
Entering the following data (for 1927 June 29 at <IMG WIDTH="49" HEIGHT="17" ALIGN="BOTTOM" BORDER="0"
SRC="img320.gif"
ALT="$5^{\rm h}\,25^{\rm m}$"> ET
and the position of Preston, UK.):
<PRE>
1927 6 29
5 25
-2 42
53 46
</PRE>
produces the following report:
<PRE>
Sun 06 28 14.03 +23 17 17.5 1887.8
Mercury 08 08 58.62 +19 20 57.3 9.3
Venus 09 38 53.64 +15 35 32.9 22.8
Moon 06 28 18.30 +23 18 37.3 1903.9
Mars 09 06 49.34 +17 52 26.7 4.0
Jupiter 00 11 12.06 -00 10 57.5 41.1
Saturn 16 01 43.34 -18 36 55.9 18.2
Uranus 00 13 33.53 +00 39 36.0 3.5
Neptune 09 49 35.75 +13 38 40.8 2.2
Pluto 07 05 29.50 +21 25 04.2 .1
</PRE>
Inspection of the Sun and Moon data reveals that
a total solar eclipse is in progress.
<P>
SLALIB also provides for the case where orbital elements (with respect
to the J2000 equinox and ecliptic)
are available. This allows predictions to be made for minor-planets and
(if you ignore non-gravitational effects)
comets. Furthermore, if major-planet elements for an epoch close to the date
in question are available, more accurate predictions can be made than
are offered by
sla_RDPLAN and
sla_PLANET.
<P>
The SLALIB planetary-prediction
routines that work with orbital elements are
sla_PLANTE (the orbital-elements equivalent of
sla_RDPLAN), which predicts the topocentric <IMG WIDTH="42" HEIGHT="29" ALIGN="MIDDLE" BORDER="0"
SRC="img3.gif"
ALT="$[\,\alpha,\delta\,]$">, and
sla_PLANEL (the orbital-elements equivalent of
sla_PLANET), which predicts the heliocentric <IMG WIDTH="106" HEIGHT="29" ALIGN="MIDDLE" BORDER="0"
SRC="img51.gif"
ALT="$[\,x,y,z,\dot{x},\dot{y},\dot{z}\,]$"> with respect to the
J2000 equinox and equator. In addition, the routine
sla_PV2EL does the inverse of
sla_PLANEL, transforming <IMG WIDTH="106" HEIGHT="29" ALIGN="MIDDLE" BORDER="0"
SRC="img51.gif"
ALT="$[\,x,y,z,\dot{x},\dot{y},\dot{z}\,]$"> into <I>osculating elements.</I>
<P>
Osculating elements describe the unperturbed 2-body orbit. This is
a good approximation to the actual orbit for a few weeks either
side of the specified epoch, outside which perturbations due to
the other bodies of the Solar System lead to
increasing errors. Given a minor planet's osculating elements for
a particular date, predictions for a date even just
100 days earlier or later
are likely to be in error by several arcseconds.
These errors can
be reduced if new elements are generated which take account of
the perturbations of the major planets, and this is what the routine
sla_PERTEL does. Once
sla_PERTEL has been called, to provide osculating elements
close to the required date, the elements can be passed to
sla_PLANEL or
sla_PLANTE in the normal way. Predictions of arcsecond accuracy
over a span of a decade or more are available using this
technique.
<P>
Three different combinations of orbital elements are
provided for, matching the usual conventions
for major planets, minor planets and
comets respectively. The choice is made through the
argument <TT>JFORM</TT>:
<BR>
<P><TABLE CELLPADDING=3 BORDER="1">
<TR VALIGN="TOP"><TD ALIGN="CENTER" NOWRAP><TT>JFORM=1</TT></TD>
<TD ALIGN="CENTER" NOWRAP><TT>JFORM=2</TT></TD>
<TD ALIGN="CENTER" NOWRAP><TT>JFORM=3</TT></TD>
</TR>
<TR VALIGN="TOP"><TD ALIGN="CENTER" NOWRAP><I>t<SUB>0</SUB></I></TD>
<TD ALIGN="CENTER" NOWRAP><I>t<SUB>0</SUB></I></TD>
<TD ALIGN="CENTER" NOWRAP><I>T</I></TD>
</TR>
<TR VALIGN="TOP"><TD ALIGN="CENTER" NOWRAP><I>i</I></TD>
<TD ALIGN="CENTER" NOWRAP><I>i</I></TD>
<TD ALIGN="CENTER" NOWRAP><I>i</I></TD>
</TR>
<TR VALIGN="TOP"><TD ALIGN="CENTER" NOWRAP><IMG WIDTH="14" HEIGHT="13" ALIGN="BOTTOM" BORDER="0"
SRC="img99.gif"
ALT="$\Omega$"></TD>
<TD ALIGN="CENTER" NOWRAP><IMG WIDTH="14" HEIGHT="13" ALIGN="BOTTOM" BORDER="0"
SRC="img99.gif"
ALT="$\Omega$"></TD>
<TD ALIGN="CENTER" NOWRAP><IMG WIDTH="14" HEIGHT="13" ALIGN="BOTTOM" BORDER="0"
SRC="img99.gif"
ALT="$\Omega$"></TD>
</TR>
<TR VALIGN="TOP"><TD ALIGN="CENTER" NOWRAP><IMG WIDTH="16" HEIGHT="14" ALIGN="BOTTOM" BORDER="0"
SRC="img100.gif"
ALT="$\varpi$"></TD>
<TD ALIGN="CENTER" NOWRAP><IMG WIDTH="13" HEIGHT="14" ALIGN="BOTTOM" BORDER="0"
SRC="img101.gif"
ALT="$\omega$"></TD>
<TD ALIGN="CENTER" NOWRAP><IMG WIDTH="13" HEIGHT="14" ALIGN="BOTTOM" BORDER="0"
SRC="img101.gif"
ALT="$\omega$"></TD>
</TR>
<TR VALIGN="TOP"><TD ALIGN="CENTER" NOWRAP><I>a</I></TD>
<TD ALIGN="CENTER" NOWRAP><I>a</I></TD>
<TD ALIGN="CENTER" NOWRAP><I>q</I></TD>
</TR>
<TR VALIGN="TOP"><TD ALIGN="CENTER" NOWRAP><I>e</I></TD>
<TD ALIGN="CENTER" NOWRAP><I>e</I></TD>
<TD ALIGN="CENTER" NOWRAP><I>e</I></TD>
</TR>
<TR VALIGN="TOP"><TD ALIGN="CENTER" NOWRAP><I>L</I></TD>
<TD ALIGN="CENTER" NOWRAP><I>M</I></TD>
<TD ALIGN="CENTER" NOWRAP> </TD>
</TR>
<TR VALIGN="TOP"><TD ALIGN="CENTER" NOWRAP><I>n</I></TD>
<TD ALIGN="CENTER" NOWRAP> </TD>
<TD ALIGN="CENTER" NOWRAP> </TD>
</TR>
</TABLE>
<BR>
<BR>
<BR>
<BR>
<BR>
<BR>
The symbols have the following meanings:
<PRE><TT>
<I>t<SUB>0</SUB></I> epoch at which the elements were correct
<I>T</I> epoch of perihelion passage
<I>i</I> inclination of the orbit
<IMG WIDTH="14" HEIGHT="13" ALIGN="BOTTOM" BORDER="0"
SRC="img99.gif"
ALT="$\Omega$"> longitude of the ascending node
<IMG WIDTH="16" HEIGHT="14" ALIGN="BOTTOM" BORDER="0"
SRC="img100.gif"
ALT="$\varpi$"> longitude of perihelion (<IMG WIDTH="81" HEIGHT="27" ALIGN="MIDDLE" BORDER="0"
SRC="img321.gif"
ALT="$\varpi = \Omega + \omega$">) <IMG WIDTH="13" HEIGHT="14" ALIGN="BOTTOM" BORDER="0"
SRC="img101.gif"
ALT="$\omega$"> argument of perihelion
<I>a</I> semi-major axis of the orbital ellipse
<I>q</I> perihelion distance
<I>e</I> orbital eccentricity
<I>L</I> mean longitude (<IMG WIDTH="87" HEIGHT="27" ALIGN="MIDDLE" BORDER="0"
SRC="img322.gif"
ALT="$L = \varpi + M$">) <I>M</I> mean anomaly
<I>n</I> mean motion
</TT></PRE>
<P>
The mean motion, <I>n</I>, tells sla_PLANEL the mass of the planet.
If it is not available, it should be claculated
from <I>n<SUP>2</SUP></I> <I>a<SUP>3</SUP></I> = <I>k<SUP>2</SUP></I> (1+<I>m</I>), where <I>k</I> = 0.01720209895 and
m is the mass of the planet (<IMG WIDTH="59" HEIGHT="27" ALIGN="MIDDLE" BORDER="0"
SRC="img323.gif"
ALT="$M_\odot = 1$">); <I>a</I> is in AU.
<P>
Conventional elements are not the only way of specifying an orbit.
The <IMG WIDTH="106" HEIGHT="29" ALIGN="MIDDLE" BORDER="0"
SRC="img51.gif"
ALT="$[\,x,y,z,\dot{x},\dot{y},\dot{z}\,]$"> state vector is an equally valid specification,
and the so-called <I>method of universal variables</I> allows
orbital calculations to be made directly, bypassing angular
quantities and avoiding Kepler's Equation. The universal-variables
approach has various advantages, including better handling of
near-parabolic cases and greater efficiency.
SLALIB uses universal variables for its internal
calculations and also offers a number of routines which
applications can call.
<P>
The universal elements are the <IMG WIDTH="106" HEIGHT="29" ALIGN="MIDDLE" BORDER="0"
SRC="img51.gif"
ALT="$[\,x,y,z,\dot{x},\dot{y},\dot{z}\,]$"> and its epoch, plus the mass
of the body. The SLALIB routines supplement these elements with
certain redundant values in order to
avoid unnecessary recomputation when the elements are next used.
<P>
The routines
sla_EL2UE and
sla_UE2EL transform conventional elements into the
universal form and <I>vice versa.</I>
The routine
sla_PV2UE takes an <IMG WIDTH="106" HEIGHT="29" ALIGN="MIDDLE" BORDER="0"
SRC="img51.gif"
ALT="$[\,x,y,z,\dot{x},\dot{y},\dot{z}\,]$"> and forms the set of universal
elements;
sla_UE2PV takes a set of universal elements and predicts the <IMG WIDTH="106" HEIGHT="29" ALIGN="MIDDLE" BORDER="0"
SRC="img51.gif"
ALT="$[\,x,y,z,\dot{x},\dot{y},\dot{z}\,]$"> for a specified epoch.
The routine
sla_PERTUE provides updated universal elements,
taking into account perturbations from the major planets.
<P>
<BR> <HR>
<A NAME="tex2html2693" HREF="node225.html">
<IMG WIDTH="37" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="next" SRC="next_motif.gif"></A>
<A NAME="tex2html2691" HREF="node197.html">
<IMG WIDTH="26" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="up" SRC="up_motif.gif"></A>
<A NAME="tex2html2685" HREF="node223.html">
<IMG WIDTH="63" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="previous" SRC="previous_motif.gif"></A> <A HREF="sun67.html#stardoccontents"><IMG ALIGN="BOTTOM" BORDER="0"
SRC="contents_motif.gif"></A>
<BR>
<B> Next:</B> <A NAME="tex2html2694" HREF="node225.html">Radial Velocity and Light-Time Corrections</A>
<BR>
<B>Up:</B> <A NAME="tex2html2692" HREF="node197.html">EXPLANATION AND EXAMPLES</A>
<BR>
<B> Previous:</B> <A NAME="tex2html2686" HREF="node223.html">Geocentric Coordinates</A>
<BR> <HR> <P>
<P><!--End of Navigation Panel-->
<ADDRESS>
<I>SLALIB --- Positional Astronomy Library<BR>Starlink User Note 67<BR>P. T. Wallace<BR>12 October 1999<BR>E-mail:ptw@star.rl.ac.uk</I>
</ADDRESS>
</BODY>
</HTML>
|