大規模凝結スキームにより, 温度と比湿を調節します.
Adjust temperature and specific humidity by large scale condensation
scheme.
subroutine LScaleCond( xyz_Temp, xyz_QVap, xy_Rain, xyz_DTempDt, xyz_DQVapDt, xyz_Press, xyr_Press, xyz_DQH2OLiqDt )
!
! 大規模凝結スキームにより, 温度と比湿を調節します.
!
! Adjust temperature and specific humidity by
! large scale condensation scheme.
!
! モジュール引用 ; USE statements
!
! 物理定数設定
! Physical constants settings
!
use constants, only: Grav, CpDry, LatentHeat
! $ L $ [J kg-1] .
! 凝結の潜熱.
! Latent heat of condensation
! 時刻管理
! Time control
!
use timeset, only: DelTime, TimeN, TimesetClockStart, TimesetClockStop
! ヒストリデータ出力
! History data output
!
use gtool_historyauto, only: HistoryAutoPut
! 飽和比湿の算出
! Evaluate saturation specific humidity
!
use saturate, only: xyz_CalcQVapSat, xyz_CalcDQVapSatDTemp
! 宣言文 ; 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(inout):: xy_Rain (0:imax-1, 1:jmax)
! 降水量.
! Precipitation
real(DP), intent(inout):: xyz_DTempDt (0:imax-1, 1:jmax, 1:kmax)
! 温度変化率.
! Temperature tendency
real(DP), intent(inout):: xyz_DQVapDt (0:imax-1, 1:jmax, 1:kmax)
! 比湿変化率.
! Specific humidity tendency
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)
real(DP), intent(out), optional :: xyz_DQH2OLiqDt(0:imax-1,1:jmax,1:kmax)
! 作業変数
! Work variables
!
real(DP):: xy_RainLsc (0:imax-1, 1:jmax)
! 降水量.
! Precipitation
real(DP):: xyz_DTempDtLsc (0:imax-1, 1:jmax, 1:kmax)
! 温度変化率.
! Temperature tendency
real(DP):: xyz_DQVapDtLsc (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):: xyz_QVapSat (0:imax-1, 1:jmax, 1:kmax)
! 飽和比湿.
! Saturation specific humidity.
real(DP):: xyz_DQVapSatDTemp(0:imax-1, 1:jmax, 1:kmax)
! $ \DD{q_{\rm{sat}}}{T} $
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
real(DP):: xyz_RainLSC(0:imax-1, 1:jmax, 1:kmax)
real(DP):: TempBefAdj
real(DP):: QVapBefAdj
logical:: xyz_FlagSaturated(0:imax-1, 1:jmax, 1:kmax)
! 実行文 ; Executable statement
!
! 計算時間計測開始
! Start measurement of computation time
!
call TimesetClockStart( module_name )
! 初期化
! Initialization
!
if ( .not. lscond_inited ) call LSCondInit
if ( .not. FlagUse ) return
! 調節前 "QVap", "Temp" の保存
! Store "QVap", "Temp" before adjustment
!
xyz_QVapB = xyz_QVap
xyz_TempB = xyz_Temp
! 調節
! Adjustment
!
! 飽和比湿計算
! Calculate saturation specific humidity
!
xyz_QVapSat = xyz_CalcQVapSat( xyz_Temp, xyz_Press )
do k = 1, kmax
do j = 1, jmax
do i = 0, imax-1
if ( ( xyz_QVap(i,j,k) / xyz_QVapSat(i,j,k) ) >= CrtlRH ) then
xyz_FlagSaturated(i,j,k) = .true.
else
xyz_FlagSaturated(i,j,k) = .false.
end if
end do
end do
end do
do itr = 1, ItrtMax
! 飽和比湿計算
! Calculate saturation specific humidity
!
xyz_QVapSat = xyz_CalcQVapSat ( xyz_Temp, xyz_Press )
xyz_DQVapSatDTemp = xyz_CalcDQVapSatDTemp( xyz_Temp, xyz_QVapSat )
do k = 1, kmax
do j = 1, jmax
do i = 0, imax-1
! 飽和していたら, 温度と比湿の変化を計算
! Calculate tendency of temperature and humidity
! if moist is saturation.
!
if ( xyz_FlagSaturated(i,j,k) ) then
! 温度の変化分をニュートン法で求める
! Calculate variation of temperature
!
DelTemp = LatentHeat * ( xyz_QVap(i,j,k) - xyz_QVapSat(i,j,k) ) / ( CpDry + LatentHeat * xyz_DQVapSatDTemp(i,j,k) )
!=========
! check routine
!---------
!!$ TempBefAdj = xyz_Temp(i,j,k)
!!$ QVapBefAdj = xyz_QVap(i,j,k)
!=========
! 温度と比湿の調節
! Adjust temperature and specific humidity
!
xyz_Temp(i,j,k) = xyz_Temp(i,j,k) + DelTemp
xyz_QVap(i,j,k) = xyz_QVapSat(i,j,k) + xyz_DQVapSatDTemp(i,j,k) * DelTemp
!=========
! check routine
!---------
!!$ write( 6, * ) '====='
!!$ write( 6, * ) 'Energy difference before and after adjustment and each energy'
!!$ write( 6, * ) &
!!$ & CpDry * TempBefAdj + LatentHeat * QVapBefAdj &
!!$ & - ( CpDry * xyz_Temp(i,j,k) + LatentHeat * xyz_QVap(i,j,k) ), &
!!$ & CpDry * TempBefAdj + LatentHeat * QVapBefAdj, &
!!$ & ( CpDry * xyz_Temp(i,j,k) + LatentHeat * xyz_QVap(i,j,k) )
!=========
end if
end do
end do
end do
end do
! 比湿変化率, 温度変化率, 降水量の算出
! Calculate specific humidity tendency, temperature tendency,
! precipitation
!
xyz_DQVapDtLsc = ( xyz_QVap - xyz_QVapB ) / ( 2.0_DP * DelTime )
xyz_DTempDtLsc = ( xyz_Temp - xyz_TempB ) / ( 2.0_DP * DelTime )
! OLD CODE TO BE DELETED
!
!!$ xy_RainLsc = 0.
!!$ do k = kmax, 1, -1
!!$ xy_RainLsc = xy_RainLsc &
!!$ & + ( xyz_Temp(:,:,k) - xyz_TempB(:,:,k) ) &
!!$ & * CpDry / ( 2.0_DP * DelTime ) &
!!$ & * ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav &
!!$ & / LatentHeat
!!$ end do
!!$ i = 0
!!$ j = jmax/2+1
!!$ write( 6, * ) xy_RainLsc(i,j)
xy_RainLsc = 0.0d0
do k = kmax, 1, -1
xy_RainLsc = xy_RainLsc + ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav * ( xyz_QVapB(:,:,k) - xyz_QVap(:,:,k) )
end do
xy_RainLsc = xy_RainLsc / ( 2.0_DP * DelTime )
!!$ i = 0
!!$ j = jmax/2+1
!!$ write( 6, * ) xy_RainLsc(i,j)
!!$ write( 6, * ) '---'
xy_Rain = xy_Rain + xy_RainLsc
xyz_DTempDt = xyz_DTempDt + xyz_DTempDtLsc
xyz_DQVapDt = xyz_DQVapDt + xyz_DQVapDtLsc
! ヒストリデータ出力
! History data output
!
call HistoryAutoPut( TimeN, 'RainLsc', xy_RainLsc * LatentHeat )
call HistoryAutoPut( TimeN, 'DTempDtLsc', xyz_DTempDtLsc )
call HistoryAutoPut( TimeN, 'DQVapDtLsc', xyz_DQVapDtLsc )
if ( present( xyz_DQH2OLiqDt ) ) then
xyz_DQH2OLiqDt = - xyz_DQVapDtLsc
end if
! 計算時間計測一時停止
! Pause measurement of computation time
!
call TimesetClockStop( module_name )
end subroutine LScaleCond