|
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) |
|
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:
- primordial_init() at the beginning (anytime after perturbations_init() and before harmonic_init())
- primordial_spectrum_at_k() at any time for computing P(k) at any k
- primordial_free() at the end
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 (from -1 to 1) instead of
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
-
ppm | Input: pointer to primordial structure containing tabulated primordial spectrum |
index_md | Input: index of mode (scalar, tensor, ...) |
mode | Input: linear or logarithmic |
input | Input: wavenumber in 1/Mpc (linear mode) or its logarithm (logarithmic mode) |
output | Output: for each pair of initial conditions, primordial spectra P(k) in (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
This routine initializes the primordial structure (in particular, it computes table of primordial spectrum values)
- Parameters
-
ppr | Input: pointer to precision structure (defines method and precision for all computations) |
ppt | Input: pointer to perturbation structure (useful for knowing k_min, k_max, etc.) |
ppm | Output: 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 '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 or
- deal with the case of external calculation of
- compute second derivative of each 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)
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)
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 (slow-roll prediction). If the found value of 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
-
ppm | Input: pointer to primordial structure |
ppr | Input: pointer to precision structure |
phi_0 | Input: field value at which we wish to find the solution |
precision | Input: tolerance on output values (if too large, an attractor will always considered to be found) |
y | Input: running vector of background variables, already allocated and initialized |
dy | Input: running vector of background derivatives, already allocated |
H_0 | Output: Hubble value at phi_0 for attractor solution |
dphidt_0 | Output: field derivative value at phi_0 for attractor solution |
- Returns
- the error status
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 (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
-
ppm | Input: pointer to primordial structure |
ppr | Input: pointer to precision structure |
y | Input/output: running vector of background variables, already allocated and initialized |
dy | Input: running vector of background derivatives, already allocated |
target | Input: whether the goal is to reach a given aH or |
stop | Input: the target value of either aH or |
check_epsilon | Input: whether we should impose inflation (epsilon>1) at each step |
direction | Input: whether we should integrate forward or backward in time |
time | Input: definition of time (proper or conformal) |
- Returns
- the error status
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
-
tau | Input: time (not used explicitly inside the routine, but requested by the generic integrator) |
y | Input/output: running vector of background variables, already allocated and initialized |
dy | Input: running vector of background derivatives, already allocated |
parameters_and_workspace | Input: all necessary input variables apart from y |
error_message | Output: error message |
- Returns
- the error status