spx_EoS.f90 Source File


Contents

Source Code


Source Code

!> author: 左志华
!> date: 2022-06-24
!>
!> 气体状态方程
module spx_eos

    use spx_kinds, only: rk
    use seakeeping_module, only: g
    use spx_configuration, only: spc_obj
    implicit none

    private
    public :: eow_obj, p_water_tait, p_water_monaghan, p_water_morris, p_ideal_gas, B

    real(rk), parameter :: rho0 = 1000.0_rk
    integer, parameter :: gamma = 7
    real(rk) :: B

    abstract interface

        !> EoS for water <br>
        !> 水的状态方程
        pure subroutine p_water_fcn(rho, p, c)
            import :: rk
            real(rk), intent(in) :: rho     !! density <br>
                                            !! 密度
            real(rk), intent(out) :: p      !! pressure <br>
                                            !! 压力
            real(rk), intent(out) :: c      !! sound speed <br>
                                            !! 声速
        end subroutine p_water_fcn

    end interface

    !> 液体状态方程
    type eos_water
        procedure(p_water_fcn), pointer, nopass :: p_water
    contains
        procedure :: init
    end type eos_water

    type(eos_water) :: eow_obj  !! 水的状态方程

contains

    !> 初始化
    subroutine init(self)
        class(eos_water), intent(inout) :: self

        select case (spc_obj%eos_water)
        case ('tait')
            self%p_water => p_water_tait
        case ('monaghan')
            self%p_water => p_water_monaghan
        case ('morris')
            self%p_water => p_water_morris
        case default
            self%p_water => p_water_tait
        end select

    end subroutine init

    !> 应用状态方程通过密度和能量计算理想气体压力和声速
    pure subroutine p_ideal_gas(rho, u, p, c)
        real(rk), intent(in) :: rho     !! 密度
        real(rk), intent(in) :: u       !! 内部能量
        real(rk), intent(out) :: p      !! 压力
        real(rk), intent(out) :: c      !! 声速
        real(rk), parameter :: gamma = 1.4_rk

        associate (tmp => (gamma - 1.0_rk)*u)
            p = rho*tmp
            c = sqrt(gamma*tmp)
        end associate

    end subroutine p_ideal_gas

    !> Suitable for free surface flow <br>
    !> Tait simplified equation of state, 适用于自由面流动 <br>
    !> 孙鹏楠, 2018, 博士毕业论文, equ(2-12)
    pure subroutine p_water_tait(rho, p, c)
        real(rk), intent(in) :: rho     !! density <br>
                                        !! 密度
        real(rk), intent(out) :: p      !! pressure <br>
                                        !! 压力
        real(rk), intent(out) :: c      !! sound speed <br>
                                        !! 声速

        c = spc_obj%c0
        p = c*c*(rho - rho0)

    end subroutine p_water_tait

    !> Not suitable for free surface flow <br>
    !> 不适用于自由面流动
    !> @note 密度决定液体压强
    pure subroutine p_water_morris(rho, p, c)
        real(rk), intent(in) :: rho     !! density <br>
                                        !! 密度
        real(rk), intent(out) :: p      !! pressure <br>
                                        !! 压力
        real(rk), intent(out) :: c      !! sound speed <br>
                                        !! 声速

        c = spc_obj%c0
        p = c*c*rho

    end subroutine p_water_morris

    !> Suitable for free surface flow <br>
    !> 自由面流动
    pure subroutine p_water_monaghan(rho, p, c)
        real(rk), intent(in) :: rho     !! density <br>
                                        !! 密度
        real(rk), intent(out) :: p      !! pressure <br>
                                        !! 压力
        real(rk), intent(out) :: c      !! sound speed <br>
                                        !! 声速

        p = B*((rho/rho0)**gamma - 1)
        c = spc_obj%c0*((rho/rho0)**3)

    end subroutine p_water_monaghan

    !> Calculate sound speed <br>
    !> 计算声速
    pure subroutine get_sound_speed(h, c, d)
        real(rk), intent(in) :: h       !! specific enthalpy <br>
                                        !! 特定热量
        real(rk), intent(out) :: c, d

        d = 2*g*h
        c = 2*sqrt(d)

    end subroutine get_sound_speed

end module spx_eos