Return speed of sound.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(eos_compressible), | intent(in) | :: | self | Equation of state. |
||
| real(kind=R8P), | intent(in) | :: | density | Density value. |
||
| real(kind=R8P), | intent(in) | :: | pressure | Pressure value. |
Speed of sound value.
elemental function density(self, energy, pressure, speed_of_sound, temperature) result(density_)
!< Return density.
class(eos_compressible), intent(in) :: self !< Equation of state.
real(R8P), intent(in), optional :: energy !< Specific internal energy value.
real(R8P), intent(in), optional :: pressure !< Pressure value.
real(R8P), intent(in), optional :: speed_of_sound !< Speed of sound value.
real(R8P), intent(in), optional :: temperature !< Temperature value.
real(R8P) :: density_ !< Density value.
density_ = 0._R8P
if (present(energy).and.present(pressure)) then
density_ = pressure / ((self%g_ - 1._R8P) * energy)
elseif (present(pressure).and.present(speed_of_sound)) then
density_ = self%g_ * pressure / (speed_of_sound * speed_of_sound)
elseif (present(pressure).and.present(temperature)) then
density_ = pressure / (self%R_ * temperature)
endif
endfunction density
pure function description(self, prefix) result(desc)
!< Return a pretty-formatted object description.
class(eos_compressible), intent(in) :: self !< Equation of state.
character(*), intent(in), optional :: prefix !< Prefixing string.
character(len=:), allocatable :: desc !< Description.
character(len=:), allocatable :: prefix_ !< Prefixing string, local variable.
character(len=1), parameter :: NL=new_line('a') !< New line character.
prefix_ = '' ; if (present(prefix)) prefix_ = prefix
desc = ''
desc = desc//prefix_//'cp = '//trim(str(n=self%cp_))//NL
desc = desc//prefix_//'cv = '//trim(str(n=self%cv_))
endfunction description
elemental function energy(self, density, pressure, temperature) result(energy_)
!< Return specific internal energy.
class(eos_compressible), intent(in) :: self !< Equation of state.
real(R8P), intent(in), optional :: density !< Density value.
real(R8P), intent(in), optional :: pressure !< Pressure value.
real(R8P), intent(in), optional :: temperature !< Temperature value.
real(R8P) :: energy_ !< Energy value.
energy_ = 0._R8P
if (present(density).and.present(pressure)) then
energy_ = pressure / ((self%g_ - 1._R8P) * density)
elseif (present(temperature)) then
energy_ = self%cv() * temperature
endif
endfunction energy
elemental function eta(self) result(eta_)
!< Return `2 * gamma / (gamma - 1)`.
class(eos_compressible), intent(in) :: self !< Equation of state.
real(R8P) :: eta_ !< `2 * gamma / (gamma - 1)` value.
eta_ = self%eta_
endfunction eta
elemental function g(self) result(g_)
!< Return specific heats ratio `gamma=cp/cv`.
class(eos_compressible), intent(in) :: self !< Equation of state.
real(R8P) :: g_ !< Specific heats ratio value.
g_ = self%g_
endfunction g
elemental function gm1(self) result(gm1_)
!< Return `gamma - 1`.
class(eos_compressible), intent(in) :: self !< Equation of state.
real(R8P) :: gm1_ !< `gamma - 1` value.
gm1_ = self%gm1_
endfunction gm1
elemental function gp1(self) result(gp1_)
!< Return `gamma + 1`.
class(eos_compressible), intent(in) :: self !< Equation of state.
real(R8P) :: gp1_ !< `gamma + 1` value.
gp1_ = self%gp1_
endfunction gp1
elemental function pressure(self, density, energy, temperature) result(pressure_)
!< Return pressure.
class(eos_compressible), intent(in) :: self !< Equation of state.
real(R8P), intent(in), optional :: density !< Density value.
real(R8P), intent(in), optional :: energy !< Specific internal energy value.
real(R8P), intent(in), optional :: temperature !< Temperature value.
real(R8P) :: pressure_ !< Pressure value.
pressure_ = 0._R8P
if (present(density).and.present(energy)) then
pressure_ = density * (self%g_ - 1._R8P) * energy
elseif (present(density).and.present(temperature)) then
pressure_ = density * self%R_ * temperature
endif
endfunction pressure
elemental function R(self) result(R_)
!< Return fluid constant `R=cp-cv`.
class(eos_compressible), intent(in) :: self !< Equation of state.
real(R8P) :: R_ !< Fluid constant value.
R_ = self%R_
endfunction R
elemental function speed_of_sound(self, density, pressure) result(speed_of_sound_)
!< Return speed of sound.
class(eos_compressible), intent(in) :: self !< Equation of state.
real(R8P), intent(in) :: density !< Density value.
real(R8P), intent(in) :: pressure !< Pressure value.
real(R8P) :: speed_of_sound_ !< Speed of sound value.
speed_of_sound_ = sqrt(self%g_ * pressure / density)
endfunction speed_of_sound