dustmap modules¶
bayestar (Green et al. 2015, 2018)¶
- class dustmaps.bayestar.BayestarQuery(map_fname=None, max_samples=None, version='bayestar2019')[source]¶
Bases:
dustmaps.map_base.DustMapQueries the Bayestar 3D dust maps (Green, Schlafly, Finkbeiner et al. 2015, 2018). The maps cover the Pan-STARRS 1 footprint (dec > -30 deg) amounting to three-quarters of the sky.
- __init__(map_fname=None, max_samples=None, version='bayestar2019')[source]¶
- Parameters
map_fname (Optional[
str]) – Filename of the Bayestar map. Defaults toNone, meaning that the default location is used.max_samples (Optional[
int]) – Maximum number of samples of the map to load. Use a lower number in order to decrease memory usage. Defaults toNone, meaning that all samples will be loaded.version (Optional[
str]) – The map version to download. Valid versions are'bayestar2019'(Green, Schlafly, Finkbeiner et al. 2019),'bayestar2017'(Green, Schlafly, Finkbeiner et al. 2018) and'bayestar2015'(Green, Schlafly, Finkbeiner et al. 2015). Defaults to'bayestar2015'.
- property distances¶
Returns the distance bin edges that the map uses. The return type is
astropy.units.Quantity, which stores unit-full quantities.
- property distmods¶
Returns the distance modulus bin edges that the map uses. The return type is
astropy.units.Quantity, with units of mags.
- query(coords, mode='random_sample', return_flags=False, pct=None)[source]¶
Returns reddening at the requested coordinates. There are several different query modes, which handle the probabilistic nature of the map differently.
- Parameters
coords (
astropy.coordinates.SkyCoord) – The coordinates to query.mode (Optional[
str]) – Seven different query modes are available: ‘random_sample’, ‘random_sample_per_pix’ ‘samples’, ‘median’, ‘mean’, ‘best’ and ‘percentile’. Themodedetermines how the output will reflect the probabilistic nature of the Bayestar dust maps.return_flags (Optional[
bool]) – IfTrue, then QA flags will be returned in a second numpy structured array. That is, the query will returnret, :obj:’flags`, whereretis the normal return value, containing reddening. Defaults toFalse.pct (Optional[
floator list/array offloat]) – If the mode ispercentile, thenpctspecifies which percentile(s) is (are) returned.
- Returns
Reddening at the specified coordinates, in magnitudes of reddening.
The conversion to E(B-V) (or other reddening units) depends on whether
version='bayestar2019'(the default),'bayestar2017'or'bayestar2015'was selected when theBayestarQueryobject was created. To convert Bayestar2019 to Pan-STARRS 1 extinctions, multiply by the coefficients given in Table 1 of Green et al. (2019). For Bayestar2017, use the coefficients given in Table 1 of Green et al. (2018). Conversion to extinction in non-PS1 passbands depends on the choice of extinction law. To convert Bayestar2015 to extinction in various passbands, multiply by the coefficients in Table 6 of Schlafly & Finkbeiner (2011). See Green et al. (2015, 2018) for more detailed discussion of how to convert the Bayestar dust maps into reddenings or extinctions in different passbands.The shape of the output depends on the
mode, and on whethercoordscontains distances.If
coordsdoes not specify distance(s), then the shape of the output begins withcoords.shape. Ifcoordsdoes specify distance(s), then the shape of the output begins withcoords.shape + ([number of distance bins],).If
modeis'random_sample', then at each coordinate/distance, a random sample of reddening is given.If
modeis'random_sample_per_pix', then the sample chosen for each angular pixel of the map will be consistent. For example, if two query coordinates lie in the same map pixel, then the same random sample will be chosen from the map for both query coordinates.If
modeis'median', then at each coordinate/distance, the median reddening is returned.If
modeis'mean', then at each coordinate/distance, the mean reddening is returned.If
modeis'best', then at each coordinate/distance, the maximum posterior density reddening is returned (the “best fit”).If
modeis'percentile', then an additional keyword argument,pct, must be specified. At each coordinate/distance, the requested percentiles (inpct) will be returned. Ifpctis a list/array, then the last axis of the output will correspond to different percentiles.Finally, if
modeis'samples', then at each coordinate/distance, all samples are returned. The last axis of the output will correspond to different samples.If
return_flagsisTrue, then in addition to reddening, a structured array containing QA flags will be returned. If the input coordinates include distances, the QA flags will be"converged"(whether or not the line-of-sight fit converged in a given pixel) and"reliable_dist"(whether or not the requested distance is within the range considered reliable, based on the inferred stellar distances). If the input coordinates do not include distances, then instead of"reliable_dist", the flags will include"min_reliable_distmod"and"max_reliable_distmod", the minimum and maximum reliable distance moduli in the given pixel.
- class dustmaps.bayestar.BayestarWebQuery(api_url=None, version='bayestar2019')[source]¶
Bases:
dustmaps.map_base.WebDustMapRemote query over the web for the Bayestar 3D dust maps (Green, Schlafly, Finkbeiner et al. 2015, 2018, 2019). The maps cover the Pan-STARRS 1 footprint (dec > -30 deg) amounting to three-quarters of the sky.
This query object does not require a local version of the data, but rather an internet connection to contact the web API. The query functions have the same inputs and outputs as their counterparts in
BayestarQuery.- __init__(api_url=None, version='bayestar2019')[source]¶
- Parameters
version (Optional[
str]) – The map version to download. Valid versions are'bayestar2019'(Green, Schlafly, Finkbeiner et al. 2019),'bayestar2017'(Green, Schlafly, Finkbeiner et al. 2018) and'bayestar2015'(Green, Schlafly, Finkbeiner et al. 2015). Defaults to'bayestar2019'.
- dustmaps.bayestar.fetch(version='bayestar2019')[source]¶
Downloads the specified version of the Bayestar dust map.
- Parameters
version (Optional[
str]) – The map version to download. Valid versions are'bayestar2019'(Green, Schlafly, Finkbeiner et al. 2019),'bayestar2017'(Green, Schlafly, Finkbeiner et al. 2018) and'bayestar2015'(Green, Schlafly, Finkbeiner et al. 2015). Defaults to'bayestar2019'.- Raises
ValueError – The requested version of the map does not exist.
DownloadError – Either no matching file was found under the given DOI, or the MD5 sum of the file was not as expected.
requests.exceptions.HTTPError – The given DOI does not exist, or there was a problem connecting to the Dataverse.
- dustmaps.bayestar.lb2pix(nside, l, b, nest=True)[source]¶
Converts Galactic (l, b) to HEALPix pixel index.
- Parameters
nside (
int) – The HEALPixnsideparameter.l (
float, or array offloat) – Galactic longitude, in degrees.b (
float, or array offloat) – Galactic latitude, in degrees.nest (Optional[
bool]) – IfTrue(the default), nested pixel ordering will be used. IfFalse, ring ordering will be used.
- Returns
The HEALPix pixel index or indices. Has the same shape as the input
landb.
bh (Burstein & Heiles 1982)¶
- class dustmaps.bh.BHQuery(bh_dir=None)[source]¶
Bases:
dustmaps.map_base.DustMapQueries the Burstein & Heiles (1982) reddening map.
- __init__(bh_dir=None)[source]¶
- Parameters
bh_dir (Optional[str]) – The directory containing the Burstein & Heiles dust map. Defaults to None, meaning that the default directory is used.
- query(coords)[source]¶
Returns E(B-V) at the specified location(s) on the sky.
- Parameters
coords (astropy.coordinates.SkyCoord) – The coordinates to query.
- Returns
A float array of reddening, in units of E(B-V), at the given coordinates. The shape of the output is the same as the shape of the coordinates stored by coords.
chen2014 (Chen et al. 2014)¶
- class dustmaps.chen2014.Chen2014Query(map_fname=None)[source]¶
Bases:
dustmaps.unstructured_map.UnstructuredDustMapThe 3D dust map of Chen et al. (2014), based on stellar photometry from the Xuyi Schmidt Telescope Photometric Survey of the Galactic Anticentre. The map covers 140 deg < l < 240 deg, -60 deg < b < 40 deg.
- __init__(map_fname=None)[source]¶
- Parameters
map_fname (Optional[
str]) – Filename at which the map is stored. Defaults toNone, meaning that the default filename is used.
- property distances¶
Returns the distance bins that the map uses. The return type is
astropy.units.Quantity, which stores unit-full quantities.
- query(coords, return_sigma=False)[source]¶
Returns r-band extinction, A_r, at the given coordinates. Can also return uncertainties.
- Parameters
coords (
astropy.coordinates.SkyCoord) – The coordinates to query.return_sigma (Optional[
bool]) – IfTrue, returns the uncertainty in extinction as well. Defaults toFalse.
- Returns
Extinction in the r-band at the specified coordinates, in mags. The shape of the output depends on whether
coordscontains distances.If
coordsdoes not specify distance(s), then the shape of the output begins withcoords.shape. Ifcoordsdoes specify distance(s), then the shape of the output begins withcoords.shape + ([number of distance bins],).
- dustmaps.chen2014.ascii2h5(dat_fname, h5_fname)[source]¶
Converts from the original ASCII format of the Chen+ (2014) 3D dust map to the HDF5 format.
- Parameters
dat_fname (
str) – Filename of the original ASCII .dat file.h5_fname (
str) – Output filename to write the resulting HDF5 file to.
- dustmaps.chen2014.fetch(clobber=False)[source]¶
Downloads the Chen et al. (2014) dust map.
- Parameters
clobber (Optional[
bool]) – IfTrue, any existing file will be overwritten, even if it appears to match. IfFalse(the default),fetch()will attempt to determine if the dataset already exists. This determination is not 100% robust against data corruption.
iphas (Sale et al. 2014)¶
- class dustmaps.iphas.IPHASQuery(map_fname=None)[source]¶
Bases:
dustmaps.unstructured_map.UnstructuredDustMapThe 3D dust map of Sale et al. (2014), based on IPHAS imaging in the Galactic plane. The map covers 30 deg < l < 115 deg, -5 deg < b < 5 deg.
- __init__(map_fname=None)[source]¶
- Parameters
map_fname (Optional[
str]) – Filename at which the map is stored. Defaults toNone, meaning that the default filename is used.
- property distances¶
Returns the distance bins that the map uses. The return type is
astropy.units.Quantity, which stores unit-full quantities.
- query(coords, mode='random_sample')[source]¶
Returns A0 at the given coordinates. There are several different query modes, which handle the probabilistic nature of the map differently.
- Parameters
coords (
astropy.coordinates.SkyCoord) – The coordinates to query.mode (Optional[
str]) – Five different query modes are available:'random_sample','random_sample_per_pix','samples','median'and'mean'. Themodedetermines how the output will reflect the probabilistic nature of the IPHAS dust map.
- Returns
Monochromatic extinction, A0, at the specified coordinates, in mags. The shape of the output depends on the
mode, and on whethercoordscontains distances.If
coordsdoes not specify distance(s), then the shape of the output begins with coords.shape. If coords does specify distance(s), then the shape of the output begins withcoords.shape + ([number of distance bins],).If
modeis'random_sample', then at each coordinate/distance, a random sample of reddening is given.If
modeis'random_sample_per_pix', then the sample chosen for each angular pixel of the map will be consistent. For example, if two query coordinates lie in the same map pixel, then the same random sample will be chosen from the map for both query coordinates.If
modeis'median', then at each coordinate/distance, the median reddening is returned.If
modeis'mean', then at each coordinate/distance, the mean reddening is returned.Finally, if
modeis'samples', then all at each coordinate/distance, all samples are returned.
- dustmaps.iphas.ascii2h5(dirname, output_fname)[source]¶
Converts from a directory of tarballed ASCII “.samp” files to a single HDF5 file. Essentially, converts from the original release format to a single HDF5 file.
- dustmaps.iphas.fetch(clobber=False)[source]¶
Downloads the IPHAS 3D dust map of Sale et al. (2014).
- Parameters
clobber (Optional[bool]) – If
True, any existing file will be overwritten, even if it appears to match. IfFalse(the default),fetch()will attempt to determine if the dataset already exists. This determination is not 100% robust against data corruption.
lenz2017 (Lenz, Hensley & Doré 2017)¶
- class dustmaps.lenz2017.Lenz2017Query(map_fname=None)[source]¶
Bases:
dustmaps.healpix_map.HEALPixFITSQueryQueries the Lenz, Hensley & Doré (2017) dust map: http://arxiv.org/abs/1706.00011
marshall (Marshall et al. 2006)¶
- class dustmaps.marshall.MarshallQuery(map_fname=None)[source]¶
Bases:
dustmaps.map_base.DustMapGalactic-plane 3D dust map of Marshall et al. (2006), based on 2MASS photometry.
- __init__(map_fname=None)[source]¶
- Parameters
map_fname (Optional[
str]) – Filename at which the map is stored. Defaults toNone, meaning that the default filename is used.
- query(coords, return_sigma=False)[source]¶
Returns 2MASS Ks-band extinction at the given coordinates.
- Parameters
coords (
astropy.coordinates.SkyCoord) – The coordinates to query. Must contain distances.return_sigma (Optional[
bool]) – IfTrue, returns the uncertainty in extinction as well. Defaults toFalse.
- Returns
Extinction at the specified coordinates, in mags of 2MASS Ks-band extinction. If
return_sigmaisTrue, then the uncertainty in reddening is also returned, so that the output is(A, sigma_A), where bothAandsigma_Ahave the same shape as the input coordinates.
- dustmaps.marshall.dat2hdf5(table_dir)[source]¶
Convert the Marshall et al. (2006) map from *.dat.gz to *.hdf5.
- dustmaps.marshall.fetch(clobber=False)[source]¶
Downloads the Marshall et al. (2006) dust map, which is based on 2MASS stellar photometry.
- Parameters
clobber (Optional[
bool]) – IfTrue, any existing file will be overwritten, even if it appears to match. IfFalse(the default),fetch()will attempt to determine if the dataset already exists. This determination is not 100% robust against data corruption.
pg2010 (Peek & Graves 2010)¶
- class dustmaps.pg2010.PG2010Query(map_dir=None, component='dust')[source]¶
Bases:
dustmaps.sfd.SFDBaseQueries the Peek & Graves (2010) correction to the SFD’98 dust reddening map.
- __init__(map_dir=None, component='dust')[source]¶
- Parameters
map_dir (Optional[
str]) – The directory containing the SFD map. Defaults toNone, which means thatdustmapswill look in its default data directory.component (Optional[
str]) –'dust'(the default) to load the correction to E(B-V), or'err'to load the uncertainty in the correction.
- query(coords, order=1)[source]¶
Returns the P&G (2010) correction to the SFD’98 E(B-V) at the specified location(s) on the sky. If component is ‘err’, then return the uncertainty in the correction.
- Parameters
coords (
astropy.coordinates.SkyCoord) – The coordinates to query.order (Optional[
int]) – Interpolation order to use. Defaults to1, for linear interpolation.
- Returns
A float array containing the P&G (2010) correction (or its uncertainty) to SFD’98 at every input coordinate. The shape of the output will be the same as the shape of the coordinates stored by
coords.
planck (Planck Collaboration 2013)¶
- class dustmaps.planck.PlanckQuery(map_fname=None, component='extragalactic')[source]¶
Bases:
dustmaps.healpix_map.HEALPixFITSQueryQueries the Planck Collaboration (2013) dust map.
- __init__(map_fname=None, component='extragalactic')[source]¶
- Parameters
map_fname (Optional[
str]) – Filename of the Planck map. Defaults to`None, meaning that the default location is used.component (Optional[str]) – Which measure of reddening to use. There are seven valid components. Three denote reddening measures:
'extragalactic','tau'and'radiance'. Four refer to dust properties:'temperature','beta','err_temp'and'err_beta'. Defaults to'extragalactic'.
- query(coords, **kwargs)[source]¶
Returns E(B-V) (or a different Planck dust inference, depending on how the class was intialized) at the specified location(s) on the sky.
- Parameters
coords (
astropy.coordinates.SkyCoord) – The coordinates to query.- Returns
A float array of the selected Planck component, at the given coordinates. The shape of the output is the same as the shape of the coordinates stored by
coords. If extragalactic E(B-V), tau_353 or radiance was chosen, then the output has units of magnitudes of E(B-V). If the selected Planck component is temperature (or temperature error), then anastropy.Quantityis returned, with units of Kelvin. If beta (or beta error) was chosen, then the output is unitless.
sfd (Schlegel, Finkbeiner & Davis 1998)¶
- class dustmaps.sfd.SFDBase(base_fname)[source]¶
Bases:
dustmaps.map_base.DustMapQueries maps stored in the same format as Schlegel, Finkbeiner & Davis (1998).
- __init__(base_fname)[source]¶
- Parameters
base_fname (str) – The map should be stored in two FITS files, named
base_fname + '_' + X + '.fits', whereXis'ngp'and'sgp'.
- query(coords, order=1)[source]¶
Returns the map value at the specified location(s) on the sky.
- Parameters
coords (astropy.coordinates.SkyCoord) – The coordinates to query.
order (Optional[int]) – Interpolation order to use. Defaults to 1, for linear interpolation.
- Returns
A float array containing the map value at every input coordinate. The shape of the output will be the same as the shape of the coordinates stored by coords.
- class dustmaps.sfd.SFDQuery(map_dir=None)[source]¶
Bases:
dustmaps.sfd.SFDBaseQueries the Schlegel, Finkbeiner & Davis (1998) dust reddening map.
- __init__(map_dir=None)[source]¶
- Parameters
map_dir (Optional[str]) – The directory containing the SFD map. Defaults to None, which means that dustmaps will look in its default data directory.
- query(coords, order=1)[source]¶
Returns E(B-V) at the specified location(s) on the sky. See Table 6 of Schlafly & Finkbeiner (2011) for instructions on how to convert this quantity to extinction in various passbands.
- Parameters
coords (astropy.coordinates.SkyCoord) – The coordinates to query.
order (Optional[int]) – Interpolation order to use. Defaults to 1, for linear interpolation.
- Returns
A float array containing the SFD E(B-V) at every input coordinate. The shape of the output will be the same as the shape of the coordinates stored by coords.
- class dustmaps.sfd.SFDWebQuery(api_url=None)[source]¶
Bases:
dustmaps.map_base.WebDustMapRemote query over the web for the Schlegel, Finkbeiner & Davis (1998) dust map.
This query object does not require a local version of the data, but rather an internet connection to contact the web API. The query functions have the same inputs and outputs as their counterparts in
SFDQuery.- __init__(api_url=None)[source]¶
Initialize the
WebDustMapobject.- Parameters
api_url (Optional[
str]) – The base URL for the API. Defaults to'http://argonaut.skymaps.info/api/v2/'.map_name (Optional[
str]) – The name of the dust map to query. For example, the Green et al. (2015) dust map is hosted athttp://argonaut.skymaps.info/api/v2/bayestar2015, so the correct specifier for that map ismap_name='bayestar2015'.
fetch_utils¶
- exception dustmaps.fetch_utils.DownloadError[source]¶
Bases:
dustmaps.fetch_utils.ErrorAn exception that occurs while trying to download a file.
- exception dustmaps.fetch_utils.Error[source]¶
Bases:
Exception- __weakref__¶
list of weak references to the object (if defined)
- class dustmaps.fetch_utils.FileTransferProgressBar(content_length)[source]¶
Bases:
progressbar.bar.ProgressBar
- dustmaps.fetch_utils.check_md5sum(fname, md5sum, chunk_size=1024)[source]¶
Checks that a file exists, and has the correct MD5 checksum.
- Parameters
fname (str) – The filename of the file.
md5sum (str) – The expected MD5 sum.
chunk_size (Optional[int]) – Process in chunks of this size (in Bytes). Defaults to 1024.
- dustmaps.fetch_utils.dataverse_download_doi(doi, local_fname=None, file_requirements={}, clobber=False)[source]¶
Downloads a file from the Dataverse, using a DOI and set of metadata parameters to locate the file.
- Parameters
doi (str) – Digital Object Identifier (DOI) containing the file.
local_fname (Optional[str]) – Local filename to download the file to. If None, then use the filename provided by the Dataverse. Defaults to None.
file_requirements (Optional[dict]) – Select the file containing the given metadata entries. If multiple files meet these requirements, only the first in downloaded. Defaults to {}, corresponding to no requirements.
- Raises
DownloadError – Either no matching file was found under the given DOI, or the MD5 sum of the file was not as expected.
requests.exceptions.HTTPError – The given DOI does not exist, or there was a problem connecting to the Dataverse.
- dustmaps.fetch_utils.dataverse_search_doi(doi)[source]¶
Fetches metadata pertaining to a Digital Object Identifier (DOI) in the Harvard Dataverse.
- Parameters
doi (str) – The Digital Object Identifier (DOI) of the entry in the Dataverse.
- Raises
requests.exceptions.HTTPError – The given DOI does not exist, or there was a problem connecting to the Dataverse.
- dustmaps.fetch_utils.download(url, fname=None)[source]¶
Downloads a file.
- Parameters
url (str) – The URL to download.
fname (Optional[str]) – The filename to store the downloaded file in. If None, take the filename from the URL. Defaults to None.
- Returns
The filename the URL was downloaded to.
- Raises
requests.exceptions.HTTPError – There was a problem connecting to the URL.
- dustmaps.fetch_utils.download_and_verify(url, md5sum, fname=None, chunk_size=1024, clobber=False, verbose=True)[source]¶
Downloads a file and verifies the MD5 sum.
- Parameters
url (str) – The URL to download.
md5sum (str) – The expected MD5 sum.
fname (Optional[str]) – The filename to store the downloaded file in. If None, infer the filename from the URL. Defaults to None.
chunk_size (Optional[int]) – Process in chunks of this size (in Bytes). Defaults to 1024.
clobber (Optional[bool]) – If True, any existing, identical file will be overwritten. If False, the MD5 sum of any existing file with the destination filename will be checked. If the MD5 sum does not match, the existing file will be overwritten. Defaults to False.
verbose (Optional[bool]) – If True (the default), then a progress bar will be shownd during downloads.
- Returns
The filename the URL was downloaded to.
- Raises
DownloadError – The MD5 sum of the downloaded file does not match md5sum.
requests.exceptions.HTTPError – There was a problem connecting to the URL.
- dustmaps.fetch_utils.get_md5sum(fname, chunk_size=1024)[source]¶
Returns the MD5 checksum of a file.
- Parameters
fname (str) – Filename
chunk_size (Optional[int]) – Size (in Bytes) of the chunks that should be read in at once. Increasing chunk size reduces the number of reads required, but increases the memory usage. Defaults to 1024.
- Returns
The MD5 checksum of the file, which is a string.
- dustmaps.fetch_utils.h5_file_exists(fname, size_guess=None, rtol=0.1, atol=1.0, dsets={})[source]¶
Returns
Trueif an HDF5 file exists, has the expected file size, and contains (at least) the given datasets, with the correct shapes.- Parameters
fname (str) – Filename to check.
size_guess (Optional[int]) – Expected size (in Bytes) of the file. If
None(the default), then filesize is not checked.rtol (Optional[float]) – Relative tolerance for filesize.
atol (Optional[float]) – Absolute tolerance (in Bytes) for filesize.
dsets (Optional[dict]) – Dictionary specifying expected datasets. Each key is the name of a dataset, while each value is the expected shape of the dataset. Defaults to
{}, meaning that no datasets are checked.
- Returns
Trueif the file matches by all given criteria.
map_base¶
- class dustmaps.map_base.DustMap[source]¶
Bases:
objectBase class for querying dust maps. For each individual dust map, a different subclass should be written, implementing the
query()function.- __call__(coords, **kwargs)[source]¶
An alias for
DustMap.query.
- query(coords, **kwargs)[source]¶
Query the map at a set of coordinates.
- Parameters
coords (
astropy.coordinates.SkyCoord) – The coordinates at which to query the map.- Raises
NotImplementedError – This function must be defined by derived classes.
- query_equ(ra, dec, d=None, frame='icrs', **kwargs)[source]¶
Query using Equatorial coordinates. By default, the ICRS frame is used, although other frames implemented by
astropy.coordinatesmay also be specified.- Parameters
ra (
float, scalar or array-like) – Galactic longitude, in degrees, or as anastropy.unit.Quantity.dec (float, scalar or array-like) – Galactic latitude, in degrees, or as an
astropy.unit.Quantity.d (Optional[
float, scalar or array-like]) – Distance from the Solar System, in kpc, or as anastropy.unit.Quantity. Defaults toNone, meaning no distance is specified.frame (Optional[
icrs]) – The coordinate system. Can be'icrs'(the default),'fk5','fk4'or'fk4noeterms'.**kwargs – Any additional keyword arguments accepted by derived classes.
- Returns
The results of the query, which must be implemented by derived classes.
- query_gal(l, b, d=None, **kwargs)[source]¶
Query using Galactic coordinates.
- Parameters
l (
float, scalar or array-like) – Galactic longitude, in degrees, or as anastropy.unit.Quantity.b (
float, scalar or array-like) – Galactic latitude, in degrees, or as anastropy.unit.Quantity.d (Optional[
float, scalar or array-like]) – Distance from the Solar System, in kpc, or as anastropy.unit.Quantity. Defaults toNone, meaning no distance is specified.**kwargs – Any additional keyword arguments accepted by derived classes.
- Returns
The results of the query, which must be implemented by derived classes.
- class dustmaps.map_base.WebDustMap(api_url=None, map_name='')[source]¶
Bases:
objectBase class for querying dust maps through a web API. For each individual dust map, a different subclass should be written, specifying the base URL.
- __call__(coords, **kwargs)[source]¶
An alias for
WebDustMap.query().
- query(coords, **kwargs)[source]¶
A web API version of
DustMap.query. See the documentation for the corresponding local query object.- Parameters
coords (
astropy.coordinates.SkyCoord) – The coordinates at which to query the map.
- query_equ(ra, dec, d=None, frame='icrs', **kwargs)[source]¶
A web API version of
DustMap.query_equ(). See the documentation for the corresponding local query object. Queries using Equatorial coordinates. By default, the ICRS frame is used, although other frames implemented byastropy.coordinatesmay also be specified.- Parameters
ra (
float, scalar or array-like) – Galactic longitude, in degrees, or as anastropy.unit.Quantity.dec (
float, scalar or array-like) – Galactic latitude, in degrees, or as anastropy.unit.Quantity.d (Optional[
float, scalar or array-like]) – Distance from the Solar System, in kpc, or as anastropy.unit.Quantity. Defaults toNone, meaning no distance is specified.frame (Optional[icrs]) – The coordinate system. Can be ‘icrs’ (the default), ‘fk5’, ‘fk4’ or ‘fk4noeterms’.
**kwargs – Any additional keyword arguments accepted by derived classes.
- Returns
The results of the query.
- query_gal(l, b, d=None, **kwargs)[source]¶
A web API version of
DustMap.query_gal(). See the documentation for the corresponding local query object. Queries using Galactic coordinates.- Parameters
l (
float, scalar or array-like) – Galactic longitude, in degrees, or as anastropy.unit.Quantity.b (
float, scalar or array-like) – Galactic latitude, in degrees, or as anastropy.unit.Quantity.d (Optional[
float, scalar or array-like]) – Distance from the Solar System, in kpc, or as anastropy.unit.Quantity. Defaults toNone, meaning no distance is specified.**kwargs – Any additional keyword arguments accepted by derived classes.
- Returns
The results of the query.
- dustmaps.map_base.coord2healpix(coords, frame, nside, nest=True)[source]¶
Calculate HEALPix indices from an astropy SkyCoord. Assume the HEALPix system is defined on the coordinate frame
frame.- Parameters
coords (
astropy.coordinates.SkyCoord) – The input coordinates.frame (
str) – The frame in which the HEALPix system is defined.nside (
int) – The HEALPix nside parameter to use. Must be a power of 2.nest (Optional[
bool]) –True(the default) if nested HEALPix ordering is desired.Falsefor ring ordering.
- Returns
An array of pixel indices (integers), with the same shape as the input SkyCoord coordinates (
coords.shape).- Raises
dustexceptions.CoordFrameError – If the specified frame is not supported.
- dustmaps.map_base.ensure_coord_type(f)[source]¶
A decorator for class methods of the form
Class.method(self, coords, **kwargs)
where
coordsis anastropy.coordinates.SkyCoordobject.The decorator raises a
TypeErrorif thecoordsthat gets passed toClass.methodis not anastropy.coordinates.SkyCoordinstance.- Parameters
f (class method) – A function with the signature
(self, coords, **kwargs), wherecoordsis aSkyCoordobject containing an array.- Returns
A function that raises a
TypeErrorifcoordsis not anastropy.coordinates.SkyCoordobject, but which otherwise behaves the same as the decorated function.
- dustmaps.map_base.ensure_flat_coords(f)[source]¶
A decorator for class methods of the form
Class.method(self, coords, **kwargs)
where
coordsis anastropy.coordinates.SkyCoordobject.The decorator ensures that the
coordsthat gets passed toClass.methodis a flat array. It also reshapes the output ofClass.methodto have the same shape (possibly scalar) as the inputcoords. If the output ofClass.methodis a tuple or list (instead of an array), each element in the output is reshaped instead.- Parameters
f (class method) – A function with the signature
(self, coords, **kwargs), wherecoordsis aSkyCoordobject containing an array.- Returns
A function that takes
SkyCoordinput with any shape (including scalar).
- dustmaps.map_base.ensure_flat_galactic(f)[source]¶
A decorator for class methods of the form
Class.method(self, coords, **kwargs)
where
coordsis anastropy.coordinates.SkyCoordobject.The decorator ensures that the
coordsthat gets passed toClass.methodis a flat array of Galactic coordinates. It also reshapes the output ofClass.methodto have the same shape (possibly scalar) as the inputcoords. If the output ofClass.methodis a tuple or list (instead of an array), each element in the output is reshaped instead.- Parameters
f (class method) – A function with the signature
(self, coords, **kwargs), wherecoordsis aSkyCoordobject containing an array.- Returns
A function that takes
SkyCoordinput with any shape (including scalar).
healpix_map¶
- class dustmaps.healpix_map.HEALPixFITSQuery(fname, coord_frame, hdu=0, field=None, dtype='f8')[source]¶
Bases:
dustmaps.healpix_map.HEALPixQueryA HEALPix map class that is initialized from a FITS file.
- __init__(fname, coord_frame, hdu=0, field=None, dtype='f8')[source]¶
- Parameters
fname (str, HDUList, TableHDU or BinTableHDU) – The filename, HDUList or HDU from which the map should be loaded.
coord_frame (str) – The coordinate system in which the HEALPix map is defined. Must be a coordinate frame which
astropyunderstands.hdu (Optional[int or str]) – Specifies which HDU to load the map from. Defaults to
0.field (Optional[int or str]) – Specifies which field (column) to load the map from. Defaults to
None, meaning thathdu.data[:]is used.dtype (Optional[str or type]) – The data will be coerced to this datatype. Can be any type specification that numpy understands. Defaults to
'f8', for IEEE754 double precision.
- class dustmaps.healpix_map.HEALPixQuery(pix_val, nest, coord_frame)[source]¶
Bases:
dustmaps.map_base.DustMapA class for querying HEALPix maps.
- __init__(pix_val, nest, coord_frame)[source]¶
- Parameters
pix_val (array) – Value of the map in every pixel. The length of the array must be of the form 12 * nside**2, where nside is a power of two.
nest (bool) – True if the map uses nested ordering. False if ring ordering is used.
coord_frame (str) – The coordinate system that the HEALPix map is in. Should be one of the frames supported by astropy.coordinates.
unstructured_map¶
- class dustmaps.unstructured_map.UnstructuredDustMap(pix_coords, max_pix_scale, metric_p=2, frame=None)[source]¶
Bases:
dustmaps.map_base.DustMapA class for querying dust maps with unstructured pixels. Sky coordinates are assigned to the nearest pixel.
- __init__(pix_coords, max_pix_scale, metric_p=2, frame=None)[source]¶
- Parameters
pix_coords (array-like
astropy.coordinates.SkyCoord) – The sky coordinates of the pixels.max_pix_scale (scalar
astropy.units.Quantity) – Maximum angular extent of a pixel. If no pixel is within this distance of a query point, NaN will be returned for that query point.metric_p (Optional[
float]) – The metric to use. Defaults to 2, which is the Euclidean metric. A value of 1 corresponds to the Manhattan metric, while a value approaching infinity yields the maximum component metric.frame (Optional[
str]) – The coordinate frame to use internally. Must be a frame understood byastropy.coordinates.SkyCoord. Defaults toNone, meaning that the frame will be inferred frompix_coords.
config¶
- class dustmaps.config.Configuration(fname)[source]¶
Bases:
objectA class that stores the package configuration.
- get(key, default=None)[source]¶
Gets a configuration option, returning a default value if the specified key isn’t set.
- save(force=False)[source]¶
Saves the configuration to a JSON, in the standard config location.
- Parameters
force (Optional[
bool]) – Continue writing, even if the original config file was not loaded properly. This is dangerous, because it could cause the previous configuration options to be lost. Defaults toFalse.- Raises
ConfigError – if the configuration file was not successfully loaded on initialization of the class, and
forceisFalse.
- dustmaps.config.config = <dustmaps.config.Configuration object>¶
The package configuration. This is the object that the user should interact with in order to change settings. For example, to set the directory where large files (e.g., dust maps) will be stored:
from dustmaps.config import config config['data_dir'] = '/path/to/data/directory'
std_paths¶
- dustmaps.std_paths.data_dir()[source]¶
Returns the directory used to store large data files (e.g., dust maps).
json_serializers¶
- class dustmaps.json_serializers.MultiJSONDecoder(*args, **kwargs)[source]¶
Bases:
json.decoder.JSONDecoder- A JSON decoder that can handle:
numpy.ndarraynumpy.dtypeastropy.units.Quantityastropy.coordinates.SkyCoord
- __init__(*args, **kwargs)[source]¶
object_hook, if specified, will be called with the result of every JSON object decoded and its return value will be used in place of the givendict. This can be used to provide custom deserializations (e.g. to support JSON-RPC class hinting).object_pairs_hook, if specified will be called with the result of every JSON object decoded with an ordered list of pairs. The return value ofobject_pairs_hookwill be used instead of thedict. This feature can be used to implement custom decoders. Ifobject_hookis also defined, theobject_pairs_hooktakes priority.parse_float, if specified, will be called with the string of every JSON float to be decoded. By default this is equivalent to float(num_str). This can be used to use another datatype or parser for JSON floats (e.g. decimal.Decimal).parse_int, if specified, will be called with the string of every JSON int to be decoded. By default this is equivalent to int(num_str). This can be used to use another datatype or parser for JSON integers (e.g. float).parse_constant, if specified, will be called with one of the following strings: -Infinity, Infinity, NaN. This can be used to raise an exception if invalid JSON numbers are encountered.If
strictis false (true is the default), then control characters will be allowed inside strings. Control characters in this context are those with character codes in the 0-31 range, including'\t'(tab),'\n','\r'and'\0'.
- dustmaps.json_serializers.deserialize_dtype(d)[source]¶
Deserializes a JSONified
numpy.dtype.- Parameters
d (
dict) – A dictionary representation of adtypeobject.- Returns
A
dtypeobject.
- dustmaps.json_serializers.deserialize_ndarray(d)[source]¶
Deserializes a JSONified
numpy.ndarray. Can handle arrays serialized using any of the methods in this module:"npy","b64","readable".- Parameters
d (dict) – A dictionary representation of an
ndarrayobject.- Returns
An
ndarrayobject.
- dustmaps.json_serializers.deserialize_ndarray_npy(d)[source]¶
Deserializes a JSONified
numpy.ndarraythat was created using numpy’ssavefunction.- Parameters
d (
dict) – A dictionary representation of anndarrayobject, created usingnumpy.save.- Returns
An
ndarrayobject.
- dustmaps.json_serializers.deserialize_quantity(d)[source]¶
Deserializes a JSONified
astropy.units.Quantity.- Parameters
d (
dict) – A dictionary representation of aQuantityobject.- Returns
A
Quantityobject.
- dustmaps.json_serializers.deserialize_skycoord(d)[source]¶
Deserializes a JSONified
astropy.coordinates.SkyCoord.- Parameters
d (
dict) – A dictionary representation of aSkyCoordobject.- Returns
A
SkyCoordobject.
- dustmaps.json_serializers.deserialize_tuple(d)[source]¶
Deserializes a JSONified tuple.
- Parameters
d (
dict) – A dictionary representation of the tuple.- Returns
A tuple.
- dustmaps.json_serializers.get_encoder(ndarray_mode='b64')[source]¶
- Returns a JSON encoder that can handle:
numpy.ndarraynumpy.floating(converted tofloat)numpy.integer(converted toint)numpy.dtypeastropy.units.Quantityastropy.coordinates.SkyCoord
- Parameters
ndarray_mode (Optional[
str]) – Which method to use to serializenumpy.ndarrayobjects. Defaults to'b64', which converts the array data to binary64 encoding (non-human-readable), and stores the datatype/shape in human-readable formats. Other options are'readable', which produces fully human-readable output, and'npy', which uses numpy’s built-insavefunction and produces completely unreadable output. Of all the methods'npy'is the most reliable, but also least human-readable.'readable'produces the most human-readable output, but is the least reliable and loses precision.- Returns
A subclass of
json.JSONEncoder.
- dustmaps.json_serializers.hint_tuples(o)[source]¶
Annotates tuples before JSON serialization, so that they can be reconstructed during deserialization. Each tuple is converted into a dictionary of the form:
{‘_type’: ‘tuple’, ‘items’: (…)}
This function acts recursively on lists, so that tuples nested inside a list (or doubly nested, triply nested, etc.) will also be annotated.
- dustmaps.json_serializers.serialize_dtype(o)[source]¶
Serializes a
numpy.dtype.- Parameters
o (
numpy.dtype) –dtypeto be serialized.- Returns
A dictionary that can be passed to
json.dumps.
- dustmaps.json_serializers.serialize_ndarray_b64(o)[source]¶
Serializes a
numpy.ndarrayin a format where the datatype and shape are human-readable, but the array data itself is binary64 encoded.- Parameters
o (
numpy.ndarray) –ndarrayto be serialized.- Returns
A dictionary that can be passed to
json.dumps.
- dustmaps.json_serializers.serialize_ndarray_npy(o)[source]¶
Serializes a
numpy.ndarrayusing numpy’s built-insavefunction. This produces totally unreadable (and very un-JSON-like) results (in “npy” format), but it’s basically guaranteed to work in 100% of cases.- Parameters
o (
numpy.ndarray) –ndarrayto be serialized.- Returns
A dictionary that can be passed to
json.dumps.
- dustmaps.json_serializers.serialize_ndarray_readable(o)[source]¶
Serializes a
numpy.ndarrayin a human-readable format.- Parameters
o (
numpy.ndarray) –ndarrayto be serialized.- Returns
A dictionary that can be passed to
json.dumps.