CLASS MANUAL
lensing.c File Reference
#include "lensing.h"
#include <time.h>
+ Include dependency graph for lensing.c:

Functions

int lensing_cl_at_l (struct lensing *ple, int l, double *cl_lensed)
 
int lensing_init (struct precision *ppr, struct perturbations *ppt, struct harmonic *phr, struct fourier *pfo, struct lensing *ple)
 
int lensing_free (struct lensing *ple)
 
int lensing_indices (struct precision *ppr, struct harmonic *phr, struct lensing *ple)
 
int lensing_lensed_cl_tt (double *ksi, double **d00, double *w8, int nmu, struct lensing *ple)
 
int lensing_addback_cl_tt (struct lensing *ple, double *cl_tt)
 
int lensing_lensed_cl_te (double *ksiX, double **d20, double *w8, int nmu, struct lensing *ple)
 
int lensing_addback_cl_te (struct lensing *ple, double *cl_te)
 
int lensing_lensed_cl_ee_bb (double *ksip, double *ksim, double **d22, double **d2m2, double *w8, int nmu, struct lensing *ple)
 
int lensing_addback_cl_ee_bb (struct lensing *ple, double *cl_ee, double *cl_bb)
 
int lensing_d00 (double *mu, int num_mu, int lmax, double **d00)
 
int lensing_d11 (double *mu, int num_mu, int lmax, double **d11)
 
int lensing_d1m1 (double *mu, int num_mu, int lmax, double **d1m1)
 
int lensing_d2m2 (double *mu, int num_mu, int lmax, double **d2m2)
 
int lensing_d22 (double *mu, int num_mu, int lmax, double **d22)
 
int lensing_d20 (double *mu, int num_mu, int lmax, double **d20)
 
int lensing_d31 (double *mu, int num_mu, int lmax, double **d31)
 
int lensing_d3m1 (double *mu, int num_mu, int lmax, double **d3m1)
 
int lensing_d3m3 (double *mu, int num_mu, int lmax, double **d3m3)
 
int lensing_d40 (double *mu, int num_mu, int lmax, double **d40)
 
int lensing_d4m2 (double *mu, int num_mu, int lmax, double **d4m2)
 
int lensing_d4m4 (double *mu, int num_mu, int lmax, double **d4m4)
 

Detailed Description

Documented lensing module

Simon Prunet and Julien Lesgourgues, 6.12.2010

This module computes the lensed temperature and polarization anisotropy power spectra $ C_l^{X}, P(k), ... $'s given the unlensed temperature, polarization and lensing potential spectra.

Follows Challinor and Lewis full-sky method, astro-ph/0502425

The following functions can be called from other modules:

  1. lensing_init() at the beginning (but after harmonic_init())
  2. lensing_cl_at_l() at any time for computing Cl_lensed at any l
  3. lensing_free() at the end

Function Documentation

◆ lensing_cl_at_l()

int lensing_cl_at_l ( struct lensing ple,
int  l,
double *  cl_lensed 
)

Anisotropy power spectra $ C_l$'s for all types, modes and initial conditions. SO FAR: ONLY SCALAR

This routine evaluates all the lensed $ C_l$'s at a given value of l by picking it 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 lensing_init() has been called before, and lensing_free() has not been called yet.

Parameters
pleInput: pointer to lensing structure
lInput: multipole number
cl_lensedOutput: lensed $ C_l$'s for all types (TT, TE, EE, etc..)
Returns
the error status

◆ lensing_init()

int lensing_init ( struct precision ppr,
struct perturbations ppt,
struct harmonic phr,
struct fourier pfo,
struct lensing ple 
)

This routine initializes the lensing structure (in particular, computes table of lensed anisotropy spectra $ C_l^{X} $)

Parameters
pprInput: pointer to precision structure
pptInput: pointer to perturbation structure (just in case, not used in current version...)
phrInput: pointer to harmonic structure
pfoInput: pointer to fourier structure
pleOutput: pointer to initialized lensing structure
Returns
the error status

Summary:

  • Define local variables
  • check that we really want to compute at least one spectrum
  • initialize indices and allocate some of the arrays in the lensing structure
  • put all precision variables hare; will be stored later in precision structure
  • Last element in $ \mu $ will be for $ \mu=1 $, needed for sigma2. The rest will be chosen as roots of a Gauss-Legendre quadrature
  • allocate array of $ \mu $ values, as well as quadrature weights
  • Compute $ d^l_{mm'} (\mu) $
  • Allocate main contiguous buffer
  • compute $ Cgl(\mu)$, $ Cgl2(\mu) $ and sigma2( $\mu$)
  • Locally store unlensed temperature $ cl_{tt}$ and potential $ cl_{pp}$ spectra
  • Compute sigma2 $(\mu)$ and Cgl2( $\mu$)
  • compute ksi, ksi+, ksi-, ksiX
  • --> ksi is for TT
  • --> ksiX is for TE
  • --> ksip, ksim for EE, BB
  • compute lensed $ C_l$'s by integration
  • spline computed $ C_l$'s in view of interpolation
  • Free lots of stuff
  • Exit

◆ lensing_free()

int lensing_free ( struct lensing ple)

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

To be called at the end of each run, only when no further calls to lensing_cl_at_l() are needed.

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

◆ lensing_indices()

int lensing_indices ( struct precision ppr,
struct harmonic phr,
struct lensing ple 
)

This routine defines indices and allocates tables in the lensing structure

Parameters
pprInput: pointer to precision structure
phrInput: pointer to harmonic structure
pleInput/output: pointer to lensing structure
Returns
the error status

◆ lensing_lensed_cl_tt()

int lensing_lensed_cl_tt ( double *  ksi,
double **  d00,
double *  w8,
int  nmu,
struct lensing ple 
)

This routine computes the lensed power spectra by Gaussian quadrature

Parameters
ksiInput: Lensed correlation function (ksi[index_mu])
d00Input: Legendre polynomials ( $ d^l_{00}$[l][index_mu])
w8Input: Legendre quadrature weights (w8[index_mu])
nmuInput: Number of quadrature points (0<=index_mu<=nmu)
pleInput/output: Pointer to the lensing structure
Returns
the error status

Integration by Gauss-Legendre quadrature.

◆ lensing_addback_cl_tt()

int lensing_addback_cl_tt ( struct lensing ple,
double *  cl_tt 
)

This routine adds back the unlensed $ cl_{tt}$ power spectrum Used in case of fast (and BB inaccurate) integration of correlation functions.

Parameters
pleInput/output: Pointer to the lensing structure
cl_ttInput: Array of unlensed power spectrum
Returns
the error status

◆ lensing_lensed_cl_te()

int lensing_lensed_cl_te ( double *  ksiX,
double **  d20,
double *  w8,
int  nmu,
struct lensing ple 
)

This routine computes the lensed power spectra by Gaussian quadrature

Parameters
ksiXInput: Lensed correlation function (ksiX[index_mu])
d20Input: Wigner d-function ( $ d^l_{20}$[l][index_mu])
w8Input: Legendre quadrature weights (w8[index_mu])
nmuInput: Number of quadrature points (0<=index_mu<=nmu)
pleInput/output: Pointer to the lensing structure
Returns
the error status

Integration by Gauss-Legendre quadrature.

◆ lensing_addback_cl_te()

int lensing_addback_cl_te ( struct lensing ple,
double *  cl_te 
)

This routine adds back the unlensed $ cl_{te}$ power spectrum Used in case of fast (and BB inaccurate) integration of correlation functions.

Parameters
pleInput/output: Pointer to the lensing structure
cl_teInput: Array of unlensed power spectrum
Returns
the error status

◆ lensing_lensed_cl_ee_bb()

int lensing_lensed_cl_ee_bb ( double *  ksip,
double *  ksim,
double **  d22,
double **  d2m2,
double *  w8,
int  nmu,
struct lensing ple 
)

This routine computes the lensed power spectra by Gaussian quadrature

Parameters
ksipInput: Lensed correlation function (ksi+[index_mu])
ksimInput: Lensed correlation function (ksi-[index_mu])
d22Input: Wigner d-function ( $ d^l_{22}$[l][index_mu])
d2m2Input: Wigner d-function ( $ d^l_{2-2}$[l][index_mu])
w8Input: Legendre quadrature weights (w8[index_mu])
nmuInput: Number of quadrature points (0<=index_mu<=nmu)
pleInput/output: Pointer to the lensing structure
Returns
the error status

Integration by Gauss-Legendre quadrature.

◆ lensing_addback_cl_ee_bb()

int lensing_addback_cl_ee_bb ( struct lensing ple,
double *  cl_ee,
double *  cl_bb 
)

This routine adds back the unlensed $ cl_{ee}$, $ cl_{bb}$ power spectra Used in case of fast (and BB inaccurate) integration of correlation functions.

Parameters
pleInput/output: Pointer to the lensing structure
cl_eeInput: Array of unlensed power spectrum
cl_bbInput: Array of unlensed power spectrum
Returns
the error status

◆ lensing_d00()

int lensing_d00 ( double *  mu,
int  num_mu,
int  lmax,
double **  d00 
)

This routine computes the d00 term

Parameters
muInput: Vector of cos(beta) values
num_muInput: Number of cos(beta) values
lmaxInput: maximum multipole
d00Input/output: Result is stored here

Wigner d-functions, computed by recurrence actual recurrence on $ \sqrt{(2l+1)/2} d^l_{mm'} $ for stability Formulae from Kostelec & Rockmore 2003

◆ lensing_d11()

int lensing_d11 ( double *  mu,
int  num_mu,
int  lmax,
double **  d11 
)

This routine computes the d11 term

Parameters
muInput: Vector of cos(beta) values
num_muInput: Number of cos(beta) values
lmaxInput: maximum multipole
d11Input/output: Result is stored here

Wigner d-functions, computed by recurrence actual recurrence on $ \sqrt{(2l+1)/2} d^l_{mm'} $ for stability Formulae from Kostelec & Rockmore 2003

◆ lensing_d1m1()

int lensing_d1m1 ( double *  mu,
int  num_mu,
int  lmax,
double **  d1m1 
)

This routine computes the d1m1 term

Parameters
muInput: Vector of cos(beta) values
num_muInput: Number of cos(beta) values
lmaxInput: maximum multipole
d1m1Input/output: Result is stored here

Wigner d-functions, computed by recurrence actual recurrence on $ \sqrt{(2l+1)/2} d^l_{mm'} $ for stability Formulae from Kostelec & Rockmore 2003

◆ lensing_d2m2()

int lensing_d2m2 ( double *  mu,
int  num_mu,
int  lmax,
double **  d2m2 
)

This routine computes the d2m2 term

Parameters
muInput: Vector of cos(beta) values
num_muInput: Number of cos(beta) values
lmaxInput: maximum multipole
d2m2Input/output: Result is stored here

Wigner d-functions, computed by recurrence actual recurrence on $ \sqrt{(2l+1)/2} d^l_{mm'} $ for stability Formulae from Kostelec & Rockmore 2003

◆ lensing_d22()

int lensing_d22 ( double *  mu,
int  num_mu,
int  lmax,
double **  d22 
)

This routine computes the d22 term

Parameters
muInput: Vector of cos(beta) values
num_muInput: Number of cos(beta) values
lmaxInput: maximum multipole
d22Input/output: Result is stored here

Wigner d-functions, computed by recurrence actual recurrence on $ \sqrt{(2l+1)/2} d^l_{mm'} $ for stability Formulae from Kostelec & Rockmore 2003

◆ lensing_d20()

int lensing_d20 ( double *  mu,
int  num_mu,
int  lmax,
double **  d20 
)

This routine computes the d20 term

Parameters
muInput: Vector of cos(beta) values
num_muInput: Number of cos(beta) values
lmaxInput: maximum multipole
d20Input/output: Result is stored here

Wigner d-functions, computed by recurrence actual recurrence on $ \sqrt{(2l+1)/2} d^l_{mm'} $ for stability Formulae from Kostelec & Rockmore 2003

◆ lensing_d31()

int lensing_d31 ( double *  mu,
int  num_mu,
int  lmax,
double **  d31 
)

This routine computes the d31 term

Parameters
muInput: Vector of cos(beta) values
num_muInput: Number of cos(beta) values
lmaxInput: maximum multipole
d31Input/output: Result is stored here

Wigner d-functions, computed by recurrence actual recurrence on $ \sqrt{(2l+1)/2} d^l_{mm'} $ for stability Formulae from Kostelec & Rockmore 2003

◆ lensing_d3m1()

int lensing_d3m1 ( double *  mu,
int  num_mu,
int  lmax,
double **  d3m1 
)

This routine computes the d3m1 term

Parameters
muInput: Vector of cos(beta) values
num_muInput: Number of cos(beta) values
lmaxInput: maximum multipole
d3m1Input/output: Result is stored here

Wigner d-functions, computed by recurrence actual recurrence on $ \sqrt{(2l+1)/2} d^l_{mm'} $ for stability Formulae from Kostelec & Rockmore 2003

◆ lensing_d3m3()

int lensing_d3m3 ( double *  mu,
int  num_mu,
int  lmax,
double **  d3m3 
)

This routine computes the d3m3 term

Parameters
muInput: Vector of cos(beta) values
num_muInput: Number of cos(beta) values
lmaxInput: maximum multipole
d3m3Input/output: Result is stored here

Wigner d-functions, computed by recurrence actual recurrence on $ \sqrt{(2l+1)/2} d^l_{mm'} $ for stability Formulae from Kostelec & Rockmore 2003

◆ lensing_d40()

int lensing_d40 ( double *  mu,
int  num_mu,
int  lmax,
double **  d40 
)

This routine computes the d40 term

Parameters
muInput: Vector of cos(beta) values
num_muInput: Number of cos(beta) values
lmaxInput: maximum multipole
d40Input/output: Result is stored here

Wigner d-functions, computed by recurrence actual recurrence on $ \sqrt{(2l+1)/2} d^l_{mm'} $ for stability Formulae from Kostelec & Rockmore 2003

◆ lensing_d4m2()

int lensing_d4m2 ( double *  mu,
int  num_mu,
int  lmax,
double **  d4m2 
)

This routine computes the d4m2 term

Parameters
muInput: Vector of cos(beta) values
num_muInput: Number of cos(beta) values
lmaxInput: maximum multipole
d4m2Input/output: Result is stored here

Wigner d-functions, computed by recurrence actual recurrence on $ \sqrt{(2l+1)/2} d^l_{mm'} $ for stability Formulae from Kostelec & Rockmore 2003

◆ lensing_d4m4()

int lensing_d4m4 ( double *  mu,
int  num_mu,
int  lmax,
double **  d4m4 
)

This routine computes the d4m4 term

Parameters
muInput: Vector of cos(beta) values
num_muInput: Number of cos(beta) values
lmaxInput: maximum multipole
d4m4Input/output: Result is stored here

Wigner d-functions, computed by recurrence actual recurrence on $ \sqrt{(2l+1)/2} d^l_{mm'} $ for stability Formulae from Kostelec & Rockmore 2003