大規模凝結スキームにより, 温度と比湿を調節します.
Adjust temperature and specific humidity by large scale condensation
scheme.
subroutine LScaleCond( xyz_Temp, xyz_QVap, xyz_Press, xyr_Press )
!
! 大規模凝結スキームにより, 温度と比湿を調節します.
!
! Adjust temperature and specific humidity by
! large scale condensation scheme.
!
! モジュール引用 ; USE statements
!
! 物理定数設定
! Physical constants settings
!
use constants, only: GasRUniv, Grav, CpDry, LatentHeat, EpsV
! $ \epsilon_v $ .
! 水蒸気分子量比.
! Molecular weight of water vapor
! 時刻管理
! Time control
!
use timeset, only: DelTime, TimeN, TimesetClockStart, TimesetClockStop
! ヒストリデータ出力
! History data output
!
use gtool_historyauto, only: HistoryAutoPut
! 宣言文 ; Declaration statements
!
implicit none
real(DP), intent(inout):: xyz_Temp (0:imax-1, 1:jmax, 1:kmax)
! $ T $ . 温度. Temperature
real(DP), intent(inout):: xyz_QVap (0:imax-1, 1:jmax, 1:kmax)
! $ q $ . 比湿. Specific humidity
real(DP), intent(in):: xyz_Press (0:imax-1, 1:jmax, 1:kmax)
! $ p $ . 気圧 (整数レベル).
! Air pressure (full level)
real(DP), intent(in):: xyr_Press (0:imax-1, 1:jmax, 0:kmax)
! $ \hat{p} $ . 気圧 (半整数レベル).
! Air pressure (half level)
! 作業変数
! Work variables
!
real(DP):: xy_Rain (0:imax-1, 1:jmax)
! 降水量.
! Precipitation
real(DP):: xyz_DTempDt (0:imax-1, 1:jmax, 1:kmax)
! 温度変化率.
! Temperature tendency
real(DP):: xyz_DQVapDt (0:imax-1, 1:jmax, 1:kmax)
! 比湿変化率.
! Specific humidity tendency
real(DP):: xyz_QVapB (0:imax-1, 1:jmax, 1:kmax)
! 調節前の比湿.
! Specific humidity before adjust.
real(DP):: xyz_TempB (0:imax-1, 1:jmax, 1:kmax)
! 調節前の温度.
! Temperature before adjust.
!
real(DP):: QVapSat
! 飽和比湿.
! Saturation specific humidity.
real(DP):: DQVapSatDTemp
! $ \DD{q_{\rm{sat}}}{T} $
real(DP):: DelQVap
! 調節による比湿の変化量.
! Specific humidity variation by adjustment
real(DP):: DelTemp
! 調節による温度変化量.
! Temperature variation by adjustment
integer:: i ! 経度方向に回る DO ループ用作業変数
! Work variables for DO loop in longitude
integer:: j ! 緯度方向に回る DO ループ用作業変数
! Work variables for DO loop in latitude
integer:: k ! 鉛直方向に回る DO ループ用作業変数
! Work variables for DO loop in vertical direction
integer:: itr ! イテレーション方向に回る DO ループ用作業変数
! Work variables for DO loop in iteration direction
! 実行文 ; Executable statement
!
! 計算時間計測開始
! Start measurement of computation time
!
call TimesetClockStart( module_name )
! 初期化
! Initialization
!
if ( .not. lscond_inited ) call LSCondInit
! 調節前 "QVap", "Temp" の保存
! Store "QVap", "Temp" before adjustment
!
xyz_QVapB = xyz_QVap
xyz_TempB = xyz_Temp
! 調節
! Adjustment
!
do k = kmax, 1, -1
do i = 0, imax-1
do j = 1, jmax
! 飽和比湿計算
! Calculate saturation specific humidity
!
QVapSat = EpsV * ( P0Nha92 / xyz_Press(i,j,k) ) * exp ( - LatHeatNha92 / ( GasRUniv * xyz_Temp(i,j,k) ) )
!!$ call CalcQVapSat( phy_lsc % phy_sat, & ! (in)
!!$ & Temp = xyz_Temp(i,j,k), & ! (in)
!!$ & Press = xyz_Press(i,j,k), & ! (in)
!!$ & QVapSat = QVapSat ) ! (out)
!!$
!!$ QVapSat = &
!!$ & EpsV * ES0 &
!!$ & * exp( LatentHeat / RVap &
!!$ & * ( 1.0_DP / 273.0_DP - 1.0_DP / xyz_Temp(i,j,k) ) ) &
!!$ & / xyz_Press(i,j,k)
! 飽和していたら, 温度と比湿の変化を計算
! Calculate tendency of temperature and humidity
! if moist is saturation.
!
if ( ( xyz_QVap(i,j,k) / QVapSat ) >= CrtlRH ) then
do itr = 1, ItrtMax
! 飽和比湿計算
! Calculate saturation specific humidity
!
QVapSat = EpsV * ( P0Nha92 / xyz_Press(i,j,k) ) * exp ( - LatHeatNha92 / ( GasRUniv * xyz_Temp(i,j,k) ) )
DQVapSatDTemp = QVapSat * ( LatHeatNha92 / ( GasRUniv * xyz_Temp(i,j,k)**2 ) )
!!$ call CalcQVapSat( phy_lsc % phy_sat, & ! (in)
!!$ & Temp = xyz_Temp(i,j,k), & ! (in)
!!$ & Press = xyz_Press(i,j,k), & ! (in)
!!$ & QVapSat = QVapSat ) ! (out)
!!$
!!$ call CalcDQVapSatDTemp( phy_lsc % phy_sat, & ! (in)
!!$ & Temp = xyz_Temp(i,j,k), & ! (in)
!!$ & Press = xyz_Press(i,j,k), & ! (in)
!!$ & DQVapSatDTemp = DQVapSatDTemp ) ! (out)
!!$
!!$ QVapSat = &
!!$ & EpsV * ES0 &
!!$ & * exp( LatentHeat / RVap &
!!$ & * ( 1.0_DP / 273.0_DP - 1.0_DP / xyz_Temp(i,j,k) ) ) &
!!$ & / xyz_Press(i,j,k)
!!$
!!$ DQVapSatDTemp = &
!!$ & LatentHeat * QVapSat / ( RVap * xyz_Temp(i,j,k) * xyz_Temp(i,j,k) )
! 温度と比湿の変化分をニュートン法で求める
! Calculate variation of temperature and specific humidity
! with Newton method
!
DelTemp = LatentHeat / CpDry * ( xyz_QVap(i,j,k) - QVapSat ) / ( 1.0_DP + LatentHeat / CpDry * DQVapSatDTemp )
DelQVap = DQVapSatDTemp * DelTemp
! 温度と比湿の調節
! Adjust temperature and specific humidity
!
xyz_Temp(i,j,k) = xyz_Temp(i,j,k) + DelTemp
xyz_QVap(i,j,k) = QVapSat + DelQVap
end do
end if
end do
end do
end do
! 比湿変化率, 温度変化率, 降水量の算出
! Calculate specific humidity tendency, temperature tendency,
! precipitation
!
xy_Rain = 0.0_DP
xyz_DTempDt = 0.0_DP
xyz_DQvapDt = 0.0_DP
xyz_DQVapDt = xyz_DQVapDt + ( xyz_QVap - xyz_QVapB ) / ( 2.0_DP * DelTime )
xyz_DTempDt = xyz_DTempDt + ( xyz_Temp - xyz_TempB ) / ( 2.0_DP * DelTime )
do k = kmax, 1, -1
xy_Rain = xy_Rain + ( xyz_Temp(:,:,k) - xyz_TempB(:,:,k) ) * CpDry / ( 2.0_DP * DelTime ) * ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav
end do
! ヒストリデータ出力
! History data output
!
call HistoryAutoPut( TimeN, 'RainLsc', xy_Rain )
call HistoryAutoPut( TimeN, 'DTempDtLsc', xyz_DTempDt )
call HistoryAutoPut( TimeN, 'DQVapDtLsc', xyz_DQVapDt )
! 計算時間計測一時停止
! Pause measurement of computation time
!
call TimesetClockStop( module_name )
end subroutine LScaleCond