| Class | lscond_LL91 |
| In: |
lscond/lscond_LL91.f90
|
Note that Japanese and English are described in parallel.
大規模凝結過程によって温度と比湿を調節します.
Adjust temperature and specific humidity by a large scale condensation process.
Le Treut, H., and Z.-X. Li, Sensitivity of an atmospheric general circulation model to prescribed SST changes: feedback effects associated with the simulation of cloud optical properties, Clim. Dyn., 5, 175-187, 1991. Manabe, S., J. Smagorinsky, R. F. Strickler, Simulated climatology of a general circulation model with a hydrologic cycle, Mon. Wea. Rev., 93, 769-798, 1965.
| LScaleCondLL91 : | 温度と比湿の調節 |
| LScaleCondLL91Init : | 初期化 |
| —————— : | ———— |
| LScaleCondLL91 : | Adjust temperature and specific humidity |
| LScaleCondLL91Init : | Initialization |
| Subroutine : | |||
| xyz_Temp(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(inout)
| ||
| xyz_QVap(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(inout)
| ||
| xyz_DTempDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(inout)
| ||
| xyz_DQVapDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(inout)
| ||
| xyz_Press(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
| xyr_Press(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
| xyz_DQH2OLiqDt(0:imax-1,1:jmax,1:kmax) : | real(DP), intent(out) |
大規模凝結スキームにより, 温度と比湿を調節します.
Adjust temperature and specific humidity by large scale condensation scheme.
subroutine LScaleCondLL91( xyz_Temp, xyz_QVap, 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):: 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) :: xyz_DQH2OLiqDt(0:imax-1,1:jmax,1:kmax)
! 作業変数
! Work variables
!
real(DP) :: xyz_CldFrac (0:imax-1, 1:jmax, 1:kmax)
!
! Cloud fraction
real(DP) :: xyz_QVapHalRan(0:imax-1, 1:jmax, 1:kmax)
!
! Half range of specific humidity
real(DP) :: xyz_QVapCler (0:imax-1, 1:jmax, 1:kmax)
!
! Specific humidity in clear region
real(DP) :: xyz_QVapCldy (0:imax-1, 1:jmax, 1:kmax)
!
! Specific humidity in cloudy region
real(DP) :: xyz_TempCler (0:imax-1, 1:jmax, 1:kmax)
!
! Temperature in clear region
real(DP) :: xyz_TempCldy (0:imax-1, 1:jmax, 1:kmax)
!
! Temperature in cloudy region
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)
logical :: xyz_FlagSaturated(0:imax-1, 1:jmax, 1:kmax)
! Variables for debug
!!$ real(DP) :: TempBefAdj
!!$ real(DP) :: QVapBefAdj
! 実行文 ; Executable statement
!
! 初期化確認
! Initialization check
!
if ( .not. lscond_LL91_inited ) then
call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
end if
!!$ if ( .not. FlagUse ) return
! 計算時間計測開始
! Start measurement of computation time
!
call TimesetClockStart( module_name )
! 調節前 "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 )
! Preparation for temperature in clear and cloudy regions
!
xyz_TempCler = xyz_Temp
xyz_TempCldy = xyz_Temp
! Calculation of cloud fraction, specific humidities in clear and cloudy regions
!
xyz_QVapHalRan = Gamma * xyz_QVap
xyz_CldFrac = ( xyz_QVap + xyz_QVapHalRan - xyz_QVapSat ) / ( 2.0_DP * xyz_QVapHalRan + 1.0e-100_DP )
xyz_QVapCler = ( xyz_QVapSat + xyz_QVap - xyz_QVapHalRan ) / 2.0_DP
xyz_QVapCldy = ( xyz_QVapSat + xyz_QVap + xyz_QVapHalRan ) / 2.0_DP
do k = 1, kmax
do j = 1, jmax
do i = 0, imax-1
if ( xyz_CldFrac(i,j,k) >= 1.0_DP ) then
xyz_QVapCler(i,j,k) = 0.0_DP
xyz_QVapCldy(i,j,k) = xyz_QVap(i,j,k)
else if ( xyz_CldFrac(i,j,k) <= 0.0_DP ) then
xyz_QVapCler(i,j,k) = xyz_QVap(i,j,k)
xyz_QVapCldy(i,j,k) = 0.0_DP
end if
end do
end do
end do
do k = 1, kmax
do j = 1, jmax
do i = 0, imax-1
if ( xyz_CldFrac(i,j,k) > 1.0_DP ) then
xyz_CldFrac(i,j,k) = 1.0_DP
else if ( xyz_CldFrac(i,j,k) < 0.0_DP ) then
xyz_CldFrac(i,j,k) = 0.0_DP
end if
end do
end do
end do
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
if ( xyz_CldFrac(i,j,k) > 0.0_DP ) 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_TempCldy, xyz_Press )
xyz_DQVapSatDTemp = xyz_CalcDQVapSatDTemp( xyz_TempCldy, 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_QVapCldy(i,j,k) - xyz_QVapSat(i,j,k) ) / ( CpDry + LatentHeat * xyz_DQVapSatDTemp(i,j,k) )
!=========
! check routine
!---------
!!$ TempBefAdj = xyz_TempCldy(i,j,k)
!!$ QVapBefAdj = xyz_QVapCldy(i,j,k)
!=========
! 温度と比湿の調節
! Adjust temperature and specific humidity
!
xyz_TempCldy(i,j,k) = xyz_TempCldy(i,j,k) + DelTemp
xyz_QVapCldy(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', xyz_CldFrac(i,j,k)
!!$ write( 6, * ) &
!!$ & ( CpDry*TempBefAdj + LatentHeat*QVapBefAdj &
!!$ & - ( CpDry*xyz_TempCldy(i,j,k) + LatentHeat*xyz_QVapCldy(i,j,k) ) ) &
!!$ & /( CpDry*TempBefAdj + LatentHeat*QVapBefAdj ), &
!!$ & CpDry*TempBefAdj + LatentHeat*QVapBefAdj, &
!!$ & ( CpDry*xyz_TempCldy(i,j,k) + LatentHeat*xyz_QVapCldy(i,j,k) )
!=========
end if
end do
end do
end do
end do
! Calculation of temperature and specific humidity averaged over a grid
!
xyz_Temp = xyz_CldFrac * xyz_TempCldy + ( 1.0_DP - xyz_CldFrac ) * xyz_TempCler
xyz_QVap = xyz_CldFrac * xyz_QVapCldy + ( 1.0_DP - xyz_CldFrac ) * xyz_QVapCler
! 比湿変化率, 温度変化率, 降水量の算出
! 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 )
xyz_DTempDt = xyz_DTempDt + xyz_DTempDtLsc
xyz_DQVapDt = xyz_DQVapDt + xyz_DQVapDtLsc
!!$ 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)
!!$ xyz_DQH2OLiqDt = - xyz_DQVapDtLsc
!!$ xy_RainLsc = 0.0d0
!!$ do k = kmax, 1, -1
!!$ xy_RainLsc = xy_RainLsc &
!!$ & + xyz_DQH2OLiqDt(:,:,k) * ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav
!!$ end do
!!$ i = 0
!!$ j = jmax/2+1
!!$ write( 6, * ) xy_RainLsc(i,j)
!!$ write( 6, * ) '---'
!!$ xy_Rain = xy_Rain + xy_RainLsc
xyz_DQH2OLiqDt = - xyz_DQVapDtLsc
! calculation for output
xy_RainLsc = 0.0_DP
do k = kmax, 1, -1
xy_RainLsc = xy_RainLsc + xyz_DQH2OLiqDt(:,:,k) * ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav
end do
! ヒストリデータ出力
! History data output
!
call HistoryAutoPut( TimeN, 'RainLsc', xy_RainLsc * LatentHeat )
call HistoryAutoPut( TimeN, 'DTempDtLsc', xyz_DTempDtLsc )
call HistoryAutoPut( TimeN, 'DQVapDtLsc', xyz_DQVapDtLsc )
! 計算時間計測一時停止
! Pause measurement of computation time
!
call TimesetClockStop( module_name )
end subroutine LScaleCondLL91
| Subroutine : | |||
| xyz_Press(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
| xyr_Press(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
| xyz_Temp(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(inout)
| ||
| xyz_QVap(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(inout)
| ||
| xyz_QLiq(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(inout)
| ||
| xyz_QSol(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(inout)
| ||
| xyz_DTempDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
| xyz_DQVapDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
| ||
| xyz_DQLiqDt(0:imax-1,1:jmax,1:kmax) : | real(DP), intent(out) | ||
| xyz_DQSolDt(0:imax-1,1:jmax,1:kmax) : | real(DP), intent(out) |
大規模凝結スキームにより, 温度と比湿を調節します.
Adjust temperature and specific humidity by large scale condensation scheme.
subroutine LScaleCondLL91Ice( xyz_Press, xyr_Press, xyz_Temp, xyz_QVap, xyz_QLiq, xyz_QSol, xyz_DTempDt, xyz_DQVapDt, xyz_DQLiqDt, xyz_DQSolDt )
!
! 大規模凝結スキームにより, 温度と比湿を調節します.
!
! Adjust temperature and specific humidity by
! large scale condensation scheme.
!
! モジュール引用 ; USE statements
!
! 物理定数設定
! Physical constants settings
!
use constants, only: Grav, CpDry, LatentHeat, LatentHeatFusion
! 時刻管理
! 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, SaturateWatFraction
! 宣言文 ; Declaration statements
!
implicit none
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(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):: xyz_QLiq (0:imax-1, 1:jmax, 1:kmax)
! $ q_l $ . Specific liquid water content
real(DP), intent(inout):: xyz_QSol (0:imax-1, 1:jmax, 1:kmax)
! $ q_i $ . Specific ice water content
real(DP), intent(out):: xyz_DTempDt (0:imax-1, 1:jmax, 1:kmax)
! 温度変化率.
! Temperature tendency
real(DP), intent(out):: xyz_DQVapDt (0:imax-1, 1:jmax, 1:kmax)
! 比湿変化率.
! Specific humidity tendency
real(DP), intent(out) :: xyz_DQLiqDt(0:imax-1,1:jmax,1:kmax)
real(DP), intent(out) :: xyz_DQSolDt(0:imax-1,1:jmax,1:kmax)
! 作業変数
! Work variables
!
real(DP) :: xyz_CldFrac (0:imax-1, 1:jmax, 1:kmax)
!
! Cloud fraction
real(DP) :: QCld
real(DP) :: WatFrac
real(DP) :: IceFrac
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_DQLiqDtLsc (0:imax-1, 1:jmax, 1:kmax)
real(DP) :: xyz_DQSolDtLsc (0:imax-1, 1:jmax, 1:kmax)
real(DP) :: xyz_TempTentative(0:imax-1, 1:jmax, 1:kmax)
real(DP) :: xyz_TempB (0:imax-1, 1:jmax, 1:kmax)
! 調節前の温度.
! Temperature before adjust.
real(DP) :: xyz_QVapB (0:imax-1, 1:jmax, 1:kmax)
! 調節前の比湿.
! Specific humidity before adjust.
real(DP) :: xyz_QLiqB (0:imax-1, 1:jmax, 1:kmax)
real(DP) :: xyz_QSolB (0:imax-1, 1:jmax, 1:kmax)
!
real(DP) :: xyz_QTot (0:imax-1, 1:jmax, 1:kmax)
real(DP) :: xyz_DelQ (0:imax-1, 1:jmax, 1:kmax)
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)
logical :: xyz_FlagSaturated(0:imax-1, 1:jmax, 1:kmax)
! Variables for debug
!!$ real(DP) :: TempBefAdj
!!$ real(DP) :: QVapBefAdj
! 実行文 ; Executable statement
!
! 初期化確認
! Initialization check
!
if ( .not. lscond_LL91_inited ) then
call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
end if
!!$ if ( .not. FlagUse ) return
! 計算時間計測開始
! Start measurement of computation time
!
call TimesetClockStart( module_name )
! 調節前 "QVap", "Temp" の保存
! Store "QVap", "Temp" before adjustment
!
xyz_TempB = xyz_Temp
xyz_QVapB = xyz_QVap
xyz_QLiqB = xyz_QLiq
xyz_QSolB = xyz_QSol
xyz_QTot = xyz_QVap + xyz_QLiq + xyz_QSol
xyz_DelQ = Gamma * xyz_QTot
! All cloud water and ice are evaporated temporarily
! After this temporal evaporation, adjustment will be done below.
!
xyz_TempTentative = xyz_Temp - LatentHeat * xyz_QLiq - ( LatentHeat + LatentHeatFusion ) * xyz_QSol
! 調節
! Adjustment
!
! 飽和比湿計算
! Calculate saturation specific humidity
!
xyz_QVapSat = xyz_CalcQVapSat( xyz_TempTentative, 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
if ( ( ( xyz_QTot(i,j,k) + xyz_DelQ(i,j,k) ) / xyz_QVapSat(i,j,k) ) >= 1.0_DP ) 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
! Liquid water and ice fractions
call SaturateWatFraction( xyz_Temp(i,j,k), WatFrac )
IceFrac = 1.0_DP - WatFrac
! 温度の変化分をニュートン法で求める
! Calculate variation of temperature
!
DelTemp = ( LatentHeat * ( xyz_QVap(i,j,k) - ( xyz_QVapSat(i,j,k) + xyz_QTot(i,j,k) - xyz_DelQ(i,j,k) ) / 2.0_DP ) - LatentHeatFusion * ( xyz_QSol(i,j,k) - IceFrac * ( xyz_QVapSat(i,j,k) + xyz_QTot(i,j,k) + xyz_DelQ(i,j,k) ) / 2.0_DP ) ) / ( CpDry + ( LatentHeat - LatentHeatFusion * IceFrac ) * xyz_DQVapSatDTemp(i,j,k) / 2.0_DP )
! 温度と比湿の調節
! Adjust temperature and specific humidity
!
xyz_QVapSat(i,j,k) = xyz_QVapSat(i,j,k) + xyz_DQVapSatDTemp(i,j,k) * DelTemp
xyz_Temp(i,j,k) = xyz_Temp(i,j,k) + DelTemp
QCld = ( xyz_QVapSat(i,j,k) + xyz_QTot(i,j,k) + xyz_DelQ(i,j,k) ) / 2.0_DP
xyz_QVap(i,j,k) = ( xyz_QVapSat(i,j,k) + xyz_QTot(i,j,k) - xyz_DelQ(i,j,k) ) / 2.0_DP
xyz_QLiq(i,j,k) = WatFrac * QCld
xyz_QSol(i,j,k) = IceFrac * QCld
end if
end do
end do
end do
xyz_CldFrac = ( xyz_QTot + xyz_DelQ - xyz_QVapSat )
xyz_CldFrac = max( min( xyz_CldFrac, 1.0_DP ), 0.0_DP )
end do
! 比湿変化率, 温度変化率, 降水量の算出
! Calculate specific humidity tendency, temperature tendency,
! precipitation
!
xyz_DTempDt = ( xyz_Temp - xyz_TempB ) / ( 2.0_DP * DelTime )
xyz_DQVapDt = ( xyz_QVap - xyz_QVapB ) / ( 2.0_DP * DelTime )
xyz_DQLiqDt = ( xyz_QLiq - xyz_QLiqB ) / ( 2.0_DP * DelTime )
xyz_DQSolDt = ( xyz_QSol - xyz_QSolB ) / ( 2.0_DP * DelTime )
! calculation for output
xy_RainLsc = 0.0_DP
do k = kmax, 1, -1
xy_RainLsc = xy_RainLsc + ( xyz_DQLiqDt(:,:,k) + xyz_DQSolDt(:,:,k) ) * ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav
end do
! ヒストリデータ出力
! History data output
!
call HistoryAutoPut( TimeN, 'RainLsc', xy_RainLsc * LatentHeat )
call HistoryAutoPut( TimeN, 'DTempDtLsc', xyz_DTempDtLsc )
call HistoryAutoPut( TimeN, 'DQVapDtLsc', xyz_DQVapDtLsc )
! 計算時間計測一時停止
! Pause measurement of computation time
!
call TimesetClockStop( module_name )
end subroutine LScaleCondLL91Ice
| Subroutine : |
lscond モジュールの初期化を行います. NAMELIST#lscond_nml の読み込みはこの手続きで行われます.
"lscond" module is initialized. "NAMELIST#lscond_nml" is loaded in this procedure.
This procedure input/output NAMELIST#lscond_nml .
subroutine LScaleCondLL91Init
!
! lscond モジュールの初期化を行います.
! NAMELIST#lscond_nml の読み込みはこの手続きで行われます.
!
! "lscond" module is initialized.
! "NAMELIST#lscond_nml" is loaded in this procedure.
!
! モジュール引用 ; USE statements
!
! NAMELIST ファイル入力に関するユーティリティ
! Utilities for NAMELIST file input
!
use namelist_util, only: namelist_filename, NmlutilMsg, NmlutilAryValid
! ファイル入出力補助
! File I/O support
!
use dc_iounit, only: FileOpen
! 種別型パラメタ
! Kind type parameter
!
use dc_types, only: STDOUT ! 標準出力の装置番号. Unit number of standard output
! 文字列操作
! Character handling
!
use dc_string, only: StoA
! ヒストリデータ出力
! History data output
!
use gtool_historyauto, only: HistoryAutoAddVariable
! 宣言文 ; Declaration statements
!
implicit none
integer:: unit_nml ! NAMELIST ファイルオープン用装置番号.
! Unit number for NAMELIST file open
integer:: iostat_nml ! NAMELIST 読み込み時の IOSTAT.
! IOSTAT of NAMELIST read
! NAMELIST 変数群
! NAMELIST group name
!
namelist /lscond_nml/ Gamma, ItrtMax
!
! デフォルト値については初期化手続 "lscond#LSCondInit"
! のソースコードを参照のこと.
!
! Refer to source codes in the initialization procedure
! "lscond#LSCondInit" for the default values.
!
! 実行文 ; Executable statement
!
if ( lscond_LL91_inited ) return
! デフォルト値の設定
! Default values settings
!
!!$ CrtlRH = 1.0_DP
Gamma = 0.2_DP
! This value is from Le Treut and Li (1991).
ItrtMax = 3
!!$ FlagUse = .true.
! NAMELIST の読み込み
! NAMELIST is input
!
if ( trim(namelist_filename) /= '' ) then
call FileOpen( unit_nml, namelist_filename, mode = 'r' ) ! (in)
rewind( unit_nml )
read( unit_nml, nml = lscond_nml, iostat = iostat_nml ) ! (out)
close( unit_nml )
call NmlutilMsg( iostat_nml, module_name ) ! (in)
if ( iostat_nml == 0 ) write( STDOUT, nml = lscond_nml )
end if
! ヒストリデータ出力のためのへの変数登録
! Register of variables for history data output
!
call HistoryAutoAddVariable( 'RainLsc', (/ 'lon ', 'lat ', 'time' /), 'precipitation by large scale condensation', 'W m-2' )
call HistoryAutoAddVariable( 'DTempDtLsc', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'large-scale condensation heating', 'K s-1' )
call HistoryAutoAddVariable( 'DQVapDtLsc', (/ 'lon ', 'lat ', 'sig ', 'time' /), 'large-scale condensation moistening', 'kg kg-1 s-1' )
! 印字 ; Print
!
call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
!!$ call MessageNotify( 'M', module_name, ' FlagUse = %b', l = (/ FlagUse /) )
!!$ call MessageNotify( 'M', module_name, ' CrtlRH = %f', d = (/ CrtlRH /) )
call MessageNotify( 'M', module_name, ' Gamma = %f', d = (/ Gamma /) )
call MessageNotify( 'M', module_name, ' ItrtMax = %d', i = (/ ItrtMax /) )
call MessageNotify( 'M', module_name, '' )
call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )
lscond_LL91_inited = .true.
end subroutine LScaleCondLL91Init
| Variable : | |||||
| Gamma : | real(DP), save
|
| Variable : | |||||
| ItrtMax : | integer, save
|
| Variable : | |||
| lscond_LL91_inited = .false. : | logical, save
|
| Constant : | |||
| version = ’$Name: $’ // ’$Id: lscond_LL91.f90,v 1.4 2015/01/29 12:02:56 yot Exp $’ : | character(*), parameter
|