def from_mpc(cls, targetids, epochs=None, location='500', **kwargs): """Load ephemerides from the `Minor Planet Center <https://minorplanetcenter.net>`_. Parameters ---------- targetids : str 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 : `~astropy.time.Time` object, or dictionary, optional Request ephemerides at these epochs. May be a single epoch or multiple epochs as `~astropy.time.Time` (see :ref:`epochs`)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 a `~sbpy.data.TimeScaleWarning` will be raised. If ``None`` (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/or ``number`` (number of epochs total) are used. Only one of ``stop`` and ``number`` may be specified at a time. ``step``, ``stop``, and ``number`` are optional. The values of of ``start`` and ``stop`` must be `~astropy.time.Time` objects. ``number`` should be an integer value; ``step`` should be a `~astropy.units.Quantity` with 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 `~sbpy.data.TimeScaleWarning` will be raised. location : various, optional Location of the observer as an IAU observatory code [OBSCODES]_ (string), a 3-element array of Earth longitude, latitude, altitude, or an `~astropy.coordinates.EarthLocation`. Longitude and latitude should be parseable by `~astropy.coordinates.Angle`, and altitude should be parsable by `~astropy.units.Quantity` (with units of length). If ``None``, then the geocenter (code 500) is used. **kwargs Additional keyword arguments are passed to `~astroquery.mpc.MPC.get_ephemerides`: ``eph_type``, ``ra_format``, ``dec_format``, ``proper_motion``, ``proper_motion_unit``, ``suppress_daytime``, ``suppress_set``, ``perturbed``, ``unc_links``, ``cache``. Returns ------- `~Ephem` object The resulting object will be populated with columns as defined in `~astroquery.mpc.get_ephemerides`; refer to that document on information on how to modify the list of queried parameters. 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 Notes ----- * All properties are provided in the J2000.0 reference system. * See `astroquery.mpc.MPC.get_ephemerides` and the Minor Planet Ephemeris Service user's guide [MPES]_ for details, including acceptable target names. References ---------- .. [MPES] Wiliams, G. The Minor Planet Ephemeris Service. https://minorplanetcenter.org/iau/info/MPES.pdf .. [OBSCODES] IAU Minor Planet Center. List of observatory codes. https://minorplanetcenter.org/iau/lists/ObsCodesF.html """ # parameter check # if targetids is a list, run separate Horizons queries and append if not isinstance(targetids, (list, ndarray, tuple)): targetids = [targetids] _epochs = None # avoid modifying epochs in-place start = None if epochs is None: _epochs = Time([Time.now()]) elif isinstance(epochs, Time): if not iterable(epochs): _epochs = Time([epochs]) else: _epochs = epochs if _epochs.scale is not 'utc': warn(('converting {} epochs to utc for use in ' 'astroquery.mpc').format(_epochs.scale), TimeScaleWarning) _epochs = _epochs.utc elif isinstance(epochs, dict): _epochs = epochs.copy() start = _epochs['start'] # required if start.scale is not 'utc': warn(('converting {} start epoch to utc for use in ' 'astroquery.mpc').format(start.scale), TimeScaleWarning) start = start.utc step = _epochs.get('step') stop = _epochs.get('stop') if stop is not None and stop.scale is not 'utc': warn(('converting {} stop epoch to utc for use in ' 'astroquery.mpc').format(stop.scale), TimeScaleWarning) stop = stop.utc number = _epochs.get('number') if step is not None and stop is None: step = u.Quantity(step) if step.unit not in (u.d, u.h, u.min, u.s): raise QueryError( 'step must have units of days, hours, minutes,' ' or seconds') if stop is not None: if step is not None and number is None: # start and stop both defined, estimate number of steps dt = (Time(stop).jd - Time(start).jd) * u.d number = int((dt / step).decompose()) + 1 elif step is None and number is not None: step = int( (stop - start).jd * 1440 / (number - 1)) * u.minute else: raise QueryError( ('epoch definition unclear; step xor number ' 'must be provided with start and stop')) else: raise ValueError('Invalid `epochs` parameter') # append ephemerides table for each targetid all_eph = None for targetid in targetids: try: # get ephemeris if start is None: eph = [] for i in range(len(_epochs)): e = MPC.get_ephemeris(targetid, location=location, start=Time(_epochs[i], scale='utc'), number=1, **kwargs) e['Date'] = e['Date'].iso # for vstack to work eph.append(e) eph = vstack(eph) eph['Date'] = Time(eph['Date'], scale='utc') else: eph = MPC.get_ephemeris(targetid, location=location, start=start, step=step, number=number, **kwargs) except InvalidQueryError as e: raise QueryError('Error raised by astroquery.mpc: {:s}'.format( str(e))) # add targetname column eph.add_column(Column([targetid] * len(eph), name='Targetname'), index=0) if all_eph is None: all_eph = eph else: all_eph = vstack([all_eph, eph]) # if ra_format or dec_format is defined, then units must be # dropped or else QTable will raise an exception because # strings cannot have units if 'ra_format' in kwargs: all_eph['RA'].unit = None if 'dec_format' in kwargs: all_eph['Dec'].unit = None return cls.from_table(all_eph)