Obs¶
- class sbpy.data.Obs[source]¶
Bases:
Ephem
Class for querying, storing, and manipulating observations
Methods Summary
from_mpc
(targetid[, id_type])Load available observations for a target from the Minor Planet Center using
get_observations
.supplement
([service, id_field, epoch_field, ...])Supplement observational data with ephemerides queried from the selected service.
Methods Documentation
- classmethod from_mpc(targetid, id_type=None, **kwargs)[source]¶
Load available observations for a target from the Minor Planet Center using
get_observations
.- Parameters:
- targetidstr
Target identifier, resolvable by the Minor Planet Ephemeris Service.
- id_typestr, optional
'asteroid number'
,'asteroid designation'
,'comet number'
,'comet designation'
, orNone
.None
attempts to automatically identify the target id type usingNames
. Default:None
- **kwargsoptional
Additional keyword arguments are passed to
query_object
- Returns:
Obs
objectThe resulting object will be populated with the same fields defined in
get_observations
.
Examples
>>> from sbpy.data import Obs >>> obs = Obs.from_mpc('12893') >>> orb[:3] number desig discovery note1 ... DEC mag band observatory ... deg mag ------ --------- --------- ----- ... ------------------- --- ---- ----------- 12893 1998 QS55 -- -- ... -15.78888888888889 0.0 -- 413 12893 1998 QS55 -- -- ... -15.788944444444445 0.0 -- 413 12893 1998 QS55 * 4 ... 5.526472222222222 0.0 -- 809
- supplement(service='jplhorizons', id_field='targetname', epoch_field='epoch', location='500', modify_fieldnames='obs', **kwargs)[source]¶
Supplement observational data with ephemerides queried from the selected service.
- Parameters:
- servicestr, optional
Service from which to acquire data:
'jplhorizons'
,'mpc'
, or'miriade'
, corresponding to the JPL Horizons system (usingfrom_horizons
), the Minor Planet Center ephemeris service (usingfrom_mpc
), and the IMCCE Miriade service (usingfrom_miriade
). Default:'jplhorizons'
- id_fieldstr, optional
Field name that corresponds to a suitable target identifier in this
Obs
object. Default:'targetname'
- epoch_fieldstr, optional
Field name that corresponds to a suitable epoch identifier in this
Obs
object. The corresponding column must be of typeTime
. Default:'epoch'
- locationstr, optional
Location of the observer for the data stored in this
Obs
object. Default:'500'
(geocenter)- modify_fieldnamesstr, optional
Defines whether field names in this
Obs
object ('obs'
) or in the supplemental data to be queried ('eph'
) will be modified by adding a suffix in case of field name collisions. Default:'obs'
- **kwargsoptional
Additional keyword arguments are passed to the corresponding ephemerides query service.
- Returns:
Notes
Not all available service are equally suited for this kind of query: only the JPL Horizons system enables quick queries for a large number of epochs. Queries using the other services may take a long time depending on the number of epochs and targets.
Examples
>>> from sbpy.data import Obs >>> obs = Obs.from_mpc('2019 AA', id_type='asteroid designation') >>> data = obs.supplement(id_field='designation') >>> data.field_names <TableColumns names=('number','desig','discovery','note1','note2','epoch','RA_obs','DEC_obs','mag','band','observatory','target','RA','DEC','delta','V','alpha','elong','RAcosD_rate','DEC_rate','delta_rate')>