| Class | Turbulence |
| In: |
physics/turbulence.f90
|
モデルの物理過程を計算するために必要となる関数群を束ねたモジュール 具体的には以下の項を計算するための関数を格納する.
* 乱流拡散項 * 乱流エネルギーの時間発展方程式に含まれる各項 * 散逸加熱項
| 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)
|
水平速度に対する乱流拡散
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 モジュールの初期化ルーチン
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)
|
鉛直速度に対する乱流拡散
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)
|
乱流エネルギーの浮力項を計算. 乾燥大気版.
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)
|
温位の散逸加熱項
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_ExnerBasicZ * CpDry * (Cm ** 2.0d0) * (MixLen ** 4.0d0) )
!値の保管
call StoreDisp( 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)
|
乱流エネルギーの消散項
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)
|
速度シアーによる乱流エネルギー生成項
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 方向に半格子ずれた点での乱流拡散
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))
end function xz_TurbScalar
| Function : | |||
| xz_TurbScalar2(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 方向に半格子ずれた点での乱流拡散
function xz_TurbScalar2(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_TurbScalar2(DimXMin:DimXMax, DimZMin:DimZMax)
!スカラー量の水平乱流拡散
xz_TurbScalar2 = 0.0d0 !初期化
xz_TurbScalar2 = 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 StoreTurb( xz_TurbScalar2 )
end function xz_TurbScalar2