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
|
/*@@
@file EOS_GP.F90
@date Mon Mar 14 16:42:15 2005
@author Ian Hawke
@desc
The functions that actually set the polytropic EOS
@enddesc
@@*/
#include "cctk.h"
#include "cctk_Parameters.h"
subroutine EOS_GP_Pressure(nelems, rho, press)
USE EOS_GP_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: rho
CCTK_REAL, dimension(nelems), intent(out) :: press
press = p_geom_factor * eos_k_cgs * &
(rho * rho_geom_factor_inv) ** eos_gamma_local
end subroutine EOS_GP_Pressure
subroutine EOS_GP_IntEn(nelems, rho, IntEn)
USE EOS_GP_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: rho
CCTK_REAL, dimension(nelems), intent(out) :: IntEn
IntEn = p_geom_factor * eos_k_cgs * &
(rho * rho_geom_factor_inv) ** eos_gamma_local / &
(rho * (eos_gamma_local - 1.d0) )
end subroutine EOS_GP_IntEn
subroutine EOS_GP_DPDRho(nelems, rho, dpdrho)
USE EOS_GP_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: rho
CCTK_REAL, dimension(nelems), intent(out) :: dpdrho
dpdrho = p_geom_factor * eos_k_cgs * &
eos_gamma_local * rho_geom_factor_inv * &
(rho * rho_geom_factor_inv) ** (eos_gamma_local - 1.d0)
end subroutine EOS_GP_DPDRho
subroutine EOS_GP_DPDIE(nelems, rho, dpdie)
USE EOS_GP_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: rho
CCTK_REAL, dimension(nelems), intent(out) :: dpdie
dpdie = 0.d0
end subroutine EOS_GP_DPDIE
subroutine EOS_GP_cs2(nelems, rho, cs2)
USE EOS_GP_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: rho
CCTK_REAL, dimension(nelems), intent(out) :: cs2
cs2 = (p_geom_factor * eos_k_cgs * eos_gamma_local * rho_geom_factor_inv * &
(rho * rho_geom_factor_inv) ** (eos_gamma_local - 1.d0)) / &
(1.d0 + p_geom_factor * eos_k_cgs * eos_gamma_local * &
rho_geom_factor_inv / (eos_gamma_local - 1.d0) * &
(rho * rho_geom_factor_inv) ** (eos_gamma_local - 1.d0))
end subroutine EOS_GP_cs2
/*@@
@routine Inverse functions
@date Mon Apr 4 11:03:04 2005
@author Ian Hawke
@desc
Give P find the results...
@enddesc
@calls
@calledby
@history
@endhistory
@@*/
subroutine EOS_GP_Inv_Rho(nelems, press, rho)
USE EOS_GP_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: press
CCTK_REAL, dimension(nelems), intent(out) :: rho
rho = rho_geom_factor * &
(press / p_geom_factor / eos_k_cgs) ** (1.d0 / eos_gamma_local)
end subroutine EOS_GP_Inv_Rho
subroutine EOS_GP_Inv_IntEn(nelems, press, IntEn)
USE EOS_GP_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: press
CCTK_REAL, dimension(nelems), intent(out) :: IntEn
IntEn = press / &
( (eos_gamma_local - 1.d0) * rho_geom_factor * &
(press / p_geom_factor / eos_k_cgs) ** (1.d0 / eos_gamma_local) )
end subroutine EOS_GP_Inv_IntEn
subroutine EOS_GP_Inv_DPDRho(nelems, press, dpdrho)
USE EOS_GP_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: press
CCTK_REAL, dimension(nelems), intent(out) :: dpdrho
dpdrho = p_geom_factor * eos_k_cgs * &
eos_gamma_local * rho_geom_factor_inv * &
(press / (p_geom_factor * eos_k_cgs) ) ** &
( (eos_gamma_local - 1.d0) / eos_gamma_local )
end subroutine EOS_GP_Inv_DPDRho
subroutine EOS_GP_Inv_DPDIE(nelems, press, dpdie)
USE EOS_GP_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: press
CCTK_REAL, dimension(nelems), intent(out) :: dpdie
dpdie = 0.d0
end subroutine EOS_GP_Inv_DPDIE
subroutine EOS_GP_Inv_cs2(nelems, press, cs2)
USE EOS_GP_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: press
CCTK_REAL, dimension(nelems), intent(out) :: cs2
cs2 = (p_geom_factor * eos_k_cgs * eos_gamma_local * rho_geom_factor_inv * &
(press / (p_geom_factor * eos_k_cgs) ) ** &
( (eos_gamma_local - 1.d0) / eos_gamma_local ) ) / &
(1.d0 + p_geom_factor * eos_k_cgs * eos_gamma_local * &
rho_geom_factor_inv / (eos_gamma_local - 1.d0) * &
(press / (p_geom_factor * eos_k_cgs) ) ** &
( (eos_gamma_local - 1.d0) / eos_gamma_local ) )
end subroutine EOS_GP_Inv_cs2
|