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
|
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
#include "cctk_Functions.h"
! eoskey:
! 1 --- polytropic EOS
! 2 --- gamma-law EOS
! 3 --- hybrid EOS
! 4 --- finite-T microphysical NSE EOS
subroutine EOS_Omni_EOS_Press_f_hrho_v2_rhoW(eoskey,keytemp,rf_precision,npoints,&
hrho,v2,rhoW,eps,temp,ye,press,keyerr,anyerr)
use EOS_Omni_Module
implicit none
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_FUNCTIONS
CCTK_INT, intent(in) :: eoskey,keytemp,npoints
CCTK_INT, intent(out) :: keyerr(npoints)
CCTK_INT, intent(out) :: anyerr
CCTK_REAL, intent(in) :: rf_precision
CCTK_REAL, intent(in) :: hrho(npoints),v2(npoints)
CCTK_REAL, intent(in) :: rhoW(npoints),ye(npoints)
CCTK_REAL, intent(inout) :: eps(npoints), temp(npoints)
CCTK_REAL, intent(out) :: press(npoints)
! local vars
integer :: i
real*8 :: gtmp1,gtmp2
character(256) :: warnstring
real*8 :: hybrid_local_gamma, hybrid_local_k_cgs, &
hybrid_p_poly, hybrid_p_th
real*8,parameter :: zero = 0.0d0
anyerr = 0
keyerr(:) = 0
select case (eoskey)
case (1)
! polytropic EOS
call CCTK_WARN(0,"Polytropic EOS not implemented for press_f_hro_v2_rhoW")
case (2)
! gamma-law EOS
if(keytemp.eq.1) then
call CCTK_WARN(0,"keytemp.eq.1 not implemented for press_f_hro_v2_rhoW")
endif
do i=1,npoints
gtmp1 = 1.0d0 - v2(i)
gtmp2 = rhoW(i)*sqrt(gtmp1)
press(i) = (gl_gamma - 1.0d0)/(gl_gamma) * &
(hrho(i)*gtmp1 - gtmp2)
eps(i) = press(i) / (gl_gamma - 1.0d0) / gtmp2
enddo
case (3)
! hybrid EOS
call CCTK_WARN(0,"Hybrid EOS not implemented for press_f_hro_v2_rhoW")
case DEFAULT
write(warnstring,*) "eoskey ",eoskey," not implemented!"
call CCTK_WARN(0,warnstring)
end select
end subroutine EOS_Omni_EOS_Press_f_hrho_v2_rhoW
subroutine EOS_Omni_EOS_dpdhrho_f_hrho_v2_rhoW(eoskey,keytemp,rf_precision,npoints,&
hrho,v2,rhoW,eps,temp,ye,dpdhrho,keyerr,anyerr)
use EOS_Omni_Module
implicit none
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_FUNCTIONS
CCTK_INT, intent(in) :: eoskey,keytemp,npoints
CCTK_INT, intent(out) :: keyerr(npoints)
CCTK_INT, intent(out) :: anyerr
CCTK_REAL, intent(in) :: rf_precision
CCTK_REAL, intent(in) :: hrho(npoints),v2(npoints)
CCTK_REAL, intent(in) :: rhoW(npoints),ye(npoints)
CCTK_REAL, intent(inout) :: eps(npoints), temp(npoints)
CCTK_REAL, intent(out) :: dpdhrho(npoints)
! local vars
integer :: i
character(256) :: warnstring
real*8 :: hybrid_local_gamma, hybrid_local_k_cgs, &
hybrid_p_poly, hybrid_p_th
real*8,parameter :: zero = 0.0d0
anyerr = 0
keyerr(:) = 0
select case (eoskey)
case (1)
! polytropic EOS
call CCTK_WARN(0,"Polytropic EOS not implemented for press_f_hro_v2_rhoW")
case (2)
! gamma-law EOS
if(keytemp.eq.1) then
call CCTK_WARN(0,"keytemp.eq.1 not implemented for press_f_hro_v2_rhoW")
endif
do i=1,npoints
dpdhrho(i) = (gl_gamma - 1.0d0) * (1.0d0 - v2(i)) / (gl_gamma)
enddo
case (3)
! hybrid EOS
call CCTK_WARN(0,"Hybrid EOS not implemented for press_f_hro_v2_rhoW")
case DEFAULT
write(warnstring,*) "eoskey ",eoskey," not implemented!"
call CCTK_WARN(0,warnstring)
end select
end subroutine EOS_Omni_EOS_dpdhrho_f_hrho_v2_rhoW
subroutine EOS_Omni_EOS_dpdv2_f_hrho_v2_rhoW(eoskey,keytemp,rf_precision,npoints,&
hrho,v2,rhoW,eps,temp,ye,dpdv2,keyerr,anyerr)
use EOS_Omni_Module
implicit none
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_FUNCTIONS
CCTK_INT, intent(in) :: eoskey,keytemp,npoints
CCTK_INT, intent(out) :: keyerr(npoints)
CCTK_INT, intent(out) :: anyerr
CCTK_REAL, intent(in) :: rf_precision
CCTK_REAL, intent(in) :: hrho(npoints),v2(npoints)
CCTK_REAL, intent(in) :: rhoW(npoints),ye(npoints)
CCTK_REAL, intent(inout) :: eps(npoints), temp(npoints)
CCTK_REAL, intent(out) :: dpdv2(npoints)
! local vars
integer :: i
character(256) :: warnstring
real*8 :: hybrid_local_gamma, hybrid_local_k_cgs, &
hybrid_p_poly, hybrid_p_th
real*8,parameter :: zero = 0.0d0
anyerr = 0
keyerr(:) = 0
select case (eoskey)
case (1)
! polytropic EOS
call CCTK_WARN(0,"Polytropic EOS not implemented for press_f_hro_v2_rhoW")
case (2)
! gamma-law EOS
if(keytemp.eq.1) then
call CCTK_WARN(0,"keytemp.eq.1 not implemented for press_f_hro_v2_rhoW")
endif
do i=1,npoints
dpdv2(i) = (gl_gamma - 1.0d0) * (0.5d0 * rhoW(i) &
/ sqrt(1.0d0 - v2(i)) - hrho(i)) / gl_gamma
enddo
case (3)
! hybrid EOS
call CCTK_WARN(0,"Hybrid EOS not implemented for press_f_hro_v2_rhoW")
case DEFAULT
write(warnstring,*) "eoskey ",eoskey," not implemented!"
call CCTK_WARN(0,warnstring)
end select
end subroutine EOS_Omni_EOS_dpdv2_f_hrho_v2_rhoW
|