Orbit¶
- class sbpy.data.Orbit[source]¶
Bases:
DataClassClass for querying, manipulating, integrating, and fitting orbital elements
Methods Summary
D_criterion(obj[, version])Evaluate orbit similarity D-criterion
from_horizons(targetids[, id_type, epochs, ...])Load target orbital elements from JPL Horizons using
astroquery.jplhorizons.HorizonsClass.elementsfrom_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
Orbitobject.oo_transform(orbittype[, ephfile])Uses pyoorb to transform this orbit object to a different orbit type definition.
tisserand([planet, epoch])Tisserand parameter with respect to a planet
Methods Documentation
- D_criterion(obj, version='sh')[source]¶
Evaluate orbit similarity D-criterion
Three different versions of D-criterion are defined and often compared to each other, includingthe Southworth-Hawkins function [SH63], Drummond function [D81], and the hybrid version [J93]. See review by [W19].
- Parameters:
- obj
Orbitobject Object(s) against which to calculate D-criterion
- version[‘sh’, ‘d’, ‘h’], optional
Select the versions of D-criterion formula. Case insensitive. ‘sh’ : Southworth-Hawkins function ‘d’ : Drummond function ‘h’ : Hybrid function See references for the details of each version.
- obj
- Returns:
- float or numpy.ndarray
References
Examples
>>> from sbpy.data import Orbit >>> comets = Orbit.from_horizons(['252P', 'P/2016 BA14'], ... id_type='designation', closest_apparition=True) ... >>> # Southworth & Hawkins function >>> D_SH = comets[0].D_criterion(comets[1]) >>> # Drummond function >>> D_D = comets[0].D_criterion(comets[1], version='d') >>> # hybrid function >>> D_H = comets[0].D_criterion(comets[1], version='h')
- 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
targetidsprovided; 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'- epochs
Timeor dict, optional Epochs of elements to be queried; requires a
Timeobject with a single or multiple epochs. A dictionary including keywordsstartandstop, as well as eitherstepornumber, can be used to generate a range of epochs.startandstophave to beTimeobjects (see Epochs and the use of astropy.time). Ifstepis provided, a range of epochs will be queries starting atstartand ending atstopin steps ofstep;stephas to be provided as aQuantityobject with integer value and a unit of either minutes, hours, days, or years. Ifnumberis provided as an integer, the interval defined bystartandstopis split intonumberequidistant intervals. IfNoneis provided, current date and time are used. All epochs should be provided in TDB; if not, they will be converted to TDB and aTimeScaleWarningwill 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:
Orbitobject
Notes
For detailed explanations of the queried fields, refer to
astroquery.jplhorizons.HorizonsClass.elementsand the JPL Horizons documentation.By default, elements are provided in the J2000.0 reference system and relative to the ecliptic and mean equinox of the reference epoch. Different settings can be chosen using additional keyword arguments as used by
astroquery.jplhorizons.HorizonsClass.elements.
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', orNoneto indicate type of identifiers provided intargetids. IfNone, automatic identification is attempted usingnames. Default:None- target_typestr, optional
'asteroid','comet', orNoneto indicate target type. IfNone, automatic identification is attempted usingnames. Default:None- **kwargsoptional
Additional keyword arguments are passed to
query_object
- Returns:
OrbitobjectThe 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
Orbitobject. 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 aueccentricity (
'e', for Keplerian or cometary orbit) or y-component of state vector ('y', for cartesian orbit) in auinclination (
'i', for Keplerian or cometary orbit) in degrees or z-component of state vector ('z', for cartesian orbit) in aulongitude of the ascending node (
'Omega', for Keplerian or cometary orbit) in degrees or x-component of velocity vector ('vx', for cartesian orbit), au/dayargument of the periapsis (
'w', for Keplerian or cometary orbit) in degrees or y-component of velocity vector ('vy', for cartesian orbit) in au/daymean 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/dayepoch (
'epoch') asTimeobjectabsolute magnitude (
'H') in magphotometric phase slope (
'G')
- Parameters:
- epochs
Timeobject Epoch to which the orbit will be propagated to. Must be a
Timeobject holding a single epoch or multiple epochs (see Epochs and the use of astropy.time). The resultingOrbitobject will have the same time scale asepochs.- 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'
- epochs
- Returns:
Orbitobject
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 aueccentricity (
'e', for Keplerian or cometary orbit) or y-component of state vector ('y', for cartesian orbit) in auinclination (
'i', for Keplerian or cometary orbit) in degrees or z-component of state vector ('z', for cartesian orbit) in aulongitude of the ascending node (
'Omega', for Keplerian or cometary orbit) in degrees or x-component of velocity vector ('vx', for cartesian orbit), au/dayargument of the periapsis (
'w', for Keplerian or cometary orbit) in degrees or y-component of velocity vector ('vy', for cartesian orbit) in au/daymean 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/dayepoch (
'epoch') asTimeobjectabsolute magnitude (
'H') in magphotometric 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:
Orbitobject
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
- tisserand(planet='599', epoch=None)[source]¶
Tisserand parameter with respect to a planet
- Parameters:
- planetstr,
Orbitobject, or sequence of str, optional Planet(s) against which the Tisserand parameter is calculated. If
str, then the orbital elements of the planet will be automaticallyl pulled from JPL Horizons. Ifselfand/orplanetcontains more than one object, thennumpybroadcasting rules apply. Default is Jupiter.- epoch
Time, optional The epoch of planet orbit if pulled from JPL Horizons. This parameter will be passed to
from_horizons. If the planet orbit is passed directly, then this parameter has no effect.
- planetstr,
- Returns:
- u.Quantity
The Tisserand parameter(s)
References
Levison, H. F., Duncan, M. J. 1997, Icarus 127, 13
Examples
>>> from sbpy.data import Orbit >>> comets = Orbit.from_horizons(['252P', 'P/2016 BA14'], ... id_type='designation', closest_apparition=True) ... >>> T_J = comets.tisserand() >>> >>> halley = Orbit.from_horizons('1P', id_type='designation', ... closest_apparition=True) >>> T = halley.tisserand(['599', '699', '799', '899'])