Module orbdetpy.conversion

Functions

def azel_to_radec(time: float, az: float, el: float, lat: float, lon: float, alt: float, frame: int) ‑> Tuple[float, float]

Convert Azimuth/Elevation to Right Ascension/Declination.

Parameters

time : Offset in TT from J2000 epoch [s]. az : Azimuth [rad]. el : Elevation [rad]. lat : Observer WGS-84 latitude [rad]. lon : Observer WGS-84 longitude [rad]. alt : Observer height above WGS-84 reference ellipsoid [m]. frame : Destination reference frame; Frame.GCRF or Frame.EME2000.

Returns

Right Ascension and Declination.

def elem_to_pv(frame: int, time: float, sma: float, ecc: float, inc: float, raan: float, argp: float, anom: float, anom_type: int) ‑> List[float]

Convert osculating orbital elements to Cartesian state vector.

Parameters

frame : Inertial reference frame; a constant from Frame. time : Offset in TT from J2000 epoch [s]. Give a list for bulk conversions. sma : Semi-major axis [m]. Give a list for bulk conversions. ecc : Eccentricity. Give a list for bulk conversions. inc : Inclination [rad]. Give a list for bulk conversions. raan : RA of ascending node [rad]. Give a list for bulk conversions. argp : Argument of perigee [rad]. Give a list for bulk conversions. anom : Anomaly angle [rad]. Give a list for bulk conversions. anom_type : Anomaly angle type; a constant from PositionAngle. Give a list for bulk conversions.

Returns

Cartesian state vector.

def get_J2000_epoch_offset(utc: str) ‑> float

Get TT offset corresponding to ISO-8601 formatted UTC string.

Parameters

utc : ISO-8601 formatted UTC string or list of strings.

Returns

Offset in TT from J2000 epoch [s] or list of offsets.

def get_UTC_string(j2000_offset: float, precision: int = 3) ‑> str

Get ISO-8601 formatted UTC string corresponding to TT offset.

Parameters

j2000_offset : Offset in TT from J2000 epoch [s] or list of offsets. precision : Desired precision. Defaults to 3, i.e. milliseconds.

Returns

ISO-8601 formatted UTC string or list of strings.

def get_epoch_difference(from_epoch: int, to_epoch: int) ‑> str

Get the constant time difference between two epochs.

Parameters

from_epoch : From epoch; a value from the Epoch enumeration. to_epoch : To epoch; a value from the Epoch enumeration.

Returns

to_epoch - from_epoch in seconds.

def get_lvlh_rotation(state: List[float]) ‑> 

Construct the inertial->LVLH rotation matrix.

Parameters

state : Inertial state vector.

Returns

Inertial->LVLH frame rotation matrix.

def lla_to_pos(time: float, lla: List[float]) ‑> List[float]

Convert WGS-84 lat/lon/alt to Cartesian position.

Parameters

time : Offset in TT from J2000 epoch [s]. Give a list for bulk conversions. lla : WGS-84 latitude, longitude, altitude. Provide a list of lists for bulk conversions.

Returns

ITRF geocentric Cartesian position [m].

def ltr_to_matrix(lower_triangle: List[float]) ‑> List[List[float]]

Construct a symmetric matrix from its lower triangle.

Parameters

lower_triangle : Lower triangle of a symmetric matrix in row-major order.

Returns

Symmetric matrix as a list of lists or None on error.

def pos_to_lla(frame: int, time: float, pos: List[float]) ‑> List[float]

Convert a Cartesian position to WGS-84 lat/lon/alt.

Parameters

frame : Position reference frame; a constant from Frame. time : Offset in TT from J2000 epoch [s]. Give a list for bulk conversions. pos : Position vector [m]. Provide a list of lists for bulk conversions.

Returns

WGS-84 latitude [rad], longitude [rad], altitude [m].

def pv_to_elem(frame: int, time: float, pv: List[float]) ‑> List[float]

Convert Cartesian state vector to osculating orbital elements.

Parameters

frame : Inertial reference frame; a constant from Frame. time : Offset in TT from J2000 epoch [s]. Give a list for bulk conversions. pv : Inertial Cartesian state vector. Provide a list of lists for bulk conversions.

Returns

SMA, eccentricity, inclination, RAAN, perigee argument, mean anomaly, true anomaly, eccentric anomaly
 
def radec_to_azel(frame: int, time: float, ra: float, dec: float, lat: float, lon: float, alt: float) ‑> Tuple[float, float]

Convert Right Ascension/Declination to Azimuth/Elevation.

Parameters

frame : Source reference frame; Frame.GCRF or Frame.EME2000. time : Offset in TT from J2000 epoch [s]. ra : Right Ascension [rad]. dec : Declination [rad]. lat : Observer WGS-84 latitude [rad]. lon : Observer WGS-84 longitude [rad]. alt : Observer height above WGS-84 reference ellipsoid [m].

Returns

Azimuth and Elevation.

def transform_frame(src_frame: int, time: float, pva: List[float], dest_frame: int) ‑> List[float]

Transform a state vector from one frame to another.

Parameters

src_frame : Source reference frame; a constant from Frame.
time : Offset in TT from J2000 epoch [s]. Give a list for bulk transforms.
pva : State vector to transform, can be pos or pos+vel or pos+vel+acc. Provide a
list of lists for bulk frame transforms.

dest_frame : Destination reference frame; a constant from Frame.

Returns

State vector transformed to the destination frame.

Classes

class AnglesInput (*args, **kwargs)

A ProtocolMessage

Ancestors

  • google.protobuf.pyext._message.CMessage
  • google.protobuf.message.Message

Class variables

var DESCRIPTOR

Instance variables

var altitude

Field AnglesInput.altitude

var angle1

Field AnglesInput.angle1

var angle2

Field AnglesInput.angle2

var frame

Field AnglesInput.frame

var latitude

Field AnglesInput.latitude

var longitude

Field AnglesInput.longitude

var time

Field AnglesInput.time

class DoubleArray (*args, **kwargs)

A ProtocolMessage

Ancestors

  • google.protobuf.pyext._message.CMessage
  • google.protobuf.message.Message

Class variables

var DESCRIPTOR

Instance variables

var array

Field DoubleArray.array

class IntegerArray (*args, **kwargs)

A ProtocolMessage

Ancestors

  • google.protobuf.pyext._message.CMessage
  • google.protobuf.message.Message

Class variables

var DESCRIPTOR

Instance variables

var array

Field IntegerArray.array

class TransformFrameInput (*args, **kwargs)

A ProtocolMessage

Ancestors

  • google.protobuf.pyext._message.CMessage
  • google.protobuf.message.Message

Class variables

var DESCRIPTOR

Instance variables

var UTC_time

Field TransformFrameInput.UTC_time

var dest_frame

Field TransformFrameInput.dest_frame

var pva

Field TransformFrameInput.pva

var src_frame

Field TransformFrameInput.src_frame

var time

Field TransformFrameInput.time