source: trunk/MagicSoft/slalib/rdplan.c@ 10146

Last change on this file since 10146 was 731, checked in by tbretz, 24 years ago
*** empty log message ***
  • Property svn:executable set to *
File size: 5.1 KB
Line 
1#include "slalib.h"
2#include "slamac.h"
3void slaRdplan ( double date, int np, double elong, double phi,
4 double *ra, double *dec, double *diam )
5/*
6** - - - - - - - - - -
7** s l a R d p l a n
8** - - - - - - - - - -
9**
10** Approximate topocentric apparent RA,Dec of a planet, and its
11** angular diameter.
12**
13** Given:
14** date double MJD of observation (JD - 2400000.5)
15** np int planet: 1 = Mercury
16** 2 = Venus
17** 3 = Moon
18** 4 = Mars
19** 5 = Jupiter
20** 6 = Saturn
21** 7 = Uranus
22** 8 = Neptune
23** 9 = Pluto
24** else = Sun
25** elong,phi double observer's east longitude and geodetic
26** latitude (radians)
27**
28** Returned:
29** ra,dec double RA, Dec (topocentric apparent, radians)
30** diam double angular diameter (equatorial, radians)
31**
32** Notes:
33**
34** 1 The date is in a dynamical timescale (TDB, formerly ET) and is
35** in the form of a Modified Julian Date (JD-2400000.5). For all
36** practical purposes, TT can be used instead of TDB, and for many
37** applications UT will do (except for the Moon).
38**
39** 2 The longitude and latitude allow correction for geocentric
40** parallax. This is a major effect for the Moon, but in the
41** context of the limited accuracy of the present routine its
42** effect on planetary positions is small (negligible for the
43** outer planets). Geocentric positions can be generated by
44** appropriate use of the routines slaDmoon and slaPlanet.
45**
46** 3 The direction accuracy (arcsec, 1000-3000AD) is of order:
47**
48** Sun 5
49** Mercury 2
50** Venus 10
51** Moon 30
52** Mars 50
53** Jupiter 90
54** Saturn 90
55** Uranus 90
56** Neptune 10
57** Pluto 1 (1885-2099AD only)
58**
59** The angular diameter accuracy is about 0.4% for the Moon,
60** and 0.01% or better for the Sun and planets.
61**
62** Called: slaGmst, slaDt, slaEpj, slaDmoon, slaPvobs, slaPrenut,
63** slaPlanet, slaDmxv, slaDcc2s, slaDranrm
64**
65** Last revision: 27 May 1997
66**
67** Copyright P.T.Wallace. All rights reserved.
68*/
69
70#define AUKM 1.49597870e8 /* AU in km */
71#define TAU 499.004782 /* Light time for unit distance (sec) */
72
73{
74 int ip, j, i;
75 double stl, vgm[6], v[6], rmat[3][3], vse[6], vsg[6], vsp[6],
76 vgo[6], dx, dy, dz, r, tl;
77
78/* Equatorial radii (km) */
79 static double eqrau[] = {
80 696000.0, /* Sun */
81 2439.7, /* Mercury */
82 6051.9, /* Venus */
83 1738.0, /* Moon */
84 3397.0, /* Mars */
85 71492.0, /* Jupiter */
86 60268.0, /* Saturn */
87 25559.0, /* Uranus */
88 24764.0, /* Neptune */
89 1151.0 /* Pluto */
90 };
91
92
93
94/* Classify NP. */
95 ip = ( np >= 1 && np <= 9 ) ? np : 0;
96
97/* Approximate local ST. */
98 stl = slaGmst ( date - slaDt ( slaEpj ( date ) ) / 86400.0 ) + elong;
99
100/* Geocentre to Moon (mean of date). */
101 slaDmoon ( date, v );
102
103/* Nutation, to true of date. */
104 slaNut ( date, rmat );
105 slaDmxv ( rmat, &v[0], &vgm[0] );
106 slaDmxv ( rmat, &v[3], &vgm[3] );
107
108/* Moon? */
109 if ( ip == 3 ) {
110
111 /* Yes: geocentre to Moon (true of date). */
112 for ( i = 0; i <= 5; i++ ) v[i] = vgm[i];
113
114 } else {
115
116 /* No: precession/nutation matrix, J2000 to date. */
117 slaPrenut ( 2000.0, date, rmat );
118
119 /* Sun to Earth-Moon Barycentre (J2000). */
120 slaPlanet ( date, 3, v, &j );
121
122 /* Precession and nutation to date. */
123 slaDmxv ( rmat, &v[0], &vse[0] );
124 slaDmxv ( rmat, &v[3], &vse[3] );
125
126 /* Sun to geocentre. */
127 for ( i = 0; i <= 5; i++ ) vsg[i] = vse[i] - 0.012150581 * vgm[i];
128
129 /* Sun? */
130 if ( ip == 0 ) {
131
132 /* Yes: geocentre to Sun. */
133 for ( i = 0; i <= 5; i++ ) v[i] = - vsg[i];
134
135 } else {
136
137 /* No: Sun to Planet. */
138 slaPlanet ( date, ip, v, &j );
139
140 /* Precession and nutation to date. */
141 slaDmxv ( rmat, &v[0], &vsp[0] );
142 slaDmxv ( rmat, &v[3], &vsp[3] );
143
144 /* Geocentre to planet. */
145 for ( i = 0; i <= 5; i++ ) v[i] = vsp[i] - vsg[i];
146 }
147 }
148
149/* Refer to origin at the observer. */
150 slaPvobs ( phi, 0.0, stl, vgo );
151 for ( i = 0; i <= 5; i++ ) v[i] -= vgo[i];
152
153/* Geometric distance (AU). */
154 dx = v[0];
155 dy = v[1];
156 dz = v[2];
157 r = sqrt ( dx * dx + dy * dy + dz * dz );
158
159/* Light time (sec). */
160 tl = TAU * r;
161
162/* Correct position for planetary aberration. */
163 for ( i = 0; i <= 2; i++ ) v[i] -= tl * v[i+3];
164
165/* To RA,Dec. */
166 slaDcc2s ( v, ra, dec );
167 *ra = slaDranrm ( *ra );
168
169/* Angular diameter (radians). */
170 *diam = 2.0 * asin ( eqrau[ip] / ( r * AUKM ) );
171}
Note: See TracBrowser for help on using the repository browser.