tle_linestrings_from_orbital_parameters#
- cysgp4.tle_linestrings_from_orbital_parameters(unicode sat_name, int sat_nr, double mjd_epoch, double inclination_deg, double raan_deg, double eccentricity, double argument_of_perigee_deg, double mean_anomaly_deg, double mean_motion_per_day) tuple #
Generate TLE strings from orbital parameters.
See TLE line strings for a description of the two-line element format. The parameters are the Orbital/Keplarian elements and Wikipedia has a really good description what each of them means.
Note: The epoch (date/time) used in TLEs has a very strange format: first two digits are the year, next three digits are the day from beginning of year, then fraction of a day is given, e.g. 20180.25 would be 2020, day 180, 6 hours (UTC).
- Parameters:
- sat_name
str
The satellite name.
- sat_nr
str
The satellite number.
- mjd_epoch
float
The epoch of the orbital parameters as MJD [day]. The function will convert this to the TLE epoch format.
- inclination_deg
float
Inclination of the orbit [deg].
- raan_deg
float
Right ascension of the ascending node of the orbit [deg].
- eccentricity
float
Eccentricity of the orbit [dimensionless]. Note that in the TLE only the decimal digits are stored.
- argument_of_perigee_deg
float
Argument of Perigee [deg].
- mean_anomaly_deg
float
Mean anomaly of the node [deg].
- mean_motion_per_day
float
Mean motion of the satellite [1 / day].
- sat_name
- Returns:
Examples
A simple use case would be like:
>>> from datetime import datetime >>> import cysgp4 >>> # Define satellite orbital parameters >>> sat_name, sat_nr = 'MYSAT', 1 >>> alt_km = 3000. # satellite altitude >>> mean_motion = cysgp4.satellite_mean_motion(alt_km) >>> inclination = 10. # deg >>> raan = 35. # deg >>> eccentricity = 0.0001 >>> argument_of_perigee = 0. # deg >>> mean_anomaly = 112. # deg >>> # assume, the parameters are valid for the following time >>> dt = datetime(2019, 11, 2, 2, 5, 14) >>> pydt = cysgp4.PyDateTime(dt) >>> pydt <PyDateTime: 2019-11-02 02:05:14.000000 UTC> >>> mjd_epoch = pydt.mjd >>> tle_tuple = cysgp4.tle_linestrings_from_orbital_parameters( ... sat_name, ... sat_nr, ... mjd_epoch, ... inclination, ... raan, ... eccentricity, ... argument_of_perigee, ... mean_anomaly, ... mean_motion, ... ) >>> tle_tuple ('MYSAT', '1 00001U 20001A 19306.08696759 .00000000 00000-0 50000-4 0 05', '2 00001 10.0000 35.0000 0001000 0.0000 112.0000 9.55934723 04')