astroquery:docs

SDSS

class astroquery.sdss.core.SDSS[source] [edit on github]

Bases: astroquery.query.BaseQuery

Attributes Summary

AVAILABLE_TEMPLATES
BASE_URL str(object=’‘) -> string
IMAGING str(object=’‘) -> string
MAXQUERIES int(x=0) -> int or long
QUERY_URL str(object=’‘) -> string
SPECTRO_1D str(object=’‘) -> string
TEMPLATES str(object=’‘) -> string

Methods Summary

get_images(*args, **kwds) Download an image from SDSS.
get_images_async(*args, **kwds) Download an image from SDSS.
get_spectra(*args, **kwds) Download spectrum from SDSS.
get_spectra_async(*args, **kwds) Download spectrum from SDSS.
get_spectral_template(*args, **kwds) Download spectral templates from SDSS DR-2, which are located here: http://www.sdss.org/dr5/algorithms/spectemplates/ There 32 spectral templates available from DR-2, from stellar spectra, to galaxies, to quasars.
get_spectral_template_async(*args, **kwds) Download spectral templates from SDSS DR-2, which are located here: http://www.sdss.org/dr5/algorithms/spectemplates/ There 32 spectral templates available from DR-2, from stellar spectra, to galaxies, to quasars.
query_region(*args, **kwds) Queries the service and returns a table object.
query_region_async(*args, **kwds) Used to query a region around given coordinates.

Attributes Documentation

AVAILABLE_TEMPLATES = {'star_L1': 15, 'star_FA': 5, 'galaxy': [24, 25, 26], 'star_G': [8, 9], 'star_F': [6, 7], 'star_OB': 1, 'galaxy_lrg': 28, 'galaxy_late': 27, 'star_M5': 13, 'star_M3': 12, 'star_M1': 11, 'galaxy_early': 23, 'qso': 29, 'star_M8': 14, 'star_K': 10, 'star_O': 0, 'star_B': 2, 'star_A': [3, 4], 'star_carbon': [17, 18, 19], 'qso_bal': [30, 31], 'star_wd': [16, 20, 21], 'qso_bright': 32, 'star_Ksubdwarf': 22}
BASE_URL = 'http://das.sdss.org'
IMAGING = 'http://das.sdss.org/www/cgi-bin/drC'
MAXQUERIES = 1
QUERY_URL = 'http://cas.sdss.org/public/en/tools/search/x_sql.asp'
SPECTRO_1D = 'http://das.sdss.org/spectro/1d_26'
TEMPLATES = 'http://www.sdss.org/dr5/algorithms/spectemplates/spDR2'

Methods Documentation

static get_images(*args, **kwds)[source] [edit on github]

Download an image from SDSS.

Querying SDSS for images will return the entire plate. For subsequent analyses of individual objects

Parameters :

crossID : dict

Dictionary that must contain the run, rerun, camcol, and field for desired image. These parameters can be passed separately as well. All are required. Most convenient to pass the result of function astroquery.sdss.crossID.

band : str, list

Could be individual band, or list of bands. Options: u, g, r, i, or z

Returns :

List of PyFITS HDUList objects. :

static get_images_async(*args, **kwds)[source] [edit on github]

Download an image from SDSS.

Querying SDSS for images will return the entire plate. For subsequent analyses of individual objects

Parameters :

crossID : dict

Dictionary that must contain the run, rerun, camcol, and field for desired image. These parameters can be passed separately as well. All are required. Most convenient to pass the result of function astroquery.sdss.crossID.

band : str, list

Could be individual band, or list of bands. Options: u, g, r, i, or z

Returns :

List of PyFITS HDUList objects. :

static get_spectra(*args, **kwds)[source] [edit on github]

Download spectrum from SDSS.

Parameters :matches : astropy.table.Table instance (result of query_region).
Returns :List of PyFITS HDUList objects. :
static get_spectra_async(*args, **kwds)[source] [edit on github]

Download spectrum from SDSS.

Parameters :matches : astropy.table.Table instance (result of query_region).
Returns :A list of context-managers that yield readable file-like objects :
static get_spectral_template(*args, **kwds)[source] [edit on github]

Download spectral templates from SDSS DR-2, which are located here:

There 32 spectral templates available from DR-2, from stellar spectra, to galaxies, to quasars. To see the available templates, do:

from astroquery.sdss import SDSS print sdss.AVAILABLE_TEMPLATES
Parameters :

kind : str, list

Which spectral template to download? Options are stored in the dictionary astroquery.sdss.SDSS.AVAILABLE_TEMPLATES

Returns :

List of PyFITS HDUList objects. :

Examples

>>> qso = SDSS.get_spectral_template(kind='qso')
>>> Astar = SDSS.get_spectral_template(kind='star_A')
>>> Fstar = SDSS.get_spectral_template(kind='star_F')
static get_spectral_template_async(*args, **kwds)[source] [edit on github]

Download spectral templates from SDSS DR-2, which are located here:

There 32 spectral templates available from DR-2, from stellar spectra, to galaxies, to quasars. To see the available templates, do:

from astroquery.sdss import SDSS print sdss.AVAILABLE_TEMPLATES
Parameters :

kind : str, list

Which spectral template to download? Options are stored in the dictionary astroquery.sdss.SDSS.AVAILABLE_TEMPLATES

Returns :

List of PyFITS HDUList objects. :

Examples

>>> qso = SDSS.get_spectral_template(kind='qso')
>>> Astar = SDSS.get_spectral_template(kind='star_A')
>>> Fstar = SDSS.get_spectral_template(kind='star_F')
static query_region(*args, **kwds) [edit on github]

Queries the service and returns a table object. Used to query a region around given coordinates. Equivalent to the object cross-ID from the web interface.

Parameters :

coordinates : str or astropy.coordinates object

The target around which to search. It may be specified as a string in which case it is resolved using online services or as the appropriate astropy.coordinates object. ICRS coordinates may also be entered as strings as specified in the astropy.coordinates module.

radius : str or astropy.units.Quantity object, optional

The string must be parsable by astropy.coordinates.Angle. The appropriate Quantity object from astropy.units may also be used. Defaults to 2 arcsec.

fields : list, optional

SDSS PhotoObj or SpecObj quantities to return. If None, defaults to quantities required to find corresponding spectra and images of matched objects (e.g. plate, fiberID, mjd, etc.).

spectro : bool, optional

Look for spectroscopic match in addition to photometric match? If True, objects will only count as a match if photometry and spectroscopy exist. If False, will look for photometric matches only.

Returns :

An `astropy.table.Table` object :

Examples

>>> from astroquery.sdss import SDSS
>>> from astropy import coordinates as coords
>>> result = SDSS.query_region(coords.ICRSCoordinates('0h8m05.63s +14d50m23.3s'))
static query_region_async(*args, **kwds)[source] [edit on github]

Used to query a region around given coordinates. Equivalent to the object cross-ID from the web interface.

Parameters :

coordinates : str or astropy.coordinates object

The target around which to search. It may be specified as a string in which case it is resolved using online services or as the appropriate astropy.coordinates object. ICRS coordinates may also be entered as strings as specified in the astropy.coordinates module.

radius : str or astropy.units.Quantity object, optional

The string must be parsable by astropy.coordinates.Angle. The appropriate Quantity object from astropy.units may also be used. Defaults to 2 arcsec.

fields : list, optional

SDSS PhotoObj or SpecObj quantities to return. If None, defaults to quantities required to find corresponding spectra and images of matched objects (e.g. plate, fiberID, mjd, etc.).

spectro : bool, optional

Look for spectroscopic match in addition to photometric match? If True, objects will only count as a match if photometry and spectroscopy exist. If False, will look for photometric matches only.

Returns :

result : astropy.table.Table

The result of the query as an astropy.table.Table object.

Examples

>>> from astroquery.sdss import SDSS
>>> from astropy import coordinates as coords
>>> result = SDSS.query_region(coords.ICRSCoordinates('0h8m05.63s +14d50m23.3s'))

Page Contents