CLASS MANUAL
background.h File Reference
#include "common.h"
#include "quadrature.h"
#include "growTable.h"
#include "arrays.h"
#include "dei_rkck.h"
#include "parser.h"
+ Include dependency graph for background.h:
+ This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Data Structures

struct  background
 
struct  background_parameters_and_workspace
 
struct  background_parameters_for_distributions
 

Enumerations

enum  spatial_curvature
 
enum  equation_of_state
 
enum  varconst_dependence
 
enum  vecback_format
 
enum  interpolation_method
 

Detailed Description

Documented includes for background module


Data Structure Documentation

◆ background

struct background

background structure containing all the background information that other modules need to know.

Once initialized by the backgound_init(), contains all necessary information on the background evolution (except thermodynamics), and in particular, a table of all background quantities as a function of time and scale factor, used for interpolation in other modules.

Data Fields
double H0

$ H_0 $: Hubble parameter (in fact, [ $H_0/c$]) in $ Mpc^{-1} $

double h

reduced Hubble parameter

double Omega0_g

$ \Omega_{0 \gamma} $: photons

double T_cmb

$ T_{cmb} $: current CMB temperature in Kelvins

double Omega0_b

$ \Omega_{0 b} $: baryons

double Omega0_ur

$ \Omega_{0 \nu r} $: ultra-relativistic neutrinos

double Omega0_cdm

$ \Omega_{0 cdm} $: cold dark matter

double Omega0_idm

$ \Omega_{0 idm} $: interacting dark matter with photons, baryons, and idr

double Omega0_idr

$ \Omega_{0 idr} $: interacting dark radiation

double T_idr

$ T_{idr} $: current temperature of interacting dark radiation in Kelvins

double Omega0_dcdmdr

$ \Omega_{0 dcdm}+\Omega_{0 dr} $: decaying cold dark matter (dcdm) decaying to dark radiation (dr)

double Omega_ini_dcdm

$ \Omega_{ini,dcdm} $: rescaled initial value for dcdm density (see 1407.2418 for definitions)

double Gamma_dcdm

$ \Gamma_{dcdm} $: decay constant for decaying cold dark matter

double tau_dcdm
int N_ncdm

Number of distinguishable ncdm species

char * ncdm_psd_files

list of filenames for tabulated p-s-d

int * got_files

list of flags for each species, set to true if p-s-d is passed through file

double * ncdm_psd_parameters

list of parameters for specifying/modifying ncdm p.s.d.'s, to be customized for given model (could be e.g. mixing angles)

double * M_ncdm

vector of masses of non-cold relic: dimensionless ratios m_ncdm/T_ncdm

double * m_ncdm_in_eV

list of ncdm masses in eV (inferred from M_ncdm and other parameters above)

double * Omega0_ncdm
double Omega0_ncdm_tot

Omega0_ncdm for each species and for the total Omega0_ncdm

double * T_ncdm
double T_ncdm_default

list of 1st parameters in p-s-d of non-cold relics: relative temperature T_ncdm1/T_gamma; and its default value

double * ksi_ncdm
double ksi_ncdm_default

list of 2nd parameters in p-s-d of non-cold relics: relative chemical potential ksi_ncdm1/T_ncdm1; and its default value

double * deg_ncdm
double deg_ncdm_default

vector of degeneracy parameters in factor of p-s-d: 1 for one family of neutrinos (= one neutrino plus its anti-neutrino, total g*=1+1=2, so deg = 0.5 g*); and its default value

int * ncdm_input_q_size

Vector of numbers of q bins

double * ncdm_qmax

Vector of maximum value of q

double Omega0_k

$ \Omega_{0_k} $: curvature contribution

double Omega0_lambda

$ \Omega_{0_\Lambda} $: cosmological constant

double Omega0_fld

$ \Omega_{0 de} $: fluid

double Omega0_scf

$ \Omega_{0 scf} $: scalar field

short use_ppf

flag switching on PPF perturbation equations instead of true fluid equations for perturbations. It could have been defined inside perturbation structure, but we leave it here in such way to have all fld parameters grouped.

double c_gamma_over_c_fld

ppf parameter defined in eq. (16) of 0808.3125 [astro-ph]

enum equation_of_state fluid_equation_of_state

parametrisation scheme for fluid equation of state

double w0_fld

$ w0_{DE} $: current fluid equation of state parameter

double wa_fld

$ wa_{DE} $: fluid equation of state parameter derivative

double cs2_fld

$ c^2_{s~DE} $: sound speed of the fluid in the frame comoving with the fluid (so, this is not [delta p/delta rho] in the synchronous or newtonian gauge!)

double Omega_EDE

$ wa_{DE} $: Early Dark Energy density parameter

double * scf_parameters

list of parameters describing the scalar field potential

short attractor_ic_scf

whether the scalar field has attractor initial conditions

int scf_tuning_index

index in scf_parameters used for tuning

double phi_ini_scf

$ \phi(t_0) $: scalar field initial value

double phi_prime_ini_scf

$ d\phi(t_0)/d\tau $: scalar field initial derivative wrt conformal time

int scf_parameters_size

size of scf_parameters

double varconst_alpha

finestructure constant for varying fundamental constants

double varconst_me

electron mass for varying fundamental constants

enum varconst_dependence varconst_dep

dependence of the varying fundamental constants as a function of time

double varconst_transition_redshift

redshift of transition between varied fundamental constants and normal fundamental constants in the 'varconst_instant' case

double age

age in Gyears

double conformal_age

conformal age in Mpc

double K

$ K $: Curvature parameter $ K=-\Omega0_k*a_{today}^2*H_0^2$;

int sgnK

K/|K|: -1, 0 or 1

double Neff

so-called "effective neutrino number", computed at earliest time in interpolation table

double Omega0_dcdm

$ \Omega_{0 dcdm} $: decaying cold dark matter

double Omega0_dr

$ \Omega_{0 dr} $: decay radiation

double Omega0_m

total non-relativistic matter today

double Omega0_r

total ultra-relativistic radiation today

double Omega0_de

total dark energy density today, currently defined as 1 - Omega0_m - Omega0_r - Omega0_k

double Omega0_nfsm

total non-free-streaming matter, that is, cdm, baryons and wdm

double a_eq

scale factor at radiation/matter equality

double H_eq

Hubble rate at radiation/matter equality [Mpc^-1]

double z_eq

redshift at radiation/matter equality

double tau_eq

conformal time at radiation/matter equality [Mpc]

int index_bg_a

scale factor (in fact (a/a_0), see normalisation conventions explained at beginning of background.c)

int index_bg_H

Hubble parameter in $Mpc^{-1}$

int index_bg_H_prime

its derivative w.r.t. conformal time

int index_bg_rho_g

photon density

int index_bg_rho_b

baryon density

int index_bg_rho_cdm

cdm density

int index_bg_rho_idm

idm density

int index_bg_rho_lambda

cosmological constant density

int index_bg_rho_fld

fluid density

int index_bg_w_fld

fluid equation of state

int index_bg_rho_idr

density of interacting dark radiation

int index_bg_rho_ur

relativistic neutrinos/relics density

int index_bg_rho_dcdm

dcdm density

int index_bg_rho_dr

dr density

int index_bg_phi_scf

scalar field value

int index_bg_phi_prime_scf

scalar field derivative wrt conformal time

int index_bg_V_scf

scalar field potential V

int index_bg_dV_scf

scalar field potential derivative V'

int index_bg_ddV_scf

scalar field potential second derivative V''

int index_bg_rho_scf

scalar field energy density

int index_bg_p_scf

scalar field pressure

int index_bg_p_prime_scf

scalar field pressure

int index_bg_rho_ncdm1

density of first ncdm species (others contiguous)

int index_bg_p_ncdm1

pressure of first ncdm species (others contiguous)

int index_bg_pseudo_p_ncdm1

another statistical momentum useful in ncdma approximation

int index_bg_rho_tot

Total density

int index_bg_p_tot

Total pressure

int index_bg_p_tot_prime

Conf. time derivative of total pressure

int index_bg_Omega_r

relativistic density fraction ( $ \Omega_{\gamma} + \Omega_{\nu r} $)

int index_bg_rho_crit

critical density

int index_bg_Omega_m

non-relativistic density fraction ( $ \Omega_b + \Omega_cdm + \Omega_{\nu nr} $)

int index_bg_conf_distance

conformal distance (from us) in Mpc

int index_bg_ang_distance

angular diameter distance in Mpc

int index_bg_lum_distance

luminosity distance in Mpc

int index_bg_time

proper (cosmological) time in Mpc

int index_bg_rs

comoving sound horizon in Mpc

int index_bg_D

scale independent growth factor D(a) for CDM perturbations

int index_bg_f

corresponding velocity growth factor [dlnD]/[dln a]

int index_bg_varc_alpha

value of fine structure constant in varying fundamental constants

int index_bg_varc_me

value of effective electron mass in varying fundamental constants

int bg_size_short

size of background vector in the "short format"

int bg_size_normal

size of background vector in the "normal format"

int bg_size

size of background vector in the "long format"

int bt_size

number of lines (i.e. time-steps) in the four following array

double * loga_table

vector loga_table[index_loga] with values of log(a) (in fact $ log(a/a0) $, logarithm of relative scale factor compared to today)

double * tau_table

vector tau_table[index_loga] with values of conformal time $ \tau $ (in fact $ a_0 c tau $, see normalisation conventions explained at beginning of background.c)

double * z_table

vector z_table[index_loga] with values of $ z $ (redshift)

double * background_table

table background_table[index_tau*pba->bg_size+pba->index_bg] with all other quantities (array of size bg_size*bt_size)

double * d2tau_dz2_table

vector d2tau_dz2_table[index_loga] with values of $ d^2 \tau / dz^2 $ (conformal time)

double * d2z_dtau2_table

vector d2z_dtau2_table[index_loga] with values of $ d^2 z / d\tau^2 $ (conformal time)

double * d2background_dloga2_table

table d2background_dtau2_table[index_loga*pba->bg_size+pba->index_bg] with values of $ d^2 b_i / d\log(a)^2 $

int index_bi_rho_dcdm

{B} dcdm density

int index_bi_rho_dr

{B} dr density

int index_bi_rho_fld

{B} fluid density

int index_bi_phi_scf

{B} scalar field value

int index_bi_phi_prime_scf

{B} scalar field derivative wrt conformal time

int index_bi_time

{C} proper (cosmological) time in Mpc

int index_bi_rs

{C} sound horizon

int index_bi_tau

{C} conformal time in Mpc

int index_bi_D

{C} scale independent growth factor D(a) for CDM perturbations.

int index_bi_D_prime

{C} D satisfies $ [D''(\tau)=-aHD'(\tau)+3/2 a^2 \rho_M D(\tau) $

int bi_B_size

Number of {B} parameters

int bi_size

Number of {B}+{C} parameters

short has_cdm

presence of cold dark matter?

short has_idm

presence of interacting dark matter with photons, baryons, and idr

short has_dcdm

presence of decaying cold dark matter?

short has_dr

presence of relativistic decay radiation?

short has_scf

presence of a scalar field?

short has_ncdm

presence of non-cold dark matter?

short has_lambda

presence of cosmological constant?

short has_fld

presence of fluid with constant w and cs2?

short has_ur

presence of ultra-relativistic neutrinos/relics?

short has_idr

presence of interacting dark radiation?

short has_curvature

presence of global spatial curvature?

short has_varconst

presence of varying fundamental constants?

int * ncdm_quadrature_strategy

Vector of integers according to quadrature strategy.

double ** q_ncdm_bg

Pointers to vectors of background sampling in q

double ** w_ncdm_bg

Pointers to vectors of corresponding quadrature weights w

double ** q_ncdm

Pointers to vectors of perturbation sampling in q

double ** w_ncdm

Pointers to vectors of corresponding quadrature weights w

double ** dlnf0_dlnq_ncdm

Pointers to vectors of logarithmic derivatives of p-s-d

int * q_size_ncdm_bg

Size of the q_ncdm_bg arrays

int * q_size_ncdm

Size of the q_ncdm arrays

double * factor_ncdm

List of normalization factors for calculating energy density etc.

short shooting_failed

flag is set to true if shooting failed.

ErrorMsg shooting_error

Error message from shooting failed.

short background_verbose

flag regulating the amount of information sent to standard output (none if set to zero)

ErrorMsg error_message

zone for writing error messages

◆ background_parameters_and_workspace

struct background_parameters_and_workspace

temporary parameters and workspace passed to the background_derivs function

◆ background_parameters_for_distributions

struct background_parameters_for_distributions

temporary parameters and workspace passed to phase space distribution function

Enumeration Type Documentation

◆ spatial_curvature

list of possible types of spatial curvature

◆ equation_of_state

list of possible parametrisations of the DE equation of state

◆ varconst_dependence

list of possible parametrizations of the varying fundamental constants

◆ vecback_format

list of formats for the vector of background quantities

◆ interpolation_method

list of interpolation methods: search location in table either by bisection (inter_normal), or step by step starting from given index (inter_closeby)