subroutine physics_lscond( xyz_Temp, xyz_Qvap, xy_LscRain, xyz_DLscTempDt, xyz_DLscQvapDt, xyz_Press, xyr_Press, DelTimePhy )
!==== Dependency
use type_mod, only: REKIND, DBKIND, INTKIND, TOKEN, STRING
use grid_3d_mod, only: im, jm, km
use constants_mod, only: Cp , EL , EpsV , ES0 , RVap , Grav ! 重力加速度
use dc_trace, only: SetDebug, BeginSub, EndSub, DbgMessage, DataDump
implicit none
!==== Input
!
real(DBKIND), intent(in) :: xyz_Press(im,jm,km) , xyr_Press(im,jm,km+1) , DelTimePhy ! (in) 2Δt
!==== Output
!
real(DBKIND), intent(out) :: xy_LscRain(im,jm) , xyz_DLscTempDt(im,jm,km), xyz_DLscQvapDt(im,jm,km) ! (out) 比湿変化率
!==== In/Out
!
real(DBKIND), intent(inout) :: xyz_Temp(im,jm,km) , xyz_Qvap(im,jm,km) ! (inout) 比湿
!----- 作業用内部変数 -----
character(STRING), parameter:: subname = "physics_lscond"
real(DBKIND) :: xyz_Qvap_b(im,jm,km) , xyz_Temp_b(im,jm,km) , QvapSat , DQvapSatDTemp , DelTemp, DelQvap
real(DBKIND), parameter :: CrtlRH = 1.0d0 ! 臨界相対湿度
integer(INTKIND), parameter :: IterationMax = 3 ! イテレーション回数
! do ループ用作業変数 (東西 i*、南北 j*、鉛直 k*、波数 l*用)
integer(INTKIND) :: i, j, k
integer(INTKIND) :: Iteration
continue
!----------------------------------------------------------------
! 開始処理
!----------------------------------------------------------------
call BeginSub(subname)
!----------------------------------------------------------------
! 大規模凝結
!----------------------------------------------------------------
!----- 調節前 Qvap の保存 -----
xyz_Qvap_b = xyz_Qvap
xyz_Temp_b = xyz_Temp
!----- 変数初期化 -----
xy_LscRain = 0.0d0
xyz_DLscTempDt = 0.0d0
xyz_DLscQvapDt = 0.0d0
!----- 調節 ------
do k = km, 1, -1
do i = 1, im
do j = 1, jm
! クラウジウスクラペイロンの式より飽和比湿計算
QvapSat = EpsV * ES0 * EXP( EL / RVap * ( 1./273. - 1./xyz_Temp(i,j,k) ) ) / xyz_Press(i,j,k)
! 飽和していたら, 温度と比湿の変化の収束計算
if ( ( xyz_Qvap(i,j,k) / QvapSat ) .GE. CrtlRH ) then
do Iteration = 1, IterationMax
! 飽和比湿計算
QvapSat = EpsV * ES0 * EXP( EL / RVap * (1./273. - 1./xyz_Temp(i,j,k))) / xyz_Press(i,j,k)
DQvapSatDTemp = EL * QvapSat / ( RVap * xyz_Temp(i,j,k) * xyz_Temp(i,j,k) )
! 温度と比湿の変化分を Newton 法で求める
DelTemp = EL / Cp * ( xyz_Qvap(i,j,k) - QvapSat ) / ( 1. + EL / Cp * DQvapSatDTemp )
DelQvap = DQvapSatDTemp * DelTemp
! 調節
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
!----- 比湿変化率, 温度変化率, 降水量の算出 -----
do k = km, 1, -1
do i = 1, im
do j = 1, jm
! 比湿変化率
xyz_DLscQvapDt(i,j,k) = xyz_DLscQvapDt(i,j,k) + ( xyz_Qvap(i,j,k) - xyz_Qvap_b(i,j,k) ) / DelTimePhy
! 温度変化率
xyz_DLscTempDt(i,j,k) = xyz_DLscTempDt(i,j,k) + ( xyz_Temp(i,j,k) - xyz_Temp_b(i,j,k) ) / DelTimePhy
! 降水量
xy_LscRain(i,j) = xy_LscRain(i,j) + ( xyz_Temp(i,j,k) - xyz_Temp_b(i,j,k) ) * CP / DelTimePhy * ( xyr_Press(i,j,k) - xyr_Press(i,j,k+1) ) /Grav
end do
end do
end do
!----------------------------------------------------------------
! 終了処理
!----------------------------------------------------------------
call EndSub(subname)
end subroutine physics_lscond