CLASS MANUAL
|
#include "harmonic.h"
Functions | |
int | harmonic_cl_at_l (struct harmonic *phr, double l, double *cl_tot, double **cl_md, double **cl_md_ic) |
int | harmonic_init (struct precision *ppr, struct background *pba, struct perturbations *ppt, struct primordial *ppm, struct fourier *pfo, struct transfer *ptr, struct harmonic *phr) |
int | harmonic_free (struct harmonic *phr) |
int | harmonic_indices (struct background *pba, struct perturbations *ppt, struct transfer *ptr, struct primordial *ppm, struct harmonic *phr) |
int | harmonic_cls (struct background *pba, struct perturbations *ppt, struct transfer *ptr, struct primordial *ppm, struct harmonic *phr) |
int | harmonic_compute_cl (struct background *pba, struct perturbations *ppt, struct transfer *ptr, struct primordial *ppm, struct harmonic *phr, int index_md, int index_ic1, int index_ic2, int index_l, int cl_integrand_num_columns, double *cl_integrand, double *primordial_pk, double *transfer_ic1, double *transfer_ic2) |
int | harmonic_pk_at_z (struct background *pba, struct harmonic *phr, enum linear_or_logarithmic mode, double z, double *output_tot, double *output_ic, double *output_cb_tot, double *output_cb_ic) |
int | harmonic_pk_at_k_and_z (struct background *pba, struct primordial *ppm, struct harmonic *phr, double k, double z, double *pk_tot, double *pk_ic, double *pk_cb_tot, double *pk_cb_ic) |
int | harmonic_pk_nl_at_z (struct background *pba, struct harmonic *phr, enum linear_or_logarithmic mode, double z, double *output_tot, double *output_cb_tot) |
int | harmonic_pk_nl_at_k_and_z (struct background *pba, struct primordial *ppm, struct harmonic *phr, double k, double z, double *pk_tot, double *pk_cb_tot) |
int | harmonic_fast_pk_at_kvec_and_zvec (struct background *pba, struct harmonic *phr, double *kvec, int kvec_size, double *zvec, int zvec_size, double *pk_tot_out, double *pk_cb_tot_out, int nonlinear) |
int | harmonic_sigma (struct background *pba, struct primordial *ppm, struct harmonic *phr, double R, double z, double *sigma) |
int | harmonic_sigma_cb (struct background *pba, struct primordial *ppm, struct harmonic *phr, double R, double z, double *sigma_cb) |
int | harmonic_tk_at_z (struct background *pba, struct harmonic *phr, double z, double *output) |
int | harmonic_tk_at_k_and_z (struct background *pba, struct harmonic *phr, double k, double z, double *output) |
Documented harmonic module
Julien Lesgourgues, 1.11.2019
This module computes the harmonic power spectra 's given the transfer functions and the primordial spectra.
The following functions can be called from other modules:
int harmonic_cl_at_l | ( | struct harmonic * | phr, |
double | l, | ||
double * | cl_tot, | ||
double ** | cl_md, | ||
double ** | cl_md_ic | ||
) |
Anisotropy power spectra 's for all types, modes and initial conditions.
This routine evaluates all the 's at a given value of l by interpolating in the pre-computed table. When relevant, it also sums over all initial conditions for each mode, and over all modes.
This function can be called from whatever module at whatever time, provided that harmonic_init() has been called before, and harmonic_free() has not been called yet.
phr | Input: pointer to harmonic structure (containing pre-computed table) |
l | Input: multipole number |
cl_tot | Output: total 's for all types (TT, TE, EE, etc..) |
cl_md | Output: 's for all types (TT, TE, EE, etc..) decomposed mode by mode (scalar, tensor, ...) when relevant |
cl_md_ic | Output: 's for all types (TT, TE, EE, etc..) decomposed by pairs of initial conditions (adiabatic, isocurvatures) for each mode (usually, only for the scalar mode) when relevant |
Summary:
int harmonic_init | ( | struct precision * | ppr, |
struct background * | pba, | ||
struct perturbations * | ppt, | ||
struct primordial * | ppm, | ||
struct fourier * | pfo, | ||
struct transfer * | ptr, | ||
struct harmonic * | phr | ||
) |
This routine initializes the harmonic structure (in particular, computes table of anisotropy and Fourier spectra )
ppr | Input: pointer to precision structure |
pba | Input: pointer to background structure (will provide H, Omega_m at redshift of interest) |
ppt | Input: pointer to perturbation structure |
ptr | Input: pointer to transfer structure |
ppm | Input: pointer to primordial structure |
pfo | Input: pointer to fourier structure |
phr | Output: pointer to initialized harmonic structure |
Summary:
int harmonic_free | ( | struct harmonic * | phr | ) |
This routine frees all the memory space allocated by harmonic_init().
To be called at the end of each run, only when no further calls to harmonic_cls_at_l(), harmonic_pk_at_z(), harmonic_pk_at_k_and_z() are needed.
phr | Input: pointer to harmonic structure (which fields must be freed) |
int harmonic_indices | ( | struct background * | pba, |
struct perturbations * | ppt, | ||
struct transfer * | ptr, | ||
struct primordial * | ppm, | ||
struct harmonic * | phr | ||
) |
This routine defines indices and allocates tables in the harmonic structure
pba | Input: pointer to background structure |
ppt | Input: pointer to perturbation structure |
ptr | Input: pointer to transfer structure |
ppm | Input: pointer to primordial structure |
phr | Input/output: pointer to harmonic structure |
int harmonic_cls | ( | struct background * | pba, |
struct perturbations * | ppt, | ||
struct transfer * | ptr, | ||
struct primordial * | ppm, | ||
struct harmonic * | phr | ||
) |
This routine computes a table of values for all harmonic spectra 's, given the transfer functions and primordial spectra.
pba | Input: pointer to background structure |
ppt | Input: pointer to perturbation structure |
ptr | Input: pointer to transfer structure |
ppm | Input: pointer to primordial structure |
phr | Input/Output: pointer to harmonic structure |
Summary:
int harmonic_compute_cl | ( | struct background * | pba, |
struct perturbations * | ppt, | ||
struct transfer * | ptr, | ||
struct primordial * | ppm, | ||
struct harmonic * | phr, | ||
int | index_md, | ||
int | index_ic1, | ||
int | index_ic2, | ||
int | index_l, | ||
int | cl_integrand_num_columns, | ||
double * | cl_integrand, | ||
double * | primordial_pk, | ||
double * | transfer_ic1, | ||
double * | transfer_ic2 | ||
) |
This routine computes the 's for a given mode, pair of initial conditions and multipole, but for all types (TT, TE...), by convolving the transfer functions with the primordial spectra.
pba | Input: pointer to background structure |
ppt | Input: pointer to perturbation structure |
ptr | Input: pointer to transfer structure |
ppm | Input: pointer to primordial structure |
phr | Input/Output: pointer to harmonic structure (result stored here) |
index_md | Input: index of mode under consideration |
index_ic1 | Input: index of first initial condition in the correlator |
index_ic2 | Input: index of second initial condition in the correlator |
index_l | Input: index of multipole under consideration |
cl_integrand_num_columns | Input: number of columns in cl_integrand |
cl_integrand | Input: an allocated workspace |
primordial_pk | Input: table of primordial spectrum values |
transfer_ic1 | Input: table of transfer function values for first initial condition |
transfer_ic2 | Input: table of transfer function values for second initial condition |
int harmonic_pk_at_z | ( | struct background * | pba, |
struct harmonic * | phr, | ||
enum linear_or_logarithmic | mode, | ||
double | z, | ||
double * | output_tot, | ||
double * | output_ic, | ||
double * | output_cb_tot, | ||
double * | output_cb_ic | ||
) |
Matter power spectrum for arbitrary redshift and for all initial conditions.
This function is deprecated since v2.8. Try using fourier_pk_at_z() instead.
pba | Input: pointer to background structure (used for converting z into tau) |
phr | Input: pointer to harmonic structure (containing pre-computed table) |
mode | Input: linear or logarithmic |
z | Input: redshift |
output_tot | Output: total matter power spectrum P(k) in (linear mode), or its logarithms (logarithmic mode) |
output_ic | Output: for each pair of initial conditions, matter power spectra P(k) in (linear mode), or their logarithms and cross-correlation angles (logarithmic mode) |
output_cb_tot | Output: CDM+baryon power spectrum P_cb(k) in (linear mode), or its logarithms (logarithmic mode) |
output_cb_ic | Output: for each pair of initial conditions, CDM+baryon power spectra P_cb(k) in (linear mode), or their logarithms and cross-correlation angles (logarithmic mode) |
int harmonic_pk_at_k_and_z | ( | struct background * | pba, |
struct primordial * | ppm, | ||
struct harmonic * | phr, | ||
double | k, | ||
double | z, | ||
double * | pk_tot, | ||
double * | pk_ic, | ||
double * | pk_cb_tot, | ||
double * | pk_cb_ic | ||
) |
Matter power spectrum for arbitrary wavenumber, redshift and initial condition.
This function is deprecated since v2.8. Try using fourier_pk_linear_at_k_and_z() instead.
pba | Input: pointer to background structure (used for converting z into tau) |
ppm | Input: pointer to primordial structure (used only in the case 0 < k < kmin) |
phr | Input: pointer to harmonic structure (containing pre-computed table) |
k | Input: wavenumber in 1/Mpc |
z | Input: redshift |
pk_tot | Output: total matter power spectrum P(k) in |
pk_ic | Output: for each pair of initial conditions, matter power spectra P(k) in |
pk_cb_tot | Output: b+CDM power spectrum P(k) in |
pk_cb_ic | Output: for each pair of initial conditions, b+CDM power spectra P(k) in |
int harmonic_pk_nl_at_z | ( | struct background * | pba, |
struct harmonic * | phr, | ||
enum linear_or_logarithmic | mode, | ||
double | z, | ||
double * | output_tot, | ||
double * | output_cb_tot | ||
) |
Non-linear total matter power spectrum for arbitrary redshift.
This function is deprecated since v2.8. Try using fourier_pk_at_z() instead.
pba | Input: pointer to background structure (used for converting z into tau) |
phr | Input: pointer to harmonic structure (containing pre-computed table) |
mode | Input: linear or logarithmic |
z | Input: redshift |
output_tot | Output: total matter power spectrum P(k) in (linear mode), or its logarithms (logarithmic mode) |
output_cb_tot | Output: b+CDM power spectrum P(k) in (linear mode), or its logarithms (logarithmic mode) |
int harmonic_pk_nl_at_k_and_z | ( | struct background * | pba, |
struct primordial * | ppm, | ||
struct harmonic * | phr, | ||
double | k, | ||
double | z, | ||
double * | pk_tot, | ||
double * | pk_cb_tot | ||
) |
Non-linear total matter power spectrum for arbitrary wavenumber and redshift.
This function is deprecated since v2.8. Try using fourier_pk_at_k_and_z() instead.
pba | Input: pointer to background structure (used for converting z into tau) |
ppm | Input: pointer to primordial structure (used only in the case 0 < k < kmin) |
phr | Input: pointer to harmonic structure (containing pre-computed table) |
k | Input: wavenumber in 1/Mpc |
z | Input: redshift |
pk_tot | Output: total matter power spectrum P(k) in |
pk_cb_tot | Output: b+CDM power spectrum P(k) in |
int harmonic_fast_pk_at_kvec_and_zvec | ( | struct background * | pba, |
struct harmonic * | phr, | ||
double * | kvec, | ||
int | kvec_size, | ||
double * | zvec, | ||
int | zvec_size, | ||
double * | pk_tot_out, | ||
double * | pk_cb_tot_out, | ||
int | nonlinear | ||
) |
Return the P(k,z) for a grid of (k_i,z_j) passed in input, for all available pk types (_m, _cb), either linear or nonlinear depending on input.
This function is deprecated since v2.8. Try using fourier_pks_at_kvec_and_zvec() instead.
pba | Input: pointer to background structure |
phr | Input: pointer to harmonic structure |
kvec | Input: array of wavenumbers in ascending order (in 1/Mpc) |
kvec_size | Input: size of array of wavenumbers |
zvec | Input: array of redshifts in arbitrary order |
zvec_size | Input: size of array of redshifts |
pk_tot_out | Output: P(k_i,z_j) for total matter (if available) in Mpc**3 |
pk_cb_tot_out | Output: P_cb(k_i,z_j) for cdm+baryons (if available) in Mpc**3 |
nonlinear | Input: TRUE or FALSE (to output nonlinear or linear P(k,z)) |
int harmonic_sigma | ( | struct background * | pba, |
struct primordial * | ppm, | ||
struct harmonic * | phr, | ||
double | R, | ||
double | z, | ||
double * | sigma | ||
) |
This routine computes sigma(R) given P(k) for total matter power spectrum (does not check that k_max is large enough)
This function is deprecated since v2.8. Try using fourier_sigmas_at_z() instead.
pba | Input: pointer to background structure |
ppm | Input: pointer to primordial structure |
phr | Input: pointer to harmonic structure |
R | Input: radius in Mpc |
z | Input: redshift |
sigma | Output: variance in a sphere of radius R (dimensionless) |
int harmonic_sigma_cb | ( | struct background * | pba, |
struct primordial * | ppm, | ||
struct harmonic * | phr, | ||
double | R, | ||
double | z, | ||
double * | sigma_cb | ||
) |
This routine computes sigma(R) given P(k) for baryon+cdm power spectrum (does not check that k_max is large enough)
This function is deprecated since v2.8. Try using fourier_sigmas_at_z() instead.
pba | Input: pointer to background structure |
ppm | Input: pointer to primordial structure |
phr | Input: pointer to harmonic structure |
R | Input: radius in Mpc |
z | Input: redshift |
sigma_cb | Output: variance in a sphere of radius R (dimensionless) |
int harmonic_tk_at_z | ( | struct background * | pba, |
struct harmonic * | phr, | ||
double | z, | ||
double * | output | ||
) |
Obsolete function, superseeded by perturbations_sources_at_tau() (at the time of the switch, this function was anyway never used anywhere)
pba | Input: pointer to background structure (used for converting z into tau) |
phr | Input: pointer to harmonic structure (containing pre-computed table) |
z | Input: redshift |
output | Output: matter transfer functions |
int harmonic_tk_at_k_and_z | ( | struct background * | pba, |
struct harmonic * | phr, | ||
double | k, | ||
double | z, | ||
double * | output | ||
) |
Obsolete function, superseeded by perturbations_sources_at_tau() (at the time of the switch, this function was anyway never used anywhere)
pba | Input: pointer to background structure (used for converting z into tau) |
phr | Input: pointer to harmonic structure (containing pre-computed table) |
k | Input: wavenumber in 1/Mpc |
z | Input: redshift |
output | Output: matter transfer functions |