flatsky – Flatsky analysis tools

rec_lens – quadratic lensing reconstruction

cmblensplus.flatsky.rec_lens.qtt(nx, ny, D, rL, fC, T1, T2, gtype='')

Reconstructing CMB lensing potential and its curl mode from the temperature quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

Temperature power spectrum on 2D grid, with bounds (nx,ny)

T1[lx,ly] (dcmplx):

2D Fourier modes of 1st inverse-variance filtered temperature, with bounds (nx,ny)

T2[lx,ly] (dcmplx):

2D Fourier modes of 2nd inverse-variance filtered temperature, with bounds (nx,ny)

Args(optional):
gtype (str):

Type of output, i.e., convergence (gtype=’k’) or lensing potential (gtype=’’, default)

Returns:
glm[lx,ly] (dcmplx):

2D Fourier modes of CMB lensing potential, with bounds (nx,ny)

clm[lx,ly] (dcmplx):

2D Fourier modes of Curl mode (pseudo lensing potential), with bounds (nx,ny)

Usage:
glm,clm = flatsky.rec_lens.qtt(nx,ny,D,rL,fC,T1,T2,gtype):

cmblensplus.flatsky.rec_lens.qte(nx, ny, D, rL, fC, T, E, gtype='')

Reconstructing CMB lensing potential and its curl mode from the suboptimal TE quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

TE cross power spectrum on 2D grid, with bounds (nx,ny)

T[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered temperature, with bounds (nx,ny)

E[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered E-mode, with bounds (nx,ny)

Args(optional):
gtype (str):

Type of output, i.e., convergence (gtype=’k’) or lensing potential (gtype=’’, default)

Returns:
glm[lx,ly] (dcmplx):

2D Fourier modes of CMB lensing potential, with bounds (nx,ny)

clm[lx,ly] (dcmplx):

2D Fourier modes of Curl mode (pseudo lensing potential), with bounds (nx,ny)

Usage:
glm,clm = flatsky.rec_lens.qte(nx,ny,D,rL,fC,T,E,gtype):

cmblensplus.flatsky.rec_lens.qtb(nx, ny, D, rL, fC, T, B, gtype='')

Reconstructing CMB lensing potential and its curl mode from the TB quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

TE cross power spectrum on 2D grid, with bounds (nx,ny)

T[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered temperature, with bounds (nx,ny)

B[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered B-mode, with bounds (nx,ny)

Args(optional):
gtype (str):

Type of output, i.e., convergence (gtype=’k’) or lensing potential (gtype=’’, default)

Returns:
glm[lx,ly] (dcmplx):

2D Fourier modes of CMB lensing potential, with bounds (nx,ny)

clm[lx,ly] (dcmplx):

2D Fourier modes of Curl mode (pseudo lensing potential), with bounds (nx,ny)

Usage:
glm,clm = flatsky.rec_lens.qtb(nx,ny,D,rL,fC,T,B,gtype):

cmblensplus.flatsky.rec_lens.qee(nx, ny, D, rL, fC, E1, E2, gtype='')

Reconstructing CMB lensing potential and its curl mode from the EE quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

EE power spectrum on 2D grid, with bounds (nx,ny)

E1[lx,ly] (dcmplx):

2D Fourier modes of 1st inverse-variance filtered E-mode, with bounds (nx,ny)

E2[lx,ly] (dcmplx):

2D Fourier modes of 2nd inverse-variance filtered E-mode, with bounds (nx,ny)

Args(optional):
gtype (str):

Type of output, i.e., convergence (gtype=’k’) or lensing potential (gtype=’’, default)

Returns:
glm[lx,ly] (dcmplx):

2D Fourier modes of CMB lensing potential, with bounds (nx,ny)

clm[lx,ly] (dcmplx):

2D Fourier modes of Curl mode (pseudo lensing potential), with bounds (nx,ny)

Usage:
glm,clm = flatsky.rec_lens.qee(nx,ny,D,rL,fC,E1,E2,gtype):

cmblensplus.flatsky.rec_lens.qeb(nx, ny, D, rL, fC, E, B, gtype='')

Reconstructing CMB lensing potential and its curl mode from the EB quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

EE power spectrum on 2D grid, with bounds (nx,ny)

E[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered E-mode, with bounds (nx,ny)

B[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered B-mode, with bounds (nx,ny)

Args(optional):
gtype (str):

Type of output, i.e., convergence (gtype=’k’) or lensing potential (gtype=’’, default)

Returns:
glm[lx,ly] (dcmplx):

2D Fourier modes of CMB lensing potential, with bounds (nx,ny)

clm[lx,ly] (dcmplx):

2D Fourier modes of Curl mode (pseudo lensing potential), with bounds (nx,ny)

Usage:
glm,clm = flatsky.rec_lens.qeb(nx,ny,D,rL,fC,E,B,gtype):

cmblensplus.flatsky.rec_lens.qbb(nx, ny, D, rL, fC, B1, B2, gtype='')

Reconstructing CMB lensing potential and its curl mode from the BB quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

BB power spectrum on 2D grid, with bounds (nx,ny)

B1[lx,ly] (dcmplx):

2D Fourier modes of 1st inverse-variance filtered B-mode, with bounds (nx,ny)

B2[lx,ly] (dcmplx):

2D Fourier modes of 2nd inverse-variance filtered B-mode, with bounds (nx,ny)

Args(optional):
gtype (str):

Type of output, i.e., convergence (gtype=’k’) or lensing potential (gtype=’’, default)

Returns:
glm[lx,ly] (dcmplx):

2D Fourier modes of CMB lensing potential, with bounds (nx,ny)

clm[lx,ly] (dcmplx):

2D Fourier modes of Curl mode (pseudo lensing potential), with bounds (nx,ny)

Usage:
glm,clm = flatsky.rec_lens.qbb(nx,ny,D,rL,fC,B1,B2,gtype):

rec_rot – quadratic pol. rot. reconstruction

cmblensplus.flatsky.rec_rot.qte(nx, ny, D, rL, fC, T, E)

Reconstructing anisotropic pol. rot. angles from TE quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

TE cross power spectrum on 2D grid, with bounds (nx,ny)

T[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered temperature, with bounds (nx,ny)

E[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered E-mode, with bounds (nx,ny)

Returns:
alm[lx,ly] (dcmplx):

2D Fourier modes of anisotropic pol. rot. angles, with bounds (nx,ny)

Usage:
alm = flatsky.rec_rot.qte(nx,ny,D,rL,fC,T,E):

cmblensplus.flatsky.rec_rot.qtb(nx, ny, D, rL, fC, T, B)

Reconstructing anisotropic pol. rot. angles from TB quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

TE cross power spectrum on 2D grid, with bounds (nx,ny)

T[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered temperature, with bounds (nx,ny)

B[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered B-mode, with bounds (nx,ny)

Returns:
alm[lx,ly] (dcmplx):

2D Fourier modes of anisotropic pol. rot. angles, with bounds (nx,ny)

Usage:
alm = flatsky.rec_rot.qtb(nx,ny,D,rL,fC,T,B):

cmblensplus.flatsky.rec_rot.qee(nx, ny, D, rL, fC, E1, E2)

Reconstructing anisotropic pol. rot. angles from EE quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

EE power spectrum on 2D grid, with bounds (nx,ny)

E1[lx,ly] (dcmplx):

2D Fourier modes of 1st inverse-variance filtered E-mode, with bounds (nx,ny)

E2[lx,ly] (dcmplx):

2D Fourier modes of 2nd inverse-variance filtered E-mode, with bounds (nx,ny)

Returns:
alm[lx,ly] (dcmplx):

2D Fourier modes of anisotropic pol. rot. angles, with bounds (nx,ny)

Usage:
alm = flatsky.rec_rot.qee(nx,ny,D,rL,fC,E1,E2):

cmblensplus.flatsky.rec_rot.qeb(nx, ny, D, rL, EE, E, B, BB=0)

Reconstructing anisotropic pol. rot. angles from EB quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

EE power spectrum on 2D grid, with bounds (nx,ny)

E[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered E-mode, with bounds (nx,ny)

B[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered B-mode, with bounds (nx,ny)

Args(Optional):
BB[lx,ly] (double):

Theory B-mode spectrum on 2D grid, with bounds (nx,ny), default to BB=0

Returns:
alm[lx,ly] (dcmplx):

2D Fourier modes of anisotropic pol. rot. angles, with bounds (nx,ny)

Usage:
alm = flatsky.rec_rot.qeb(nx,ny,D,rL,EE,E,B,BB):

rec_tau – quadratic patchy tau reconstruction

cmblensplus.flatsky.rec_tau.qtt(nx, ny, D, rL, fC, T1, T2)

Reconstructing patchy tau from the temperature quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

Temperature power spectrum on 2D grid, with bounds (nx,ny)

T1[lx,ly] (dcmplx):

2D Fourier modes of 1st inverse-variance filtered temperature, with bounds (nx,ny)

T2[lx,ly] (dcmplx):

2D Fourier modes of 2nd inverse-variance filtered temperature, with bounds (nx,ny)

Returns:
tlm[lx,ly] (dcmplx):

2D Fourier modes of patchy tau, with bounds (nx,ny)

Usage:
tlm = flatsky.rec_tau.qtt(nx,ny,D,rL,fC,T1,T2):

cmblensplus.flatsky.rec_tau.qte(nx, ny, D, rL, fC, T, E)

Reconstructing patchy tau from the suboptimal TE quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

TE cross power spectrum on 2D grid, with bounds (nx,ny)

T[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered temperature, with bounds (nx,ny)

E[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered E-mode, with bounds (nx,ny)

Returns:
tlm[lx,ly] (dcmplx):

2D Fourier modes of patchy tau, with bounds (nx,ny)

Usage:
tlm = flatsky.rec_tau.qte(nx,ny,D,rL,fC,T,E):

cmblensplus.flatsky.rec_tau.qtb(nx, ny, D, rL, fC, T, B)

Reconstructing patchy tau from the TB quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

TE cross power spectrum on 2D grid, with bounds (nx,ny)

T[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered temperature, with bounds (nx,ny)

B[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered B-mode, with bounds (nx,ny)

Returns:
tlm[lx,ly] (dcmplx):

2D Fourier modes of patchy tau, with bounds (nx,ny)

Usage:
tlm = flatsky.rec_tau.qtb(nx,ny,D,rL,fC,T,B):

cmblensplus.flatsky.rec_tau.qee(nx, ny, D, rL, fC, E1, E2)

Reconstructing patchy tau from the EE quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

EE power spectrum on 2D grid, with bounds (nx,ny)

E1[lx,ly] (dcmplx):

2D Fourier modes of 1st inverse-variance filtered E-mode, with bounds (nx,ny)

E2[lx,ly] (dcmplx):

2D Fourier modes of 2nd inverse-variance filtered E-mode, with bounds (nx,ny)

Returns:
tlm[lx,ly] (dcmplx):

2D Fourier modes of patchy tau, with bounds (nx,ny)

Usage:
tlm = flatsky.rec_tau.qee(nx,ny,D,rL,fC,E1,E2):

cmblensplus.flatsky.rec_tau.qeb(nx, ny, D, rL, fE, fB, E, B)

Reconstructing patchy tau from the EB quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fE[lx,ly] (double):

EE power spectrum on 2D grid, with bounds (nx,ny)

fB[lx,ly] (double):

BB power spectrum on 2D grid, with bounds (nx,ny)

E[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered E-mode, with bounds (nx,ny)

B[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered B-mode, with bounds (nx,ny)

Returns:
tlm[lx,ly] (dcmplx):

2D Fourier modes of patchy tau, with bounds (nx,ny)

Usage:
tlm = flatsky.rec_tau.qeb(nx,ny,D,rL,fE,fB,E,B):

rec_src – quadratic point-src reconstruction

cmblensplus.flatsky.rec_src.qtt(nx, ny, D, rL, T1, T2)

Reconstructing point source fields from the temperature quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

Temperature power spectrum on 2D grid, with bounds (nx,ny)

T1[lx,ly] (dcmplx):

2D Fourier modes of 1st inverse-variance filtered temperature, with bounds (nx,ny)

T2[lx,ly] (dcmplx):

2D Fourier modes of 2nd inverse-variance filtered temperature, with bounds (nx,ny)

Returns:
slm[lx,ly] (dcmplx):

2D Fourier modes of point source fields, with bounds (nx,ny)

Usage:
slm = flatsky.rec_src.qtt(nx,ny,D,rL,T1,T2):

cmblensplus.flatsky.rec_src.qte(nx, ny, D, rL, T, E)

Reconstructing point source fields from the suboptimal TE quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

TE cross power spectrum on 2D grid, with bounds (nx,ny)

T[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered temperature, with bounds (nx,ny)

E[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered E-mode, with bounds (nx,ny)

Returns:
slm[lx,ly] (dcmplx):

2D Fourier modes of point source fields, with bounds (nx,ny)

Usage:
slm = flatsky.rec_src.qte(nx,ny,D,rL,T,E):

cmblensplus.flatsky.rec_src.qtb(nx, ny, D, rL, T, B)

Reconstructing point source fields from the TB quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

TE cross power spectrum on 2D grid, with bounds (nx,ny)

T[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered temperature, with bounds (nx,ny)

B[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered B-mode, with bounds (nx,ny)

Returns:
slm[lx,ly] (dcmplx):

2D Fourier modes of point source fields, with bounds (nx,ny)

Usage:
slm = flatsky.rec_src.qtb(nx,ny,D,rL,T,B):

cmblensplus.flatsky.rec_src.qee(nx, ny, D, rL, E1, E2)

Reconstructing point source fields from the EE quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

EE power spectrum on 2D grid, with bounds (nx,ny)

E1[lx,ly] (dcmplx):

2D Fourier modes of 1st inverse-variance filtered E-mode, with bounds (nx,ny)

E2[lx,ly] (dcmplx):

2D Fourier modes of 2nd inverse-variance filtered E-mode, with bounds (nx,ny)

Returns:
slm[lx,ly] (dcmplx):

2D Fourier modes of point source fields, with bounds (nx,ny)

Usage:
slm = flatsky.rec_src.qee(nx,ny,D,rL,E1,E2):

cmblensplus.flatsky.rec_src.qeb(nx, ny, D, rL, E, B)

Reconstructing point source fields from the EB quadratic estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

fC[lx,ly] (double):

EE power spectrum on 2D grid, with bounds (nx,ny)

E[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered E-mode, with bounds (nx,ny)

B[lx,ly] (dcmplx):

2D Fourier modes of inverse-variance filtered B-mode, with bounds (nx,ny)

Returns:
slm[lx,ly] (dcmplx):

2D Fourier modes of point source fields, with bounds (nx,ny)

Usage:
slm = flatsky.rec_src.qeb(nx,ny,D,rL,E,B):

norm_lens – normalization of quadratic lensing reconstruction

cmblensplus.flatsky.norm_lens.qtt(nx, ny, D, rL, OT, TT, eL)

Normalization of the temperature quadratic estimator for CMB lensing potential and its curl mode

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

OT[lx,ly] (double):

Inverse of Observed temperature power spectrum on 2D grid, with bounds (nx,ny)

TT[lx,ly] (double):

Theory temperature power spectrum on 2D grid, with bounds (nx,ny)

eL[2] (int):

Minimum and maximum multipole of output normalization spectrum, with bounds (2)

Returns:
Ag[lx,ly] (dcmplx):

Normalization of CMB lensing potential on 2D grid, with bounds (nx,ny)

Ac[lx,ly] (dcmplx):

Normalization of Curl mode (pseudo lensing potential) on 2D grid, with bounds (nx,ny)

Usage:
Ag,Ac = flatsky.norm_lens.qtt(nx,ny,D,rL,OT,TT,eL):

cmblensplus.flatsky.norm_lens.qte(nx, ny, D, rL, OT, OE, TE, eL)

Normalization of the TE quadratic estimator for CMB lensing potential and its curl mode

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

OT[lx,ly] (double):

Inverse of Observed temperature power spectrum on 2D grid, with bounds (nx,ny)

OE[lx,ly] (double):

Inverse of Observed E-mode power spectrum on 2D grid, with bounds (nx,ny)

TE[lx,ly] (double):

Theory TE cross spectrum on 2D grid, with bounds (nx,ny)

eL[2] (int):

Minimum and maximum multipole of output normalization spectrum, with bounds (2)

Returns:
Ag[lx,ly] (dcmplx):

Normalization of CMB lensing potential on 2D grid, with bounds (nx,ny)

Ac[lx,ly] (dcmplx):

Normalization of Curl mode (pseudo lensing potential) on 2D grid, with bounds (nx,ny)

Usage:
Ag,Ac = flatsky.norm_lens.qte(nx,ny,D,rL,OT,OE,TE,eL):

cmblensplus.flatsky.norm_lens.qtb(nx, ny, D, OT, OB, TE, rL, eL)

Normalization of the TB quadratic estimator for CMB lensing potential and its curl mode

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

OT[lx,ly] (double):

Inverse of Observed temperature power spectrum on 2D grid, with bounds (nx,ny)

OB[lx,ly] (double):

Inverse of Observed B-mode power spectrum on 2D grid, with bounds (nx,ny)

TE[lx,ly] (double):

Theory TE cross spectrum on 2D grid, with bounds (nx,ny)

eL[2] (int):

Minimum and maximum multipole of output normalization spectrum, with bounds (2)

Returns:
Ag[lx,ly] (dcmplx):

Normalization of CMB lensing potential on 2D grid, with bounds (nx,ny)

Ac[lx,ly] (dcmplx):

Normalization of Curl mode (pseudo lensing potential) on 2D grid, with bounds (nx,ny)

Usage:
Ag,Ac = flatsky.norm_lens.qtb(nx,ny,D,OT,OB,TE,rL,eL):

cmblensplus.flatsky.norm_lens.qee(nx, ny, D, OE, EE, rL, eL)

Normalization of the EE quadratic estimator for CMB lensing potential and its curl mode

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

OE[lx,ly] (double):

Inverse of Observed E-mode power spectrum on 2D grid, with bounds (nx,ny)

EE[lx,ly] (double):

Theory E-mode spectrum on 2D grid, with bounds (nx,ny)

eL[2] (int):

Minimum and maximum multipole of output normalization spectrum, with bounds (2)

Returns:
Ag[lx,ly] (dcmplx):

Normalization of CMB lensing potential on 2D grid, with bounds (nx,ny)

Ac[lx,ly] (dcmplx):

Normalization of Curl mode (pseudo lensing potential) on 2D grid, with bounds (nx,ny)

Usage:
Ag,Ac = flatsky.norm_lens.qee(nx,ny,D,OE,EE,rL,eL):

cmblensplus.flatsky.norm_lens.qeb(nx, ny, D, OE, OB, EE, rL, eL)

Normalization of the EB quadratic estimator for CMB lensing potential and its curl mode

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

OE[lx,ly] (double):

Inverse of Observed E-mode power spectrum on 2D grid, with bounds (nx,ny)

OB[lx,ly] (double):

Inverse of Observed B-mode power spectrum on 2D grid, with bounds (nx,ny)

EE[lx,ly] (double):

Theory E-mode spectrum on 2D grid, with bounds (nx,ny)

eL[2] (int):

Minimum and maximum multipole of output normalization spectrum, with bounds (2)

Returns:
Ag[lx,ly] (dcmplx):

Normalization of CMB lensing potential on 2D grid, with bounds (nx,ny)

Ac[lx,ly] (dcmplx):

Normalization of Curl mode (pseudo lensing potential) on 2D grid, with bounds (nx,ny)

Usage:
Ag,Ac = flatsky.norm_lens.qeb(nx,ny,D,OE,OB,EE,rL,eL):

norm_rot – normalization of quadratic pol. rot. reconstruction

cmblensplus.flatsky.norm_rot.qeb(nx, ny, D, rL, IE, IB, EE, eL, BB=0)

Normalization of the EB quadratic estimator for anisotropic pol. rot. angles

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

IE[lx,ly] (double):

Inverse of observed E-mode power spectrum on 2D grid, with bounds (nx,ny)

IB[lx,ly] (double):

Inverse of observed B-mode power spectrum on 2D grid, with bounds (nx,ny)

EE[lx,ly] (double):

Theory E-mode spectrum on 2D grid, with bounds (nx,ny)

eL[2] (int):

Minimum and maximum multipole of output normalization spectrum, with bounds (2)

Args(Optional):
BB[lx,ly] (double):

Theory B-mode spectrum on 2D grid, with bounds (nx,ny), default to BB=0

Returns:
Aa[lx,ly] (dcmplx):

Normalization of anisotropic pol. rot. angles on 2D grid, with bounds (nx,ny)

Usage:
Aa = flatsky.norm_rot.qeb(nx,ny,D,rL,IE,IB,EE,eL,BB):

norm_tau – normalization of quadratic patchy tay reconstruction

cmblensplus.flatsky.norm_tau.qtt(nx, ny, D, rL, OT, TT, eL)

Normalization of the temperature quadratic estimator for patchy tau

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

OT[lx,ly] (double):

Inverse of Observed temperature power spectrum on 2D grid, with bounds (nx,ny)

TT[lx,ly] (double):

Theory temperature power spectrum on 2D grid, with bounds (nx,ny)

eL[2] (int):

Minimum and maximum multipole of output normalization spectrum, with bounds (2)

Returns:
At[lx,ly] (dcmplx):

Normalization of patchy tau on 2D grid, with bounds (nx,ny)

Usage:
At = flatsky.norm_tau.qtt(nx,ny,D,rL,OT,TT,eL):

cmblensplus.flatsky.norm_tau.qeb(nx, ny, D, rL, IE, IB, EE, eL)

Normalization of the EB quadratic estimator for patchy tau

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

rL[2] (int):

Minimum and maximum multipole of CMB for reconstruction

IE[lx,ly] (double):

Inverse of observed E-mode power spectrum on 2D grid, with bounds (nx,ny)

IB[lx,ly] (double):

Inverse of observed B-mode power spectrum on 2D grid, with bounds (nx,ny)

EE[lx,ly] (double):

Theory E-mode spectrum on 2D grid, with bounds (nx,ny)

eL[2] (int):

Minimum and maximum multipole of output normalization spectrum, with bounds (2)

Returns:
At[lx,ly] (dcmplx):

Normalization of patchy tau on 2D grid, with bounds (nx,ny)

Usage:
At = flatsky.norm_tau.qeb(nx,ny,D,rL,IE,IB,EE,eL):

bispec – bispectrum tools

cmblensplus.flatsky.bispec.bispec_norm(nx, ny, D, bp, dbin_max=-1, bn=1)
Usage:
norm = flatsky.bispec.bispec_norm(nx,ny,D,bp,dbin_max,bn):

cmblensplus.flatsky.bispec.bispec_bin(kmap, bp, kn=1, bn=1, nx=0, ny=0, dbin_max=-1)
Usage:
bispec = flatsky.bispec.bispec_bin(kn,bn,nx,ny,kmap,bp,dbin_max):

cmblensplus.flatsky.bispec.binfilter(nx, ny, D, bp, bn=1)

The multipole bin binary-mask given on the 2D Fourier grids so that M = 1 inside the multipole bin and 0 otherwise

Args:
nx, ny (int):

Number of Lx and Ly grids

D[2] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

bp[bn+1] (double):

Multipole bin edges

Returns:
bf[bn,nx,ny] (double):

The multipole bin binary-mask for each multipol bin

Usage:
bf = flatsky.bispec.binfilter(nx,ny,D,bp,bn):

cmblensplus.flatsky.bispec.bispec_norm_1d(nx, ny, D, bfs, bn=1)

Normalization of the 1D binned bispectrum estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[2] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

bfs[3,bn,nx,ny] (double):

Multipole bin mask on 2D grids obtained from the binfilter function

Returns:
bnorm[bn] (double):

Normalization of 1D binned bispectrum at each multipole bin

Usage:
bnorm = flatsky.bispec.bispec_norm_1d(nx,ny,D,bfs,bn):

cmblensplus.flatsky.bispec.bispec_bin_1d(nx, ny, D, bfs, bnorm, alm, bn=1)

1D binned bispectrum estimator

Args:
nx, ny (int):

Number of Lx and Ly grids

D[2] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

bfs[3,bn,nx,ny] (double):

Multipole bin mask on 2D grids obtained from the binfilter function

bnorm[bn] (double):

Normalization of 1D binned bispectrum at each multipole bin

alm[3,nx,ny] (dcmplx):

Fourier modes for each leg

Returns:
bispec[bn] (double):

1D binned bispectrum at each multipole bin

Usage:
bispec = flatsky.bispec.bispec_bin_1d(nx,ny,D,bfs,bnorm,alm,bn):

utils – other tools

cmblensplus.flatsky.utils.map2alm(nx, ny, D, map)

DFT for 2D array.

Args:
nx, ny (int):

Number of x and y grids

D[2] (double):

Side length (x and y) of map

map[x,y] (double):

Map on 2D grid with bounds (nx,ny)

Returns:
alm[x,y] (dcmplx):

Fourier modes on 2D grid, with bounds (nx,ny)

Usage:
alm = flatsky.utils.map2alm(nx,ny,D,map):

cmblensplus.flatsky.utils.alm2map(nx, ny, D, alm)

DFT for 2D array.

Args:
nx, ny (int):

Number of Lx and Ly grids

D[2] (double):

Side length (x and y) of map

alm[x,y] (dcmplx):

Fourier modes on 2D grid to be transformed, with bounds (nx,ny)

Returns:
map[x,y] (double):

Map on 2D grid, with bounds (nx,ny)

Usage:
map = flatsky.utils.alm2map(nx,ny,D,alm):

cmblensplus.flatsky.utils.el2d(nx, ny, D)

Return absolute value of multipole in 2D grids

Args:
nx, ny (int):

number of Lx and Ly grids

D[xy] (double):

map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

Returns:
els[nx,ny] (double):

absolute value of Fourier mode, (Lx**2+Ly**2)**0.5, with bounds (nx,ny)

Usage:
els = flatsky.utils.el2d(nx,ny,D):

cmblensplus.flatsky.utils.elarrays(nx, ny, D)

Return Lx, Ly, absolute value of multipole, and its inverse in 2D grids

Args:
nx, ny (int):

number of Lx and Ly grids

D[xy] (double):

map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

Returns:
elx[nx,ny] (double):

Lx, with bounds (nx,ny)

ely[nx,ny] (double):

Ly, with bounds (nx,ny)

els[nx,ny] (double):

absolute value of Fourier mode, (Lx**2+Ly**2)**0.5, with bounds (nx,ny)

eli[nx,ny] (double):

inverse of els, with bounds (nx,ny)

Usage:
elx,ely,els,eli = flatsky.utils.elarrays(nx,ny,D):

cmblensplus.flatsky.utils.elmask(nx, ny, D, lmin=0, lmax=1000, lxcut=0, lycut=0)

Return mask in 2D Fourier space. The lmask is unity at lmin<=|L|<=lmax, |Lx|>=lxcut, |Ly|>=lycut, and otherwize zero.

Args:
nx, ny (int):

number of Lx and Ly grids

D[xy] (double):

map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

Args(optional):
lmin/lmax (int):

Minimum/Maximum of multipoles

lxcut/lycut (int):

Remove |Lx|<lxcut / |Ly|<lycut cut of multipoles

Returns:
lmask[nx,ny] (double):

Mask, with bounds (nx,ny)

Usage:
lmask = flatsky.utils.elmask(nx,ny,D,lmin,lmax,lxcut,lycut):

cmblensplus.flatsky.utils.alm2bcl(bn, oL, nx, ny, D, alm1, alm2=None, spc='')

Compute angular power spectrum from Fourier modes, with multipole binning

Args:
bn (int):

number of multipole bin

oL[2] (int):

minimum and maximum multipoles of the output cl

nx, ny (int):

number of Lx and Ly grids

D[xy] (double):

map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

alm1[nx,ny] (dcmplx):

Fourier mode, with bounds (nx,ny)

Args(optional):
alm2[nx,ny] (dcmplx):

Fourier mode, with bounds (nx,ny), default to None

spc (str):

type of multipole binning, i.e., linear spacing (spc=’’, default), or log spacing (spc=’log’)

Returns:
Cb[bin] (double):

angular power spectrum with multipole binning, with bounds (bn)

Usage:
Cb = flatsky.utils.alm2bcl(bn,oL,nx,ny,D,alm1,alm2,spc):

cmblensplus.flatsky.utils.c2d2bcl(nx, ny, D, c2d, bn, oL, spc='')

Return 1D angular power spectrum with multipole binning from a 2D power spectrum

Args:
nx, ny (int):

number of Lx and Ly grids

D[xy] (double):

map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

c2d[nx,ny] (double):

2D power spectrum, with bounds (nx,ny)

bn (int):

number of multipole bin

oL[2] (int):

minimum and maximum multipoles of the output cl

Args(optional):
spc (str):

type of multipole binning, i.e., linear spacing (spc=’’, default), or log spacing (spc=’log’)

Returns:
Cb[bin] (double):

angular power spectrum with multipole binning, with bounds (bn)

Usage:
Cb = flatsky.utils.c2d2bcl(nx,ny,D,c2d,bn,oL,spc):

cmblensplus.flatsky.utils.cl2c2d(nx, ny, D, lmin, lmax, Cl, method='linear')

Assign values of 1D angular power spectrum on to 2D grid with linear interpolation

Args:
nx, ny (int):

number of Lx and Ly grids

D[xy] (double):

map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

lmin (int):

minimum multipole of cl to be interpolated

lmax (int):

maximum multipole of cl to be interpolated

Cl[l] (double):

1D power spectrum, with bounds (0:lmax)

Args(optional):
method (str):

type of interpolation method, i.e., linear interpolation (method=’linear’, default), or step (method=’step’)

Returns:
c2d[nx,ny] (double):

2D power spectrum, with bounds (nx,ny)

Usage:
c2d = flatsky.utils.cl2c2d(nx,ny,D,lmin,lmax,Cl,method):

cmblensplus.flatsky.utils.cb2c2d(bn, bc, nx, ny, D, lmin, lmax, Cb, method='')

Assign values of 1D angular power spectrum on to 2D grid with linear interpolation

Args:
bn (int):

number of multipole bins

bc[bin] (double):

multipole bin center, with bounds (bn)

nx, ny (int):

number of Lx and Ly grids

D[xy] (double):

map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

lmin (int):

minimum multipole of cl to be interpolated

lmax (int):

maximum multipole of cl to be interpolated

Cb[bin] (double):

1D power spectrum with multipole binning, with bounds (bn)

Args(optional):
method (str):

interpolation method from binned to unbinned angular spectrum, i.e., spline (‘’, default), or linear (‘linear’) interpolation

Returns:
c2d[nx,ny] (double):

2D power spectrum, with bounds (nx,ny)

Usage:
C2d = flatsky.utils.cb2c2d(bn,bc,nx,ny,D,lmin,lmax,Cb,method):

cmblensplus.flatsky.utils.gauss1alm(nx, ny, D, lmin, lmax, Cl)

Generate random gaussian fields in 2D Fourier space for a given isotropic spectrum, satisfying a^*_l = a_{-l}

Args:
nx, ny (int):

number of Lx and Ly grids

D[xy] (double):

map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

lmin (int):

minimum multipole of cl to be interpolated

lmax (int):

maximum multipole of cl to be interpolated

Cl[l] (double):

1D power spectrum, with bounds (0:lmax)

Returns:
alm[lx,ly] (dcmplx):

random gaussian fields on 2D Fourier plane, with bounds (nx,ny)

Usage:
alm = flatsky.utils.gauss1alm(nx,ny,D,lmin,lmax,Cl):

cmblensplus.flatsky.utils.gauss2alm(nx, ny, D, lmin, lmax, TT, TE, EE)

Generate two correlated random gaussian fields in 2D Fourier space for a given isotropic spectrum

Args:
nx, ny (int):

number of Lx and Ly grids

D[xy] (double):

map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

lmin (int):

minimum multipole of cl to be interpolated

lmax (int):

maximum multipole of cl to be interpolated

TT[l] (double):

the 1st 1D power spectrum, with bounds (0:lmax)

TE[l] (double):

the cross 1D power spectrum, with bounds (0:lmax)

EE[l] (double):

the 2nd 1D power spectrum, with bounds (0:lmax)

Returns:
tlm[lx,ly] (dcmplx):

the 1st random gaussian fields on 2D Fourier plane, with bounds (nx,ny)

elm[lx,ly] (dcmplx):

the 2nd random gaussian fields on 2D Fourier plane, with bounds (nx,ny)

Usage:
tlm,elm = flatsky.utils.gauss2alm(nx,ny,D,lmin,lmax,TT,TE,EE):

cmblensplus.flatsky.utils.window_sin(nx, ny, D, ap=1, cut=1)

Return a sin window function.

Args:
nx, ny (int):

Number of Lx and Ly grids

D[xy] (double):

Map side length, or equivalent to dLx/2pi, dLy/2pi, with bounds (2)

Args(Optional):
ap (double):

Apodization parameter defined by apodized-range = (1-ap) x (cut)mapsize, from 0 (full apodization) to 1 (no apodization). Default to 1.

cut (double):

Map cut scale defined by cutmapsize = cut x mapsize, from 0 (full cut) to 1 (no cut). Default to 1.

Return:
W[x,y] (double):

Window function, with bounds (nx,ny)

Usage:
W = flatsky.utils.window_sin(nx,ny,D,ap,cut):

cmblensplus.flatsky.utils.window_norm(nx, ny, wind, num)
Usage:
wn = flatsky.utils.window_norm(nx,ny,wind,num):

cmblensplus.flatsky.utils.window_norm_x(nx, ny, W1, W2, num)
Usage:
Wn = flatsky.utils.window_norm_x(nx,ny,W1,W2,num):

cmblensplus.flatsky.utils.rotation(nx, ny, rot, QU, rtype)
Usage:
rQU = flatsky.utils.rotation(nx,ny,rot,QU,rtype):

cmblensplus.flatsky.utils.get_angle(nx, ny, D)
Usage:
theta,phi = flatsky.utils.get_angle(nx,ny,D):

cmblensplus.flatsky.utils.cutmap(ox, oy, cx, cy, omap)
Usage:
cmap = flatsky.utils.cutmap(ox,oy,cx,cy,omap):