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
|
/*@@
@file EOS_Polytrope.F
@date Wed Mar 20 14:56:35 2002
@author Ian Hawke
@desc
Routines to calculate a polytropic EOS. This is a faked
2d version that ALWAYS assumes you want the pressure from
the density. As such the specific internal energy is not set.
@enddesc
@@*/
#include "cctk.h"
#include "cctk_Parameters.h"
/*@@
@routine EOS_Polytrope
@date Mon Aug 26 15:12:39 2002
@author Ian Hawke
@desc
@enddesc
@calls
@calledby
@history
@endhistory
@@*/
CCTK_REAL function EOS_Polytrope_Pressure(rho, eps)
USE EOS_Polytrope_Scalars
implicit none
DECLARE_CCTK_PARAMETERS
CCTK_REAL rho, eps
EOS_Polytrope_Pressure = p_geom_factor * eos_k_cgs *
. (rho * rho_geom_factor_inv) ** eos_gamma
end function EOS_Polytrope_Pressure
/*@@
@routine EOS_Polytrope
@date Mon Aug 26 15:12:39 2002
@author Ian Hawke
@desc
@enddesc
@calls
@calledby
@history
@endhistory
@@*/
CCTK_REAL function EOS_Polytrope_SpecificIE(rho, press)
USE EOS_Polytrope_Scalars
implicit none
DECLARE_CCTK_PARAMETERS
CCTK_REAL rho, press, EOS_Polytrope_Pressure
EOS_Polytrope_SpecificIE = EOS_Polytrope_Pressure(rho, press)
. / ( rho * (eos_gamma - 1.d0) )
end function EOS_Polytrope_SpecificIE
/*@@
@routine EOS_Polytrope
@date Mon Aug 26 15:12:39 2002
@author Ian Hawke
@desc
@enddesc
@calls
@calledby
@history
@endhistory
@@*/
CCTK_REAL function EOS_Polytrope_RestMassDens(eps, press)
USE EOS_Polytrope_Scalars
implicit none
DECLARE_CCTK_PARAMETERS
CCTK_REAL eps, press
EOS_Polytrope_RestMassDens = press / ((eos_gamma - 1.d0) * eps)
end function EOS_Polytrope_RestMassDens
/*@@
@routine EOS_Polytrope
@date Mon Aug 26 15:12:39 2002
@author Ian Hawke
@desc
@enddesc
@calls
@calledby
@history
@endhistory
@@*/
CCTK_REAL function EOS_Polytrope_DPressByDRho(rho, eps)
USE EOS_Polytrope_Scalars
implicit none
DECLARE_CCTK_PARAMETERS
CCTK_REAL rho, eps
EOS_Polytrope_DPressByDRho = p_geom_factor * eos_k_cgs *
. eos_gamma * rho_geom_factor_inv *
. (rho*rho_geom_factor_inv) ** (eos_gamma - 1.d0)
end function EOS_Polytrope_DPressByDRho
/*@@
@routine EOS_Polytrope
@date Mon Aug 26 15:12:39 2002
@author Ian Hawke
@desc
@enddesc
@calls
@calledby
@history
@endhistory
@@*/
CCTK_REAL function EOS_Polytrope_DPressByDEps(rho, eps)
USE EOS_Polytrope_Scalars
implicit none
DECLARE_CCTK_PARAMETERS
CCTK_REAL rho, eps
EOS_Polytrope_DPressByDEps = 0.d0
end function EOS_Polytrope_DPressByDEps
|