CLASS MANUAL
background.c File Reference
#include "background.h"
+ Include dependency graph for background.c:

Functions

int background_at_z (struct background *pba, double z, enum vecback_format return_format, enum interpolation_method inter_mode, int *last_index, double *pvecback)
 
int background_at_tau (struct background *pba, double tau, enum vecback_format return_format, enum interpolation_method inter_mode, int *last_index, double *pvecback)
 
int background_tau_of_z (struct background *pba, double z, double *tau)
 
int background_z_of_tau (struct background *pba, double tau, double *z)
 
int background_functions (struct background *pba, double a, double *pvecback_B, enum vecback_format return_format, double *pvecback)
 
int background_w_fld (struct background *pba, double a, double *w_fld, double *dw_over_da_fld, double *integral_fld)
 
int background_varconst_of_z (struct background *pba, double z, double *alpha, double *me)
 
int background_init (struct precision *ppr, struct background *pba)
 
int background_free (struct background *pba)
 
int background_free_noinput (struct background *pba)
 
int background_free_input (struct background *pba)
 
int background_indices (struct background *pba)
 
int background_ncdm_distribution (void *pbadist, double q, double *f0)
 
int background_ncdm_test_function (void *pbadist, double q, double *test)
 
int background_ncdm_init (struct precision *ppr, struct background *pba)
 
int background_ncdm_momenta (double *qvec, double *wvec, int qsize, double M, double factor, double z, double *n, double *rho, double *p, double *drho_dM, double *pseudo_p)
 
int background_ncdm_M_from_Omega (struct precision *ppr, struct background *pba, int n_ncdm)
 
int background_checks (struct precision *ppr, struct background *pba)
 
int background_solve (struct precision *ppr, struct background *pba)
 
int background_initial_conditions (struct precision *ppr, struct background *pba, double *pvecback, double *pvecback_integration, double *loga_ini)
 
int background_find_equality (struct precision *ppr, struct background *pba)
 
int background_output_titles (struct background *pba, char titles[_MAXTITLESTRINGLENGTH_])
 
int background_output_data (struct background *pba, int number_of_titles, double *data)
 
int background_derivs (double loga, double *y, double *dy, void *parameters_and_workspace, ErrorMsg error_message)
 
int background_sources (double loga, double *y, double *dy, int index_loga, void *parameters_and_workspace, ErrorMsg error_message)
 
int background_timescale (double loga, void *parameters_and_workspace, double *timescale, ErrorMsg error_message)
 
int background_output_budget (struct background *pba)
 
double V_e_scf (struct background *pba, double phi)
 
double V_p_scf (struct background *pba, double phi)
 
double V_scf (struct background *pba, double phi)
 

Detailed Description

Documented background module

  • Julien Lesgourgues, 17.04.2011
  • routines related to ncdm written by T. Tram in 2011
  • new integration scheme written by N. Schoeneberg in 2020

Deals with the cosmological background evolution. This module has two purposes:

  • at the beginning, to initialize the background, i.e. to integrate the background equations, and store all background quantities as a function of conformal time inside an interpolation table.
  • to provide routines which allow other modules to evaluate any background quantity for a given value of the conformal time (by interpolating within the interpolation table), or to find the correspondence between redshift and conformal time.

The overall logic in this module is the following:

  1. most background parameters that we will call {A} (e.g. rho_gamma, ..) can be expressed as simple analytical functions of the scale factor 'a' plus a few variables that we will call {B} (e.g. (phi, phidot) for quintessence, or some temperature for exotic particles, etc...). [Side note: for simplicity, all variables {B} are declared redundently inside {A}.]
  2. in turn, quantities {B} can be found as a function of the the scale factor [or rather (a/a_0)] by integrating the background equations. Thus {B} also includes the density of species which energy conservation equation must be integrated explicitely, like the density of fluids or of decaying dark matter.
  3. some other quantities that we will call {C} (like e.g. proper and conformal time, the sound horizon, the analytic scale-invariant growth factor) also require an explicit integration with respect to (a/a_0) [or rather log(a/a_p)], since they cannot be inferred analytically from (a/a_0) and parameters {B}. The difference between {B} and {C} parameters is that {C} parameters do not need to be known in order to get {A}.

So, we define the following routines:

  • background_functions() returns all background quantities {A} as a function of (a/a_0) and of quantities {B}.
  • background_solve() integrates the quantities {B} and {C} with respect to log(a/a_0); this integration requires many calls to background_functions().
  • the result is stored in the form of a big table in the background structure. There is one column for the scale factor, and one for each quantity {A} or {C} [Side note: we don;t include {B} here because the {B} variables are already decalred redundently also as {A} quantitites.]

Later in the code:

  • If we know the variables (a/a_0) + {B} and need some quantity {A} (but not {C}), the quickest and most precise way is to call directly background_functions() (for instance, in simple models, if we want H at a given value of the scale factor).
  • If we know 'tau' and want any other quantity, we can call background_at_tau(), which interpolates in the table and returns all values.
  • If we know 'z' but not the {B} variables, or if we know 'z' and we want {C} variables, we need to call background_at_z(), which interpolates in the table and returns all values.
  • Finally, it can be useful to get 'tau' for a given redshift 'z' or vice-versa: this can be done with background_tau_of_z() or background_z_of_tau().

In order to save time, background_at_tau() ans background_at_z() can be called in three modes: short_info, normal_info, long_info (returning only essential quantities, or useful quantities, or rarely useful quantities). Each line in the interpolation table is a vector whose first few elements correspond to the short_info format; a larger fraction contribute to the normal format; and the full vector corresponds to the long format. The guideline is that short_info returns only geometric quantities like a, H, H'; normal format returns quantities strictly needed at each step in the integration of perturbations; long_info returns quantities needed only occasionally.

In summary, the following functions can be called from other modules:

  1. background_init() at the beginning background_at_tau(),
  2. background_at_z(), background_tau_of_z(), background_z_of_tau() at any later time
  3. background_free() at the end, when no more calls to the previous functions are needed

For units and normalisation conventions, there are two guiding principles:

1) All quantities are expressed in natural units in which everything is in powers of Mpc, e.g.:

  • t stands for (cosmological or proper time)*c in Mpc
  • tau stands for (conformal time)*c in Mpc
  • H stands for (Hubble parameter)/c in $ Mpc^{-1} $
  • etc.

2) New since v3.0: all quantities that should normally scale with some power of a_0^n are renormalised by a_0^{-n}, in order to be independent of a_0, e.g.

  • a in the code stands for $ a/a_0 $ in reality
  • tau in the code stands for $ a_0 \tau c $ in Mpc
  • any prime in the code stands for $ (1/a_0) d/d\tau $
  • r stands for any comoving radius times a_0
  • etc.

Function Documentation

◆ background_at_z()

int background_at_z ( struct background pba,
double  z,
enum vecback_format  return_format,
enum interpolation_method  inter_mode,
int *  last_index,
double *  pvecback 
)

Background quantities at given redshift z.

Evaluates all background quantities at a given value of redshift by reading the pre-computed table and interpolating.

Parameters
pbaInput: pointer to background structure (containing pre-computed table)
zInput: redshift
return_formatInput: format of output vector (short_info, normal_info, long_info)
inter_modeInput: interpolation mode (normal or closeby)
last_indexInput/Output: index of the previous/current point in the interpolation array (input only for closeby mode, output for both)
pvecbackOutput: vector (assumed to be already allocated)
Returns
the error status

Summary:

  • define local variables
  • check that log(a) = log(1/(1+z)) = -log(1+z) is in the pre-computed range
  • deduce length of returned vector from format mode
  • interpolate from pre-computed table with array_interpolate() or array_interpolate_growing_closeby() (depending on interpolation mode)

◆ background_at_tau()

int background_at_tau ( struct background pba,
double  tau,
enum vecback_format  return_format,
enum interpolation_method  inter_mode,
int *  last_index,
double *  pvecback 
)

Background quantities at given conformal time tau.

Evaluates all background quantities at a given value of conformal time by reading the pre-computed table and interpolating.

Parameters
pbaInput: pointer to background structure (containing pre-computed table)
tauInput: value of conformal time
return_formatInput: format of output vector (short_info, normal_info, long_info)
inter_modeInput: interpolation mode (normal or closeby)
last_indexInput/Output: index of the previous/current point in the interpolation array (input only for closeby mode, output for both)
pvecbackOutput: vector (assumed to be already allocated)
Returns
the error status

Summary:

  • define local variables
  • Get current redshift
  • Get background at corresponding redshift

◆ background_tau_of_z()

int background_tau_of_z ( struct background pba,
double  z,
double *  tau 
)

Conformal time at given redshift.

Returns tau(z) by interpolation from pre-computed table.

Parameters
pbaInput: pointer to background structure
zInput: redshift
tauOutput: conformal time
Returns
the error status

Summary:

  • define local variables
  • check that $ z $ is in the pre-computed range
  • interpolate from pre-computed table with array_interpolate()

◆ background_z_of_tau()

int background_z_of_tau ( struct background pba,
double  tau,
double *  z 
)

Redshift at given conformal time.

Returns z(tau) by interpolation from pre-computed table.

Parameters
pbaInput: pointer to background structure
tauInput: conformal time
zOutput: redshift
Returns
the error status

Summary:

  • define local variables
  • check that $ tau $ is in the pre-computed range
  • interpolate from pre-computed table with array_interpolate()

◆ background_functions()

int background_functions ( struct background pba,
double  a,
double *  pvecback_B,
enum vecback_format  return_format,
double *  pvecback 
)

Function evaluating all background quantities which can be computed analytically as a function of a and of {B} quantities (see discussion at the beginning of this file).

Parameters
pbaInput: pointer to background structure
aInput: scale factor (in fact, with our normalisation conventions, this is (a/a_0) )
pvecback_BInput: vector containing all {B} quantities
return_formatInput: format of output vector
pvecbackOutput: vector of background quantities (assumed to be already allocated)
Returns
the error status

Summary:

  • define local variables
  • initialize local variables
  • pass value of $ a$ to output
  • compute each component's density and pressure

<– This depends on a_prime_over_a, so we cannot add it now!

See e.g. Eq. A6 in 1811.00904.

  • compute expansion rate H from Friedmann equation: this is the only place where the Friedmann equation is assumed. Remember that densities are all expressed in units of $ [3c^2/8\pi G] $, ie $ \rho_{class} = [8 \pi G \rho_{physical} / 3 c^2]$
  • compute derivative of H with respect to conformal time

The contribution of scf was not added to dp_dloga, add p_scf_prime here:

  • compute critical density
  • compute relativistic density to total density ratio
  • compute other quantities in the exhaustive, redundant format
  • store critical density
  • compute Omega_m
  • cosmological time
  • comoving sound horizon
  • growth factor
  • velocity growth factor
  • Varying fundamental constants

◆ background_w_fld()

int background_w_fld ( struct background pba,
double  a,
double *  w_fld,
double *  dw_over_da_fld,
double *  integral_fld 
)

Single place where the fluid equation of state is defined. Parameters of the function are passed through the background structure. Generalisation to arbitrary functions should be simple.

Parameters
pbaInput: pointer to background structure
aInput: current value of scale factor (in fact, with our conventions, of (a/a_0))
w_fldOutput: equation of state parameter w_fld(a)
dw_over_da_fldOutput: function dw_fld/da
integral_fldOutput: function $ \int_{a}^{a_0} da 3(1+w_{fld})/a $
Returns
the error status
  • first, define the function w(a)
  • then, give the corresponding analytic derivative dw/da (used by perturbation equations; we could compute it numerically, but with a loss of precision; as long as there is a simple analytic expression of the derivative of the previous function, let's use it!
  • finally, give the analytic solution of the following integral: $ \int_{a}^{a0} da 3(1+w_{fld})/a $. This is used in only one place, in the initial conditions for the background, and with a=a_ini. If your w(a) does not lead to a simple analytic solution of this integral, no worry: instead of writing something here, the best would then be to leave it equal to zero, and then in background_initial_conditions() you should implement a numerical calculation of this integral only for a=a_ini, using for instance Romberg integration. It should be fast, simple, and accurate enough.

note: of course you can generalise these formulas to anything, defining new parameters pba->w..._fld. Just remember that so far, HyRec explicitely assumes that w(a)= w0 + wa (1-a/a0); but Recfast does not assume anything

◆ background_varconst_of_z()

int background_varconst_of_z ( struct background pba,
double  z,
double *  alpha,
double *  me 
)

Single place where the variation of fundamental constants is defined. Parameters of the function are passed through the background structure. Generalisation to arbitrary functions should be simple.

Parameters
pbaInput: pointer to background structure
zInput: current value of redhsift
alphaOutput: fine structure constant relative to its current value
meOutput: effective electron mass relative to its current value
Returns
the error status

◆ background_init()

int background_init ( struct precision ppr,
struct background pba 
)

Initialize the background structure, and in particular the background interpolation table.

Parameters
pprInput: pointer to precision structure
pbaInput/Output: pointer to initialized background structure
Returns
the error status

Summary:

  • write class version
  • if shooting failed during input, catch the error here
  • assign values to all indices in vectors of background quantities
  • check that input parameters make sense and write additional information about them
  • integrate the background over log(a), allocate and fill the background table
  • find and store a few derived parameters at radiation-matter equality

◆ background_free()

int background_free ( struct background pba)

Free all memory space allocated by background_init() and by input_read_parameters().

Parameters
pbaInput: pointer to background structure (to be freed)
Returns
the error status

◆ background_free_noinput()

int background_free_noinput ( struct background pba)

Free only the memory space NOT allocated through input_read_parameters(), but through background_init()

Parameters
pbaInput: pointer to background structure (to be freed)
Returns
the error status

◆ background_free_input()

int background_free_input ( struct background pba)

Free pointers inside background structure which were allocated in input_read_parameters()

Parameters
pbaInput: pointer to background structure
Returns
the error status

◆ background_indices()

int background_indices ( struct background pba)

Assign value to each relevant index in vectors of background quantities.

Parameters
pbaInput: pointer to background structure
Returns
the error status

Summary:

  • define local variables
  • initialize all flags: which species are present?
  • initialize all indices

◆ background_ncdm_distribution()

int background_ncdm_distribution ( void *  pbadist,
double  q,
double *  f0 
)

This is the routine where the distribution function f0(q) of each ncdm species is specified (it is the only place to modify if you need a partlar f0(q))

Parameters
pbadistInput: structure containing all parameters defining f0(q)
qInput: momentum
f0Output: phase-space distribution
  • extract from the input structure pbadist all the relevant information
  • shall we interpolate in file, or shall we use analytical formula below?
  • a) deal first with the case of interpolating in files
  • b) deal now with case of reading analytical function

Next enter your analytic expression(s) for the p.s.d.'s. If you need different p.s.d.'s for different species, put each p.s.d inside a condition, like for instance: if (n_ncdm==2) {*f0=...}. Remember that n_ncdm = 0 refers to the first species.

This form is only appropriate for approximate studies, since in reality the chemical potentials are associated with flavor eigenstates, not mass eigenstates. It is easy to take this into account by introducing the mixing angles. In the later part (not read by the code) we illustrate how to do this.

◆ background_ncdm_test_function()

int background_ncdm_test_function ( void *  pbadist,
double  q,
double *  test 
)

This function is only used for the purpose of finding optimal quadrature weights. The logic is: if we can accurately convolve f0(q) with this function, then we can convolve it accurately with any other relevant function.

Parameters
pbadistInput: structure containing all background parameters
qInput: momentum
testOutput: value of the test function test(q)

Using a + bq creates problems for otherwise acceptable distributions which diverges as $ 1/r $ or $ 1/r^2 $ for $ r\to 0 $

◆ background_ncdm_init()

int background_ncdm_init ( struct precision ppr,
struct background pba 
)

This function finds optimal quadrature weights for each ncdm species

Parameters
pprInput: precision structure
pbaInput/Output: background structure

Automatic q-sampling for this species

  • in verbose mode, inform user of number of sampled momenta for background quantities

Manual q-sampling for this species. Same sampling used for both perturbation and background sampling, since this will usually be a high precision setting anyway

  • in verbose mode, inform user of number of sampled momenta for background quantities

◆ background_ncdm_momenta()

int background_ncdm_momenta ( double *  qvec,
double *  wvec,
int  qsize,
double  M,
double  factor,
double  z,
double *  n,
double *  rho,
double *  p,
double *  drho_dM,
double *  pseudo_p 
)

For a given ncdm species: given the quadrature weights, the mass and the redshift, find background quantities by a quick weighted sum over. Input parameters passed as NULL pointers are not evaluated for speed-up

Parameters
qvecInput: sampled momenta
wvecInput: quadrature weights
qsizeInput: number of momenta/weights
MInput: mass
factorInput: normalization factor for the p.s.d.
zInput: redshift
nOutput: number density
rhoOutput: energy density
pOutput: pressure
drho_dMOutput: derivative used in next function
pseudo_pOutput: pseudo-pressure used in perturbation module for fluid approx

Summary:

  • rescale normalization at given redshift
  • initialize quantities
  • loop over momenta
  • adjust normalization
+ Here is the caller graph for this function:

◆ background_ncdm_M_from_Omega()

int background_ncdm_M_from_Omega ( struct precision ppr,
struct background pba,
int  n_ncdm 
)

When the user passed the density fraction Omega_ncdm or omega_ncdm in input but not the mass, infer the mass with Newton iteration method.

Parameters
pprInput: precision structure
pbaInput/Output: background structure
n_ncdmInput: index of ncdm species
+ Here is the call graph for this function:

◆ background_checks()

int background_checks ( struct precision ppr,
struct background pba 
)

Perform some check on the input background quantities, and send to standard output some information about them

Parameters
pprInput: pointer to precision structure
pbaInput: pointer to initialized background structure
Returns
the error status
  • define local variables
  • control that we have photons and baryons in the problem
  • control that cosmological parameter values make sense, otherwise inform user
  • in verbose mode, send to standard output some additional information on non-obvious background parameters

◆ background_solve()

int background_solve ( struct precision ppr,
struct background pba 
)

This function integrates the background over time, allocates and fills the background table

Parameters
pprInput: precision structure
pbaInput/Output: background structure

Summary:

  • define local variables
  • setup background workspace
  • allocate vector of quantities to be integrated
  • impose initial conditions with background_initial_conditions()
  • Determine output vector
  • allocate background tables
  • define values of loga at which results will be stored
  • choose the right evolver
  • perform the integration
  • recover some quantities today
  • In a loop over lines, fill rest of background table for quantities that depend on numbers like "conformal_age" or "D_today" that were calculated just before
  • fill tables of second derivatives (in view of spline interpolation)
  • compute remaining "related parameters"
  • so-called "effective neutrino number", computed at earliest time in interpolation table. This should be seen as a definition: Neff is the equivalent number of instantaneously-decoupled neutrinos accounting for the radiation density, beyond photons
  • send information to standard output
  • store information in the background structure

◆ background_initial_conditions()

int background_initial_conditions ( struct precision ppr,
struct background pba,
double *  pvecback,
double *  pvecback_integration,
double *  loga_ini 
)

Assign initial values to background integrated variables.

Parameters
pprInput: pointer to precision structure
pbaInput: pointer to background structure
pvecbackInput: vector of background quantities used as workspace
pvecback_integrationOutput: vector of background quantities to be integrated, returned with proper initial values
loga_iniOutput: value of loga (in fact with our conventions log(a/a_0)) at initial time
Returns
the error status

Summary:

  • define local variables
  • fix initial value of $ a $

If we have ncdm species, perhaps we need to start earlier than the standard value for the species to be relativistic. This could happen for some WDM models.

  • We must add the relativistic contribution from NCDM species
  • f is the critical density fraction of DR. The exact solution is:

f = -Omega_rad+pow(pow(Omega_rad,3./2.)+0.5*pow(a,6)*pvecback_integration[pba->index_bi_rho_dcdm]*pba->Gamma_dcdm/pow(pba->H0,3),2./3.);

but it is not numerically stable for very small f which is always the case. Instead we use the Taylor expansion of this equation, which is equivalent to ignoring f(a) in the Hubble rate.

There is also a space reserved for a future case where dr is not sourced by dcdm

  • Fix initial value of $ \phi, \phi' $ set directly in the radiation attractor => fixes the units in terms of rho_ur

TODO:

  • There seems to be some small oscillation when it starts.
  • Check equations and signs. Sign of phi_prime?
  • is rho_ur all there is early on?
  • --> If there is no attractor solution for scf_lambda, assign some value. Otherwise would give a nan.
  • --> If no attractor initial conditions are assigned, gets the provided ones.
  • compute initial proper time, assuming radiation-dominated universe since Big Bang and therefore $ t=1/(2H) $ (good approximation for most purposes)
  • compute initial conformal time, assuming radiation-dominated universe since Big Bang and therefore $ \tau=1/(aH) $ (good approximation for most purposes)
  • compute initial sound horizon, assuming $ c_s=1/\sqrt{3} $ initially
  • set initial value of D and D' in RD. D and D' need only be set up to an overall constant, since they will later be re-normalized. From Ma&Bertschinger, one can derive D ~ (ktau)^2 at early times, from which one finds D'/D = 2 aH (assuming aH=1/tau during RD)
  • return the value finally chosen for the initial log(a)

◆ background_find_equality()

int background_find_equality ( struct precision ppr,
struct background pba 
)

Find the time of radiation/matter equality and store characteristic quantitites at that time in the background structure..

Parameters
pprInput: pointer to precision structure
pbaInput/Output: pointer to background structure
Returns
the error status

◆ background_output_titles()

int background_output_titles ( struct background pba,
char  titles[_MAXTITLESTRINGLENGTH_] 
)

Subroutine for formatting background output

Parameters
pbaInput: pointer to background structure
titlesOuput: name of columns when printing the background table
Returns
the error status
  • Length of the column title should be less than OUTPUTPRECISION+6 to be indented correctly, but it can be as long as .

◆ background_output_data()

int background_output_data ( struct background pba,
int  number_of_titles,
double *  data 
)

Subroutine for writing the background output

Parameters
pbaInput: pointer to background structure
number_of_titlesInput: number of background quantities to print at each time step
dataOuput: 1d array storing all the background table
Returns
the error status

Stores quantities

◆ background_derivs()

int background_derivs ( double  loga,
double *  y,
double *  dy,
void *  parameters_and_workspace,
ErrorMsg  error_message 
)

Subroutine evaluating the derivative with respect to loga of quantities which are integrated (tau, t, etc).

This is one of the few functions in the code which is passed to the generic_integrator() routine. Since generic_integrator() should work with functions passed from various modules, the format of the arguments is a bit special:

  • fixed input parameters and workspaces are passed through a generic pointer. Here, this is just a pointer to the background structure and to a background vector, but generic_integrator() doesn't know its fine structure.
  • the error management is a bit special: errors are not written as usual to pba->error_message, but to a generic error_message passed in the list of arguments.
Parameters
logaInput: current value of log(a)
yInput: vector of variable
dyOutput: its derivative (already allocated)
parameters_and_workspaceInput: pointer to fixed parameters (e.g. indices)
error_messageOutput: error message

Summary:

  • define local variables
  • scale factor a (in fact, given our normalisation conventions, this stands for a/a_0)
  • calculate functions of $ a $ with background_functions()
  • Short hand notation for Hubble
  • calculate derivative of cosmological time $ dt/dloga = 1/H $
  • calculate derivative of conformal time $ d\tau/dloga = 1/aH $
  • calculate detivative of sound horizon $ drs/dloga = drs/dtau * dtau/dloga = c_s/aH $
  • solve second order growth equation $ [D''(\tau)=-aHD'(\tau)+3/2 a^2 \rho_M D(\tau) $ written as $ dD/dloga = D' / (aH) $ and $ dD'/dloga = -D' + (3/2) (a/H) \rho_M D $
  • compute dcdm density $ d\rho/dloga = -3 \rho - \Gamma/H \rho $
  • Compute dr density $ d\rho/dloga = -4\rho - \Gamma/H \rho $
  • Compute fld density $ d\rho/dloga = -3 (1+w_{fld}(a)) \rho $
  • Scalar field equation: $ \phi'' + 2 a H \phi' + a^2 dV = 0 $ (note H is wrt cosmological time) written as $ d\phi/dlna = phi' / (aH) $ and $ d\phi'/dlna = -2*phi' - (a/H) dV $

◆ background_sources()

int background_sources ( double  loga,
double *  y,
double *  dy,
int  index_loga,
void *  parameters_and_workspace,
ErrorMsg  error_message 
)

At some step during the integraton of the background equations, this function extracts the qantities that we want to keep memory of, and stores them in a row of the background table (as well as extra tables: z_table, tau_table).

This is one of the few functions in the code which is passed to the generic_integrator() routine. Since generic_integrator() should work with functions passed from various modules, the format of the arguments is a bit special:

  • fixed parameters and workspaces are passed through a generic pointer. generic_integrator() doesn't know the content of this pointer.
  • the error management is a bit special: errors are not written as usual to pba->error_message, but to a generic error_message passed in the list of arguments.
Parameters
logaInput: current value of log(a)
yInput: current vector of integrated quantities (with index_bi)
dyInput: current derivative of y w.r.t log(a)
index_logaInput: index of the log(a) value within the background_table
parameters_and_workspaceInput/output: fixed parameters (e.g. indices), workspace, background structure where the output is written...
error_messageOutput: error message
  • localize the row inside background_table where the current values must be stored
  • scale factor a (in fact, given our normalisation conventions, this stands for a/a_0)
  • corresponding redhsift 1/a-1
  • corresponding conformal time

-> compute all other quantities depending only on a + {B} variables and get them stored in one row of background_table The value of {B} variables in pData are also copied to pvecback.

◆ background_timescale()

int background_timescale ( double  loga,
void *  parameters_and_workspace,
double *  timescale,
ErrorMsg  error_message 
)

Evalute the typical timescale for the integration of he background over loga=log(a/a_0). This is only required for rkck, but not for the ndf15 evolver.

The evolver will take steps equal to this value times ppr->background_integration_stepsize. Since our variable of integration is loga, and the time steps are (delta a)/a, the reference timescale is precisely one, i.e., the code will take some steps such that (delta a)/a = ppr->background_integration_stepsize.

The argument list is predetermined by the format of generic_evolver; however in this particular case, they are never used.

This is one of the few functions in the code which is passed to the generic_integrator() routine. Since generic_integrator() should work with functions passed from various modules, the format of the arguments is a bit special:

  • fixed parameters and workspaces are passed through a generic pointer (void *). generic_integrator() doesn't know the content of this pointer.
  • the error management is a bit special: errors are not written as usual to pba->error_message, but to a generic error_message passed in the list of arguments.
Parameters
logaInput: current value of log(a/a_0)
parameters_and_workspaceInput: fixed parameters (e.g. indices), workspace, approximation used, etc.
timescaleOutput: perturbation variation timescale
error_messageOutput: error message

◆ background_output_budget()

int background_output_budget ( struct background pba)

Function outputting the fractions Omega of the total critical density today, and also the reduced fractions omega=Omega*h*h

It also prints the total budgets of non-relativistic, relativistic, and other contents, and of the total

Parameters
pbaInput: Pointer to background structure
Returns
the error status

◆ V_e_scf()

double V_e_scf ( struct background pba,
double  phi 
)

Scalar field potential and its derivatives with respect to the field _scf For Albrecht & Skordis model: 9908085

  • $ V = V_{p_{scf}}*V_{e_{scf}} $
  • $ V_e = \exp(-\lambda \phi) $ (exponential)
  • $ V_p = (\phi - B)^\alpha + A $ (polynomial bump)

TODO:

  • Add some functionality to include different models/potentials (tuning would be difficult, though)
  • Generalize to Kessence/Horndeski/PPF and/or couplings
  • A default module to numerically compute the derivatives when no analytic functions are given should be added.
  • Numerical derivatives may further serve as a consistency check.

    The units of phi, tau in the derivatives and the potential V are the following:

    • phi is given in units of the reduced Planck mass $ m_{pl} = (8 \pi G)^{(-1/2)}$
    • tau in the derivative is given in units of Mpc.
    • the potential $ V(\phi) $ is given in units of $ m_{pl}^2/Mpc^2 $. With this convention, we have $ \rho^{class} = (8 \pi G)/3 \rho^{physical} = 1/(3 m_{pl}^2) \rho^{physical} = 1/3 * [ 1/(2a^2) (\phi')^2 + V(\phi) ] $ and $ \rho^{class} $ has the proper dimension $ Mpc^-2 $.
+ Here is the caller graph for this function:

◆ V_p_scf()

double V_p_scf ( struct background pba,
double  phi 
)

parameters and functions for the polynomial coefficient $ V_p = (\phi - B)^\alpha + A $(polynomial bump)

double scf_alpha = 2;

double scf_B = 34.8;

double scf_A = 0.01; (values for their Figure 2)

+ Here is the caller graph for this function:

◆ V_scf()

double V_scf ( struct background pba,
double  phi 
)

Fianlly we can obtain the overall potential $ V = V_p*V_e $

+ Here is the call graph for this function: