basic
– Some basic tools¶
delens
– analytic calculation for delensing¶
-
basic.delens.
lintemplate
(lmax, elmin, elmax, klmin, klmax, CE, Cm, WE, Wm, gtype='p')¶ Estimate of lensing template B-mode power spectrum (Wiener filters as inputs)
- Args:
lmax (int): Maximum multipole of output spectrum elmin (int): Minimum multipole of E elmax (int): Maximum multipole of E klmin (int): Minimum multipole of lensing mass klmax (int): Maximum multipole of lensing mass CE[l] (double): Power spectrum of E-mode, with bounds (0:dlmax) Cp[l] (double): Power spectrum of lensing pontential, with bounds (0:dlmax) WE[l] (double): Wiener filter of E-mode, with bounds (0:dlmax) Wp[l] (double): Wiener filter of lensing potential, with bountd (0:dlmax) - Args(optional):
gtype (str): specify type of the input Cp, p (default) or k. - Returns:
CB[l] (double): Lensing B-mode power spectrum, with bounds (0:lmax) - Usage:
CB = basic.delens.lintemplate(lmax,elmin,elmax,klmin,klmax,CE,Cm,WE,Wm,gtype):
-
basic.delens.
lensingbb
(lmax, dlmin, dlmax, CE, Cp)¶ Lensing B-mode power spectrum as a convolution of ClEE and Clpp
- Args:
lmax (int): Maximum multipole of output spectrum dlmin (int): Minimum multipole of E and lensing for delensing dlmax (int): Maximum multipole of E and lensing for delensing CE[l] (double): Power spectrum of E-mode, with bounds (0:dlmax) Cp[l] (double): Power spectrum of lensing pontential, with bounds (0:dlmax) - Returns:
CB[l] (double): Lensing B-mode power spectrum, with bounds (0:lmax) - Usage:
CB = basic.delens.lensingbb(lmax,dlmin,dlmax,CE,Cp):
-
basic.delens.
delensbias_dom
(lmax, dlmin, dlmax, CE, CB, Cp, NP1, NP2, Ag)¶ Dominant term of the delensing bias in the B-mode internal delensing
- Args:
lmax (int): Maximum multipole of output spectrum dlmin (int): Minimum multipole of E and lensing for delensing dlmax (int): Maximum multipole of E and lensing for delensing CE[l] (double): Power spectrum of E-mode, with bounds (0:dlmax) CB[l] (double): Power spectrum of B-mode, with bounds (0:dlmax) Cp[l] (double): Power spectrum of lensing pontential, with bounds (0:dlmax) NP1[l] (double): Pol. noise spectrum for lensing reconstruction, with bounds (0:dlmax) NP2[l] (double): Pol. noise spectrum for B-mode to be delensed, with bounds (0:dlmax) Ag[l] (double): Lensing reconstruction noise, with bounds (0:dlmax) - Returns:
DB[l] (double): Lensing B-mode power spectrum, with bounds (0:lmax) - Usage:
DB = basic.delens.delensbias_dom(lmax,dlmin,dlmax,CE,CB,Cp,NP1,NP2,Ag):
bispec
– analytic calculation of bispectrum¶
-
basic.bispec.
cl_flat
(cpmodel, z, dz, zs, lmax, k, pk0, zn=0, kn=0, pktype='T12', cltype='kk', dNdz=None, wdel=None)¶ Compute flat-sky lensing powerspectrum analytically
- Args:
cpmodel (str): cosmological parameter model (model0, modelw, or modelp) z[zn] (double): redshift points for the z-integral dz[zn] (double): interval of z zs[2] (double): two source redshifts lmin/lmax (int): minimum/maximum multipoles of the bispectrum k[kn] (double): k for the matter power spectrum in unit of [h/Mpc] pk0[kn] (double): the linear matter power spectrum at z=0 in unit of [Mpc^3/h^3] - Args(optional):
pktype (str): fitting formula for the matter power spectrum (Lin, S02 or T12) cltype (str): kk, gk, or gg dNdz[zn] (double): redshift distribution of galaxy, only used when cltype includes g wdel[zn,l] (double): modified chi-kernel function for z-cleaning at l=0 to lmax - Returns:
cl[l] (double): power spectrum from LSS contributions at [lmin,lmax] - Usage:
cl = basic.bispec.cl_flat(cpmodel,z,dz,zs,lmax,k,pk0,zn,kn,pktype,cltype,dNdz,wdel):
-
basic.bispec.
bispeclens
(shap, cpmodel, model, z, dz, zs, lmin, lmax, k, pk0, zn=0, kn=0, lan=0.0, kan=0.0, pktype='T12', ltype='', btype='kkk', dNdz=None, wdel=None)¶ Compute lensing bispectrum analytically
- Args:
shap (str): shape of the bispectrum (equi, fold, sque, or isos) cpmodel (str): cosmological parameter model (model0, modelw, or modelp) model (str): fitting formula of the matter bispectrum (LN=linear, SC=SC03, GM=Gil-Marin+12, 3B=3-shape-bispectrum, or RT=Takahashi+19) z[zn] (double): redshift points for the z-integral dz[zn] (double): interval of z zs[3] (double): source redshifts lmin/lmax (int): minimum/maximum multipoles of the bispectrum k[kn] (double): k for the matter power spectrum in unit of [h/Mpc] pk0[kn] (double): the linear matter power spectrum at z=0 in unit of [Mpc^3/h^3] - Args(optional):
lan, kan (double): parameters for the modified gravity extension, default to lan=kan=1 (GR) pktype (str): fitting formula for the matter power spectrum (Lin=Linear, S02=Smith et al. 2002 or T12=Takahashi et al. 2012), default to T12 ltype (str): fullsky correction (full) or not btype (str): bispectrum type, i.e., kkk (lens-lens-lens), gkk (density-lens-lens), ggk (density-density-lens), or ggg (density-density-density) dNdz[zn] (double): redshift distribution of galaxy, only used when btype includes g wdel[zn,l] (double): modified chi-kernel function by z-cleaning at l=0 to lmax - Returns:
bl0[l] (double): lensing bispectrum from LSS contributions at [lmin,lmax] bl1[l] (double): lensing bispectrum from post-Born contributions at [lmin,lmax] - Usage:
bl0,bl1 = basic.bispec.bispeclens(shap,cpmodel,model,z,dz,zs,lmin,lmax,k,pk0,lan,kan,zn,kn,pktype,ltype,btype,dNdz,wdel):
-
basic.bispec.
bispeclens_bin
(shap, cpmodel, model, z, dz, zs, lmin, lmax, bn, k, pk0, zn=0, kn=0, lan=0.0, kan=0.0, pktype='T12', btype='kkk', dNdz=None, wdel=None)¶ Compute binned lensing bispectrum analytically
- Args:
shap (str): shape of the bispectrum (equi, fold, sque, or isos) cpmodel (str): cosmological parameter model (model0, modelw, or modelp) model (str): fitting formula of the matter bispectrum (LN=linear, SC=SC03, GM=Gil-Marin+12, 3B=3-shape-bispectrum, or RT=Takahashi+19) z[zn] (double): redshift points for the z-integral dz[zn] (double): interval of z zs[3] (double): source redshifts lmin/lmax (int): minimum/maximum multipoles of the bispectrum bn (int): number of multipole bins k[kn] (double): k for the matter power spectrum in unit of [h/Mpc] pk0[kn] (double): the linear matter power spectrum at z=0 in unit of [Mpc^3/h^3] - Args(optional):
lan, kan (double): parameters for the modified gravity extension, default to lan=kan=1 (GR) pktype (str): fitting formula for the matter power spectrum (Lin, S02 or T12) btype (str): bispectrum type, i.e., kkk (lens-lens-lens), gkk (density-lens-lens), ggk (density-density-lens), or ggg (density-density-density) dNdz[zn] (double): redshift distribution of galaxy, only used when btype includes g wdel[zn,l] (double): modified chi-kernel function by z-cleaning at l=0 to lmax - Returns:
bc[bn] (double): multipole bin centers bl0[bn] (double): binned lensing bispectrum from LSS contributions bl1[bn] (double): binned lensing bispectrum from post-Born contributions - Usage:
bc,bl0,bl1 = basic.bispec.bispeclens_bin(shap,cpmodel,model,z,dz,zn,zs,lmin,lmax,bn,k,pk0,kn,lan,kan,pktype,btype,dNdz,wdel):
-
basic.bispec.
bispeclens_snr
(cpmodel, model, z, dz, zs, lmin, lmax, cl, k, pk0, zn=0, kn=0, pktype='T12', btype='kkk', dNdz=None, wdel=None, cgg=None, ro=100)¶ Compute SNR of lensing bispectrum analytically
- Args:
cpmodel (str): cosmological parameter model (model0, modelw, or modelp) model (str): fitting formula of the matter bispectrum (LN=linear, SC=SC03, GM=Gil-Marin+12, 3B=3-shape-bispectrum, or RT=Takahashi+19) z[zn] (double): redshift points for the z-integral dz[zn] (double): interval of z zs[3] (double): source redshifts lmin/lmax (int): minimum/maximum multipoles of the bispectrum cl[l] (int): observed angular power spectrum at 0<=l<=lmax k[kn] (double): k for the matter power spectrum in unit of [h/Mpc] pk0[kn] (double): the linear matter power spectrum at z=0 in unit of [Mpc^3/h^3] - Args(optional):
pktype (str): fitting formula for the matter power spectrum (Lin, S02 or T12) btype (str): bispectrum type, i.e., kkk (lens-lens-lens), gkk (density-lens-lens), ggk (density-density-lens), or ggg (density-density-density) dNdz[zn] (double): redshift distribution of galaxy, only used when btype includes g wdel[zn,l] (double): modified chi-kernel function by z-cleaning at l=0 to lmax cgg[l] (double): observed galaxy spectrum ro (int): output progress for every “ro” multipoles (ro=100, default) - Returns:
snr[2] (double): total SNR amd LSS-only SNR - Usage:
snr = basic.bispec.bispeclens_snr(cpmodel,model,z,dz,zn,zs,lmin,lmax,cl,k,pk0,kn,pktype,btype,dNdz,cgg,ro,wdel):
-
basic.bispec.
bispeclens_gauss_bin
(shap, bn, lmin, lmax, cl)¶ Compute binned bispectrum analytically for the quadratic gaussian model
- Args:
shap (str): shape of the bispectrum (equi, fold, sque, or isos) bn (int): number of multipole bins lmin/lmax (int): minimum/maximum multipoles of the bispectrum cl[l] (double): the power spectrum at [0:lmax+1] - Returns:
bc[bn] (double): multipole bin centers bl[bn] (double): binned bispectrum - Usage:
bc,bl = basic.bispec.bispeclens_gauss_bin(shap,bn,lmin,lmax,cl):
-
basic.bispec.
zpoints
(zmin, zmax, zn, zspace=1)¶ Precomputing interpolation points for z
- Args:
zmin/zmax (double): minimum/maximum redshifts zn (int): number of redshifts - Args(optional):
zspace (int): type of spacing - Returns:
z[zn] (double): redshifts dz[zn] (double): redshift intervals - Usage:
z,dz = basic.bispec.zpoints(zmin,zmax,zn,zspace):
-
basic.bispec.
skewspeclens
(cpmodel, model, zmin, zmax, zn, zs, ols, lmin, lmax, k, pk0, pktype='T12', pb=True, theta=0.0, Om=0.3, H0=70.0, w0=-1.0, wa=0.0, mnu=0.06, ns=0.965, bn=0, kn=0, verbose=True, wdel=None)¶ Compute skew spectrum using a matter bispectrum fitting formula
- Args:
cpmodel (str): cosmological parameter model (model0, modelw, modelp, or input) model (str): fitting formula of the matter bispectrum (LN=linear, SC=SC03, GM=Gil-Marin+12, 3B=3-shape-bispectrum, or RT=Takahashi+19) zmin/zmax (double): minimum/maximum z for z-integral zn (int): number of redshifts for the z-integral zs[2] (double): source redshifts where zs[2] is used for the squared map lmin/lmax (int): minimum/maximum multipoles of alms included in the skew spectrum ols[bn] (int): output multipoles to be computed for skew spectrum k[kn] (double): k for the matter power spectrum [h/Mpc] pk0 (double): the linear matter power spectrum at z=0 [Mpc^3/h^3] - Args(optional):
pktype (str): fitting formula for the matter power spectrum (Lin, S02 or T12) theta (double): kappa map resolution in arcmin pb (bool): with post-Born correction or not (default=True) verbose (bool): output messages wdel[zn,l] (double): modified chi-kernel function by z-cleaning at l=0 to lmax - Returns:
skew (double): skew-spectrum - Usage:
skew = basic.bispec.skewspeclens(cpmodel,model,zmin,zmax,zn,zs,bn,ols,lmin,lmax,k,pk0,kn,theta,pktype,pb,Om,H0,w0,wa,mnu,ns,verbose,wdel):
galaxy
– tools for galaxy z distribution¶
-
basic.galaxy.
dndz_sf
(z, a, b, zm=0, z0=0, zn=None)¶ A model of galaxy z distribution
- Args:
z[zn] (double): redshifts at which dNdz is returned a, b (double): shape parameters of Schechter-like galaxy distribution - Args(optional):
zm (double): mean redshift, default to 0 z0 (double): a parameter relevant to zm, default to 0. Either zm or z0 has to be specified. - Returns:
dndz[zn] (double): galaxy z distribution - Usage:
dndz = basic.galaxy.dndz_sf(zn,z,a,b,zm,z0):
-
basic.galaxy.
photoz_error
(z, zi, zn=None, sigma=0.03, zbias=0.0)¶ Photo-z error on z distribution which is multiplied to original galaxy z distribution. See Eq.(13) of arXiv:1103.1118 for details.
- Args:
z[zn] (double): redshifts at which photoz error function is returned zi[2] (double): z-bin edges sigma (double): a parameter of photo-z error which is given by, sigma x (1+z) zbias (double): photo-z mean bias - Returns:
pz[zn] (double): photoz error function - Usage:
pz = basic.galaxy.photoz_error(zn,z,zi,sigma,zbias):
-
basic.galaxy.
zbin
(zn, a, b, zm=0, z0=0, verbose=False)¶ Computing z-interval of z-bin so that number of galaxies at each z-bin is equal
- Args:
zn (int): number of z-bins a, b (double): shape parameters of Schechter-like galaxy distribution - Args(optional):
zm (double): mean redshift, default to 0 z0 (double): a parameter relevant to zm, default to 0. Either zm or z0 has to be specified. verbose (bool): output messages (default to True) - Returns:
zb[zn+1] (double): z-intervals - Usage:
zb = basic.galaxy.zbin(zn,a,b,zm,z0,verbose):
-
basic.galaxy.
frac
(zn, zb, a, b, zm, verbose=False, zbias=0.0, sigma=0.0)¶ Computing z-interval of z-bin so that number of galaxies at each z-bin is equal
- Args:
zn (int): number of z-bins a, b (double): shape parameters of Schechter-like galaxy distribution zm (double): mean redshift zb[zn+1] (double): z-intervals - Args(optional):
verbose (bool): output messages (default to True) zbias (double): constant bias to true photo-z sigma (double): uncertaines of photo-z - Returns:
nfrac[zn] (double): fraction of galaxy number at each bin, defined by int_zi^zi+1 dz N(z)/int dz N(z) - Usage:
nfrac = basic.galaxy.frac(zn,zb,a,b,zm,zbias,sigma,verbose):
cosmofuncs
– tools for cosmology related functions¶
-
basic.cosmofuncs.
hubble
(z, H0=70.0, Ov=0.7, Om=0.3, w0=-1.0, wa=0.0, zn=0, divc=False)¶ Compute the expansion rate in unit of 1/Mpc, H/c, or in unit of km/s/Mpc, H.
- Args:
z[zn] (double): Redshifts at which H is computed - Args(optional):
H0 (double): The current value of hubble parameter in km/s/Mpc, default to 70 km/s/Mpc Om (double): The current value of Omega_matter, default to 0.3 Ov (double): The current value of Omega_Dark-energy, default to 0.7 w0, wa (double): The EoS of Dark Energy, default to w0=-1 and wa=0. divc (bool): Divide H by c or not, default to False. - Returns:
Hz[zn] (double): The expansion rate, H(z), divided by c or not. - Usage:
Hz = basic.cosmofuncs.hubble(z,H0,Om,Ov,w0,wa,zn,divc):
-
basic.cosmofuncs.
dhubble_dz
(z, H0=70.0, Ov=0.7, Om=0.3, w0=-1.0, wa=0.0, zn=0)¶ - Usage:
dHdz = basic.cosmofuncs.dhubble_dz(z,H0,Om,Ov,w0,wa,zn):
-
basic.cosmofuncs.
dist2z
(rz, H0=70.0, Ov=0.7, Om=0.3, w0=-1.0, wa=0.0, zn=0)¶ - Usage:
z = basic.cosmofuncs.dist2z(rz,H0,Om,Ov,w0,wa,zn):
-
basic.cosmofuncs.
dist_comoving
(z, H0=70.0, Ov=0.7, Om=0.3, w0=-1.0, wa=0.0, zn=0)¶ - Usage:
rz = basic.cosmofuncs.dist_comoving(z,H0,Om,Ov,w0,wa,zn):
-
basic.cosmofuncs.
dist_luminosity
(z, H0=70.0, Ov=0.7, Om=0.3, w0=-1.0, wa=0.0, zn=0)¶ - Usage:
DLz = basic.cosmofuncs.dist_luminosity(z,H0,Om,Ov,w0,wa,zn):
-
basic.cosmofuncs.
growth_factor
(z, normed=False, H0=70.0, Ov=0.7, Om=0.3, w0=-1.0, wa=0.0, zn=0)¶ - Usage:
Dz = basic.cosmofuncs.growth_factor(z,H0,Om,Ov,w0,wa,zn,normed):
-
basic.cosmofuncs.
growth_rate
(z, H0=70.0, Ov=0.7, Om=0.3, w0=-1.0, wa=0.0, zn=0)¶ - Usage:
fz = basic.cosmofuncs.growth_rate(z,H0,Om,Ov,w0,wa,zn):
-
basic.cosmofuncs.
nz_gw
(z, Cz, Hz, ntype='CH06', dotn0=1e-06, Tobs=3.0)¶ Distribution function of NS binary sources per redshift (dN/dz) at z
- Args:
z (double): redshift Cz (double): comoving distance Hz (double): expansion rate - Args(optional):
ntype (str): type of dotn functional form, i.e, CH06 (default) or none. dotn0 (double): current merger-rate Tobs (double): total observation time - Returns:
nz (double): distribution function at z - Usage:
nz = basic.cosmofuncs.nz_gw(z,Cz,Hz,ntype,dotn0,Tobs):
-
basic.cosmofuncs.
drate_dz
(z, ntype='CH06', zn=0)¶ - Usage:
dndz = basic.cosmofuncs.drate_dz(z,zn,ntype):
flat
– cross-check tools for flat-sky normalization¶
-
basic.flat.
alxy
(qest, qtype, lmax, rlmin, rlmax, fC, W1, W2, gln=100, lxcut=0, gle=1e-14)¶ Compute flat-sky quadratic estimator normalization CAUTION: This code interpolates the input Cl at the non-integer multipole by simply Cl(int(ell)) which leads to a small discrepancy in the normalization computed from the FFT-based method (which uses linear interpolation) and from this code. It is desireble to use the FFT-based normalization if you want to normalize the simulation results.
- Args:
qest (str): estimator combination (TT, TE, TB, EE, EB, or BB) qtype (str): estimator type (lensing, patchytau) lmax (double): output maximum multipole rlmax/rlmin (double): input CMB multipole range for reconstruction fC[rlmax] (double): power spectrum in the numerator :W1/W2[rlmax] : inverse of the observed power spectrum
- Args(optional):
gln (int): number of the GL integration points lxcut (int): multipole cut in x-direction, |l_x| < lx gle (double): convergence parameter for the GL integration - Returns:
Ag/Ac[l] (double): normalization for even and odd estimator pairs - Usage:
Ag,Ac = basic.flat.alxy(qest,qtype,lmax,rlmin,rlmax,fC,W1,W2,gln,gle,lxcut):
-
basic.flat.
alxy_asym
(qest, qtype, lmax, rlmin, rlmax, fC, AA, BB, AB, gln=100, lxcut=0, gle=1e-14)¶ - Usage:
Ag,Ac = basic.flat.alxy_asym(qest,qtype,lmax,rlmin,rlmax,fC,AA,BB,AB,gln,gle,lxcut):
-
basic.flat.
bbxy
(lmax, rlmin, rlmax, XX, YY, weight='lensing', gln=100, gle=1e-14)¶ - Usage:
BB = basic.flat.bbxy(lmax,rlmin,rlmax,XX,YY,weight,gln,gle):
wigner_funcs
– Wigner 3j symbols¶
-
basic.wigner_funcs.
wigner_3j
(l2, l3, m2, m3)¶ Compute wigner3j for all possible l1 where w3j = (j1,j2,j3/m1,m2,m3)
- Args:
shap (str): shape of the bispectrum (equi, fold, sque, or isos) cpmodel (str): cosmological parameter model (model0, modelw, or modelp) model (str): fitting formula of the matter bispectrum (LN=linear, SC=SC03, GM=Gil-Marin+12, 3B=3-shape-bispectrum, or RT=Takahashi+19) z[zn] (double): redshift points for the z-integral zn (int): number of redshifts for the z-integral dz[zn] (double): interval of z zs[3] (double): source redshifts lmin/lmax (int): minimum/maximum multipoles of the bispectrum k[kn] (double): k for the matter power spectrum pk0 (double): the linear matter power spectrum at z=0 kn (int): size of k - Returns:
bl0[l] (double): lensing bispectrum from LSS contributions at [lmin,lmax] bl1[l] (double): lensing bispectrum from post-Born contributions at [lmin,lmax] - Usage:
w3j = basic.wigner_funcs.wigner_3j(l2,l3,m2,m3):