Commit 58be8805 authored by Syl's avatar Syl
Browse files

few corrections

parent f1d0cc3b
......@@ -3,6 +3,7 @@ Test script for xQML
Author: Vanneste
"""
from __future__ import division
import timeit
......@@ -11,107 +12,135 @@ import healpy as hp
from pylab import *
import astropy.io.fits as fits
from xqml import xQML
from libcov import compute_ds_dcb
import xqml
import libcov
import simulation
from xqml_utils import progress_bar
from simulation import getstokes, muKarcmin2var, GetBinningMatrix
from simulation import muKarcmin2var, GetBinningMatrix
from simulation import extrapolpixwin
if __name__ == "__main__":
nside = 4
lmax = 3 * nside - 1
Slmax = 3 * nside - 1
deltal = 1
nsimu = 500
clth = np.array(hp.read_cl('planck_base_planck_2015_TTlowP.fits'))
lth = arange(2, lmax+1)
ellbins = arange(2, lmax + 2, deltal)
ellbins[-1] = lmax + 1
P, Q, ell, ellval = GetBinningMatrix(ellbins, lmax)
nbins = len(ellbins) - 1
# Create mask
t, p = hp.pix2ang(nside, range(hp.nside2npix(nside)))
mask = np.ones(hp.nside2npix(nside), bool)
mask[abs(90 - rad2deg(t)) < 10] = False
npix = sum(mask)
fwhm = 0.5
bl = hp.gauss_beam(deg2rad(fwhm), lmax=Slmax + 1)
allStoke, der, ind = getstokes(polar=True, temp=False, EBTB=False)
nder = len(der)
muKarcmin = 0.1
pixvar = muKarcmin2var(muKarcmin, nside)
varmap = ones((2 * npix)) * pixvar
NoiseVar = np.diag(varmap)
cmb = np.array(hp.synfast(
clth, nside, fwhm=deg2rad(fwhm), pixwin=True, new=True, verbose=False))
noise = (randn(len(varmap)) * varmap**0.5).reshape(2, -1)
# dm = cmb[1:, mask] + noise
dm = cmb[1:][:, mask] + noise
# ############## Compute ds_dcb ###############
ip = arange(hp.nside2npix(nside))
ipok = ip[mask]
Pl, S = compute_ds_dcb(
ellbins, nside, ipok, bl, clth, Slmax,
polar=True, temp=False, EBTB=False,
pixwining=True, timing=True, MC=False)
# Pl = Pl.reshape((nder)*(np.shape(Pl)[1]), 2 * npix, 2 * npix)
# ############## Compute spectra ###############
esti = xQML(mask, ellbins, clth, Pl=Pl, S=S, fwhm=fwhm)
esti.construct_esti(NoiseVar, NoiseVar)
cl = esti.get_spectra(dm, dm)
V = esti.get_covariance()
# ############## Construct MC ###############
allcl = []
esti = xQML(mask, ellbins, clth, Pl=Pl, fwhm=fwhm)
esti.construct_esti(NoiseVar, NoiseVar)
fpixw = extrapolpixwin(nside, lmax+2, pixwining=True)
start = timeit.default_timer()
for n in np.arange(nsimu):
ion()
# if __name__ == "__main__":
nside = 4
lmax = 2 * nside - 1
Slmax = 2 * nside - 1
deltal = 1
nsimu = 10000
clth = np.array(hp.read_cl('planck_base_planck_2015_TTlowP.fits'))
Clthshape = zeros(((6,)+shape(clth)[1:]))
Clthshape[:4] = clth
clth = Clthshape
EB = 0.5
clth[4] = EB*sqrt(clth[1]*clth[2])
TB = 0.5
clth[5] = TB*sqrt(clth[0]*clth[2])
lth = arange(2, lmax+1)
spec = ['EB', 'TE', 'TB']
temp = True
polar = True
corr = False
pixwin = False
ellbins = arange(2, lmax + 2, deltal)
ellbins[-1] = lmax + 1
P, Q, ell, ellval = GetBinningMatrix(ellbins, lmax)
nbins = len(ellbins) - 1
# Create mask
t, p = hp.pix2ang(nside, range(hp.nside2npix(nside)))
mask = np.ones(hp.nside2npix(nside), bool)
mask[abs(90 - rad2deg(t)) < 60] = False
npix = sum(mask)
fwhm = 0.5
bl = hp.gauss_beam(deg2rad(fwhm), lmax=Slmax + 1)
stokes, spec, istokes, ispecs = simulation.getstokes(
spec=spec, polar=polar, temp=temp, corr=corr)
print(stokes, spec, istokes, ispecs)
nspec = len(spec)
nstoke = len(stokes)
# ############## Compute ds_dcb ###############
ip = arange(hp.nside2npix(nside))
ipok = ip[mask]
Pl, S = libcov.compute_ds_dcb(ellbins, nside, ipok, bl, clth, Slmax, spec=spec,
pixwin=pixwin, timing=True, MC=False)
# ############## Compute spectra ###############
muKarcmin = 1.0
pixvar = muKarcmin2var(muKarcmin, nside)
varmap = ones((nstoke * npix)) * pixvar
NoiseVar = np.diag(varmap)
cmb = np.array(hp.synfast(clth, nside, fwhm=deg2rad(fwhm), pixwin=pixwin,
new=True, verbose=False, lmax=Slmax))
noise = (randn(len(varmap)) * varmap**0.5).reshape(nstoke, -1)
dm = cmb[istokes][:, mask] + noise
esti = xqml.xQML(mask, ellbins, clth, Pl=Pl, S=S, fwhm=fwhm,
spec=spec, temp=temp, polar=polar, corr=corr)
esti.construct_esti(NoiseVar, NoiseVar)
cl = esti.get_spectra(dm, dm)
V = esti.get_covariance()
# ############## Construct MC ###############
allcl = []
# allcmb = []
# esti = xqml.xQML(mask, ellbins, clth, Pl=Pl, fwhm=fwhm, spec=spec, temp=temp,
# polar=polar, corr=corr)
esti.construct_esti(NoiseVar, NoiseVar)
fpixw = extrapolpixwin(nside, lmax+2, pixwin=pixwin)
start = timeit.default_timer()
for n in np.arange(nsimu):
progress_bar(n, nsimu, timeit.default_timer() - start)
cmb = np.array(hp.synfast(clth[:, :len(fpixw)]*fpixw**2, nside,
fwhm=deg2rad(fwhm), new=True, verbose=False))
dm = cmb[1:, mask] + (randn(2 * npix) * sqrt(varmap)).reshape(2, npix)
allcl.append(esti.get_spectra(dm, dm))
figure()
subplot(3, 1, 1)
plot(lth, clth.transpose()[lth, 1: 3], '--k')
hcl = mean(allcl, 0).transpose()
scl = std(allcl, 0).transpose()
plot(ellval, hcl, 'b.')
plot(ellval, hcl + scl, 'r--', label=r"$\pm 1\sigma$")
plot(ellval, hcl - scl, 'r--')
ylabel(r"$C_\ell$")
semilogy()
legend(loc=4)
subplot(3, 1, 2)
cosmic = sqrt(2./(2 * lth + 1)) / mean(mask) * clth[1: 3, lth]
plot(lth, cosmic.transpose(), '--k')
plot(ellval, scl, 'r-', label=r"$\sigma_{\rm MC}$")
plot(ellval, sqrt(diag(V)).reshape(nder, -1).transpose(), 'b.')
ylabel(r"$\sigma(C_\ell)$")
semilogy()
legend(loc=4)
subplot(3, 1, 3)
plot(ellval, (hcl-clth[1: 3, lth].T)/(scl/sqrt(nsimu)), 'o')
ylabel(r"$R[C_\ell]$")
xlabel(r"$\ell$")
ylim(-3, 3)
show()
cmb = np.array(hp.synfast(clth[:, :len(fpixw)]*(fpixw*bl)**2, nside,
lmax=Slmax, fwhm=deg2rad(fwhm), new=True, verbose=False))
cmbm = cmb[istokes][:, mask]
dmA = cmbm + (randn(nstoke * npix) * sqrt(varmap)).reshape(nstoke, npix)
dmB = cmbm + (randn(nstoke * npix) * sqrt(varmap)).reshape(nstoke, npix)
# allcmb.append(cmbm)
allcl.append(esti.get_spectra(dmA, dmB))
figure(figsize=[10, 8])
clf()
subplot(3, 1, 1)
plot(lth, clth[ispecs][:, lth].T, '--k')
hcl = mean(allcl, 0)
scl = std(allcl, 0)
[plot(ellval, hcl[s], 'o', color='C%i' % s, label=r"$%s$" % spec[s])
for s in np.arange(nspec)]
[fill_between(ellval, (hcl - scl/sqrt(nsimu))[s], (hcl + scl/sqrt(nsimu))[s],
color='C%i' % s, alpha=0.2) for s in np.arange(nspec)]
ylabel(r"$C_\ell$")
semilogy()
legend(loc=4)
subplot(3, 1, 2)
cosmic = sqrt(2./(2 * lth + 1)) / mean(mask) * clth[ispecs][:, lth]
# plot(lth, cosmic.transpose(), '--k')
[plot(ellval, scl[s], '--', color='C%i' % s, label=r"$\sigma^{%s}_{\rm MC}$" %
spec[s]) for s in np.arange(nspec)]
[plot(ellval, sqrt(diag(V)).reshape(nspec, -1)[s], 'o', color='C%i' % s)
for s in np.arange(nspec)]
ylabel(r"$\sigma(C_\ell)$")
semilogy()
legend(loc=4)
subplot(3, 1, 3)
[plot(ellval, (hcl[s]-clth[ispecs[s]][lth])/(scl[s]/sqrt(nsimu)), '--o',
color='C%i' % s) for s in np.arange(nspec)]
ylabel(r"$R[C_\ell]$")
xlabel(r"$\ell$")
ylim(-3, 3)
grid()
show()
# savefig("../Plots/Git/"+"test0.png")
if __name__ == "__main__":
"""
......
......@@ -42,7 +42,7 @@ def Pl(ds_dcb):
return np.copy(ds_dcb).reshape(2 * (np.shape(ds_dcb)[1]), nnpix, nnpix)
def CorrelationMatrix(Clth, Pl, ellbins, polar=True, temp=False, EBTB=False):
def CorrelationMatrix(Clth, Pl, ellbins, polar=True, temp=False, corr=False):
"""
Compute correlation matrix S = sum_l Pl*Cl
......@@ -58,7 +58,7 @@ def CorrelationMatrix(Clth, Pl, ellbins, polar=True, temp=False, EBTB=False):
If True, get Stokes parameters for polar (default: True)
temp : bool
If True, get Stokes parameters for temperature (default: False)
EBTB : bool
corr : bool
If True, get Stokes parameters for EB and TB (default: False)
Returns
......@@ -76,7 +76,7 @@ def CorrelationMatrix(Clth, Pl, ellbins, polar=True, temp=False, EBTB=False):
[[ 0 280 560 840 1120]
[1400 1680 1960 2240 2520]]
"""
if EBTB:
if corr:
xx = ['TT', 'EE', 'BB', 'TE', 'EB', 'TB']
ind = [0, 1, 2, 3, 4, 5]
else:
......@@ -85,7 +85,7 @@ def CorrelationMatrix(Clth, Pl, ellbins, polar=True, temp=False, EBTB=False):
if not temp:
allStoke = ['Q', 'U']
if EBTB:
if corr:
xx = ['EE', 'BB', 'EB']
ind = [1, 2, 5]
else:
......
"""
Set of routines to ...
Author: Vanneste
"""
from __future__ import division
import numpy as np
def polrotangle(ri, rj):
"""
Computes cosine and sine of twice the angle between pixels i and j.
Parameters
----------
ri : 3D array of floats
Coordinates of vector corresponding to input pixels i following
healpy.pix2vec(nside,ipix) output
rj : 3D array of floats
Coordinates of vector corresponding to input pixels j following
healpy.pix2vec(nside,jpix) output
Returns
----------
cos2a : 1D array of floats
Cosine of twice the angle between pixels i and j
sin2a : 1D array of floats
Sine of twice the angle between pixels i and j
Example
----------
>>> cos2a, sin2a = polrotangle([0.1,0.2,0.3], [0.4,0.5,0.6])
>>> print(round(cos2a,5),round(sin2a,5))
(0.06667, 0.37333)
"""
z = np.array([0.0, 0.0, 1.0])
# Compute ri^rj : unit vector for the great circle connecting i and j
rij = np.cross(ri, rj)
norm = np.sqrt(np.dot(rij, np.transpose(rij)))
# case where pixels are identical or diametrically opposed on the sky
if norm <= 1e-15:
cos2a = 1.0
sin2a = 0.0
return cos2a, sin2a
rij = rij / norm
# Compute z^ri : unit vector for the meridian passing through pixel i
ris = np.cross(z, ri)
norm = np.sqrt(np.dot(ris, np.transpose(ris)))
# case where pixels is at the pole
if norm <= 1e-15:
cos2a = 1.0
sin2a = 0.0
return cos2a, sin2a
ris = ris / norm
# Now, the angle we want is that
# between these two great circles: defined by
cosa = np.dot(rij, np.transpose(ris))
# the sign is more subtle : see tegmark et de oliveira costa 2000 eq. A6
rijris = np.cross(rij, ris)
sina = np.dot(rijris, np.transpose(ri))
# so now we have directly cos2a and sin2a
cos2a = 2.0 * cosa * cosa - 1.0
sin2a = 2.0 * cosa * sina
return cos2a, sin2a
if __name__ == "__main__":
"""
Run the doctest using
python simulation.py
If the tests are OK, the script should exit gracefuly, otherwise the
failure(s) will be printed out.
"""
import doctest
if np.__version__ >= "1.14.0":
np.set_printoptions(legacy="1.13")
doctest.testmod()
......@@ -9,19 +9,20 @@ import timeit
import numpy as np
import healpy as hp
import math
from scipy import special
from spin_functions import dlss, pl0
from spin_functions import F1l2, F2l2
from libangles import polrotangle
# from spin_functions import dlss, pl0
# from spin_functions import F1l2, F2l2
# from libangles import polrotangle
from simulation import getstokes
from simulation import extrapolpixwin
from xqml_utils import progress_bar
def compute_ds_dcb(
ellbins, nside, ipok, bl, clth, Slmax,
polar=True, temp=True, EBTB=False,
pixwining=False, timing=False, MC=0, Sonly=False):
ellbins, nside, ipok, bl, clth, Slmax, spec,
pixwin=False, timing=False, MC=0, Sonly=False):
"""
Compute the Pl = dS/dCl matrices.
......@@ -40,13 +41,9 @@ def compute_ds_dcb(
Fiducial power spectra
Slmax : int
Maximum lmax computed for the pixel covariance pixel matrix
polar : bool
If True, get Stokes parameters for polar. Default: True
temp : bool
If True, get Stokes parameters for temperature. Default: False
EBTB : bool
If True, get Stokes parameters for EB and TB. Default: False
pixwining : bool
spec : 1D array of string
Spectra list
pixwin : bool
If True, multiplies the beam window function by the pixel
window function. Default: True
timing : bool
......@@ -60,55 +57,33 @@ def compute_ds_dcb(
Returns
----------
dcov : ndarray of floats
Pl : ndarray of floats
Normalize Legendre polynomials dS/dCl
Smatrix : 2D array of floats
S : 2D array of floats
Pixel signal covariance matrix S
Example
----------
>>> Pl, S = compute_ds_dcb(
... np.array([2,4,5,8]), 2, np.array([0,1,4,10,11]), np.arange(10),
... clth=np.arange(30).reshape(3,-1), Slmax=8,
... polar=True, temp=False, EBTB=False,
... pixwining=False, timing=False, MC=0, Sonly=False)
>>> print(round(np.sum(Pl),5), round(np.sum(S),5))
(1158.35868, 38277.82933)
>>> Pl, S = compute_ds_dcb(
... np.array([2,4,7,8]), 2, np.array([0,3,4,10,11]), np.arange(10),
... clth=np.arange(60).reshape(6,-1), Slmax=8,
... polar=True, temp=False, EBTB=True,
... pixwining=True, timing=False, MC=0, Sonly=False)
... spec=["TT", "EE", "BB", "TE", "EB", "TB"],
... pixwin=True, timing=False, MC=0, Sonly=False)
>>> print(round(np.sum(Pl),5), round(np.sum(S),5))
(570.84071, 15859.29629)
(1128.61613, 29501.83265)
>>> Pl, S = compute_ds_dcb(
... np.array([2,4,7,8]), 2, np.array([0,3,4,10,11]), np.arange(10),
... np.array([2,4,5,8]), 2, np.array([0,1,4,10,11]), np.arange(10),
... clth=np.arange(60).reshape(6,-1), Slmax=8,
... polar=False, temp=True, EBTB=False,
... pixwining=True, timing=False, MC=0, Sonly=False)
>>> print(round(np.sum(Pl),5), round(np.sum(S),5))
(109.73253, 679.79218)
>>> import pylab
>>> pylab.seed(0)
>>> Pl, S = compute_ds_dcb(
... np.array([2,3,4]), 4, np.array([0,1,3]), np.arange(13),
... clth=np.arange(4*13).reshape(4,-1), Slmax=11,
... polar=True, temp=False, EBTB=False,
... pixwining=True, timing=False, MC=100, Sonly=False)
... spec=["TT", "EE", "BB", "TE", "EB", "TB"],
... pixwin=False, timing=False, MC=0, Sonly=False)
>>> print(round(np.sum(Pl),5), round(np.sum(S),5))
(12.68135, 406990.12056)
(2590.55388, 65797.6979)
"""
if Slmax < ellbins[-1]-1:
print("WARNING : Slmax < lmax")
allStoke, der, ind = getstokes(polar=polar, temp=temp, EBTB=EBTB)
nder = len(der)
# ### define pixels
rpix = np.array(hp.pix2vec(nside, ipok))
allcosang = np.dot(np.transpose(rpix), rpix)
......@@ -117,34 +92,33 @@ def compute_ds_dcb(
if Sonly:
if MC:
Smatrix = S_bins_MC(
S = S_bins_MC(
ellbins, nside, ipok, allcosang, bl, clth, Slmax, MC,
polar=polar, temp=temp, EBTB=EBTB,
pixwining=pixwining, timing=timing)
polar=polar, temp=temp, corr=corr,
pixwin=pixwin, timing=timing)
else:
Smatrix = S_bins_fast(
S = compute_S(
ellbins, nside, ipok, allcosang, bl, clth, Slmax,
polar=polar, temp=temp, EBTB=EBTB,
pixwining=pixwining, timing=timing)
return Smatrix
polar=polar, temp=temp, corr=corr,
pixwin=pixwin, timing=timing)
return S
if MC:
dcov, Smatrix = covth_bins_MC(
Pl, S = covth_bins_MC(
ellbins, nside, ipok, allcosang, bl, clth, Slmax, MC,
polar=polar, temp=temp, EBTB=EBTB,
pixwining=pixwining, timing=timing)
polar=polar, temp=temp, corr=corr,
pixwin=pixwin, timing=timing)
else:
dcov, Smatrix = covth_bins_fast(
Pl, S = compute_PlS(
ellbins, nside, ipok, allcosang, bl, clth, Slmax,
polar=polar, temp=temp, EBTB=EBTB,
pixwining=pixwining, timing=timing)
spec=spec, pixwin=pixwin, timing=timing)
return (dcov, Smatrix)
return Pl, S
def covth_bins_fast(
def compute_PlS(
ellbins, nside, ipok, allcosang, bl, clth, Slmax,
polar=True, temp=True, EBTB=False, pixwining=False, timing=False):
spec, pixwin=True, timing=False):
"""
Computes Legendre polynomes Pl = dS/dCb and signal matrix S.
......@@ -163,13 +137,9 @@ def covth_bins_fast(
Fiducial power spectra
Slmax : int
Maximum lmax computed for the pixel covariance pixel matrix
polar : bool
If True, get Stokes parameters for polar. Default: True
temp : bool
If True, get Stokes parameters for temperature. Default: False
EBTB : bool
If True, get Stokes parameters for EB and TB. Default: False
pixwining : bool
spec : 1D array of string
Spectra list
pixwin : bool
If True, multiplies the beam window function by the pixel
window function. Default: True
timing : bool
......@@ -177,36 +147,28 @@ def covth_bins_fast(
Returns
----------
dcov : ndarray of floats
Pl : ndarray of floats
Normalize Legendre polynomials dS/dCl
Smatrix : 2D array of floats
S : 2D array of floats
Pixel signal covariance matrix S
Example
----------
>>> dcov, S = covth_bins_fast(np.array([2,3,4,7]),4,ipok=np.array([0,1,3]),
... allcosang=np.linspace(0,1,15).reshape(3,-1), bl=np.arange(13),
... clth=np.arange(4*13).reshape(4,-1), Slmax=11,
... polar=True, temp=False, EBTB=False, pixwining=True, timing=False)
>>> print(round(np.sum(dcov),5), round(np.sum(S),5))
(-22.04764, -2340.05934)
>>> dcov, S = covth_bins_fast(np.array([2,3,4,7]),4,ipok=np.array([0,1,3]),
>>> Pl, S = compute_PlS(np.array([2,3,4,7]),4,ipok=np.array([0,1,3]),
... allcosang=np.linspace(0,1,15).reshape(3,-1), bl=np.arange(13),
... clth=np.arange(6*13).reshape(6,-1), Slmax=11,
... polar=True, temp=False, EBTB=True, pixwining=True, timing=False)
>>> print(round(np.sum(dcov),5), round(np.sum(S),5))
(56.1748, 8597.38767)
... spec=["TT", "EE", "BB", "TE", "EB", "TB"], pixwin=True, timing=False)
>>> print(round(np.sum(Pl),5), round(np.sum(S),5))
(14.38023, -4125.5419)
>>> dcov, S = covth_bins_fast(np.array([2,3,4,7]),4,ipok=np.array([0,1,3]),
>>> Pl, S = compute_PlS(np.array([2,3,4,7]),4,ipok=np.array([0,1,3]),
... allcosang=np.linspace(0,1,15).reshape(3,-1), bl=np.arange(13),
... clth=np.arange(4*13).reshape(4,-1), Slmax=11,
... polar=False, temp=True, EBTB=False, pixwining=True, timing=False)
>>> print(round(np.sum(dcov),5), round(np.sum(S),5))
(-4.98779, -1270.80223)
... clth=np.arange(6*13).reshape(6,-1), Slmax=11,
... spec=["TT", "EE", "BB", "TE", "EB", "TB"], pixwin=False, timing=False)
>>> print(round(np.sum(Pl),5), round(np.sum(S),5))
(17.44139, 3086.13152)
"""
lmax = ellbins[-1]
......@@ -216,216 +178,207 @@ def covth_bins_fast(
maxell = np.array(ellbins[1:nbins+1]) - 1
ellval = (minell + maxell) * 0.5
allStoke, der, ind = getstokes(polar=polar, temp=temp, EBTB=EBTB)
nder = len(der)
npi = ipok.size
stokes, spec, istokes, ispecs = getstokes(spec)
nspec = len(spec)
nsto = len(stokes)
temp = "TT" in spec
polar = "EE" in spec or "BB" in spec
TE = 'TE' in spec
EB = 'EB' in spec
TB = 'TB' in spec
te = spec.index('TE') if TE else 0
tb = spec.index('TB') if TB else 0
eb = spec.index('EB') if EB else 0
ponbins = nbins*temp
ponpi = npi*temp
tenbins = te*nbins
tbnbins = tb*nbins
ebnbins = eb*nbins
rpix = np.array(hp.pix2vec(nside, ipok))
ll = np.arange(Slmax+2)
fpixwin = extrapolpixwin(nside, Slmax+2, pixwining)
norm = (2*ll[2:]+1)/(4.*np.pi)*(fpixwin[2:]**2)*(bl[2:Slmax+2]**2)
ll = np.arange(Slmax+1)
fpixwin = extrapolpixwin(nside, Slmax+1, pixwin)
norm = (2*ll[2:]+1)/(4.*np.pi)*(fpixwin[2:]**2)*(bl[2:Slmax+1]**2)
clthn = clth[:, 2: Slmax+1]