Ephem¶
- class sbpy.data.Ephem[source]¶
Bases:
DataClassClass for querying, manipulating, and calculating ephemerides
Methods Summary
from_horizons(targetids[, id_type, epochs, ...])Load target ephemerides from JPL Horizons using
astroquery.jplhorizons.HorizonsClass.ephemeridesfrom_miriade(targetids[, objtype, epochs, ...])Load target ephemerides from IMCCE Miriade using
astroquery.imcce.MiriadeClass.get_ephemeridesfrom_mpc(targetids[, epochs, location, ...])Load ephemerides from the Minor Planet Center.
from_oo(orbit[, epochs, location, scope, ...])Uses pyoorb to derive ephemerides from an
Orbitobject.Methods Documentation
- classmethod from_horizons(targetids, id_type='smallbody', epochs=None, location='500', **kwargs)[source]¶
Load target ephemerides from JPL Horizons using
astroquery.jplhorizons.HorizonsClass.ephemerides- 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
Timeobject, or dictionary, optional Epochs of elements to be queried;
Timeobjects support single and 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 UTC; if not, they will be converted to UTC and aTimeScaleWarningwill be raised. Default:None- locationstr or
EarthLocation, optional Location of the observer using IAU observatory codes (see IAU observatory codes) or as
EarthLocation. Default:'500'(geocentric)- **kwargsoptional
Arguments that will be provided to
astroquery.jplhorizons.HorizonsClass.ephemerides.
- Returns:
EphemobjectThe resulting object will be populated with columns as defined in
ephemerides; refer to that document on information on how to modify the list of queried parameters.
Notes
For detailed explanations of the queried fields, refer to
astroquery.jplhorizons.HorizonsClass.ephemeridesand the JPL Horizons documentation.By default, all properties are provided in the J2000.0 reference system. Different settings can be chosen using additional keyword arguments as used by
astroquery.jplhorizons.HorizonsClass.ephemerides.
Examples
>>> from sbpy.data import Ephem >>> from astropy.time import Time >>> epoch = Time('2018-05-14', scale='utc') >>> eph = Ephem.from_horizons('ceres', epochs=epoch) ...
- classmethod from_miriade(targetids, objtype='asteroid', epochs=None, location='500', **kwargs)[source]¶
Load target ephemerides from IMCCE Miriade using
astroquery.imcce.MiriadeClass.get_ephemerides- Parameters:
- targetidsstr or iterable of str
Target identifier, i.e., a number, name, designation, or JPL Horizons record number, for one or more targets.
- objtypestr, optional
The nature of
targetidsprovided; possible values are'asteroid','comet','dwarf planet','planet', or'satellite'. Default:'asteroid'- epochs
Timeobject, or dictionary, optional Epochs of elements to be queried;
Timeobjects support single and 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 queried starting atstartand ending atstopin steps ofstep;stephas to be provided as aQuantityobject with integer value and a unit of either seconds, minutes, hours, or days. 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 UTC; if not, they will be converted to UTC and aTimeScaleWarningwill be raised. Default:None- locationstr or
EarthLocation, optional Location of the observer using IAU observatory codes (see IAU observatory codes) or as
EarthLocation. Default:'500'(geocentric)- **kwargsoptional
Arguments that will be provided to
astroquery.imcce.MiriadeClass.get_ephemerides.
- Returns:
EphemobjectThe resulting object will be populated with columns as defined in
get_ephemerides; refer to that document on information on how to modify the list of queried parameters.
Notes
For detailed explanations of the queried fields, refer to
astroquery.imcce.MiriadeClass.get_ephemeridesand the Miriade documentation.By default, all properties are provided in the J2000.0 reference system. Different settings can be chosen using additional keyword arguments as used by
astroquery.imcce.MiriadeClass.get_ephemerides.
Examples
>>> from sbpy.data import Ephem >>> from astropy.time import Time >>> epoch = Time('2018-05-14', scale='utc') >>> eph = Ephem.from_horizons('ceres', epochs=epoch) ...
- classmethod from_mpc(targetids, epochs=None, location='500', ra_format=None, dec_format=None, **kwargs)[source]¶
Load ephemerides from the Minor Planet Center.
- Parameters:
- targetidsstr or iterable of str
Target identifier, resolvable by the Minor Planet Ephemeris Service [MPES], e.g., 2P, C/1995 O1, P/Encke, (1), 3200, Ceres, and packed designations, for one or more targets.
- epochs
Timeobject, or dictionary, optional Request ephemerides at these epochs. May be a single epoch or multiple epochs as
Time(see Epochs and the use of astropy.time)or a dictionary describing a linearly-spaced array of epochs. All epochs should be provided in UTC; if not, they will be converted to UTC and aTimeScaleWarningwill be raised. IfNone(default), the current date and time will be used.For the dictionary format, the keys
start(start epoch),step(step size),stop(end epoch), and/ornumber(number of epochs total) are used. Only one ofstopandnumbermay be specified at a time.step,stop, andnumberare optional. The values of ofstartandstopmust beTimeobjects.numbershould be an integer value;stepshould be aQuantitywith an integer value and units of seconds, minutes, hours, or days.All epochs should be provided in UTC; if not, they will be converted to UTC and a
TimeScaleWarningwill be raised.- locationvarious, optional
Location of the observer as an IAU observatory code [OBSCODES] (string), a 3-element array of Earth longitude, latitude, altitude, or an
EarthLocation. Longitude and latitude should be parseable byAngle, and altitude should be parsable byQuantity(with units of length). IfNone, then the geocenter (code 500) is used.- ra_formatdict, optional
Format the RA column with
to_stringusing these keyword arguments, e.g.,{'sep': ':', 'unit': 'hourangle', 'precision': 1}.- dec_formatdict, optional
Format the Dec column with
to_stringusing these keyword arguments, e.g.,{'sep': ':', 'precision': 0}.- **kwargs
Additional keyword arguments are passed to
get_ephemerides:eph_type,proper_motion,proper_motion_unit,suppress_daytime,suppress_set,perturbed,unc_links,cache.
- Returns:
EphemobjectThe resulting object will be populated with columns as defined in
get_ephemerides; refer to that document on information on how to modify the list of queried parameters.
Notes
All properties are provided in the J2000.0 reference system.
See
astroquery.mpc.MPC.get_ephemeridesand the Minor Planet Ephemeris Service user’s guide [MPES] for details, including acceptable target names.
References
[MPES] (1,2)Wiliams, G. The Minor Planet Ephemeris Service. https://minorplanetcenter.net/iau/info/MPES.pdf
[OBSCODES]IAU Minor Planet Center. List of observatory codes. https://minorplanetcenter.net/iau/lists/ObsCodesF.html
Examples
Query a single set of ephemerides of Ceres as observed from Maunakea: >>> from sbpy.data import Ephem >>> from astropy.time import Time >>> epoch = Time(‘2018-05-14’, scale=’utc’) >>> eph = Ephem.from_mpc(‘ceres’, epoch, 568) # doctest: +REMOTE_DATA
Query a range of ephemerides of comet 2P/Encke as observed from Maunakea: >>> epochs = {‘start’: Time(‘2019-01-01’), … ‘step’: 1*u.d, ‘number’: 365} >>> eph = Ephem.from_mpc(‘2P’, epochs, 568) # doctest: +REMOTE_DATA
- classmethod from_oo(orbit, epochs=None, location='500', scope='full', dynmodel='N', ephfile='de430')[source]¶
Uses pyoorb to derive ephemerides from an
Orbitobject. For a list of output parameters, please read the pyoorb documentation.- Parameters:
- orbit
Orbitobject Can contain any number of orbits, ephemerides will be calculated for each orbit. 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') asTimeabsolute magnitude (
'H') in magphotometric phase slope (
'G')
- epochs
Timeobject, optional Epochs of elements to be queried; must be a
Timeobject holding a single or multiple epochs (see Epochs and the use of astropy.time). IfNoneis provided, current date and time are used. The same time scale that is used inepochswill be applied to the results. Default:None- locationstr, optional, default
'500'(geocentric) Location of the observer.
- scopestr
Scope of data to be determined:
'full'obtains all available properties,'basic'obtains only a limited amount of data. Default:'full'- 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'
- orbit
- Returns:
Ephemobject
Examples
Compute ephemerides for Ceres as seen from the Discovery Channel Telescope for the next 10 days at 1hr intervals:
>>> import numpy as np >>> from sbpy.data import Orbit, Ephem >>> from astropy.time import Time >>> epochs = Time(Time.now().jd + np.arange(0, 10, 1/24), format='jd') >>> ceres = Orbit.from_horizons('1') >>> eph = Ephem.from_oo(ceres, epochs, 'G37') >>> eph <QTable length=240> targetname epoch ... obsz trueanom d ... AU deg str7 float64 ... float64 float64 ---------- ------------------ ... ----------------------- ----------------- 1 Ceres 2458519.316966272 ... 3.2083678848104924e-06 68.0863831954328 1 Ceres 2458519.3586329385 ... 2.7022422510736277e-07 68.09589266358881 1 Ceres 2458519.4002996054 ... -3.111046209036683e-06 68.10540191585879 1 Ceres 2458519.441966272 ... -6.700369254264427e-06 68.11491095202307 1 Ceres 2458519.4836329385 ... -1.0248419404668141e-05 68.12441977218093 1 Ceres 2458519.5252996054 ... -1.3508703580356052e-05 68.13392837643161 ... ... ... ... ... 1 Ceres 2458529.066966272 ... 1.2522500440509399e-05 70.30569661787204 1 Ceres 2458529.1086329385 ... 1.4101698473351076e-05 70.31515536712485 1 Ceres 2458529.1502996054 ... 1.4771304981564537e-05 70.3246138990413 1 Ceres 2458529.191966272 ... 1.448582020449618e-05 70.33407221340468 1 Ceres 2458529.2336329385 ... 1.326517587380005e-05 70.34353031031534 1 Ceres 2458529.2752996054 ... 1.1193369555934085e-05 70.35298818987367