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
|
/*@@
@file EOS_GIF.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_GIF_Pressure(nelems, rho, eps, press)
USE EOS_GIF_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: rho, eps
CCTK_REAL, dimension(nelems), intent(out) :: press
press = (eos_if_gamma_local - 1.d0) * rho * eps
end subroutine EOS_GIF_Pressure
subroutine EOS_GIF_DPDRho(nelems, rho, eps, dpdrho)
USE EOS_GIF_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: rho, eps
CCTK_REAL, dimension(nelems), intent(out) :: dpdrho
dpdrho = (eos_if_gamma_local - 1.d0) * eps
end subroutine EOS_GIF_DPDRho
subroutine EOS_GIF_DPDIE(nelems, rho, eps, dpdie)
USE EOS_GIF_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: rho, eps
CCTK_REAL, dimension(nelems), intent(out) :: dpdie
dpdie = (eos_if_gamma_local - 1.d0) * rho
end subroutine EOS_GIF_DPDIE
subroutine EOS_GIF_cs2(nelems, rho, eps, cs2)
USE EOS_GIF_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: rho, eps
CCTK_REAL, dimension(nelems), intent(out) :: cs2
cs2 = eos_if_gamma_local * (eos_if_gamma_local - 1.d0) * eps / &
(1.d0 + eos_if_gamma_local * eps)
end subroutine EOS_GIF_cs2
subroutine EOS_GIF_Temp(nelems, rho, eps, temperature)
USE EOS_GIF_Scalars
implicit none
CCTK_INT, intent(in) :: nelems
CCTK_REAL, dimension(nelems), intent(in) :: rho, eps
CCTK_REAL, dimension(nelems), intent(out) :: temperature
temperature = (eos_if_gamma_local - 1.d0) * eps * &
mean_molecular_weight_local / &
k_boltzmann
end subroutine EOS_GIF_Temp
|