void conics_c ( ConstSpiceDouble elts[8],
SpiceDouble et,
SpiceDouble state[6] )
Determine the state (position, velocity) of an orbiting body
from a set of elliptic, hyperbolic, or parabolic orbital
elements.
None.
CONIC
EPHEMERIS
VARIABLE I/O DESCRIPTION
-------- --- --------------------------------------------------
elts I Conic elements.
et I Input time.
state O State of orbiting body at et.
elts are conic osculating elements describing the orbit of a
body around a primary. The elements are, in order:
RP Perifocal distance.
ECC Eccentricity.
INC Inclination.
LNODE Longitude of the ascending node.
ARGP Argument of periapse.
M0 Mean anomaly at epoch.
T0 Epoch.
MU Gravitational parameter.
Units are km, rad, rad/sec, km**3/sec**2.
The epoch T0 is given in ephemeris seconds past J2000.
T0 is the instant at which the state of the body is
specified by the elements.
The same elements are used to describe all three types
(elliptic, hyperbolic, and parabolic) of conic orbit.
et is the time at which the state of the orbiting body
is to be determined, in ephemeris seconds J2000.
state is the state (position and velocity) of the body at
time `et'. Components are x, y, z, dx/dt, dy/dt, dz/dt.
None.
1) If the eccentricity supplied is less than 0, the error
SPICE(BADECCENTRICITY) is signaled.
2) If a non-positive periapse distance is supplied, the error
SPICE(BADPERIAPSEVALUE) is signaled.
3) If a non-positive value for the attracting mass is supplied,
the error SPICE(BADGM), is signaled.
4) Errors such as an out of bounds value for `et' are diagnosed
by routines in the call tree of this routine.
None.
None.
Let vinit contain the initial state of a spacecraft relative to the
center of a planet at epoch `et', and let `gm' be the gravitation
parameter of the planet. The call
oscelt_c ( vinit, et, gm, elts );
produces a set of osculating elements describing the nominal
orbit that the spacecraft would follow in the absence of all
other bodies in the solar system and non-gravitational forces
on the spacecraft.
Now let `state' contain the state of the same spacecraft at some
other, later epoch. The difference between this state and the
state predicted by the nominal orbit at the same epoch can be
computed as follows.
conics_c ( elts, later, nominal );
vsubg_c ( nominal, state, 6, diff );
printf( "Perturbation in x, dx/dt = %f %f", diff[0], diff[3] );
printf( " y, dy/dt = %f %f", diff[1], diff[4] );
printf( " z, dz/dt = %f %f", diff[2], diff[5] );
None.
[1] Roger Bate, Fundamentals of Astrodynamics, Dover, 1971.
N.J. Bachman (JPL)
I.M. Underwood (JPL)
W.L. Taber (JPL)
E.D. Wright (JPL)
-CSPICE Version 1.1.1, 29-JUL-2003 (NJB)
Various header corrections were made.
-CSPICE Version 1.1.0, 24-JUL-2001 (NJB)
Changed protoype: input elts is now type (ConstSpiceDouble *).
Implemented interface macro for casting input array to const.
-CSPICE Version 1.0.1, 08-FEB-1998 (EDW)
Corrected and clarified header entries.
-CSPICE Version 1.0.0, 10-NOV-1997 (EDW)
state from conic elements
Link to routine conics_c source file conics_c.c
|