Class Turbulence
In: physics/turbulence.f90

モデルの物理過程を計算するために必要となる関数群を束ねたモジュール 具体的には以下の項を計算するための関数を格納する.

 * 乱流拡散項
 * 乱流エネルギーの時間発展方程式に含まれる各項
 * 散逸加熱項

Methods

Included Modules

gridset basicset average differentiate_center2 StorePotTemp StoreMixRt

Public Instance methods

Subroutine :
pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
xz_PotTemp(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
xz_Km(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
xz_KmA(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)

乱流エネルギーを定常とした場合の渦拡散係数, 渦粘性係数を求める. この乱流パラメタリゼーションは, Mellor and Yamada (1974) の Level 1 Closure に対応するが, Level1 Closure に存在する bar{theta^2} は無視されている.

[Source]

  subroutine EddyViscosity(pz_VelX, xr_VelZ, xz_PotTemp, xz_Km, xz_KmA)
    !
    ! 乱流エネルギーを定常とした場合の渦拡散係数, 渦粘性係数を求める. 
    ! この乱流パラメタリゼーションは, Mellor and Yamada (1974) の
    ! Level 1 Closure に対応するが, Level1 Closure に存在する 
    ! \bar{\theta^2} は無視されている. 

    !--- 暗黙の型宣言禁止
    implicit none

    real(8), intent(in)  :: pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8), intent(in)  :: xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8), intent(in)  :: xz_PotTemp(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(in)  :: xz_Km(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8), intent(out) :: xz_KmA(DimXMin:DimXMax, DimZMin:DimZMax) 

    xz_KmA = sqrt ( max( 0.0d0, ( xz_ShearKm(xz_Km, pz_VelX, xr_VelZ) + xz_BuoyKm(xz_PotTemp) ) * MixLen ** 2.0d0 ) )

  end subroutine EddyViscosity
Function :
pz_TurbVelX(DimXMin:DimXMax, DimZMin:DimZMax) :real(8)
: スカラー量の水平乱流拡散
xz_Km(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 乱流拡散係数
pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 水平速度
xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 鉛直速度

水平速度に対する乱流拡散

[Source]

  function pz_TurbVelX(xz_Km, pz_VelX, xr_VelZ)
    !
    ! 水平速度に対する乱流拡散
    !
    
    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8), intent(in) :: pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !水平速度
    real(8), intent(in) :: xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !鉛直速度
    real(8), intent(in) :: xz_Km(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !乱流拡散係数
    real(8)             :: pz_TurbVelX(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !スカラー量の水平乱流拡散
  
!    pz_TurbVelX = 0.0d0  !初期化  
    pz_TurbVelX = 2.0d0 * pz_dx_xz( xz_Km * xz_dx_pz( pz_VelX ) ) + pz_dz_pr( pr_avr_xz( xz_Km ) * pr_dx_xr( xr_VelZ ) + pr_avr_xz( xz_Km ) * pr_dz_pz( pz_VelX ) ) - 2.0d0 * pz_dx_xz( ( xz_Km ** 2.0d0 ) ) / ( 3.0d0 * ( Cm ** 2.0d0 ) * ( MixLen ** 2.0d0 ) )
    
  end function pz_TurbVelX
Subroutine :

Turbulence モジュールの初期化ルーチン

[Source]

  subroutine turbulence_init()
    !
    ! Turbulence モジュールの初期化ルーチン
    ! 
    
    !暗黙の型宣言禁止
    implicit none

    !混合距離
!    MixLen = sqrt(DelX * DelZ) 
    MixLen = min(DelX,  DelZ) 

  end subroutine turbulence_init
Function :
xr_TurbVelZ(DimXMin:DimXMax, DimZMin:DimZMax) :real(8)
: スカラー量の水平乱流拡散
xz_Km(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 乱流拡散係数
pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 水平速度
xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 鉛直速度

鉛直速度に対する乱流拡散

[Source]

  function xr_TurbVelZ(xz_Km, pz_VelX, xr_VelZ)
    !
    !鉛直速度に対する乱流拡散
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8), intent(in) :: pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !水平速度
    real(8), intent(in) :: xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !鉛直速度
    real(8), intent(in) :: xz_Km(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !乱流拡散係数
    real(8)             :: xr_TurbVelZ(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !スカラー量の水平乱流拡散
  
!    xr_TurbVelZ = 0.0d0  !初期化
    xr_TurbVelZ = + 2.0d0 * xr_dz_xz( xz_Km * xz_dz_xr( xr_VelZ ) ) + xr_dx_pr( pr_avr_xz( xz_Km ) * pr_dz_pz( pz_VelX ) + pr_avr_xz( xz_Km ) * pr_dx_xr( xr_VelZ ) ) - 2.0d0 * xr_dz_xz( ( xz_Km ** 2.0d0 ) ) / ( 3.0d0 * ( Cm ** 2.0d0 ) * ( MixLen ** 2.0d0 ) )
    
  end function xr_TurbVelZ
Function :
xz_BuoyKm(DimXMin:DimXMax, DimZMin:DimZMax) :real(8)
: 浮力項
xz_PotTemp(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 温位擾乱

乱流エネルギーの浮力項を計算. 乾燥大気版.

[Source]

  function xz_BuoyKm(xz_PotTemp)
    !
    !乱流エネルギーの浮力項を計算. 乾燥大気版. 
    !
    
    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8), intent(in)  :: xz_PotTemp(DimXMin:DimXMax, DimZMin:DimZMax)
                                               !温位擾乱
    real(8)              :: xz_BuoyKm(DimXMin:DimXMax, DimZMin:DimZMax)
                                               !浮力項

    !浮力項の計算
!    xz_BuoyKm = 0.0d0
    xz_BuoyKm = - 3.0d0 * Grav * ( Cm ** 2.0d0 ) * ( MixLen ** 2.0d0 ) * xz_avr_xr( xr_dz_xz( xz_PotTemp + xz_PotTempBasicZ ) ) / ( 2.0d0 * xz_PotTempBasicZ )
    
  end function xz_BuoyKm
Function :
xz_DispHeat(DimXMin:DimXMax, DimZMin:DimZMax) :real(8)
: 乱流エネルギーの消散
xz_Km(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 渦粘性係数

温位の散逸加熱項

[Source]

  function xz_DispHeat(xz_Km)
    !
    !温位の散逸加熱項
    !
    
    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)             :: xz_DispHeat(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !乱流エネルギーの消散
    real(8), intent(in) :: xz_Km(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !渦粘性係数

!    xz_DispHeat = 0.0d0
    xz_DispHeat = (xz_Km ** 3.0d0) * xz_EffMolWtBasicZ / (xz_ExnerBasicZ * CpDry * (Cm ** 2.0d0) * (MixLen ** 4.0d0))

    !値の保管
    call StorePotTempDisp( xz_DispHeat )
    
  end function xz_DispHeat
Function :
xz_DispKm(DimXMin:DimXMax, DimZMin:DimZMax) :real(8)
: 乱流エネルギーの消散
xz_Km(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 渦粘性係数

乱流エネルギーの消散項

[Source]

  function xz_DispKm(xz_Km)
    !
    !乱流エネルギーの消散項
    !
    
    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)             :: xz_DispKm(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !乱流エネルギーの消散
    real(8), intent(in) :: xz_Km(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !渦粘性係数

!    xz_DispKm = 0.0d0
    xz_DispKm = - (xz_Km ** 2.0d0) * 5.0d-1 / (MixLen ** 2.0d0)
    
  end function xz_DispKm
Function :
xz_ShearKm(DimXMin:DimXMax, DimZMin:DimZMax) :real(8)
: 渦粘性係数の移流
xz_Km(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 渦粘性係数
pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 水平速度
xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 鉛直速度

速度シアーによる乱流エネルギー生成項

[Source]

  function xz_ShearKm(xz_Km, pz_VelX, xr_VelZ)
    !
    !速度シアーによる乱流エネルギー生成項
    !
    
    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)             :: xz_ShearKm(DimXMin:DimXMax, DimZMin:DimZMax)
                                                        !渦粘性係数の移流
    real(8), intent(in) :: xz_Km(DimXMin:DimXMax, DimZMin:DimZMax)
                                                        !渦粘性係数
    real(8), intent(in) :: pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax)
                                                        !水平速度
    real(8), intent(in) :: xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax)
                                                        !鉛直速度
  
!    xz_ShearKm = 0.0d0
    xz_ShearKm = ( Cm ** 2.0d0 ) * ( MixLen ** 2.0d0 ) * ( (xz_dx_pz(pz_VelX)) ** 2.0d0 + (xz_dz_xr(xr_VelZ)) ** 2.0d0 + 5.0d-1 * ( ( xz_avr_pr(pr_dz_pz(pz_VelX)) + xz_avr_pr(pr_dx_xr(xr_VelZ)) ) ** 2.0d0 ) ) - xz_Km * (xz_dx_pz(pz_VelX) + xz_dz_xr(xr_VelZ)) / 3.0d0
    
  end function xz_ShearKm
Function :
xz_TurbScalar(DimXMin:DimXMax, DimZMin:DimZMax) :real(8)
: スカラー量の水平乱流拡散
xz_Var(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: スカラー量
xz_Kh(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 乱流拡散係数

x, z 方向に半格子ずれた点での乱流拡散

[Source]

  function xz_TurbScalar(xz_Var, xz_Kh)
    !
    ! x, z 方向に半格子ずれた点での乱流拡散
    !
    
    !暗黙の型宣言禁止
    implicit none

    !変数定義
    real(8), intent(in) :: xz_Var(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !スカラー量
    real(8), intent(in) :: xz_Kh(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !乱流拡散係数
    real(8)             :: xz_TurbScalar(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !スカラー量の水平乱流拡散
  
!    xz_TurbScalar = 0.0d0  !初期化
    xz_TurbScalar = xz_dx_pz(pz_avr_xz(xz_Kh) * pz_dx_xz(xz_Var)) + xz_dz_xr(xr_avr_xz(xz_Kh) * xr_dz_xz(xz_Var))

    call StorePotTempTurb( xz_TurbScalar )    
    
  end function xz_TurbScalar
Function :
xza_TurbScalar(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum) :real(8)
: スカラー量の水平乱流拡散
xza_Var(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum) :real(8), intent(in)
: スカラー量
xz_Kh(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
: 乱流拡散係数

x, z 方向に半格子ずれた点での乱流拡散

[Source]

  function xza_TurbScalar(xza_Var, xz_Kh)
    !
    ! x, z 方向に半格子ずれた点での乱流拡散
    !
    
    !暗黙の型宣言禁止
    implicit none

    !変数定義
    real(8), intent(in) :: xza_Var(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
                                                    !スカラー量
    real(8), intent(in) :: xz_Kh(DimXMin:DimXMax, DimZMin:DimZMax)
                                                    !乱流拡散係数
    real(8)             :: xza_TurbScalar(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
                                                    !スカラー量の水平乱流拡散
    integer             :: s
  
    do s = 1, SpcNum
      xza_TurbScalar(:,:,s) = xz_dx_pz(pz_avr_xz(xz_Kh) * pz_dx_xz(xza_Var(:,:,s))) + xz_dz_xr(xr_avr_xz(xz_Kh) * xr_dz_xz(xza_Var(:,:,s)))
    end do

    call StoreMixRtTurb( xza_TurbScalar )    
    
  end function xza_TurbScalar

[Validate]