clusters module¶
Galaxy Cluster Ensemble Calculations.
Cluster mass-richness and mass-concentration scaling relations, and NFW halo profiles for weak lensing shear and magnification, including the effects of cluster miscentering offsets.
This framework calculates properties and profiles for every individual cluster, storing the data in tabular form. This is useful for fitting measured stacked weak lensing profiles, e.g. when you want to account for the full redshift, mass, and/or centroid offset distributions, and avoid fitting a single average mass at a single effective redshift.
-
class
clusterlensing.clusters.
ClusterEnsemble
(redshifts, cosmology=FlatLambdaCDM(name="Planck13", H0=67.8 km / (Mpc s), Om0=0.307, Tcmb0=2.725 K, Neff=3.05, m_nu=[ 0. 0. 0.06] eV, Ob0=0.0483), cm='DuttonMaccio')¶ Bases:
object
Ensemble of galaxy clusters and their properties.
The ClusterEnsemble object contains parameters and calculated values for every individual cluster in a sample. Initializing with a collection of redshifts z will fix the number of clusters described by the object. Setting the n200 values will then populate the object with the full set of available attributes (except for the NFW profiles). The values of z, n200, massrich_norm, and massrich_slope can be altered and these changes will propogate to the other (dependent) attributes.
In order to generate the attributes containing the NFW halo profiles, sigma_nfw and deltasigma_nfw, pass the desired radial bins (and additional optional parameters) to the calc_nfw method.
Parameters: z : array_like
Redshifts for each cluster in the sample. Should be 1D.
Other Parameters: cosmology : astropy.cosmology instance, optional
Cosmology to use in calculations, default astropy.cosmology.Planck13. Must be an instance of astropy.cosmology with ‘h’ and ‘Om0’ attributes.
cm : str, optional
Concentration-mass relation to use, default ‘DuttonMaccio’. Other choices are ‘Prada’ and ‘Duffy’.
Attributes
z
Cluster redshifts. n200
Cluster richness values. m200
Cluster masses. c200
Cluster concentration parameters. delta_c
Characteristic overdensities of the cluster halos. r200
Cluster Radii. rs
Cluster scale radii. Dang_l
Angular diameter distances to clusters. dataframe
Pandas DataFrame of cluster properties. massrich_norm
Normalization of Mass-Richness relation: massrich_slope
Slope of Mass-Richness relation: describe str Short description of the ClusterEnsemble object. number int Number of clusters in the sample. sigma_nfw Quantity Surface mass density of every cluster, (1D ndarray) with astropy.units of Msun/pc/pc. Generated for each element of rbins by running calc_nfw(rbins) method. deltasigma_nfw Quantity Differential surface mass density of every cluster, (1D ndarray) with astropy.units of Msun/pc/pc. Generated for each element of rbins by running calc_nfw(rbins) method. Methods
-
Dang_l
¶ Angular diameter distances to clusters.
Property: Returns distances in Mpc Type: Quantity (1D ndarray, with astropy.units of Mpc)
-
c200
¶ Cluster concentration parameters.
c200 is calculated from m200 and z using the mass-concentration relation specified when ClusterEnsemble object was created (default is relation from Dutton & Maccio 2015). Note that c200 = r200/rs.
Property: Returns c200 Type: ndarray
-
calc_nfw
(rbins, offsets=None, numTh=200, numRoff=200, numRinner=20, factorRouter=3)¶ Calculates Sigma and DeltaSigma profiles.
Generates the surface mass density (sigma_nfw attribute of parent object) and differential surface mass density (deltasigma_nfw attribute of parent object) profiles of each cluster, assuming a spherical NFW model. Optionally includes the effect of cluster miscentering offsets.
Parameters: rbins : array_like
Radial bins (in Mpc) for calculating cluster profiles. Should be 1D, optionally with astropy.units of Mpc.
offsets : array_like, optional
Parameter describing the width (in Mpc) of the Gaussian distribution of miscentering offsets. Should be 1D, optionally with astropy.units of Mpc.
Other Parameters: numTh : int, optional
Parameter to pass to SurfaceMassDensity(). Number of bins to use for integration over theta, for calculating offset profiles (no effect for offsets=None). Default 200.
numRoff : int, optional
Parameter to pass to SurfaceMassDensity(). Number of bins to use for integration over R_off, for calculating offset profiles (no effect for offsets=None). Default 200.
numRinner : int, optional
Parameter to pass to SurfaceMassDensity(). Number of bins at r < min(rbins) to use for integration over Sigma(<r), for calculating DeltaSigma (no effect for Sigma ever, and no effect for DeltaSigma if offsets=None). Default 20.
factorRouter : int, optional
Parameter to pass to SurfaceMassDensity(). Factor increase over number of rbins, at min(r) < r < max(r), of bins that will be used at for integration over Sigma(<r), for calculating DeltaSigma (no effect for Sigma, and no effect for DeltaSigma if offsets=None). Default 3.
-
dataframe
¶ Pandas DataFrame of cluster properties.
Property: Returns DataFrame Type: pandas.core.frame.DataFrame
-
delta_c
¶ Characteristic overdensities of the cluster halos.
Property: Returns characteristic overdensity Type: ndarray
-
m200
¶ Cluster masses.
Mass interior to a sphere of radius r200. If m200 is set directly, then richness n200 is calculated from m200 using the mass-richness scaling relation specified by the parameters massrich_norm and massrich_slope. If n200 is set directly, then m200 is calculated from n200 using the same scaling relation. Changes to m200 will propagate to all mass-dependant variables.
Property: Returns cluster masses in Msun Property type: Quantity 1D ndarray, with astropy.units of Msun. Setter: Sets cluster mass values in Msun Setter type: array_like Should be 1D array or list, optionally with units.
-
massrich_norm
¶ Normalization of Mass-Richness relation:
M200 = norm * (N200 / 20) ^ slope.
Changes to massrich_norm will propagate to all mass-dependant variables. (This will take current n200 values and convert them to m200; in order to retain original values of m200, save them in a temporary variable and reset them after this change).
Property: Returns normalization in Msun Property type: Quantity float, with astropy.units of Msun. Default is 2.7e+13 Msun. Setter: Sets normalization in Msun Setter type: float (optionally in astropy.units of Msun)
-
massrich_parameters
()¶ Print values of M200-N200 scaling relation parameters.
-
massrich_slope
¶ Slope of Mass-Richness relation:
M200 = norm * (N200 / 20) ^ slope.
Changes to massrich_slope will propagate to all mass-dependant variables. (This will take current n200 values and convert them to m200; in order to retain original values of m200, save them in a temporary variable and reset them after this change).
Property: Returns slope Property type: float Default value is 1.4. Setter: Sets slope Setter type: float
-
n200
¶ Cluster richness values.
If n200 is set directly, then mass m200 is calculated from n200 using the mass-richness scaling relation specified by the parameters massrich_norm and massrich_slope. If m200 is set directly, then n200 is calculated from m200 using the same scaling relation. Changes to n200 will propagate to all mass-dependant variables.
Property: Returns cluster richness values Property type: ndarray Setter: Sets cluster richness values Setter type: array_like
-
r200
¶ Cluster Radii.
r200 is the cluster radius within which the mean density is 200 times the critical energy density of the universe at that z.
Property: Returns r200 in Mpc Type: Quantity (1D ndarray, in astropy.units of Mpc)
-
rs
¶ Cluster scale radii.
Property: Returns scale radius in Mpc Type: Quantity (1D ndarray, in astropy.units of Mpc)
-
show
(notebook=True)¶ Display cluster properties and scaling relation parameters.
-
z
¶ Cluster redshifts.
Property: Returns cluster redshifts Property type: ndarray Setter: Sets cluster redshifts Setter type: array_like
-
-
clusterlensing.clusters.
calc_delta_c
(c200)¶ Calculate characteristic overdensity from concentration.
Parameters: c200 : ndarray or float
Cluster concentration parameter.
Returns: ndarray or float :
Cluster characteristic overdensity, of same type as c200.
-
clusterlensing.clusters.
mass_to_richness
(mass, norm=27000000000000.0, slope=1.4)¶ Calculate richness from mass.
Mass-richness relation assumed is: mass = norm * (richness / 20) ^ slope.
Parameters: mass : ndarray or float
Cluster mass value(s) in units of solar masses.
norm : float, optional
Normalization of mass-richness relation in units of solar masses, defaults to 2.7e13.
slope : float, optional
Slope of mass-richness relation in units of solar masses, defaults to 1.4.
Returns: ndarray or float :
Cluster richness(es), of same type as mass.
See also
richness_to_mass
- The inverse of this function.
-
clusterlensing.clusters.
richness_to_mass
(richness, norm=27000000000000.0, slope=1.4)¶ Calculate mass from richness.
Mass-richness relation assumed is: mass = norm * (richness / 20) ^ slope.
Parameters: richness : ndarray or float
Cluster richness value(s).
norm : float, optional
Normalization of mass-richness relation in units of solar masses, defaults to 2.7e13.
slope : float, optional
Slope of mass-richness relation in units of solar masses, defaults to 1.4.
Returns: ndarray or float :
Cluster mass(es) in units of solar masses, of same type as richness.
See also
mass_to_richness
- The inverse of this function.