deepszsim package
Documentation for DeepSZSim.
To create a simulation of two SZ clusters as quickly as possible, run
tc = deepszsim.simulate_clusters(num_halos=2).get_T_maps(), which will return a dictionary with two entries, accessed by tc.clusters. Each of these entries includes a params subdirectory and a maps subdirectory.
Submodules
deepszsim.dm_halo_dist module
creates a mass and redshift distribution of halos
- deepszsim.dm_halo_dist.flatdist_halo(zmin, zmax, m200min_SM, m200max_SM, size, seed=None)
Creates a random uniform distribution of redshifts and masses for use in creating simulations
Parameters:
- zmin: float
minimum value of the redshift distribution
- zmax: float
maximum value of the redshift distribution
- m200min_SM: float
minimum value of the mass distribution, in units of solar masses
- m200max_SM: float
maximum value of the mass distribution, in units of solar masses
- size: int
size of the distribution
- seed: int
seed for random number generation
Returns:
- zdist: float array
distribution of random uniform redshifts starting at zmin ending at zmax with size size
- mdist: float array
distribution of random uniform redshifts starting at m500min_SM ending at m500max_SM with size size
deepszsim.filters module
submap filtering tools to analyze SZ submaps
- deepszsim.filters.get_tSZ_signal_aperture_photometry(dT_map, radmax_arcmin, pixel_scale, fmax=np.float64(1.4142135623730951))
Retrieve aperture photometry of the tSZ signal
Parameters:
- dT_map: float array
map in uK
- radmax_arcmin: float
radius of the inner aperture in arcmin
- pixel_scale: float
arcmin per pixel for the current settings
- fmax: float
fmax * radmax_arcmin is the radius of the outer aperture in arcmin
Returns:
- disk_mean: float
average value within an annulus of inner radius r
- ring_mean: float
average value within an annulus of outer radius sqrt(2)*r
- tSZ signal: float
thermal SZ effect signal
deepszsim.io_params module
CAMB and yaml-handling tools
- class deepszsim.io_params.config_obj(user_config='/home/docs/checkouts/readthedocs.org/user_builds/deepszsim/checkouts/latest/deepszsim/settings/user_config.yaml')
Bases:
objectConfiguration object that is used to obtain power spectra
Attributes:
- CAMBparams: CAMBparams instance
CAMB object; returns results from CAMB
- UserParams: dict
dictionary of values that the user has specified (smaller than the corresponding dictionary that would be necessary to fully specify a CAMBparams instance)
- dict_iterables: dict
dictionary of all of the iterables that the user has specified, which will be made available to loop over in camb_power_spectrum.CAMBPowerSpectrum
- cosmology_param(ref)
Set cosmology parameters
Parameters:
- ref: dict
dictionary containing cosmology parameters
Returns:
- (cosmo, sigma8, ns): tuple
tuple of (cosmology in astropy form, sigma8 constant, ns constant)
deepszsim.load_vars module
read-in functions for input parameter yaml and h5 files
- deepszsim.load_vars.load_vars(file_name='config_simACTDR5.yaml', file_dir='/home/docs/checkouts/readthedocs.org/user_builds/deepszsim/checkouts/latest/deepszsim/Settings', survey_num: int = None, survey_name: str = None, survey_freq_val: int = None, cosmo_name: str = None, enforce_odd_pix: bool = True)
Load variables from YAML file.
Parameters:
- file_name: str
name of a yaml file that specifies the variables that you want to simulate
- file_dir: str
path to the yaml file named in file_name
- survey_num: None or int
optional – if there are multiple surveys, choose one by its position in the list
- survey_name: None or str
optional – if there are multiple surveys, choose one by its name
- survey_freq_val: None or int
optional – if a survey has multiple frequencies, choose a frequency by its value (in GHz)
- cosmo_name: None or str
optional – if there are multiple cosmologies, choose one by its name
- enforce_odd_pix: bool
whether or not to enforce an odd number of pixels in the simulated image (defaults to true so that the central pixel has the maximal y-value
Returns:
- dict: dict
a dictionary that will enable you to simulate SZ clusters with a given cosmology and set of observational properties.
- deepszsim.load_vars.readh5(fname, fdir=None)
Read data from h5 file.
Parameters:
- fname: str
name of the file you want to read (must include the suffix, eg ‘.h5’ or ‘.hd5’ or ‘.hdf5’)
- fdir: None or str
path to the directory that contains the file you want to read
Returns:
- fmdict: dict
dictionary that represents a single cluster (if fname represents a single cluster) or a nested set of clusters. Each cluster that is returned will itself have a dictionary with ‘maps’ and ‘params’ as keys, containing the output maps and the parameters that were used to make them, respectively
deepszsim.make_sz_cluster module
pressure profile, Compton-y, R200, C200 and temperature submap generating functions based on halo redshift and mass information
- deepszsim.make_sz_cluster.P200_Battaglia2012(M200_SM, redshift_z, load_vars_dict, R200_Mpc=None)
Calculates the P200 pressure of a cluster, as defined in Battaglia 2012.
Parameters:
- M200_SM: float
mass contained within R200, in units of solar masses
- redshift_z: float
redshift of the cluster (unitless)
- load_vars_dict: instance of load_vars.load_vars()
dictionary that includes background cosmology, ns, and sigma8 (necessary for the calculation of R200 and c200)
- R200_Mpc: None or float
if None, will calculate the radius that corresponds to the mass M200, the redshift redshift_z, and the cosmology contained in load_vars_dict
Returns:
- P200_kevcm3: instance
thermal pressure of the shell defined by R200 in units of keV/cm**3
- deepszsim.make_sz_cluster.Pe_to_y(profile, radii_mpc, M200_SM, redshift_z, load_vars_dict, alpha=1.0, gamma=-0.3, R200_Mpc=None)
Converts from an electron pressure profile to a compton-y profile; integrates over line of sight from -1 to 1 Mpc relative to center.
Parameters:
- profile: method
Method to get thermal pressure profile, accepts radius, M200, redshift_z, cosmo, and two additional parameters alpha and gamma that are usually fixed
- radii_mpc: array
array of radii corresponding to the profile in Mpc
- M200_SM: float
mass contained within R200, in units of solar masses
- redshift_z: float
redshift of the cluster (unitless)
- load_vars_dict: instance of load_vars.load_vars()
dictionary that includes background cosmology, ns, and sigma8 (necessary for the calculation of R200 and c200)
- alpha: float
variable fixed by Battaglia et al 2012 to 1.0
- gamma: float
variable fixed by Battaglia et al 2012 to -0.3
- R200_Mpc: None or float
if None, will calculate the radius that corresponds to the mass M200, the redshift redshift_z, and the cosmology contained in load_vars_dict
Returns:
- y_pro: array
Compton-y profile corresponding to the radii.
- deepszsim.make_sz_cluster.Pth_Battaglia2012(radius_mpc, M200_SM, redshift_z, load_vars_dict=None, alpha=1.0, gamma=-0.3, R200_Mpc=None)
Calculates the Pth profile using the Battaglia profile, Battaglia 2012, Equation 10. Pth is unitless. It is normalized by P200.
Parameters:
- radius_mpc: float or float array
radius for the pressure to be calculated at, in units of Mpc
- M200_SM: float
mass contained within R200, in units of solar masses
- redshift_z: float
redshift of the cluster (unitless)
- load_vars_dict: instance of load_vars.load_vars()
dictionary that includes background cosmology, ns, and sigma8 (necessary for the calculation of R200 and c200)
- alpha: float
fixed by Battaglia et al 2012 to 1.0
- gamma: float
fixed by Battaglia et al 2012 to -0.3
- R200_Mpc: None or float
if None, will calculate the radius that corresponds to the mass M200, the redshift redshift_z, and the cosmology contained in load_vars_dict
Returns:
- Pth: float or float array
thermal pressure profile normalized by P200, which has units keV/cm**3
- deepszsim.make_sz_cluster.generate_y_submap(M200_SM, redshift_z, profile='Battaglia2012', image_size_pixels=None, pixel_size_arcmin=None, load_vars_dict=None, alpha=1.0, gamma=-0.3, R200_Mpc=None)
Converts from an electron pressure profile to a compton-y profile, integrates over line of sight from -1 to 1 Mpc relative to center.
Parameters:
- M200_SM: float
the mass contained within R200 in solar masses
- redshift_z: float
the redshift of the cluster (unitless)
- profile: str
name of profile, currently only supports “Battaglia2012”
- image_size_pixels: float
num pixels to each side of center; end shape of submap will be (image_size_pixels, image_size_pixels)
- pixel_size_arcmin: float
size of each pixel in arcmin
- load_vars_dict: dict
result of running the load_vars() function, which includes a dictionary of cosmological and experimental parameters
- alpha: float
variable fixed by Battaglia et al 2012 to 1.0
- gamma: float
variable fixed by Battaglia et al 2012 to -0.3
- R200_Mpc: None or float
if None, will calculate the radius that corresponds to the mass M200, the redshift redshift_z, and the cosmology contained in load_vars_dict
Returns:
- y_map: float array
Compton-y submap with dimension (image_size_pixels, image_size_pixels)
- deepszsim.make_sz_cluster.get_r200_angsize_and_c200(M200_SM, redshift_z, load_vars_dict, angsize_density=None)
Get radius r200 and concentration c200
Parameters:
- M200_SM: float
mass contained within R200, in units of solar masses
- redshift_z: float
redshift of the cluster (unitless)
- load_vars_dict: dict
must contain ‘cosmo’ (a FlatLambaCDM instance describing the background cosmology), ‘sigma8’ (float, around 0.8), and ‘ns’ (float, around 0.96)
- angsize_density: None or str
density measure at which to calculate the angular size, if desired. If None, will not calculate an angular size. Otherwise, use a valid choice as specified in colossus.halo.mass_adv such as 500c
Returns:
- M200_SM: float
mass contained within R200, in units of solar masses
- R200_Mpc: float
radius of the cluster at 200 times the critical density of the universe in Mpc
- c200: float
concentration parameter
deepszsim.noise module
functions to generate simulated CMB map noise
- deepszsim.noise.generate_noise_map(image_size, noise_level, pix_size)
Generates a white noise map based on the noise level and beam size.
Parameters:
- image_size: int
Size of the noise map with dimensions N x N
- noise_level: float
Noise level of the survey
- pix_size: int
size of pixels in arcminutes
Returns:
- noise_map: float array
2D map of the white noise
deepszsim.read_yaml module
parsing a yaml file and returning a directory
deepszsim.simclusters module
- class deepszsim.simclusters.simulate_clusters(M200=None, redshift_z=None, num_halos=None, halo_params_dict=None, R200_Mpc=None, profile='Battaglia2012', image_size_pixels=None, image_size_arcmin=None, pixel_size_arcmin=None, alpha=1.0, gamma=-0.3, load_vars_yaml='/home/docs/checkouts/readthedocs.org/user_builds/deepszsim/checkouts/latest/deepszsim/Settings/config_simACTDR5.yaml', seed=None, tqverb=False)
Bases:
objectClass for simulating a distribution of clusters
Parameters:
- M200: float or array-like of float
mass contained within R200 in solar masses (same length as z_dist).
- redshift_z: float or array-like of float
redshift of the cluster (unitless) (same length as M200_dist).
- num_halos: None or int
number of halos to simulate if none supplied.
- halo_params_dict: None or dict
parameters from which to sample halos if num_halos specified, must contain zmin, zmax, m200min_SM, m200max_SM.
- R200_Mpc: None or float or np.ndarray(float)
if None, will calculate the R200 values corresponding to a given set of M200 and redshift_z values for the specified cosmology.
- profile: str
Name of profile, currently only supports “Battaglia2012”.
- image_size_pixels: None or int
image size in pixels (should be odd; if even, will return images whose sides are image_size_pixels+1 in length).
- image_size_arcmin: None or float
image size in arcmin.
- pixel_size_arcmin: None or float
pixel size in arcmin.
- alpha: float
fixed to equal 1.0 in Battaglia 2012.
- gamma: float
fixed to equal -0.3 in Battaglia 2012.
- load_vars_yaml: None or str
path to yaml file with params; if None, must explicitly include image specifications.
- seed: None or int
random seed value to sample with.
- tqverb: bool
whether or not to display tqdm progress bar while making T maps.
- get_T_maps(add_CMB=True, returnval=False)
Get ‘T’ maps from object
Parameters:
- add_CMB: bool
whether or not to include the CMB contribution to the final map.
- returnval: bool
whether or not to return the T maps themselves or simply update internal attribute.
Returns:
- self.clusters: float array
self._size many maps of the sky in units of uK, each of which is image_size_pixels x image_size_pixels in size.
- get_dT_maps()
get ‘dT’ maps from object
Parameters:
none
Returns:
- self.dT_maps: float array
self._size many maps of the dT values in units of uK, each of which is image_size_pixels x image_size_pixels in size.
- get_y_maps()
Get ‘y’ maps from object
Parameters:
none
Returns:
- self.y_maps: float array
self._size many maps of the Compton y value, each of which is image_size_pixels x image_size_pixels in size.
- ith_T_map(i, add_CMB=True)
Get ith ‘T’ map from object
Parameters:
- i: int
index of the returned map
- add_CMB: bool
whether or not to include the CMB contribution to the final map
Returns:
- self.clusters: array(float)
the ith map of the sky in units of uK, which is image_size_pixels x image_size_pixels in size
- save_map(i=None, nest_h5=True, nest_name=None, savedir='/home/docs/checkouts/readthedocs.org/user_builds/deepszsim/checkouts/latest/deepszsim/outfiles')
Save map to file
Parameters:
- i: None or int
map you want to save, if you only want to save a single map
- nest_h5: bool
Whether or not to nest the clusters into a single h5 file, assuming that you are saving all of the clusters that you have calculated
- nest_name: None or str
Name for the overall file if you are nesting them (otherwise, it will name it with the number of clusters plus the date plus a random string)
- savedir: str
Path to the directory you want to save into
Returns:
None
deepszsim.simtools module
beam convolution, SZ function and CMB simulation tools
- deepszsim.simtools.add_cmb_map_and_convolve(dT_map_uK, ps, pix_size_arcmin, beam_size_fwhp_arcmin)
Add CMB to the dT map and convolve with beam
Parameters:
- dT_map_uK: array
map to add to the CMB, units of -uK
- ps: array
power spectrum with shape (3, 3, lmax); clTT spectrum at ps[0][0]
- pix_size_arcmin: float
size of each pixel in arcmin
- beam_size_fwhp_arcmin: float or None
beam size in arcmin
Returns:
- dT submap: array
dT submap with same shape as dT_map, in units -uK
- deepszsim.simtools.convolve_map_with_gaussian_beam(pix_size_arcmin, beam_size_fwhm_arcmin, map_to_convolve)
Convolve the map with a Gaussian beam
Parameters:
- pix_size_arcmin: float
size of each pixel in arcmin
- beam_size_fwhm_arcmin: float
beam size in arcmin
- map_to_convolve: array
image to apply beam convolution to
Returns:
- convolved_map: array
map that has been convolved with beam
- deepszsim.simtools.f_sz(freq_ghz, T_CMB_K)
leading order correction to blackbody from Compton scattering see Eq 3.31 of https://background.uchicago.edu/~whu/thesis/chap3.pdf
Parameters:
- freq_ghz: float
Observation frequency f, in units of GHz
- T_CMB_K: instance of temperature spectrum
Temperature of CMB in K
Returns:
- fsz: float
radiation frequency
- deepszsim.simtools.get_cls(ns, cosmo, lmax=2000)
Makes a cmb temperature map based on the given power spectrum
Parameters:
- ns: float
scalar spectral index of the power-spectrum
- cosmo: FlatLambaCDM instance
background cosmology
Returns:
- ps: array
power spectrum as can be used in szluster.make_cmb_map
- deepszsim.simtools.make_cmb_map(shape, pix_size_arcmin, ps, seed=None)
Makes a cmb temperature map based on the given power spectrum
Parameters:
- shape: tuple
shape of submap in arcmin
- pix_size_arcmin: float
size of each pixel in arcmin
- ps: array
power spectrum with shape (3, 3, lmax); clTT spectrum at ps[0][0]
Returns:
- omap: array
CMB temperature map
deepszsim.utils module
unit conversions and data saving functions
- deepszsim.utils.Mpc_to_arcmin(r_arcmin, redshift_z, cosmo)
Changes the units of r from Mpc to arcmin
Parameters:
- r_arcmin: float or array
distance r, in units of arcmin
- redshift_z: float
redshift (unitless)
- cosmo: FlatLambaCDM instance
background cosmology for density calculation
Returns:
- r_arcmin / Mpc_per_arcmin: float or array (same type as r_arcmin)
distance r in units of Mpc; r_arcmin converted to Mpc
- deepszsim.utils.arcmin_to_Mpc(r_Mpc, redshift_z, cosmo)
Changes the units of r from arcmin to Mpc
Parameters:
- r_Mpc: float or array
distance r, in units of Mpc
- redshift_z: float
redshift (unitless)
- cosmo: FlatLambaCDM instance
background cosmology for density calculation
Returns:
- r_Mpc / Mpc_per_arcmin: float or array (same type as r_Mpc)
distance r in units of Mpc; r_Mpc converted to Mpc r_arcmin
- deepszsim.utils.gaussian_kernal(pix_size_arcmin, beam_size_fwhp_arcmin)
Create a Gaussian kernel for the beam
Parameters:
- pix_size_arcmin: float
size of each pixel in arcmin
- beam_size_fwhp_arcmin: float
beam size in arcmin
Returns:
- gaussian: array
2d gaussian kernal
- deepszsim.utils.save_sim_to_h5(file, name, data, attributes={}, overwrite=False)
Save data to h5 file
Parameters:
- file: h5py.File
The HDF5 file where the data will be saved
- name: str
The name under which the data will be stored in the HDF5 file
- data: dict
A dictionary where keys are dataset names (strings) and values are the corresponding dataset arrays or values
Returns:
none
deepszsim.visualization module
plotting functions
- deepszsim.visualization.plot_graphs(image, title=None, xlabel=None, ylabel=None, cbarlabel=None, width=None, specs=None, extend=False, logNorm=False)
Plotting tool function for our 2D submaps and CMB maps
Parameters:
- image: float array
graph that is plotted
- title: str
title of the graph
- xlabel: str
label of the x-axis
- ylabel: str
label of the y-axis
- cbarlabel: str
label of the color bar
- width: int
half rounded down of the width of output plot in pixels (eg, image size = 2*width+1)
- specs: dict
optional dictionary to pass title, xlabel, ylabel, cbarlabel, and width
- logNorm: bool
if true, uses a logarithmic normalization for the plot (using SymLogNorm in case values are negative)
Returns:
none
- deepszsim.visualization.plotting_specs(cluster)
Specifications for plot formatting
Parameters:
- cluster: array
dictionary representing a cluster instance, which contains a parameters dictionary which itself has ‘M200’, ‘redshift’, and ‘image_size_pixels’ keys
Returns:
- out: dict
dictionary that can be passed as the value of the specs kwarg in plot_graphs
Module contents
package for fast simulations of the Sunyaev-Zel’dovich effect from galaxy clusters