so_cov - a module for covariance matrix estimation#

Tools for analytical covariance matrix estimation. For more details on computation of the matrix see https://pspy.readthedocs.io/en/latest/scientific_doc.pdf.

pspy.so_cov.bin_mat(mat, binning_file, lmax, speclist=['TT'])#

Take a matrix and bin it Mbb’= Pbl Pb’l’ Mll’ with Pbl =1/Nb sum_(l in b)

Parameters:
  • binning_file (data file) – a binning file with format bin low, bin high, bin mean

  • lmax (int) – the maximum multipole to consider

pspy.so_cov.block_diagonal_mult(slices, mbb_inv_ab_list, mbb_inv_cd_list, analytic_cov)#

Suggestion by adrien to do an operation of the type A M B, where A and B are block diagonal matrix

Parameters:
  • slices (list of tuple of integer) – give the min and max indices of each block

  • mbb_inv_ab_list (list of 2d array) – list of the inverse of the mode coupling matrix

  • mbb_inv_cd_list (list of 2d array) – list of the transpose of the inverse of the mode coupling matrix

  • analytic_cov (2d array) – analytic covariance matrix

pspy.so_cov.chi(alpha, gamma, beta, eta, ns, Dl, DNl, id='TTTT')#

doc not ready yet

pspy.so_cov.corr2cov(corr, var)#

Go from correlation and variance to covariance matrix

Parameters:
  • corr (2d array) – the correlation matrix

  • var (1d array) – vector of variance of the random variables

pspy.so_cov.cov2corr(cov, remove_diag=False)#

Go from covariance to correlation matrix, also setting the diagonal to zero

Parameters:

cov (2d array) – the covariance matrix

pspy.so_cov.cov_coupling_spin0(win, lmax, niter=3, save_file=None, l_exact=None, l_band=None, l_toep=None)#

Compute the coupling kernels corresponding to the T only covariance matrix

Parameters:
  • win (so_map or dictionnary of so_map) – the window functions, can be a so_map or a dictionnary containing so_map if the later, the entry of the dictionnary should be Ta,Tb,Tc,Td.

  • lmax (integer) – the maximum multipole to consider

  • niter (int) – the number of iteration performed while computing the alm

  • save_file (string) – the name of the file in which the coupling kernel will be saved (npy format)

  • l_toep (int) – parameter for the toeplitz approximation (see arXiv:2010.14344)

  • l_band (int) – parameter for the toeplitz approximation (see arXiv:2010.14344)

  • l_exact (int) – parameter for the toeplitz approximation (see arXiv:2010.14344)

pspy.so_cov.cov_coupling_spin0and2_simple(win, lmax, niter=3, save_file=None, planck=False, l_exact=None, l_band=None, l_toep=None)#

Compute the coupling kernels corresponding to the T and E covariance matrix

Parameters:
  • win (so_map or dictionnary of so_map) – the window functions, can be a so_map or a dictionnary containing so_map if the later, the entry of the dictionnary should be Ta,Tb,Tc,Td,Pa,Pb,Pc,Pd

  • lmax (integer) – the maximum multipole to consider

  • niter (int) – the number of iteration performed while computing the alm

  • save_file (string) – the name of the file in which the coupling kernel will be saved (npy format)

  • l_toep (int) – parameter for the toeplitz approximation (see arXiv:2010.14344)

  • l_band (int) – parameter for the toeplitz approximation (see arXiv:2010.14344)

  • l_exact (int) – parameter for the toeplitz approximation (see arXiv:2010.14344)

pspy.so_cov.cov_spin0(Clth_dict, coupling_dict, binning_file, lmax, mbb_inv_ab, mbb_inv_cd, binned_mcm=True)#

From the two point functions and the coupling kernel construct the spin0 analytical covariance matrix of <(C_ab- Clth)(C_cd-Clth)>

Parameters:
  • Clth_dict (dictionnary) – A dictionnary of theoretical power spectrum (auto and cross) for the different split combinaison (‘TaTb’ etc)

  • coupling_dict (dictionnary) – a dictionnary containing the coupling kernel

  • binning_file (data file) – a binning file with format bin low, bin high, bin mean

  • lmax (int) – the maximum multipole to consider

  • mbb_inv_ab (2d array) – the inverse mode coupling matrix for the ‘TaTb’ power spectrum

  • mbb_inv_cd (2d array) – the inverse mode coupling matrix for the ‘TcTd’ power spectrum

  • binned_mcm (boolean) – specify if the mode coupling matrices are binned or not

pspy.so_cov.cov_spin0_aniso1(Clth, Rl, couplings, binning_file, lmax, mbb_inv_ab, mbb_inv_cd, binned_mcm=True)#
Estimate the analytic covariance in the case of anisotropic noise pixel variance,

using the theoretical Cls, the ratios of noise to white noise power spectra, the coupling matrices estimated from the masks and variance maps, and the mode-coupling matrices.

Parameters:
  • Clth (dictionary) – A dictionary of theoretical power spectrum (auto and cross) for the different split combinations (‘TaTb’ etc)

  • Rl (dictionary) – a dictionary containing the ratios of noise power spectra to white noise

  • couplings (list of arrays) – a list containing the eight coupling matrices involved in this covariance

  • binning_file (data file) – a binning file with format bin low, bin high, bin mean

  • lmax (int) – the maximum multipole to consider

  • mbb_inv_ab (2d array) – the inverse mode coupling matrix for the ‘TaTb’ power spectrum

  • mbb_inv_cd (2d array) – the inverse mode coupling matrix for the ‘TcTd’ power spectrum

  • binned_mcm (boolean) – specify if the mode coupling matrices are binned or not

pspy.so_cov.cov_spin0and2(Clth_dict, coupling_dict, binning_file, lmax, mbb_inv_ab, mbb_inv_cd, binned_mcm=True, cov_T_E_only=True, dtype=<class 'numpy.float64'>, old_way=False)#

From the two point functions and the coupling kernel construct the T and E analytical covariance matrix of <(C_ab- Clth)(C_cd-Clth)>

Parameters:
  • Clth_dict (dictionnary) – A dictionnary of theoretical power spectrum (auto and cross) for the different split combinaison (‘XaYb’ etc)

  • coupling_dict (dictionnary) – a dictionnary containing the coupling kernel

  • binning_file (data file) – a binning file with format bin low, bin high, bin mean

  • lmax (int) – the maximum multipole to consider

  • mbb_inv_ab (2d array) – the inverse mode coupling matrix for the ‘XaYb’ power spectrum

  • mbb_inv_cd (2d array) – the inverse mode coupling matrix for the ‘XcYd’ power spectrum

  • binned_mcm (boolean) – specify if the mode coupling matrices are binned or not

  • binned_mcm – specify if the mode coupling matrices are binned or not

  • cov_T_E_only (boolean) – if true don’t do B

pspy.so_cov.covariance_element_beam(id_element, ps_all, norm_beam_cov, binning_file, lmax, cov_T_E_only=True)#

This routine compute the contribution from beam errors to the analytical covariance of the power spectra We want to compute the beam covariance between the two spectra C1 = Wa * Xb, C2 = Yc * Zd Here W, X, Y, Z can be either T or E and a,b,c,d will be an index corresponding to the survey and array we consider so for example a = dr6&pa5_150 or a = dr6&pa4_090 The formula for the analytic covariance of C1, C2 is given by let’s denote the normalised beam covariance <BB>_ac = < delta B_a delta B_c >/np.outer(B_a, B_c)

Cov(Wa * Xb, Yc * Zd) = Dl^{WaXb} Dl^{YcZd} ( <BB>_ac + <BB>_ad + <BB>_bc + <BB>_bd )

Parameters:
  • id_element (list) – a list of the form [a,b,c,d] where a = dr6_pa4_090, etc, this identify which pair of power spectrum we want the covariance of

  • ps_all (dict) – this dict contains the theoretical best power spectra, convolve with the beam for example ps[“dr6&pa5_150”, “dr6&pa4_150”, “TT”] = (Dl^{CMB}_TT + fg_TT)

  • norm_beam_cov (dict) – this dict contains the normalized beam covariance for each survey and array

  • binning_file (str) – a binning file with three columns bin low, bin high, bin mean

  • lmax (int) – the maximum multipole to consider

  • cov_T_E_only (boolean) – if true don’t do B

pspy.so_cov.delta2(a, b)#

Simple delta function

pspy.so_cov.delta3(a, b, c)#

Delta function (3 variables)

pspy.so_cov.delta4(a, b, c, d)#

Delta function (4 variables)

pspy.so_cov.extract_EEEBBB_mbb(mbb_inv)#

this routine extract the E and B part of the mode coupling matrix

Parameters:

mbb_inv (2d array) – the inverse spin0 and 2 mode coupling matrix

pspy.so_cov.extract_mbb(mbb_inv, cov_T_E_only=True, dtype=<class 'numpy.float64'>)#

The mode coupling marix is computed for T,E,B but for if cov_T_E_only=True we only construct analytical covariance matrix for T and E

Parameters:
  • mbb_inv (dict of 2d arrays) – the inverse spin0 and 2 mode coupling matrix

  • cov_T_E_only (boolean) – if true don’t do B

pspy.so_cov.extract_mbb_list(mbb_inv, cov_T_E_only=True, dtype=<class 'numpy.float64'>, transpose=False)#

extract a list of inverse mode coupling matrix, you can optionnaly choose the dtype, if you want only the TT-TE-ET-EE part, and if you want to transpose it

Parameters:
  • mbb_inv (dict of 2d arrays) – the inverse spin0 and 2 mode coupling matrix

  • cov_T_E_only (boolean) – if true don’t do B

  • dtype (np dtype) – choose the type to save memory

  • transpose (boolean) – wether to transpose it or not

pspy.so_cov.f(a, b, c, d, ns)#

f combination factor in the covariance computation

pspy.so_cov.fast_cov_coupling_spin0and2(sq_win_alms_dir, id_element, lmax, l_exact=None, l_band=None, l_toep=None)#

Compute the coupling kernels corresponding to the T and E covariance matrix This routine assume that you have already precomputed the alms of the square of the windows and that the window in temperature and polarisation are the same. Both conditions lead to a very important speed up, that explains the name of the routine

Parameters:
  • sq_win_alms_dir (string) – the folder with precomputed alms corresponding to the sq of the window function this is what enters the analytic computation

  • id_element (list) – a list of the form [a,b,c,d] where a = dr6_pa4_090, etc, this identify which pair of power spectrum we want the covariance of

  • lmax (integer) – the maximum multipole to consider

  • l_toep (int) – parameter for the toeplitz approximation (see arXiv:2010.14344)

  • l_band (int) – parameter for the toeplitz approximation (see arXiv:2010.14344)

  • l_exact (int) – parameter for the toeplitz approximation (see arXiv:2010.14344)

pspy.so_cov.g(a, b, c, d, ns)#

g combination factor in the covariance computation

pspy.so_cov.generalized_cov_spin0and2(coupling_dict, id_element, ns, ps_all, nl_all, lmax, binning_file, mbb_inv_ab, mbb_inv_cd, binned_mcm=True, return_full_cov=False, cov_T_E_only=True, dtype=<class 'numpy.float64'>, old_way=False)#

This routine deserves some explanation We want to compute the covariance between two power spectra C1 = Wa * Xb, C2 = Yc * Zd Here W, X, Y, Z can be either T or E and a,b,c,d will be an index corresponding to the survey and array we consider so for example a = s17&pa5_150 or a = dr6&pa4_090 The formula for the analytic covariance of C1, C2 is given by Cov( Wa * Xb, Yc * Zd) = < Wa Yc> <Xb Zd> + < Wa Zd> <Xb Yc> (this is just from the wick theorem) In practice we need to include the effect of the mask (so we have to introduce the coupling dict D) and we need to take into account that we use as spectra an average of cross power spectra, that is why we use the chi function (and what make it different from cov_spin0and2) Cov( Wa * Xb, Yc * Zd) = D(Wa*Yc,Xb Zd) chi(Wa,Yc,Xb Zd) + D(Wa*Zd,Xb*Yc) chi(Wa,Zd,Xb,Yc)

Parameters:
  • coupling_dict (dictionnary) – a dictionnary that countains the coupling terms arising from the window functions

  • id_element (list) – a list of the form [a,b,c,d] where a = dr6_pa4_090, etc, this identify which pair of power spectrum we want the covariance of

  • ns (dict) – this dictionnary contains the number of split we consider for each of the survey

  • ps_all (dict) – this dict contains the theoretical best power spectra, convolve with the beam for example ps[“dr6&pa5_150”, “dr6&pa4_150”, “TT”] = bl_dr6_pa5_150 * bl_dr6_pa4_150 * (Dl^{CMB}_TT + fg_TT)

  • nl_all (dict) – this dict contains the estimated noise power spectra, note that it correspond to the noise power spectrum per split e.g nl[“dr6&pa5_150”, “dr6&pa4_150”, “TT”]

  • binning_file – a binning file with three columns bin low, bin high, bin mean

  • mbb_inv_cd (mbb_inv_ab and) – the inverse mode coupling matrices corresponding to the C1 = Wa * Xb and C2 = Yc * Zd power spectra

  • binned_mcm (boolean) – specify if the mode coupling matrices are binned or not

  • return_full_cov (boolean) – an option to return the lbyl cov (if binned_mcm=False) mostly used for debugging

  • cov_T_E_only (boolean) – if true don’t do B

pspy.so_cov.generate_aniso_couplings(survey_name, win, var, lmax, niter=0)#

survey_name : list of names of the four splits win : dict with windows with keys ‘Ta’, “Tb’, … var : dict with variance maps, with keys ‘Ta’, ‘Tb’

pspy.so_cov.get_sigma(cov, spectra_order, n_bins, spectrum)#

get the error of the spectrum for the given cov mat

Parameters:
  • cov (2d array) – the covariance matrix

  • spectra_order (list of strings) – the arangement of the different block

  • n_bins (int) – the number of bins for each block

  • spectrum (string) – the spectrum we consider

pspy.so_cov.mc_cov_from_spectra_list(spec_list_a, spec_list_b, spectra=None)#

Return the montecarlo covariance of two spectra list estimated from simulations This routine is inefficient and need some more work

Parameters:
  • spec_list_a (first list of spectra (or list of dict if spectra is not None)) – a list of ps computed from simulation (each entry of the list is a ps)

  • spec_list_b (second list of spectra (or list of dict if spectra is not None)) – a list of ps computed from simulation (each entry of the list is a ps)

  • spectra (list of strings) – needed for spin0 and spin2 cross correlation, the arrangement of the spectra

pspy.so_cov.measure_white_noise_level(var, mask, dtype=<class 'numpy.float64'>)#

take enmaps of map variance and mask, outputs white noise power spectrum amplitude :param var: variance map :type var: enmap :param mask: window map :type mask: enmap

pspy.so_cov.plot_cov_matrix(mat, color_range=None, color='pwhite', file_name=None)#

plot the covariance matrix at full resolution using pixell plotting routines

Parameters:
  • mat (2d array) – the covariance matrix

  • color_range (float) – the range of the plot

  • color (pixell colormap) – the colormap for the plot (have to be pixell compatible)

  • file_name (string) – file_name is the name of the png file that will be created, if None the plot will be displayed.

pspy.so_cov.read_coupling(file, spectra=None)#

Read a precomputed coupling kernels the code use the size of the array to infer what type of survey it corresponds to

Parameters:

file (string) – the name of the npy file

pspy.so_cov.selectblock(cov, spectra_order, n_bins, block='TTTT')#

Select a block in a spin0 and 2 covariance matrix

Parameters:
  • cov (2d array) – the covariance matrix

  • spectra_order (list of strings) – the arangement of the different block

  • n_bins (int) – the number of bins for each block

  • block (string) – the block you want to look at

pspy.so_cov.symmetrize(Clth, mode='arithm')#

Take a power spectrum Cl and return a symmetric array C_l1l2=f(Cl)

Parameters:
  • Clth (1d array) – the power spectrum to be made symmetric

  • mode (string) – geometric or arithmetic mean if geo return C_l1l2 = sqrt( |Cl1 Cl2 |) if arithm return C_l1l2 = (Cl1 + Cl2)/2