libnova  v 0.15.0
elliptic_motion.h
1 /*
2  * This library is free software; you can redistribute it and/or
3  * modify it under the terms of the GNU Lesser General Public
4  * License as published by the Free Software Foundation; either
5  * version 2 of the License, or (at your option) any later version.
6  *
7  * This library is distributed in the hope that it will be useful,
8  * but WITHOUT ANY WARRANTY; without even the implied warranty of
9  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
10  * Lesser General Public License for more details.
11  *
12  * You should have received a copy of the GNU General Public License
13  * along with this program; if not, write to the Free Software
14  * Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
15  *
16  * Copyright (C) 2000 - 2005 Liam Girdwood
17  */
18 
19 #ifndef _LN_ELLIPTIC_MOTION_H
20 #define _LN_ELLIPTIC_MOTION_H
21 
22 #include <libnova/ln_types.h>
23 
24 #ifdef __cplusplus
25 extern "C" {
26 #endif
27 
39 double LIBNOVA_EXPORT ln_solve_kepler (double e, double M);
40 
45 double LIBNOVA_EXPORT ln_get_ell_mean_anomaly (double n, double delta_JD);
46 
51 double LIBNOVA_EXPORT ln_get_ell_true_anomaly (double e, double E);
52 
57 double LIBNOVA_EXPORT ln_get_ell_radius_vector (double a, double e, double E);
58 
63 double LIBNOVA_EXPORT ln_get_ell_smajor_diam (double e, double q);
64 
69 double LIBNOVA_EXPORT ln_get_ell_sminor_diam (double e, double a);
70 
75 double LIBNOVA_EXPORT ln_get_ell_mean_motion (double a);
76 
81 void LIBNOVA_EXPORT ln_get_ell_geo_rect_posn (struct ln_ell_orbit* orbit, double JD, struct ln_rect_posn* posn);
82 
87 void LIBNOVA_EXPORT ln_get_ell_helio_rect_posn (struct ln_ell_orbit* orbit, double JD, struct ln_rect_posn* posn);
88 
93 double LIBNOVA_EXPORT ln_get_ell_orbit_len (struct ln_ell_orbit * orbit);
94 
99 double LIBNOVA_EXPORT ln_get_ell_orbit_vel (double JD, struct ln_ell_orbit * orbit);
100 
105 double LIBNOVA_EXPORT ln_get_ell_orbit_pvel (struct ln_ell_orbit * orbit);
106 
111 double LIBNOVA_EXPORT ln_get_ell_orbit_avel (struct ln_ell_orbit * orbit);
112 
117 double LIBNOVA_EXPORT ln_get_ell_body_phase_angle (double JD, struct ln_ell_orbit * orbit);
118 
123 double LIBNOVA_EXPORT ln_get_ell_body_elong (double JD, struct ln_ell_orbit * orbit);
124 
130 double LIBNOVA_EXPORT ln_get_ell_body_solar_dist (double JD, struct ln_ell_orbit * orbit);
131 
137 double LIBNOVA_EXPORT ln_get_ell_body_earth_dist (double JD, struct ln_ell_orbit * orbit);
138 
144 void LIBNOVA_EXPORT ln_get_ell_body_equ_coords (double JD, struct ln_ell_orbit * orbit, struct ln_equ_posn * posn);
145 
150 int LIBNOVA_EXPORT ln_get_ell_body_rst (double JD, struct ln_lnlat_posn * observer, struct ln_ell_orbit * orbit, struct ln_rst_time * rst);
151 
156 int LIBNOVA_EXPORT ln_get_ell_body_rst_horizon (double JD, struct ln_lnlat_posn * observer, struct ln_ell_orbit * orbit, double horizon, struct ln_rst_time * rst);
157 
162 int LIBNOVA_EXPORT ln_get_ell_body_next_rst (double JD, struct ln_lnlat_posn * observer, struct ln_ell_orbit * orbit, struct ln_rst_time * rst);
163 
168 int LIBNOVA_EXPORT ln_get_ell_body_next_rst_horizon (double JD, struct ln_lnlat_posn * observer, struct ln_ell_orbit * orbit, double horizon, struct ln_rst_time * rst);
169 
174 int LIBNOVA_EXPORT ln_get_ell_body_next_rst_horizon_future (double JD, struct ln_lnlat_posn * observer, struct ln_ell_orbit * orbit, double horizon, int day_limit, struct ln_rst_time * rst);
175 
180 double LIBNOVA_EXPORT ln_get_ell_last_perihelion (double epoch_JD, double M, double n);
181 
182 #ifdef __cplusplus
183 };
184 #endif
185 
186 #endif
int LIBNOVA_EXPORT ln_get_ell_body_rst_horizon(double JD, struct ln_lnlat_posn *observer, struct ln_ell_orbit *orbit, double horizon, struct ln_rst_time *rst)
Calculate the time of rise, set and transit for a body with an elliptic orbit.
Definition: elliptic_motion.c:517
double LIBNOVA_EXPORT ln_get_ell_radius_vector(double a, double e, double E)
Calculate the radius vector.
Definition: elliptic_motion.c:129
Elliptic Orbital elements.
Definition: ln_types.h:265
Rectangular coordinates.
Definition: ln_types.h:237
void LIBNOVA_EXPORT ln_get_ell_geo_rect_posn(struct ln_ell_orbit *orbit, double JD, struct ln_rect_posn *posn)
Calculate the objects rectangular geocentric position.
Definition: elliptic_motion.c:243
double LIBNOVA_EXPORT ln_get_ell_body_solar_dist(double JD, struct ln_ell_orbit *orbit)
Calculate the distance between a body and the Sun.
Definition: elliptic_motion.c:378
double LIBNOVA_EXPORT ln_get_ell_mean_motion(double a)
Calculate the mean daily motion (degrees/day).
Definition: elliptic_motion.c:165
double LIBNOVA_EXPORT ln_get_ell_body_elong(double JD, struct ln_ell_orbit *orbit)
Calculate the bodies elongation to the Sun..
Definition: elliptic_motion.c:456
double LIBNOVA_EXPORT ln_get_ell_orbit_vel(double JD, struct ln_ell_orbit *orbit)
Calculate orbital velocity in km/s.
Definition: elliptic_motion.c:329
double LIBNOVA_EXPORT ln_get_ell_last_perihelion(double epoch_JD, double M, double n)
Calculate the julian day of the last perihelion.
Definition: elliptic_motion.c:592
int LIBNOVA_EXPORT ln_get_ell_body_next_rst_horizon(double JD, struct ln_lnlat_posn *observer, struct ln_ell_orbit *orbit, double horizon, struct ln_rst_time *rst)
Calculate the time of rise, set and transit for a body with an elliptic orbit.
Definition: elliptic_motion.c:558
double LIBNOVA_EXPORT ln_get_ell_mean_anomaly(double n, double delta_JD)
Calculate the mean anomaly.
Definition: elliptic_motion.c:96
Equatorial Coordinates.
Definition: ln_types.h:170
double LIBNOVA_EXPORT ln_get_ell_true_anomaly(double e, double E)
Calculate the true anomaly.
Definition: elliptic_motion.c:109
void LIBNOVA_EXPORT ln_get_ell_helio_rect_posn(struct ln_ell_orbit *orbit, double JD, struct ln_rect_posn *posn)
Calculate the objects rectangular heliocentric position.
Definition: elliptic_motion.c:179
Ecliptical (or celestial) Longitude and Latitude.
Definition: ln_types.h:200
double LIBNOVA_EXPORT ln_get_ell_smajor_diam(double e, double q)
Calculate the semi major diameter.
Definition: elliptic_motion.c:142
void LIBNOVA_EXPORT ln_get_ell_body_equ_coords(double JD, struct ln_ell_orbit *orbit, struct ln_equ_posn *posn)
Calculate a bodies equatorial coords.
Definition: elliptic_motion.c:269
Rise, Set and Transit times.
Definition: ln_types.h:317
int LIBNOVA_EXPORT ln_get_ell_body_rst(double JD, struct ln_lnlat_posn *observer, struct ln_ell_orbit *orbit, struct ln_rst_time *rst)
Calculate the time of rise, set and transit for a body with an elliptic orbit.
Definition: elliptic_motion.c:498
double LIBNOVA_EXPORT ln_get_ell_orbit_pvel(struct ln_ell_orbit *orbit)
Calculate orbital velocity at perihelion in km/s.
Definition: elliptic_motion.c:346
int LIBNOVA_EXPORT ln_get_ell_body_next_rst(double JD, struct ln_lnlat_posn *observer, struct ln_ell_orbit *orbit, struct ln_rst_time *rst)
Calculate the time of rise, set and transit for a body with an elliptic orbit.
Definition: elliptic_motion.c:537
double LIBNOVA_EXPORT ln_get_ell_orbit_avel(struct ln_ell_orbit *orbit)
Calculate the orbital velocity at aphelion in km/s.
Definition: elliptic_motion.c:361
double LIBNOVA_EXPORT ln_solve_kepler(double e, double M)
Calculate the eccentric anomaly.
Definition: elliptic_motion.c:53
int LIBNOVA_EXPORT ln_get_ell_body_next_rst_horizon_future(double JD, struct ln_lnlat_posn *observer, struct ln_ell_orbit *orbit, double horizon, int day_limit, struct ln_rst_time *rst)
Calculate the time of rise, set and transit for a body with an elliptic orbit.
Definition: elliptic_motion.c:580
double LIBNOVA_EXPORT ln_get_ell_body_phase_angle(double JD, struct ln_ell_orbit *orbit)
Calculate the phase angle of the body. The angle Sun - body - Earth.
Definition: elliptic_motion.c:423
double LIBNOVA_EXPORT ln_get_ell_body_earth_dist(double JD, struct ln_ell_orbit *orbit)
Calculate the distance between a body and the Earth.
Definition: elliptic_motion.c:401
double LIBNOVA_EXPORT ln_get_ell_sminor_diam(double e, double a)
Calculate the semi minor diameter.
Definition: elliptic_motion.c:154
double LIBNOVA_EXPORT ln_get_ell_orbit_len(struct ln_ell_orbit *orbit)
Calculate the orbital length in AU.
Definition: elliptic_motion.c:308