source: trunk/FACT++/pal/palPlante.c@ 18347

Last change on this file since 18347 was 18347, checked in by tbretz, 9 years ago
File size: 9.9 KB
Line 
1/*
2*+
3* Name:
4* palPlante
5
6* Purpose:
7* Topocentric RA,Dec of a Solar-System object from heliocentric orbital elements
8
9* Language:
10* Starlink ANSI C
11
12* Type of Module:
13* Library routine
14
15* Invocation:
16* void palPlante ( double date, double elong, double phi, int jform,
17* double epoch, double orbinc, double anode, double perih,
18* double aorq, double e, double aorl, double dm,
19* double *ra, double *dec, double *r, int *jstat );
20
21
22* Description:
23* Topocentric apparent RA,Dec of a Solar-System object whose
24* heliocentric orbital elements are known.
25
26* Arguments:
27* date = double (Given)
28* TT MJD of observation (JD-2400000.5)
29* elong = double (Given)
30* Observer's east longitude (radians)
31* phi = double (Given)
32* Observer's geodetic latitude (radians)
33* jform = int (Given)
34* Element set actually returned (1-3; Note 6)
35* epoch = double (Given)
36* Epoch of elements (TT MJD)
37* orbinc = double (Given)
38* inclination (radians)
39* anode = double (Given)
40* longitude of the ascending node (radians)
41* perih = double (Given)
42* longitude or argument of perihelion (radians)
43* aorq = double (Given)
44* mean distance or perihelion distance (AU)
45* e = double (Given)
46* eccentricity
47* aorl = double (Given)
48* mean anomaly or longitude (radians, JFORM=1,2 only)
49* dm = double (Given)
50* daily motion (radians, JFORM=1 only)
51* ra = double * (Returned)
52* Topocentric apparent RA (radians)
53* dec = double * (Returned)
54* Topocentric apparent Dec (radians)
55* r = double * (Returned)
56* Distance from observer (AU)
57* jstat = int * (Returned)
58* status: 0 = OK
59* - -1 = illegal jform
60* - -2 = illegal e
61* - -3 = illegal aorq
62* - -4 = illegal dm
63* - -5 = numerical error
64
65* Authors:
66* PTW: Pat Wallace (STFC)
67* TIMJ: Tim Jenness (JAC, Hawaii)
68* {enter_new_authors_here}
69
70* Notes:
71* - DATE is the instant for which the prediction is required. It is
72* in the TT timescale (formerly Ephemeris Time, ET) and is a
73* Modified Julian Date (JD-2400000.5).
74* - The longitude and latitude allow correction for geocentric
75* parallax. This is usually a small effect, but can become
76* important for near-Earth asteroids. Geocentric positions can be
77* generated by appropriate use of routines palEpv (or palEvp) and
78* palUe2pv.
79* - The elements are with respect to the J2000 ecliptic and equinox.
80* - A choice of three different element-set options is available:
81*
82* Option JFORM = 1, suitable for the major planets:
83*
84* EPOCH = epoch of elements (TT MJD)
85* ORBINC = inclination i (radians)
86* ANODE = longitude of the ascending node, big omega (radians)
87* PERIH = longitude of perihelion, curly pi (radians)
88* AORQ = mean distance, a (AU)
89* E = eccentricity, e (range 0 to <1)
90* AORL = mean longitude L (radians)
91* DM = daily motion (radians)
92*
93* Option JFORM = 2, suitable for minor planets:
94*
95* EPOCH = epoch of elements (TT MJD)
96* ORBINC = inclination i (radians)
97* ANODE = longitude of the ascending node, big omega (radians)
98* PERIH = argument of perihelion, little omega (radians)
99* AORQ = mean distance, a (AU)
100* E = eccentricity, e (range 0 to <1)
101* AORL = mean anomaly M (radians)
102*
103* Option JFORM = 3, suitable for comets:
104*
105* EPOCH = epoch of elements and perihelion (TT MJD)
106* ORBINC = inclination i (radians)
107* ANODE = longitude of the ascending node, big omega (radians)
108* PERIH = argument of perihelion, little omega (radians)
109* AORQ = perihelion distance, q (AU)
110* E = eccentricity, e (range 0 to 10)
111*
112* Unused arguments (DM for JFORM=2, AORL and DM for JFORM=3) are not
113* accessed.
114* - Each of the three element sets defines an unperturbed heliocentric
115* orbit. For a given epoch of observation, the position of the body
116* in its orbit can be predicted from these elements, which are
117* called "osculating elements", using standard two-body analytical
118* solutions. However, due to planetary perturbations, a given set
119* of osculating elements remains usable for only as long as the
120* unperturbed orbit that it describes is an adequate approximation
121* to reality. Attached to such a set of elements is a date called
122* the "osculating epoch", at which the elements are, momentarily,
123* a perfect representation of the instantaneous position and
124* velocity of the body.
125*
126* Therefore, for any given problem there are up to three different
127* epochs in play, and it is vital to distinguish clearly between
128* them:
129*
130* . The epoch of observation: the moment in time for which the
131* position of the body is to be predicted.
132*
133* . The epoch defining the position of the body: the moment in time
134* at which, in the absence of purturbations, the specified
135* position (mean longitude, mean anomaly, or perihelion) is
136* reached.
137*
138* . The osculating epoch: the moment in time at which the given
139* elements are correct.
140*
141* For the major-planet and minor-planet cases it is usual to make
142* the epoch that defines the position of the body the same as the
143* epoch of osculation. Thus, only two different epochs are
144* involved: the epoch of the elements and the epoch of observation.
145*
146* For comets, the epoch of perihelion fixes the position in the
147* orbit and in general a different epoch of osculation will be
148* chosen. Thus, all three types of epoch are involved.
149*
150* For the present routine:
151*
152* . The epoch of observation is the argument DATE.
153*
154* . The epoch defining the position of the body is the argument
155* EPOCH.
156*
157* . The osculating epoch is not used and is assumed to be close
158* enough to the epoch of observation to deliver adequate accuracy.
159* If not, a preliminary call to palPertel may be used to update
160* the element-set (and its associated osculating epoch) by
161* applying planetary perturbations.
162* - Two important sources for orbital elements are Horizons, operated
163* by the Jet Propulsion Laboratory, Pasadena, and the Minor Planet
164* Center, operated by the Center for Astrophysics, Harvard.
165*
166* The JPL Horizons elements (heliocentric, J2000 ecliptic and
167* equinox) correspond to PAL/SLALIB arguments as follows.
168*
169* Major planets:
170*
171* JFORM = 1
172* EPOCH = JDCT-2400000.5
173* ORBINC = IN (in radians)
174* ANODE = OM (in radians)
175* PERIH = OM+W (in radians)
176* AORQ = A
177* E = EC
178* AORL = MA+OM+W (in radians)
179* DM = N (in radians)
180*
181* Epoch of osculation = JDCT-2400000.5
182*
183* Minor planets:
184*
185* JFORM = 2
186* EPOCH = JDCT-2400000.5
187* ORBINC = IN (in radians)
188* ANODE = OM (in radians)
189* PERIH = W (in radians)
190* AORQ = A
191* E = EC
192* AORL = MA (in radians)
193*
194* Epoch of osculation = JDCT-2400000.5
195*
196* Comets:
197*
198* JFORM = 3
199* EPOCH = Tp-2400000.5
200* ORBINC = IN (in radians)
201* ANODE = OM (in radians)
202* PERIH = W (in radians)
203* AORQ = QR
204* E = EC
205*
206* Epoch of osculation = JDCT-2400000.5
207*
208* The MPC elements correspond to SLALIB arguments as follows.
209*
210* Minor planets:
211*
212* JFORM = 2
213* EPOCH = Epoch-2400000.5
214* ORBINC = Incl. (in radians)
215* ANODE = Node (in radians)
216* PERIH = Perih. (in radians)
217* AORQ = a
218* E = e
219* AORL = M (in radians)
220*
221* Epoch of osculation = Epoch-2400000.5
222*
223* Comets:
224*
225* JFORM = 3
226* EPOCH = T-2400000.5
227* ORBINC = Incl. (in radians)
228* ANODE = Node. (in radians)
229* PERIH = Perih. (in radians)
230* AORQ = q
231* E = e
232*
233* Epoch of osculation = Epoch-2400000.5
234
235* History:
236* 2012-03-12 (TIMJ):
237* Initial version direct conversion of SLA/F.
238* Adapted with permission from the Fortran SLALIB library.
239* {enter_further_changes_here}
240
241* Copyright:
242* Copyright (C) 2004 Patrick T. Wallace
243* Copyright (C) 2012 Science and Technology Facilities Council.
244* All Rights Reserved.
245
246* Licence:
247* This program is free software; you can redistribute it and/or
248* modify it under the terms of the GNU General Public License as
249* published by the Free Software Foundation; either version 3 of
250* the License, or (at your option) any later version.
251*
252* This program is distributed in the hope that it will be
253* useful, but WITHOUT ANY WARRANTY; without even the implied
254* warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
255* PURPOSE. See the GNU General Public License for more details.
256*
257* You should have received a copy of the GNU General Public License
258* along with this program; if not, write to the Free Software
259* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
260* MA 02110-1301, USA.
261
262* Bugs:
263* {note_any_bugs_here}
264*-
265*/
266
267#include "pal.h"
268
269void palPlante ( double date, double elong, double phi, int jform,
270 double epoch, double orbinc, double anode, double perih,
271 double aorq, double e, double aorl, double dm,
272 double *ra, double *dec, double *r, int *jstat ) {
273
274 double u[13];
275
276 /* Transform conventional elements to universal elements */
277 palEl2ue( date, jform, epoch, orbinc, anode, perih, aorq, e, aorl,
278 dm, u, jstat );
279
280 /* If succcessful, make the prediction */
281 if (*jstat == 0) palPlantu( date, elong, phi, u, ra, dec, r, jstat );
282
283}
284
285
286
Note: See TracBrowser for help on using the repository browser.