-
Notifications
You must be signed in to change notification settings - Fork 159
/
Copy pathDarkEnergyFluid.f90
293 lines (238 loc) · 11.3 KB
/
DarkEnergyFluid.f90
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
module DarkEnergyFluid
use DarkEnergyInterface
use results
use constants
use classes
implicit none
type, extends(TDarkEnergyEqnOfState) :: TDarkEnergyFluid
!comoving sound speed is always exactly 1 for quintessence
!(otherwise assumed constant, though this is almost certainly unrealistic)
contains
procedure :: ReadParams => TDarkEnergyFluid_ReadParams
procedure, nopass :: PythonClass => TDarkEnergyFluid_PythonClass
procedure, nopass :: SelfPointer => TDarkEnergyFluid_SelfPointer
procedure :: Init =>TDarkEnergyFluid_Init
procedure :: PerturbedStressEnergy => TDarkEnergyFluid_PerturbedStressEnergy
procedure :: PerturbationEvolve => TDarkEnergyFluid_PerturbationEvolve
end type TDarkEnergyFluid
!Example implementation of fluid model using specific analytic form
!(approximate effective axion fluid model from arXiv:1806.10608, with c_s^2=1 if n=infinity (w_n=1))
!This is an example, it's not supposed to be a rigorous model! (not very well tested)
type, extends(TDarkEnergyModel) :: TAxionEffectiveFluid
real(dl) :: w_n = 1._dl !Effective equation of state when oscillating
real(dl) :: fde_zc = 0._dl ! energy density fraction at a_c (not the same as peak dark energy fraction)
real(dl) :: zc !transition redshift (scale factor a_c)
real(dl) :: theta_i = const_pi/2 !Initial value
!om is Omega of the early DE component today (assumed to be negligible compared to omega_lambda)
!omL is the lambda component of the total dark energy omega
real(dl), private :: a_c, pow, om, omL, acpow, freq, n !cached internally
contains
procedure :: ReadParams => TAxionEffectiveFluid_ReadParams
procedure, nopass :: PythonClass => TAxionEffectiveFluid_PythonClass
procedure, nopass :: SelfPointer => TAxionEffectiveFluid_SelfPointer
procedure :: Init => TAxionEffectiveFluid_Init
procedure :: w_de => TAxionEffectiveFluid_w_de
procedure :: grho_de => TAxionEffectiveFluid_grho_de
procedure :: PerturbedStressEnergy => TAxionEffectiveFluid_PerturbedStressEnergy
procedure :: PerturbationEvolve => TAxionEffectiveFluid_PerturbationEvolve
end type TAxionEffectiveFluid
contains
subroutine TDarkEnergyFluid_ReadParams(this, Ini)
use IniObjects
class(TDarkEnergyFluid) :: this
class(TIniFile), intent(in) :: Ini
call this%TDarkEnergyEqnOfState%ReadParams(Ini)
this%cs2_lam = Ini%Read_Double('cs2_lam', 1.d0)
end subroutine TDarkEnergyFluid_ReadParams
function TDarkEnergyFluid_PythonClass()
character(LEN=:), allocatable :: TDarkEnergyFluid_PythonClass
TDarkEnergyFluid_PythonClass = 'DarkEnergyFluid'
end function TDarkEnergyFluid_PythonClass
subroutine TDarkEnergyFluid_SelfPointer(cptr,P)
use iso_c_binding
Type(c_ptr) :: cptr
Type (TDarkEnergyFluid), pointer :: PType
class (TPythonInterfacedClass), pointer :: P
call c_f_pointer(cptr, PType)
P => PType
end subroutine TDarkEnergyFluid_SelfPointer
subroutine TDarkEnergyFluid_Init(this, State)
use classes
class(TDarkEnergyFluid), intent(inout) :: this
class(TCAMBdata), intent(in), target :: State
call this%TDarkEnergyEqnOfState%Init(State)
if (this%is_cosmological_constant) then
this%num_perturb_equations = 0
else
if (this%use_tabulated_w) then
if (any(this%equation_of_state%F<-1) .and. any(this%equation_of_state%F>-1)) &
error stop 'Fluid dark energy model does not allow w crossing -1'
elseif (this%wa/=0 .and. &
((1+this%w_lam < -1.e-6_dl) .or. 1+this%w_lam + this%wa < -1.e-6_dl)) then
error stop 'Fluid dark energy model does not allow w crossing -1'
end if
this%num_perturb_equations = 2
end if
end subroutine TDarkEnergyFluid_Init
subroutine TDarkEnergyFluid_PerturbedStressEnergy(this, dgrhoe, dgqe, &
a, dgq, dgrho, grho, grhov_t, w, gpres_noDE, etak, adotoa, k, kf1, ay, ayprime, w_ix)
class(TDarkEnergyFluid), intent(inout) :: this
real(dl), intent(out) :: dgrhoe, dgqe
real(dl), intent(in) :: a, dgq, dgrho, grho, grhov_t, w, gpres_noDE, etak, adotoa, k, kf1
real(dl), intent(in) :: ay(*)
real(dl), intent(inout) :: ayprime(*)
integer, intent(in) :: w_ix
if (this%no_perturbations) then
dgrhoe=0
dgqe=0
else
dgrhoe = ay(w_ix) * grhov_t
dgqe = ay(w_ix + 1) * grhov_t * (1 + w)
end if
end subroutine TDarkEnergyFluid_PerturbedStressEnergy
subroutine TDarkEnergyFluid_PerturbationEvolve(this, ayprime, w, w_ix, &
a, adotoa, k, z, y)
class(TDarkEnergyFluid), intent(in) :: this
real(dl), intent(inout) :: ayprime(:)
real(dl), intent(in) :: a, adotoa, w, k, z, y(:)
integer, intent(in) :: w_ix
real(dl) Hv3_over_k, loga
Hv3_over_k = 3*adotoa* y(w_ix + 1) / k
!density perturbation
ayprime(w_ix) = -3 * adotoa * (this%cs2_lam - w) * (y(w_ix) + (1 + w) * Hv3_over_k) &
- (1 + w) * k * y(w_ix + 1) - (1 + w) * k * z
if (this%use_tabulated_w) then
!account for derivatives of w
loga = log(a)
if (loga > this%equation_of_state%Xmin_interp .and. loga < this%equation_of_state%Xmax_interp) then
ayprime(w_ix) = ayprime(w_ix) - adotoa*this%equation_of_state%Derivative(loga)* Hv3_over_k
end if
elseif (this%wa/=0) then
ayprime(w_ix) = ayprime(w_ix) + Hv3_over_k*this%wa*adotoa*a
end if
!velocity
if (abs(w+1) > 1e-6) then
ayprime(w_ix + 1) = -adotoa * (1 - 3 * this%cs2_lam) * y(w_ix + 1) + &
k * this%cs2_lam * y(w_ix) / (1 + w)
else
ayprime(w_ix + 1) = 0
end if
end subroutine TDarkEnergyFluid_PerturbationEvolve
subroutine TAxionEffectiveFluid_ReadParams(this, Ini)
use IniObjects
class(TAxionEffectiveFluid) :: this
class(TIniFile), intent(in) :: Ini
call this%TDarkEnergyModel%ReadParams(Ini)
if (Ini%HasKey('AxionEffectiveFluid_a_c')) then
error stop 'AxionEffectiveFluid inputs changed to AxionEffectiveFluid_fde_zc and AxionEffectiveFluid_zc'
end if
this%w_n = Ini%Read_Double('AxionEffectiveFluid_w_n')
this%fde_zc = Ini%Read_Double('AxionEffectiveFluid_fde_zc')
this%zc = Ini%Read_Double('AxionEffectiveFluid_zc')
call Ini%Read('AxionEffectiveFluid_theta_i', this%theta_i)
end subroutine TAxionEffectiveFluid_ReadParams
function TAxionEffectiveFluid_PythonClass()
character(LEN=:), allocatable :: TAxionEffectiveFluid_PythonClass
TAxionEffectiveFluid_PythonClass = 'AxionEffectiveFluid'
end function TAxionEffectiveFluid_PythonClass
subroutine TAxionEffectiveFluid_SelfPointer(cptr,P)
use iso_c_binding
Type(c_ptr) :: cptr
Type (TAxionEffectiveFluid), pointer :: PType
class (TPythonInterfacedClass), pointer :: P
call c_f_pointer(cptr, PType)
P => PType
end subroutine TAxionEffectiveFluid_SelfPointer
subroutine TAxionEffectiveFluid_Init(this, State)
use classes
class(TAxionEffectiveFluid), intent(inout) :: this
class(TCAMBdata), intent(in), target :: State
real(dl) :: grho_rad, F, p, mu, xc, n
select type(State)
class is (CAMBdata)
this%is_cosmological_constant = this%fde_zc==0
this%pow = 3*(1+this%w_n)
this%a_c = 1/(1+this%zc)
this%acpow = this%a_c**this%pow
!Omega in early de at z=0
this%om = 2*this%fde_zc/(1-this%fde_zc)*&
(State%grho_no_de(this%a_c)/this%a_c**4/State%grhocrit + State%Omega_de)/(1 + 1/this%acpow)
this%omL = State%Omega_de - this%om !Omega_de is total dark energy density today
this%num_perturb_equations = 2
if (this%w_n < 0.9999) then
! n <> infinity
!get (very) approximate result for sound speed parameter; arXiv:1806.10608 Eq 30 (but mu may not exactly agree with what they used)
n = nint((1+this%w_n)/(1-this%w_n))
!Assume radiation domination, standard neutrino model; H0 factors cancel
grho_rad = (kappa/c**2*4*sigma_boltz/c**3*State%CP%tcmb**4*Mpc**2*(1+default_nnu*7._dl/8*(4._dl/11)**(4._dl/3)))
xc = this%a_c**2/2/sqrt(grho_rad/3)
F=7./8
p=1./2
mu = 1/xc*(1-cos(this%theta_i))**((1-n)/2.)*sqrt((1-F)*(6*p+2)*this%theta_i/n/sin(this%theta_i))
this%freq = mu*(1-cos(this%theta_i))**((n-1)/2.)* &
sqrt(const_pi)*Gamma((n+1)/(2.*n))/Gamma(1+0.5/n)*2.**(-(n**2+1)/(2.*n))*3.**((1./n-1)/2)*this%a_c**(-6./(n+1)+3) &
*( this%a_c**(6*n/(n+1.))+1)**(0.5*(1./n-1))
this%n = n
end if
end select
end subroutine TAxionEffectiveFluid_Init
function TAxionEffectiveFluid_w_de(this, a)
class(TAxionEffectiveFluid) :: this
real(dl) :: TAxionEffectiveFluid_w_de
real(dl), intent(IN) :: a
real(dl) :: rho, apow, acpow
apow = a**this%pow
acpow = this%acpow
rho = this%omL+ this%om*(1+acpow)/(apow+acpow)
TAxionEffectiveFluid_w_de = this%om*(1+acpow)/(apow+acpow)**2*(1+this%w_n)*apow/rho - 1
end function TAxionEffectiveFluid_w_de
function TAxionEffectiveFluid_grho_de(this, a) !relative density (8 pi G a^4 rho_de /grhov)
class(TAxionEffectiveFluid) :: this
real(dl) :: TAxionEffectiveFluid_grho_de, apow
real(dl), intent(IN) :: a
if(a == 0.d0)then
TAxionEffectiveFluid_grho_de = 0.d0
else
apow = a**this%pow
TAxionEffectiveFluid_grho_de = (this%omL*(apow+this%acpow)+this%om*(1+this%acpow))*a**4 &
/((apow+this%acpow)*(this%omL+this%om))
endif
end function TAxionEffectiveFluid_grho_de
subroutine TAxionEffectiveFluid_PerturbationEvolve(this, ayprime, w, w_ix, &
a, adotoa, k, z, y)
class(TAxionEffectiveFluid), intent(in) :: this
real(dl), intent(inout) :: ayprime(:)
real(dl), intent(in) :: a, adotoa, w, k, z, y(:)
integer, intent(in) :: w_ix
real(dl) Hv3_over_k, deriv, apow, acpow, cs2, fac
if (this%w_n < 0.9999) then
fac = 2*a**(2-6*this%w_n)*this%freq**2
cs2 = (fac*(this%n-1) + k**2)/(fac*(this%n+1) + k**2)
else
cs2 = 1
end if
apow = a**this%pow
acpow = this%acpow
Hv3_over_k = 3*adotoa* y(w_ix + 1) / k
! dw/dlog a/(1+w)
deriv = (acpow**2*(this%om+this%omL)+this%om*acpow-apow**2*this%omL)*this%pow &
/((apow+acpow)*(this%omL*(apow+acpow)+this%om*(1+acpow)))
!density perturbation
ayprime(w_ix) = -3 * adotoa * (cs2 - w) * (y(w_ix) + Hv3_over_k) &
- k * y(w_ix + 1) - (1 + w) * k * z - adotoa*deriv* Hv3_over_k
!(1+w)v
ayprime(w_ix + 1) = -adotoa * (1 - 3 * cs2 - deriv) * y(w_ix + 1) + &
k * cs2 * y(w_ix)
end subroutine TAxionEffectiveFluid_PerturbationEvolve
subroutine TAxionEffectiveFluid_PerturbedStressEnergy(this, dgrhoe, dgqe, &
a, dgq, dgrho, grho, grhov_t, w, gpres_noDE, etak, adotoa, k, kf1, ay, ayprime, w_ix)
class(TAxionEffectiveFluid), intent(inout) :: this
real(dl), intent(out) :: dgrhoe, dgqe
real(dl), intent(in) :: a, dgq, dgrho, grho, grhov_t, w, gpres_noDE, etak, adotoa, k, kf1
real(dl), intent(in) :: ay(*)
real(dl), intent(inout) :: ayprime(*)
integer, intent(in) :: w_ix
dgrhoe = ay(w_ix) * grhov_t
dgqe = ay(w_ix + 1) * grhov_t
end subroutine TAxionEffectiveFluid_PerturbedStressEnergy
end module DarkEnergyFluid