-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathEH98_funcs.py
More file actions
171 lines (148 loc) · 6.46 KB
/
EH98_funcs.py
File metadata and controls
171 lines (148 loc) · 6.46 KB
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
import numpy as np
import scipy.constants as conts
def EH98(kvector, redshift, scaling_factor, cosmo=None):
#This is the link to the paper: https://arxiv.org/pdf/astro-ph/9710252.pdf
#The input kvector should be in unit of h/Mpc after rescaling by rd.
cdict = cosmo.get_current_derived_parameters(['z_d'])
h = cosmo.h()
H_at_z = cosmo.Hubble(redshift) * conts.c /1000. /(100.*h)
Omm = cosmo.Omega_m()
Omb = cosmo.Omega_b()
#Cannot find the following function.
# Omc = cosmo.omegach2()/h**2.
Omc = cosmo.Omega0_cdm()
Omm_at_z = Omm*(1.+redshift)**3./H_at_z**2.
OmLambda_at_z = 1.-Omm_at_z
ns = cosmo.n_s()
rs = cosmo.rs_drag()*h/scaling_factor
Omnu = Omm-Omb-Omc
fnu = Omnu/Omm
fb = Omb/Omm
fnub = (Omb+Omnu)/Omm
fc = Omc/Omm
fcb = (Omc+Omb)/Omm
pc = 1./4.*(5-np.sqrt(1+24*fc))
pcb = 1./4.*(5-np.sqrt(1+24*fcb))
Neff = cosmo.Neff()
Omg = cosmo.Omega_g()
Omr = Omg * (1. + Neff * (7./8.)*(4./11.)**(4./3.))
aeq = Omr/(Omb+Omc)/(1-fnu)
zeq = 1./aeq -1.
Heq = cosmo.Hubble(zeq)/h
keq = aeq*Heq*scaling_factor
zd = cdict['z_d']
yd = (1.+zeq)/(1.+zd)
growth = cosmo.scale_independent_growth_factor(redshift)
if (fnu==0):
Nnu = 0.
else:
Nnu = 1.
#alpha_gamma = 1 - 0.328*np.log(431*Omm*h**2)*Omb/Omm + 0.38*np.log(22.3*Omm*h**2)*(Omb/Omm)**2
#There seems to be a mistake in this equation.
# alpha_nu = fc/fcb * (5 - 2*(pc+pcb))/(5-4*pcb) * (1-0.553*fnub + 0.126*fnub**3.) / (1 - 0.193*np.sqrt(fnu)*Nnu**0.2 + 0.169*fnu) \
# *(1.+yd)**(pcb-pc) * (1+(pc-pcb)/2*(1+1./(3-4*pc)/(7-4*pcb))*(1.+yd)**(-1.))
alpha_nu = fc/fcb * (5 - 2*(pc+pcb))/(5-4*pcb) * (1-0.553*fnub + 0.126*fnub**3.) / (1 - 0.193*np.sqrt(fnu*Nnu) + 0.169*fnu*(Nnu)**0.2) \
*(1.+yd)**(pcb-pc) * (1+(pc-pcb)/2*(1+1./(3-4*pc)/(7-4*pcb))*(1.+yd)**(-1.))
#eff_shape = (alpha_gamma + (1.-alpha_gamma)/(1+(0.43*kvector*rs)**4.))
eff_shape = (np.sqrt(alpha_nu) + (1.-np.sqrt(alpha_nu))/(1+(0.43*kvector*rs)**4.))
#add in q_test.
q_test = kvector/keq*7.46e-2
q0 = kvector/(keq/7.46e-2)/eff_shape
betac = (1.-0.949*fnub)**(-1.)
L0 = np.log(np.exp(1.)+1.84*np.sqrt(alpha_nu)*betac*q0)
# C0 = 14.4 + 325./(1+60.5*q0**1.08)
C0 = 14.4 + 325./(1+60.5*q0**1.11)
T0 = L0/(L0+C0*q0**2.)
if (fnu==0):
yfs=0.
qnu=3.92*q_test
else:
yfs = 17.2*fnu*(1.+0.488*fnu**(-7./6.))*(Nnu*q_test/fnu)**2.
qnu = 3.92*q_test*np.sqrt(Nnu/fnu)
D1 = (1.+zeq)/(1.+redshift)*5*Omm_at_z/2./(Omm_at_z**(4./7.) - OmLambda_at_z + (1.+Omm_at_z/2.)*(1.+OmLambda_at_z/70.))
Dcbnu = (fcb**(0.7/pcb)+(D1/(1.+yfs))**0.7)**(pcb/0.7) * D1**(1.-pcb)
Bk = 1. + 1.24*fnu**(0.64)*Nnu**(0.3+0.6*fnu)/(qnu**(-1.6)+qnu**0.8)
#
Tcbnu = T0*Dcbnu/D1*Bk
deltah = 1.94e-5 * Omm**(-0.785-0.05*np.log(Omm))*np.exp(-0.95*(ns-1)-0.169*(ns-1)**2.)
#The output power spectrum will be in the unit of (Mpc/h)^3.
Pk = 2*np.pi**2. * deltah**2. * (kvector)**(ns) * Tcbnu**2. * growth**2. /cosmo.Hubble(0)**(3.+ns)
return Pk
def EH98_transfer(kvector, redshift, scaling_factor, cosmo=None):
#This is the link to the paper: https://arxiv.org/pdf/astro-ph/9710252.pdf
#The input kvector should be in unit of h/Mpc after rescaling by rd.
cdict = cosmo.get_current_derived_parameters(['z_d'])
h = cosmo.h()
H_at_z = cosmo.Hubble(redshift) * conts.c /1000. /(100.*h)
Omm = cosmo.Omega_m()
Omb = cosmo.Omega_b()
#Cannot find the following function.
# Omc = cosmo.omegach2()/h**2.
Omc = cosmo.Omega0_cdm()
Omm_at_z = Omm*(1.+redshift)**3./H_at_z**2.
OmLambda_at_z = 1.-Omm_at_z
ns = cosmo.n_s()
rs = cosmo.rs_drag()*h/scaling_factor
Omnu = Omm-Omb-Omc
fnu = Omnu/Omm
fb = Omb/Omm
fnub = (Omb+Omnu)/Omm
fc = Omc/Omm
fcb = (Omc+Omb)/Omm
pc = 1./4.*(5-np.sqrt(1+24*fc))
pcb = 1./4.*(5-np.sqrt(1+24*fcb))
Neff = cosmo.Neff()
Omg = cosmo.Omega_g()
Omr = Omg * (1. + Neff * (7./8.)*(4./11.)**(4./3.))
aeq = Omr/(Omb+Omc)/(1-fnu)
zeq = 1./aeq -1.
Heq = cosmo.Hubble(zeq)/h
keq = aeq*Heq*scaling_factor
zd = cdict['z_d']
yd = (1.+zeq)/(1.+zd)
growth = cosmo.scale_independent_growth_factor(redshift)
if (fnu==0):
Nnu = 0.
else:
Nnu = 1.
#alpha_gamma = 1 - 0.328*np.log(431*Omm*h**2)*Omb/Omm + 0.38*np.log(22.3*Omm*h**2)*(Omb/Omm)**2
#There seems to be a mistake in this equation.
# alpha_nu = fc/fcb * (5 - 2*(pc+pcb))/(5-4*pcb) * (1-0.553*fnub + 0.126*fnub**3.) / (1 - 0.193*np.sqrt(fnu)*Nnu**0.2 + 0.169*fnu) \
# *(1.+yd)**(pcb-pc) * (1+(pc-pcb)/2*(1+1./(3-4*pc)/(7-4*pcb))*(1.+yd)**(-1.))
alpha_nu = fc/fcb * (5 - 2*(pc+pcb))/(5-4*pcb) * (1-0.553*fnub + 0.126*fnub**3.) / (1 - 0.193*np.sqrt(fnu*Nnu) + 0.169*fnu*(Nnu)**0.2) \
*(1.+yd)**(pcb-pc) * (1+(pc-pcb)/2.0*(1+1./(3-4*pc)/(7-4*pcb))*(1.+yd)**(-1.))
#eff_shape = (alpha_gamma + (1.-alpha_gamma)/(1+(0.43*kvector*rs)**4.))
eff_shape = (np.sqrt(alpha_nu) + (1.-np.sqrt(alpha_nu))/(1+(0.43*kvector*rs)**4.))
#add in q_test.
q_test = kvector/keq*7.46e-2
q0 = kvector/(keq/7.46e-2)/eff_shape
betac = (1.-0.949*fnub)**(-1.)
L0 = np.log(np.exp(1.)+1.84*np.sqrt(alpha_nu)*betac*q0)
# C0 = 14.4 + 325./(1+60.5*q0**1.08)
C0 = 14.4 + 325./(1+60.5*q0**1.11)
T0 = L0/(L0+C0*q0**2.)
if (fnu==0):
yfs=0.
qnu=3.92*q_test
# qnu=3.92*q0
else:
yfs = 17.2*fnu*(1.+0.488*fnu**(-7./6.))*(Nnu*q_test/fnu)**2.
qnu = 3.92*q_test*np.sqrt(Nnu/fnu)
# yfs = 17.2*fnu*(1.+0.488*fnu**(-7./6.))*(Nnu*q0/fnu)**2.
# qnu = 3.92*q0*np.sqrt(Nnu/fnu)
#The original code seems to be missing a factor of 5.
D1 = (1.+zeq)/(1.+redshift)*5*Omm_at_z/2./(Omm_at_z**(4./7.) - OmLambda_at_z + (1.+Omm_at_z/2.)*(1.+OmLambda_at_z/70.))
Dcbnu = (fcb**(0.7/pcb)+(D1/(1.+yfs))**0.7)**(pcb/0.7) * D1**(1.-pcb)
Bk = 1. + 1.24*fnu**(0.64)*Nnu**(0.3+0.6*fnu)/(qnu**(-1.6)+qnu**0.8)
#
Tcbnu = T0*Dcbnu/D1*Bk
return Tcbnu
def Primordial(k, A_s, n_s, k_p = 0.05):
#k_p is in the unit of 1/Mpc, so the input power spectrum also need to be in 1/Mpc.
P_R = A_s*(k/k_p)**(n_s - 1.0)
return P_R
def slope_at_x(xvector,yvector):
#find the slope
diff = np.diff(yvector)/np.diff(xvector)
diff = np.append(diff,diff[-1])
return diff