ittc_tpwes_waveamplitude Subroutine

public pure subroutine ittc_tpwes_waveamplitude(omega, domega, t1, hs, wave_amplitude)

12th-ITTC (1969) 双参数波能谱模型计算频段有义波幅 ITTC two-parameter wave energy spectrum

Arguments

Type IntentOptional Attributes Name
real(kind=rk), intent(in) :: omega

波浪角频率, rad/s

real(kind=rk), intent(in) :: domega

频率间隔, rad/s

real(kind=rk), intent(in) :: t1

平均周期, s, T1 = 2pim0/m1

real(kind=rk), intent(in) :: hs

有义波高, m

real(kind=rk), intent(out) :: wave_amplitude

频段有义波幅, m


Contents


Source Code

    pure subroutine ittc_tpwes_waveamplitude(omega, domega, t1, hs, wave_amplitude)
        real(kind=rk), intent(in) :: omega  !! 波浪角频率, rad/s
        real(kind=rk), intent(in) :: domega !! 频率间隔, rad/s
        real(kind=rk), intent(in) :: t1     !! 平均周期, s, T1 = 2*pi*m0/m1
        real(kind=rk), intent(in) :: hs     !! 有义波高, m
        real(kind=rk), intent(out) :: wave_amplitude !! 频段有义波幅, m
        real(kind=rk) :: s

        associate (a => 173*hs**2/t1**4, b => 691/t1**4)
            s = a*exp(-b/omega**4)/omega**5
            wave_amplitude = sqrt(2*s*domega)  ! A = sqrt(2*S_\omega*d\omega)
        end associate

    end subroutine ittc_tpwes_waveamplitude