dustmaps documentation¶
dustmaps
provides a unified interface for several 2D and 3D
maps of interstellar dust reddening and extinction.
To get started, take a look at Installation and Examples. To see a list of all available maps, take a look at Available Dust Maps. For a complete reference to the API, see dustmap modules.
If you make use of dustmaps
in your research, please cite
Green (2018):
@ARTICLE{2018JOSS....3..695M,
author = {{Green}, {Gregory M.}},
title = "{dustmaps: A Python interface for maps of interstellar dust}",
journal = {The Journal of Open Source Software},
year = "2018",
month = "Jun",
volume = {3},
number = {26},
pages = {695},
doi = {10.21105/joss.00695},
adsurl = {https://ui.adsabs.harvard.edu/abs/2018JOSS....3..695G},
adsnote = {Provided by the SAO/NASA Astrophysics Data System}
}
Contents¶
Installation¶
There are two ways to install dustmaps
.
1. Using pip
¶
From the commandline, run
pip install dustmaps
You may have to use sudo
.
Next, we’ll configure the package and download the dust maps we’ll want to use. Start up a python interpreter and type:
from dustmaps.config import config
config['data_dir'] = '/path/to/store/maps/in'
import dustmaps.sfd
dustmaps.sfd.fetch()
import dustmaps.csfd
dustmaps.csfd.fetch()
import dustmaps.planck
dustmaps.planck.fetch()
import dustmaps.planck
dustmaps.planck.fetch(which='GNILC')
import dustmaps.bayestar
dustmaps.bayestar.fetch()
import dustmaps.iphas
dustmaps.iphas.fetch()
import dustmaps.marshall
dustmaps.marshall.fetch()
import dustmaps.chen2014
dustmaps.chen2014.fetch()
import dustmaps.lenz2017
dustmaps.lenz2017.fetch()
import dustmaps.pg2010
dustmaps.pg2010.fetch()
import dustmaps.leike_ensslin_2019
dustmaps.leike_ensslin_2019.fetch()
import dustmaps.leike2020
dustmaps.leike2020.fetch()
import dustmaps.edenhofer2023
dustmaps.edenhofer2023.fetch()
import dustmaps.gaia_tge
dustmaps.gaia_tge.fetch()
All the dust maps should now be in the path you gave to
config['data_dir']
. Note that these dust maps can be very large - some
are several Gigabytes! Only download those you think you’ll need.
Note that there are two versions of the Bayestar dust map. By default,
dustmaps.bayestar.fetch()
will download Bayestar19 (Green et al. 2019).
In order to download earlier version of the map (Green et al. 2015, 2018), you can
provide the keyword argument version='bayestar2017'
(Green et al. 2018) or
version='bayestar2015'
(Green et al. 2015).
2. Using setup.py
¶
An alternative way to download dustmaps
, if you don’t want to use
pip
, is to download or clone the respository from
https://github.com/gregreen/dustmaps.
In this case, you will have to manually make sure that the dependencies are satisfied:
numpy
scipy
astropy
h5py
healpy
requests
six
progressbar2
These packages can typically be installed using the Python package manager,
pip
.
Once these dependencies are installed, run the following command from the root
directory of the dustmaps
package:
python setup.py install --large-data-dir=/path/to/store/maps/in
Then, fetch the maps you’d like to use. Depending on which dust maps you choose to download, this step can take up several Gigabytes of disk space. Be careful to only download those you think you’ll need:
python setup.py fetch --map-name=sfd
python setup.py fetch --map-name=csfd
python setup.py fetch --map-name=planck
python setup.py fetch --map-name=planckGNILC
python setup.py fetch --map-name=bayestar
python setup.py fetch --map-name=iphas
python setup.py fetch --map-name=marshall
python setup.py fetch --map-name=chen2014
python setup.py fetch --map-name=lenz2017
python setup.py fetch --map-name=leikeensslin2019
python setup.py fetch --map-name=leike2020
python setup.py fetch --map-name=edenhofer2023
python setup.py fetch --map-name=gaia_tge
That’s it!
Note that the above code will download the latest version of the Bayestar dust map (the 2019 version). If you want to download the 2015 and 2017 versions, you can enter the commands
python setup.py fetch --map-name=bayestar2015
python setup.py fetch --map-name=bayestar2017
3. Custom configuration file location (Optional)¶
By default, a configuration file is stored in ~/.dustmapsrc
. This
file might look like the following:
{"data_dir": "/path/to/store/maps/in"}
If you would like dustmaps
to use a different configuration file,
then you can set the environmental variable DUSTMAPS_CONFIG_FNAME
.
For example, in a bash
terminal,
export DUSTMAPS_CONFIG_FNAME=/path/to/custom/config/file.json
python script_using_dustmaps.py
The paths listed in the configuration file can also include environmental
variables, which will be expanded when dustmaps
is loaded. For example,
the configuration file could contain the following:
{"data_dir": "/path/with/${VARIABLE}/included"}
If the environmental variable VARIABLE
is set to "foo"
,
for example, then dustmaps
will expand data_dir
to
"/path/with/foo/included"
.
Examples¶
Getting Started¶
Here, we’ll look up the reddening at a number of different locations on the sky. We specify coordinates on the sky using astropy.coordinates.SkyCoord objects. This allows us a great deal of flexibility in how we specify sky coordinates. We can use different coordinate frames (e.g., Galactic, equatorial, ecliptic), different units (e.g., degrees, radians, hour angles), and either scalar or vector input.
For our first example, let’s load the Schlegel, Finkbeiner & Davis (1998) – or “SFD” – dust reddening map, and then query the reddening at one location on the sky:
from __future__ import print_function
from astropy.coordinates import SkyCoord
from dustmaps.sfd import SFDQuery
coords = SkyCoord('12h30m25.3s', '15d15m58.1s', frame='icrs')
sfd = SFDQuery()
ebv = sfd(coords)
coords = SkyCoord('12h30m25.3s', '15d15m58.1s', frame='icrs')
print('E(B-V) = {:.3f} mag'.format(ebv))
>>> E(B-V) = 0.030 mag
A couple of things to note here:
- In this example, we used
from __future__ import print_function
in order to ensure compatibility with both Python 2 and 3. - Above, we used the
ICRS coordinate system,
by specifying
frame=’icrs’
. SFDQuery
returns reddening in a unit that is similar to magnitudes of E(B-V). However, care should be taken: a unit of SFD reddening is not quite equivalent to a magnitude of E(B-V). The way to correctly convert SFD units to extinction in various broadband filters is to use the conversions in Table 6 of Schlafly & Finkbeiner (2011).
We can query the other maps in the dustmaps
package with only minor
modification to the above code. For example, here’s how we would query the
Planck Collaboration (2013) dust map:
from __future__ import print_function
from astropy.coordinates import SkyCoord
from dustmaps.planck import PlanckQuery
coords = SkyCoord('12h30m25.3s', '15d15m58.1s', frame='icrs')
planck = PlanckQuery()
ebv = planck(coords)
print('E(B-V) = {:.3f} mag'.format(ebv))
>>> E(B-V) = 0.035 mag
Querying Reddening at an Array of Coordinates¶
We can also query an array of coordinates, as follows:
from __future__ import print_function
import numpy as np
from astropy.coordinates import SkyCoord
from dustmaps.planck import PlanckQuery
from dustmaps.sfd import SFDQuery
l = np.array([0., 90., 180.])
b = np.array([15., 0., -15.])
coords = SkyCoord(l, b, unit='deg', frame='galactic')
planck = PlanckQuery()
planck(coords)
>>> array([ 0.50170666, 1.62469053, 0.29259142])
sfd = SFDQuery()
sfd(coords)
>>> array([ 0.55669367, 2.60569382, 0.37351534], dtype=float32)
The input need not be a flat array. It can have any shape – the shape of the output will match the shape of the input:
from __future__ import print_function
import numpy as np
from astropy.coordinates import SkyCoord
from dustmaps.planck import PlanckQuery
l = np.linspace(0., 180., 12)
b = np.zeros(12, dtype='f8')
l.shape = (3, 4)
b.shape = (3, 4)
coords = SkyCoord(l, b, unit='deg', frame='galactic')
planck = PlanckQuery()
ebv = planck(coords)
print(ebv)
>>> [[ 315.52438354 28.11778831 23.53047562 20.72829247]
[ 2.20861101 15.68559361 1.46233201 1.70338535]
[ 0.94013882 1.11140835 0.38023439 0.81017196]]
print(ebv.shape)
>>> (3, 4)
Querying 3D Reddening Maps¶
When querying a 3D dust map, there are two slight complications:
- There is an extra axis – distance – to care about.
- Many 3D dust maps are probabilistic, so we need to specify whether we want the median reddening, mean reddening, a random sample of the reddening, etc.
Let’s see how this works out with the “Bayestar” dust map of Green, Schlafly & Finkbeiner (2015).
How Distances are Handled¶
If we don’t provide distances in our input, dustmaps
will assume we want dust
reddening along the entire line of sight.
from __future__ import print_function
from astropy.coordinates import SkyCoord
from dustmaps.bayestar import BayestarQuery
coords = SkyCoord(180., 0., unit='deg', frame='galactic')
# Note that below, we could use version='bayestar2017' to get the newer
# version of the map. Note, however, that the reddening units are not
# identical in the two versions of the map. See Green et al. (2018) for
# an explanation of the units.
bayestar = BayestarQuery(max_samples=2, version='bayestar2015')
ebv = bayestar(coords, mode='random_sample')
print(ebv)
>>> [ 0.00476 0.00616 0.0073 0.00773 0.00796 0.07453
0.07473 0.0748 0.07807 0.07831 0.18957999 0.2013
0.20448001 0.20734 0.21008 0.73733997 0.75415999 0.93702
0.93956 1.09001005 1.09141004 1.11407995 1.11925006 1.12212002
1.12284994 1.12289 1.12296999 1.12305999 1.12308002 1.12309003
1.12311995]
Here, the Bayestar map has given us a single random sample of the cumulative dust reddening along the entire line of sight – that is, to a set of distances. To see what those distances are, we can call:
bayestar.distances
>>> <Quantity [ 0.06309573, 0.07943282, 0.1 , 0.12589255,
0.15848933, 0.19952621, 0.25118864, 0.31622776,
0.3981072 , 0.50118726, 0.63095725, 0.79432821,
1. , 1.2589252 , 1.58489335, 1.99526215,
2.51188707, 3.1622777 , 3.98107076, 5.01187277,
6.3095727 , 7.94328403, 10. , 12.58925152,
15.84893322, 19.95262146, 25.11886978, 31.62277603,
39.81070709, 50.11872864, 63.09572601] kpc>
The return type is an astropy.unit.Quantity instance, which keeps track of units.
If we provide Bayestar with distances, then it will do the distance interpolation for us, returning the cumulative dust reddening out to specific distances:
import astropy.units as units
coords = SkyCoord(180.*units.deg, 0.*units.deg,
distance=500.*units.pc, frame='galactic')
ebv = bayestar(coords, mode='median')
print(ebv)
>>> 0.10705789
Because we have explicitly told Bayestar what distance to evaluate the map at, it returns only a single value.
How Probability is Handled¶
The Bayestar 3D dust map is probabilistic, meaning that it stores random samples
of how dust reddening could increase along each sightline. Sometimes we might be
interested in the median reddening to a given point in space, or we might want
to have all the samples of reddening out to that point. We specify how we want
to deal with the probabilistic nature of the map by providing the keyword
argument mode
to dustmaps.bayestar.BayestarQuery.__call__
.
For example, if we want all the reddening samples, we invoke:
l = np.array([30., 60., 90.]) * units.deg
b = np.array([10., -10., 15.]) * units.deg
d = np.array([1.5, 0.3, 4.0]) * units.kpc
coords = SkyCoord(l, b, distance=d, frame='galactic')
ebv = bayestar(coords, mode='samples')
print(ebv.shape) # (# of coordinates, # of samples)
>>> (3, 2)
print(ebv)
>>> [[ 0.24641787 0.27142054] # Two samples at the first coordinate
[ 0.01696703 0.0149225 ] # Two samples at the second coordinate
[ 0.08348 0.11068 ]] # Two samples at the third coordinate
If we instead ask for the mean reddening, the shape of the output is different:
ebv = bayestar(coords, mode='mean')
print(ebv.shape) # (# of coordinates)
>>> (3,)
print(ebv)
>>> [ 0.25891921 0.09121627 0.09708 ]
The only axis is for the different coordinates, because we have reduced the samples axis by taking the mean.
In general, the shape of the output from the Bayestar map is:
(coordinate, distance, sample)
where any of the axes can be missing (e.g., if only one coordinate was specified, if distances were provided, or if the median reddening was requested).
Percentiles are handled in much the same way as samples. In the following query, we request the 16th, 50th and 84th percentiles of reddening at each coordinate, using the same coordinates as we generated in the previous example:
ebv = bayestar(coords, mode='percentile', pct=[16., 50., 84.])
print(ebv)
>>> [[ 0.24789949 0.25583497 0.26986977] # Percentiles at 1st coordinate
[ 0.01505572 0.01814967 0.02750403] # Percentiles at 2nd coordinate
[ 0.0860716 0.09787634 0.10787529]] # Percentiles at 3rd coordinate
We can also pass a single percentile:
ebv = bayestar(coords, mode='percentile', pct=25.)
print(ebv)
>>> [ 0.24930404 0.01524667 0.08961 ] # 25th percentile at 3 coordinates
Getting Quality Assurance Flags from the Bayestar Dust Maps¶
For the Bayestar dust maps, one can retrieve QA flags by providing the keyword
argument return_flags=True
:
ebv, flags = bayestar(coords, mode='median', return_flags=True)
print(flags.dtype)
>>> [('converged', '?'), ('reliable_dist', '?')]
print(flags['converged']) # Whether or not fit converged in each pixel
>>> [ True True True]
# Whether or not map is reliable at requested distances
print(flags['reliable_dist'])
>>> [ True False True]
If the coordinates do not include distances, then instead of
'reliable_dist'
, the query will return the minimum and maxmimum reliable
distance moduli of the map in each requested coordinate:
l = np.array([30., 60., 90.]) * units.deg
b = np.array([10., -10., 15.]) * units.deg
coords = SkyCoord(l, b, frame='galactic')
ebv, flags = bayestar(coords, mode='median', return_flags=True)
print(flags['min_reliable_distmod'])
>>> [ 7.875 8.24800014 6.87300014]
print(flags['max_reliable_distmod'])
>>> [ 15.18599987 15.25500011 15.00699997]
We can see from the above that in the previous example, the reason the second coordinate was labeled unreliable was because the requested distance (300 pc) was closer than a distance modulus of 8.248 (corresponding to ~450 pc).
Plotting the Dust Maps¶
We’ll finish by plotting a comparison of the SFD, Planck Collaboration and Bayestar Dust maps. First, we’ll import the necessary modules:
from __future__ import print_function
import matplotlib
import matplotlib.pyplot as plt
import numpy as np
import astropy.units as units
from astropy.coordinates import SkyCoord
from dustmaps.sfd import SFDQuery
from dustmaps.planck import PlanckQuery
from dustmaps.bayestar import BayestarQuery
Next, we’ll set up a grid of coordinates to plot, centered on the Aquila South cloud:
l0, b0 = (37., -16.)
l = np.arange(l0 - 5., l0 + 5., 0.05)
b = np.arange(b0 - 5., b0 + 5., 0.05)
l, b = np.meshgrid(l, b)
coords = SkyCoord(l*units.deg, b*units.deg,
distance=1.*units.kpc, frame='galactic')
Then, we’ll load up and query three different dust maps:
sfd = SFDQuery()
Av_sfd = 2.742 * sfd(coords)
planck = PlanckQuery()
Av_planck = 3.1 * planck(coords)
bayestar = BayestarQuery(max_samples=1)
Av_bayestar = 2.742 * bayestar(coords)
We’ve assumed \(R_V = 3.1\), and used the coefficient from Table 6 of Schlafly & Finkbeiner (2011) to convert SFD and Bayestar reddenings to magnitudes of \(A_V\).
Finally, we create the figure using matplotlib
:
fig = plt.figure(figsize=(12,4), dpi=150)
for k,(Av,title) in enumerate([(Av_sfd, 'SFD'),
(Av_planck, 'Planck'),
(Av_bayestar, 'Bayestar')]):
ax = fig.add_subplot(1,3,k+1)
ax.imshow(
np.sqrt(Av)[::,::-1],
vmin=0.,
vmax=2.,
origin='lower',
interpolation='nearest',
cmap='binary',
aspect='equal'
)
ax.axis('off')
ax.set_title(title)
fig.subplots_adjust(wspace=0., hspace=0.)
plt.savefig('comparison.png', dpi=150)
Here’s the result:

Querying the web server¶
Some of the maps included in this package are large, and can take up a lot of memory, or be slow to load. To make it easier to work with these maps, some of them are available to query over the internet. As of now, the following maps can be queried remotely:
- Bayestar (all versions)
- SFD
The API for querying these maps remotely is almost identical to the API for local queries. For example, the following code queries SFD remotely:
from __future__ import print_function
from astropy.coordinates import SkyCoord
from dustmaps.sfd import SFDWebQuery
l = [180., 160.]
b = [30., 45.]
coords = SkyCoord(l, b, unit='deg', frame='galactic')
sfd = SFDWebQuery()
ebv = sfd(coords)
print(ebv)
>>> [0.04704102 0.02022794]
The following example queries the Bayestar2019 dust map remotely. The web interface takes the same arguments as the local interface:
import astropy.units as u
from dustmaps.bayestar import BayestarWebQuery
l = [90., 150., 35.] * u.deg
b = [10., 12., -25.] * u.deg
d = [500., 3500., 1000.] * u.pc
coords = SkyCoord(l, b, distance=d, frame='galactic')
q = BayestarWebQuery(version='bayestar2019')
E = q(coords, mode='median')
print(E)
>>> [0.13 0.63 0.09999999]
The query_gal()
and query_equ()
convenience functions also
work with web queries. Continuing from the previous example,
E = q.query_gal([120., 125.], [-5., -10.],
d=[1.5, 1.3],
mode='random_sample')
print(E)
>>> [0.32 0.24]
Please take it easy on our web server. If you want to query multiple coordinates, then bundle them up into one query. If you want to query a very large number of coordinates, consider downloading the maps and querying them locally instead.
Available Dust Maps¶
Two-Dimensional Dust Maps¶
SFD¶
A two-dimensional map of dust reddening across the entire sky. The “SFD” dust map is based on far-infrared emission of dust. The authors model the temperature and optical depth of the dust, and then calibrate a relationship between the dust’s far-infrared optical depth and optical reddening. This calibration was updated by Schlafly & Finkbeiner (2011).
In order to convert SFD values of E(B-V) to extinction, one should use the conversions provided in Table 6 of Schlafly & Finkbeiner (2011).
- Reference: Schlegel, Finkbeiner & Davis (1998)
- Recalibration: Schlafly & Finkbeiner (2011)
CSFD (Chiang 2023)¶
“Corrected SFD,” a 2D dust map based on a reanalysis of SFD, using tomographically constrained templates from WISE galaxy density fields to remove extragalactic contamination from the cosmic infrared background (CIB).
- Reference: Chiang (2023)
- Website: Project description
- Data: Zenodo
Gaia Total Galactic Extinction (2022)¶
A two-dimensional map of A0, the monochromatic extinction at 541.4 nm. The map is based on extinction estimates for giants beyond 300 pc. The individual exitnction estimates estimates were obtained by fitting Gaia BP/RP spectra, parallaxes and G-band apparent magnitudes.
The map comes in multiple HEALPix levels (6 to 9). By default, an “optimum” map is loaded, with an adaptive HEALPix level, based on the local number of stars (at least 3 stars are required per pixel).
- Reference: Delchambre et al. (2022).
Lenz, Hensley & Doré (2017)¶
A two-dimensional map of dust reddening, covering 40% of the sky with a 16.1’ resolution. This map is derived from emission from low-velocity (l.o.s. velocity < 90 km/s) HI, which is found to correlate much more strongly with E(B-V) than emission from high-velocity HI. The underlying data comes from the HI4PI Survey. This map reports E(B-V) in magnitudes.
- Reference: Lenz, Hensley & Doré (2017).
- See also: GitHub page.
Planck Collaboration (2013)¶
Two-dimensional maps of dust reddening across the entire sky. The Planck Collaboration (2013) fits a modified blackbody dust emission model to the Planck and IRAS far-infrared maps, and provides three different conversions to dust reddening.
The three maps provided by Planck Collaboration (2013) are based on:
- τ353: dust optical depth at 353 GHz.
- ℛ: thermal dust radiance.
- A recommended extragalactic reddening estimate, based on thermal dust radiance, but with point sources removed.
- Reference: Planck Collaboration (2013)
- Website: Planck Explanatory Supplement
Planck Collaboration (2016; “GNILC”)¶
Two-dimensional maps of dust reddening across the entire sky, using the generalized needlet internal linear combination (GNILC) method to separate out Galactic dust emission from CIB anisotropies.
This map contains both reddening estimates and estimated uncertainties.
- Reference: Planck Collaboration (2016)
- Website: Planck Explanatory Supplement
Peek & Graves (2010)¶
A correction to the SFD’98 dust map, based on color excess measurements of “standard crayons” – spectroscopically selected passively evolving galaxies. The maps have an angular resolution of 4.5°, and have a 1σ uncertainty of 1.5 mmag in E(B-V). Subtract this map from SFD’98 to obtain the corrected E(B-V) reddening.
- Reference: Peek & Graves (2010)
Burstein & Heiles¶
Primarily of historical interest, the Burstein & Heiles (1982) dust reddening maps are derived from HI column density and galaxy counts.
- Reference: Burstein & Heiles (1982)
Three-Dimensional Dust Maps¶
Bayestar¶
A three-dimensional map of Milky Way dust reddening, covering the three quarters of the sky north of a declination of -30°. The map is probabilistic. containing samples of the reddening along each line of sight. The “Bayestar” dust map is inferred from stellar photometry of 800 million stars observed by Pan-STARRS 1, and 2MASS photometry for a quarter of the stars. The latest version of Bayestar also makes use of Gaia DR2 parallaxes.
There are three versions of Bayestar, called Bayestar19, Bayestar17 and
Bayestar15 here. By default, dustmaps
will use the latest version,
Bayestar19, although the earlier versions of the map can be selected by providing
the keyword argument version='bayestar2017'
or version='bayestar2015'
in routines such as dustmaps.bayestar.fetch
,
dustmaps.bayestar.BayestarQuery
and dustmaps.bayestar.BayestarWebQuery
.
If you want to make sure that your code will always use the same version of the
map, even as new versions of Bayestar are released, then set the version
keyword explicitly.
The units of reddening used by each map are slightly different:
- Bayestar19 reports reddening in an arbitrary unit that can be converted to extinction in different bands using the coefficients given in Table 1 of Green, Schlafly, Finkbeiner et al. (2019).
- Bayestar17 reports reddening in an arbitrary unit that can be converted to extinction in different bands using the coefficients given in Table 1 of Green, Schlafly, Finkbeiner et al. (2018).
- Bayestar15 reports reddening in the same units as those used by SFD. Therefore, in order to convert Bayestar15 reddenings to extinction in different bands, one should use the conversions provided in Table 6 of Schlafly & Finkbeiner (2011).
Chen et al. (2014)¶
A three-dimensional map of dust extinction in the Galactic anticenter. The map covers about 6000 deg2, from 140° < ℓ < 240° and -60° < b < 40°, and is based on stellar photometry from the Xuyi Schmidt Telescope Photometric Survey of the Galactic Anticentre (XSTPS-GAC), 6MASS and WISE. The map has an angular resolution of 3 to 9 arcminutes, and reports r-band extinction, along with Gaussian error estimates.
- Reference: Chen et al. (2014)
- Website: http://lamost973.pku.edu.cn
Edenhofer et al. (2023)¶
A three-dimensional map of Milky Way dust extinction, with a Gaussian process prior on the logarithm of the dust extinction density. The prior is implemented on a spherical grid. The map starts at 69 pc and extends out to 1.25 kpc in distance from the Sun. It has an angular resolution of 14’ and a maximum distance voxalization of 0.4 pc at 69 pc and a minimum distance voxalization of 7 pc at 1250 pc. The map is based on the stellar distance and extinction estimates of Zhang, Green & Rix (2023), and therefore reports extinctions in their units. Accompanying the main reconstruciton is an additional map that uses less data but extends out to 2 kpc from the Sun.
- Reference: Edenhofer et al. (2023)
- Data: Zenodo
IPHAS¶
A three-dimensional map of Milky Way dust extinction, covering a 10°-thick strip of the Galactic plane, between 30° < ℓ < 120°. The map is probabilistic, containing samples of the cumulative extinction along each line of sight. The map is based on IPHAS imaging of stars. The map returns A0, the monochromatic extinction.
- Reference: Sale et al. (2014)
- Website: www.iphas.org/extinction
Leike & Enßlin (2019)¶
A three-dimensional map of Milky Way dust extinction, incorporating a Gaussian process prior on the log of the dust extinction density. The map is based on the Gaia DR2 catalog parallaxes and G-band extinctions, and spans a (600 pc)³ box centered on the Sun.
- Reference: Leike & Enßlin (2019)
- Data: Zenodo
Leike, Glatzle & Enßlin (2020)¶
A three-dimensional map of Milky Way dust extinction, incorporating a Gaussian process prior on the log of the dust extinction density, similar to Leike & Enßlin (2019). The map is based on data from Gaia, 2MASS, Pan-STARRS 1 and ALLWISE, and is calculated on a Cartesian grid spanning a (740 pc)×(740 pc)×(540 pc) box (in Galactic x, y and z, respectively) centered on the Sun.
- References: Leike, Glatzle & Enßlin (2020)
- Data: Zenodo
Marshall et al. (2006)¶
A three-dimensional map of Milky Way dust extinction, covering a 20°-thick strip of the Galactic plane, between -100° < ℓ < 100°. The map is contains 2MASS Ks-band extinctions with a Gaussian uncertainty estimates. The map is based on a comparison of 2MASS colors of stars with expectations from the Besançon model of the Galaxy.
- Reference: Marshall et al. (2006)
- Website: http://cds.u-strasbg.fr/
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.DustMap
Queries 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'
.
- map_fname (Optional[
-
distances
¶ Returns the distance bin edges that the map uses. The return type is
astropy.units.Quantity
, which stores unit-full quantities.
-
distmods
¶ Returns the distance modulus bin edges that the map uses. The return type is
astropy.units.Quantity
, with units of mags.
-
query
(coords, **kwargs)[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’. Themode
determines 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`, whereret
is the normal return value, containing reddening. Defaults toFalse
. - pct (Optional[
float
or list/array offloat
]) – If the mode ispercentile
, thenpct
specifies 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 theBayestarQuery
object 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 whethercoords
contains distances.If
coords
does not specify distance(s), then the shape of the output begins withcoords.shape
. Ifcoords
does specify distance(s), then the shape of the output begins withcoords.shape + ([number of distance bins],)
.If
mode
is'random_sample'
, then at each coordinate/distance, a random sample of reddening is given.If
mode
is'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
mode
is'median'
, then at each coordinate/distance, the median reddening is returned.If
mode
is'mean'
, then at each coordinate/distance, the mean reddening is returned.If
mode
is'best'
, then at each coordinate/distance, the maximum posterior density reddening is returned (the “best fit”).If
mode
is'percentile'
, then an additional keyword argument,pct
, must be specified. At each coordinate/distance, the requested percentiles (inpct
) will be returned. Ifpct
is a list/array, then the last axis of the output will correspond to different percentiles.Finally, if
mode
is'samples'
, then at each coordinate/distance, all samples are returned. The last axis of the output will correspond to different samples.If
return_flags
isTrue
, 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.- coords (
-
-
class
dustmaps.bayestar.
BayestarWebQuery
(api_url=None, version='bayestar2019')[source]¶ Bases:
dustmaps.map_base.WebDustMap
Remote 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 HEALPixnside
parameter. - 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
l
andb
.- nside (
bh (Burstein & Heiles 1982)¶
-
class
dustmaps.bh.
BHQuery
(bh_dir=None)[source]¶ Bases:
dustmaps.map_base.DustMap
Queries 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, **kwargs)[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.UnstructuredDustMap
The 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.
-
distances
¶ Returns the distance bins that the map uses. The return type is
astropy.units.Quantity
, which stores unit-full quantities.
-
query
(coords, **kwargs)[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
coords
contains distances.If
coords
does not specify distance(s), then the shape of the output begins withcoords.shape
. Ifcoords
does specify distance(s), then the shape of the output begins withcoords.shape + ([number of distance bins],)
.- coords (
-
-
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.
- dat_fname (
-
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.
csfd (Chiang 2023)¶
-
class
dustmaps.csfd.
CSFDQuery
(map_fname=None, mask_fname=None)[source]¶ Bases:
dustmaps.healpix_map.HEALPixQuery
Queries the Corrected SFD dust map of Chiang (2023). This map is based on SFD, but contains a correction to remove contamination from large-scale structure (i.e., external galaxies).
-
__init__
(map_fname=None, mask_fname=None)[source]¶ Parameters: - map_fname (Optional[
str
]) – Filename of the CSFD EBV map. Defaults to`None
, meaning that the default location is used. - mask_fname (Optional[
str
]) – Filename of the CSFD mask map. Defaults to`None
, meaning that the default location is used.
- map_fname (Optional[
-
query
(coords, **kwargs)[source]¶ Returns CSFD reddening on the same scale as SFD (similar to E(B-V)) at the specified location(s) on the sky. Also optionally returns a bit mask, where the bits (ordered from least to most significant) have the following meanings:
Bit 0: 'LSS_corr' - This bit is set in the footprint within which the LSS is reconstructed, and CSFD = SFD - LSS (otherwise CSFD = SFD). Bit 1: 'no_IRAS' - Set in the area with no IRAS data (DIRBE data filled in SFD); LSS removal in CSFD is done using a 1 deg smoothed LSS. Bit 2: 'cosmology' - Set in the area where both the LSS and CSFD are most reliable for precision cosmology analyses.
Parameters: - coords (
astropy.coordinates.SkyCoord
) – The coordinates to query. - return_flags (Optional[
bool
]) – IfTrue
, then a bit mask is returned as well, indicating where CSFD has been corrected for large-scale structure, where IRAS data was used, and where the map is suitable for cosmology. See above description of bits. Defaults toFalse
.
Returns: A float array of the reddening, at the given coordinates. The shape of the output is the same as the shape of the input coordinate array,
coords
. Ifreturn_flags
isTrue
, a second array (a bit mask) of the same shape is returned. See above description of the meaning of each bit.- coords (
-
-
dustmaps.csfd.
fetch
(clobber=False)[source]¶ Downloads the Corrected SFD dust map of Chiang (2023).
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.
gaia_tge (Delchambre et al. 2022)¶
-
class
dustmaps.gaia_tge.
GaiaTGEQuery
(map_fname=None, healpix_level='optimum')[source]¶ Bases:
dustmaps.healpix_map.HEALPixQuery
Queries the Gaia Total Galactic Extinction (Delchambre 2022) dust map, which contains estimates of monochromatic extinction, A0, in mags.
-
__init__
(map_fname=None, healpix_level='optimum')[source]¶ Parameters: - map_fname (Optional[str]) – Filename of the Gaia TGE map.
Defaults to
None
, meaning that the default location is used. - healpix_level (Optional[int or str]) – Which HEALPix level to load into the map. If “optimum” (the default), loads the optimum HEALPix level available at each location. If an int, instead loads the specified HEALPix level.
- map_fname (Optional[str]) – Filename of the Gaia TGE map.
Defaults to
-
query
(coords, **kwargs)[source]¶ Returns a numpy array containing A0 at the specified location(s) on the sky. Optionally, returns a 2nd array containing flags at the same location(s).
Parameters: - coords (astropy.coordinates.SkyCoord) – The coordinates to query.
- return_flags (Optional[bool]) – If True, returns a 2nd array containing flags at each coordinate. Defaults to False.
Returns: A numpy array containing A0 at the specified coordinates. The shape of the output is the same as the shape of the input coordinate array,
coords
. If return_flags is True, a 2nd record array containing flags at each coordinate is also returned.
-
iphas (Sale et al. 2014)¶
-
class
dustmaps.iphas.
IPHASQuery
(map_fname=None)[source]¶ Bases:
dustmaps.unstructured_map.UnstructuredDustMap
The 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.
-
distances
¶ Returns the distance bins that the map uses. The return type is
astropy.units.Quantity
, which stores unit-full quantities.
-
query
(coords, **kwargs)[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'
. Themode
determines 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 whethercoords
contains distances.If
coords
does 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
mode
is'random_sample'
, then at each coordinate/distance, a random sample of reddening is given.If
mode
is'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
mode
is'median'
, then at each coordinate/distance, the median reddening is returned.If
mode
is'mean'
, then at each coordinate/distance, the mean reddening is returned.Finally, if
mode
is'samples'
, then all at each coordinate/distance, all samples are returned.- coords (
-
-
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.
leike_ensslin_2019 (Leike & Enßlin 2019)¶
-
class
dustmaps.leike_ensslin_2019.
LeikeEnsslin2019Query
(map_fname=None)[source]¶ Bases:
dustmaps.map_base.DustMap
A class for querying the Leike & Ensslin (2019) dust map.
-
__init__
(map_fname=None)[source]¶ Parameters: map_fname (Optional[str]) – Filename of the map. Defaults to None
, meaning that the default location is used.
-
query
(coords, **kwargs)[source]¶ Returns the extinction density (in e-foldings / kpc, in Gaia G-band) at the given coordinates.
Parameters: - coords (
astropy.coordinates.SkyCoord
) – Coordinates at which to query the extinction. Must be 3D (i.e., include distance information). - component (str) – Which component to return. Allowable values are ‘mean’ (for the mean extinction density) and ‘std’ (for the standard deviation of extinction density). Defaults to ‘mean’.
Returns: The extinction density, in units of e-foldings / kpc, as either a numpy array or float, with the same shape as the input
coords
.- coords (
-
-
dustmaps.leike_ensslin_2019.
fetch
(clobber=False)[source]¶ Downloads the 3D dust map of Leike & Ensslin (2019).
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.
leike2020 (Leike, Glatzle & Enßlin 2020)¶
-
class
dustmaps.leike2020.
Leike2020Query
(map_fname=None)[source]¶ Bases:
dustmaps.map_base.DustMap
A class for querying the Leike, Glatzle & Ensslin (2020) dust map.
For details on how to use this map, see the original paper: https://ui.adsabs.harvard.edu/abs/2020A%26A…639A.138L/abstract.
The data is deposited at Zenodo: https://doi.org/10.5281/zenodo.3993082.
-
__init__
(map_fname=None)[source]¶ Parameters: map_fname (Optional[str]) – Filename of the map. Defaults to None
, meaning that the default location is used.
-
query
(coords, **kwargs)[source]¶ Returns the extinction density (in e-foldings / kpc, in Gaia G-band) at the given coordinates.
Parameters: - coords (
astropy.coordinates.SkyCoord
) – Coordinates at which to query the extinction. Must be 3D (i.e., include distance information). - component (str) – Which component to return. Allowable values are ‘mean’ (for the mean extinction density) and ‘std’ (for the standard deviation of extinction density). Defaults to ‘mean’.
Returns: The extinction density, in units of e-foldings / kpc, as either a numpy array or float, with the same shape as the input
coords
.- coords (
-
-
dustmaps.leike2020.
fetch
(clobber=False, fetch_samples=False)[source]¶ Downloads the 3D dust map of Leike & Ensslin (2020).
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. - fetch_samples (Optional[bool]) – If
True
, the samples will also be downloaded. IfFalse
(the default), only the mean and standard deviation will be downloaded. The samples take up 14 GB, which is why the default is not to download them.
- clobber (Optional[bool]) – If
edenhofer2023 (Edenhofer et al. 2023)¶
-
class
dustmaps.edenhofer2023.
Edenhofer2023Query
(map_fname=None, load_samples=False, integrated=False, flavor='main', seed=None)[source]¶ Bases:
dustmaps.map_base.DustMap
A class for querying the Edenhofer et al. (2023) 3D dust map.
The map is in units of E of Zhang, Green, and Rix (2023) per parsec but can be translated to an extinction at any given wavelength by using the extinction curve published at https://doi.org/10.5281/zenodo.6674521 . For further details on the map, see the original paper.
If you use this map in your work, please cite Edenhofer et al. (2023).
The data is deposited at https://doi.org/10.5281/zenodo.8187943.
-
__init__
(map_fname=None, load_samples=False, integrated=False, flavor='main', seed=None)[source]¶ Parameters: - map_fname (Optional[str]) – Filename of the map. Defaults
to
None
, meaning that the default location is used. - load_samples (Optional[bool]) – Whether to load the posterior samples
of the extinction density. The samples give more accurate
interpolation resoluts and are required for standard deviations
of integrated extinctions. Defaults to
False
. - integrated (Optional[bool]) – Whether to return integrated extinction
density. In this case, the units of the map are E of Zhang,
Green, and Rix (2023). The pre-processing for efficient access
to the integrated extinction map can take a couple of minutes
but subsequent queries will be fast. Defaults to
False
. - flavor (Optional[str]) – Flavor of the map to use. Must be in (‘main’, ‘less_data_but_2kpc’).
- seed (Optional[int]) – A random seed, used when drawing random samples from the map. Set this seed in order to make the pseudo-random draws reproducible.
- map_fname (Optional[str]) – Filename of the map. Defaults
to
-
distance_bounds
¶ Returns the distance bin edges that the map uses. The return type is
astropy.units.Quantity
, which stores unit-full quantities.
-
distances
¶ Returns the distance bin centers that the map uses. The return type is
astropy.units.Quantity
, which stores unit-full quantities.
-
flavor
¶ Returns the flavor of the map. This is a string that is either
"main"
or"less_data_but_2kpc"
.
-
integrated
¶ Returns
True
if the map contains integrated extinction, orFalse
if it contains extinction density.
-
n_samples
¶ Returns the number of samples stored in the map. If no samples have been loaded, then
None
is returned.
-
query
(coords, **kwargs)[source]¶ Returns the 3D dust extinction density or integrated extinction (depending on how the class was initialized) from Edenhofer et al. (2023) at the given coordinates. The map is in units of E of Zhang, Green, and Rix (2023) per parsec or if integrated simply in units of E.
Parameters: - coords (
astropy.coordinates.SkyCoord
) – Coordinates at which to query the extinction. Must be 3D (i.e., include distance information). - mode (str) – Which mode to return. Allowed values are ‘mean’ (for the mean), ‘std’ (for the standard deviation), ‘samples’ (for all posterior samples), or ‘random_sample’ (for a single sample). Defaults to ‘mean’.
Notes
To query integrated extinction, set integrated=True during initizalization.
Returns: Depending on how the class was initialized, either extinction density or extinction are returned. The output is either a numpy array or float, with the same shape as the input coords
.- coords (
-
-
dustmaps.edenhofer2023.
fetch
(clobber=False, fetch_samples=False, fetch_2kpc=False)[source]¶ Downloads the 3D dust map of Edenhofer et al. (2023).
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. - fetch_samples (Optional[bool]) – If
True
, the samples will also be downloaded. IfFalse
(the default), only the mean and standard deviation will be downloaded taking up about 3.2 GB in size while the samples take up 19 GB. - fetch_2kpc (Optional[bool]) – If
True
, the validation run using less data though which extends out to 2kpc in distance will also be downloaded.
- clobber (Optional[bool]) – If
lenz2017 (Lenz, Hensley & Doré 2017)¶
-
class
dustmaps.lenz2017.
Lenz2017Query
(map_fname=None)[source]¶ Bases:
dustmaps.healpix_map.HEALPixFITSQuery
Queries 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.DustMap
Galactic-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, **kwargs)[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_sigma
isTrue
, then the uncertainty in reddening is also returned, so that the output is(A, sigma_A)
, where bothA
andsigma_A
have the same shape as the input coordinates.- coords (
-
-
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.SFDBase
Queries 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 thatdustmaps
will 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.
- map_dir (Optional[
-
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
.- coords (
-
planck (Planck Collaboration 2013, 2016)¶
-
class
dustmaps.planck.
PlanckGNILCQuery
(map_fname=None, load_errors=False)[source]¶ Bases:
dustmaps.healpix_map.HEALPixFITSQuery
Queries the Planck Collaboration (2016) GNILC dust map.
-
__init__
(map_fname=None, load_errors=False)[source]¶ Parameters: - map_fname (Optional[
str
]) – Filename of the Planck map. Defaults to`None
, meaning that the default location is used. - load_errors (Optional[
str
]) – IfTrue
, then the error estimates will be loaded as well, and returned with any query. IfFalse
(the default), then queries will only return the the reddening estimate, without any error estimate.
- map_fname (Optional[
-
query
(coords, **kwargs)[source]¶ Returns E(B-V) at the specified location(s) on the sky.
Parameters: coords ( astropy.coordinates.SkyCoord
) – The coordinates to query.Returns: If the error estimates have been loaded, then a structured array containing 'EBV'
and'EBV_err'
is returned. Otherwise, returns a float array of E(B-V), at the given coordinates. The shape of the output is the same as the shape of the input coordinate array,coords
.
-
-
class
dustmaps.planck.
PlanckQuery
(map_fname=None, component='extragalactic')[source]¶ Bases:
dustmaps.healpix_map.HEALPixFITSQuery
Queries 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'
.
- map_fname (Optional[
-
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 input coordinate array, 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.Quantity
is returned, with units of Kelvin. If beta (or beta error) was chosen, then the output is unitless.
-
-
dustmaps.planck.
fetch
(which='2013')[source]¶ Downloads the Planck dust maps, placing them in the default
dustmaps
directory. There are two different Planck dust maps that can be downloaded: the Planck Collaboration (2013) map and the “GNILC” (Planck Collaboration 2016) map.Parameters: which (Optional[ str
]) – The name of the dust map to download. Should be either2013
(the default) orGNILC
.
sfd (Schlegel, Finkbeiner & Davis 1998)¶
-
class
dustmaps.sfd.
SFDBase
(base_fname)[source]¶ Bases:
dustmaps.map_base.DustMap
Queries 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'
, whereX
is'ngp'
and'sgp'
.
-
query
(coords, **kwargs)[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.SFDBase
Queries 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.WebDustMap
Remote 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
WebDustMap
object.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'
.
- api_url (Optional[
-
fetch_utils¶
-
exception
dustmaps.fetch_utils.
DownloadError
[source]¶ Bases:
dustmaps.fetch_utils.Error
An exception that occurs while trying to download a file.
-
exception
dustmaps.fetch_utils.
Error
[source]¶ Bases:
exceptions.Exception
-
__weakref__
¶ list of weak references to the object (if defined)
-
-
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
True
if 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: True
if the file matches by all given criteria.
map_base¶
-
class
dustmaps.map_base.
DustMap
[source]¶ Bases:
object
Base 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.coordinates
may 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.
- ra (
-
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.
- l (
-
-
class
dustmaps.map_base.
WebDustMap
(api_url=None, map_name='')[source]¶ Bases:
object
Base 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
(*args, **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.coordinates
may 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.
- ra (
-
query_gal
(*args, **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.
- l (
-
-
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.False
for 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.
- coords (
-
dustmaps.map_base.
ensure_coord_type
(f)[source]¶ A decorator for class methods of the form
Class.method(self, coords, **kwargs)
where
coords
is anastropy.coordinates.SkyCoord
object.The decorator raises a
TypeError
if thecoords
that gets passed toClass.method
is not anastropy.coordinates.SkyCoord
instance.Parameters: f (class method) – A function with the signature (self, coords, **kwargs)
, wherecoords
is aSkyCoord
object containing an array.Returns: A function that raises a TypeError
ifcoords
is not anastropy.coordinates.SkyCoord
object, 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
coords
is anastropy.coordinates.SkyCoord
object.The decorator ensures that the
coords
that gets passed toClass.method
is a flat array. It also reshapes the output ofClass.method
to have the same shape (possibly scalar) as the inputcoords
. If the output ofClass.method
is 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)
, wherecoords
is aSkyCoord
object containing an array.Returns: A function that takes SkyCoord
input 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
coords
is anastropy.coordinates.SkyCoord
object.The decorator ensures that the
coords
that gets passed toClass.method
is a flat array of Galactic coordinates. It also reshapes the output ofClass.method
to have the same shape (possibly scalar) as the inputcoords
. If the output ofClass.method
is 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)
, wherecoords
is aSkyCoord
object containing an array.Returns: A function that takes SkyCoord
input with any shape (including scalar).
healpix_map¶
-
class
dustmaps.healpix_map.
HEALPixFITSQuery
(fname, coord_frame, hdu=0, field=None, dtype='f8', scale=None)[source]¶ Bases:
dustmaps.healpix_map.HEALPixQuery
A HEALPix map class that is initialized from a FITS file.
-
__init__
(fname, coord_frame, hdu=0, field=None, dtype='f8', scale=None)[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
astropy
understands. - 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,
including a structured datatype, if multiple fields are to be
loaded. Defaults to
'f8'
, for IEEE754 double precision. - scale (Optional[
float
]) – Scale factor to be multiplied into the data.
-
-
class
dustmaps.healpix_map.
HEALPixQuery
(pix_val, nest, coord_frame, flags=None)[source]¶ Bases:
dustmaps.map_base.DustMap
A class for querying HEALPix maps.
-
__init__
(pix_val, nest, coord_frame, flags=None)[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.
-
query
(coords, return_flags=False)[source]¶ Parameters: - coords (astropy.coordinates.SkyCoord) – The coordinates to query.
- return_flags ([Optional[
bool
]) – If True, return flags at each pixel. Only possible if flags were provided during initialization.
Returns: A float array of the value of the map at the given coordinates. The shape of the output is the same as the shape of the coordinates stored by coords. If return_flags is True, then a second array, containing flags at each pixel, is also returned.
-
unstructured_map¶
-
class
dustmaps.unstructured_map.
UnstructuredDustMap
(pix_coords, max_pix_scale, metric_p=2, frame=None)[source]¶ Bases:
dustmaps.map_base.DustMap
A 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
.
- pix_coords (array-like
-
config¶
-
class
dustmaps.config.
Configuration
(fname)[source]¶ Bases:
object
A class that stores the package configuration.
By default, the configuration is loaded from
~/.dustmapsrcThis can be overridden by setting the environmental variable
DUSTMAPS_CONFIG_FNAME
.Paths stored in the configuration file (such as the data directory,
data_dir
, can include environmental variables, which will be expanded.-
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 force
isFalse
.
-
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.ndarray
numpy.dtype
astropy.units.Quantity
astropy.coordinates.SkyCoord
-
__init__
(*args, **kwargs)[source]¶ encoding
determines the encoding used to interpret anystr
objects decoded by this instance (utf-8 by default). It has no effect when decodingunicode
objects.Note that currently only encodings that are a superset of ASCII work, strings of other encodings should be passed in as
unicode
.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_hook
will be used instead of thedict
. This feature can be used to implement custom decoders that rely on the order that the key and value pairs are decoded (for example, collections.OrderedDict will remember the order of insertion). Ifobject_hook
is also defined, theobject_pairs_hook
takes 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
strict
is 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 adtype
object.Returns: A dtype
object.
-
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 ndarray
object.Returns: An ndarray
object.
-
dustmaps.json_serializers.
deserialize_ndarray_npy
(d)[source]¶ Deserializes a JSONified
numpy.ndarray
that was created using numpy’ssave
function.Parameters: d ( dict
) – A dictionary representation of anndarray
object, created usingnumpy.save
.Returns: An ndarray
object.
-
dustmaps.json_serializers.
deserialize_quantity
(d)[source]¶ Deserializes a JSONified
astropy.units.Quantity
.Parameters: d ( dict
) – A dictionary representation of aQuantity
object.Returns: A Quantity
object.
-
dustmaps.json_serializers.
deserialize_skycoord
(d)[source]¶ Deserializes a JSONified
astropy.coordinates.SkyCoord
.Parameters: d ( dict
) – A dictionary representation of aSkyCoord
object.Returns: A SkyCoord
object.
-
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.ndarray
numpy.floating
(converted tofloat
)numpy.integer
(converted toint
)numpy.dtype
astropy.units.Quantity
astropy.coordinates.SkyCoord
Parameters: ndarray_mode (Optional[ str
]) – Which method to use to serializenumpy.ndarray
objects. 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-insave
function 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
) –dtype
to be serialized.Returns: A dictionary that can be passed to json.dumps
.
-
dustmaps.json_serializers.
serialize_ndarray_b64
(o)[source]¶ Serializes a
numpy.ndarray
in a format where the datatype and shape are human-readable, but the array data itself is binary64 encoded.Parameters: o ( numpy.ndarray
) –ndarray
to be serialized.Returns: A dictionary that can be passed to json.dumps
.
-
dustmaps.json_serializers.
serialize_ndarray_npy
(o)[source]¶ Serializes a
numpy.ndarray
using numpy’s built-insave
function. 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
) –ndarray
to be serialized.Returns: A dictionary that can be passed to json.dumps
.
-
dustmaps.json_serializers.
serialize_ndarray_readable
(o)[source]¶ Serializes a
numpy.ndarray
in a human-readable format.Parameters: o ( numpy.ndarray
) –ndarray
to be serialized.Returns: A dictionary that can be passed to json.dumps
.
License¶
The dustmaps
documentation is covered by the MIT License, as given below.
The MIT License (MIT)¶
Copyright (c) 2016 Gregory M. Green
Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the “Software”), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:
The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.
THE SOFTWARE IS PROVIDED “AS IS”, WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.