MiriadeClass¶
- class astroquery.solarsystem.MiriadeClass[source]¶
Bases:
BaseQueryA class for querying the IMCCE/Miriade service.
Attributes Summary
URI used in query to service.
Methods Summary
get_ephemerides(*args, **kwargs)Queries the service and returns a table object.
get_ephemerides_async(targetname, *[, ...])Query the IMCCE Miriade ephemcc service.
Attributes Documentation
- uri¶
URI used in query to service.
Methods Documentation
- get_ephemerides(*args, **kwargs)¶
Queries the service and returns a table object.
Query the IMCCE Miriade ephemcc service.
- Parameters:
- targetnamestr
Name of the target to be queried.
- objtypestr, optional
Type of the object to be queried. Available are:
'asteroid','comet','dwarf planet','planet','satellite'. Default:'asteroid'- epoch
Timeobject, float, str,``None``, optional Start epoch of the query. If a float is provided, it is expected to be a Julian Date; if a str is provided, it is expected to be an iso date of the form
'YYYY-MM-DD HH-MM-SS'. IfNoneis provided, the current date and time are used as epoch. Default:None- epoch_stepstr, optional
Step size for ephemerides calculation. Must consist of a decimal number followed by a single character: (d)ays, (h)ours, (m)inutes or (s)econds. Default:
'1d'- epoch_nstepsint, optional
Number of increments of
epoch_stepstarting fromepochfor which ephemerides are calculated. Maximum number of steps is 5000. Default: 1- locationstr, optional
Location of the observer on Earth as a code or a set of coordinates. See the Miriade manual for details. Default: geocentric location (
'500')- coordtypeint, optional
Type of coordinates to be calculated:
1: spherical,2: rectangular,3: local coordinates (azimuth and elevation),4: hour angle coordinates,5: dedicated to observation,6: dedicated to AO observation. Default:1- timescalestr, optional
The time scale used in the computation of the ephemerides:
'UTC'or'TT'. Default:'UTC'- planetary_theorystr, optional
Planetary ephemerides set to be utilized in the calculations:
'INPOP','DE405','DE406'. Default:'INPOP'- ephtypeint, optional
Type of ephemerides to be calculated:
1: astrometric J2000,2: apparent of the date,3: mean of the date,4: mean J2000, Default:1- refplanestr, optional
Reference plane:
'equator'or'ecliptic'. Default:'equator'- elementsstr, optional
Set of osculating elements to be used in the calculations:
'ASTORB'or'MPCORB'. Default:'ASTORB'- radial_velocitybool, optional
Calculate additional information on the target’s radial velocity. Default:
False- get_query_payloadbool, optional
When set to
Truethe method returns the HTTP request parameters as a dict, default:False- get_raw_responsebool, optional
Return raw data as obtained by Miriade without parsing the data into a table, default:
False- cachebool, optional
If
Truethe query will be cached. Default:True
- Returns:
- tableA
Tableobject.
- tableA
Notes
The following parameters can be queried using this function. Note that different
coordtypesetting provide different sets of parameters; number in parentheses denote whichcoordtypesettings include the parameters.Column Name
Definition
targetTarget name (str, 1, 2, 3, 4, 5, 6 )
epochEphemerides epoch (JD, float, 1, 2, 3, 4, 5, 6)
RATarget RA at
ephtype(deg, float, 1)DECTarget declination at
ephtype(deg, float, 1, 4, 5)RAJ2000Target RA at J2000 (deg, float, 5, 6)
DECJ2000Target declination at J2000 (deg, float, 5, 6)
AZTarget azimuth (deg, float, 3, 5)
ELTarget elevation (deg, float, 3, 5)
deltaDistance from observer (au, float, 1, 2, 3, 4, 5, 6)
delta_rateRate in observer distance (km/s, float, 1, 5, 6)
VApparent visual magnitude (mag, float, 1, 2, 3, 4, 5, 6)
alphaSolar phase angle (deg, 1, 2, 3, 4, 5, 6)
elongSolar elongation angle (deg, 1, 2, 3, 4, 5, 6)
RAcosD_rateRate of motion in RA * cos(DEC) (arcsec/min, float, 1, 5, 6)
DEC_rateRate of motion in DEC (arcsec/min, float, 1, 5, 6)
xX position state vector (au, float, 2)
yY position state vector (au, float, 2)
zZ position state vector (au, float, 2)
vxX velocity state vector (au/d, float, 2)
vyY velocity state vector (au/d, float, 2)
vzZ velocity state vector (au/d, float, 2)
rvRadial velocity (km/s, float, 2)
heldistTarget heliocentric distance (au, float, 2, 5, 6)
x_hX heliocentric position vector (au, float, 2)
y_hY heliocentric position vector (au, float, 2)
z_hZ heliocentric position vector (au, float, 2)
vx_hX heliocentric vel. vector (au/d, float, 2)
vy_hY heliocentric vel. vector (au/d, float, 2)
vz_hZ heliocentric vel. vector (au/d, float, 2)
hourangleTarget hour angle (deg, float, 4, 5)
siderealtimeLocal sidereal time (hr, float, 5, 6)
refractionAtmospheric refraction (arcsec, float, 5, 6)
airmassTarget airmass (float, 5, 6)
posuncPositional uncertainty (arcsec, float, 5, 6)
Examples
>>> from astroquery.imcce import Miriade >>> from astropy.time import Time >>> epoch = Time('2019-01-01', format='iso') >>> Miriade.get_ephemerides('3552', epoch=epoch) <Table masked=True length=1> target epoch RA ... DEC_rate delta_rate d deg ... arcs / min km / s bytes20 float64 float64 ... float64 float64 ----------- -------------------- ------------------ ... ---------- ------------ Don Quixote 2458484.5 16.105294999999998 ... -0.25244 31.4752734
- get_ephemerides_async(targetname, *, objtype='asteroid', epoch=None, epoch_step='1d', epoch_nsteps=1, location=500, coordtype=1, timescale='UTC', planetary_theory='INPOP', ephtype=1, refplane='equator', elements='ASTORB', radial_velocity=False, get_query_payload=False, get_raw_response=False, cache=True)[source]¶
Query the IMCCE Miriade ephemcc service.
- Parameters:
- targetnamestr
Name of the target to be queried.
- objtypestr, optional
Type of the object to be queried. Available are:
'asteroid','comet','dwarf planet','planet','satellite'. Default:'asteroid'- epoch
Timeobject, float, str,``None``, optional Start epoch of the query. If a float is provided, it is expected to be a Julian Date; if a str is provided, it is expected to be an iso date of the form
'YYYY-MM-DD HH-MM-SS'. IfNoneis provided, the current date and time are used as epoch. Default:None- epoch_stepstr, optional
Step size for ephemerides calculation. Must consist of a decimal number followed by a single character: (d)ays, (h)ours, (m)inutes or (s)econds. Default:
'1d'- epoch_nstepsint, optional
Number of increments of
epoch_stepstarting fromepochfor which ephemerides are calculated. Maximum number of steps is 5000. Default: 1- locationstr, optional
Location of the observer on Earth as a code or a set of coordinates. See the Miriade manual for details. Default: geocentric location (
'500')- coordtypeint, optional
Type of coordinates to be calculated:
1: spherical,2: rectangular,3: local coordinates (azimuth and elevation),4: hour angle coordinates,5: dedicated to observation,6: dedicated to AO observation. Default:1- timescalestr, optional
The time scale used in the computation of the ephemerides:
'UTC'or'TT'. Default:'UTC'- planetary_theorystr, optional
Planetary ephemerides set to be utilized in the calculations:
'INPOP','DE405','DE406'. Default:'INPOP'- ephtypeint, optional
Type of ephemerides to be calculated:
1: astrometric J2000,2: apparent of the date,3: mean of the date,4: mean J2000, Default:1- refplanestr, optional
Reference plane:
'equator'or'ecliptic'. Default:'equator'- elementsstr, optional
Set of osculating elements to be used in the calculations:
'ASTORB'or'MPCORB'. Default:'ASTORB'- radial_velocitybool, optional
Calculate additional information on the target’s radial velocity. Default:
False- get_query_payloadbool, optional
When set to
Truethe method returns the HTTP request parameters as a dict, default:False- get_raw_responsebool, optional
Return raw data as obtained by Miriade without parsing the data into a table, default:
False- cachebool, optional
If
Truethe query will be cached. Default:True
Notes
The following parameters can be queried using this function. Note that different
coordtypesetting provide different sets of parameters; number in parentheses denote whichcoordtypesettings include the parameters.Column Name
Definition
targetTarget name (str, 1, 2, 3, 4, 5, 6 )
epochEphemerides epoch (JD, float, 1, 2, 3, 4, 5, 6)
RATarget RA at
ephtype(deg, float, 1)DECTarget declination at
ephtype(deg, float, 1, 4, 5)RAJ2000Target RA at J2000 (deg, float, 5, 6)
DECJ2000Target declination at J2000 (deg, float, 5, 6)
AZTarget azimuth (deg, float, 3, 5)
ELTarget elevation (deg, float, 3, 5)
deltaDistance from observer (au, float, 1, 2, 3, 4, 5, 6)
delta_rateRate in observer distance (km/s, float, 1, 5, 6)
VApparent visual magnitude (mag, float, 1, 2, 3, 4, 5, 6)
alphaSolar phase angle (deg, 1, 2, 3, 4, 5, 6)
elongSolar elongation angle (deg, 1, 2, 3, 4, 5, 6)
RAcosD_rateRate of motion in RA * cos(DEC) (arcsec/min, float, 1, 5, 6)
DEC_rateRate of motion in DEC (arcsec/min, float, 1, 5, 6)
xX position state vector (au, float, 2)
yY position state vector (au, float, 2)
zZ position state vector (au, float, 2)
vxX velocity state vector (au/d, float, 2)
vyY velocity state vector (au/d, float, 2)
vzZ velocity state vector (au/d, float, 2)
rvRadial velocity (km/s, float, 2)
heldistTarget heliocentric distance (au, float, 2, 5, 6)
x_hX heliocentric position vector (au, float, 2)
y_hY heliocentric position vector (au, float, 2)
z_hZ heliocentric position vector (au, float, 2)
vx_hX heliocentric vel. vector (au/d, float, 2)
vy_hY heliocentric vel. vector (au/d, float, 2)
vz_hZ heliocentric vel. vector (au/d, float, 2)
hourangleTarget hour angle (deg, float, 4, 5)
siderealtimeLocal sidereal time (hr, float, 5, 6)
refractionAtmospheric refraction (arcsec, float, 5, 6)
airmassTarget airmass (float, 5, 6)
posuncPositional uncertainty (arcsec, float, 5, 6)
Examples
>>> from astroquery.imcce import Miriade >>> from astropy.time import Time >>> epoch = Time('2019-01-01', format='iso') >>> Miriade.get_ephemerides('3552', epoch=epoch) <Table masked=True length=1> target epoch RA ... DEC_rate delta_rate d deg ... arcs / min km / s bytes20 float64 float64 ... float64 float64 ----------- -------------------- ------------------ ... ---------- ------------ Don Quixote 2458484.5 16.105294999999998 ... -0.25244 31.4752734