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

Functions

int primordial_spectrum_at_k (struct primordial *ppm, int index_md, enum linear_or_logarithmic mode, double input, double *output)
 
int primordial_init (struct precision *ppr, struct perturbations *ppt, struct primordial *ppm)
 
int primordial_free (struct primordial *ppm)
 
int primordial_indices (struct perturbations *ppt, struct primordial *ppm)
 
int primordial_get_lnk_list (struct primordial *ppm, double kmin, double kmax, double k_per_decade)
 
int primordial_analytic_spectrum_init (struct perturbations *ppt, struct primordial *ppm)
 
int primordial_analytic_spectrum (struct primordial *ppm, int index_md, int index_ic1_ic2, double k, double *pk)
 
int primordial_inflation_potential (struct primordial *ppm, double phi, double *V, double *dV, double *ddV)
 
int primordial_inflation_hubble (struct primordial *ppm, double phi, double *H, double *dH, double *ddH, double *dddH)
 
int primordial_inflation_indices (struct primordial *ppm)
 
int primordial_inflation_solve_inflation (struct perturbations *ppt, struct primordial *ppm, struct precision *ppr)
 
int primordial_inflation_analytic_spectra (struct perturbations *ppt, struct primordial *ppm, struct precision *ppr, double *y_ini)
 
int primordial_inflation_spectra (struct perturbations *ppt, struct primordial *ppm, struct precision *ppr, double *y_ini)
 
int primordial_inflation_one_wavenumber (struct perturbations *ppt, struct primordial *ppm, struct precision *ppr, double *y_ini, int index_k)
 
int primordial_inflation_one_k (struct primordial *ppm, struct precision *ppr, double k, double *y, double *dy, double *curvature, double *tensor)
 
int primordial_inflation_find_attractor (struct primordial *ppm, struct precision *ppr, double phi_0, double precision, double *y, double *dy, double *H_0, double *dphidt_0)
 
int primordial_inflation_evolve_background (struct primordial *ppm, struct precision *ppr, double *y, double *dy, enum target_quantity target, double stop, short check_epsilon, enum integration_direction direction, enum time_definition time)
 
int primordial_inflation_check_potential (struct primordial *ppm, double phi, double *V, double *dV, double *ddV)
 
int primordial_inflation_check_hubble (struct primordial *ppm, double phi, double *H, double *dH, double *ddH, double *dddH)
 
int primordial_inflation_get_epsilon (struct primordial *ppm, double phi, double *epsilon)
 
int primordial_inflation_find_phi_pivot (struct primordial *ppm, struct precision *ppr, double *y, double *dy)
 
int primordial_inflation_derivs (double tau, double *y, double *dy, void *parameters_and_workspace, ErrorMsg error_message)
 
int primordial_external_spectrum_init (struct perturbations *ppt, struct primordial *ppm)
 

Detailed Description

Documented primordial module.

Julien Lesgourgues, 24.08.2010

This module computes the primordial spectra. It can be used in different modes: simple parametric form, evolving inflaton perturbations, etc. So far only the mode corresponding to a simple analytic form in terms of amplitudes, tilts and runnings has been developed.

The following functions can be called from other modules:

  1. primordial_init() at the beginning (anytime after perturbations_init() and before harmonic_init())
  2. primordial_spectrum_at_k() at any time for computing P(k) at any k
  3. primordial_free() at the end

Function Documentation

◆ primordial_spectrum_at_k()

int primordial_spectrum_at_k ( struct primordial ppm,
int  index_md,
enum linear_or_logarithmic  mode,
double  input,
double *  output 
)

Primordial spectra for arbitrary argument and for all initial conditions.

This routine evaluates the primordial spectrum at a given value of k by interpolating in the pre-computed table.

When k is not in the pre-computed range but the spectrum can be found analytically, it finds it. Otherwise returns an error.

Can be called in two modes; linear or logarithmic:

  • linear: takes k, returns P(k)
  • logarithmic: takes ln(k), return ln(P(k))

One little subtlety: in case of several correlated initial conditions, the cross-correlation spectrum can be negative. Then, in logarithmic mode, the non-diagonal elements contain the cross-correlation angle $ P_{12}/\sqrt{P_{11} P_{22}}$ (from -1 to 1) instead of $\ln{P_{12}}$

This function can be called from whatever module at whatever time, provided that primordial_init() has been called before, and primordial_free() has not been called yet.

Parameters
ppmInput: pointer to primordial structure containing tabulated primordial spectrum
index_mdInput: index of mode (scalar, tensor, ...)
modeInput: linear or logarithmic
inputInput: wavenumber in 1/Mpc (linear mode) or its logarithm (logarithmic mode)
outputOutput: for each pair of initial conditions, primordial spectra P(k) in $Mpc^3$ (linear mode), or their logarithms and cross-correlation angles (logarithmic mode)
Returns
the error status

Summary:

  • define local variables
  • infer ln(k) from input. In linear mode, reject negative value of input k value.
  • if ln(k) is not in the interpolation range, return an error, unless we are in the case of a analytic spectrum, for which a direct computation is possible
  • otherwise, interpolate in the pre-computed table

◆ primordial_init()

int primordial_init ( struct precision ppr,
struct perturbations ppt,
struct primordial ppm 
)

This routine initializes the primordial structure (in particular, it computes table of primordial spectrum values)

Parameters
pprInput: pointer to precision structure (defines method and precision for all computations)
pptInput: pointer to perturbation structure (useful for knowing k_min, k_max, etc.)
ppmOutput: pointer to initialized primordial structure
Returns
the error status

Summary:

  • define local variables
  • check that we really need to compute the primordial spectra
  • get kmin and kmax from perturbation structure. Test that they make sense.
  • allocate and fill values of $ \ln{k}$'s
  • define indices and allocate tables in primordial structure
  • deal with case of analytic primordial spectra (with amplitudes, tilts, runnings, etc.)
  • deal with case of inflation with given $V(\phi)$ or $H(\phi)$
  • deal with the case of external calculation of $ P_k $
  • compute second derivative of each $ \ln{P_k} $ versus lnk with spline, in view of interpolation
  • derive spectral parameters from numerically computed spectra (not used by the rest of the code, but useful to keep in memory for several types of investigation)
  • expression for alpha_s comes from:

ns_2 = (lnpk_plus-lnpk_pivot)/(dlnk)+1

ns_1 = (lnpk_pivot-lnpk_minus)/(dlnk)+1

alpha_s = dns/dlnk = (ns_2-ns_1)/dlnk = (lnpk_plus-lnpk_pivot-lnpk_pivot+lnpk_minus)/(dlnk)/(dlnk)

  • expression for beta_s:

ppm->beta_s = (alpha_plus-alpha_minus)/dlnk = (lnpk_plusplus-2.*lnpk_plus+lnpk_pivot - (lnpk_pivot-2.*lnpk_minus+lnpk_minusminus)/pow(dlnk,3)

◆ primordial_free()

int primordial_free ( struct primordial ppm)

This routine frees all the memory space allocated by primordial_init().

To be called at the end of each run.

Parameters
ppmInput: pointer to primordial structure (which fields must be freed)
Returns
the error status

◆ primordial_indices()

int primordial_indices ( struct perturbations ppt,
struct primordial ppm 
)

This routine defines indices and allocates tables in the primordial structure

Parameters
pptInput: pointer to perturbation structure
ppmInput/output: pointer to primordial structure
Returns
the error status

◆ primordial_get_lnk_list()

int primordial_get_lnk_list ( struct primordial ppm,
double  kmin,
double  kmax,
double  k_per_decade 
)

This routine allocates and fills the list of wavenumbers k

Parameters
ppmInput/output: pointer to primordial structure
kminInput: first value
kmaxInput: last value that we should encompass
k_per_decadeInput: number of k per decade
Returns
the error status

◆ primordial_analytic_spectrum_init()

int primordial_analytic_spectrum_init ( struct perturbations ppt,
struct primordial ppm 
)

This routine interprets and stores in a condensed form the input parameters in the case of a simple analytic spectra with amplitudes, tilts, runnings, in such way that later on, the spectrum can be obtained by a quick call to the routine primordial_analytic_spectrum(()

Parameters
pptInput: pointer to perturbation structure
ppmInput/output: pointer to primordial structure
Returns
the error status

◆ primordial_analytic_spectrum()

int primordial_analytic_spectrum ( struct primordial ppm,
int  index_md,
int  index_ic1_ic2,
double  k,
double *  pk 
)

This routine returns the primordial spectrum in the simple analytic case with amplitudes, tilts, runnings, for each mode (scalar/tensor...), pair of initial conditions, and wavenumber.

Parameters
ppmInput/output: pointer to primordial structure
index_mdInput: index of mode (scalar, tensor, ...)
index_ic1_ic2Input: pair of initial conditions (ic1, ic2)
kInput: wavenumber in same units as pivot scale, i.e. in 1/Mpc
pkOutput: primordial power spectrum A (k/k_pivot)^(n+...)
Returns
the error status

◆ primordial_inflation_potential()

int primordial_inflation_potential ( struct primordial ppm,
double  phi,
double *  V,
double *  dV,
double *  ddV 
)

This routine encodes the inflaton scalar potential

Parameters
ppmInput: pointer to primordial structure
phiInput: background inflaton field value in units of Mp
VOutput: inflaton potential in units of $ Mp^4$
dVOutput: first derivative of inflaton potential wrt the field
ddVOutput: second derivative of inflaton potential wrt the field
Returns
the error status

◆ primordial_inflation_hubble()

int primordial_inflation_hubble ( struct primordial ppm,
double  phi,
double *  H,
double *  dH,
double *  ddH,
double *  dddH 
)

This routine encodes the function $ H(\phi)$

Parameters
ppmInput: pointer to primordial structure
phiInput: background inflaton field value in units of Mp
HOutput: Hubble parameters in units of Mp
dHOutput: $ dH / d\phi $
ddHOutput: $ d^2H / d\phi^2 $
dddHOutput: $ d^3H / d\phi^3 $
Returns
the error status

◆ primordial_inflation_indices()

int primordial_inflation_indices ( struct primordial ppm)

This routine defines indices used by the inflation simulator

Parameters
ppmInput/output: pointer to primordial structure
Returns
the error status

◆ primordial_inflation_solve_inflation()

int primordial_inflation_solve_inflation ( struct perturbations ppt,
struct primordial ppm,
struct precision ppr 
)

Main routine of inflation simulator. Its goal is to check the background evolution before and after the pivot value phi=phi_pivot, and then, if this evolution is suitable, to call the routine primordial_inflation_spectra().

Parameters
pptInput: pointer to perturbation structure
ppmInput/output: pointer to primordial structure
pprInput: pointer to precision structure
Returns
the error status

Summary:

  • define local variables
  • allocate vectors for background/perturbed quantities
  • eventually, needs first to find phi_pivot
  • compute H_pivot at phi_pivot
  • check positivity and negative slope of potential in field pivot value, and find value of phi_dot and H for field's pivot value, assuming slow-roll attractor solution has been reached. If no solution, code will stop there.
  • check positivity and negative slope of $ H(\phi)$ in field pivot value, and get H_pivot
  • find a_pivot, value of scale factor when k_pivot crosses horizon while phi=phi_pivot
  • integrate background solution starting from phi_pivot and until k_max>>aH. This ensures that the inflationary model considered here is valid and that the primordial spectrum can be computed. Otherwise, if slow-roll brakes too early, model is not suitable and run stops.
  • starting from this time, i.e. from y_ini[ ], we run the routine which takes care of computing the primordial spectrum.
  • before ending, we want to compute and store the values of $ \phi $ corresponding to k=aH for k_min and k_max
  • finally, we can de-allocate

◆ primordial_inflation_analytic_spectra()

int primordial_inflation_analytic_spectra ( struct perturbations ppt,
struct primordial ppm,
struct precision ppr,
double *  y_ini 
)

Routine for the computation of an analytic apporoximation to the the primordial spectrum. In general, should be used only for comparing with exact numerical computation performed by primordial_inflation_spectra().

Parameters
pptInput: pointer to perturbation structure
ppmInput/output: pointer to primordial structure
pprInput: pointer to precision structure
y_iniInput: initial conditions for the vector of background/perturbations, already allocated and filled
Returns
the error status

Summary

  • allocate vectors for background/perturbed quantities
  • initialize the background part of the running vector
  • loop over Fourier wavenumbers
  • read value of phi at time when k=aH
  • get potential (and its derivatives) at this value
  • calculate the analytic slow-roll formula for the spectra
  • store the obtained result for curvature and tensor perturbations

◆ primordial_inflation_spectra()

int primordial_inflation_spectra ( struct perturbations ppt,
struct primordial ppm,
struct precision ppr,
double *  y_ini 
)

Routine with a loop over wavenumbers for the computation of the primordial spectrum. For each wavenumber it calls primordial_inflation_one_wavenumber()

Parameters
pptInput: pointer to perturbation structure
ppmInput/output: pointer to primordial structure
pprInput: pointer to precision structure
y_iniInput: initial conditions for the vector of background/perturbations, already allocated and filled
Returns
the error status

◆ primordial_inflation_one_wavenumber()

int primordial_inflation_one_wavenumber ( struct perturbations ppt,
struct primordial ppm,
struct precision ppr,
double *  y_ini,
int  index_k 
)

Routine coordinating the computation of the primordial spectrum for one wavenumber. It calls primordial_inflation_one_k() to integrate the perturbation equations, and then it stores the result for the scalar/tensor spectra.

Parameters
pptInput: pointer to perturbation structure
ppmInput/output: pointer to primordial structure
pprInput: pointer to precision structure
y_iniInput: initial conditions for the vector of background/perturbations, already allocated and filled
index_kInput: index of wavenumber to be considered
Returns
the error status

Summary

  • allocate vectors for background/perturbed quantities
  • initialize the background part of the running vector
  • evolve the background until the relevant initial time for integrating perturbations
  • evolve the background/perturbation equations from this time and until some time after Horizon crossing
  • store the obtained result for curvature and tensor perturbations

◆ primordial_inflation_one_k()

int primordial_inflation_one_k ( struct primordial ppm,
struct precision ppr,
double  k,
double *  y,
double *  dy,
double *  curvature,
double *  tensor 
)

Routine integrating the background plus perturbation equations for each wavenumber, and returning the scalar and tensor spectrum.

Parameters
ppmInput: pointer to primordial structure
pprInput: pointer to precision structure
kInput: Fourier wavenumber
yInput: running vector of background/perturbations, already allocated and initialized
dyInput: running vector of background/perturbation derivatives, already allocated
curvatureOutput: curvature perturbation
tensorOutput: tensor perturbation
Returns
the error status

Summary:

  • define local variables
  • initialize the generic integrator (same integrator already used in background, thermodynamics and perturbation modules)
  • initialize variable used for deciding when to stop the calculation (= when the curvature remains stable)
  • initialize conformal time to arbitrary value (here, only variations of tau matter: the equations that we integrate do not depend explicitly on time)
  • compute derivative of initial vector and infer first value of adaptive time-step
  • loop over time
  • clean the generic integrator
  • store final value of curvature for this wavenumber
  • store final value of tensor perturbation for this wavenumber

◆ primordial_inflation_find_attractor()

int primordial_inflation_find_attractor ( struct primordial ppm,
struct precision ppr,
double  phi_0,
double  precision,
double *  y,
double *  dy,
double *  H_0,
double *  dphidt_0 
)

Routine searching for the inflationary attractor solution at a given phi_0, by iterations, with a given tolerance. If no solution found within tolerance, returns error message. The principle is the following. The code starts integrating the background equations from various values of phi, corresponding to earlier and earlier value before phi_0, and separated by a small arbitrary step size, corresponding roughly to 1 e-fold of inflation. Each time, the integration starts with the initial condition $ \phi=-V'/3H$ (slow-roll prediction). If the found value of $\phi'$ in phi_0 is stable (up to the parameter "precision"), the code considers that there is an attractor, and stops iterating. If this process does not converge, it returns an error message.

Parameters
ppmInput: pointer to primordial structure
pprInput: pointer to precision structure
phi_0Input: field value at which we wish to find the solution
precisionInput: tolerance on output values (if too large, an attractor will always considered to be found)
yInput: running vector of background variables, already allocated and initialized
dyInput: running vector of background derivatives, already allocated
H_0Output: Hubble value at phi_0 for attractor solution
dphidt_0Output: field derivative value at phi_0 for attractor solution
Returns
the error status

◆ primordial_inflation_evolve_background()

int primordial_inflation_evolve_background ( struct primordial ppm,
struct precision ppr,
double *  y,
double *  dy,
enum target_quantity  target,
double  stop,
short  check_epsilon,
enum integration_direction  direction,
enum time_definition  time 
)

Routine integrating background equations only, from initial values stored in y, to a final value (if target = aH, until aH = aH_stop; if target = phi, till phi = phi_stop; if target = end_inflation, until $ d^2a/dt^2 = 0$ (here t = proper time)). In output, y contains the final background values. In addition, if check_epsilon is true, the routine controls at each step that the expansion is accelerated and that inflation holds (wepsilon>1), otherwise it returns an error. Thanks to the last argument, it is also possible to specify whether the integration should be carried forward or backward in time. For the inflation_H case, only a 1st order differential equation is involved, so the forward and backward case can be done exactly without problems. For the inflation_V case, the equation of motion is 2nd order. What the module will do in the backward case is to search for an approximate solution, corresponding to the (first-order) attractor inflationary solution. This approximate backward solution is used in order to estimate some initial times, but the approximation made here will never impact the final result: the module is written in such a way that after using this approximation, the code always computes (and relies on) the exact forward solution.

Parameters
ppmInput: pointer to primordial structure
pprInput: pointer to precision structure
yInput/output: running vector of background variables, already allocated and initialized
dyInput: running vector of background derivatives, already allocated
targetInput: whether the goal is to reach a given aH or $ \phi $
stopInput: the target value of either aH or $ \phi $
check_epsilonInput: whether we should impose inflation (epsilon>1) at each step
directionInput: whether we should integrate forward or backward in time
timeInput: definition of time (proper or conformal)
Returns
the error status

◆ primordial_inflation_check_potential()

int primordial_inflation_check_potential ( struct primordial ppm,
double  phi,
double *  V,
double *  dV,
double *  ddV 
)

Routine checking positivity and negative slope of potential. The negative slope is an arbitrary choice. Currently the code can only deal with monotonic variations of the inflaton during inflation. So the slope had to be always negative or always positive... we took the first option.

Parameters
ppmInput: pointer to primordial structure
phiInput: field value where to perform the check
VOutput: inflaton potential in units of $ Mp^4$
dVOutput: first derivative of inflaton potential wrt the field
ddVOutput: second derivative of inflaton potential wrt the field
Returns
the error status

◆ primordial_inflation_check_hubble()

int primordial_inflation_check_hubble ( struct primordial ppm,
double  phi,
double *  H,
double *  dH,
double *  ddH,
double *  dddH 
)

Routine checking positivity and negative slope of $ H(\phi)$. The negative slope is an arbitrary choice. Currently the code can only deal with monotonic variations of the inflaton during inflation. And H can only decrease with time. So the slope $ dH/d\phi$ has to be always negative or always positive... we took the first option: phi increases, H decreases.

Parameters
ppmInput: pointer to primordial structure
phiInput: field value where to perform the check
HOutput: Hubble parameters in units of Mp
dHOutput: $ dH / d\phi $
ddHOutput: $ d^2H / d\phi^2 $
dddHOutput: $ d^3H / d\phi^3 $
Returns
the error status

◆ primordial_inflation_get_epsilon()

int primordial_inflation_get_epsilon ( struct primordial ppm,
double  phi,
double *  epsilon 
)

Routine computing the first slow-roll parameter epsilon

Parameters
ppmInput: pointer to primordial structure
phiInput: field value where to compute epsilon
epsilonOutput: result
Returns
the error status

◆ primordial_inflation_find_phi_pivot()

int primordial_inflation_find_phi_pivot ( struct primordial ppm,
struct precision ppr,
double *  y,
double *  dy 
)

Routine searching phi_pivot when a given amount of inflation is requested.

Parameters
ppmInput/output: pointer to primordial structure
pprInput: pointer to precision structure
yInput: running vector of background variables, already allocated and initialized
dyInput: running vector of background derivatives, already allocated
Returns
the error status

Summary:

  • define local variables
  • check whether in vicinity of phi_end, inflation is still ongoing
  • case in which epsilon>1: hence we must find the value phi_stop < phi_end where inflation ends up naturally
  • --> find latest value of the field such that epsilon = primordial_inflation_small_epsilon (default: 0.1)
  • --> bracketing right-hand value is phi_end (but the potential will not be evaluated exactly there, only closeby
  • --> bracketing left-hand value is found by iterating with logarithmic step until epsilon < primordial_inflation_small_epsilon
  • --> find value such that epsilon = primordial_inflation_small_epsilon by bisection
  • --> value found and stored as phi_small_epsilon
  • --> find inflationary attractor in phi_small_epsilon (should exist since epsilon<<1 there)
  • --> compute amount of inflation between this phi_small_epsilon and the end of inflation
  • --> by starting from phi_small_epsilon and integrating an approximate solution backward in time, try to estimate roughly a value close to phi_pivot but a bit smaller. This is done by trying to reach an amount of inflation equal to the requested one, minus the amount after phi_small_epsilon, and plus primordial_inflation_extra_efolds efolds (default: two). Note that it is not aggressive to require two extra e-folds of inflation before the pivot, since the calculation of the spectrum in the observable range will require even more.
  • --> find attractor in phi_try
  • --> check the total amount of inflation between phi_try and the end of inflation
  • --> go back to phi_try, and now find phi_pivot such that the amount of inflation between phi_pivot and the end of inflation is exactly the one requested.
  • case in which epsilon<1:
  • --> find inflationary attractor in phi_small_epsilon (should exist since epsilon<1 there)
  • --> by starting from phi_end and integrating an approximate solution backward in time, try to estimate roughly a value close to phi_pivot but a bit smaller. This is done by trying to reach an amount of inflation equal to the requested one, minus the amount after phi_small_epsilon, and plus primordial_inflation_extra_efolds efolds (default: two). Note that it is not aggressive to require two extra e-folds of inflation before the pivot, since the calculation of the spectrum in the observable range will require even more.
  • --> we now have a value phi_try believed to be close to and slightly smaller than phi_pivot
  • --> find attractor in phi_try
  • --> check the total amount of inflation between phi_try and the end of inflation
  • --> go back to phi_try, and now find phi_pivot such that the amount of inflation between phi_pivot and the end of inflation is exactly the one requested.
  • --> In verbose mode, check that phi_pivot is correct. Done by restarting from phi_pivot and going again till the end of inflation.

◆ primordial_inflation_derivs()

int primordial_inflation_derivs ( double  tau,
double *  y,
double *  dy,
void *  parameters_and_workspace,
ErrorMsg  error_message 
)

Routine returning derivative of system of background/perturbation variables. Like other routines used by the generic integrator (background_derivs, thermodynamics_derivs, perturbations_derivs), this routine has a generic list of arguments, and a slightly different error management, with the error message returned directly in an ErrMsg field.

Parameters
tauInput: time (not used explicitly inside the routine, but requested by the generic integrator)
yInput/output: running vector of background variables, already allocated and initialized
dyInput: running vector of background derivatives, already allocated
parameters_and_workspaceInput: all necessary input variables apart from y
error_messageOutput: error message
Returns
the error status

◆ primordial_external_spectrum_init()

int primordial_external_spectrum_init ( struct perturbations ppt,
struct primordial ppm 
)

This routine reads the primordial spectrum from an external command, and stores the tabulated values. The sampling of the k's given by the external command is preserved.

Author: Jesus Torrado (torra.nosp@m.doca.nosp@m.cho@l.nosp@m.oren.nosp@m.tz.le.nosp@m.iden.nosp@m.univ..nosp@m.nl) Date: 2013-12-20

Parameters
pptInput/output: pointer to perturbation structure
ppmInput/output: pointer to primordial structure
Returns
the error status

Summary:

  • Initialization
  • Launch the command and retrieve the output
  • Store the read results into CLASS structures
  • Make room
  • Store values
  • Release the memory used locally
  • Tell CLASS that there are scalar (and tensor) modes