-
Notifications
You must be signed in to change notification settings - Fork 18
/
lensing.py
771 lines (593 loc) · 27.4 KB
/
lensing.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
from __future__ import print_function
import numpy as np
np.seterr(divide='ignore', invalid='ignore')
from orphics import maps
from pixell import enmap, utils
from scipy.special import factorial
try:
from pixell import lensing as enlensing
except:
print("WARNING: Couldn't load pixell lensing. Some features will be unavailable.")
from scipy.integrate import simps
from scipy.interpolate import splrep,splev
from scipy.fftpack import fftshift,ifftshift,fftfreq
from scipy.interpolate import interp1d
from pixell.fft import fft,ifft
from orphics.stats import bin2D
import time
from six.moves import cPickle as pickle
from orphics import stats
import os,sys
from pyfisher import get_lensing_nl as get_nl
def validate_geometry(shape,wcs,verbose=False):
area = enmap.area(shape,wcs)*(180./np.pi)**2.
if verbose: print("Geometry area : ", area, " sq.deg.")
if area>41252.:
print("WARNING: Geometry has area larger than full-sky.")
print(shape,wcs)
if area<(1./60./60.):
print("WARNING: Geometry has area less than 1 arcmin^2.")
print(shape,wcs)
res = np.rad2deg(maps.resolution(shape,wcs))
if verbose: print("Geometry pixel width : ", res*60., " arcmin.")
if res>30.0:
print("WARNING: Geometry has pixel larger than 30 degrees.")
print(shape,wcs)
if res<(1./60./60.):
print("WARNING: Geometry has pixel smaller than 1 arcsecond.")
print(shape,wcs)
def binned_nfw(mass,z,conc,cc,shape,wcs,bin_edges_arcmin,lmax=None,lmin=None,overdensity=200.,
critical=False,at_cluster_z=True,kmask=None,
sigma_mis=None,improved=True,hm=None,exclude_2h=False):
# mass in msolar/h
# cc Cosmology object
modrmap = enmap.modrmap(shape,wcs)
binner = bin2D(modrmap,bin_edges_arcmin*np.pi/180./60.)
if improved:
thetas = np.linspace(bin_edges_arcmin.min(),bin_edges_arcmin.max(),101) * utils.arcmin
Ms = [mass]
concs = [conc]
zsource = 1100
sig = sigma_mis*utils.arcmin if sigma_mis is not None else None
k1h = hm.kappa_1h_profiles(thetas,Ms,concs,zsource,sig_theta=sig,delta=overdensity,rho='critical' if critical else 'mean',rho_at_z=at_cluster_z)
k2h = hm.kappa_2h_profiles(thetas,Ms,zsource,delta=overdensity,rho='critical' if critical else 'mean',rho_at_z=at_cluster_z,lmin=2,lmax=10000) if not(exclude_2h) else np.asarray(k1h).T*0
k1h[~np.isfinite(k1h)] = 0
k1h = np.asarray(k1h[0])
k2h = k2h[:,0]
k = enmap.enmap(maps.interp(thetas,k1h+k2h)(modrmap),wcs)
else:
k = nfw_kappa(mass,modrmap,cc,zL=z,concentration=conc,overdensity=overdensity,critical=critical,atClusterZ=at_cluster_z)
if kmask is None: kmask = maps.mask_kspace(shape,wcs,lmin=lmin,lmax=lmax)
kf = maps.filter_map(k,kmask)
cents,k1d = binner.bin(kf)
return cents,k1d
def fit_nfw_profile(profile_data,profile_cov,masses,z,conc,cc,shape,wcs,bin_edges_arcmin,lmax,lmin=None,
overdensity=200.,critical=False,at_cluster_z=True,
mass_guess=2e14,sigma_guess=2e13,kmask=None,sigma_mis=None,improved=True):
"""
Returns
lnlikes - actual lnlike as function of masses
like_fit - gaussian fit as function of masses
fit_mass - fit mass
mass_err - fit mass err
"""
from orphics.stats import fit_gauss
cinv = np.linalg.inv(profile_cov)
lnlikes = []
fprofiles = []
if improved:
import hmvec
ks = np.geomspace(1e-2,10,200)
ms = np.geomspace(1e8,4e15,102)
zs = [z]
hm = hmvec.HaloModel(zs,ks,ms=ms,params={},mass_function="tinker",
halofit=None,mdef='mean',nfw_numeric=False,skip_nfw=True)
else:
hm = None
for mass in masses:
cents,profile_theory = binned_nfw(mass,z,conc,cc,shape,wcs,bin_edges_arcmin,lmax,lmin,
overdensity,critical,at_cluster_z,kmask=kmask,
sigma_mis=sigma_mis,improved=improved,hm=hm)
diff = profile_data - profile_theory
fprofiles.append(profile_theory)
lnlike = -0.5 * np.dot(np.dot(diff,cinv),diff)
lnlikes.append(lnlike)
fit_mass,mass_err,_,_ = fit_gauss(masses,np.exp(lnlikes),mu_guess=mass_guess,sigma_guess=sigma_guess)
gaussian = lambda t,mu,sigma: np.exp(-(t-mu)**2./2./sigma**2.)/np.sqrt(2.*np.pi*sigma**2.)
like_fit = gaussian(masses,fit_mass,mass_err)
cents,fit_profile = binned_nfw(fit_mass,z,conc,cc,shape,wcs,bin_edges_arcmin,lmax,lmin,
overdensity,critical,at_cluster_z,kmask=kmask,
sigma_mis=sigma_mis,improved=improved,hm=hm)
return np.array(lnlikes),np.array(like_fit),fit_mass,mass_err,fprofiles,fit_profile
def mass_estimate(kappa_recon,kappa_noise_2d,mass_guess,concentration,z):
"""Given a cutout kappa map centered on a cluster and a redshift,
returns a mass estimate and variance of the estimate by applying a matched filter.
Imagine that reliable richness estimates and redshifts exist for each cluster.
We split the sample into n richness bins.
We go through each bin. This bin has a mean richness and a mean redshift. We convert this to a fiducial mean mass and mean concentration.
We choose a template with this mean mass and mean concentration.
We do a reconstruction on each cluster for each array. We now have a 2D kappa measurement. We apply MF to this with the template.
We get a relative amplitude for each cluster, which we convert to a mass for each cluster.
We want a mean mass versus mean richness relationship.
Step 1
Loop through each cluster.
Cut out patches from each array split-coadd and split.
Estimate noise spectrum from splits for each array.
Use noise spectra to get optimal coadd of all arrays and splits of coadds.
Use coadd for reconstruction, with coadd noise from splits in weights.
For gradient, use Planck.
Save reconstructions and 2d kappa noise to disk.
Repeat above for 100x random locations and save only mean to disk.
Postprocess by loading each reconstruction, subtract meanfield, crop out region where taper**2 < 1 with some threshold.
Verify above by performing on simulations to check that cluster profile is recovered.
Step 2
For each richness bin, use fiducial mass and concentration and mean redshift to choose template.
For each cluster in each richness bin, apply MF and get masses. Find mean mass in bin. Iterate above on mass until converged.
This provides a mean mass, concentration for each richness bin.
"""
shape,wcs = kappa_recon.shape,kappa_recon.wcs
mf = maps.MatchedFilter(shape,wcs,template,kappa_noise_2d)
mf.apply(kappa_recon,kmask=kmask)
def flat_taylens(phi,imap,taylor_order = 5):
"""
Lens a map imap with lensing potential phi
using the Taylens algorithm up to taylor_order.
The original routine is from Thibaut Louis.
It has been modified here to work with pixell.
"""
f = lambda x: enmap.fft(x,normalize='phys')
invf = lambda x: enmap.ifft(x,normalize='phys')
def binomial(n,k):
"Compute n factorial by a direct multiplicative method"
if k > n-k: k = n-k # Use symmetry of Pascal's triangle
accum = 1
for i in range(1,k+1):
accum *= (n - (k - i))
accum /= i
return accum
kmap = f(phi)
Ny,Nx = phi.shape
ly,lx = enmap.laxes(phi.shape,phi.wcs)
ly_array,lx_array = phi.lmap()
alphaX=np.real(invf(1j*lx_array*kmap))
alphaY=np.real(invf(1j*ly_array*kmap))
iy,ix = np.mgrid[0:Ny,0:Nx]
py,px = enmap.pixshape(phi.shape,phi.wcs)
alphaX0 = np.array(np.round(alphaX/ px),dtype='int64')
alphaY0 = np.array(np.round(alphaY/ py),dtype='int64')
delta_alphaX=alphaX-alphaX0*px
delta_alphaY=alphaY-alphaY0*py
lensed_T_Map = imap.copy()
cont = imap.copy()
lensed_T_Map = imap[(iy+alphaY0)%Ny, (ix+alphaX0)%Nx]
kmap=f(imap)
for n in range(1,taylor_order):
cont[:]=0
for k in range(n+1):
fac=1j**n*binomial(n,k)*lx_array**(n-k)*ly_array**k/(factorial(n))
T_add=np.real(invf(fac*kmap))[(iy+alphaY0)%Ny, (ix+alphaX0)%Nx]*delta_alphaX**(n-k)*delta_alphaY**k
lensed_T_Map[:] += T_add
cont[:] += T_add
return lensed_T_Map
def alpha_from_kappa(kappa=None,posmap=None,phi=None,grad=True):
if phi is None:
phi,_ = kappa_to_phi(kappa,kappa.modlmap(),return_fphi=True)
shape,wcs = phi.shape,phi.wcs
else:
shape,wcs = phi.shape,phi.wcs
grad_phi = enmap.grad(phi)
if grad: return grad_phi
if posmap is None: posmap = enmap.posmap(shape,wcs)
pos = posmap + grad_phi
alpha_pix = enmap.sky2pix(shape,wcs,pos, safe=False)
return enmap.enmap(alpha_pix,wcs)
class FlatLensingSims(object):
def __init__(self,shape,wcs,theory,beam_arcmin,noise_uk_arcmin,noise_e_uk_arcmin=None,noise_b_uk_arcmin=None,pol=False,fixed_lens_kappa=None):
# assumes theory in uK^2
from orphics import cosmology
if len(shape)<3 and pol: shape = (3,)+shape
if noise_e_uk_arcmin is None: noise_e_uk_arcmin = np.sqrt(2.)*noise_uk_arcmin
if noise_b_uk_arcmin is None: noise_b_uk_arcmin = noise_e_uk_arcmin
self.modlmap = enmap.modlmap(shape,wcs)
Ny,Nx = shape[-2:]
lmax = self.modlmap.max()
ells = np.arange(0,lmax,1)
ps_cmb = cosmology.power_from_theory(ells,theory,lensed=False,pol=pol)
self.mgen = maps.MapGen(shape,wcs,ps_cmb)
if fixed_lens_kappa is not None:
self._fixed = True
self.kappa = fixed_lens_kappa
self.alpha = alpha_from_kappa(self.kappa)
else:
self._fixed = False
ps_kk = theory.gCl('kk',self.modlmap).reshape((1,1,Ny,Nx))
self.kgen = maps.MapGen(shape[-2:],wcs,ps_kk)
self.posmap = enmap.posmap(shape[-2:],wcs)
self.ps_kk = ps_kk
self.kbeam = maps.gauss_beam(self.modlmap,beam_arcmin)
ncomp = 3 if pol else 1
ps_noise = np.zeros((ncomp,ncomp,Ny,Nx))
ps_noise[0,0] = (noise_uk_arcmin*np.pi/180./60.)**2.
if pol:
ps_noise[1,1] = (noise_e_uk_arcmin*np.pi/180./60.)**2.
ps_noise[2,2] = (noise_b_uk_arcmin*np.pi/180./60.)**2.
self.ngen = maps.MapGen(shape,wcs,ps_noise)
self.ps_noise = ps_noise
def update_kappa(self,kappa):
self.kappa = kappa
self.alpha = alpha_from_kappa(self.kappa)
def get_unlensed(self,seed=None):
return self.mgen.get_map(seed=seed)
def get_kappa(self,seed=None):
return self.kgen.get_map(seed=seed)
def get_sim(self,seed_cmb=None,seed_kappa=None,seed_noise=None,lens_order=5,return_intermediate=False,skip_lensing=False,cfrac=None):
unlensed = self.get_unlensed(seed_cmb)
if skip_lensing:
lensed = unlensed
kappa = enmap.samewcs(lensed.copy()[0]*0,lensed)
else:
if not(self._fixed):
kappa = self.get_kappa(seed_kappa)
self.kappa = kappa
self.alpha = alpha_from_kappa(kappa,posmap=self.posmap)
else:
kappa = None
assert seed_kappa is None
lensed = enlensing.displace_map(unlensed, self.alpha, order=lens_order)
beamed = maps.filter_map(lensed,self.kbeam)
noise_map = self.ngen.get_map(seed=seed_noise)
observed = beamed + noise_map
if return_intermediate:
return [ maps.get_central(x,cfrac) for x in [unlensed,kappa,lensed,beamed,noise_map,observed] ]
else:
return maps.get_central(observed,cfrac)
def lens_cov_pol(shape,wcs,iucov,alpha_pix,lens_order=5,kbeam=None,npixout=None,comm=None):
"""Given the pix-pix covariance matrix for the unlensed CMB,
returns the lensed covmat for a given pixel displacement model.
ucov -- (ncomp,ncomp,Npix,Npix) array where Npix = Ny*Nx
alpha_pix -- (2,Ny,Nx) array of lensing displacements in pixel units
kbeam -- (Ny,Nx) array of 2d beam wavenumbers
"""
from pixell import lensing as enlensing
assert iucov.ndim==4
ncomp = iucov.shape[0]
assert ncomp==iucov.shape[1]
assert 1 <= ncomp <= 3
if len(shape)==2: shape = (1,)+shape
n = shape[-2]
assert n==shape[-1]
ucov = iucov.copy()
ucov = np.transpose(ucov,(0,2,1,3))
ucov = ucov.reshape((ncomp*n**2,ncomp*n**2))
npix = ncomp*n**2
if comm is None:
from orphics import mpi
comm = mpi.MPI.COMM_WORLD
def efunc(vec):
unlensed = enmap.enmap(vec.reshape(shape),wcs)
lensed = enlensing.displace_map(unlensed, alpha_pix, order=lens_order)
if kbeam is not None: lensed = maps.filter_map(lensed,kbeam) # TODO: replace with convolution
# because for ~(60x60) arrays, it is probably much faster. >1 threads means worse performance
# with FFTs for these array sizes.
return np.asarray(lensed).reshape(-1)
Scov = np.zeros(ucov.shape,dtype=ucov.dtype)
for i in range(comm.rank, npix, comm.size):
Scov[i,:] = efunc(ucov[i,:])
Scov2 = utils.allreduce(Scov, comm)
Scov = np.zeros(ucov.shape,dtype=ucov.dtype)
for i in range(comm.rank, npix, comm.size):
Scov[:,i] = efunc(Scov2[:,i])
Scov = utils.allreduce(Scov, comm)
Scov = Scov.reshape((ncomp,n*n,ncomp,n*n))
if (npixout is not None) and (npixout!=n):
Scov = Scov.reshape((ncomp,n,n,ncomp,n,n))
s = n//2-npixout//2
e = s + npixout
Scov = Scov[:,s:e,s:e,:,s:e,s:e].reshape((ncomp,npixout**2,ncomp,npixout**2))
Scov = np.transpose(Scov,(0,2,1,3))
return Scov
def lens_cov(shape,wcs,ucov,alpha_pix,lens_order=5,kbeam=None,bshape=None):
"""Given the pix-pix covariance matrix for the unlensed CMB,
returns the lensed covmat for a given pixel displacement model.
ucov -- (Npix,Npix) array where Npix = Ny*Nx
alpha_pix -- (2,Ny,Nx) array of lensing displacements in pixel units
kbeam -- (Ny,Nx) array of 2d beam wavenumbers
"""
from pixell import lensing as enlensing
Scov = ucov.copy()
for i in range(ucov.shape[0]):
unlensed = enmap.enmap(Scov[i,:].copy().reshape(shape),wcs)
lensed = enlensing.displace_map(unlensed, alpha_pix, order=lens_order)
if kbeam is not None: lensed = maps.filter_map(lensed,kbeam)
Scov[i,:] = lensed.ravel()
for j in range(ucov.shape[1]):
unlensed = enmap.enmap(Scov[:,j].copy().reshape(shape),wcs)
lensed = enlensing.displace_map(unlensed, alpha_pix, order=lens_order)
if kbeam is not None: lensed = maps.filter_map(lensed,kbeam)
Scov[:,j] = lensed.ravel()
if (bshape is not None) and (bshape!=shape):
ny,nx = shape
Scov = Scov.reshape((ny,nx,ny,nx))
bny,bnx = bshape
sy = ny//2-bny//2
ey = sy + bny
sx = nx//2-bnx//2
ex = sx + bnx
Scov = Scov[sy:ey,sx:ex,sy:ey,sx:ex].reshape((np.prod(bshape),np.prod(bshape)))
return Scov
def beam_cov(ucov,kbeam):
"""Given the pix-pix covariance matrix for the lensed CMB,
returns the beamed covmat. The beam can be a ratio of beams to
readjust the beam in a given matrix.
ucov -- (Npix,Npix) array where Npix = Ny*Nx
kbeam -- (Ny,Nx) array of 2d beam wavenumbers
"""
Scov = ucov.copy()
wcs = ucov.wcs
shape = kbeam.shape[-2:]
for i in range(Scov.shape[0]):
lensed = enmap.enmap(Scov[i,:].copy().reshape(shape) ,wcs)
lensed = maps.filter_map(lensed,kbeam)
Scov[i,:] = lensed.ravel()
for j in range(Scov.shape[1]):
lensed = enmap.enmap(Scov[:,j].copy().reshape(shape),wcs)
lensed = maps.filter_map(lensed,kbeam)
Scov[:,j] = lensed.ravel()
return Scov
def kappa_to_phi(kappa,modlmap,return_fphi=False):
fphi = enmap.samewcs(kappa_to_fphi(kappa,modlmap),kappa)
phi = enmap.samewcs(enmap.ifft(fphi,normalize='phys').real, kappa)
if return_fphi:
return phi, fphi
else:
return phi
def kappa_to_fphi(kappa,modlmap):
return fkappa_to_fphi(enmap.fft(kappa,normalize='phys'),modlmap)
def fkappa_to_fphi(fkappa,modlmap):
kmap = np.nan_to_num(2.*fkappa/modlmap/(modlmap+1.))
kmap[modlmap<2.] = 0.
return kmap
def fillLowEll(ells,cls,ellmin):
# Fill low ells with the same value
low_index = np.where(ells>ellmin)[0][0]
lowest_ell = ells[low_index]
lowest_val = cls[low_index]
fill_ells = np.arange(2,lowest_ell,1)
new_ells = np.append(fill_ells,ells[low_index:])
fill_cls = np.array([lowest_val]*len(fill_ells))
new_cls = np.append(fill_cls,cls[low_index:])
return new_ells,new_cls
def sanitizePower(Nlbinned):
Nlbinned[Nlbinned<0.] = np.nan
# fill nans with interp
ok = ~np.isnan(Nlbinned)
xp = ok.ravel().nonzero()[0]
fp = Nlbinned[~np.isnan(Nlbinned)]
x = np.isnan(Nlbinned).ravel().nonzero()[0]
Nlbinned[np.isnan(Nlbinned)] = np.interp(x, xp, fp)
return Nlbinned
## HALOS
# g(x) = g(theta/thetaS) HuDeDeoVale 2007
gnfw = lambda x: np.piecewise(x, [x>1., x<1., x==1.], \
[lambda y: (1./(y*y - 1.)) * \
( 1. - ( (2./np.sqrt(y*y - 1.)) * np.arctan(np.sqrt((y-1.)/(y+1.))) ) ), \
lambda y: (1./(y*y - 1.)) * \
( 1. - ( (2./np.sqrt(-(y*y - 1.))) * np.arctanh(np.sqrt(-((y-1.)/(y+1.)))) ) ), \
lambda y: (1./3.)])
f_c = lambda c: np.log(1.+c) - (c/(1.+c))
def nfw_kappa(massOverh,modrmap_radians,cc,zL=0.7,concentration=3.2,overdensity=180.,critical=False,atClusterZ=False):
sgn = 1. if massOverh>0. else -1.
comS = cc.results.comoving_radial_distance(cc.cmbZ)*cc.h
comL = cc.results.comoving_radial_distance(zL)*cc.h
winAtLens = (comS-comL)/comS
kappa,r500 = NFWkappa(cc,np.abs(massOverh),concentration,zL,modrmap_radians* 180.*60./np.pi,winAtLens,
overdensity=overdensity,critical=critical,atClusterZ=atClusterZ)
return sgn*kappa
def NFWkappa(cc,massOverh,concentration,zL,thetaArc,winAtLens,overdensity=500.,critical=True,atClusterZ=True):
comL = (cc.results.comoving_radial_distance(zL) )*cc.h
c = concentration
M = massOverh
zdensity = 0.
if atClusterZ: zdensity = zL
if critical:
r500 = cc.rdel_c(M,zdensity,overdensity).flatten()[0] # R500 in Mpc/h
else:
r500 = cc.rdel_m(M,zdensity,overdensity) # R500 in Mpc/h
conv=np.pi/(180.*60.)
theta = thetaArc*conv # theta in radians
rS = r500/c
thetaS = rS/ comL
const12 = 9.571e-20 # 2G/c^2 in Mpc / solar mass
fc = np.log(1.+c) - (c/(1.+c))
#const3 = comL * comLS * (1.+zL) / comS # Mpc
const3 = comL * (1.+zL) *winAtLens # Mpc
const4 = M / (rS*rS) #solar mass / MPc^2
const5 = 1./fc
kappaU = gnfw(theta/thetaS)+theta*0. # added for compatibility with enmap
consts = const12 * const3 * const4 * const5
kappa = consts * kappaU
if thetaArc.shape[0]%2==1 and thetaArc.shape[1]%2==1:
Ny,Nx = thetaArc.shape
cx = int(Nx/2.)
cy = int(Ny/2.)
kappa[cy,cx] = kappa[cy-1,cx]
assert np.all(np.isfinite(kappa))
return kappa, r500
def NFWMatchedFilterSN(clusterCosmology,log10Moverh,c,z,ells,Nls,kellmax,overdensity=500.,critical=True,atClusterZ=True,arcStamp=100.,pxStamp=0.05,saveId=None,verbose=False,rayleighSigmaArcmin=None,returnKappa=False,winAtLens=None):
if rayleighSigmaArcmin is not None: assert rayleighSigmaArcmin>=pxStamp
M = 10.**log10Moverh
shape,wcs = maps.rect_geometry(width_deg=arcStamp/60.,px_res_arcmin=pxStamp)
kellmin = 2.*np.pi/arcStamp*np.pi/60./180.
modLMap = enmap.modlmap(shape,wcs)
xMap,yMap,modRMap,xx,yy = maps.get_real_attributes(shape,wcs)
cc = clusterCosmology
cmb = False
if winAtLens is None:
cmb = True
comS = cc.results.comoving_radial_distance(cc.cmbZ)*cc.h
comL = cc.results.comoving_radial_distance(z)*cc.h
winAtLens = (comS-comL)/comS
kappaReal, r500 = NFWkappa(cc,M,c,z,modRMap*180.*60./np.pi,winAtLens,overdensity=overdensity,critical=critical,atClusterZ=atClusterZ)
dAz = cc.results.angular_diameter_distance(z) * cc.h
# print ("daz " , dAz , " mpc")
# print ("r500 " , r500 , " mpc")
th500 = r500/dAz
#fiveth500 = 10.*np.pi/180./60. #5.*th500
fiveth500 = 5.*th500
# print ("5theta500 " , fiveth500*180.*60./np.pi , " arcminutes")
# print ("maximum theta " , modRMap.max()*180.*60./np.pi, " arcminutes")
kInt = kappaReal.copy()
kInt[modRMap>fiveth500] = 0.
# print "mean kappa inside theta500 " , kInt[modRMap<fiveth500].mean()
# print "area of th500 disc " , np.pi*fiveth500**2.*(180.*60./np.pi)**2.
# print "estimated integral " , kInt[modRMap<fiveth500].mean()*np.pi*fiveth500**2.
k500 = simps(simps(kInt, yy), xx)
if verbose: print(("integral of kappa inside disc ",k500))
kappaReal[modRMap>fiveth500] = 0. #### !!!!!!!!! Might not be necessary!
# if cmb: print z,fiveth500*180.*60./np.pi
Ukappa = kappaReal/k500
# pl = Plotter()
# pl.plot2d(Ukappa)
# pl.done("output/kappa.png")
ellmax = kellmax
ellmin = kellmin
Uft = fft(Ukappa,axes=[-2,-1])
if rayleighSigmaArcmin is not None:
Prayleigh = rayleigh(modRMap*180.*60./np.pi,rayleighSigmaArcmin)
outDir = "/gpfs01/astro/www/msyriac/plots/"
# io.quickPlot2d(Prayleigh,outDir+"rayleigh.png")
rayK = fft(ifftshift(Prayleigh),axes=[-2,-1])
rayK /= rayK[modLMap<1.e-3]
Uft = Uft.copy()*rayK
Upower = np.real(Uft*Uft.conjugate())
# pl = Plotter()
# pl.plot2d(fftshift(Upower))
# pl.done("output/upower.png")
Nls[Nls<0.]=0.
s = splrep(ells,Nls,k=3)
Nl2d = splev(modLMap,s)
Nl2d[modLMap<ellmin]=np.inf
Nl2d[modLMap>ellmax] = np.inf
Ny,Nx = shape
pixScaleY,pixScaleX = enmap.pixshape(shape,wcs)
area = Nx*Ny*pixScaleX*pixScaleY
Upower = Upower *area / (Nx*Ny)**2
filter = np.nan_to_num(Upower/Nl2d)
#filter = np.nan_to_num(1./Nl2d)
filter[modLMap>ellmax] = 0.
filter[modLMap<ellmin] = 0.
# pl = Plotter()
# pl.plot2d(fftshift(filter))
# pl.done("output/filter.png")
# if (cmb): print Upower.sum()
# if not(cmb) and z>2.5:
# bin_edges = np.arange(500,ellmax,100)
# binner = bin2D(modLMap, bin_edges)
# centers, nl2dells = binner.bin(Nl2d)
# centers, upowerells = binner.bin(np.nan_to_num(Upower))
# centers, filterells = binner.bin(filter)
# from orphics.tools.io import Plotter
# pl = Plotter(scaleY='log')
# pl.add(centers,upowerells,label="upower")
# pl.add(centers,nl2dells,label="noise")
# pl.add(centers,filterells,label="filter")
# pl.add(ells,Nls,ls="--")
# pl.legendOn(loc='upper right')
# #pl._ax.set_ylim(0,1e-8)
# pl.done("output/filterells.png")
# sys.exit()
varinv = filter.sum()
std = np.sqrt(1./varinv)
sn = k500/std
if verbose: print(sn)
if saveId is not None:
np.savetxt("data/"+saveId+"_m"+str(log10Moverh)+"_z"+str(z)+".txt",np.array([log10Moverh,z,1./sn]))
if returnKappa:
return sn,ifft(Uft,axes=[-2,-1],normalize=True).real*k500
return sn, k500, std
def rayleigh(theta,sigma):
sigmasq = sigma*sigma
#return np.exp(-0.5*theta*theta/sigmasq)
return theta/sigmasq*np.exp(-0.5*theta*theta/sigmasq)
# NFW dimensionless form
fnfw = lambda x: 1./(x*((1.+x)**2.))
Gval = 4.517e-48 # Newton G in Mpc,seconds,Msun units
cval = 9.716e-15 # speed of light in Mpc,second units
# NFW density (M/L^3) as a function of distance from center of cluster
def rho_nfw(M,c,R):
return lambda r: 1./(4.*np.pi)*((c/R)**3.)*M/f_c(c)*fnfw(c*r/R)
# NFW projected along line of sight (M/L^2) as a function of angle on the sky in radians
def proj_rho_nfw(theta,comL,M,c,R):
thetaS = R/c/comL
return 1./(4.*np.pi)*((c/R)**2.)*M/f_c(c)*(2.*gnfw(theta/thetaS))
# Generic profile projected along line of sight (M/L^2) as a function of angle on the sky in radians
# rhoFunc is density (M/L^3) as a function of distance from center of cluster
def projected_rho(thetas,comL,rhoFunc,pmaxN=2000,numps=500000):
# default integration times are good to 0.01% for z=0.1 to 3
# increase numps for lower z/theta and pmaxN for higher z/theta
# g(x) = \int dl rho(sqrt(l**2+x**2)) = g(theta/thetaS)
pzrange = np.linspace(-pmaxN,pmaxN,numps)
g = np.array([np.trapz(rhoFunc(np.sqrt(pzrange**2.+(theta*comL)**2.)),pzrange) for theta in thetas])
return g
def kappa_nfw_generic(theta,z,comLMpcOverh,M,c,R,windowAtLens):
return 4.*np.pi*Gval*(1+z)*comLMpcOverh*windowAtLens*proj_rho_nfw(theta,comLMpcOverh,M,c,R)/cval**2.
def kappa_generic(theta,z,comLMpcOverh,rhoFunc,windowAtLens,pmaxN=2000,numps=500000):
# default integration times are good to 0.01% for z=0.1 to 3
# increase numps for lower z/theta and pmaxN for higher z/theta
return 4.*np.pi*Gval*(1+z)*comLMpcOverh*windowAtLens*projected_rho(theta,comLMpcOverh,rhoFunc,pmaxN,numps)/cval**2.
def kappa_from_rhofunc(M,c,R,theta,cc,z,rhoFunc=None):
if rhoFunc is None: rhoFunc = rho_nfw(M,c,R)
sgn = 1. if M>0. else -1.
comS = cc.results.comoving_radial_distance(cc.cmbZ)*cc.h
comL = cc.results.comoving_radial_distance(z)*cc.h
winAtLens = (comS-comL)/comS
kappa = kappa_generic(theta,z,comL,rhoFunc,winAtLens)
return sgn*kappa
def kappa_nfw(M,c,R,theta,cc,z):
sgn = 1. if M>0. else -1.
comS = cc.results.comoving_radial_distance(cc.cmbZ)*cc.h
comL = cc.results.comoving_radial_distance(z)*cc.h
winAtLens = (comS-comL)/comS
kappa = kappa_nfw_generic(theta,z,comL,np.abs(M),c,R,winAtLens)
return sgn*kappa
class SplitLensing(object):
def __init__(self,shape,wcs,qest,XY="TT"):
# PS calculator
self.fc = maps.FourierCalc(shape,wcs)
self.qest = qest
self.est = XY
def qpower(self,k1,k2):
# PS func
return self.fc.f2power(k1,k2)
def qfrag(self,a,b):
# kappa func (accepts fts, returns ft)
if self.est=='TT':
k1 = self.qest.kappa_from_map(self.est,T2DData=a.copy(),T2DDataY=b.copy(),alreadyFTed=True,returnFt=True)
elif self.est=='EE': # wrong!
k1 = self.qest.kappa_from_map(self.est,T2DData=a.copy(),E2DData=a.copy(),B2DData=a.copy(),
T2DDataY=b.copy(),E2DDataY=b.copy(),B2DDataY=b.copy(),alreadyFTed=True,returnFt=True)
return k1
def cross_estimator(self,ksplits):
# 4pt from splits
splits = ksplits
splits = np.asanyarray(ksplits)
insplits = splits.shape[0]
nsplits = float(insplits)
s = np.mean(splits,axis=0)
k = self.qfrag(s,s)
kiisum = 0.
psum = 0.
psum2 = 0.
for i in range(insplits):
mi = splits[i]
ki = (self.qfrag(mi,s)+self.qfrag(s,mi))/2.
kii = self.qfrag(mi,mi)
kiisum += kii
kic = ki - (1./nsplits)*kii
psum += self.qpower(kic,kic)
for j in range(i+1,int(insplits)):
mj = splits[j]
kij = (self.qfrag(mi,mj)+self.qfrag(mj,mi))/2.
psum2 += self.qpower(kij,kij)
kc = k - (1./nsplits**2.)*kiisum
return (nsplits**4.*self.qpower(kc,kc)-4.*nsplits**2.*psum+4.*psum2)/nsplits/(nsplits-1.)/(nsplits-2.)/(nsplits-3.)