pub struct TLE {Show 19 fields
pub name: String,
pub intl_desig: String,
pub sat_num: i32,
pub desig_year: i32,
pub desig_launch: i32,
pub desig_piece: String,
pub epoch: AstroTime,
pub mean_motion_dot: f64,
pub mean_motion_dot_dot: f64,
pub bstar: f64,
pub ephem_type: u8,
pub element_num: i32,
pub inclination: f64,
pub raan: f64,
pub eccen: f64,
pub arg_of_perigee: f64,
pub mean_anomaly: f64,
pub mean_motion: f64,
pub rev_num: i32,
/* private fields */
}Expand description
Stucture representing a Two-Line Element Set (TLE), a satellite ephemeris format from the 1970s that is still somehow in use today and can be used to calculate satellite position and velcocity in the “TEME” frame (not-quite GCRF) using the “Simplified General Perturbations-4” (SGP-4) mathemematical model that is also included in this package.
For details, see: https://en.wikipedia.org/wiki/Two-line_element_set
The TLE format is still commonly used to represent satellite ephemerides, and satellite ephemerides catalogs in this format are publicly availalble at www.space-track.org (registration required)
TLEs sometimes have a “line 0” that includes the name of the satellite
§Example Usage:
use satkit::TLE;
use satkit::AstroTime;
use satkit::sgp4::sgp4;
use satkit::frametransform;
use satkit::itrfcoord::ITRFCoord;
let lines = vec!["0 INTELSAT 902",
"1 26900U 01039A 06106.74503247 .00000045 00000-0 10000-3 0 8290",
"2 26900 0.0164 266.5378 0003319 86.1794 182.2590 1.00273847 16981 9300."];
let mut tle = TLE::load_3line(lines[0], lines[1], lines[2]).unwrap();
let tm = AstroTime::from_datetime(2006, 5, 1, 11, 0, 0.0);
// Use SGP4 to get position,
let (pTEME, vTEME, errs) = sgp4(&mut tle, &[tm]);
println!("pTEME = {}", pTEME);
// Rotate the position to the ITRF frame (Earth-fixed)
// Since pTEME is a 3xN array where N is the number of times
// (we are just using a single time)
// we need to convert to a fixed matrix to rotate
let pITRF = frametransform::qteme2itrf(&tm) * pTEME.fixed_view::<3,1>(0,0);
println!("pITRF = {}", pITRF);
// Convert to an "ITRF Coordinate" and print geodetic position
let itrf = ITRFCoord::from_slice(&pTEME.fixed_view::<3,1>(0,0).as_slice()).unwrap();
println!("latitude = {} deg", itrf.latitude_deg());
println!("longitude = {} deg", itrf.longitude_deg());
println!("altitude = {} m", itrf.hae());
Fields§
§name: StringName of satellite
intl_desig: StringString describing launch
sat_num: i32Satellite NORAD number
desig_year: i32Launch year
desig_launch: i32Numbered launch of year
desig_piece: StringPiece of launch
epoch: AstroTimeTLE epoch
mean_motion_dot: f64One half of 1st derivative of mean motion wrt time, in revs/day^2
mean_motion_dot_dot: f64One sixth of 2nd derivative of mean motion wrt tim, in revs/day^3
bstar: f64Starred ballistic coefficient, in units of inverse Earth radii
ephem_type: u8Usually 0
element_num: i32Bulliten number
inclination: f64Inclination, degrees
raan: f64Right ascension of ascending node, degrees
eccen: f64Eccentricity
arg_of_perigee: f64Argument of perigee, degrees
mean_anomaly: f64Mean anomaly, degrees
mean_motion: f64Mean motion, revs / day
rev_num: i32Revolution number
Implementations§
source§impl TLE
impl TLE
pub fn from_lines(lines: &Vec<String>) -> SKResult<Vec<TLE>>
sourcepub fn load_3line(line0: &str, line1: &str, line2: &str) -> Result<TLE, String>
pub fn load_3line(line0: &str, line1: &str, line2: &str) -> Result<TLE, String>
Load 3 lines as string into a structure representing a Two-Line Element Set (TLE)
The TLE can then be used to compute satellite position and velocity as a function of time.
For details, see here
§Arguments:
line0- the “0“th line of the TLE, which sometimes contains the satellite nameline1- the 1st line of TLEline2- the 2nd line of the TLE
§Returns:
- A TLE object or string indicating error condition
§Example
use satkit::TLE;
let line0: &str = "0 INTELSAT 902";
let line1: &str = "1 26900U 01039A 06106.74503247 .00000045 00000-0 10000-3 0 8290";
let line2: &str = "2 26900 0.0164 266.5378 0003319 86.1794 182.2590 1.00273847 16981 9300.";
let tle = TLE::load_3line(&line0.to_string(),
&line1.to_string(),
&line2.to_string()
).unwrap();
sourcepub fn load_2line(line1: &str, line2: &str) -> Result<TLE, String>
pub fn load_2line(line1: &str, line2: &str) -> Result<TLE, String>
Load 2 lines as strings into a structure representing a Two-Line Element Set (TLE)
The TLE can then be used to compute satellite position and velocity as a function of time.
For details, see here
§Arguments:
line1- the 1st line of TLEline2- the 2nd line of the TLE
§Returns:
- A TLE object or string indicating error condition
§Example
use satkit::TLE;
let line1: &str = "1 26900U 01039A 06106.74503247 .00000045 00000-0 10000-3 0 8290";
let line2: &str = "2 26900 0.0164 266.5378 0003319 86.1794 182.2590 1.00273847 16981 9300.";
let tle = TLE::load_2line(
&line1.to_string(),
&line2.to_string()
).unwrap();
pub fn to_pretty_string(&self) -> String
Trait Implementations§
source§impl PartialOrd for TLE
impl PartialOrd for TLE
1.0.0 · source§fn le(&self, other: &Rhs) -> bool
fn le(&self, other: &Rhs) -> bool
self and other) and is used by the <=
operator. Read moreimpl StructuralPartialEq for TLE
Auto Trait Implementations§
impl Freeze for TLE
impl RefUnwindSafe for TLE
impl Send for TLE
impl Sync for TLE
impl Unpin for TLE
impl UnwindSafe for TLE
Blanket Implementations§
source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
source§impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
source§fn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
self from the equivalent element of its
superset. Read moresource§fn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
self is actually part of its subset T (and can be converted to it).source§fn to_subset_unchecked(&self) -> SS
fn to_subset_unchecked(&self) -> SS
self.to_subset but without any property checks. Always succeeds.source§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
self to the equivalent element of its superset.