source: Sophya/trunk/SophyaExt/XephemAstroLib/twobody.c@ 2937

Last change on this file since 2937 was 2818, checked in by cmv, 20 years ago

Update de Xephem 3.7 cmv 21/08/2005

File size: 6.5 KB
Line 
1/*
2 *
3 * TWOBODY.C
4 *
5 * Computation of planetary position, two-body computation
6 *
7 * Paul Schlyter, 1987-06-15
8 *
9 * Decreased EPSILON from 2E-4 to 3E-8, 1988-12-05
10 *
11 * 1990-01-01: Bug fix in almost parabolic orbits: now the routine
12 * doesn't bomb there (an if block was too large)
13 *
14 * 2000-12-06: Donated to Elwood Downey if he wants to use it in XEphem
15 */
16
17
18#include <stdio.h>
19#include <stdlib.h>
20#include <math.h>
21
22
23/* Constants used when solving Kepler's equation */
24#undef EPSILON
25#define EPSILON 3E-8
26#undef INFINITY
27#define INFINITY 1E+10
28
29/* Math constants */
30#undef PI
31#define PI 3.14159265358979323846
32#define RADEG ( 180.0 / PI )
33#define DEGRAD ( PI / 180.0 )
34
35/* Trig functions in degrees */
36#define sind(x) sin(x*DEGRAD)
37#define cosd(x) cos(x*DEGRAD)
38#define atand(x) (RADEG*atan(x))
39#define atan2d(y,x) (RADEG*atan2(y,x))
40
41/* Gauss' grav.-konstant */
42#define K 1.720209895E-2
43#define KD ( K * 180.0 / PI )
44#define K2 ( K / 2.0 )
45
46
47
48
49static double cubroot( double x )
50/* Cubic root */
51{
52 double a,b;
53
54 if ( x == 0.0 )
55 return 0.0;
56 else
57 {
58 a = fabs(x);
59 b = exp( log(a) / 3.0 );
60 return x > 0.0 ? b : -b;
61 }
62} /* cubroot */
63
64
65static double rev180( double ang )
66/* Normalize angle to between +180 and -180 degrees */
67{
68 return ang - 360.0 * floor(ang*(1.0/360.0) + 0.5);
69} /* rev180 */
70
71
72
73static double kepler( double m, double ex )
74/*
75 * Solves Kepler's equation
76 * m = mean anomaly
77 * ex = eccentricity
78 * kepler = eccentric anomaly
79 */
80{
81 double m1, sinm, cosm, exd, exan, dexan, lim1, adko, adk, denom;
82 int converged;
83
84 m1 = rev180(m);
85 sinm = sind(m1);
86 cosm = cosd(m1);
87 /* 1st approximation: */
88 exan = atan2d(sinm,cosm-ex);
89 if ( ex > 0.008 )
90 { /* Iteration formula: */
91 exd = ex * RADEG;
92 lim1 = 1E-3 / ex;
93 adko = INFINITY;
94 denom = 1.0 - ex * cosd(exan);
95 do
96 {
97 dexan = (m1 + exd * sind(exan) - exan) / denom;
98 exan = exan + dexan;
99 adk = fabs(dexan);
100 converged = adk < EPSILON || adk >= adko ;
101 adko = adk;
102 if ( !converged && adk > lim1 )
103 denom = 1.0 - ex * cosd(exan);
104 } while ( !converged );
105 }
106 return exan;
107} /* kepler */
108
109
110static void vr( double *v, double *r, double m, double e, double a )
111/*
112 * Elliptic orbits only:
113 * computes: v = true anomaly (degrees)
114 * r = radius vector (a.u.)
115 * from: m = mean anomaly (degrees)
116 * e = eccentricity
117 * a = semimajor axis (a.u.)
118 */
119{
120 double ean, x, y;
121
122 ean = kepler(m,e);
123 x = a*(cosd(ean)-e);
124 y = a*sqrt(1.-e*e)*sind(ean);
125 *r = sqrt(x*x+y*y);
126 *v = atan2d(y,x);
127} /* vr */
128
129
130/* return 0 if ok, else -1 */
131int vrc( double *v, double *r, double tp, double e, double q )
132/*
133 * Elliptic, hyperbolic and near-parabolic orbits:
134 * computes: v = true anomaly (degrees)
135 * r = radius vector (a.u.)
136 * from: tp = time from perihelion (days)
137 * e = eccentricity
138 * q = perihelion distance (a.u.)
139 */
140{
141
142 double lambda;
143
144 double a, b, w, w2, w4, c, c1, c2, c3, c5, a0, a1, a2,
145 a3, m, n, g, adgg, adgg2, gs, dg;
146
147 if ( tp == 0.0 ) /* In perihelion */
148 {
149 *v = 0.0;
150 *r = q;
151 return 0;
152 }
153
154
155 lambda = (1.0-e) / (1.0+e);
156
157 if ( fabs(lambda) < 0.01 )
158 { /* Near-parabolic orbits */
159 a = K2 * sqrt((1.0+e)/(q*q*q)) * tp;
160 b = sqrt( 1.0 + 2.25*a*a );
161 w = cubroot( b + 1.5*a ) - cubroot( b - 1.5*a );
162
163 /* Test if it's accuate enough to compute this as a near-parabolic orbit */
164 if ( fabs(w*w*lambda) > 0.2 )
165 {
166 if ( fabs(lambda) < 0.0002 )
167 {
168 /* Sorry, but we cannot compute this at all -- we must give up!
169 *
170 * This happens very rarely, in orbits having an eccentricity
171 * some 2% away from 1.0 AND if the body is very very far from
172 * perihelion. E.g. a Kreutz sun-grazing comet having
173 * eccentricity near 0.98 or 1.02, and being outside
174 * the orbit of Pluto. For any reasonable orbit this will
175 * never happen in practice.
176 *
177 * You might want to code a more graceful error exit here though.
178 *
179 */
180 printf( "\nNear-parabolic orbit: inaccurate result."
181 "\n e = %f, lambda = %f, w = %f", e, lambda, w );
182 return -1;
183 }
184 else
185 {
186 /* We cannot compute this as a near-parabolic orbit, so let's
187 compute it as an elliptic or hyperbolic orbit instead. */
188 goto ellipse_hyperbola;
189 }
190 }
191
192 /* Go ahead computing the near-parabolic case */
193 c = 1.0 + 1.0 / (w*w);
194 c1 = 1.0 / c;
195 c2 = c1*c1;
196 c3 = c1*c2;
197 c5 = c3*c2;
198 w2 = w*w;
199 w4 = w2*w2;
200 a0 = w;
201 a1 = 2.0 * w * (0.33333333 + 0.2*w2) * c1;
202 a2 = 0.2 * w * (7.0 + 0.14285714 * (33.0*w2+7.4*w4)) * c3;
203 a3 = 0.022857143 * (108.0 + 37.177777*w2 + 5.1111111*w4) * c5;
204 w = (( lambda * a3 + a2 ) * lambda + a1 ) * lambda + a0;
205 w2 = w*w;
206 *v = 2.0 * atand(w);
207 *r = q * (1+w2) / ( 1.0 + w2*lambda );
208 return 0; /* Near-parabolic orbit */
209 }
210
211
212ellipse_hyperbola:
213
214 if ( lambda > 0.0 )
215 { /* Elliptic orbit: */
216 a = q / (1.0-e); /* Semi-major axis */
217 m = KD * tp / sqrt(a*a*a); /* Mean Anomaly */
218 vr( v, r, m, e, a ); /* Solve Kepler's equation, etc */
219 }
220 else
221 { /* Hyperbolic orbit: */
222 a = q / (e-1.0); /* Semi-major axis */
223 n = K * tp / sqrt(a*a*a); /* "Daily motion" */
224 g = n/e;
225 adgg = INFINITY;
226 do
227 {
228 adgg2 = adgg;
229 gs = sqrt(g*g+1.0);
230 dg = -( e*g - log(g+gs) - n ) / ( e - 1.0/gs );
231 g = g + dg;
232 adgg = fabs(dg/g);
233 } while ( adgg < adgg2 && adgg > 1E-5 );
234 gs = sqrt(g*g+1.0);
235 *v = 2.0 * atand( sqrt( (e+1.0)/(e-1.0) ) * g / (gs+1.0) );
236 *r = q * (1.0+e) / ( 1.0 + e*cosd(*v) );
237 }
238 return 0;
239
240} /* vrc */
241
242/* For RCS Only -- Do Not Edit */
243static char *rcsid[2] = {(char *)rcsid, "@(#) $RCSfile: twobody.c,v $ $Date: 2005-08-21 10:02:40 $ $Revision: 1.4 $ $Name: not supported by cvs2svn $"};
Note: See TracBrowser for help on using the repository browser.