pymgp.tleorb.tle2vec#
- pymgp.tleorb.tle2vec(tle, t, satid, propagation='J2', verbose=0)[source]#
Compute satellite position and velocity in ECI from NORAD Two Line Elements.
Compute satellite position and velocity in ECI coordinates from NORAD two line elements, for selected satellites, at given times.
- Parameters:
- tlelist with named tuples
List with two line elements as named tuples as return by tleread.
- tarray_like, float, with shape (n,)
Array with serial date numbers in UT1 or date range as used by tledatenum. Valid inputs are:
three element list with start date, end date (or duration in minutes) and data interval in minutes,
string with the date or a list with date strings, or,
ndarray with serial datenumbers, datestrings, or datetime64
- satidstr, int, array_like with int or str
String or array with the satellite names, index or index array. For strings, satellites are selected that start with satid. Uses tlefind to resolve the satellite identifiers.
- propagation{‘J2’, ‘NOJ2’, ‘SGP4’}, optional
Propagation method to for the orbital elements:
J2 Include secular effects of J2, but nothing else (default). This gives acceptable results for plotting purposes and planning computations NOJ2 Ignores effect of J2 on orbit propagation. Should only be used for educational purposes SGP4 SGP4 orbit propagator of NORAD. Not implemented yet.
- verboseint, optional
Verbosity level. If ‘verbose=1’ (default) print a message for the found TLE(s), if ‘verbose=0’ the function is quiet.
- Returns:
- xsat, vsatndarray with shape (n,3) each
The satellite position (m) and velocity (m/s) in an ECI reference frame. The second to last axis is time t, the last axis are the coordinates X, Y and Z, respectively velocities VX, VY and VZ.
- tarray_like, float, with shape (n,)
Array with serial date numbers in UT1.
- satidsarray_like, str
Array with the satellite names.
See also
Notes
The default propagation method is J2 which includes the secular effects of J2. This gives acceptable results for plotting purposes and planning computations. Formally, TLE are only compatible with the SGP4 orbit propagator of NORAD, but this is not (yet) supported by this function. The other available option is NOJ2, which ignores the effect of J2 on the orbit propagation and should only be used for educational purposes.
Examples
Read Earth resource satellite TLE’s and get position and velocity of RADARSAT-2 at a specific time
>>> tleERS = tleread('resource-20240802.tle', verbose=0) >>> xsat, vsat, *_ = tle2vec(tleERS,'2024-08-02 05:58:10','RADARSAT-2') >>> print(xsat, vsat) [[5458404.60763595 4591181.8357549 -732126.12944127]] [[ 144.29836955 -1339.50959494 -7333.2677589 ]]
Compute RADARSAT-2 position and velocities for Aug 2, 2024, with one minute intervals
>>> xsat, vsat, t, *_ = tle2vec(tleERS,['2024-08-02 0:00', 24*60, 1],'RADARSAT-2') >>> print(xsat.shape, vsat.shape, t.shape) (1441, 3) (1441, 3) (1441,)
>>> t = tledatenum(['2024-08-02 0:00', 24*60, 1]) # array with date numbers ... >>> xsat, vsat, t, *_ = tle2vec(tleERS,t,'RADARSAT-2') >>> print(xsat.shape, vsat.shape, t.shape) (1441, 3) (1441, 3) (1441,)
Position and velocity, on Aug 2, 2024, for every SENTINEL satellite, with one minute interval
>>> xsat, vsat, t, satids = tle2vec(tleERS,['2024-08-02 0:00', 24*60, 1],'SENTINEL') >>> print(xsat.shape, vsat.shape, t.shape, satids.shape) >>> print(satids) (7, 1441, 3) (7, 1441, 3) (1441,) (7,) ['SENTINEL-1A' 'SENTINEL-2A' 'SENTINEL-3A' 'SENTINEL-2B' 'SENTINEL-5P' 'SENTINEL-3B' 'SENTINEL-6']