-
Notifications
You must be signed in to change notification settings - Fork 12
/
palEl2ue.c
351 lines (307 loc) · 10.7 KB
/
palEl2ue.c
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
/*
*+
* Name:
* palEl2ue
* Purpose:
* Transform conventional elements into "universal" form
* Language:
* Starlink ANSI C
* Type of Module:
* Library routine
* Invocation:
* void palEl2ue ( double date, int jform, double epoch, double orbinc,
* double anode, double perih, double aorq, double e,
* double aorl, double dm, double u[13], int *jstat );
* Arguments:
* date = double (Given)
* Epoch (TT MJD) of osculation (Note 3)
* jform = int (Given)
* Element set actually returned (1-3; Note 6)
* epoch = double (Given)
* Epoch of elements (TT MJD)
* orbinc = double (Given)
* inclination (radians)
* anode = double (Given)
* longitude of the ascending node (radians)
* perih = double (Given)
* longitude or argument of perihelion (radians)
* aorq = double (Given)
* mean distance or perihelion distance (AU)
* e = double (Given)
* eccentricity
* aorl = double (Given)
* mean anomaly or longitude (radians, JFORM=1,2 only)
* dm = double (Given)
* daily motion (radians, JFORM=1 only)
* u = double [13] (Returned)
* Universal orbital elements (Note 1)
* - (0) combined mass (M+m)
* - (1) total energy of the orbit (alpha)
* - (2) reference (osculating) epoch (t0)
* - (3-5) position at reference epoch (r0)
* - (6-8) velocity at reference epoch (v0)
* - (9) heliocentric distance at reference epoch
* - (10) r0.v0
* - (11) date (t)
* - (12) universal eccentric anomaly (psi) of date, approx
* jstat = int * (Returned)
* status: 0 = OK
* - -1 = illegal JFORM
* - -2 = illegal E
* - -3 = illegal AORQ
* - -4 = illegal DM
* - -5 = numerical error
* Description:
* Transform conventional osculating elements into "universal" form.
* Authors:
* PTW: Pat Wallace (STFC)
* TIMJ: Tim Jenness (JAC, Hawaii)
* {enter_new_authors_here}
* Notes:
* - The "universal" elements are those which define the orbit for the
* purposes of the method of universal variables (see reference).
* They consist of the combined mass of the two bodies, an epoch,
* and the position and velocity vectors (arbitrary reference frame)
* at that epoch. The parameter set used here includes also various
* quantities that can, in fact, be derived from the other
* information. This approach is taken to avoiding unnecessary
* computation and loss of accuracy. The supplementary quantities
* are (i) alpha, which is proportional to the total energy of the
* orbit, (ii) the heliocentric distance at epoch, (iii) the
* outwards component of the velocity at the given epoch, (iv) an
* estimate of psi, the "universal eccentric anomaly" at a given
* date and (v) that date.
* - The companion routine is palUe2pv. This takes the set of numbers
* that the present routine outputs and uses them to derive the
* object's position and velocity. A single prediction requires one
* call to the present routine followed by one call to palUe2pv;
* for convenience, the two calls are packaged as the routine
* palPlanel. Multiple predictions may be made by again calling the
* present routine once, but then calling palUe2pv multiple times,
* which is faster than multiple calls to palPlanel.
* - DATE is the epoch of osculation. It is in the TT timescale
* (formerly Ephemeris Time, ET) and is a Modified Julian Date
* (JD-2400000.5).
* - The supplied orbital elements are with respect to the J2000
* ecliptic and equinox. The position and velocity parameters
* returned in the array U are with respect to the mean equator and
* equinox of epoch J2000, and are for the perihelion prior to the
* specified epoch.
* - The universal elements returned in the array U are in canonical
* units (solar masses, AU and canonical days).
* - Three different element-format options are available:
*
* Option JFORM=1, suitable for the major planets:
*
* EPOCH = epoch of elements (TT MJD)
* ORBINC = inclination i (radians)
* ANODE = longitude of the ascending node, big omega (radians)
* PERIH = longitude of perihelion, curly pi (radians)
* AORQ = mean distance, a (AU)
* E = eccentricity, e (range 0 to <1)
* AORL = mean longitude L (radians)
* DM = daily motion (radians)
*
* Option JFORM=2, suitable for minor planets:
*
* EPOCH = epoch of elements (TT MJD)
* ORBINC = inclination i (radians)
* ANODE = longitude of the ascending node, big omega (radians)
* PERIH = argument of perihelion, little omega (radians)
* AORQ = mean distance, a (AU)
* E = eccentricity, e (range 0 to <1)
* AORL = mean anomaly M (radians)
*
* Option JFORM=3, suitable for comets:
*
* EPOCH = epoch of perihelion (TT MJD)
* ORBINC = inclination i (radians)
* ANODE = longitude of the ascending node, big omega (radians)
* PERIH = argument of perihelion, little omega (radians)
* AORQ = perihelion distance, q (AU)
* E = eccentricity, e (range 0 to 10)
*
* - Unused elements (DM for JFORM=2, AORL and DM for JFORM=3) are
* not accessed.
* - The algorithm was originally adapted from the EPHSLA program of
* D.H.P.Jones (private communication, 1996). The method is based
* on Stumpff's Universal Variables.
*
* See Also:
* Everhart & Pitkin, Am.J.Phys. 51, 712 (1983).
* History:
* 2012-03-12 (TIMJ):
* Initial version taken directly from SLA/F.
* Adapted with permission from the Fortran SLALIB library.
* {enter_further_changes_here}
* Copyright:
* Copyright (C) 2005 Patrick T. Wallace
* Copyright (C) 2012 Science and Technology Facilities Council.
* All Rights Reserved.
* Licence:
* This program is free software; you can redistribute it and/or
* modify it under the terms of the GNU General Public License as
* published by the Free Software Foundation; either version 3 of
* the License, or (at your option) any later version.
*
* This program is distributed in the hope that it will be
* useful, but WITHOUT ANY WARRANTY; without even the implied
* warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
* PURPOSE. See the GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
* MA 02110-1301, USA.
* Bugs:
* {note_any_bugs_here}
*-
*/
#include <math.h>
#include "pal.h"
#include "palmac.h"
void palEl2ue ( double date, int jform, double epoch, double orbinc,
double anode, double perih, double aorq, double e,
double aorl, double dm, double u[13], int *jstat ) {
/* Sin and cos of J2000 mean obliquity (IAU 1976) */
const double SE=0.3977771559319137;
const double CE=0.9174820620691818;
int J;
double PHT,ARGPH,Q,W,CM,ALPHA,PHS,SW,CW,SI,CI,SO,CO,
X,Y,Z,PX,PY,PZ,VX,VY,VZ,DT,FC,FP,PSI,
UL[13],PV[6];
/* Validate arguments. */
if (jform < 1 || jform > 3) {
*jstat = -1;
return;
}
if (e < 0.0 || e > 10.0 || (e >= 1.0 && jform != 3)) {
*jstat = -2;
return;
}
if (aorq <= 0.0) {
*jstat = -3;
return;
}
if (jform == 1 && dm <= 0.0) {
*jstat = -4;
return;
}
/*
* Transform elements into standard form:
*
* PHT = epoch of perihelion passage
* ARGPH = argument of perihelion (little omega)
* Q = perihelion distance (q)
* CM = combined mass, M+m (mu)
*/
if (jform == 1) {
/* Major planet. */
PHT = epoch-(aorl-perih)/dm;
ARGPH = perih-anode;
Q = aorq*(1.0-e);
W = dm/PAL__GCON;
CM = W*W*aorq*aorq*aorq;
} else if (jform == 2) {
/* Minor planet. */
PHT = epoch-aorl*sqrt(aorq*aorq*aorq)/PAL__GCON;
ARGPH = perih;
Q = aorq*(1.0-e);
CM = 1.0;
} else {
/* Comet. */
PHT = epoch;
ARGPH = perih;
Q = aorq;
CM = 1.0;
}
/* The universal variable alpha. This is proportional to the total
* energy of the orbit: -ve for an ellipse, zero for a parabola,
* +ve for a hyperbola. */
ALPHA = CM*(e-1.0)/Q;
/* Speed at perihelion. */
PHS = sqrt(ALPHA+2.0*CM/Q);
/* In a Cartesian coordinate system which has the x-axis pointing
* to perihelion and the z-axis normal to the orbit (such that the
* object orbits counter-clockwise as seen from +ve z), the
* perihelion position and velocity vectors are:
*
* position [Q,0,0]
* velocity [0,PHS,0]
*
* To express the results in J2000 equatorial coordinates we make a
* series of four rotations of the Cartesian axes:
*
* axis Euler angle
*
* 1 z argument of perihelion (little omega)
* 2 x inclination (i)
* 3 z longitude of the ascending node (big omega)
* 4 x J2000 obliquity (epsilon)
*
* In each case the rotation is clockwise as seen from the +ve end of
* the axis concerned.
*/
/* Functions of the Euler angles. */
SW = sin(ARGPH);
CW = cos(ARGPH);
SI = sin(orbinc);
CI = cos(orbinc);
SO = sin(anode);
CO = cos(anode);
/* Position at perihelion (AU). */
X = Q*CW;
Y = Q*SW;
Z = Y*SI;
Y = Y*CI;
PX = X*CO-Y*SO;
Y = X*SO+Y*CO;
PY = Y*CE-Z*SE;
PZ = Y*SE+Z*CE;
/* Velocity at perihelion (AU per canonical day). */
X = -PHS*SW;
Y = PHS*CW;
Z = Y*SI;
Y = Y*CI;
VX = X*CO-Y*SO;
Y = X*SO+Y*CO;
VY = Y*CE-Z*SE;
VZ = Y*SE+Z*CE;
/* Time from perihelion to date (in Canonical Days: a canonical day
* is 58.1324409... days, defined as 1/PAL__GCON). */
DT = (date-PHT)*PAL__GCON;
/* First approximation to the Universal Eccentric Anomaly, PSI,
* based on the circle (FC) and parabola (FP) values. */
FC = DT/Q;
W = pow(3.0*DT+sqrt(9.0*DT*DT+8.0*Q*Q*Q), 1.0/3.0);
FP = W-2.0*Q/W;
PSI = (1.0-e)*FC+e*FP;
/* Assemble local copy of element set. */
UL[0] = CM;
UL[1] = ALPHA;
UL[2] = PHT;
UL[3] = PX;
UL[4] = PY;
UL[5] = PZ;
UL[6] = VX;
UL[7] = VY;
UL[8] = VZ;
UL[9] = Q;
UL[10] = 0.0;
UL[11] = date;
UL[12] = PSI;
/* Predict position+velocity at epoch of osculation. */
palUe2pv( date, UL, PV, &J );
if (J != 0) {
*jstat = -5;
return;
}
/* Convert back to universal elements. */
palPv2ue( PV, date, CM-1.0, u, &J );
if (J != 0) {
*jstat = -5;
return;
}
/* OK exit. */
*jstat = 0;
}