Calculation results

class camb.results.CAMBdata(*args, **kwargs)[source]

An object for storing calculational data, parameters and transfer functions. Results for a set of parameters (given in a CAMBparams instance) are returned by the camb.get_background(), camb.get_transfer_functions() or camb.get_results() functions. Exactly which quantities are already calculated depends on which of these functions you use and the input parameters.

To quickly make a fully calculated CAMBdata instance for a set of parameters you can call camb.get_results().

Variables:
  • Paramscamb.model.CAMBparams

  • ThermoDerivedParams – (float64 array) array of derived parameters, see get_derived_params() to get as a dictionary

  • flat – (boolean) flat universe

  • closed – (boolean) closed universe

  • grhocrit – (float64) kappa*a^2*rho_c(0)/c^2 with units of Mpc**(-2)

  • grhog – (float64) kappa/c^2*4*sigma_B/c^3 T_CMB^4

  • grhor – (float64) 7/8*(4/11)^(4/3)*grhog (per massless neutrino species)

  • grhob – (float64) baryon contribution

  • grhoc – (float64) CDM contribution

  • grhov – (float64) Dark energy contribution

  • grhornomass – (float64) grhor*number of massless neutrino species

  • grhok – (float64) curvature contribution to critical density

  • taurst – (float64) time at start of recombination

  • dtaurec – (float64) time step in recombination

  • taurend – (float64) time at end of recombination

  • tau_maxvis – (float64) time at peak visibility

  • adotrad – (float64) da/d tau in early radiation-dominated era

  • omega_de – (float64) Omega for dark energy today

  • curv – (float64) curvature K

  • curvature_radius – (float64) \(1/\sqrt{|K|}\)

  • Ksign – (float64) Ksign = 1,0 or -1

  • tau0 – (float64) conformal time today

  • chi0 – (float64) comoving angular diameter distance of big bang; rofChi(tau0/curvature_radius)

  • scale – (float64) relative to flat. e.g. for scaling L sampling

  • akthom – (float64) sigma_T * (number density of protons now)

  • fHe – (float64) n_He_tot / n_H_tot

  • Nnow – (float64) number density today

  • z_eq – (float64) matter-radiation equality redshift assuming all neutrinos relativistic

  • grhormass – (float64 array)

  • nu_masses – (float64 array)

  • num_transfer_redshifts – (integer) Number of calculated redshift outputs for the matter transfer (including those for CMB lensing)

  • transfer_redshifts – (float64 array) Calculated output redshifts

  • PK_redshifts_index – (integer array) Indices of the requested PK_redshifts

  • OnlyTransfers – (boolean) Only calculating transfer functions, not power spectra

  • HasScalarTimeSources – (boolean) calculate and save time source functions, not power spectra

angular_diameter_distance(z)[source]

Get (non-comoving) angular diameter distance to redshift z.

Must have called calc_background(), calc_background_no_thermo() or calculated transfer functions or power spectra.

Parameters:

z – redshift or array of redshifts

Returns:

angular diameter distances, matching rank of z

angular_diameter_distance2(z1, z2)[source]

Get angular diameter distance between two redshifts \(\frac{r}{1+z_2}\text{sin}_K\left(\frac{\chi(z_2) - \chi(z_1)}{r}\right)\) where \(r\) is curvature radius and \(\chi\) is the comoving radial distance. If z_1 >= z_2 returns zero.

Must have called calc_background(), calc_background_no_thermo() or calculated transfer functions or power spectra.

Parameters:
  • z1 – redshift 1, or orray of redshifts

  • z2 – redshift 2, or orray of redshifts

Returns:

result (scalar or array of distances between pairs of z1, z2)

calc_background(params)[source]

Calculate the background evolution and thermal history. e.g. call this if you want to get derived parameters and call background functions

Parameters:

paramsCAMBparams instance to use

calc_background_no_thermo(params, do_reion=False)[source]

Calculate the background evolution without calculating thermal or ionization history. e.g. call this if you want to just use angular_diameter_distance() and similar background functions

Parameters:
  • paramsCAMBparams instance to use

  • do_reion – whether to initialize the reionization model

calc_power_spectra(params=None)[source]

Calculates transfer functions and power spectra.

Parameters:

params – optional CAMBparams instance with parameters to use

calc_transfers(params, only_transfers=True, only_time_sources=False)[source]

Calculate the transfer functions (for CMB and matter power, as determined by params.WantCls, params.WantTransfer).

Parameters:
  • paramsCAMBparams instance with parameters to use

  • only_transfers – only calculate transfer functions, no power spectra

  • only_time_sources – only calculate time transfer functions, no (p,l,k) transfer functions or non-linear scaling

Returns:

non-zero if error, zero if OK

comoving_radial_distance(z, tol=0.0001)[source]

Get comoving radial distance from us to redshift z in Mpc. This is efficient for arrays.

Must have called calc_background(), calc_background_no_thermo() or calculated transfer functions or power spectra.

Parameters:
  • z – redshift

  • tol – numerical tolerance parameter

Returns:

comoving radial distance (Mpc)

conformal_time(z, presorted=None, tol=None)[source]

Conformal time from hot big bang to redshift z in Megaparsec.

Parameters:
  • z – redshift or array of redshifts

  • presorted – if True, redshifts already sorted to be monotonically increasing, if False decreasing, or if None unsorted. If presorted is True or False no checks are done.

  • tol – integration tolerance

Returns:

eta(z)/Mpc

conformal_time_a1_a2(a1, a2)[source]

Get conformal time between two scale factors (=comoving radial distance travelled by light on light cone)

Parameters:
  • a1 – scale factor 1

  • a2 – scale factor 2

Returns:

eta(a2)-eta(a1) = chi(a1)-chi(a2) in Megaparsec

copy()

Make an independent copy of this object.

Returns:

a deep copy of self

cosmomc_theta()[source]

Get \(\theta_{\rm MC}\), an approximation of the ratio of the sound horizon to the angular diameter distance at recombination.

Returns:

\(\theta_{\rm MC}\)

classmethod dict(state)

Make an instance of the class from a dictionary of field values (used to restore from repr)

Parameters:

state – dictionary of values

Returns:

new instance

get_BAO(redshifts, params)[source]

Get BAO parameters at given redshifts, using parameters in params

Parameters:
  • redshifts – list of redshifts

  • params – optional CAMBparams instance to use

Returns:

array of rs/DV, H, DA, F_AP for each redshift as 2D array

get_Omega(var, z=0)[source]

Get density relative to critical density of variables var

Parameters:
  • var – one of ‘K’, ‘cdm’, ‘baryon’, ‘photon’, ‘neutrino’ (massless), ‘nu’ (massive neutrinos), ‘de’

  • z – redshift

Returns:

\(\Omega_i(a)\)

get_background_densities(a, vars=['tot', 'K', 'cdm', 'baryon', 'photon', 'neutrino', 'nu', 'de'], format='dict')[source]

Get the individual densities as a function of scale factor. Returns \(8\pi G a^4 \rho_i\) in Mpc units. \(\Omega_i\) can be simply obtained by taking the ratio of the components to tot.

Parameters:
  • a – scale factor or array of scale factors

  • vars – list of variables to output (default all)

  • format – ‘dict’ or ‘array’, for either dict of 1D arrays indexed by name, or 2D array

Returns:

n_a x len(vars) 2D numpy array or dict of 1D arrays of \(8\pi G a^4 \rho_i\) in Mpc units.

get_background_outputs()[source]

Get BAO values for redshifts set in Params.z_outputs

Returns:

rs/DV, H, DA, F_AP for each requested redshift (as 2D array)

get_background_redshift_evolution(z, vars=['x_e', 'opacity', 'visibility', 'cs2b', 'T_b', 'dopacity', 'ddopacity', 'dvisibility', 'ddvisibility'], format='dict')[source]

Get the evolution of background variables a function of redshift. For the moment a and H are rather perversely only available via get_time_evolution()

Parameters:
  • z – array of requested redshifts to output

  • vars – list of variable names to output

  • format – ‘dict’ or ‘array’, for either dict of 1D arrays indexed by name, or 2D array

Returns:

n_eta x len(vars) 2D numpy array of outputs or dict of 1D arrays

get_background_time_evolution(eta, vars=['x_e', 'opacity', 'visibility', 'cs2b', 'T_b', 'dopacity', 'ddopacity', 'dvisibility', 'ddvisibility'], format='dict')[source]

Get the evolution of background variables a function of conformal time. For the moment a and H are rather perversely only available via get_time_evolution()

Parameters:
  • eta – array of requested conformal times to output

  • vars – list of variable names to output

  • format – ‘dict’ or ‘array’, for either dict of 1D arrays indexed by name, or 2D array

Returns:

n_eta x len(vars) 2D numpy array of outputs or dict of 1D arrays

get_cmb_correlation_functions(params=None, lmax=None, spectrum='lensed_scalar', xvals=None, sampling_factor=1)[source]

Get the CMB correlation functions from the power spectra. By default evaluated at points \(\cos(\theta)\) = xvals that are roots of Legendre polynomials, for accurate back integration with correlations.corr2cl(). If xvals is explicitly given, instead calculates correlations at provided \(\cos(\theta)\) values.

Parameters:
  • params – optional CAMBparams instance with parameters to use. If None, must have previously set parameters and called calc_power_spectra() (e.g. if you got this instance using camb.get_results()),

  • lmax – optional maximum L to use from the cls arrays

  • spectrum – type of CMB power spectrum to get; default ‘lensed_scalar’, one of [‘total’, ‘unlensed_scalar’, ‘unlensed_total’, ‘lensed_scalar’, ‘tensor’]

  • xvals – optional array of \(\cos(\theta)\) values at which to calculate correlation function.

  • sampling_factor – multiple of lmax for the Gauss-Legendre order if xvals not given (default 1)

Returns:

if xvals not given: corrs, xvals, weights; if xvals specified, just corrs. corrs is 2D array corrs[i, ix], where ix=0,1,2,3 are T, Q+U, Q-U and cross, and i indexes xvals

get_cmb_power_spectra(params=None, lmax=None, spectra=('total', 'unlensed_scalar', 'unlensed_total', 'lensed_scalar', 'tensor', 'lens_potential'), CMB_unit=None, raw_cl=False)[source]

Get CMB power spectra, as requested by the ‘spectra’ argument. All power spectra are \(\ell(\ell+1)C_\ell/2\pi\) self-owned numpy arrays (0..lmax, 0..3), where 0..3 index are TT, EE, BB, TE, unless raw_cl is True in which case return just \(C_\ell\). For the lens_potential the power spectrum returned is that of the deflection.

Note that even if lmax is None, all spectra a returned to the same lmax, appropriate for lensed spectra. Use the individual functions instead if you want to the full unlensed and lensing potential power spectra to the higher lmax actually computed.

Parameters:
  • params – optional CAMBparams instance with parameters to use. If None, must have previously set parameters and called calc_power_spectra (e.g. if you got this instance using camb.get_results()),

  • lmax – maximum L

  • spectra – list of names of spectra to get

  • CMB_unit – scale results from dimensionless. Use ‘muK’ for \(\mu K^2\) units for CMB \(C_\ell\) and \(\mu K\) units for lensing cross.

  • raw_cl – return \(C_\ell\) rather than \(\ell(\ell+1)C_\ell/2\pi\)

Returns:

dictionary of power spectrum arrays, indexed by names of requested spectra

get_cmb_transfer_data(tp='scalar')[source]

Get \(C_\ell\) transfer functions

Returns:

ClTransferData instance holding output arrays (copies, not pointers)

get_cmb_unlensed_scalar_array_dict(params=None, lmax=None, CMB_unit=None, raw_cl=False)[source]

Get all unlensed auto and cross power spectra, including any custom source functions set using model.CAMBparams.set_custom_scalar_sources().

Parameters:
  • params – optional CAMBparams instance with parameters to use. If None, must have previously set parameters and called calc_power_spectra() (e.g. if you got this instance using camb.get_results()),

  • lmax – maximum \(\ell\)

  • CMB_unit – scale results from dimensionless. Use ‘muK’ for \(\mu K^2\) units for CMB \(C_\ell\) and \(\mu K\) units for lensing cross.

  • raw_cl – return \(C_\ell\) rather than \(\ell(\ell+1)C_\ell/2\pi\)

Returns:

dictionary of power spectrum arrays, index as TxT, TxE, PxW1, W1xW2, custom_name_1xT… etc. Note that P is the lensing deflection, lensing windows Wx give convergence.

get_dark_energy_rho_w(a)[source]

Get dark energy density in units of the dark energy density today, and equation of state parameter \(w\equiv P/\rho\)

Parameters:

a – scalar factor or array of scale factors

Returns:

rho, w arrays at redshifts \(1/a-1\) [or scalars if \(a\) is scalar]

get_derived_params()[source]
Returns:

dictionary of derived parameter values, indexed by name (‘kd’, ‘age’, etc..)

get_fsigma8()[source]

Get \(f\sigma_8\) growth values (must previously have calculated power spectra). For general models \(f\sigma_8\) is defined as in the Planck 2015 parameter paper in terms of the velocity-density correlation: \(\sigma^2_{vd}/\sigma_{dd}\) for \(8 h^{-1} {\rm Mpc}\) spheres.

Returns:

array of f*sigma_8 values, in order of increasing time (decreasing redshift)

get_lens_potential_cls(lmax=None, CMB_unit=None, raw_cl=False)[source]

Get lensing deflection angle potential power spectrum, and cross-correlation with T and E. Must have already calculated power spectra. Power spectra are \([L(L+1)]^2C_L^{\phi\phi}/2\pi\) and corresponding deflection cross-correlations.

Parameters:
  • lmax – lmax to output to

  • CMB_unit – scale results from dimensionless. Use ‘muK’ for \(\mu K\) units for lensing cross.

  • raw_cl – return lensing potential \(C_L\) rather than \([L(L+1)]^2C_L/2\pi\)

Returns:

numpy array CL[0:lmax+1,0:3], where 0..2 indexes PP, PT, PE.

get_lensed_cls_with_spectrum(clpp, lmax=None, CMB_unit=None, raw_cl=False)[source]

Get lensed CMB power spectra using curved-sky correlation function method, using cpp as the lensing spectrum. Useful for e.g. getting partially-delensed spectra.

Parameters:
  • clpp – array of \([L(L+1)]^2 C_L^{\phi\phi}/2\pi\) lensing potential power spectrum (zero based)

  • lmax – lmax to output to

  • CMB_unit – scale results from dimensionless. Use ‘muK’ for \(\mu K^2\) units for CMB \(C_\ell\)

  • raw_cl – return \(C_\ell\) rather than \(\ell(\ell+1)C_\ell/2\pi\)

Returns:

numpy array CL[0:lmax+1,0:4], where 0..3 indexes TT, EE, BB, TE.

get_lensed_gradient_cls(lmax=None, CMB_unit=None, raw_cl=False, clpp=None)[source]

Get lensed gradient scalar CMB power spectra in flat sky approximation (arXiv:1101.2234). Note that lmax used to calculate results may need to be substantially larger than the lmax output from this function (there is no extrapolation as in the main lensing routines). Lensed power spectra must be already calculated.

Parameters:
  • lmax – lmax to output to

  • CMB_unit – scale results from dimensionless. Use ‘muK’ for \(\mu K^2\) units for CMB \(C_\ell\)

  • raw_cl – return \(C_\ell\) rather than \(\ell(\ell+1)C_\ell/2\pi\)

  • clpp – custom array of \([L(L+1)]^2 C_L^{\phi\phi}/2\pi\) lensing potential power spectrum to use (zero based), rather than calculated specturm from this model

Returns:

numpy array CL[0:lmax+1,0:8], where CL[:,i] are \(T\nabla T\), \(E\nabla E\), \(B\nabla B\), \(PP_\perp\), \(T\nabla E\), \(TP_\perp\), \((\nabla T)^2\), \(\nabla T\nabla T\) where the first six are as defined in appendix C of 1101.2234.

get_lensed_scalar_cls(lmax=None, CMB_unit=None, raw_cl=False)[source]

Get lensed scalar CMB power spectra. Must have already calculated power spectra.

Parameters:
  • lmax – lmax to output to

  • CMB_unit – scale results from dimensionless. Use ‘muK’ for \(\mu K^2\) units for CMB \(C_\ell\)

  • raw_cl – return \(C_\ell\) rather than \(\ell(\ell+1)C_\ell/2\pi\)

Returns:

numpy array CL[0:lmax+1,0:4], where 0..3 indexes TT, EE, BB, TE.

get_linear_matter_power_spectrum(var1=None, var2=None, hubble_units=True, k_hunit=True, have_power_spectra=True, params=None, nonlinear=False)[source]

Calculates \(P_{xy}(k)\), where x, y are one of model.Transfer_cdm, model.Transfer_xx etc. The output k values are not regularly spaced, and not interpolated. They are either k or k/h depending on the value of k_hunit (default True gives k/h).

For a description of outputs for different var1, var2 see Matter power spectrum and matter transfer function variables.

Parameters:
  • var1 – variable i (index, or name of variable; default delta_tot)

  • var2 – variable j (index, or name of variable; default delta_tot)

  • hubble_units – if true, output power spectrum in (Mpc/h) units, otherwise Mpc

  • k_hunit – if true, matter power is a function of k/h, if false, just k (both \({\rm Mpc}^{-1}\) units)

  • have_power_spectra – set to False if not already computed power spectra

  • params – if have_power_spectra=False, optional CAMBparams instance to specify new parameters

  • nonlinear – include non-linear correction from halo model

Returns:

k/h or k, z, PK, where kz and z are arrays of k/h or k and z respectively, and PK[i,j] is the value at z[i], k[j]/h or k[j]

get_matter_power_interpolator(nonlinear=True, var1=None, var2=None, hubble_units=True, k_hunit=True, return_z_k=False, log_interp=True, extrap_kmax=None, silent=False)[source]

Assuming transfers have been calculated, return a 2D spline interpolation object to evaluate matter power spectrum as function of z and k/h (or k). Uses self.Params.Transfer.PK_redshifts as the spline node points in z. If fewer than four redshift points are used the interpolator uses a reduced order spline in z (so results at intermediate z may be innaccurate), otherwise it uses bicubic. Usage example:

PK = results.get_matter_power_interpolator();
print('Power spectrum at z=0.5, k/h=0.1 is %s (Mpc/h)^3 '%(PK.P(0.5, 0.1)))

For a description of outputs for different var1, var2 see Matter power spectrum and matter transfer function variables.

Parameters:
  • nonlinear – include non-linear correction from halo model

  • var1 – variable i (index, or name of variable; default delta_tot)

  • var2 – variable j (index, or name of variable; default delta_tot)

  • hubble_units – if true, output power spectrum in \(({\rm Mpc}/h)^{3}\) units, otherwise \({\rm Mpc}^{3}\)

  • k_hunit – if true, matter power is a function of k/h, if false, just k (both \({\rm Mpc}^{-1}\) units)

  • return_z_k – if true, return interpolator, z, k where z, k are the grid used

  • log_interp – if true, interpolate log of power spectrum (unless any values cross zero in which case ignored)

  • extrap_kmax – if set, use power law extrapolation beyond kmax to extrap_kmax (useful for tails of integrals)

  • silent – Set True to silence warnings

Returns:

An object PK based on RectBivariateSpline, that can be called with PK.P(z,kh) or PK(z,log(kh)) to get log matter power values. If return_z_k=True, instead return interpolator, z, k where z, k are the grid used.

get_matter_power_spectrum(minkh=0.0001, maxkh=1.0, npoints=100, var1=None, var2=None, have_power_spectra=False, params=None)[source]

Calculates \(P_{xy}(k/h)\), where x, y are one of Transfer_cdm, Transfer_xx etc. The output k values are regularly log spaced and interpolated. If NonLinear is set, the result is non-linear.

For a description of outputs for different var1, var2 see Matter power spectrum and matter transfer function variables.

Parameters:
  • minkh – minimum value of k/h for output grid (very low values < 1e-4 may not be calculated)

  • maxkh – maximum value of k/h (check consistent with input params.Transfer.kmax)

  • npoints – number of points equally spaced in log k

  • var1 – variable i (index, or name of variable; default delta_tot)

  • var2 – variable j (index, or name of variable; default delta_tot)

  • have_power_spectra – set to True if already computed power spectra

  • params – if have_power_spectra=False and want to specify new parameters, a CAMBparams instance

Returns:

kh, z, PK, where kz and z are arrays of k/h and z respectively, and PK[i,j] is value at z[i], k/h[j]

get_matter_transfer_data() MatterTransferData[source]

Get matter transfer function data and sigma8 for calculated results.

Returns:

MatterTransferData instance holding output arrays (copies, not pointers)

get_nonlinear_matter_power_spectrum(var1=None, var2=None, hubble_units=True, k_hunit=True, have_power_spectra=True, params=None)[source]

Calculates \(P_{xy}(k/h)\), where x, y are one of model.Transfer_cdm, model.Transfer_xx etc. The output k values are not regularly spaced, and not interpolated.

For a description of outputs for different var1, var2 see Matter power spectrum and matter transfer function variables.

Parameters:
  • var1 – variable i (index, or name of variable; default delta_tot)

  • var2 – variable j (index, or name of variable; default delta_tot)

  • hubble_units – if true, output power spectrum in \(({\rm Mpc}/h)^{3}\) units, otherwise \({\rm Mpc}^{3}\)

  • k_hunit – if true, matter power is a function of k/h, if false, just k (both \({\rm Mpc}^{-1}\) units)

  • have_power_spectra – set to False if not already computed power spectra

  • params – if have_power_spectra=False, optional CAMBparams instance to specify new parameters

Returns:

k/h or k, z, PK, where kz and z are arrays of k/h or k and z respectively, and PK[i,j] is the value at z[i], k[j]/h or k[j]

get_partially_lensed_cls(Alens, lmax=None, CMB_unit=None, raw_cl=False)[source]

Get lensed CMB power spectra using curved-sky correlation function method, using true lensing spectrum scaled by Alens. Alens can be an array in L for realistic delensing estimates. Note that if Params.Alens is also set, the result is scaled by the product of both

Parameters:
  • Alens – scaling of the lensing relative to true, with Alens=1 being the standard result. Can be a scalar in which case all L are scaled, or a zero-based array with the L by L scaling (with L larger than the size of the array having Alens_L=1).

  • lmax – lmax to output to

  • CMB_unit – scale results from dimensionless. Use ‘muK’ for \(\mu K^2\) units for CMB \(C_\ell\)

  • raw_cl – return \(C_\ell\) rather than \(\ell(\ell+1)C_\ell/2\pi\)

Returns:

numpy array CL[0:lmax+1,0:4], where 0..3 indexes TT, EE, BB, TE.

get_redshift_evolution(q, z, vars=['k/h', 'delta_cdm', 'delta_baryon', 'delta_photon', 'delta_neutrino', 'delta_nu', 'delta_tot', 'delta_nonu', 'delta_tot_de', 'Weyl', 'v_newtonian_cdm', 'v_newtonian_baryon', 'v_baryon_cdm', 'a', 'etak', 'H', 'growth', 'v_photon', 'pi_photon', 'E_2', 'v_neutrino', 'T_source', 'E_source', 'lens_potential_source'], lAccuracyBoost=4)[source]

Get the mode evolution as a function of redshift for some k values.

Parameters:
  • q – wavenumber values to calculate (or array of k values)

  • z – array of redshifts to output

  • vars – list of variable names or camb.symbolic sympy expressions to output

  • lAccuracyBoost – boost factor for ell accuracy (e.g. to get nice smooth curves for plotting)

Returns:

nd array, A_{qti}, size(q) x size(times) x len(vars), or 2d array if q is scalar

get_sigma8()[source]

Get \(\sigma_8\) values at Params.PK_redshifts (must previously have calculated power spectra)

Returns:

array of \(\sigma_8\) values, in order of increasing time (decreasing redshift)

get_sigma8_0()[source]

Get \(\sigma_8\) value today (must previously have calculated power spectra)

Returns:

\(\sigma_8\) today

get_sigmaR(R, z_indices=None, var1=None, var2=None, hubble_units=True, return_R_z=False)[source]

Calculate \(\sigma_R\) values, the RMS linear matter fluctuation in spheres of radius R in linear theory. Accuracy depends on the sampling with which the matter transfer functions are computed.

For a description of outputs for different var1, var2 see Matter power spectrum and matter transfer function variables. Note that numerical errors are slightly different to get_sigma8 for R=8 Mpc/h.

Parameters:
  • R – radius in Mpc or h^{-1} Mpc units, scalar or array

  • z_indices – indices of redshifts in Params.Transfer.PK_redshifts to calculate (default None gives all computed in order of increasing time (decreasing redshift); -1 gives redshift 0; list gives all listed indices)

  • var1 – variable i (index, or name of variable; default delta_tot)

  • var2 – variable j (index, or name of variable; default delta_tot)

  • hubble_units – if true, R is in h^{-1} Mpc, otherwise Mpc

  • return_R_z – if true, return tuple of R, z, sigmaR (where R always Mpc units not h^{-1}Mpc and R, z are arrays)

Returns:

array of \(\sigma_R\) values, or 2D array indexed by (redshift, R)

get_source_cls_dict(params=None, lmax=None, raw_cl=False)[source]

Get all source window function and CMB lensing and cross power spectra. Does not include CMB spectra. Note that P is the deflection angle, but lensing windows return the kappa power.

Parameters:
  • params – optional CAMBparams instance with parameters to use. If None, must have previously set parameters and called calc_power_spectra() (e.g. if you got this instance using camb.get_results()),

  • lmax – maximum \(\ell\)

  • raw_cl – return \(C_\ell\) rather than \(\ell(\ell+1)C_\ell/2\pi\)

Returns:

dictionary of power spectrum arrays, index as PXP, PxW1, W1xW2, … etc.

get_tensor_cls(lmax=None, CMB_unit=None, raw_cl=False)[source]

Get tensor CMB power spectra. Must have already calculated power spectra.

Parameters:
  • lmax – lmax to output to

  • CMB_unit – scale results from dimensionless. Use ‘muK’ for \(\mu K^2\) units for CMB \(C_\ell\)

  • raw_cl – return \(C_\ell\) rather than \(\ell(\ell+1)C_\ell/2\pi\)

Returns:

numpy array CL[0:lmax+1,0:4], where 0..3 indexes TT, EE, BB, TE

get_time_evolution(q, eta, vars=['k/h', 'delta_cdm', 'delta_baryon', 'delta_photon', 'delta_neutrino', 'delta_nu', 'delta_tot', 'delta_nonu', 'delta_tot_de', 'Weyl', 'v_newtonian_cdm', 'v_newtonian_baryon', 'v_baryon_cdm', 'a', 'etak', 'H', 'growth', 'v_photon', 'pi_photon', 'E_2', 'v_neutrino', 'T_source', 'E_source', 'lens_potential_source'], lAccuracyBoost=4, frame='CDM')[source]

Get the mode evolution as a function of conformal time for some k values.

Note that gravitational potentials (e.g. Weyl) are not integrated in the code and are calculated as derived parameters; they may be numerically unstable far outside the horizon. (use the series expansion result if needed far outside the horizon)

Parameters:
  • q – wavenumber values to calculate (or array of k values)

  • eta – array of requested conformal times to output

  • vars – list of variable names or sympy symbolic expressions to output (using camb.symbolic)

  • lAccuracyBoost – factor by which to increase l_max in hierarchies compared to default - often needed to get nice smooth curves of acoustic oscillations for plotting.

  • frame – for symbolic expressions, can specify frame name if the variable is not gauge invariant. e.g. specifying Delta_g and frame=’Newtonian’ would give the Newtonian gauge photon density perturbation.

Returns:

nd array, A_{qti}, size(q) x size(times) x len(vars), or 2d array if q is scalar

get_total_cls(lmax=None, CMB_unit=None, raw_cl=False)[source]

Get lensed-scalar + tensor CMB power spectra. Must have already calculated power spectra.

Parameters:
  • lmax – lmax to output to

  • CMB_unit – scale results from dimensionless. Use ‘muK’ for \(\mu K^2\) units for CMB \(C_\ell\)

  • raw_cl – return \(C_\ell\) rather than \(\ell(\ell+1)C_\ell/2\pi\)

Returns:

numpy array CL[0:lmax+1,0:4], where 0..3 indexes TT, EE, BB, TE

get_unlensed_scalar_array_cls(lmax=None)[source]

Get array of all cross power spectra. Must have already calculated power spectra. Results are dimensionless, and not scaled by custom_scaled_ell_fac.

Parameters:

lmax – lmax to output to

Returns:

numpy array CL[0:, 0:,0:lmax+1], where 0.. index T, E, lensing potential, source window functions

get_unlensed_scalar_cls(lmax=None, CMB_unit=None, raw_cl=False)[source]

Get unlensed scalar CMB power spectra. Must have already calculated power spectra.

Parameters:
  • lmax – lmax to output to

  • CMB_unit – scale results from dimensionless. Use ‘muK’ for \(\mu K^2\) units for CMB \(C_\ell\)

  • raw_cl – return \(C_\ell\) rather than \(\ell(\ell+1)C_\ell/2\pi\)

Returns:

numpy array CL[0:lmax+1,0:4], where 0..3 indexes TT, EE, BB, TE. CL[:,2] will be zero.

get_unlensed_total_cls(lmax=None, CMB_unit=None, raw_cl=False)[source]

Get unlensed CMB power spectra, including tensors if relevant. Must have already calculated power spectra.

Parameters:
  • lmax – lmax to output to

  • CMB_unit – scale results from dimensionless. Use ‘muK’ for \(\mu K^2\) units for CMB \(C_\ell\)

  • raw_cl – return \(C_\ell\) rather than \(\ell(\ell+1)C_\ell/2\pi\)

Returns:

numpy array CL[0:lmax+1,0:4], where 0..3 indexes TT, EE, BB, TE.

h_of_z(z)[source]

Get Hubble rate at redshift z, in \({\rm Mpc}^{-1}\) units, scalar or array

Must have called calc_background(), calc_background_no_thermo() or calculated transfer functions or power spectra.

Use hubble_parameter instead if you want in [km/s/Mpc] units.

Parameters:

z – redshift

Returns:

H(z)

hubble_parameter(z)[source]

Get Hubble rate at redshift z, in km/s/Mpc units. Scalar or array.

Must have called calc_background(), calc_background_no_thermo() or calculated transfer functions or power spectra.

Parameters:

z – redshift

Returns:

H(z)/[km/s/Mpc]

luminosity_distance(z)[source]

Get luminosity distance from to redshift z.

Must have called calc_background(), calc_background_no_thermo() or calculated transfer functions or power spectra.

Parameters:

z – redshift or array of redshifts

Returns:

luminosity distance (matches rank of z)

physical_time(z)[source]

Get physical time from hot big bang to redshift z in Julian Gigayears.

Parameters:

z – redshift

Returns:

t(z)/Gigayear

physical_time_a1_a2(a1, a2)[source]

Get physical time between two scalar factors in Julian Gigayears

Must have called calc_background(), calc_background_no_thermo() or calculated transfer functions or power spectra.

Parameters:
  • a1 – scale factor 1

  • a2 – scale factor 2

Returns:

(age(a2)-age(a1))/Gigayear

power_spectra_from_transfer(initial_power_params=None, silent=False)[source]

Assuming calc_transfers() or calc_power_spectra() have already been used, re-calculate the power spectra using a new set of initial power spectrum parameters with otherwise the same cosmology. This is typically much faster that re-calculating everything, as the transfer functions can be re-used. NOTE: if non-linear lensing is on, the transfer functions have the non-linear correction included when they are calculated, so using this function with a different initial power spectrum will not give quite the same results as doing a full recalculation.

Parameters:
redshift_at_comoving_radial_distance(chi)[source]

Convert comoving radial distance array to redshift array.

Parameters:

chi – comoving radial distance (in Mpc), scalar or array

Returns:

redshift at chi, scalar or array

redshift_at_conformal_time(eta)[source]

Convert conformal time array to redshift array. Note that this function requires the transfers or background to have been calculated with no_thermo=False (the default).

Parameters:

eta – conformal time from bing bang (in Mpc), scalar or array

Returns:

redshift at eta, scalar or array

replace(instance)

Replace the content of this class with another instance, doing a deep copy (in Fortran)

Parameters:

instance – instance of the same class to replace this instance with

save_cmb_power_spectra(filename, lmax=None, CMB_unit='muK')[source]

Save CMB power to a plain text file. Output is lensed total \(\ell(\ell+1)C_\ell/2\pi\) then lensing potential and cross: L TT EE BB TE PP PT PE.

Parameters:
  • filename – filename to save

  • lmax – lmax to save

  • CMB_unit – scale results from dimensionless. Use ‘muK’ for \(\mu K^2\) units for CMB \(C_\ell\) and \(\mu K\) units for lensing cross.

set_params(params)[source]

Set parameters from params. Note that this does not recompute anything; you will need to call calc_transfers() if you change any parameters affecting the background cosmology or the transfer function settings.

Parameters:

params – a CAMBparams instance

sound_horizon(z)[source]

Get comoving sound horizon as function of redshift in Megaparsecs, the integral of the sound speed up to given redshift.

Parameters:

z – redshift or array of redshifts

Returns:

r_s(z)

class camb.results.MatterTransferData[source]

MatterTransferData is the base class for storing matter power transfer function data for various q values. In a flat universe q=k, in a closed universe q is quantized.

To get an instance of this data, call results.CAMBdata.get_matter_transfer_data().

For a description of the different Transfer_xxx outputs (and 21cm case) see Matter power spectrum and matter transfer function variables; the array is indexed by index+1 gven by:

  • Transfer_kh = 1 (k/h)

  • Transfer_cdm = 2 (cdm)

  • Transfer_b = 3 (baryons)

  • Transfer_g = 4 (photons)

  • Transfer_r = 5 (massless neutrinos)

  • Transfer_nu = 6 (massive neutrinos)

  • Transfer_tot = 7 (total matter)

  • Transfer_nonu = 8 (total matter excluding neutrinos)

  • Transfer_tot_de = 9 (total including dark energy perturbations)

  • Transfer_Weyl = 10 (Weyl potential)

  • Transfer_Newt_vel_cdm = 11 (Newtonian CDM velocity)

  • Transfer_Newt_vel_baryon = 12 (Newtonian baryon velocity)

  • Transfer_vel_baryon_cdm = 13 (relative baryon-cdm velocity)

Variables:
  • nq – number of q modes calculated

  • q – array of q values calculated

  • sigma_8 – array of \(\sigma_8\) values for each redshift

  • sigma2_vdelta_8 – array of v-delta8 correlation, so sigma2_vdelta_8/sigma_8 can define growth

  • transfer_data – numpy array T[entry, q_index, z_index] storing transfer functions for each redshift and q; entry+1 can be one of the Transfer_xxx variables above.

transfer_z(name, z_index=0)[source]

Get transfer function (function of q, for each q in self.q_trans) by name for given redshift index

Parameters:
  • name – parameter name

  • z_index – which redshift

Returns:

array of transfer function values for each calculated k

class camb.results.ClTransferData[source]

ClTransferData is the base class for storing CMB power transfer functions, as a function of q and \(\ell\). To get an instance of this data, call results.CAMBdata.get_cmb_transfer_data()

Variables:
  • NumSources – number of sources calculated (size of p index)

  • q – array of q values calculated (=k in flat universe)

  • L – int array of \(\ell\) values calculated

  • delta_p_l_k – transfer functions, indexed by source, L, q

get_transfer(source=0)[source]

Return \(C_\ell\) trasfer functions as a function of \(\ell\) and \(q\) (\(= k\) in a flat universe).

Parameters:

source – index of source: e.g. 0 for temperature, 1 for E polarization, 2 for lensing potential

Returns:

array of computed L, array of computed q, transfer functions T(L,q)