時間変化率の計算を行います.
Calculate tendencies.
subroutine PhyImplTendency( xyr_UFlux, xyr_VFlux, xyr_TempFlux, xyr_QVapFlux, xyr_RadSFlux, xyr_RadLFlux, xy_GroundTempFlux, xy_SurfTemp, xy_SurfHumidCoef, xy_SurfCond, xy_SurfHeatCapacity, xyra_DelRadLFlux, xyr_Press, xyz_Exner, xyr_Exner, xyr_VelTransCoef, xyr_TempTransCoef, xyr_QVapTransCoef, xy_SurfVelTransCoef, xy_SurfTempTransCoef, xy_SurfQVapTransCoef, xyz_DUDt, xyz_DVDt, xyz_DTempDt, xyz_DQVapDt, xy_DSurfTempDt )
!
! 時間変化率の計算を行います.
!
! Calculate tendencies.
!
! モジュール引用 ; USE statements
!
! 物理定数設定
! Physical constants settings
!
use constants, only: Grav, CpDry, LatentHeat, GasRDry
! $ R $ [J kg-1 K-1].
! 乾燥大気の気体定数.
! Gas constant of air
! 時刻管理
! Time control
!
use timeset, only: DelTime, TimeN, TimesetClockStart, TimesetClockStop
! 時刻管理
! Time control
!
use timeset, only: TimeN ! ステップ $ t $ の時刻. Time of step $ t $.
! ヒストリデータ出力
! History data output
!
use gtool_historyauto, only: HistoryAutoPut
! 飽和比湿計算
! Evaluate saturation specific humidity
!
#ifdef LIB_SATURATE_NHA1992
use saturate_nha1992, only: CalcQVapSat, CalcDQVapSatDTemp
#elif LIB_SATURATE_T1930
use saturate_t1930, only: CalcQVapSat, CalcDQVapSatDTemp
#else
use saturate_t1930, only: CalcQVapSat, CalcDQVapSatDTemp
#endif
! 宣言文 ; Declaration statements
!
implicit none
real(DP), intent(in):: xyr_UFlux (0:imax-1, 1:jmax, 0:kmax)
! 東西風速フラックス.
! Eastward wind flux
real(DP), intent(in):: xyr_VFlux (0:imax-1, 1:jmax, 0:kmax)
! 南北風速フラックス.
! Northward wind flux
real(DP), intent(in):: xyr_TempFlux (0:imax-1, 1:jmax, 0:kmax)
! 温度フラックス.
! Temperature flux
real(DP), intent(in):: xyr_QVapFlux (0:imax-1, 1:jmax, 0:kmax)
! 比湿フラックス.
! Specific humidity flux
real(DP), intent(in):: xyr_RadSFlux (0:imax-1, 1:jmax, 0:kmax)
! 短波 (日射) フラックス.
! Shortwave (insolation) flux
real(DP), intent(in):: xyr_RadLFlux (0:imax-1, 1:jmax, 0:kmax)
! 長波フラックス.
! Longwave flux
real(DP), intent(in):: xy_GroundTempFlux (0:imax-1, 1:jmax)
! 地中熱フラックス.
! Ground temperature flux
real(DP), intent(in):: xy_SurfTemp (0:imax-1, 1:jmax)
! 地表面温度.
! Surface temperature
real(DP), intent(in):: xy_SurfHumidCoef (0:imax-1, 1:jmax)
! 地表湿潤度.
! Surface humidity coefficient
integer, intent(in):: xy_SurfCond (0:imax-1, 1:jmax)
! 地表状態.
! Surface condition
real(DP), intent(in):: xy_SurfHeatCapacity (0:imax-1, 1:jmax)
! 地表熱容量.
! Surface heat capacity
real(DP), intent(in):: xyra_DelRadLFlux (0:imax-1, 1:jmax, 0:kmax, 0:1)
! 長波地表温度変化.
! Surface temperature tendency with longwave
real(DP), intent(in):: xyr_Press (0:imax-1, 1:jmax, 0:kmax)
! $ \hat{p} $ . 気圧 (半整数レベル).
! Air pressure (half level)
real(DP), intent(in):: xyz_Exner (0:imax-1, 1:jmax, 1:kmax)
! Exner 関数 (整数レベル).
! Exner function (full level)
real(DP), intent(in):: xyr_Exner (0:imax-1, 1:jmax, 0:kmax)
! Exner 関数 (半整数レベル).
! Exner function (half level)
real(DP), intent(in):: xyr_VelTransCoef (0:imax-1, 1:jmax, 0:kmax)
! 輸送係数:運動量.
! Transfer coefficient: velocity
real(DP), intent(in):: xyr_TempTransCoef (0:imax-1, 1:jmax, 0:kmax)
! 輸送係数:温度.
! Transfer coefficient: temperature
real(DP), intent(in):: xyr_QVapTransCoef (0:imax-1, 1:jmax, 0:kmax)
! 輸送係数:比湿.
! Transfer coefficient: specific humidity
real(DP), intent(in):: xy_SurfVelTransCoef (0:imax-1, 1:jmax)
! 輸送係数:運動量.
! Diffusion coefficient: velocity
real(DP), intent(in):: xy_SurfTempTransCoef (0:imax-1, 1:jmax)
! 輸送係数:温度.
! Transfer coefficient: temperature
real(DP), intent(in):: xy_SurfQVapTransCoef (0:imax-1, 1:jmax)
! 輸送係数:比湿.
! Transfer coefficient: specific humidity
real(DP), intent(out):: xyz_DUDt (0:imax-1, 1:jmax, 1:kmax)
! $ \DP{u}{t} $ . 東西風速変化.
! Eastward wind tendency
real(DP), intent(out):: xyz_DVDt (0:imax-1, 1:jmax, 1:kmax)
! $ \DP{v}{t} $ . 南北風速変化.
! Northward wind tendency
real(DP), intent(out):: xyz_DTempDt (0:imax-1, 1:jmax, 1:kmax)
! $ \DP{T}{t} $ . 温度変化.
! Temperature tendency
real(DP), intent(out):: xyz_DQVapDt (0:imax-1, 1:jmax, 1:kmax)
! $ \DP{q}{t} $ . 比湿変化.
! Temperature tendency
real(DP), intent(out):: xy_DSurfTempDt (0:imax-1, 1:jmax)
! 地表面温度変化率.
! Surface temperature tendency
! 作業変数
! Work variables
!
real(DP):: xyza_UVMtx (0:imax-1, 1:jmax, 1:kmax, -1:1)
! 速度陰解行列.
! Implicit matrix about velocity
real(DP):: xyra_TempMtx (0:imax-1, 1:jmax, 0:kmax, -1:1)
! 温度陰解行列.
! Implicit matrix about temperature
real(DP):: xyza_QVapMtx (0:imax-1, 1:jmax, 1:kmax, -1:1)
! 比湿陰解行列.
! Implicit matrix about specific humidity
real(DP):: xy_SurfUVMtx (0:imax-1, 1:jmax)
! 速度陰解行列: 地表.
! Implicit matrix about velocity: surface
real(DP):: xyaa_SurfTempMtx (0:imax-1, 1:jmax, 0:1, -1:1)
! 温度陰解行列: 地表.
! Implicit matrix about temperature: surface
real(DP):: xyaa_SurfQVapMtx (0:imax-1, 1:jmax, 0:1, -1:1)
! 比湿陰解行列: 地表.
! Implicit matrix about specific humidity: surface
real(DP):: xya_SurfRadLMtx (0:imax-1, 1:jmax, -1:1)
! $ T $ . 陰解行列: 地表.
! implicit matrix: surface
real(DP):: xyz_DelTempQVap (0:imax-1, 1:jmax, -kmax:kmax)
! $ T q $ の時間変化.
! Tendency of $ T q $
real(DP):: xyza_TempQVapLUMtx (0:imax-1, 1:jmax, -kmax:kmax, -1:1)
! LU 行列.
! LU matrix
real(DP):: xyza_UVLUMtx (0:imax-1, 1:jmax, 1:kmax,-1:1)
! LU 行列.
! LU matrix
real(DP):: xy_SurfExner (0:imax-1, 1:jmax)
! Exner 関数.
! Exner function
real(DP):: xy_SurfQVapSat (0:imax-1, 1:jmax)
! 地表飽和比湿.
! Saturated specific humidity on surface
real(DP):: xy_SurfDQVapSatDTemp (0:imax-1, 1:jmax)
! 地表飽和比湿変化.
! Saturated specific humidity tendency on surface
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:: l ! 行列用 DO ループ用作業変数
! Work variables for DO loop of matrices
! 実行文 ; Executable statement
!
! 計算時間計測開始
! Start measurement of computation time
!
call TimesetClockStart( module_name )
! 初期化
! Initialization
!
if ( .not. phy_implicit_inited ) call PhyImplInit
! 長波陰解用行列の算出
! Calculate long wave implicit matrix
!
xya_SurfRadLMtx(:,:,0) = xyra_DelRadLFlux(:,:,0,0)
xya_SurfRadLMtx(:,:,1) = xyra_DelRadLFlux(:,:,0,1)
xya_SurfRadLMtx(:,:,-1) = 0.
! 陰解法のための行列作成
! Create matrices for implicit scheme
!
do k = 1, kmax
xyza_UVMtx(:,:,k,0) = ( xyr_Press(:,:,k-1) - xyr_Press(:,:,k) ) / Grav / ( 2. * DelTime )
xyra_TempMtx(:,:,k,0) = xyza_UVMtx(:,:,k,0) * CpDry
xyza_QVapMtx(:,:,k,0) = xyza_UVMtx(:,:,k,0) * CpDry
end do
xyra_TempMtx(:,:,0,0) = xy_SurfHeatCapacity(:,:) / ( 2. * DelTime )
! 鉛直拡散スキームの輸送係数から陰解行列の計算 (速度)
! Calculate implicit matrices from transfer coefficient of vertical diffusion scheme (velocity)
!
xyza_UVMtx (:,:,1,-1) = 0.
do k = 2, kmax
xyza_UVMtx(:,:,k, 0) = xyza_UVMtx(:,:,k,0) + xyr_VelTransCoef(:,:,k-1)
xyza_UVMtx(:,:,k,-1) = - xyr_VelTransCoef(:,:,k-1)
end do
do k = 1, kmax-1
xyza_UVMtx(:,:,k, 0) = xyza_UVMtx(:,:,k,0) + xyr_VelTransCoef(:,:,k)
xyza_UVMtx(:,:,k, 1) = - xyr_VelTransCoef(:,:,k)
end do
xyza_UVMtx (:,:,kmax, 1) = 0.
! 鉛直拡散スキームの輸送係数から陰解行列の計算 (温度)
! Calculate implicit matrices from transfer coefficient of vertical diffusion scheme (temperature)
!
xyra_TempMtx (:,:,0:1,-1) = 0.
do k = 2, kmax
xyra_TempMtx(:,:, k, 0) = xyra_TempMtx(:,:,k,0) + CpDry * xyr_TempTransCoef(:,:,k-1) * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k)
xyra_TempMtx(:,:, k,-1) = - CpDry * xyr_TempTransCoef(:,:,k-1) * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k-1)
end do
xyra_TempMtx (:,:, 0, 1) = 0.
do k = 1, kmax-1
xyra_TempMtx(:,:, k, 0) = xyra_TempMtx(:,:,k,0) + CpDry * xyr_TempTransCoef(:,:,k) * xyr_Exner(:,:,k) / xyz_Exner(:,:,k)
xyra_TempMtx(:,:, k, 1) = - CpDry * xyr_TempTransCoef(:,:,k) * xyr_Exner(:,:,k) / xyz_Exner(:,:,k+1)
end do
xyra_TempMtx(:,:, kmax, 1) = 0.
! 鉛直拡散スキームの輸送係数から陰解行列の計算 (比湿)
! Calculate implicit matrices from transfer coefficient of vertical diffusion scheme (specific humidity)
!
xyza_QVapMtx (:,:,1,-1) = 0.
do k = 2, kmax
xyza_QVapMtx(:,:,k, 0) = xyza_QVapMtx(:,:,k,0) + CpDry * xyr_QVapTransCoef(:,:,k-1)
xyza_QVapMtx(:,:,k,-1) = - CpDry * xyr_QVapTransCoef(:,:,k-1)
end do
do k = 1, kmax-1
xyza_QVapMtx(:,:,k, 0) = xyza_QVapMtx(:,:,k,0) + CpDry * xyr_QVapTransCoef(:,:,k)
xyza_QVapMtx(:,:,k, 1) = - CpDry * xyr_QVapTransCoef(:,:,k)
end do
xyza_QVapMtx(:,:,kmax,1) = 0.
! Exner 関数算出
! Calculate Exner functions
!
xy_SurfExner = xyz_Exner(:,:,1) / xyr_Exner(:,:,0)
! 飽和比湿の計算
! Calculate saturated specific humidity
!
do i = 0, imax-1
do j = 1, jmax
xy_SurfQVapSat(i,j) = CalcQVapSat( xy_SurfTemp(i,j), xyr_Press(i,j,0) )
end do
end do
do i = 0, imax-1
do j = 1, jmax
xy_SurfDQVapSatDTemp(i,j) = CalcDQVapSatDTemp( xy_SurfTemp(i,j), xy_SurfQVapSat(i,j) )
end do
end do
! 地表面過程の輸送係数から陰解行列の計算
! Calculate implicit matrices from transfer coefficient of surface process
!
xyaa_SurfTempMtx(:,:,0,-1) = 0.
xyaa_SurfTempMtx(:,:,1, 1) = 0.
xyaa_SurfQVapMtx(:,:,0,-1) = 0.
xyaa_SurfQVapMtx(:,:,1, 1) = 0.
xy_SurfUVMtx = xy_SurfVelTransCoef
xyaa_SurfTempMtx(:,:,1, 0) = CpDry * xy_SurfTempTransCoef / xy_SurfExner
xyaa_SurfTempMtx(:,:,0, 1) = - CpDry * xy_SurfTempTransCoef / xy_SurfExner
xyaa_SurfQVapMtx(:,:,1, 0) = CpDry * xy_SurfQVapTransCoef * xy_SurfHumidCoef
xyaa_SurfQVapMtx(:,:,0, 1) = - CpDry * xy_SurfQVapTransCoef * xy_SurfHumidCoef
do i = 0, imax-1
do j = 1, jmax
if ( xy_SurfCond(i,j) >= 1 ) then
xyaa_SurfTempMtx(i,j,1,-1) = - CpDry * xy_SurfTempTransCoef(i,j)
xyaa_SurfTempMtx(i,j,0, 0) = CpDry * xy_SurfTempTransCoef(i,j)
xyaa_SurfQVapMtx(i,j,1,-1) = - LatentHeat * xy_SurfQVapTransCoef(i,j) * xy_SurfHumidCoef(i,j) * xy_SurfDQVapSatDTemp(i,j)
xyaa_SurfQVapMtx(i,j,0, 0) = LatentHeat * xy_SurfQVapTransCoef(i,j) * xy_SurfHumidCoef(i,j) * xy_SurfDQVapSatDTemp(i,j)
else
xyaa_SurfTempMtx(i,j,1,-1) = 0.
xyaa_SurfTempMtx(i,j,0, 0) = 0.
xyaa_SurfQVapMtx(i,j,1,-1) = 0.
xyaa_SurfQVapMtx(i,j,0, 0) = 0.
end if
end do
end do
! 東西風速, 南北風速の計算
! Calculate eastward and northward wind
!
xyza_UVLUMtx = xyza_UVMtx
xyza_UVLUMtx(:,:,1,0) = xyza_UVLUMtx(:,:,1,0) + xy_SurfUVMtx
call PhyImplLUDecomp3( xyza_UVLUMtx, imax * jmax, kmax ) ! (in)
do k = 1, kmax
xyz_DUDt(:,:,k) = xyr_UFlux(:,:,k-1) - xyr_UFlux(:,:,k)
xyz_DVDt(:,:,k) = xyr_VFlux(:,:,k-1) - xyr_VFlux(:,:,k)
end do
call PhyImplLUSolve3( xyz_DUDt, xyza_UVLUMtx, 1, imax * jmax, kmax ) ! (in)
call PhyImplLUSolve3( xyz_DVDt, xyza_UVLUMtx, 1, imax * jmax, kmax ) ! (in)
! 温度と比湿の計算
! Calculate temperature and specific humidity
!
do l = -1, 1
do k = 1, kmax
xyza_TempQVapLUMtx(:,:,k,l) = xyra_TempMtx(:,:,k,l)
xyza_TempQVapLUMtx(:,:,-k,-l) = xyza_QVapMtx(:,:,k,l)
end do
xyza_TempQVapLUMtx(:,:, 1, l) = xyra_TempMtx(:,:,1,l) + xyaa_SurfTempMtx(:,:,1,l)
xyza_TempQVapLUMtx(:,:,-1,-l) = xyza_QVapMtx(:,:,1,l) + xyaa_SurfQVapMtx(:,:,1,l)
end do
do j = 1, jmax
do i = 0, imax-1
if ( xy_SurfCond(i,j) >= 1 ) then
xyza_TempQVapLUMtx(i,j,0, 0) = xyra_TempMtx(i,j,0,0) + xyaa_SurfTempMtx(i,j,0,0) + xyaa_SurfQVapMtx(i,j,0,0) + xya_SurfRadLMtx(i,j,0)
xyza_TempQVapLUMtx(i,j,0, 1) = xyaa_SurfTempMtx(i,j,0,1) + xya_SurfRadLMtx(i,j,1)
xyza_TempQVapLUMtx(i,j,0,-1) = xyaa_SurfQVapMtx(i,j,0,1)
else
xyza_TempQVapLUMtx(i,j,0, 0) = 1.
xyza_TempQVapLUMtx(i,j,0, 1) = 0.
xyza_TempQVapLUMtx(i,j,0,-1) = 0.
end if
end do
end do
call PhyImplLUDecomp3( xyza_TempQVapLUMtx, imax * jmax, (2 * kmax) + 1 ) ! (in)
do k = 1, kmax
xyz_DelTempQVap(:,:, k) = xyr_TempFlux(:,:,k-1) - xyr_TempFlux(:,:,k)
xyz_DelTempQVap(:,:,-k) = xyr_QVapFlux(:,:,k-1) - xyr_QVapFlux(:,:,k)
end do
do j = 1, jmax
do i = 0, imax-1
if ( xy_SurfCond(i,j) >= 1 ) then
xyz_DelTempQVap(i,j,0) = - xyr_RadSFlux(i,j,0) - xyr_RadLFlux(i,j,0) - xyr_TempFlux(i,j,0) - xyr_QVapFlux(i,j,0) + xy_GroundTempFlux(i,j)
else
xyz_DelTempQVap(i,j,0) = 0.
end if
end do
end do
call PhyImplLUSolve3( xyz_DelTempQVap, xyza_TempQVapLUMtx, 1, imax * jmax , (2 * kmax) + 1 ) ! (in)
! 時間変化率の計算
! Calculate tendency
!
do k = 1, kmax
xyz_DUDt(:,:,k) = xyz_DUDt(:,:,k) / ( 2. * DelTime )
xyz_DVDt(:,:,k) = xyz_DVDt(:,:,k) / ( 2. * DelTime )
xyz_DTempDt(:,:,k) = xyz_DelTempQVap(:,:, k) / ( 2. * DelTime )
xyz_DQVapDt(:,:,k) = xyz_DelTempQVap(:,:,-k) / ( 2. * DelTime ) / LatentHeat * CpDry
end do
do j = 1, jmax
do i = 0, imax-1
if ( xy_SurfCond(i,j) >= 1 ) then
xy_DSurfTempDt(i,j) = xyz_DelTempQVap(i,j,0) / ( 2. * DelTime )
else
xy_DSurfTempDt(i,j) = 0.
end if
end do
end do
! 計算時間計測一時停止
! Pause measurement of computation time
!
call TimesetClockStop( module_name )
end subroutine PhyImplTendency