Generate mock data for a cluster ================================ In this example we generate mock data with a variety of systematic effects including photometric redshifts, source galaxy distributions, and shape noise. We then populate a galaxy cluster object. This notebooks is organised as follows: - Imports and configuration setup - Generate mock data with different source galaxy options - Generate mock data with different field-of-view options - Generate mock data with different galaxy cluster options (only available with the Numcosmo and/or CCL backends). Use the ``os.environ['CLMM_MODELING_BACKEND']`` line below to select your backend. .. code:: ipython3 import os ## Uncomment the following line if you want to use a specific modeling backend among 'ct' (cluster-toolkit), 'ccl' (CCL) or 'nc' (Numcosmo). Default is 'ccl' # os.environ['CLMM_MODELING_BACKEND'] = 'nc' .. code:: ipython3 import numpy as np import matplotlib.pyplot as plt %matplotlib inline .. code:: ipython3 import clmm Make sure we know which version we’re using .. code:: ipython3 clmm.__version__ .. parsed-literal:: '1.12.0' Import mock data module and setup the configuration --------------------------------------------------- .. code:: ipython3 from clmm.support import mock_data as mock from clmm import Cosmology Mock data generation requires a defined cosmology .. code:: ipython3 mock_cosmo = Cosmology(H0=70.0, Omega_dm0=0.27 - 0.045, Omega_b0=0.045, Omega_k0=0.0) Mock data generation requires some cluster information. The default is to work with the NFW profile, using the “200,mean” mass definition. The Numcosmo and CCL backends allow for more flexibility (see last section of this notebook) .. code:: ipython3 cosmo = mock_cosmo cluster_id = "Awesome_cluster" cluster_m = 1.0e15 # M200,m cluster_z = 0.3 src_z = 0.8 concentration = 4 ngals = 1000 # number of source galaxies # Cluster centre coordinates cluster_ra = 50.0 cluster_dec = 87.0 Generate the mock catalog with different source galaxy options -------------------------------------------------------------- - Clean data: no noise, all galaxies at the same redshift .. code:: ipython3 zsrc_min = cluster_z + 0.1 .. code:: ipython3 ideal_data = mock.generate_galaxy_catalog( cluster_m, cluster_z, concentration, cosmo, src_z, ngals=ngals, cluster_ra=cluster_ra, cluster_dec=cluster_dec, ) .. code:: ipython3 # let's put all these quantities in a single dictionary to facilitate clarity cluster_kwargs = { "cluster_m": cluster_m, "cluster_z": cluster_z, "cluster_ra": cluster_ra, "cluster_dec": cluster_dec, "cluster_c": concentration, "cosmo": cosmo, } .. code:: ipython3 ideal_data = mock.generate_galaxy_catalog(**cluster_kwargs, zsrc=src_z, ngals=ngals) - Noisy data: shape noise, all galaxies at the same redshift .. code:: ipython3 noisy_data_src_z = mock.generate_galaxy_catalog( **cluster_kwargs, zsrc=src_z, shapenoise=0.05, ngals=ngals ) - Noisy data: shape noise plus measurement error, all galaxies at the same redshift .. code:: ipython3 noisy_data_src_z_e_err = mock.generate_galaxy_catalog( **cluster_kwargs, zsrc=src_z, shapenoise=0.05, mean_e_err=0.05, ngals=ngals ) .. container:: alert alert-warning **WARNING:** Experimental feature. Uncertainties are created by simply drawing random numbers near the value specified by ``mean_e_err``. Use at your own risk. This will be improved in future releases. - Noisy data: photo-z errors (and pdfs!), all galaxies at the same redshift. At present, the pdfs are generated by drawing a random value from a Normal distribution with mean ``ztrue`` and width given by ``dz = (1+z)*photoz_sigma_unscaled``, and the pdf is this Normal distribution centered around ``z`` instead of ``ztrue``. .. code:: ipython3 noisy_data_photoz = mock.generate_galaxy_catalog( **cluster_kwargs, zsrc=src_z, shapenoise=0.05, photoz_sigma_unscaled=0.05, pz_bins=101, ngals=ngals ) In the default PDF storing (``pzpdf_type='shared_bins'``), the values of the PDF are added to the ``pzpdf`` column and the binning scheme is stored in the ``pzpdf_info`` attribute: .. code:: ipython3 noisy_data_photoz.pzpdf_info .. parsed-literal:: {'type': 'shared_bins', 'zbins': array([0. , 0.01959865, 0.0391973 , 0.05879595, 0.0783946 , 0.09799325, 0.1175919 , 0.13719055, 0.15678921, 0.17638786, 0.19598651, 0.21558516, 0.23518381, 0.25478246, 0.27438111, 0.29397976, 0.31357841, 0.33317706, 0.35277571, 0.37237436, 0.39197301, 0.41157166, 0.43117031, 0.45076896, 0.47036762, 0.48996627, 0.50956492, 0.52916357, 0.54876222, 0.56836087, 0.58795952, 0.60755817, 0.62715682, 0.64675547, 0.66635412, 0.68595277, 0.70555142, 0.72515007, 0.74474872, 0.76434737, 0.78394603, 0.80354468, 0.82314333, 0.84274198, 0.86234063, 0.88193928, 0.90153793, 0.92113658, 0.94073523, 0.96033388, 0.97993253, 0.99953118, 1.01912983, 1.03872848, 1.05832713, 1.07792578, 1.09752444, 1.11712309, 1.13672174, 1.15632039, 1.17591904, 1.19551769, 1.21511634, 1.23471499, 1.25431364, 1.27391229, 1.29351094, 1.31310959, 1.33270824, 1.35230689, 1.37190554, 1.3915042 , 1.41110285, 1.4307015 , 1.45030015, 1.4698988 , 1.48949745, 1.5090961 , 1.52869475, 1.5482934 , 1.56789205, 1.5874907 , 1.60708935, 1.626688 , 1.64628665, 1.6658853 , 1.68548395, 1.70508261, 1.72468126, 1.74427991, 1.76387856, 1.78347721, 1.80307586, 1.82267451, 1.84227316, 1.86187181, 1.88147046, 1.90106911, 1.92066776, 1.94026641, 1.95986506])} .. code:: ipython3 noisy_data_photoz[:5] .. raw:: html GCData
defined by: cosmo=None
with columns: ra, dec, e1, e2, z, ztrue, pzpdf, id
and pzpdf: shared_bins [0. 0.02 0.04 0.06 0.08 ... 1.88 1.9 1.92 1.94 1.96]
5 objects
radece1e2zztruepzpdfid
float64float64float64float64float64float64float64[101]int64
46.2572880823355286.89223054243050.08236207803522719-0.054118540767148880.69871031819970640.83.622089304592751e-13 .. 1.0180322952717169e-420
50.0777871562318987.10740150490655-0.014752252971336390.038828008441195870.73929621321758590.89.870794333079776e-15 .. 5.1047370592912794e-401
48.4080987988511986.87728907171557-1.2575884722892804e-050.020962676367362040.90093032742937560.87.709520397996465e-22 .. 3.849218150728941e-302
53.01436891183186686.9440239760556-0.011911350886457556-0.058930111442774110.84479698293503340.83.2663820149084725e-19 .. 2.0599039175494822e-333
53.43523018790920587.09017457767816-0.019995431173391790.023543677736634350.61913835656937030.82.345198551898936e-10 .. 2.86718907682299e-484
It is also possible to generate individual binning of the PDF with ``pzpdf_type='individual_bins'``. In this case, both the bin values and the PDF are added to the ``pzpdf`` column: .. code:: ipython3 mock.generate_galaxy_catalog( **cluster_kwargs, zsrc=src_z, shapenoise=0.05, photoz_sigma_unscaled=0.05, pz_bins=101, ngals=ngals, pzpdf_type="individual_bins" )[:5] .. raw:: html GCData
defined by: cosmo=None
with columns: ra, dec, e1, e2, z, ztrue, pzbins, pzpdf, id
and pzpdf: individual_bins
5 objects
radece1e2zztruepzbinspzpdfid
float64float64float64float64float64float64float64[101]float64[101]int64
51.7180958895545586.794184998981240.045171775281032730.0079544195667603690.71151929310340460.80.0 .. 1.61151929310340461.1876947502964896e-13 .. 8.54955402967386e-220
46.70800163801112686.88819976612274-0.012879815467979831-0.081754480509692020.69456815514609630.80.0 .. 1.59456815514609645.172194687003747e-13 .. 8.549554029673739e-221
49.3945199868530587.20677909110817-0.00396285839757514950.01749945061672090.83158050787762080.80.0 .. 1.7315805078776211.2824084572051592e-18 .. 8.549554029673739e-222
46.86063450494019486.93654387908563-0.0462597921376629650.0053228967998943670.69178144860663620.80.0 .. 1.59178144860663646.565150731239997e-13 .. 8.549554029673677e-223
51.26562279760640586.86646628375920.0316686284294420.0507252738466855150.75637302584879640.80.0 .. 1.65637302584879652.0400190276050062e-15 .. 8.549554029673739e-224
Note that if ``pzpdf_type=None``, the pdf is not stored: .. code:: ipython3 mock.generate_galaxy_catalog( **cluster_kwargs, zsrc=src_z, shapenoise=0.05, photoz_sigma_unscaled=0.05, pz_bins=101, ngals=ngals, pzpdf_type=None )[:5] .. raw:: html GCData
defined by: cosmo=None
with columns: ra, dec, e1, e2, z, ztrue, id
5 objects
radece1e2zztrueid
float64float64float64float64float64float64int64
49.48124527921151486.88497246906974-0.046399017054751070.0168106436341245240.65627993127820620.80
46.7658044620298987.117408915092980.008967332934760024-0.047367075980020530.75660219592537010.81
47.6768345136793887.09593810890625-1.6856807502580002e-06-0.0387854437794893440.74275551099609040.82
45.2589448651097687.08199536314510.010890284110664075-0.0021940007340567420.92925007307238240.83
53.4700617097647286.89370333944125-0.120533900351914110.069700428605588020.76726137341205340.84
- Clean data: source galaxy redshifts drawn from a redshift distribution instead of fixed ``src_z`` value. Options are ``chang13`` for Chang et al. 2013 or ``desc_srd`` for the distribution given in the DESC Science Requirement Document. No shape noise or photoz errors. .. code:: ipython3 ideal_with_src_dist = mock.generate_galaxy_catalog( **cluster_kwargs, zsrc="chang13", zsrc_min=zsrc_min, zsrc_max=7.0, ngals=ngals ) - Noisy data: galaxies following redshift distribution, redshift error, shape noise .. code:: ipython3 pzbins = np.linspace(0, 10, 1001) .. code:: ipython3 allsystematics = mock.generate_galaxy_catalog( **cluster_kwargs, zsrc="chang13", zsrc_min=zsrc_min, photoz_sigma_unscaled=0.05, ngals=ngals, pz_bins=pzbins, ) .. code:: ipython3 allsystematics2 = mock.generate_galaxy_catalog( **cluster_kwargs, zsrc="desc_srd", zsrc_min=zsrc_min, zsrc_max=7.0, photoz_sigma_unscaled=0.05, ngals=ngals, pz_bins=pzbins, shapenoise=0.05, ) Sanity check: checking that no galaxies were originally drawn below zsrc_min, before photoz errors are applied (when relevant) .. code:: ipython3 print( f"""Number of galaxies below zsrc_min: ideal_data : {np.sum(ideal_data['ztrue']GCData length=3
idxradece1e2zztrueid
05.278e+018.690e+01-4.495e-03-1.083e-028.000e-018.000e-010
15.095e+018.702e+01-2.997e-022.582e-028.000e-018.000e-011
24.637e+018.703e+01-1.019e-02-3.995e-038.000e-018.000e-012
- With photo-z errors .. code:: ipython3 for n in noisy_data_photoz.colnames: if n != "id": noisy_data_photoz[n].format = "%6.3e" noisy_data_photoz[0:3].show_in_notebook() .. raw:: html GCData length=3
idxradece1e2zztruepzpdfid
04.626e+018.689e+018.236e-02-5.412e-026.987e-018.000e-013.622e-13 .. 1.018e-420
15.008e+018.711e+01-1.475e-023.883e-027.393e-018.000e-019.871e-15 .. 5.105e-401
24.841e+018.688e+01-1.258e-052.096e-029.009e-018.000e-017.710e-22 .. 3.849e-302
- Histogram of the redshift distribution of background galaxies, for the true (originally drawn) redshift and the redshift once photoz errors have been added, and the stacked pdf. By construction no true redshift occurs below zsrc_min, but some ‘observed’ redshifts (i.e. including photoz errors) might be. .. code:: ipython3 plt.hist( allsystematics["z"], bins=50, alpha=0.3, density=True, label="measured z (i.e. including photoz error)", ) plt.hist(allsystematics["ztrue"], bins=50, alpha=0.3, density=True, label="true z") stacked_pdf = np.mean(allsystematics["pzpdf"], axis=0) plt.plot(allsystematics.pzpdf_info["zbins"], stacked_pdf, "C3", label="stacked pdf") plt.axvline(zsrc_min, color="k", label="requested zmin") plt.xlabel("Source Redshift") plt.ylabel("n(z)") plt.legend() plt.xlim(0, 5) pass .. image:: demo_mock_cluster_files/demo_mock_cluster_47_0.png .. code:: ipython3 plt.hist(allsystematics["ztrue"], bins=50, alpha=0.3, label="true z") plt.hist(allsystematics2["ztrue"], bins=50, alpha=0.3, label="true z"); .. image:: demo_mock_cluster_files/demo_mock_cluster_48_0.png Populate a galaxy cluster object -------------------------------- .. code:: ipython3 gc_object = clmm.GalaxyCluster(cluster_id, cluster_ra, cluster_dec, cluster_z, allsystematics) From a ``GalaxyCluster`` object that has photoz information, ``draw_gal_z_from_pdz`` allows to generate ``nobj`` random redshifts of each galaxy in ``galcat``, from its photoz pdf, and store the result in a new ``zcol_out`` column. .. code:: ipython3 z_random = gc_object.draw_gal_z_from_pdz(zcol_out="z_random", overwrite=False, nobj=1) The plot below shows the “observed photoz pdf” (blue), centered on the “observed z” (red), the true redshift from which the shear where computed (green) and a random redshift (orange) computed from the pdf .. code:: ipython3 # p(z) for one of the galaxies in the catalog, galid = 0 plt.fill(gc_object.galcat.pzpdf_info["zbins"], allsystematics["pzpdf"][galid], alpha=0.3) plt.plot(gc_object.galcat.pzpdf_info["zbins"], gc_object.galcat["pzpdf"][galid], label="Photoz pdf") plt.axvline(gc_object.galcat["z"][galid], label="Observed z", color="red") plt.axvline(gc_object.galcat["ztrue"][galid], label="True z", color="g") plt.axvline(gc_object.galcat["z_random"][galid], label="Random z from pdf", color="orange") plt.xlabel("Redshift") plt.ylabel("Photo-z Probability Distribution") plt.legend(loc=1) plt.xlim(gc_object.galcat["z"][galid] - 0.5, gc_object.galcat["z"][galid] + 0.5) .. parsed-literal:: (1.328538487022893, 2.328538487022893) .. image:: demo_mock_cluster_files/demo_mock_cluster_54_1.png Plot source galaxy ellipticities .. code:: ipython3 plt.scatter(gc_object.galcat["e1"], gc_object.galcat["e2"]) plt.xlim(-0.2, 0.2) plt.ylim(-0.2, 0.2) plt.xlabel("Ellipticity 1", fontsize="x-large") plt.ylabel("Ellipticity 2", fontsize="x-large") .. parsed-literal:: Text(0, 0.5, 'Ellipticity 2') .. image:: demo_mock_cluster_files/demo_mock_cluster_56_1.png Generate the mock data catalog with different field-of-view options ------------------------------------------------------------------- In the examples above, ``ngals=1000`` galaxies were simulated in a field corresponding to a 8 Mpc x 8 Mpc (proper distance) square box at the cluster redshift (this is the default). The user may however vary the field size and/or provide a galaxy density (instead of a number of galaxies). This is examplified below, using the ``allsystematics`` example. - ``ngals = 1000`` in a 4 x 4 Mpc box. Asking for the same number of galaxies in a smaller field of view yields high galaxy density .. code:: ipython3 allsystematics2 = mock.generate_galaxy_catalog( **cluster_kwargs, zsrc="chang13", zsrc_min=zsrc_min, zsrc_max=7.0, shapenoise=0.05, photoz_sigma_unscaled=0.05, field_size=4, ngals=ngals ) .. code:: ipython3 plt.scatter(allsystematics["ra"], allsystematics["dec"], marker=".", label="default 8 x 8 Mpc FoV") plt.scatter(allsystematics2["ra"], allsystematics2["dec"], marker=".", label="user-defined FoV") plt.legend() .. parsed-literal:: .. image:: demo_mock_cluster_files/demo_mock_cluster_61_1.png - Alternatively, the user may provide a galaxy density (here ~1 gal/arcmin2 to roughly match 1000 galaxies, given the configuration) and the number of galaxies to draw will automatically be adjusted to the box size. .. code:: ipython3 allsystematics3 = mock.generate_galaxy_catalog( **cluster_kwargs, zsrc="chang13", zsrc_min=zsrc_min, zsrc_max=7.0, shapenoise=0.05, photoz_sigma_unscaled=0.05, ngal_density=1.3, ) print(f"Number of drawn galaxies = {len(allsystematics3)}") .. parsed-literal:: Number of drawn galaxies = 989 .. code:: ipython3 allsystematics4 = mock.generate_galaxy_catalog( **cluster_kwargs, zsrc="desc_srd", zsrc_min=zsrc_min, zsrc_max=7.0, shapenoise=0.05, photoz_sigma_unscaled=0.05, ngal_density=1.3, ) print(f"Number of drawn galaxies = {len(allsystematics4)}") .. parsed-literal:: Number of drawn galaxies = 1020 .. code:: ipython3 plt.scatter(allsystematics["ra"], allsystematics["dec"], marker=".", label="ngals = 1000") plt.scatter( allsystematics3["ra"], allsystematics3["dec"], marker=".", label="ngal_density = 1 gal / arcmin2", ) plt.legend() .. parsed-literal:: .. image:: demo_mock_cluster_files/demo_mock_cluster_65_1.png Generate mock data with different galaxy cluster options -------------------------------------------------------- WARNING: Available options depend on the modeling backend: - Cluster-toolkit allows for other values of the overdensity parameter, but is retricted to working with the mean mass definition - Both CCL and Numcosmo allow for different values of the overdensity parameter, but work with both the mean and critical mass definition - Numcosmo further allows for the Einasto or Burkert density profiles to be used instead of the NFW profile Changing the overdensity parameter (all backend) - ``delta_so`` keyword (default = 200) ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ .. code:: ipython3 allsystematics_500mean = mock.generate_galaxy_catalog( **cluster_kwargs, zsrc="chang13", delta_so=500, zsrc_min=zsrc_min, zsrc_max=7.0, shapenoise=0.05, photoz_sigma_unscaled=0.05, ngals=ngals ) Using the critical mass definition (Numcosmo and CCL only) - ``massdef`` keyword (default = ‘mean’) ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ WARNING: error will be raised if using the cluster-toolkit backend .. code:: ipython3 allsystematics_200critical = mock.generate_galaxy_catalog( **cluster_kwargs, zsrc="chang13", massdef="critical", zsrc_min=zsrc_min, zsrc_max=7.0, shapenoise=0.05, photoz_sigma_unscaled=0.05, ngals=ngals ) Changing the halo density profile (Numcosmo and CCL only) - ``halo_profile_model`` keyword (default = ‘nfw’) ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ WARNING: error will be raised if using the cluster-toolkit or CCL backends .. code:: ipython3 allsystematics_200m_einasto = mock.generate_galaxy_catalog( **cluster_kwargs, zsrc="chang13", halo_profile_model="einasto", zsrc_min=zsrc_min, zsrc_max=7.0, shapenoise=0.05, photoz_sigma_unscaled=0.05, ngals=ngals )