{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Calculate cosmological distances with CCL\n", "In this example, we will calculate various cosmological distances for an example cosmology." ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "import numpy as np\n", "import pylab as plt\n", "import pyccl as ccl\n", "%matplotlib inline" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Set up a Cosmology object\n", "`Cosmology` objects contain the parameters and metadata needed as inputs to most functions. Each `Cosmology` object has a set of cosmological parameters attached to it. In this example, we will only use the parameters of a vanilla LCDM model, but simple extensions (like curvature, neutrino mass, and w0/wa) are also supported.\n", "\n", "`Cosmology` objects also contain precomputed data (e.g. splines) to help speed-up certain calculations. As such, `Cosmology` objects are supposed to be immutable; you should create a new `Cosmology` object when you want to change the values of any cosmological parameters." ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "cosmo = ccl.Cosmology(Omega_c=0.27, Omega_b=0.045, h=0.67, A_s=2.1e-9, n_s=0.96)" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "\n", "\tA_s = 2.1e-09\n", "\tNeff = 3.044\n", "\tOmega_b = 0.045\n", "\tOmega_c = 0.27\n", "\th = 0.67\n", "\tn_s = 0.96\n", "\textra_parameters =\n", "\t\temu = {'neutrinos': 'strict'}\n", "\tHASH_ACCURACY_PARAMS = 0xa0af9c5b75d6de86\n" ] } ], "source": [ "print(cosmo)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Parameter values can be accessed from the `Cosmology` object contains, like so:" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "0.27\n" ] } ], "source": [ "print(cosmo['Omega_c'])" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "nan\n" ] } ], "source": [ "print(cosmo['sigma8'])" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Notice $\\sigma_8$ is nan here because it has not yet been computed. Instead, the normalization of the power spectrum was defined via $A_s$ above." ] }, { "cell_type": "markdown", "metadata": { "collapsed": true }, "source": [ "### Cosmological Distances\n", "\n", "With a cosmology in hand, we can begin performing some calculations. We can start with the most basic measure, the comoving radial distance. " ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "1962.9390685778556" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "z = 0.5 \n", "ccl.comoving_radial_distance(cosmo, 1/(1+z)) # Mpc" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note that all distance function calls require scale factors, not redshifts. This function can take a `numpy` array of values as well." ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "array([ 0. , 436.6985067 , 851.3944111 , 1243.78925379,\n", " 1614.09910482, 1962.93906858, 2291.20721282, 2599.98198852,\n", " 2890.43959663, 3163.79072723])" ] }, "execution_count": 7, "metadata": {}, "output_type": "execute_result" } ], "source": [ "zs = np.arange(0, 1, 0.1)\n", "ccl.comoving_radial_distance(cosmo, 1/(1+zs))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "CCL also supports calculation of the comoving angular distance. In flat spacetime (like the cosmology we have here) it is the same as the radial distance. " ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "1962.9390685778556" ] }, "execution_count": 8, "metadata": {}, "output_type": "execute_result" } ], "source": [ "ccl.comoving_angular_distance(cosmo, 1/(1+z))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "If we create a cosmology with curvature, we'll get a different result. " ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Radial Dist. = 1992.53 Mpc \t Angular Dist. = 1999.12 Mpc\n" ] } ], "source": [ "curved_cosmo = ccl.Cosmology(Omega_k = 0.1, Omega_c=0.17, Omega_b=0.045, h=0.67, A_s=2.1e-9, n_s=0.96)\n", "\n", "chi_rad = ccl.comoving_radial_distance(curved_cosmo, 1/(1+z))\n", "chi_curved = ccl.comoving_angular_distance(curved_cosmo, 1/(1+z))\n", "\n", "print ('Radial Dist. = %.2f Mpc \\t Angular Dist. = %.2f Mpc'%(chi_rad, chi_curved))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "CCL explictly supports the calculation of the luminosity distance and the distance modulus too:" ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Luminosity Dist = 2944.41 Mpc \t Distance Modulus = 42.34 \n" ] } ], "source": [ "chi_lum = ccl.luminosity_distance(cosmo, 1/(1+z))\n", "DM = ccl.distance_modulus(cosmo, 1/(1+z))\n", "print('Luminosity Dist = %.2f Mpc \\t Distance Modulus = %.2f ' % (chi_lum, DM))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Finally, CCL supports an inverse operation, which calculates the scale factor for a given comoving distance:" ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "0.6666639215879805" ] }, "execution_count": 11, "metadata": {}, "output_type": "execute_result" } ], "source": [ "ccl.scale_factor_of_chi(cosmo, 1962.96)" ] } ], "metadata": { "kernelspec": { "display_name": "Python [conda env:ccl-dev-v3] *", "language": "python", "name": "conda-env-ccl-dev-v3-py" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.9.16" } }, "nbformat": 4, "nbformat_minor": 2 }