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

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)
 

Detailed Description

Documented harmonic module

Julien Lesgourgues, 1.11.2019

This module computes the harmonic power spectra $ C_l^{X} $'s given the transfer functions and the primordial spectra.

The following functions can be called from other modules:

  1. harmonic_init() at the beginning (but after transfer_init())
  2. harmonic_cl_at_l() at any time for computing individual $ C_l $'s at any l
  3. harmonic_free() at the end

Function Documentation

◆ harmonic_cl_at_l()

int harmonic_cl_at_l ( struct harmonic phr,
double  l,
double *  cl_tot,
double **  cl_md,
double **  cl_md_ic 
)

Anisotropy power spectra $ C_l$'s for all types, modes and initial conditions.

This routine evaluates all the $C_l$'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.

Parameters
phrInput: pointer to harmonic structure (containing pre-computed table)
lInput: multipole number
cl_totOutput: total $C_l$'s for all types (TT, TE, EE, etc..)
cl_mdOutput: $C_l$'s for all types (TT, TE, EE, etc..) decomposed mode by mode (scalar, tensor, ...) when relevant
cl_md_icOutput: $C_l$'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
Returns
the error status

Summary:

  • define local variables
  • (a) treat case in which there is only one mode and one initial condition. Then, only cl_tot needs to be filled.
  • (b) treat case in which there is only one mode with several initial condition. Fill cl_md_ic[index_md=0] and sum it to get cl_tot.
  • (c) loop over modes
  • --> (c.1.) treat case in which the mode under consideration has only one initial condition. Fill cl_md[index_md].
  • --> (c.2.) treat case in which the mode under consideration has several initial conditions. Fill cl_md_ic[index_md] and sum it to get cl_md[index_md]
  • --> (c.3.) add contribution of cl_md[index_md] to cl_tot

◆ harmonic_init()

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 $ C_l^{X}, P(k), ... $)

Parameters
pprInput: pointer to precision structure
pbaInput: pointer to background structure (will provide H, Omega_m at redshift of interest)
pptInput: pointer to perturbation structure
ptrInput: pointer to transfer structure
ppmInput: pointer to primordial structure
pfoInput: pointer to fourier structure
phrOutput: pointer to initialized harmonic structure
Returns
the error status

Summary:

  • check that we really want to compute at least one spectrum
  • initialize indices and allocate some of the arrays in the harmonic structure
  • deal with $ C_l$'s, if any
  • a pointer to the fourier structure is stored in the spectra structure. This odd, unusual and unelegant feature has been introduced in v2.8 in order to keep in use some deprecated functions harmonic_pk_...() that are now pointing at new function fourier_pk_...(). In the future, if the deprecated functions are removed, it will be possible to remove also this pointer.

◆ harmonic_free()

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.

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

◆ harmonic_indices()

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

Parameters
pbaInput: pointer to background structure
pptInput: pointer to perturbation structure
ptrInput: pointer to transfer structure
ppmInput: pointer to primordial structure
phrInput/output: pointer to harmonic structure
Returns
the error status

◆ harmonic_cls()

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 $ C_l $'s, given the transfer functions and primordial spectra.

Parameters
pbaInput: pointer to background structure
pptInput: pointer to perturbation structure
ptrInput: pointer to transfer structure
ppmInput: pointer to primordial structure
phrInput/Output: pointer to harmonic structure
Returns
the error status

Summary:

  • define local variables
  • allocate pointers to arrays where results will be stored
  • store values of l
  • loop over modes (scalar, tensors, etc). For each mode:
  • --> (a) store number of l values for this mode
  • --> (b) allocate arrays where results will be stored
  • --> (c) loop over initial conditions
  • —> loop over l values defined in the transfer module. For each l, compute the $ C_l$'s for all types (TT, TE, ...) by convolving primordial spectra with transfer functions. This elementary task is assigned to harmonic_compute_cl()
  • --> (d) now that for a given mode, all possible $ C_l$'s have been computed, compute second derivative of the array in which they are stored, in view of spline interpolation.

◆ harmonic_compute_cl()

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 $ C_l$'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.

Parameters
pbaInput: pointer to background structure
pptInput: pointer to perturbation structure
ptrInput: pointer to transfer structure
ppmInput: pointer to primordial structure
phrInput/Output: pointer to harmonic structure (result stored here)
index_mdInput: index of mode under consideration
index_ic1Input: index of first initial condition in the correlator
index_ic2Input: index of second initial condition in the correlator
index_lInput: index of multipole under consideration
cl_integrand_num_columnsInput: number of columns in cl_integrand
cl_integrandInput: an allocated workspace
primordial_pkInput: table of primordial spectrum values
transfer_ic1Input: table of transfer function values for first initial condition
transfer_ic2Input: table of transfer function values for second initial condition
Returns
the error status

◆ harmonic_pk_at_z()

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.

Parameters
pbaInput: pointer to background structure (used for converting z into tau)
phrInput: pointer to harmonic structure (containing pre-computed table)
modeInput: linear or logarithmic
zInput: redshift
output_totOutput: total matter power spectrum P(k) in $ Mpc^3 $ (linear mode), or its logarithms (logarithmic mode)
output_icOutput: for each pair of initial conditions, matter power spectra P(k) in $ Mpc^3 $ (linear mode), or their logarithms and cross-correlation angles (logarithmic mode)
output_cb_totOutput: CDM+baryon power spectrum P_cb(k) in $ Mpc^3 $ (linear mode), or its logarithms (logarithmic mode)
output_cb_icOutput: for each pair of initial conditions, CDM+baryon power spectra P_cb(k) in $ Mpc^3 $ (linear mode), or their logarithms and cross-correlation angles (logarithmic mode)
Returns
the error status

◆ harmonic_pk_at_k_and_z()

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.

Parameters
pbaInput: pointer to background structure (used for converting z into tau)
ppmInput: pointer to primordial structure (used only in the case 0 < k < kmin)
phrInput: pointer to harmonic structure (containing pre-computed table)
kInput: wavenumber in 1/Mpc
zInput: redshift
pk_totOutput: total matter power spectrum P(k) in $ Mpc^3 $
pk_icOutput: for each pair of initial conditions, matter power spectra P(k) in $ Mpc^3$
pk_cb_totOutput: b+CDM power spectrum P(k) in $ Mpc^3 $
pk_cb_icOutput: for each pair of initial conditions, b+CDM power spectra P(k) in $ Mpc^3$
Returns
the error status

◆ harmonic_pk_nl_at_z()

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.

Parameters
pbaInput: pointer to background structure (used for converting z into tau)
phrInput: pointer to harmonic structure (containing pre-computed table)
modeInput: linear or logarithmic
zInput: redshift
output_totOutput: total matter power spectrum P(k) in $ Mpc^3$ (linear mode), or its logarithms (logarithmic mode)
output_cb_totOutput: b+CDM power spectrum P(k) in $ Mpc^3$ (linear mode), or its logarithms (logarithmic mode)
Returns
the error status

◆ harmonic_pk_nl_at_k_and_z()

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.

Parameters
pbaInput: pointer to background structure (used for converting z into tau)
ppmInput: pointer to primordial structure (used only in the case 0 < k < kmin)
phrInput: pointer to harmonic structure (containing pre-computed table)
kInput: wavenumber in 1/Mpc
zInput: redshift
pk_totOutput: total matter power spectrum P(k) in $ Mpc^3$
pk_cb_totOutput: b+CDM power spectrum P(k) in $ Mpc^3$
Returns
the error status

◆ harmonic_fast_pk_at_kvec_and_zvec()

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.

Parameters
pbaInput: pointer to background structure
phrInput: pointer to harmonic structure
kvecInput: array of wavenumbers in ascending order (in 1/Mpc)
kvec_sizeInput: size of array of wavenumbers
zvecInput: array of redshifts in arbitrary order
zvec_sizeInput: size of array of redshifts
pk_tot_outOutput: P(k_i,z_j) for total matter (if available) in Mpc**3
pk_cb_tot_outOutput: P_cb(k_i,z_j) for cdm+baryons (if available) in Mpc**3
nonlinearInput: TRUE or FALSE (to output nonlinear or linear P(k,z))
Returns
the error status

◆ harmonic_sigma()

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.

Parameters
pbaInput: pointer to background structure
ppmInput: pointer to primordial structure
phrInput: pointer to harmonic structure
RInput: radius in Mpc
zInput: redshift
sigmaOutput: variance in a sphere of radius R (dimensionless)
Returns
the error status

◆ harmonic_sigma_cb()

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.

Parameters
pbaInput: pointer to background structure
ppmInput: pointer to primordial structure
phrInput: pointer to harmonic structure
RInput: radius in Mpc
zInput: redshift
sigma_cbOutput: variance in a sphere of radius R (dimensionless)
Returns
the error status

◆ harmonic_tk_at_z()

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)

Parameters
pbaInput: pointer to background structure (used for converting z into tau)
phrInput: pointer to harmonic structure (containing pre-computed table)
zInput: redshift
outputOutput: matter transfer functions
Returns
the error status

◆ harmonic_tk_at_k_and_z()

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)

Parameters
pbaInput: pointer to background structure (used for converting z into tau)
phrInput: pointer to harmonic structure (containing pre-computed table)
kInput: wavenumber in 1/Mpc
zInput: redshift
outputOutput: matter transfer functions
Returns
the error status