curvedsky
– Curvedsky analysis tools¶
rec_lens
– quadratic lensing reconstruction¶
-
curvedsky.rec_lens.
qtt
(lmax, rlmin, rlmax, fC, Tlm1, Tlm2, nside_t=0, gtype='', verbose=False)¶ Reconstructing CMB lensing potential and its curl mode from the temperature quadratic estimator
- Args:
lmax (int): Maximum multipole of output lensing potential alms rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fC [l] (double): TT spectrum, with bounds (0:rlmax) Tlm1 [l,m] (dcmplx): 1st inverse-variance filtered temperature alm, with bounds (0:rlmax,0:rlmax) Tlm2 [l,m] (dcmplx): 2nd inverse-variance filtered temperature alm, with bounds (0:rlmax,0:rlmax) - Args(optional):
nside_t (int): Nside for the convolution calculation gtype (str): Type of output, i.e., convergence (gtype=’k’) or lensing potential (gtype=’‘, default) verbose (bool): Output messages, default to False - Returns:
glm [l,m] (dcmplx): CMB lensing potential alm, with bounds (0:lmax,0:lmax) clm [l,m] (dcmplx): Curl mode (pseudo lensing potential) alm, with bounds (0:lmax,0:lmax) - Usage:
glm,clm = curvedsky.rec_lens.qtt(lmax,rlmin,rlmax,fC,Tlm1,Tlm2,nside_t,gtype,verbose):
-
curvedsky.rec_lens.
qte
(lmax, rlmin, rlmax, fC, Tlm, Elm, nside_t=0, gtype='', verbose=False)¶ Reconstructing CMB lensing potential and its curl mode from the TE quadratic estimator
- Args:
lmax (int): Maximum multipole of output lensing potential alms rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fC [l] (double): TE spectrum, with bounds (0:rlmax) Tlm [l,m] (dcmplx): Inverse-variance filtered temperature alm, with bounds (0:rlmax,0:rlmax) Elm [l,m] (dcmplx): Inverse-variance filtered E-mode alm, with bounds (0:rlmax,0:rlmax) - Args(optional):
nside_t (int): Nside for the convolution calculation gtype (str): Type of output, i.e., convergence (gtype=’k’) or lensing potential (gtype=’‘, default) verbose (bool): Output messages, default to False - Returns:
glm [l,m] (dcmplx): CMB lensing potential, with bounds (0:lmax,0:lmax) clm [l,m] (dcmplx): Curl mode (pseudo lensing potential), with bounds (0:lmax,0:lmax) - Usage:
glm,clm = curvedsky.rec_lens.qte(lmax,rlmin,rlmax,fC,Tlm,Elm,nside_t,gtype,verbose):
-
curvedsky.rec_lens.
qtb
(lmax, rlmin, rlmax, fC, Tlm, Blm, nside_t=0, gtype='', verbose=False)¶ Reconstructing CMB lensing potential and its curl mode from the TB quadratic estimator
- Args:
lmax (int): Maximum multipole of output lensing potential alms rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fC [l] (double): TE spectrum, with bounds (0:rlmax) Tlm [l,m] (dcmplx): Inverse-variance filtered temperature alm, with bounds (0:rlmax,0:rlmax) Blm [l,m] (dcmplx): Inverse-variance filtered B-mode alm, with bounds (0:rlmax,0:rlmax) - Args(optional):
nside_t (int): Nside for the convolution calculation gtype (str): Type of output, i.e., convergence (gtype=’k’) or lensing potential (gtype=’‘, default) verbose (bool): Output messages, default to False - Returns:
glm [l,m] (dcmplx): CMB lensing potential, with bounds (0:lmax,0:lmax) clm [l,m] (dcmplx): Curl mode (pseudo lensing potential), with bounds (0:lmax,0:lmax) - Usage:
glm,clm = curvedsky.rec_lens.qtb(lmax,rlmin,rlmax,fC,Tlm,Blm,nside_t,gtype,verbose):
-
curvedsky.rec_lens.
qee
(lmax, rlmin, rlmax, fC, Elm1, Elm2, nside_t=0, gtype='', verbose=False)¶ Reconstructing CMB lensing potential and its curl mode from the EE quadratic estimator
- Args:
lmax (int): Maximum multipole of output lensing potential alms rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fC [l] (double): EE spectrum, with bounds (0:rlmax) Elm1 [l,m] (dcmplx): 1st inverse-variance filtered E-mode alm, with bounds (0:rlmax,0:rlmax) Elm2 [l,m] (dcmplx): 2nd inverse-variance filtered E-mode alm, with bounds (0:rlmax,0:rlmax) - Args(optional):
nside_t (int): Nside for the convolution calculation gtype (str): Type of output, i.e., convergence (gtype=’k’) or lensing potential (gtype=’‘, default) verbose (bool): Output messages, default to False - Returns:
glm [l,m] (dcmplx): CMB lensing potential, with bounds (0:lmax,0:lmax) clm [l,m] (dcmplx): Curl mode (pseudo lensing potential), with bounds (0:lmax,0:lmax) - Usage:
glm,clm = curvedsky.rec_lens.qee(lmax,rlmin,rlmax,fC,Elm1,Elm2,nside_t,gtype,verbose):
-
curvedsky.rec_lens.
qeb
(lmax, rlmin, rlmax, fC, Elm, Blm, nside_t=0, gtype='', verbose=False)¶ Reconstructing CMB lensing potential and its curl mode from the EB quadratic estimator
- Args:
lmax (int): Maximum multipole of output lensing potential alms rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fC [l] (double): EE spectrum, with bounds (0:rlmax) Elm [l,m] (dcmplx): Inverse-variance filtered E-mode alm, with bounds (0:rlmax,0:rlmax) Blm [l,m] (dcmplx): Inverse-variance filtered B-mode alm, with bounds (0:rlmax,0:rlmax) - Args(optional):
nside_t (int): Nside for the convolution calculation gtype (str): Type of output, i.e., convergence (gtype=’k’) or lensing potential (gtype=’‘, default) verbose (bool): Output messages, default to False - Returns:
glm [l,m] (dcmplx): CMB lensing potential, with bounds (0:lmax,0:lmax) clm [l,m] (dcmplx): Curl mode (pseudo lensing potential), with bounds (0:lmax,0:lmax) - Usage:
glm,clm = curvedsky.rec_lens.qeb(lmax,rlmin,rlmax,fC,Elm,Blm,nside_t,gtype,verbose):
-
curvedsky.rec_lens.
qbb
(lmax, rlmin, rlmax, fC, Blm1, Blm2, nside_t=0, gtype='', verbose=False)¶ Reconstructing CMB lensing potential and its curl mode from the BB quadratic estimator
- Args:
lmax (int): Maximum multipole of output lensing potential alms rlmin/rlmax (int): Minimum/Maximum multipoles of CMB for reconstruction fC [l] (double): BB spectrum, with bounds (0:rlmax) Blm1 [l,m] (dcmplx): 1st inverse-variance filtered B-mode alm, with bounds (0:rlmax,0:rlmax) Blm2 [l,m] (dcmplx): 2nd inverse-variance filtered B-mode alm, with bounds (0:rlmax,0:rlmax) - Args(optional):
nside_t (int): Nside for the convolution calculation gtype (str): Type of output, i.e., convergence (gtype=’k’) or lensing potential (gtype=’‘, default) verbose (bool): Output messages, default to False - Returns:
glm [l,m] (dcmplx): CMB lensing potential, with bounds (0:lmax,0:lmax) clm [l,m] (dcmplx): Curl mode (pseudo lensing potential), with bounds (0:lmax,0:lmax) - Usage:
glm,clm = curvedsky.rec_lens.qbb(lmax,rlmin,rlmax,fC,Blm1,Blm2,nside_t,gtype,verbose):
rec_rot
– quadratic rotation reconstruction¶
-
curvedsky.rec_rot.
qeb
(lmax, rlmin, rlmax, fCE, Elm, Blm, nside_t=0, verbose=False)¶ Reconstructing pol rotation angle from the EB quadratic estimator
- Args:
lmax (int): Maximum multipole of output lensing potential alms rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fCE [l] (double): EE spectrum, with bounds (0:rlmax) Elm [l,m] (dcmplx): Inverse-variance filtered E-mode alm, with bounds (0:rlmax,0:rlmax) Blm [l,m] (dcmplx): Inverse-variance filtered B-mode alm, with bounds (0:rlmax,0:rlmax) - Args(optional):
nside_t (int): Nside for the convolution calculation verbose (bool): Output messages, default to False - Returns:
alm [l,m] (dcmplx): Rotation angle alm, with bounds (0:lmax,0:lmax) - Usage:
alm = curvedsky.rec_rot.qeb(lmax,rlmin,rlmax,fCE,Elm,Blm,nside_t,verbose):
rec_tau
– quadratic tau reconstruction¶
-
curvedsky.rec_tau.
qtt
(lmax, rlmin, rlmax, fC, Tlm1, Tlm2, nside_t=0, verbose=False)¶ Reconstructing inhomogeneous tau from the temperature quadratic estimator
- Args:
lmax (int): Maximum multipole of output tau alms rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fC [l] (double): TT spectrum, with bounds (1:rlmax) Tlm1 [l,m] (dcmplx): 1st inverse-variance filtered temperature alm, with bounds (0:rlmax,0:rlmax) Tlm2 [l,m] (dcmplx): 2nd inverse-variance filtered temperature alm, with bounds (0:rlmax,0:rlmax) - Args(optional):
nside_t (int): Nside for the convolution calculation verbose (bool): Output messages, default to False - Returns:
alm [l,m] (dcmplx): Amplitude modulation alm, with bounds (0:lmax,0:lmax) - Usage:
alm = curvedsky.rec_tau.qtt(lmax,rlmin,rlmax,fC,Tlm1,Tlm2,nside_t,verbose):
-
curvedsky.rec_tau.
qeb
(lmax, rlmin, rlmax, fCE, Elm, Blm, nside_t=0, verbose=False)¶ Reconstructing amplitude modulation from the EB quadratic estimator
- Args:
lmax (int): Maximum multipole of output lensing potential alms rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fCE [l] (double): EE spectrum, with bounds (0:rlmax) Elm [l,m] (dcmplx): Inverse-variance filtered E-mode alm, with bounds (0:rlmax,0:rlmax) Blm [l,m] (dcmplx): Inverse-variance filtered B-mode alm, with bounds (0:rlmax,0:rlmax) - Args(optional):
nside_t (int): Nside for the convolution calculation verbose (bool): Output messages, default to False - Returns:
alm [l,m] (dcmplx): Amplitude modulation alm, with bounds (0:lmax,0:lmax) - Usage:
alm = curvedsky.rec_tau.qeb(lmax,rlmin,rlmax,fCE,Elm,Blm,nside_t,verbose):
rec_src
– quadratic point-soruce reconstruction¶
-
curvedsky.rec_src.
qtt
(lmax, rlmin, rlmax, Tlm1, Tlm2, nside_t=0, verbose=False)¶ Reconstructing point sources from the temperature quadratic estimator
- Args:
lmax (int): Maximum multipole of output point-source alms rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction Tlm1 [l,m] (dcmplx): 1st inverse-variance filtered temperature alm, with bounds (0:rlmax,0:rlmax) Tlm2 [l,m] (dcmplx): 2nd inverse-variance filtered temperature alm, with bounds (0:rlmax,0:rlmax) - Args(optional):
nside_t (int): Nside for the convolution calculation verbose (bool): Output messages, default to False - Returns:
slm [l,m] (dcmplx): Point-source alm, with bounds (0:lmax,0:lmax) - Usage:
slm = curvedsky.rec_src.qtt(lmax,rlmin,rlmax,Tlm1,Tlm2,nside_t,verbose):
norm_quad
– normalization of quadratic estimator reconstruction¶
-
curvedsky.norm_quad.
qtt
(est, lmax, rlmin, rlmax, TT, OCT, lfac='')¶ Normalization of reconstructed fields from the temperature quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction TT [l] (double): Theory TT spectrum, with bounds (0:rlmax) OCT [l] (double): Observed TT spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Multiplying square of L(L+1)/2, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Al [2,l] (double): Normalizations (1 is dummy except lens = 0 and curl = 1), with bounds (0:lmax) - Usage:
Al = curvedsky.norm_quad.qtt(est,lmax,rlmin,rlmax,TT,OCT,lfac):
-
curvedsky.norm_quad.
qte
(est, lmax, rlmin, rlmax, TE, OCT, OCE, lfac='')¶ Normalization of reconstructed fields from the TE quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction TE [l] (double): Theory TE spectrum, with bounds (0:rlmax) OCT [l] (double): Observed TT spectrum, with bounds (0:rlmax) OCE [l] (double): Observed EE spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Multiplying square of L(L+1)/2, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Al [2,l] (double): Normalizations (1 is dummy except lens = 0 and curl = 1), with bounds (0:lmax) - Usage:
Al = curvedsky.norm_quad.qte(est,lmax,rlmin,rlmax,TE,OCT,OCE,lfac):
-
curvedsky.norm_quad.
qtb
(est, lmax, rlmin, rlmax, TE, OCT, OCB, lfac='')¶ Normalization of reconstructed fields from the TB quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction TE [l] (double): Theory TE spectrum, with bounds (0:rlmax) OCT [l] (double): Observed TT spectrum, with bounds (0:rlmax) OCB [l] (double): Observed BB spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Multiplying square of L(L+1)/2, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Al [2,l] (double): Normalizations (1 is dummy except lens = 0 and curl = 1), with bounds (0:lmax) - Usage:
Al = curvedsky.norm_quad.qtb(est,lmax,rlmin,rlmax,TE,OCT,OCB,lfac):
-
curvedsky.norm_quad.
qee
(est, lmax, rlmin, rlmax, EE, OCE, lfac='')¶ Normalization of reconstructed amplitude modulation from the EE quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction EE [l] (double): Theory EE spectrum, with bounds (0:rlmax) OCE [l] (double): Observed EE spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Multiplying square of L(L+1)/2, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Al [2,l] (double): Normalization, with bounds (0:lmax) - Usage:
Al = curvedsky.norm_quad.qee(est,lmax,rlmin,rlmax,EE,OCE,lfac):
-
curvedsky.norm_quad.
qeb
(est, lmax, rlmin, rlmax, EE, OCE, OCB, lfac='', BB=0)¶ Normalization of reconstructed fields from the EB quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction EE [l] (double): Theory EE spectrum, with bounds (0:rlmax) OCE [l] (double): Observed EE spectrum, with bounds (0:rlmax) OCB [l] (double): Observed BB spectrum, with bounds (0:rlmax) - Args(optionals):
BB [l] (double): Theory BB spectrum, with bounds (0:rlmax) lfac (str): Multiplying square of L(L+1)/2, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Al [2,l] (double): Normalization, with bounds (0:lmax) - Usage:
Al = curvedsky.norm_quad.qeb(est,lmax,rlmin,rlmax,EE,OCE,OCB,BB,lfac):
-
curvedsky.norm_quad.
qbb
(est, lmax, rlmin, rlmax, BB, OCB, lfac='')¶ Normalization of reconstructed amplitude modulation from the BB quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction BB [l] (double): Theory BB spectrum, with bounds (0:rlmax) OCB [l] (double): Observed BB spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Multiplying square of L(L+1)/2, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Al [2,l] (double): Normalization, with bounds (0:lmax) - Usage:
Al = curvedsky.norm_quad.qbb(est,lmax,rlmin,rlmax,BB,OCB,lfac):
-
curvedsky.norm_quad.
qttte
(est, lmax, rlmin, rlmax, fCTT, fCTE, OCT, OCE, OCTE, lfac='')¶ Correlation between unnormalized TT and TE quadratic estimators
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fCTT [l] (double): Theory TT spectrum, with bounds (0:rlmax) fCTE [l] (double): Theory TE spectrum, with bounds (0:rlmax) OCT [l] (double): Observed TT spectrum, with bounds (0:rlmax) OCE [l] (double): Observed EE spectrum, with bounds (0:rlmax) OCTE [l] (double): Observed TE spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Multiplying square of L(L+1)/2, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Ig [l] (double): Correlation between lensing potential estimators, with bounds (0:lmax) Ic [l] (double): Correlation between curl mode estimators, with bounds (0:lmax) - Usage:
Ig,Ic = curvedsky.norm_quad.qttte(est,lmax,rlmin,rlmax,fCTT,fCTE,OCT,OCE,OCTE,lfac):
-
curvedsky.norm_quad.
qttee
(est, lmax, rlmin, rlmax, fCTT, fCEE, OCT, OCE, OCTE, lfac='')¶ Correlation between unnormalized TT and EE quadratic estimators
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fCTT [l] (double): Theory TT spectrum, with bounds (0:rlmax) fCEE [l] (double): Theory EE spectrum, with bounds (0:rlmax) OCT [l] (double): Observed TT spectrum, with bounds (0:rlmax) OCE [l] (double): Observed EE spectrum, with bounds (0:rlmax) OCTE [l] (double): Observed TE spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Multiplying square of L(L+1)/2, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Ig [l] (double): Correlation between lensing potential estimators, with bounds (0:lmax) Ic [l] (double): Correlation between curl mode estimators, with bounds (0:lmax) - Usage:
Ig,Ic = curvedsky.norm_quad.qttee(est,lmax,rlmin,rlmax,fCTT,fCEE,OCT,OCE,OCTE,lfac):
-
curvedsky.norm_quad.
qteee
(est, lmax, rlmin, rlmax, fCEE, fCTE, OCT, OCE, OCTE, lfac='')¶ Correlation between unnormalized TE and EE quadratic estimators
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fCEE [l] (double): Theory EE spectrum, with bounds (0:rlmax) fCTE [l] (double): Theory TE spectrum, with bounds (0:rlmax) OCT [l] (double): Observed TT spectrum, with bounds (0:rlmax) OCE [l] (double): Observed EE spectrum, with bounds (0:rlmax) OCTE [l] (double): Observed TE spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Multiplying square of L(L+1)/2, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Ig [l] (double): Correlation between lensing potential estimators, with bounds (0:lmax) Ic [l] (double): Correlation between curl mode estimators, with bounds (0:lmax) - Usage:
Ig,Ic = curvedsky.norm_quad.qteee(est,lmax,rlmin,rlmax,fCEE,fCTE,OCT,OCE,OCTE,lfac):
-
curvedsky.norm_quad.
qtbeb
(est, lmax, rlmin, rlmax, fCEE, fCBB, fCTE, OCT, OCE, OCB, OCTE, lfac='')¶ Correlation between unnormalized TB and EB quadratic estimators
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fCEE [l] (double): Theory EE spectrum, with bounds (0:rlmax) fCBB [l] (double): Theory BB spectrum, with bounds (0:rlmax) OCT [l] (double): Observed TT spectrum, with bounds (0:rlmax) OCE [l] (double): Observed EE spectrum, with bounds (0:rlmax) OCB [l] (double): Observed BB spectrum, with bounds (0:rlmax) OCTE [l] (double): Observed TE spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Multiplying square of L(L+1)/2, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Ig [l] (double): Correlation between lensing potential estimators, with bounds (0:lmax) Ic [l] (double): Correlation between curl mode estimators, with bounds (0:lmax) - Usage:
Ig,Ic = curvedsky.norm_quad.qtbeb(est,lmax,rlmin,rlmax,fCEE,fCBB,fCTE,OCT,OCE,OCB,OCTE,lfac):
-
curvedsky.norm_quad.
qmv
(lmax, QDO, Al, Il)¶ Compute MV estimator normalization. Currently BB is ignored.
- Args:
lmax (int): Maximum multipole of the output power spectra QDO[6] (bool): Specifying which estimators to be combined for the minimum variance estimator, with size (6). The oder is TT, TE, EE, TB, EB and BB. Currently, BB is always False Al [5,l] (double): Normalizations of each estimator (TT, TE, EE, TB, EB). Il [4,l] (double): Correlation between different estimators (TTxTE, TTxEE, TExEE, TBxEB). - Returns:
MV [l] (double): Normalization of the MV estimator, with bounds (0:lmax) Nl [6,l] (double): Weights for each estimator (TT, TE, EE, TB, EB, BB=0), with bounds (0:lmax) - Usage:
MV,Nl = curvedsky.norm_quad.qmv(lmax,QDO,Al,Il):
-
curvedsky.norm_quad.
qall
(est, QDO, lmax, rlmin, rlmax, fC, OC, lfac='')¶ Compute MV estimator normalization. Currently BB is ignored.
- Args:
QDO[6] (bool): Specifying which estimators to be combined for the minimum variance estimator, with size (6). The oder is TT, TE, EE, TB, EB and BB. lmax (int): Maximum multipole of the output power spectra rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fC/OC [l] (double): Theory/Observed CMB angular power spectra (TT, EE, BB, TE), with bounds (0:rlmax) - Args(optional):
lfac (str): Multiplying square of L(L+1)/2, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Ag [6,l] (double): Normalization of the TT, TE, EE, TB, EB, and MV estimators for lensing potential, with bounds (6,0:lmax) Ac [6,l] (double): Same as Ag but for curl mode Nlg [6,l] (double): Weights for TT, TE, EE, TB, EB, and BB (=0) estimators for lensing potential, with bounds (6,0:lmax) Nlc [6,l] (double): Same as Nlg but for curl mode - Usage:
Ag,Ac,Nlg,Nlc = curvedsky.norm_quad.qall(est,QDO,lmax,rlmin,rlmax,fC,OC,lfac):
-
curvedsky.norm_quad.
qeb_iter
(lmax, elmax, rlmin, rlmax, dlmin, dlmax, CE, OCE, OCB, Cpp, iter=1, conv=1e-05)¶ Normalization of reconstructed CMB lensing potential and its curl mode from the EB quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization elmax (int): Maximum multipole of input EE spectra, CE and OCE rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction dlmin/dlmax (int): Minimum/Maximum multipole of E mode and lensing potential for delensing CE [l] (double): Theory EE angular power spectrum, with bounds (0:elmax) OCE [l] (double): Observed EE spectrum, with bounds (0:elmax) OCB [l] (double): Observed BB spectrum, with bounds (0:rlmax) Cpp [l] (double): Theory lensing potential spectrum, with bounds (0:dlmax) - Args(optional):
iter (int): number of iteration, default to 1 (no iteration) conv (double): a parameter for convergence the iteration, default to 0.00001 - Returns:
Ag [l] (double): CMB lensing potential normalization, with bounds (0:lmax) Ac [l] (double): Curl mode (pseudo lensing potential) normalization, with bounds (0:lmax) - Usage:
Ag,Ac = curvedsky.norm_quad.qeb_iter(lmax,elmax,rlmin,rlmax,dlmin,dlmax,CE,OCE,OCB,Cpp,iter,conv):
-
curvedsky.norm_quad.
xtt
(est, lmax, rlmin, rlmax, fC, OCT, lfac='')¶ Unnormalized response between lensing potential and amplitude modulation from the temperature quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fC [l] (double): Theory TT spectrum, with bounds (0:rlmax) OCT [l] (double): Observed TT spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Multiplying square of L(L+1)/2, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Ag [l] (double): Cross normalization, with bounds (0:lmax) - Usage:
Ag = curvedsky.norm_quad.xtt(est,lmax,rlmin,rlmax,fC,OCT,lfac):
-
curvedsky.norm_quad.
xeb
(est, lmax, rlmin, rlmax, EE, EB, OCE, OCB, BB=None)¶ Response of reconstructed field to other in the EB quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction EE [l] (double): Theory EE spectrum, with bounds (0:rlmax) EB [l] (double): Theory EB spectrum, with bounds (0:rlmax) OCE [l] (double): Observed EE spectrum, with bounds (0:rlmax) OCB [l] (double): Observed BB spectrum, with bounds (0:rlmax) - Args(optionals):
BB [l] (double): Theory BB spectrum, with bounds (0:rlmax) - Returns:
Aa [l] (double): Pol. rot. angle normalization, with bounds (0:lmax) - Usage:
Aa = curvedsky.norm_quad.xeb(est,lmax,rlmin,rlmax,EE,EB,OCE,OCB,BB):
norm_imag
– normalization of imaginary quadratic estimator reconstruction¶
-
curvedsky.norm_imag.
qte
(est, lmax, rlmin, rlmax, TB, OCT, OCE, lfac='')¶ Normalization of reconstructed imaginary CMB lensing potential and its curl mode from the TE quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction TB [l] (double): Theory TB spectrum, with bounds (0:rlmax) OCT [l] (double): Observed TT spectrum, with bounds (0:rlmax) OCE [l] (double): Observed EE spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Type of output, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Al [2,l] (double): Normalization, with bounds (0:lmax) - Usage:
Al = curvedsky.norm_imag.qte(est,lmax,rlmin,rlmax,TB,OCT,OCE,lfac):
-
curvedsky.norm_imag.
qtb
(est, lmax, rlmin, rlmax, TB, OCT, OCB, lfac='')¶ Normalization of reconstructed imaginary CMB lensing potential and its curl mode from the TB quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction TB [l] (double): Theory TE spectrum, with bounds (0:rlmax) OCT [l] (double): Observed TT spectrum, with bounds (0:rlmax) OCB [l] (double): Observed BB spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Type of output, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Al [2,l] (double): Normalization, with bounds (0:lmax) - Usage:
Al = curvedsky.norm_imag.qtb(est,lmax,rlmin,rlmax,TB,OCT,OCB,lfac):
-
curvedsky.norm_imag.
qee
(est, lmax, rlmin, rlmax, fC, OCE, lfac='')¶ Normalization of reconstructed imaginary CMB lensing potential and its curl mode from the E-mode quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fC [l] (double): Theory EE spectrum, with bounds (0:rlmax) OCE [l] (double): Observed EE spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Type of output, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Al [2,l] (double): Normalization, with bounds (0:lmax) - Usage:
Al = curvedsky.norm_imag.qee(est,lmax,rlmin,rlmax,fC,OCE,lfac):
-
curvedsky.norm_imag.
qeb
(est, lmax, rlmin, rlmax, fC, OCE, OCB, lfac='')¶ Normalization of reconstructed imaginary CMB lensing potential and its curl mode from the EB quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fC [l] (double): Theory EB spectrum, with bounds (0:rlmax) OCE [l] (double): Observed EE spectrum, with bounds (0:rlmax) OCB [l] (double): Observed BB spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Type of output, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Al [2,l] (double): Normalization, with bounds (0:lmax) - Usage:
Al = curvedsky.norm_imag.qeb(est,lmax,rlmin,rlmax,fC,OCE,OCB,lfac):
-
curvedsky.norm_imag.
qbb
(est, lmax, rlmin, rlmax, fC, OCB, lfac='')¶ Normalization of reconstructed imaginary CMB lensing potential and its curl mode from the B-mode quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fC [l] (double): Theory BB spectrum, with bounds (0:rlmax) OCB [l] (double): Observed BB spectrum, with bounds (0:rlmax) - Args(optional):
lfac (str): Type of output, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Al [2,l] (double): Normalization, with bounds (0:lmax) - Usage:
Al = curvedsky.norm_imag.qbb(est,lmax,rlmin,rlmax,fC,OCB,lfac):
-
curvedsky.norm_imag.
qbb_asym
(est, lmax, rlmin, rlmax, fC, OCB1, OCB2, lfac='')¶ Normalization of reconstructed imaginary CMB lensing potential and its curl mode from the B-mode quadratic estimator
- Args:
lmax (int): Maximum multipole of output normalization spectrum rlmin/rlmax (int): Minimum/Maximum multipole of CMB for reconstruction fC [l] (double): Theory BB spectrum, with bounds (0:rlmax) OCB1/2 [l] (double): Observed BB spectrum for experiment 1 and 2, with bounds (0:rlmax) - Args(optional):
lfac (str): Type of output, i.e., convergence (lfac=’k’) or lensing potential (lfac=’‘, default) - Returns:
Al [2,l] (double): Normalization, with bounds (0:lmax) - Usage:
Al = curvedsky.norm_imag.qbb_asym(est,lmax,rlmin,rlmax,fC,OCB1,OCB2,lfac):
delens
– delensing for simulation¶
-
curvedsky.delens.
lensingb
(lmax, elmin, elmax, plmin, plmax, wElm, wplm, nside_t=0, gtype='p')¶ Computing lensing B mode as a convolution of wiener-filtered E-mode and lensing potential
- Args:
lmax (int): Maximum multipole of output lensing B-mode alm elmin (int): Minimum multipole of wiener-filtered E-mode alm elmax (int): Maximum multipole of wiener-filtered E-mode alm plmin (int): Minimum multipole of wiener-filtered lensing potential alm plmax (int): Maximum multipole of wiener-filtered lensing potential alm wElm [l,m] (dcmplx): Wiener-filtered E-mode alm, with bounds (0:elmax,0:elmax) wplm [l,m] (dcmplx): Wiener-filtered lensing potential (or kappa) alm, with bounds (0:plmax,0:plmax) - Args(optional):
nside_t (int): Nside for the convolution calculation gtype (str): Type of input wplm (‘p’=phi or ‘k’=kappa), default to ‘p’ (phi). - Returns:
lBlm [l,m] (dcmplx): Lensing B-mode alm, with bounds (0:lmax,0:lmax) - Usage:
lBlm = curvedsky.delens.lensingb(lmax,elmin,elmax,plmin,plmax,wElm,wplm,nside_t,gtype):
-
curvedsky.delens.
shiftvec
(nside, lmax, plm, nremap)¶ Return the anti deflection vector, beta, at the Healpix pixel for the delensing where
beta(n) + alphaiw(n+beta(n)) = 0and alphaw is the filtered lensing deflection vector (see arXiv:1701.01712).
- Args:
nside (int): Nside of output shift vector lmax (int): Maximum multipole of the input plm plm [l,m] (dcmplx): Wiener-filtered lensing potential alm, with bounds (0:lmax,0:lmax) - Args(optional):
nremap (int): Number of iteration for computing the shift vector - Returns:
beta [pix,2] (double): 2D shift vector, with bounds (0:npix-1,1:2) - Usage:
beta = curvedsky.delens.shiftvec(nside,lmax,plm,nremap):
-
curvedsky.delens.
phi2grad
(nside, lmax, plm)¶ Return the deflection vector, grad, at the Healpix pixel
- Args:
nside (int): Nside of output deflection vector lmax (int): Maximum multipole of the input plm/clm plm [l,m] (dcmplx): Lensing potential alm, with bounds (0:lmax,0:lmax) - Returns:
grad [pix,2] (double): 2D deflection vector, with bounds (0:npix-1,1:2) - Usage:
grad = curvedsky.delens.phi2grad(nside,lmax,plm):
-
curvedsky.delens.
remap_tp
(nside, lmax, beta, alm_in)¶ Remapping CMB temperaure and polarization with a given shift vector based on a simple implementation of LensPix This function returs X(n+beta) where n is the original direction and beta is the shift vector The output is given by alms where alm[0,l,m] is temperature, alm[1,l,m] is E mode, and alm[2,l,m] is B mode.
- Args:
nside (int): Nside of input shift vector lmax (int): Maximum multipole of the input plm beta [pix,2] (double): 2D shift vector, with bounds (0:npix-1,1:2) alm_in [TEB,l,m] (dcmplx): Input T/E/B alms to be remapped, with bounds (1:3,0:lmax,0:lmax). - Returns:
alm_re [TEB,l,m] (dcmplx): Remapped T/E/B alms, with bounds (1:3,0:lmax,0:lmax). - Usage:
alm_re = curvedsky.delens.remap_tp(nside,lmax,beta,alm_in):
-
curvedsky.delens.
remap_tp_map
(nside, lmax, beta, alm_in)¶ Remapping CMB temperaure and polarization with a given shift vector based on a simple implementation of LensPix This function returs X(n+beta) where n is the original direction and beta is the shift vector The output is given by alms where alm[0,l,m] is temperature, alm[1,l,m] is E mode, and alm[2,l,m] is B mode.
- Args:
nside (int): Nside of input shift vector lmax (int): Maximum multipole of the input plm beta [pix,2] (double): 2D shift vector, with bounds (0:npix-1,1:2) alm_in [TEB,l,m] (dcmplx): Input T/E/B alms to be remapped, with bounds (1:3,0:lmax,0:lmax). - Returns:
map_re [pix,TQU] (double): Remapped T/Q/U alms, with bounds (0:lmax,0:lmax,1:3). - Usage:
map_re = curvedsky.delens.remap_tp_map(nside,lmax,beta,alm_in):
bispec
– binned reduced bispectrum¶
-
curvedsky.bispec.
make_quad_gauss
(lmax, alm)¶ Return a non-Gaussian alm. The corresponding non-Gaussian field is defined as
delta^NL(n) = delta^L(n) + delta^L(n)**2where delta^L(n) is a gaussian field obtained from the input alm.
- Args:
lmax (int): Maximum multipole of alm alm [l,m] (dcmplx): Input harmonic coefficients, with bounds (0:lmax,0:lmax). - Returns:
qlm [l,m] (dcmplx): Output harmonic coefficients of the non-Gaussian fields, with bounds (0:lmax,0:lmax). - Usage:
qlm = curvedsky.bispec.make_quad_gauss(lmax,alm):
-
curvedsky.bispec.
bispec_norm
(bn, bp, bstype='equi', bst=2, sL=[0, 0])¶ Return normalization of the binned reduced bispectrum for a given multipole bin
- Args:
bn (int): Number of multipole bins bp [edge] (double): Bin edges, with bounds (bn+1) - Args(optional):
bstype (str): Configuration of the bispectrum, default to equilateral bst (int): A parameter, bst=nside/lmax, which controls the accuracy of the calculation, default to 2. More accurate for a larger value. sL[2] (int): The fixed bin for the squeezed configuration, b[sL,eL,eL], default to the lowest multipole bin - Returns:
norm [bin] (double): Normalization of the binned reduced bispectrum at each bin, with bounds (bn) - Usage:
norm = curvedsky.bispec.bispec_norm(bn,bp,bstype,bst,sL):
-
curvedsky.bispec.
bispec_bin
(bn, bp, lmax, alm, bst=2, bstype='equi', sL=[0, 0])¶ Return the unnormalized binned reduced bispectrum for a given multipole bin
- Args:
bn (int): Number of multipole bins bp [edge] (double): Bin edges, with bounds (bn+1) lmax (int): Maximum multipole of the input alm alm [l,m] (dcmplx): Input harmonic coefficients, with bounds (0:lmax,0:lmax) - Args(optional):
bstype (str): Configuration of the bispectrum, default to equilateral bst (int): A parameter, bst=nside/lmax, which controls the accuracy of the calculation, default to 2. More accurate for a larger value. sL[2] (int): The fixed bin for the squeezed configuration, b[sL,eL,eL], default to the lowest multipole bin - Returns:
bis [bin] (double): The unnormalized binned reduced bispectrum at each bin, with bounds (bn) - Usage:
bis = curvedsky.bispec.bispec_bin(bn,bp,lmax,alm,bstype,bst,sL):
-
curvedsky.bispec.
xbispec_bin
(bn, bp, lmax, n, alm, bst=2, bstype='equi', sL=[0, 0])¶ Return the unnormalized binned reduced cross-bispectrum for a given multipole bin
- Args:
bn (int): Number of multipole bins bp [edge] (double): Bin edges, with bounds (bn+1) lmax (int): Maximum multipole of the input alm n (int): Number of alms alm [n,l,m] (dcmplx): Input harmonic coefficients, with bounds (0:n-1,0:lmax,0:lmax) - Args(optional):
bstype (str): Configuration of the bispectrum, default to equilateral bst (int): A parameter, bst=nside/lmax, which controls the accuracy of the calculation, default to 2. More accurate for a larger value. sL[2] (int): The fixed bin for the squeezed configuration, b[sL,eL,eL], default to the lowest multipole bin - Returns:
bis [bin] (double): The unnormalized binned reduced bispectrum at each bin, with bounds (bn) - Usage:
bis = curvedsky.bispec.xbispec_bin(bn,bp,lmax,n,alm,bstype,bst,sL):
cninv
– Filtering Anisotropies¶
-
curvedsky.cninv.
cnfilter_freq
(n, mn, nside, lmax, cl, bl, iNcov, maps, chn=1, lmaxs=[0], nsides=[0], itns=[1], eps=[1e-06], filter='W', verbose=False, ro=50, stat='', inl=None)¶ Combining multiple frequency CMB maps optimally. The filtering would work if the noise variance is not significantly varied with scale (multipole). Please make sure that your input maps are beam-convolved. This code deconvolves the beam during filtering and the output are the filtered alms after the beam-deconvolution.
- Args:
n (int): Number of maps, i.e., temperature only (n=1), polarization only (n=2) or both (n=3) mn (int): Number of frequencies nside (int): Nside of input map lmax (int): Maximum multipole of the input cl cl[n,l] (double): Theory signal power spectrum, with bounds (0:n-1,0:lmax) bl[mn,l] (double): Beam spectrum, with bounds (0:mn-1,0:lmax) iNcov[n,mn,pix] (double): Inverse of the noise variance at each pixel, with bounds (0:n-1,0:mn-1,0:npix-1) maps[n,mn,pix] (double): Beam-convolved T, Q, U maps, with bouds (0:n-1,0:mn-1,0:npix-1) - Args(optional):
chn (int): Number of grids for preconsitioner (chn=1 for diagonal preconditioner, default) lmaxs[chain] (int): Maximum multipole(s) at each preconditioning and lmaxs[0] is the input maximum multipole of cl nsides[chain] (int): Nside(s) of preconditoner and nsides[0] should be consistent with the input map’s nside. eps[chain] (double): Parameter to finish the iteration (i.e. terminate if the residul fraction becomes smaller than eps). Default to 1e-6. itns[chain] (int): Number of interations. filter (str): C-inverse (‘’) or Wiener filter (W), default to the Wiener filter. inl[n,mn,l] (double): Inverse noise spectrum (0 for white noise case, default). verbose (bool): Output messages, default to False ro (int): the residual fraction is output for every ro iteration (e.g. ro=2 means 1 output per 2 iterations). Default to 50. Useful for convergence speed. stat (str): Realtime status filename which contains the residual fraction, default to no output file - Returns:
xlm[n,l,m] (dcmplx): C-inverse / Wiener filtered multipoles, with bounds (0:n-1,0:lmax,0:lmax) - Usage:
xlm = curvedsky.cninv.cnfilter_freq(n,mn,nside,lmax,cl,bl,iNcov,maps,chn,lmaxs,nsides,itns,eps,filter,inl,verbose,ro,stat):
-
curvedsky.cninv.
cnfilter_kappa
(n, nside, lmax, cov, iNcov, maps, chn=1, lmaxs=[0], nsides=[0], itns=[1], eps=[1e-06], verbose=False, ro=50, stat='', inl=None)¶ Computing the inverse-variance weighted multipole, (C+N)^-1 x kappa, for multiple mass-tracers of kappa maps.
- Args:
n (int): Number of input kappa maps to be combined nside (int): Nside of input maps lmax (int): Maximum multipole of the input cl cov[n,n,l] (double): Signal covariance matrix for each multipole, with bounds (0:n-1,0:n-1,0:lmax) iNcov[n,pix] (double): Inverse of the noise variance at each pixel, with bounds (0:n-1,0:npix-1) maps[n,pix] (double): Input kappa maps, with bouds (0:n-1,0:npix-1) - Args(optional):
chn (int): Number of grids for preconsitioner (chn=1 for diagonal preconditioner, default) lmaxs[chain] (int): Maximum multipole(s) at each preconditioning and lmaxs[0] is the input maximum multipole of cl nsides[chain] (int): Nside(s) of preconditoner and nsides[0] should be consistent with the input map’s nside. eps[chain] (double): Parameter to finish the iteration (i.e. terminate if the residul fraction becomes smaller than eps). Default to 1e-6. itns[chain] (int): Number of interations. inl[n,l] (double): Inverse noise spectrum for each mass map (0 for white noise case, default). verbose (bool): Output messages, default to False ro (int): the residual fraction is output for every ro iteration (e.g. ro=2 means 1 output per 2 iterations). Default to 50. Useful for convergence speed. stat (str): Realtime status filename which contains the residual fraction, default to no output file - Returns:
xlm[n,l,m] (dcmplx): Wiener filtered multipoles, with bounds (n,0:lmax,0:lmax) - Usage:
xlm = curvedsky.cninv.cnfilter_kappa(n,nside,lmax,cov,iNcov,maps,chn,lmaxs,nsides,itns,eps,inl,verbose,ro,stat):
-
curvedsky.cninv.
cnfilter_freq_nside
(n, mn0, mn1, nside0, nside1, lmax, cl, bl0, bl1, iNcov0, iNcov1, maps0, maps1, chn=1, lmaxs=[0], nsides0=[0], nsides1=[0], itns=[1], eps=[1e-06], filter='W', verbose=False, reducmn=0, ro=50, stat='', inl=None)¶ Same as cnfilter_freq but for the maps with two different Nsides. Please make sure that your input maps are beam-convolved. This code deconvolves the beam during filtering and the output are the filtered alms after the beam-deconvolution.
- Args:
n (int): Number of maps, i.e., temperature only (n=1), polarization only (n=2) or both (n=3) mn0/1 (int): Number of frequencies nside0/1 (int): Nsides of input maps lmax (int): Maximum multipole of the input cl cl[n,l] (double): Theory signal power spectrum, with bounds (0:n-1,0:lmax) bl0/1[mn,l] (double): Beam function, with bounds (0:n-1,0:lmax) iNcov0/1[n,mn,pix] (double): Inverse of the noise variance at each pixel, with bounds (0:n-1,0:npix-1) maps0/1[n,mn,pix] (double): Beam-convolved T, Q, U maps, with bouds (0:n-1,0:npix-1) - Args(optional):
chn (int): Number of grids for preconsitioner (chn=1 for diagonal preconditioner, default) lmaxs[chain] (int): Maximum multipole(s) at each preconditioning and lmaxs[0] is the input maximum multipole of cl nsides0/1[chain] (int): Nside(s) of preconditoner and nsides[0] should be consistent with the input map’s nside. eps[chain] (double): Parameter to finish the iteration (i.e. terminate if the residul fraction becomes smaller than eps). Default to 1e-6. itns[chain] (int): Number of interations. filter (str): C-inverse (‘’) or Wiener filter (W), default to the Wiener filter. inl[n,mn,l] (double): Inverse noise spectrum, 0 for white noise case. verbose (bool): Output messages, default to False ro (int): the residual fraction is output for every ro iteration (e.g. ro=2 means 1 output per 2 iterations). Default to 50. Useful for convergence speed. stat (str): Realtime status filename which contains the residual fraction, default to no output file reducmn (int): Reducing number of maps per chain (1,2) or not (0, default). If 1, the maps are combined for the same nside inside the multigrid chain. If 2, in addition to the procedure of 1, the each nside maprs are further combined into a single map inside the second chain (chain>=3). - Returns:
xlm[n,l,m] (dcmplx): C-inverse or Wiener filtered multipoles, with bounds (0:n-1,0:lmax,0:lmax)
mgc%cv(1,mi)%nij = iNcov0(:,mi,:)maps0(:,mi,:) mgc%cv(1,mi)%nij = iNcov1(:,mi-mn0,:)maps1(:,mi-mn0,:)
- Usage:
xlm = curvedsky.cninv.cnfilter_freq_nside(n,mn0,mn1,nside0,nside1,lmax,cl,bl0,bl1,iNcov0,iNcov1,maps0,maps1,chn,lmaxs,nsides0,nsides1,itns,eps,filter,inl,verbose,reducmn,ro,stat):
utils
– other tools¶
-
curvedsky.utils.
gauss1alm
(lmax, Cl)¶ Generating alm as a random Gaussian field whose power spectrum is cl. The output alm is given by a 2D array.
- Args:
lmax (int): Maximum multipole of the output alm Cl [l] (double): Angular power spectrum, with bounds (0:lmax) - Returns:
alm [l,m] (dcmplx): Random Gaussian alm, with bounds (0:lmax,0:lmax) - Usage:
alm = curvedsky.utils.gauss1alm(lmax,Cl):
-
curvedsky.utils.
gauss2alm
(lmax, cl1, cl2, xl)¶ Generating two alms as random Gaussian fields whose power spectra are cl1, cl2 and the cross spectrum is xl.
- Args:
lmax (int): Maximum multipole of the output alm cl1 [l] (double): Angular power spectrum of the 1st alm, with bounds (0:lmax) cl2 [l] (double): Angular power spectrum of the 2nd alm, with bounds (0:lmax) xl [l] (double): Cross-power spectrum between alm1 and alm2, with bounds (0:lmax) - Returns:
alm [2,l,m] (dcmplx): Random Gaussian alms, with bounds (2,0:lmax,0:lmax) - Usage:
alm = curvedsky.utils.gauss2alm(lmax,cl1,cl2,xl):
-
curvedsky.utils.
gauss2alm_const
(lmax, cl1, cl2, xl, flm)¶ Generating two alms as random Gaussian fields whose power spectra are cl1, cl2 and the cross spectrum is xl.
- Args:
lmax (int): Maximum multipole of the output alm cl1 [l] (double): Angular power spectrum of the 1st alm, with bounds (0:lmax) cl2 [l] (double): Angular power spectrum of the 2nd alm, with bounds (0:lmax) xl [l] (double): Cross-power spectrum between alm1 and alm2, with bounds (0:lmax) flm [l,m] (dcmplx): Constrained realiation of alms whose power is cl1 - Returns:
alm [2,l,m] (dcmplx): Random Gaussian alms, with bounds (2,0:lmax,0:lmax) - Usage:
alm = curvedsky.utils.gauss2alm_const(lmax,cl1,cl2,xl,flm):
-
curvedsky.utils.
gaussTEB
(lmax, TT, EE, BB, TE)¶ Generating T/E/B alms as random Gaussian fields whose power spectra are TT, EE, BB and the cross spectrum is TE.
- Args:
lmax (int): Maximum multipole of the output alms TT [l] (double): Angular power spectrum of temperature, with bounds (0:lmax) EE [l] (double): Angular power spectrum of E mode, with bounds (0:lmax) BB [l] (double): Angular power spectrum of B mode, with bounds (0:lmax) TE [l] (double): TE cross-power spectrum, with bounds (0:lmax) - Returns:
alm [3,l,m] (dcmplx): Random Gaussian T/E/B alms, with bounds (3,0:lmax,0:lmax) - Usage:
alm = curvedsky.utils.gaussTEB(lmax,TT,EE,BB,TE):
-
curvedsky.utils.
gaussalm
(cl, n=None, lmax=None, ilm=None)¶ Generating alms as random Gaussian fields whose covariance is given by cl[i,j].
- Args:
cl [i,j,l] (double): Covariance between the gaussian fields, with bounds (n,n,0:lmax) - Args:
ilm [l,m] (dcmplx): Input alm for the cl[0,0] element (default to None). The other alms are generated to be correlated with ilm. - Returns:
alm [i,l,m] (dcmplx): Random Gaussian alms, with bounds (n,0:lmax,0:lmax) - Usage:
alm = curvedsky.utils.gaussalm(n,lmax,cl,ilm):
-
curvedsky.utils.
cov_diag_tri
(n, cl)¶ - Usage:
A = curvedsky.utils.cov_diag_tri(n,cl):
-
curvedsky.utils.
get_baseline
(npix, nside_subpatch, QU)¶ Calculate baseline of each subpatch. The subpatches have the same size. Written by Ryo Nagata.
- Args:
npix (int): pixel number of the full map nside_subpatch (int): Nside of sub patch QU [pix,2] (double): Q/U maps, with bounds (0:npix-1,2) - Returns:
blmap [pix,2] (double): baseline maps, with bounds (0:npix-1,2) - Usage:
blmap = curvedsky.utils.get_baseline(npix,nside_subpatch,QU):
-
curvedsky.utils.
get_winmap
(nside_large, nside_small, ipix_pix, apod)¶ Return apodization window for subpatch. Written by Ryo Nagata.
- Args:
nside_large (int): Nside of sub patch nside_small (int): full Nside ipix_pix (int): pixel index of full map apod (double): apodization length - Returns:
wind_out (double): aporization window at ipix_pix - Usage:
win_out = curvedsky.utils.get_winmap(nside_large,nside_small,ipix_pix,apod):
-
curvedsky.utils.
get_apod_window
(s, a)¶ A sine apodization window
- Args:
s (double): Distance from the center of the window a (double): Apodization length, nothing (a=1) to all (a=0) - Returns:
w (double): Aporization window - Usage:
w = curvedsky.utils.get_apod_window(s,a):
-
curvedsky.utils.
eb_separate
(nside, lmax, W, Q, U)¶ E/B mode seperation based on the chi-field estimator. See e.g. Sec.III.2 of arXiv:1305.7441 for numerical implimentation.
- Args:
npix (int): Pixel number of the desired map lmax (int): Maximum multipole used for the harmonic transform internally W[pix] (double): Window function satisfying zero 1st and 2nd derivatives at the boundary, with bounds (0:npix-1) Q/U[pix] (double): Input Q/U map already multiplied by W, with bounds (0:npix-1) - Returns:
Elm/Blm[l,m] (dcmplx): Seperated E/B modes, with bounds (0:lmax,0:lmax) - Usage:
Elm,Blm = curvedsky.utils.eb_separate(nside,lmax,W,Q,U):
-
curvedsky.utils.
alm2cl
(lmax, alm1, alm2=None)¶ From alm to angular power spectrum
- Args:
lmax (int): Maximum multipole of the input alm alm1 [l,m] (dcmplx): 1st harmonic coefficient, with bounds (0:lmax,0:lmax) - Args(optional):
alm2 [l,m] (dcmplx): 2nd harmonic coefficient, with bounds (0:lmax,0:lmax), default to alm1 - Returns:
cl [l] (double): Auto or cross angular power spectrum, with bounds (0:lmax) - Usage:
cl = curvedsky.utils.alm2cl(lmax,alm1,alm2):
-
curvedsky.utils.
alm2bcl
(bn, lmax, alm1, spc='', alm2=None)¶ From alm to angular power spectrum with multipole binning
- Args:
bn (int): Number of multipole bins lmax (int): maximum multipole of the input alm alm1 [l,m] (dcmplx): 1st harmonic coefficient, with bounds (0:lmax,0:lmax) - Args(optional):
alm2 [l,m] (dcmplx): 2nd harmonic coefficient, with bounds (0:lmax,0:lmax), default to alm1 spc (str): Specify bin spacing, ‘’ = linear (default), ‘log’ = log spacing, ‘log10’ = log10 spacing, ‘p2’ = power of 2 spacing, ‘p3’ = power of 3 spacing - Returns:
cb [bin] (double): Auto or cross angular power spectrum with multipole binning, with bounds (0:bn-1) - Usage:
cb = curvedsky.utils.alm2bcl(bn,lmax,alm1,alm2,spc):
-
curvedsky.utils.
alm2rho
(lmax, alm1, alm2)¶ Compute correlation coefficients between two alms
- Args:
lmax (int): Maximum multipole of the input alm alm1 [l,m] (dcmplx): 1st harmonic coefficient, with bounds (0:lmax,0:lmax) alm2 [l,m] (dcmplx): 2nd harmonic coefficient, with bounds (0:lmax,0:lmax) - Returns:
rho [l] (double): Auto or cross angular power spectrum, with bounds (0:lmax) - Usage:
rho = curvedsky.utils.alm2rho(lmax,alm1,alm2):
-
curvedsky.utils.
alm2cov
(alm, n=0, lmax=0)¶ Compute correlation coefficients between two alms
- Args:
lmax (int): Maximum multipole of the input alms alm [n,l,m] (dcmplx): 1st harmonic coefficient, with bounds (0:lmax,0:lmax) - Returns:
cov [n,n,l] (double): Auto and cross angular power spectra between alm[i] and alm[j] - Usage:
cov = curvedsky.utils.alm2cov(n,lmax,alm):
-
curvedsky.utils.
apodize
(nside, rmask, ascale, order=1, holeminsize=0)¶ Compute apodized window function. Partially use Healpix’s process_mask code.
- Args:
nside (int): Nside of the input map rmask[pix] (double): Input window function, with bounds (0:pix-1). Pixels at rmask=0 is considered as masked pixels. ascale (double): Apodization length [deg] from the closest masked pixel - Args(optional):
order (int): Pixel order, 1 for RING (default), otherwize NESTED holeminsize (double): Minimum hole size [arcmin] (i.e., holes within this size is filled), default to 0 - Returns:
amask[pix] (double): Apodization window, with bounds (0:npix-1), using the same ordering as input - Usage:
amask = curvedsky.utils.apodize(nside,rmask,ascale,order,holeminsize):
-
curvedsky.utils.
almxfl
(lmax, mmax, alm, cl)¶ Calculate xlm[l,m] = alm[l,m] x cl[l]
- Args:
lmax (int): Maximum multipole of the input alm mmax (int): Maximum m of the input alm alm [l,m] (dcmplx): Harmonic coefficient to be transformed to a map, with bounds (0:lmax,0:lmax) cl [l] (double): 1D spectrum to be multiplied to alm, with bounds (0:lmax) - Returns:
xlm [l,m] (dcmplx): Modified alm, with bounds (0:lmax,0:lmax) - Usage:
xlm = curvedsky.utils.almxfl(lmax,mmax,alm,cl):
-
curvedsky.utils.
hp_alm2map
(nside, lmax, mmax, alm)¶ Ylm transform of the map to alm with the healpix (l,m) order
- Args:
nside (int): Nside of the input map lmax (int): Maximum multipole of the input alm mmax (int): Maximum m of the input alm alm [l,m] (dcmplx): Harmonic coefficient to be transformed to a map, with bounds (0:lmax,0:mmax) - Returns:
map [pix] (double): Transformed map, with bounds (0:npix-1) - Usage:
map = curvedsky.utils.hp_alm2map(nside,lmax,mmax,alm):
-
curvedsky.utils.
hp_alm2map_spin
(nside, lmax, mmax, spin, elm, blm)¶ Ylm transform of the map to alm with the healpix (l,m) order
- Args:
nside (int): Nside of the input map lmax (int): Maximum multipole of the input alm mmax (int): Maximum m of the input alm spin (int): Spin of the transform elm [l,m] (dcmplx): Spin-s E-like harmonic coefficient to be transformed to a map, with bounds (0:lmax,0:mmax) blm [l,m] (dcmplx): Spin-s B-like harmonic coefficient to be transformed to a map, with bounds (0:lmax,0:mmax) - Returns:
map0 [pix] (double): Real part of the transformed map (Q-like map), with bounds (0:npix-1) map1 [pix] (double): Imaginary part of the transformed map (U-like map), with bounds (0:npix-1) - Usage:
map0,map1 = curvedsky.utils.hp_alm2map_spin(nside,lmax,mmax,spin,elm,blm):
-
curvedsky.utils.
hp_map2alm
(nside, lmax, mmax, map)¶ Ylm transform of the map to alm with the healpix (l,m) order
- Args:
nside (int): Nside of the input map lmax (int): Maximum multipole of the input alm mmax (int): Maximum m of the input alm map [pix] (double): Input map, with bounds (0:npix-1) - Returns:
alm [l,m] (dcmplx): Harmonic coefficient obtained from the input map, with bounds (0:lmax,0:mmax) - Usage:
alm = curvedsky.utils.hp_map2alm(nside,lmax,mmax,map):
-
curvedsky.utils.
hp_map2alm_spin
(nside, lmax, mmax, spin, map0, map1)¶ Spin Ylm transform of the map ( = map0 + i map1 ) to alm with the healpix (l,m) order. For example, if map0=Q, map1=U and spin=2, the alm contains E-mode and B-mode.
- Args:
nside (int): Nside of the input map lmax (int): Maximum multipole of the input alm mmax (int): Maximum m of the input alm spin (int): Spin of the transform map0 [pix] (double): Real part of the input map, with bounds (0:npix-1) map1 [pix] (double): Imaginary part of the input map, with bounds (0:npix-1) - Returns:
alm [2,l,m] (dcmplx): Parity-eve/odd harmonic coefficients obtained from the input map, with bounds (2,0:lmax,0:mmax) - Usage:
alm = curvedsky.utils.hp_map2alm_spin(nside,lmax,mmax,spin,map0,map1):
-
curvedsky.utils.
map_mul_lfunc
(nside, imap, lmax, lfunc)¶ Convert map to alm, multiply a function to alm and convert back again to map
- Args:
nside (int): Nside of input map imap [pix] (double): Input map, with bounds (0:12*nside**2-1) lmax (int): Maximum multipole of the input alm lfunc [l] (double): 1D spectrum to be multiplied to alm, with bounds (0:lmax) - Returns:
omap [pix] (double): Output map, with bounds (0:12*nside**2-1) - Usage:
omap = curvedsky.utils.map_mul_lfunc(nside,imap,lmax,lfunc):
-
curvedsky.utils.
mulwin
(alm, win, nside=0, lmax=0, mmax=0)¶ Multiply window to a map obtained from alm
- Args:
alm [l,m] (dcmplx): Harmonic coefficient to be multiplied at window, with bounds (0:lmax,0:mmax) win [pix] (double): Transformed map, with bounds (0:npix-1) - Returns:
wlm [l,m] (dcmplx): Harmonic coefficient of the window-multiplied map, with bounds (0:lmax,0:mmax) - Usage:
wlm = curvedsky.utils.mulwin(nside,lmax,mmax,alm,win):
-
curvedsky.utils.
mulwin_spin
(elm, blm, win, nside=0, lmax=0, mmax=0, spin=2)¶ Multiply window to a map obtained from alm
- Args:
elm [l,m] (dcmplx): Spin-s E-like harmonic coefficient to be transformed to a map, with bounds (0:lmax,0:mmax) blm [l,m] (dcmplx): Spin-s B-like harmonic coefficient to be transformed to a map, with bounds (0:lmax,0:mmax) win [pix] (double): Transformed map, with bounds (0:npix-1) - Args (Optional):
spin (int): Spin of the transform, default to 2 - Returns:
wlm [2,l,m] (dcmplx): Parity-eve/odd harmonic coefficients obtained from the window-multiplied map, with bounds (2,0:lmax,0:mmax) - Usage:
wlm = curvedsky.utils.mulwin_spin(nside,lmax,mmax,spin,elm,blm,win):
-
curvedsky.utils.
lmpy2lmax
(lmpy)¶ - Usage:
lmax = curvedsky.utils.lmpy2lmax(lmpy):
-
curvedsky.utils.
lm_healpy2healpix
(almpy, lmax, lmpy=0)¶ Transform healpy alm to healpix alm
- Args:
lmax (int): Maximum multipole of the input/output alm satisfying 2 x lmpy = (lmax+1) x (lmax+2) almpy[index] (dcmplx): Healpy alm, with bounds (0:lmpy-1) - Returns:
almpix [l,m] (dcmplx): Healpix alm, with bounds (0:lmax,0:lmax) - Usage:
almpix = curvedsky.utils.lm_healpy2healpix(lmpy,almpy,lmax):
-
curvedsky.utils.
cosin_healpix
(nside)¶ Return cos(theta) as a function of the Healpix pixel index
- Args:
nside (int): Nside of the desired map - Returns:
cosin[pix] (double): cosin(theta), with bounds (0:npix-1) - Usage:
cosin = curvedsky.utils.cosin_healpix(nside):
-
curvedsky.utils.
load_defangle_takahashi
(fname, npix, verbose=False)¶ Read theta and phi coordinates at source plane obtained by Takahashi et al. (2017)
- Args:
fname (str): file name npix (int): Number of pixels of theta and phi data - Args (optional):
verbose (bool): output messages, default to False - Returns:
theta[pix] (double): theta, with bounds (0:npix-1) phi[pix] (double): phi, with bounds (0:npix-1) - Usage:
theta,phi = curvedsky.utils.load_defangle_takahashi(fname,npix,verbose):
-
curvedsky.utils.
polcoord2angle
(npix, theta, phi, verbose=False)¶ Converting theta and phi coordinates at source plane to deflection angle. The algorithm is provided by Takashi Hamana and Ryuichi Takahashi.
- Args:
npix (int): Number of pixels of theta and phi data theta[pix] (double): theta, with bounds (0:npix-1) phi[pix] (double): phi, with bounds (0:npix-1) - Args (optional):
verbose (bool): output messages, default to False - Returns:
angle[pix,2] (double): deflection angle vector containing two components, with bounds (0:npix-1,1:2) - Usage:
angle = curvedsky.utils.polcoord2angle(npix,theta,phi,verbose):
-
curvedsky.utils.
polcoord2angle_alm
(nside, lmax, theta, phi, verbose=False)¶ Converting theta and phi coordinates at source plane to deflection angle.
- Args:
npix (int): Number of pixels of theta and phi data lmax (int): Maximum multipole of alms for gradient and curl modes theta[pix] (double): theta, with bounds (0:npix-1) phi[pix] (double): phi, with bounds (0:npix-1) - Args (optional):
verbose (bool): output messages, default to False - Returns:
glm[l,m] (dcmplx): gradient mode, with bounds (0:lmax,0:lmax) clm[l,m] (dcmplx): curl mode, with bounds (0:lmax,0:lmax) - Usage:
glm,clm = curvedsky.utils.polcoord2angle_alm(nside,lmax,theta,phi,verbose):
-
curvedsky.utils.
calc_mfs
(bn, nu, lmax, walm, nside=0)¶ Compute 2D Minkowski functionals from a given alm
- Args:
bn (int): Number of nu bins nu [bin] (double): Nu bins, with bounds (bn) lmax (int): Maximum multipole of the input walm walm [l,m] (dcmplx): Alm with filtering, possibly divided by the map variance, with bounds (0:lmax,0:lmax) - Args(optional):
nside (int): Nside of the intermediate map, default to lmax - Returns:
V [bin,type] (double): The three Minkowski functionals, V0, V1 and V2, at each nu bin, with bounds (bn,0:2) - Usage:
V = curvedsky.utils.calc_mfs(bn,nu,lmax,walm,nside):
-
curvedsky.utils.
mock_galaxy_takahashi
(fname, zn, ngz, zi, a=2.0, b=1.0, zm=1.0, sz=0.0, zbias=0.0, b0=1.0, btype='sqrtz')¶ Compute galaxy overdensity map from dark matter density map of Takahashi et al. (2017). The galaxy z distribution is assumed to have the functional form given by Eq.(7) of https://arxiv.org/abs/1810.03346
- Args:
fname (str): Filename of density map zn (int): Number of tomographic bins ngz[zn] (double): Total number of galaxies at each z-bin zi[zn+1] (double): Redshift intervals of the tomographic bins - Args(optional):
a (double): galaxy distribution shape parameter b (double): galaxy distribution shape parameter zm (double): mean redshift of galaxy distribution b0 (double): constant galaxy bias at z=0 - Returns:
gmap [pix,zbin] (double): The galaxy number density map at each zbin - Usage:
gmap = curvedsky.utils.mock_galaxy_takahashi(fname,zn,ngz,zi,b0,btype,a,b,zm,sz,zbias):