aboutsummaryrefslogtreecommitdiff
path: root/src/EOS_Ideal_Fluid.F
blob: 36e4addcdf35ba8b80e6deda0a25a3afce0146dc (plain)
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
c/*@@
c  @file      EOS_Ideal_Fluid.F
c  @date      December 1999
c  @author    Mark Miller
c  @desc 
c     Routines to calculate Ideal Fluid EOS through EOS_Base
c  @enddesc 
c@@*/

#include "cctk.h"
#include "cctk_Parameters.h"

c/*@@
c  @routine    EOS_Ideal_Fluid_Pressure
c  @date       December 1999
c  @author     Mark Miller
c  @desc 
c    Calculate the pressure  P=(gamma-1)*eps*rho
c  @enddesc 
c  @calls     
c  @calledby   
c  @history 
c
c  @endhistory 
c
c@@*/
      CCTK_REAL FUNCTION EOS_Ideal_Fluid_Pressure(rest_mass_density, 
     .                                            specific_internal_energy)

      IMPLICIT NONE
      DECLARE_CCTK_PARAMETERS

      CCTK_REAL rest_mass_density
      CCTK_REAL specific_internal_energy

      EOS_Ideal_Fluid_Pressure = (eos_ideal_fluid_gamma - 1.0d0) *
     .   rest_mass_density * specific_internal_energy

      END FUNCTION EOS_Ideal_Fluid_Pressure

c/*@@
c  @routine    EOS_Ideal_Fluid_SpecificIE
c  @date       December 1999
c  @author     Mark Miller
c  @desc 
c    Calculate the specific internal energy  
c       eps = press / ( (gamma-1)*rho)
c  @enddesc 
c  @calls     
c  @calledby   
c  @history 
c
c  @endhistory 
c
c@@*/
      CCTK_REAL FUNCTION EOS_Ideal_Fluid_SpecificIE(rest_mass_density,pressure)

      IMPLICIT NONE
      DECLARE_CCTK_PARAMETERS

      CCTK_REAL rest_mass_density
      CCTK_REAL pressure

      EOS_Ideal_Fluid_SpecificIE = pressure / ((eos_ideal_fluid_gamma - 1.0d0)*
     .   rest_mass_density)

      END FUNCTION EOS_Ideal_Fluid_SpecificIE

c/*@@
c  @routine    EOS_Ideal_Fluid_RestMassDens
c  @date       December 1999
c  @author     Mark Miller
c  @desc 
c    Calculate the rest mass density
c       rho = pressure / ((gamma-1)*eps)
c  @enddesc 
c  @calls     
c  @calledby   
c  @history 
c
c  @endhistory 
c
c@@*/
      CCTK_REAL FUNCTION EOS_Ideal_Fluid_RestMassDens(specific_internal_energy,
     .    pressure)

      IMPLICIT NONE
      DECLARE_CCTK_PARAMETERS

      CCTK_REAL specific_internal_energy
      CCTK_REAL pressure

      EOS_Ideal_Fluid_RestMassDens = pressure / 
     .   ((eos_ideal_fluid_gamma - 1.0d0)*specific_internal_energy)

      END FUNCTION EOS_Ideal_Fluid_RestMassDens

c/*@@
c  @routine    EOS_Ideal_Fluid_DPressByDRho
c  @date       December 1999
c  @author     Mark Miller
c  @desc 
c    Calculate d(pressure)/d(rho), keeping eps fixed:
c       dp/drho = (gamma-1)*eps
c  @enddesc 
c  @calls     
c  @calledby   
c  @history 
c
c  @endhistory 
c
c@@*/
      CCTK_REAL FUNCTION EOS_Ideal_Fluid_DPressByDRho(rest_mass_density, 
     .      specific_internal_energy)

      IMPLICIT NONE
      DECLARE_CCTK_PARAMETERS

c     use CCTK_DECLARE to silence compiler warnings
      CCTK_DECLARE(CCTK_REAL, rest_mass_density, )
      CCTK_REAL specific_internal_energy

      EOS_Ideal_Fluid_DPressByDRho = (eos_ideal_fluid_gamma - 1.0d0) *
     .   specific_internal_energy

      END FUNCTION EOS_Ideal_Fluid_DPressByDRho

c/*@@
c  @routine    EOS_Ideal_Fluid_DPressByDEps
c  @date       December 1999
c  @author     Mark Miller
c  @desc 
c    Calculate d(pressure)/d(eps), keeping rho fixed
c      dp/deps = (gamma-1)*rho
c  @enddesc 
c  @calls     
c  @calledby   
c  @history 
c
c  @endhistory 
c
c@@*/
      CCTK_REAL FUNCTION EOS_Ideal_Fluid_DPressByDEps(rest_mass_density, 
     .   specific_internal_energy)

      IMPLICIT NONE
      DECLARE_CCTK_PARAMETERS

      CCTK_REAL rest_mass_density
c     use CCTK_DECLARE to silence compiler warnings
      CCTK_DECLARE(CCTK_REAL, specific_internal_energy, )

      EOS_Ideal_Fluid_DPressByDEps = (eos_ideal_fluid_gamma - 1.0d0) *
     .   rest_mass_density

      END FUNCTION EOS_Ideal_Fluid_DPressByDEps