Orbit

class sbpy.data.Orbit[source]

Bases: sbpy.data.core.DataClass

Class for querying, manipulating, integrating, and fitting orbital elements

Methods Summary

from_horizons(targetids[, id_type, epochs, ...])

Load target orbital elements from JPL Horizons using astroquery.jplhorizons.HorizonsClass.elements

from_mpc(targetids[, id_type, target_type])

Load latest orbital elements from the Minor Planet Center.

oo_propagate(epochs[, dynmodel, ephfile])

Uses pyoorb to propagate this Orbit object.

oo_transform(orbittype[, ephfile])

Uses pyoorb to transform this orbit object to a different orbit type definition.

Methods Documentation

classmethod from_horizons(targetids, id_type='smallbody', epochs=None, center='500@10', **kwargs)[source]

Load target orbital elements from JPL Horizons using astroquery.jplhorizons.HorizonsClass.elements

Parameters
targetidsstr or iterable of str

Target identifier, i.e., a number, name, designation, or JPL Horizons record number, for one or more targets.

id_typestr, optional

The nature of targetids provided; possible values are 'smallbody' (asteroid or comet), 'majorbody' (planet or satellite), 'designation' (asteroid or comet designation), 'name' (asteroid or comet name), 'asteroid_name', 'comet_name', 'id' (Horizons id). Default: 'smallbody'

epochsTime or dict, optional

Epochs of elements to be queried; requires a Time object with a single or multiple epochs. A dictionary including keywords start and stop, as well as either step or number, can be used to generate a range of epochs. start and stop have to be Time objects (see Epochs and the use of astropy.time). If step is provided, a range of epochs will be queries starting at start and ending at stop in steps of step; step has to be provided as a Quantity object with integer value and a unit of either minutes, hours, days, or years. If number is provided as an integer, the interval defined by start and stop is split into number equidistant intervals. If None is provided, current date and time are used. All epochs should be provided in TDB; if not, they will be converted to TDB and a TimeScaleWarning will be raised. Default: None

centerstr, optional, default '500@10' (center of the Sun)

Elements will be provided relative to this position.

**kwargsoptional

Arguments that will be provided to astroquery.jplhorizons.HorizonsClass.elements.

Returns
Orbit object

Notes

Examples

>>> from sbpy.data import Orbit
>>> from astropy.time import Time
>>> epoch = Time('2018-05-14', scale='tdb')  
>>> eph = Orbit.from_horizons('Ceres', epochs=epoch)  
classmethod from_mpc(targetids, id_type=None, target_type=None, **kwargs)[source]

Load latest orbital elements from the Minor Planet Center.

Parameters
targetidsstr or iterable of str

Target identifier, resolvable by the Minor Planet Ephemeris Service. If multiple targetids are provided in a list, the same format (number, name, or designation) must be used.

id_typestr, optional

'name', 'number', 'designation', or None to indicate type of identifiers provided in targetids. If None, automatic identification is attempted using names. Default: None

target_typestr, optional

'asteroid', 'comet', or None to indicate target type. If None, automatic identification is attempted using names. Default: None

**kwargsoptional

Additional keyword arguments are passed to query_object

Returns
Orbit object

The resulting object will be populated with most columns as defined in query_object; refer to that document on information on how to modify the list of queried parameters.

Examples

>>> from sbpy.data import Orbit
>>> orb = Orbit.from_mpc('Ceres') 
>>> orb  
<QTable length=1>
 absmag    Q      arc      w     ...     a        Tj   moid_uranus moid_venus
  mag      AU      d      deg    ...     AU                 AU         AU
float64 float64 float64 float64  ...  float64  float64   float64    float64
------- ------- ------- -------- ... --------- ------- ----------- ----------
   3.34    2.98 79653.0 73.59764 ... 2.7691652     3.3     15.6642    1.84632
oo_propagate(epochs, dynmodel='N', ephfile='de430')[source]

Uses pyoorb to propagate this Orbit object. Required fields are:

  • target identifier ('targetname')

  • semi-major axis ('a', for Keplerian orbit) or perihelion distance ('q', for cometary orbit), typically in au or or x-component of state vector ('x', for cartesian orbit), typically in au

  • eccentricity ('e', for Keplerian or cometary orbit) or y-component of state vector ('y', for cartesian orbit) in au

  • inclination ('i', for Keplerian or cometary orbit) in degrees or z-component of state vector ('z', for cartesian orbit) in au

  • longitude of the ascending node ('Omega', for Keplerian or cometary orbit) in degrees or x-component of velocity vector ('vx', for cartesian orbit), au/day

  • argument of the periapsis ('w', for Keplerian or cometary orbit) in degrees or y-component of velocity vector ('vy', for cartesian orbit) in au/day

  • mean anomaly ('M', for Keplerian orbits) in degrees or perihelion epoch ('Tp', for cometary orbits) as astropy Time object or z-component of velocity vector ('vz', for cartesian orbit) in au/day

  • epoch ('epoch') as Time object

  • absolute magnitude ('H') in mag

  • photometric phase slope ('G')

Parameters
epochsTime object

Epoch to which the orbit will be propagated to. Must be a Time object holding a single epoch or multiple epochs (see Epochs and the use of astropy.time). The resulting Orbit object will have the same time scale as epochs.

dynmodelstr, optional

The dynamical model to be used in the propagation: 'N' for n-body simulation or '2' for a 2-body simulation. Default: 'N'

ephfilestr, optional

Planet and Lunar ephemeris file version as provided by JPL to be used in the propagation. Default: 'de430'

Returns
Orbit object

Examples

Propagate the orbit of Ceres 100 days into the future:

>>> from sbpy.data import Orbit
>>> from astropy.time import Time
>>> epoch = Time(Time.now().jd + 100, format='jd')
>>> ceres = Orbit.from_horizons('Ceres')      
>>> future_ceres = ceres.oo_propagate(epoch)  
>>> print(future_ceres)  
<QTable length=1>
   id           a                  e          ...    H       G    timescale
                AU                            ...   mag
  str7       float64            float64       ... float64 float64    str3
------- ----------------- ------------------- ... ------- ------- ---------
1 Ceres 2.769331727251861 0.07605371361208543 ...    3.34    0.12       UTC
oo_transform(orbittype, ephfile='de430')[source]

Uses pyoorb to transform this orbit object to a different orbit type definition. Required fields are:

  • target identifier ('targetname')

  • semi-major axis ('a', for Keplerian orbit) or perihelion distance ('q', for cometary orbit), typically in au or or x-component of state vector ('x', for cartesian orbit), typically in au

  • eccentricity ('e', for Keplerian or cometary orbit) or y-component of state vector ('y', for cartesian orbit) in au

  • inclination ('i', for Keplerian or cometary orbit) in degrees or z-component of state vector ('z', for cartesian orbit) in au

  • longitude of the ascending node ('Omega', for Keplerian or cometary orbit) in degrees or x-component of velocity vector ('vx', for cartesian orbit), au/day

  • argument of the periapsis ('w', for Keplerian or cometary orbit) in degrees or y-component of velocity vector ('vy', for cartesian orbit) in au/day

  • mean anomaly ('M', for Keplerian orbits) in degrees or perihelion epoch ('Tp', for cometary orbits) as astropy Time object or z-component of velocity vector ('vz', for cartesian orbit) in au/day

  • epoch ('epoch') as Time object

  • absolute magnitude ('H') in mag

  • photometric phase slope ('G')

Parameters
orbittypestr

Orbit definition to be transformed to; available orbit definitions are KEP (Keplerian elements), CART (cartesian elements), COM (cometary elements).

ephfilestr, optional

Planet and Lunar ephemeris file version as provided by JPL to be used in the propagation. Default: 'de430'

Returns
Orbit object

Examples

Obtain the current state vector (cartesian definition, CART) for asteroid Ceres.

>>> from sbpy.data import Orbit
>>> ceres = Orbit.from_horizons('Ceres')  
>>> statevec = ceres.oo_transform('CART') 
>>> print(statevec)  
<QTable length=1>
   id           x                   y          ...    H       G    timescale
                AU                  AU         ...   mag
  str7       float64             float64       ... float64 float64    str2
------- ------------------ ------------------- ... ------- ------- ---------
1 Ceres -1.967176101061908 -1.7891189971612211 ...    3.34    0.12        TT