CLASS MANUAL
fourier.h File Reference
#include "primordial.h"
#include "trigonometric_integrals.h"
+ Include dependency graph for fourier.h:
+ This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Data Structures

struct  fourier
 
struct  fourier_workspace
 

Macros

#define _M_EV_TOO_BIG_FOR_HALOFIT_   10.
 
#define _M_SUN_   1.98847e30
 

Detailed Description

Documented includes for trg module


Data Structure Documentation

◆ fourier

struct fourier

Structure containing all information on non-linear spectra.

Once initialized by fourier_init(), contains a table for all two points correlation functions and for all the ai,bj functions (containing the three points correlation functions), for each time and wave-number.

Data Fields
enum non_linear_method method

method for computing non-linear corrections (none, Halogit, etc.)

enum source_extrapolation extrapolation_method

method for analytical extrapolation of sources beyond pre-computed range

enum hmcode_baryonic_feedback_model feedback
double c_min

to choose between different baryonic feedback models in hmcode (dmonly, gas cooling, Agn or supernova feedback)

double eta_0

for HMcode: minimum concentration in Bullock 2001 mass-concentration relation

double z_infinity

for HMcode: halo bloating parameter

short has_pk_eq

for HMcode: z value at which Dark Energy correction is evaluated needs to be at early times (default flag: in case wa_fld is defined and non-zero, should we use the pk_eq method?

int index_md_scalars

set equal to phr->index_md_scalars (useful since this module only deals with scalars)

int ic_size

for a given mode, ic_size[index_md] = number of initial conditions included in computation

int ic_ic_size

for a given mode, ic_ic_size[index_md] = number of pairs of (index_ic1, index_ic2) with index_ic2 >= index_ic1; this number is just N(N+1)/2 where N = ic_size[index_md]

short * is_non_zero

for a given mode, is_non_zero[index_md][index_ic1_ic2] is set to true if the pair of initial conditions (index_ic1, index_ic2) are statistically correlated, or to false if they are uncorrelated

short has_pk_m

do we want spectra for total matter?

short has_pk_cb

do we want spectra for cdm+baryons?

int index_pk_m

index of pk for matter (defined only when has_pk_m is TRUE)

int index_pk_cb

index of pk for cold dark matter plus baryons (defined only when has_pk_cb is TRUE

int index_pk_total

always equal to index_pk_m (always defined, useful e.g. for weak lensing spectrum)

int index_pk_cluster

equal to index_pk_cb if it exists, otherwise to index_pk_m (always defined, useful e.g. for galaxy clustering spectrum)

int pk_size

k_size = total number of pk

short has_pk_matter

do we need matter Fourier spectrum?

int k_size

k_size = total number of k values

int k_size_pk

k_size = number of k values for P(k,z) and T(k,z) output)

double * k

k[index_k] = list of k values

double * ln_k

ln_k[index_k] = list of log(k) values

double * ln_tau

log(tau) array, only needed if user wants some output at z>0, instead of only z=0. This array only covers late times, used for the output of P(k) or T(k), and matching the condition z(tau) < z_max_pk

int ln_tau_size

total number of values in this array

int index_ln_tau_pk

first index relevant for output of P(k,z) and T(k,z)

double ** ln_pk_ic_l

Matter power spectrum (linear). Depends on indices index_pk, index_ic1_ic2, index_k, index_tau as: ln_pk_ic_l[index_pk][(index_tau * pfo->k_size + index_k)* pfo->ic_ic_size + index_ic1_ic2] where index-pk labels P(k) types (m = total matter, cb = baryons+CDM), while index_ic1_ic2 labels ordered pairs (index_ic1, index_ic2) (since the primordial spectrum is symmetric in (index_ic1, index_ic2)).

  • for diagonal elements (index_ic1 = index_ic2) this arrays contains ln[P(k)] where P(k) is positive by construction.
  • for non-diagonal elements this arrays contains the k-dependent cosine of the correlation angle, namely P(k)_(index_ic1, index_ic2)/sqrt[P(k)_index_ic1 P(k)_index_ic2] This choice is convenient since the sign of the non-diagonal cross-correlation could be negative. For fully correlated or anti-correlated initial conditions, this non-diagonal element is independent on k, and equal to +1 or -1.
double ** ddln_pk_ic_l

second derivative of above array with respect to log(tau), for spline interpolation. So:

  • for index_ic1 = index_ic, we spline ln[P(k)] vs. ln(k), which is good since this function is usually smooth.
  • for non-diagonal coefficients, we spline P(k)_(index_ic1, index_ic2)/sqrt[P(k)_index_ic1 P(k)_index_ic2] vs. ln(k), which is fine since this quantity is often assumed to be constant (e.g for fully correlated/anticorrelated initial conditions) or nearly constant, and with arbitrary sign.
double ** ln_pk_l

Total matter power spectrum summed over initial conditions (linear). Only depends on indices index_pk,index_k, index_tau as: ln_pk[index_pk][index_tau * pfo->k_size + index_k]

double ** ddln_pk_l

second derivative of above array with respect to log(tau), for spline interpolation.

double ** ln_pk_nl

Total matter power spectrum summed over initial conditions (nonlinear). Only depends on indices index_pk,index_k, index_tau as: ln_pk[index_pk][index_tau * pfo->k_size + index_k]

double ** ddln_pk_nl

second derivative of above array with respect to log(tau), for spline interpolation.

double * sigma8

sigma8[index_pk]

int k_size_extra
int tau_size

total number of k values of extrapolated k array (high k) tau_size = number of values

double * tau

tau[index_tau] = list of time values, covering all the values of the perturbation module

double ** nl_corr_density

nl_corr_density[index_pk][index_tau * ppt->k_size + index_k]

double ** k_nl

wavenumber at which non-linear corrections become important, defined differently by different non_linear_method's

int index_tau_min_nl

index of smallest value of tau at which nonlinear corrections have been computed (so, for tau<tau_min_nl, the array nl_corr_density only contains some factors 1

int index_pk_eq_w

index of w in table pk_eq_w_and_Omega

int index_pk_eq_Omega_m

index of Omega_m in table pk_eq_w_and_Omega

int pk_eq_size

number of indices in table pk_eq_w_and_Omega

int pk_eq_tau_size

number of times (and raws in table pk_eq_w_and_Omega)

double * pk_eq_tau

table of time values

double * pk_eq_w_and_Omega

table of background quantites

double * pk_eq_ddw_and_ddOmega

table of second derivatives

short fourier_verbose

amount of information written in standard output

ErrorMsg error_message

zone for writing error messages

◆ fourier_workspace

struct fourier_workspace

Structure containing variables used only internally in fourier module by various functions.

Data Fields
double * rtab
double * stab

List of R values

double * ddstab

List of Sigma Values

double * growtable

Splined sigma

double * ztable
double * tautable
double ** sigma_8
double ** sigma_disp
double ** sigma_disp_100
double ** sigma_prime
double dark_energy_correction

Macro Definition Documentation

◆ _M_EV_TOO_BIG_FOR_HALOFIT_

#define _M_EV_TOO_BIG_FOR_HALOFIT_   10.

above which value of non-CDM mass (in eV) do we stop trusting halofit?

◆ _M_SUN_

#define _M_SUN_   1.98847e30

Solar mass in Kg