-
Notifications
You must be signed in to change notification settings - Fork 4
/
cmb_modules.py
332 lines (272 loc) · 12 KB
/
cmb_modules.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
import numpy as np
import matplotlib
import sys
import matplotlib.cm as cm
import matplotlib.mlab as mlab
import matplotlib.pyplot as plt
import astropy.io.fits as fits
def make_CMB_T_map(N,pix_size,ell,DlTT):
"makes a realization of a simulated CMB sky map given an input DlTT as a function of ell,"
"the pixel size (pix_size) required and the number N of pixels in the linear dimension."
#np.random.seed(100)
# convert Dl to Cl
ClTT = DlTT * 2 * np.pi / (ell*(ell+1.))
ClTT[0] = 0. # set the monopole and the dipole of the Cl spectrum to zero
ClTT[1] = 0.
# make a 2D real space coordinate system
onesvec = np.ones(N)
inds = (np.arange(N)+.5 - N/2.) /(N-1.) # create an array of size N between -0.5 and +0.5
# compute the outer product matrix: X[i, j] = onesvec[i] * inds[j] for i,j
# in range(N), which is just N rows copies of inds - for the x dimension
X = np.outer(onesvec,inds)
# compute the transpose for the y dimension
Y = np.transpose(X)
# radial component R
R = np.sqrt(X**2. + Y**2.)
# now make a 2D CMB power spectrum
pix_to_rad = (pix_size/60. * np.pi/180.) # going from pix_size in arcmins to degrees and then degrees to radians
ell_scale_factor = 2. * np.pi /pix_to_rad # now relating the angular size in radians to multipoles
ell2d = R * ell_scale_factor # making a fourier space analogue to the real space R vector
ClTT_expanded = np.zeros(int(ell2d.max())+1)
# making an expanded Cl spectrum (of zeros) that goes all the way to the size of the 2D ell vector
ClTT_expanded[0:(ClTT.size)] = ClTT # fill in the Cls until the max of the ClTT vector
# the 2D Cl spectrum is defined on the multiple vector set by the pixel scale
CLTT2d = ClTT_expanded[ell2d.astype(int)]
#plt.imshow(np.log(CLTT2d))
# now make a realization of the CMB with the given power spectrum in real space
random_array_for_T = np.random.normal(0,1,(N,N))
FT_random_array_for_T = np.fft.fft2(random_array_for_T) # take FFT since we are in Fourier space
FT_2d = np.sqrt(CLTT2d) * FT_random_array_for_T # we take the sqrt since the power spectrum is T^2
#plt.imshow(np.real(FT_2d))
## make a plot of the 2D cmb simulated map in Fourier space, note the x and y axis labels need to be fixed
#Plot_CMB_Map(np.real(np.conj(FT_2d)*FT_2d*ell2d * (ell2d+1)/2/np.pi),0,np.max(np.conj(FT_2d)*FT_2d*ell2d * (ell2d+1)/2/np.pi),ell2d.max(),ell2d.max()) ###
# move back from ell space to real space
CMB_T = np.fft.ifft2(np.fft.fftshift(FT_2d))
# move back to pixel space for the map
CMB_T = CMB_T/(pix_size /60.* np.pi/180.)
# we only want to plot the real component
CMB_T = np.real(CMB_T)
## return the map
return(CMB_T)
###############################
def Plot_CMB_Map(Map_to_Plot,c_min,c_max,X_width,Y_width):
from mpl_toolkits.axes_grid1 import make_axes_locatable
print("map mean:",np.mean(Map_to_Plot),"map rms:",np.std(Map_to_Plot))
plt.gcf().set_size_inches(10, 10)
im = plt.imshow(Map_to_Plot, interpolation='bilinear', origin='lower',cmap=cm.RdBu_r)
im.set_clim(c_min,c_max)
ax=plt.gca()
divider = make_axes_locatable(ax)
cax = divider.append_axes("right", size="5%", pad=0.05)
cbar = plt.colorbar(im, cax=cax)
#cbar = plt.colorbar()
im.set_extent([0,X_width,0,Y_width])
plt.ylabel('angle $[^\circ]$')
plt.xlabel('angle $[^\circ]$')
cbar.set_label('tempearture [uK]', rotation=270)
plt.show()
return(0)
###############################
def Poisson_source_component(N,pix_size,Number_of_Sources,Amplitude_of_Sources):
"makes a realization of a naive Poisson-distributed point source map"
PSMap = np.zeros([np.int(N),np.int(N)])
i = 0
print('Number of sources required: ', Number_of_Sources)
while (i < int(Number_of_Sources)):
pix_x = np.int(N*np.random.rand())
pix_y = np.int(N*np.random.rand())
PSMap[pix_x,pix_y] += np.random.poisson(Amplitude_of_Sources)
i = i + 1
return(PSMap)
###############################
def Exponential_source_component(N,pix_size,Number_of_Sources_EX,Amplitude_of_Sources_EX):
N=int(N)
"makes a realization of a naive exponentially-distributed point source map"
PSMap = np.zeros([N,N])
i = 0
while (i < Number_of_Sources_EX):
pix_x = int(N*np.random.rand() )
pix_y = int(N*np.random.rand())
PSMap[pix_x,pix_y] += np.random.exponential(Amplitude_of_Sources_EX)
i = i + 1
return(PSMap)
###############################
def SZ_source_component(N,pix_size,Number_of_SZ_Clusters,Mean_Amplitude_of_SZ_Clusters,SZ_beta,SZ_Theta_core,do_plots):
"makes a realization of a nieve SZ map"
N=int(N)
SZMap = np.zeros([N,N])
SZcat = np.zeros([3,Number_of_SZ_Clusters]) ## catalogue of SZ sources, X, Y, amplitude
# make a distribution of point sources with varying amplitude
i = 0
while (i < Number_of_SZ_Clusters):
pix_x = np.int(N*np.random.rand())
pix_y = np.int(N*np.random.rand() )
pix_amplitude = np.random.exponential(Mean_Amplitude_of_SZ_Clusters)*(-1.)
SZcat[0,i] = pix_x
SZcat[1,i] = pix_y
SZcat[2,i] = pix_amplitude
SZMap[pix_x,pix_y] += pix_amplitude
i = i + 1
if (do_plots):
hist,bin_edges = np.histogram(SZMap,bins = 50,range=[SZMap.min(),-10])
plt.semilogy(bin_edges[0:-1],hist)
plt.xlabel('source amplitude [$\mu$K]')
plt.ylabel('number or pixels')
plt.show()
# make a beta function
beta = beta_function(N,pix_size,SZ_beta,SZ_Theta_core)
# convovle the beta funciton with the point source amplitude to get the SZ map
FT_beta = np.fft.fft2(np.fft.fftshift(beta))
FT_SZMap = np.fft.fft2(np.fft.fftshift(SZMap))
SZMap = np.fft.fftshift(np.real(np.fft.ifft2(FT_beta*FT_SZMap)))
# return the SZ map
return(SZMap,SZcat)
###############################
def beta_function(N,pix_size,SZ_beta,SZ_Theta_core):
# make a beta function
N=int(N)
ones = np.ones(N)
inds = (np.arange(N)+.5 - N/2.) * pix_size
X = np.outer(ones,inds)
Y = np.transpose(X)
R = np.sqrt(X**2. + Y**2.)
beta = (1 + (R/SZ_Theta_core)**2.)**((1-3.*SZ_beta)/2.)
# return the beta function map
return(beta)
###############################
def convolve_map_with_gaussian_beam(N,pix_size,beam_size_fwhp,Map):
"convolves a map with a gaussian beam pattern. NOTE: pix_size and beam_size_fwhp need to be in the same units"
# make a 2d gaussian
gaussian = make_2d_gaussian_beam(N,pix_size,beam_size_fwhp)
# do the convolution
FT_gaussian = np.fft.fft2(np.fft.fftshift(gaussian))
FT_Map = np.fft.fft2(np.fft.fftshift(Map))
convolved_map = np.fft.fftshift(np.real(np.fft.ifft2(FT_gaussian*FT_Map)))
# return the convolved map
return(convolved_map)
###############################
def make_2d_gaussian_beam(N,pix_size,beam_size_fwhp):
# make a 2d coordinate system
N=int(N)
ones = np.ones(N)
inds = (np.arange(N)+.5 - N/2.) * pix_size
X = np.outer(ones,inds)
Y = np.transpose(X)
R = np.sqrt(X**2. + Y**2.)
# make a 2d gaussian
beam_sigma = beam_size_fwhp / np.sqrt(8.*np.log(2))
gaussian = np.exp(-.5 *(R/beam_sigma)**2.)
gaussian = gaussian / np.sum(gaussian)
# return the gaussian
return(gaussian)
###############################
def make_noise_map(N,pix_size,white_noise_level,atmospheric_noise_level,one_over_f_noise_level):
"makes a realization of instrument noise, atmosphere and 1/f noise level set at 1 degrees"
## make a white noise map
N=int(N)
white_noise = np.random.normal(0,1,(N,N)) * white_noise_level/pix_size
## make an atmosperhic noise map
atmospheric_noise = 0.
if (atmospheric_noise_level != 0):
ones = np.ones(N)
inds = (np.arange(N)+.5 - N/2.)
X = np.outer(ones,inds)
Y = np.transpose(X)
R = np.sqrt(X**2. + Y**2.) * pix_size /60. ## angles relative to 1 degrees
mag_k = 2 * np.pi/(R+.01) ## 0.01 is a regularizaiton factor
atmospheric_noise = np.fft.fft2(np.random.normal(0,1,(N,N)))
atmospheric_noise = np.fft.ifft2(atmospheric_noise * np.fft.fftshift(mag_k**(5/3.)))* atmospheric_noise_level/pix_size
## make a 1/f map, along a single direction to illustrate striping
oneoverf_noise = 0.
if (one_over_f_noise_level != 0):
ones = np.ones(N)
inds = (np.arange(N)+.5 - N/2.)
X = np.outer(ones,inds) * pix_size /60. ## angles relative to 1 degrees
kx = 2 * np.pi/(X+.01) ## 0.01 is a regularizaiton factor
oneoverf_noise = np.fft.fft2(np.random.normal(0,1,(N,N)))
oneoverf_noise = np.fft.ifft2(oneoverf_noise * np.fft.fftshift(kx))* one_over_f_noise_level/pix_size
## return the noise map
noise_map = np.real(white_noise + atmospheric_noise + oneoverf_noise)
return(noise_map)
###############################
def Filter_Map(Map,N,N_mask):
N=int(N)
## set up a x, y, and r coordinates for mask generation
ones = np.ones(N)
inds = (np.arange(N)+.5 - N/2.)
X = np.outer(ones,inds)
Y = np.transpose(X)
R = np.sqrt(X**2. + Y**2.) ## angles relative to 1 degrees
## make a mask
mask = np.ones([N,N])
mask[np.where(np.abs(X) < N_mask)] = 0
return apply_filter(Map,mask)
def apply_filter(Map,filter2d):
## apply the filter in fourier space
FMap = np.fft.fftshift(np.fft.ifft2(np.fft.fftshift(Map)))
FMap_filtered = FMap * filter2d
Map_filtered = np.real(np.fft.fftshift(np.fft.fft2(FMap_filtered)))
## return the output
return(Map_filtered)
def cosine_window(N):
"makes a cosine window for apodizing to avoid edges effects in the 2d FFT"
# make a 2d coordinate system
ones = np.ones(N)
inds = (np.arange(N)+.5 - N/2.)/N *np.pi ## eg runs from -pi/2 to pi/2
X = np.outer(ones,inds)
Y = np.transpose(X)
# make a window map
window_map = np.cos(X) * np.cos(Y)
# return the window map
return(window_map)
###############################
def average_N_spectra(spectra,N_spectra,N_ells):
avgSpectra = np.zeros(N_ells)
rmsSpectra = np.zeros(N_ells)
# calcuate the average spectrum
i = 0
while (i < N_spectra):
avgSpectra = avgSpectra + spectra[i,:]
i = i + 1
avgSpectra = avgSpectra/(1. * N_spectra)
#calculate the rms of the spectrum
i =0
while (i < N_spectra):
rmsSpectra = rmsSpectra + (spectra[i,:] - avgSpectra)**2
i = i + 1
rmsSpectra = np.sqrt(rmsSpectra/(1. * N_spectra))
return(avgSpectra,rmsSpectra)
def calculate_2d_spectrum(Map,delta_ell,ell_max,pix_size,N,Map2=None):
"calculates the power spectrum of a 2d map by FFTing, squaring, and azimuthally averaging"
import matplotlib.pyplot as plt
# make a 2d ell coordinate system
N=int(N)
ones = np.ones(N)
inds = (np.arange(N)+.5 - N/2.) /(N-1.)
kX = np.outer(ones,inds) / (pix_size/60. * np.pi/180.)
kY = np.transpose(kX)
K = np.sqrt(kX**2. + kY**2.)
ell_scale_factor = 2. * np.pi
ell2d = K * ell_scale_factor
# make an array to hold the power spectrum results
N_bins = int(ell_max/delta_ell)
ell_array = np.arange(N_bins)
CL_array = np.zeros(N_bins)
# get the 2d fourier transform of the map
FMap = np.fft.ifft2(np.fft.fftshift(Map))
if Map2 is None: FMap2 = FMap.copy()
else: FMap2 = np.fft.ifft2(np.fft.fftshift(Map2))
# print FMap
PSMap = np.fft.fftshift(np.real(np.conj(FMap) * FMap2))
# print PSMap
# fill out the spectra
i = 0
while (i < N_bins):
ell_array[i] = (i + 0.5) * delta_ell
inds_in_bin = ((ell2d >= (i* delta_ell)) * (ell2d < ((i+1)* delta_ell))).nonzero()
CL_array[i] = np.mean(PSMap[inds_in_bin])
i = i + 1
CL_array_new = CL_array[~np.isnan(CL_array)]
ell_array_new = ell_array[~np.isnan(CL_array)]
# return the power spectrum and ell bins
return(ell_array_new,CL_array_new*np.sqrt(pix_size /60.* np.pi/180.)*2.)