Input parameter model (main modifications are here)

class isitgr.model.CAMBparams(**kwargs)[source]

Object storing the parameters for a CAMB calculation, including cosmological parameters and settings for what to calculate. When a new object is instantiated, default parameters are set automatically.

To add a new parameter, add it to the CAMBparams type in model.f90, then edit the _fields_ list in the CAMBparams class in model.py to add the new parameter in the corresponding location of the member list. After rebuilding the python version you can then access the parameter by using params.new_parameter_name where params is a CAMBparams instance. You could also modify the wrapper functions to set the field value less directly.

You can view the set of underlying parameters used by the Fortran code by printing the CAMBparams instance. In python, to set cosmology parameters it is usually best to use set_cosmology() and equivalent methods for most other parameters. Alternatively the convenience function isitgr.set_params() can construct a complete instance from a dictionary of relevant parameters.

Variables:
  • WantCls – (boolean) Calculate C_L
  • WantTransfer – (boolean) Calculate matter transfer functions and matter power spectrum
  • WantScalars – (boolean) Calculates scalar modes
  • WantTensors – (boolean) Calculate tensor modes
  • WantVectors – (boolean) Calculate vector modes
  • WantDerivedParameters – (boolean) Calculate derived parameters
  • Want_cl_2D_array – (boolean) For the C_L, include NxN matrix of all possible cross-spectra between sources
  • Want_CMB – (boolean) Calculate the temperature and polarization power spectra
  • Want_CMB_lensing – (boolean) Calculate the lensing potential power spectrum
  • DoLensing – (boolean) Include CMB lensing
  • NonLinear – (integer/string, one of: NonLinear_none, NonLinear_pk, NonLinear_lens, NonLinear_both)
  • Transferisitgr.model.TransferParams
  • want_zstar – (boolean)
  • want_zdrag – (boolean)
  • min_l – (integer) l_min for the scalar C_L (1 or 2, L=1 dipoles are Newtonian Gauge)
  • max_l – (integer) l_max for the scalar C_L
  • max_l_tensor – (integer) l_max for the tensor C_L
  • max_eta_k – (float64) Maximum k*eta_0 for scalar C_L, where eta_0 is the conformal time today
  • max_eta_k_tensor – (float64) Maximum k*eta_0 for tensor C_L, where eta_0 is the conformal time today
  • ombh2 – (float64) Omega_baryon h^2
  • omch2 – (float64) Omega_cdm h^2
  • omk – (float64) Omega_K
  • omnuh2 – (float64) Omega_massive_neutrino h^2
  • H0 – (float64) Hubble parameter is km/s/Mpc units
  • TCMB – (float64) CMB temperature today in Kelvin
  • YHe – (float64) Helium mass fraction
  • num_nu_massless – (float64) Effective number of massless neutrinos
  • num_nu_massive – (integer) Total physical (integer) number of massive neutrino species
  • nu_mass_eigenstates – (integer) Number of non-degenerate mass eigenstates
  • share_delta_neff – (boolean) Share the non-integer part of num_nu_massless between the eigenstates
  • nu_mass_degeneracies – (float64 array) Degeneracy of each distinct eigenstate
  • nu_mass_fractions – (float64 array) Mass fraction in each distinct eigenstate
  • nu_mass_numbers – (integer array) Number of physical neutrinos per distinct eigenstate
  • InitPowerisitgr.initialpower.InitialPower
  • Recombisitgr.recombination.RecombinationModel
  • Reionisitgr.reionization.ReionizationModel
  • DarkEnergyisitgr.dark_energy.DarkEnergyModel
  • NonLinearModelisitgr.nonlinear.NonLinearModel
  • Accuracyisitgr.model.AccuracyParams
  • SourceTermsisitgr.model.SourceTermParams
  • z_outputs – (float64 array) redshifts to always calculate BAO output parameters
  • scalar_initial_condition – (integer/string, one of: initial_vector, initial_adiabatic, initial_iso_CDM, initial_iso_baryon, initial_iso_neutrino, initial_iso_neutrino_vel)
  • InitialConditionVector – (float64 array) if scalar_initial_condition is initial_vector, the vector of initial condition amplitudes
  • OutputNormalization – (integer) If non-zero, multipole to normalize the C_L at
  • Alens – (float64) non-physical scaling amplitude for the CMB lensing spectrum power
  • MassiveNuMethod – (integer/string, one of: Nu_int, Nu_trunc, Nu_approx, Nu_best)
  • DoLateRadTruncation – (boolean) If true, use smooth approx to radiation perturbations after decoupling on small scales, saving evolution of irrelevant oscillatory multipole equations
  • Evolve_baryon_cs – (boolean) Evolve a separate equation for the baryon sound speed rather than using background approximation
  • Evolve_delta_xe – (boolean) Evolve ionization fraction perturbations
  • Evolve_delta_Ts – (boolean) Evolve the spin temperature perturbation (for 21cm)
  • Do21cm – (boolean) 21cm is not yet implemented via the python wrapper
  • transfer_21cm_cl – (boolean) Get 21cm C_L at a given fixed redshift
  • Log_lvalues – (boolean) Use log spacing for sampling in L
  • use_cl_spline_template – (boolean) When interpolating use a fiducial spectrum shape to define ratio to spline
  • SourceWindows – array of isitgr.sources.SourceWindow
  • CustomSourcesisitgr.model.CustomSources
  • Q0 – (float64) MG parameter for (Q,D) and (Q,R) parametrizations as used in arXiv:1002.4197v4, or also as described in Table VI in arXiv:1908.00290
  • Qinf – (float64) MG parameter for (Q,D) and (Q,R) parametrizations as used in arXiv:1002.4197v4, or also as described in Table VI in arXiv:1908.00290
  • DR0 – (float64) MG parameter for (Q,D) and (Q,R) parametrizations as used in arXiv:1002.4197v4, or also as described in Table VI in arXiv:1908.00290
  • DRinf – (float64) MG parameter for (Q,D) and (Q,R) parametrizations as used in arXiv:1002.4197v4, or also as described in Table VI in arXiv:1908.00290
  • s – (float64) MG parameter for a^s dependence
  • k_c – (float64) parameter that separates scale-bins
  • Q1 – (float64) Bin parameter for (Q,D) parameterization
  • Q2 – (float64) Bin parameter for (Q,D) parameterization
  • Q3 – (float64) Bin parameter for (Q,D) parameterization
  • Q4 – (float64) Bin parameter for (Q,D) parameterization
  • D1 – (float64) Bin parameter for (Q,D) parameterization
  • D2 – (float64) Bin parameter for (Q,D) parameterization
  • D3 – (float64) Bin parameter for (Q,D) parameterization
  • D4 – (float64) Bin parameter for (Q,D) parameterization
  • z_div – (float64) redshift at which the transition between the two bins occurs
  • z_TGR – (float64) redshift below which GR is to be tested
  • z_tw – (float64) transition width for the hyperbolic tangent function.
  • k_tw – (float64) transition width between k bins
  • t_dep – (boolean) flag for a^s dependence for (Q,D) and (Q,R)
  • R_func – (boolean) flag to swith (Q,D) for (Q,R) parametrization
  • true_bin – (boolean) flag to use traditional binning instead of hybrid binning
  • ISiTGR_BIN – (boolean) flag to use binning method instead of functional form
  • E11 – (float64) MG parameter for (mu,eta) parametrization
  • E22 – (float64) MG parameter for (mu,eta) parametrization
  • mu0 – (float64) MG parameter for (mu,Sigma) parametrization
  • Sigma0 – (float64) MG parameter for (mu,Sigma) parametrization
  • c1 – (float64) MG parameter for scale-dependence
  • c2 – (float64) MG parameter for scale-dependence
  • Lambda – (float64) MG parameter for scale-dependence
  • mu1 – (float64) Bin parameter for (mu,eta) or (mu,Sigma) parameterization
  • mu2 – (float64) Bin parameter for (mu,eta) or (mu,Sigma) parameterization
  • mu3 – (float64) Bin parameter for (mu,eta) or (mu,Sigma) parameterization
  • mu4 – (float64) Bin parameter for (mu,eta) or (mu,Sigma) parameterization
  • eta1 – (float64) Bin parameter for (mu,eta) parameterization
  • eta2 – (float64) Bin parameter for (mu,eta) parameterization
  • eta3 – (float64) Bin parameter for (mu,eta) parameterization
  • eta4 – (float64) Bin parameter for (mu,eta) parameterization
  • Sigma1 – (float64) Bin parameter for (mu,Sigma) parameterization
  • Sigma2 – (float64) Bin parameter for (mu,Sigma) parameterization
  • Sigma3 – (float64) Bin parameter for (mu,Sigma) parameterization
  • Sigma4 – (float64) Bin parameter for (mu,Sigma) parameterization
  • ISiTGR_mueta – (boolean) flag to use (mu,eta) parametrization for functional form
  • ISiTGR_muSigma – (boolean) flag to use (mu,Sigma) parametrization for functional form
  • ISiTGR_QDR – (boolean) flag to use (Q,D) or (Q,R) parametrizations for functional form (flag ‘R_func’ decides if D or R is used)
  • ISiTGR_BIN_mueta – (boolean) flag to use (mu,eta) parametrization for binning method
  • ISiTGR_BIN_muSigma – (boolean) flag to use (mu,Sigma) parametrization for binning method
  • GR – (integer) GR switch on/off
N_eff
Returns:Effective number of degrees of freedom in relativistic species at early times.
copy()

Make independent copy of this object.

Returns:deep copy of self
diff(params)[source]

Print differences between this set of parameters and params

Parameters:params – another CAMBparams instance
get_DH(ombh2=None, delta_neff=None)[source]

Get deuterium ration D/H by intepolation using the bbn.BBNPredictor instance passed to set_cosmology() (or the default one, if Y_He has not been set).

Parameters:
  • ombh2\(\Omega_b h^2\) (default: value passed to set_cosmology())
  • delta_neff – additional \(N_{\rm eff}\) relative to standard value (of 3.046) (default: from values passed to set_cosmology())
Returns:

BBN helium nucleon fraction D/H

get_Y_p(ombh2=None, delta_neff=None)[source]

Get BBN helium nucleon fraction (NOT the same as the mass fraction Y_He) by intepolation using the bbn.BBNPredictor instance passed to set_cosmology() (or the default one, if Y_He has not been set).

Parameters:
  • ombh2\(\Omega_b h^2\) (default: value passed to set_cosmology())
  • delta_neff – additional \(N_{\rm eff}\) relative to standard value (of 3.046) (default: from values passed to set_cosmology())
Returns:

\(Y_p^{\rm BBN}\) helium nucleon fraction predicted by BBN.

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
scalar_power(k)[source]

Get the primordial scalar curvature power spectrum at \(k\)

Parameters:k – wavenumber \(k\) (in \({\rm Mpc}^{-1}\) units)
Returns:power spectrum at \(k\)
set_H0_for_theta(theta, cosmomc_approx=False, theta_H0_range=(10, 100), est_H0=67.0, iteration_threshold=8)[source]

Set H0 to give a specified value of the acoustic angular scale parameter theta.

Parameters:
  • theta – value of \(r_s/D_M\) at redshift \(z_\star\)
  • cosmomc_approx – if true, use approximate fitting formula for \(z_\star\), if false do full numerical calculation
  • theta_H0_range – min, max iterval to search for H0 (in km/s/Mpc)
  • est_H0 – an initial guess for H0 in km/s/Mpc, used in the case comsomc_approx=False.
  • iteration_threshold – differnce in H0 from est_H0 for which to iterate, used for cosmomc_approx=False
set_accuracy(AccuracyBoost=1.0, lSampleBoost=1.0, lAccuracyBoost=1.0, DoLateRadTruncation=True)[source]

Set parameters determining overall calculation accuracy (large values may give big slow down). For finer control you can set individual accuracy parameters by changing CAMBParams.Accuracy (model.AccuracyParams) .

Parameters:
  • AccuracyBoost – increase AccuracyBoost to decrease integration step size, increase density of k sampling, etc.
  • lSampleBoost – increase lSampleBoost to increase density of L sampling for CMB
  • lAccuracyBoost – increase lAccuracyBoost to increase the maximum L included in the Boltzmann hierarchies
  • DoLateRadTruncation – If True, use approximation to radiation perturbation evolution at late times
Returns:

self

set_classes(dark_energy_model=None, initial_power_model=None, non_linear_model=None, recombination_model=None)[source]

Change the classes used to implement parts of the model.

Parameters:
  • dark_energy_model – ‘fluid’, ‘ppf’, or name of a DarkEnergyModel class
  • initial_power_model – name of an InitialPower class
  • non_linear_model – name of a NonLinearModel class
  • recombination_model – name of recombination_model class
set_cosmology(H0=None, ombh2=0.022, omch2=0.12, omk=0.0, cosmomc_theta=None, thetastar=None, neutrino_hierarchy='degenerate', num_massive_neutrinos=1, mnu=0.06, nnu=3.046, YHe=None, meffsterile=0.0, standard_neutrino_neff=3.046, TCMB=2.7255, tau=None, deltazrei=None, Alens=1.0, bbn_predictor=None, theta_H0_range=(10, 100), parameterization=None, Q0=0, Qinf=0, D0=0, R0=0, Dinf=0, Rinf=0, s=0, k_c=0, binning=None, E11=0, E22=0, mu0=0, Sigma0=0, c1=1, c2=1, Lambda=0, z_div=0, z_TGR=0, z_tw=0, k_tw=0, Q1=0, Q2=0, Q3=0, Q4=0, D1=0, D2=0, D3=0, D4=0, mu1=0, mu2=0, mu3=0, mu4=0, eta1=0, eta2=0, eta3=0, eta4=0, Sigma1=0, Sigma2=0, Sigma3=0, Sigma4=0)[source]

Sets cosmological parameters in terms of physical densities and parameters (e.g. as used in Planck analyses). Default settings give a single distinct neutrino mass eigenstate, by default one neutrino with mnu = 0.06eV. Set the neutrino_hierarchy parameter to normal or inverted to use a two-eigenstate model that is a good approximation to the known mass splittings seen in oscillation measurements. For more fine-grained control can set the neutrino parameters directly rather than using this function.

Instead of setting the Hubble parameter directly, you can instead set the acoustic scale parameter (cosmomc_theta, which is based on a fitting forumula for simple models, or thetastar, which is numerically calculated more generally). Note that you must have already set the dark energy model, you can’t use set_cosmology with theta and then change the background evolution (which would change theta at the calculated H0 value).Likewise the dark energy model cannot depend explicitly on H0.

Parameters:
  • H0 – Hubble parameter today in km/s/Mpc. Can leave unset and instead set thetastar or cosmomc_theta (which solves for the required H0).
  • ombh2 – physical density in baryons
  • omch2 – physical density in cold dark matter
  • omk – Omega_K curvature parameter
  • cosmomc_theta – The approximate CosmoMC theta parameter \(\theta_{\rm MC}\). The angular diamter distance is calculated numerically, but the redshift \(z_\star\) is calculated using an approximate (quite accurate but non-general) fitting formula. Leave unset to use H0 or thetastar.
  • thetastar – The angular acoustic scale parameter \(\theta_\star = r_s(z_*)/D_M(z_*)\), defined as the ratio of the photon-baryon sound horizon \(r_s\) to the angular diameter distance \(D_M\), where both quantities are evaluated at \(z_*\), the redshift at which the optical depth (excluding reionization) is unity. Leave unset to use H0 or cosmomc_theta.
  • neutrino_hierarchy – ‘degenerate’, ‘normal’, or ‘inverted’ (1 or 2 eigenstate approximation)
  • num_massive_neutrinos – number of massive neutrinos
  • mnu – sum of neutrino masses (in eV). Omega_nu is calculated approximately from this assuming neutrinos non-relativistic today; i.e. here is defined as a direct proxy for Omega_nu. Internally the actual physical mass is calculated from the Omega_nu accounting for small mass-dependent velocity corrections but neglecting spectral distortions to the neutrino distribution. Set the neutrino field values directly if you need finer control or more complex neutrino models.
  • nnu – N_eff, effective relativistic degrees of freedom
  • YHe – Helium mass fraction. If None, set from BBN consistency.
  • meffsterile – effective mass of sterile neutrinos
  • standard_neutrino_neff – default value for N_eff in standard cosmology (non-integer to allow for partial heating of neutrinos at electron-positron annihilation and QED effects)
  • TCMB – CMB temperature (in Kelvin)
  • tau – optical depth; if None, current Reion settings are not changed
  • deltazrei – redshift width of reionization; if None, uses default
  • Alens – (non-physical) scaling of the lensing potential compared to prediction
  • bbn_predictorbbn.BBNPredictor instance used to get YHe from BBN consistency if YHe is None, or name of a BBN predictor class, or file name of an interpolation table
  • theta_H0_range – if thetastar or cosmomc_theta is specified, the min, max interval of H0 values to map to; if H0 is outside this range it will raise an exception.
  • parameterization – set which MG parametrization to be used, None, “mueta”, “muSigma”, “QD”, or “QR”.
  • Q0 – MG parameter for (Q,D) and (Q,R) parametrizations as shown in arXiv:1002.4197v4. Also, it can be used as in eq 26 in arXiv:1908.00290.
  • Qinf – MG parameter for (Q,D) and (Q,R) parametrizations as shown in arXiv:1002.4197v4.
  • D0 – MG parameter for (Q,D) parametrization as in eq 28 in arXiv:1908.00290.
  • R0 – MG parameter for (Q,R) parametrization as in eq 27 in arXiv:1908.00290.
  • Dinf – MG parameter for (Q,D) and (Q,R) parametrizations as shown in arXiv:1002.4197v4.
  • Rinf – MG parameter for (Q,D) and (Q,R) parametrizations as shown in arXiv:1002.4197v4.
  • s – MG parameter for exponent of scale factor (a^s). See Table VI in arXiv:1908.00290.
  • k_c – sets scale factor at which transition between bins in scale occurs.
  • binning – allows ISiTGR to bin MG parameters. Set as “traditional”, “hybrid” or None for functional form.
  • E11 – MG parameter for functional (mu,eta) parametrization. Defined in arXiv:1502.01590v2.
  • E22 – MG parameter for functional (mu,eta) parametrization. Defined in arXiv:1502.01590v2.
  • mu0 – MG parameter for functional (mu,Sigma) parametrization as defined in 1212.3339v2 originally.
  • Sigma0 – MG parameter for functional (mu,Sigma) parametrization as defined in 1212.3339v2 originally.
  • c1 – Scale dependence parameter. See Table VI in arXiv:1908.00290.
  • c2 – Scale dependence parameter. See Table VI in arXiv:1908.00290.
  • Lambda – Scale dependence parameter. See Table VI in arXiv:1908.00290.
  • z_div – transition redshift between bins.
  • z_TGR – redshift below which GR is to be tested.
  • z_tw – transition width redshift.
  • k_tw – transition width for scale.
  • Q1 – MG parameter for binning method, using (Q,D) parametrization.
  • Q2 – MG parameter for binning method, using (Q,D) parametrization.
  • Q3 – MG parameter for binning method, using (Q,D) parametrization.
  • Q4 – MG parameter for binning method, using (Q,D) parametrization.
  • D1 – MG parameter for binning method, using (Q,D) parametrization.
  • D2 – MG parameter for binning method, using (Q,D) parametrization.
  • D3 – MG parameter for binning method, using (Q,D) parametrization.
  • D4 – MG parameter for binning method, using (Q,D) parametrization.
  • mu1 – MG parameter for binning method, using (mu,eta) and (mu,Sigma) parametrization.
  • mu2 – MG parameter for binning method, using (mu,eta) and (mu,Sigma) parametrization.
  • mu3 – MG parameter for binning method, using (mu,eta) and (mu,Sigma) parametrization.
  • mu4 – MG parameter for binning method, using (mu,eta) and (mu,Sigma) parametrization.
  • eta1 – MG parameter for binning method, using (mu,eta) parametrization.
  • eta2 – MG parameter for binning method, using (mu,eta) parametrization.
  • eta3 – MG parameter for binning method, using (mu,eta) parametrization.
  • eta4 – MG parameter for binning method, using (mu,eta) parametrization.
  • Sigma1 – MG parameter for binning method, using (mu,Sigma) parametrization.
  • Sigma2 – MG parameter for binning method, using (mu,Sigma) parametrization.
  • Sigma3 – MG parameter for binning method, using (mu,Sigma) parametrization.
  • Sigma4 – MG parameter for binning method, using (mu,Sigma) parametrization.
set_custom_scalar_sources(custom_sources, source_names=None, source_ell_scales=None, frame='CDM', code_path=None)[source]

Set custom sources for angular power spectrum using isitgr.symbolic sympy expressions.

Parameters:
  • custom_sources – list of sympy expressions for the angular power spectrum sources
  • source_names – optional list of string naes for the sources
  • source_ell_scales – list or dictionary of scalings for each source name, where for integer entry n, the source for multipole \(\ell\) is scalled by \(\sqrt{(\ell+n)!/(\ell-n)!}\), i.e. \(n=2\) for a new polarization-like source.
  • frame – if the source is not gauge invariant, frame in which to interpret result
  • code_path – optional path for output of source code for CAMB f90 source function
set_dark_energy(w=-1.0, cs2=1.0, wa=0, dark_energy_model='fluid')[source]

Set dark energy parameters (use set_dark_energy_w_a to set w(a) from numerical table instead) To use a custom dark energy model, assign the class instance to the DarkEnergy field instead.

Parameters:
  • w\(w\equiv p_{\rm de}/\rho_{\rm de}\), assumed constant
  • wa – evolution of w (for dark_energy_model=ppf)
  • cs2 – rest-frame sound speed squared of dark energy fluid
  • dark_energy_model – model to use (‘fluid’ or ‘ppf’), default is ‘fluid’
Returns:

self

set_dark_energy_w_a(a, w, dark_energy_model='fluid')[source]

Set the dark energy equation of state from tabulated values (which are cubic spline interpolated).

Parameters:
  • a – array of sampled a = 1/(1+z) values
  • w – array of w(a)
  • dark_energy_model – model to use (‘fluid’ or ‘ppf’), default is ‘fluid’
Returns:

self

set_for_lmax(lmax, max_eta_k=None, lens_potential_accuracy=0, lens_margin=150, k_eta_fac=2.5, lens_k_eta_reference=18000.0)[source]

Set parameters to get CMB power spectra accurate to specific a l_lmax. Note this does not fix the actual output L range, spectra may be calculated above l_max (but may not be accurate there). To fix the l_max for output arrays use the optional input argument to results.CAMBdata.get_cmb_power_spectra() etc.

Parameters:
  • lmax\(\ell_{\rm max}\) you want
  • max_eta_k – maximum value of \(k \eta_0\approx k\chi_*\) to use, which indirectly sets k_max. If None, sensible value set automatically.
  • lens_potential_accuracy – Set to 1 or higher if you want to get the lensing potential accurate (1 is only Planck-level accuracy)
  • lens_margin – the \(\Delta \ell_{\rm max}\) to use to ensure lensed \(C_\ell\) are correct at \(\ell_{\rm max}\)
  • k_eta_fac – k_eta_fac default factor for setting max_eta_k = k_eta_fac*lmax if max_eta_k=None
  • lens_k_eta_reference – value of max_eta_k to use when lens_potential_accuracy>0; use k_eta_max = lens_k_eta_reference*lens_potential_accuracy
Returns:

self

set_initial_power(initial_power_params)[source]

Set the InitialPower primordial power spectrum parameters

Parameters:initial_power_paramsinitialpower.InitialPowerLaw or initialpower.SplinedInitialPower instance
Returns:self
set_initial_power_function(P_scalar, P_tensor=None, kmin=1e-06, kmax=100.0, N_min=200, rtol=5e-05, effective_ns_for_nonlinear=None, args=())[source]

Set the initial power spectrum from a function P_scalar(k, *args), and optionally also the tensor spectrum. The function is called to make a pre-computed array which is then interpolated inside CAMB. The sampling in k is set automatically so that the spline is accurate, but you may also need to increase other accuracy parameters.

Parameters:
  • P_scalar – function returning normalized initial scalar curvature power as function of k (in Mpc^{-1})
  • P_tensor – optional function returning normalized initial tensor power spectrum
  • kmin – minimum wavenumber to compute
  • kmax – maximum wavenumber to compute
  • N_min – minimum number of spline points for the pre-computation
  • rtol – relative tolerance for deciding how many points are enough
  • effective_ns_for_nonlinear – an effective n_s for use with approximate non-linear corrections
  • args – optional list of arguments passed to P_scalar (and P_tensor)
Returns:

self

set_initial_power_table(k, pk=None, pk_tensor=None, effective_ns_for_nonlinear=None)[source]

Set a general intial power spectrum from tabulated values. It’s up to you to ensure the sampling of the k values is high enough that it can be interpolated accurately.

Parameters:
  • k – array of k values (Mpc^{-1})
  • pk – array of primordial curvature perturbation power spectrum values P(k_i)
  • pk_tensor – array of tensor spectrum values
  • effective_ns_for_nonlinear – an effective n_s for use with approximate non-linear corrections
set_matter_power(redshifts=(0.0, ), kmax=1.2, k_per_logint=None, nonlinear=None, accurate_massive_neutrino_transfers=False, silent=False)[source]

Set parameters for calculating matter power spectra and transfer functions.

Parameters:
  • redshifts – array of redshifts to calculate
  • kmax – maximum k to calculate (where k is just k, not k/h)
  • k_per_logint – minimum number of k steps per log k. Set to zero to use default optimized spacing.
  • nonlinear – if None, uses existing setting, otherwise boolean for whether to use non-linear matter power.
  • accurate_massive_neutrino_transfers – if you want the massive neutrino transfers accurately
  • silent – if True, don’t give warnings about sort order
Returns:

self

set_nonlinear_lensing(nonlinear)[source]

Settings for whether or not to use non-linear corrections for the CMB lensing potential. Note that set_for_lmax also sets lensing to be non-linear if lens_potential_accuracy>0

Parameters:nonlinear – true to use non-linear corrections
tensor_power(k)[source]

Get the primordial tensor curvature power spectrum at \(k\)

Parameters:k – wavenumber \(k\) (in \({\rm Mpc}^{-1}\) units)
Returns:tensor power spectrum at \(k\)
validate()[source]

Do some quick tests for sanity

Returns:True if OK
class isitgr.model.AccuracyParams[source]

Structure with parameters governing numerical accuracy. AccuracyBoost will also scale almost all the other parameters except for lSampleBoost (which is specific to the output interpolation) and lAccuracyBoost (which is specific to the multipole hierarchy evolution), e.g setting AccuracyBoost=2, IntTolBoost=1.5, means that internally the k sampling for integration will be boosed by AccuracyBoost*IntTolBoost = 3.

Variables:
  • AccuracyBoost – (float64) general accuracy setting effecting everything related to step sizes etc. (including separate settings below except the next two)
  • lSampleBoost – (float64) accuracy for sampling in ell for interpolation for the C_l (if >=50, all ell are calculated)
  • lAccuracyBoost – (float64) Boosts number of multipoles integrated in Boltzman heirarchy
  • AccuratePolarization – (boolean) Do you care about the accuracy of the polarization Cls?
  • AccurateBB – (boolean) Do you care about BB accuracy (e.g. in lensing)
  • AccurateReionization – (boolean) Do you care about pecent level accuracy on EE signal from reionization?
  • TimeStepBoost – (float64) Sampling timesteps
  • BackgroundTimeStepBoost – (float64) Number of time steps for background thermal history and source window interpolation
  • IntTolBoost – (float64) Tolerances for integrating differential equations
  • SourcekAccuracyBoost – (float64) Accuracy of k sampling for source time integration
  • IntkAccuracyBoost – (float64) Accuracy of k sampling for integration
  • TransferkBoost – (float64) Accuracy of k sampling for transfer functions
  • NonFlatIntAccuracyBoost – (float64) Accuracy of non-flat time integration
  • BessIntBoost – (float64) Accuracy of bessel integration truncation
  • LensingBoost – (float64) Accuracy of the lensing of CMB power spectra
  • NonlinSourceBoost – (float64) Accuracy of steps and kmax used for the non-linear correction
  • BesselBoost – (float64) Accuracy of bessel pre-computation sampling
  • LimberBoost – (float64) Accuracy of Limber approximation use
  • SourceLimberBoost – (float64) Scales when to switch to Limber for source windows
  • KmaxBoost – (float64) Boost max k for source window functions
  • neutrino_q_boost – (float64) Number of momenta integrated for neutrino perturbations
class isitgr.model.TransferParams[source]

Object storing parameters for the matter power spectrum calculation.

Variables:
  • high_precision – (boolean) True for more accuracy
  • accurate_massive_neutrinos – (boolean) True if you want neutrino transfer functions accurate (false by default)
  • kmax – (float64) k_max to output (no h in units)
  • k_per_logint – (integer) number of points per log k interval. If zero, set an irregular optimized spacing
  • PK_num_redshifts – (integer) number of redshifts to calculate
  • PK_redshifts – (float64 array) redshifts to output for the matter transfer and power
class isitgr.model.SourceTermParams[source]

Structure with parameters determining how galaxy/lensing/21cm power spectra and transfer functions are calculated.

Variables:
  • limber_windows – (boolean) Use Limber approximation where appropriate. CMB lensing uses Limber even if limber_window is false, but method is changed to be consistent with other sources if limber_windows is true
  • limber_phi_lmin – (integer) When limber_windows=True, the minimum L to use Limber approximation for the lensing potential and other sources (which may use higher but not lower)
  • counts_density – (boolean) Include the density perturbation source
  • counts_redshift – (boolean) Include redshift distortions
  • counts_lensing – (boolean) Include magnification bias for number counts
  • counts_velocity – (boolean) Non-redshift distortion velocity terms
  • counts_radial – (boolean) Radial displacement velocity term; does not include time delay; subset of counts_velocity, just 1 / (chi * H) term
  • counts_timedelay – (boolean) Include time delay terms * 1 / (H * chi)
  • counts_ISW – (boolean) Include tiny ISW terms
  • counts_potential – (boolean) Include tiny terms in potentials at source
  • counts_evolve – (boolean) Accout for source evolution
  • line_phot_dipole – (boolean) Dipole sources for 21cm
  • line_phot_quadrupole – (boolean) Quadrupole sources for 21cm
  • line_basic – (boolean) Include main 21cm monopole density/spin temerature sources
  • line_distortions – (boolean) Redshift distortions for 21cm
  • line_extra – (boolean) Include other sources
  • line_reionization – (boolean) Replace the E modes with 21cm polarization
  • use_21cm_mK – (boolean) Use mK units for 21cm
class isitgr.model.CustomSources[source]

Structure containing symoblic-compiled custom CMB angular power spectrum source functions. Don’t change this directly, instead call model.CAMBparams.set_custom_scalar_sources().

Variables:
  • num_custom_sources – (integer) number of sources set
  • c_source_func – (pointer) Don’t directly change this
  • custom_source_ell_scales – (integer array) scaling in L for outputs