時間変化率の計算を行います.
Calculate tendencies.
subroutine PhyImplTendency( xyr_MomFluxX, xyr_MomFluxY, xyr_HeatFlux, xyrf_QMixFlux, xyr_RadSFlux, xyr_RadLFlux, xy_GroundTempFlux, xy_SurfTemp, xy_SurfHumidCoef, xy_SurfCond, xy_SurfHeatCapacity, xyra_DelRadLFlux, xyr_Press, xyz_Exner, xyr_Exner, xyr_VirTemp, xyz_Height, xyr_VelDiffCoef, xyr_TempDiffCoef, xyr_QMixDiffCoef, xy_SurfVelTransCoef, xy_SurfTempTransCoef, xy_SurfQVapTransCoef, xyz_DUDt, xyz_DVDt, xyz_DTempDt, xyzf_DQMixDt, 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
! ヒストリデータ出力
! History data output
!
use gtool_historyauto, only: HistoryAutoPut
! 飽和比湿の算出
! Evaluate saturation specific humidity
!
use saturate, only: xy_CalcQVapSat, xy_CalcDQVapSatDTemp
! 陰解法による時間積分のためのルーチン
! Routines for time integration with implicit scheme
!
use phy_implicit_utils, only : PhyImplLUDecomp3, PhyImplLUSolve3
! 宣言文 ; Declaration statements
!
implicit none
real(DP), intent(in):: xyr_MomFluxX (0:imax-1, 1:jmax, 0:kmax)
! 東西方向運動量フラックス.
! Eastward momentum flux
real(DP), intent(in):: xyr_MomFluxY (0:imax-1, 1:jmax, 0:kmax)
! 南北方向運動量フラックス.
! Northward momentum flux
real(DP), intent(in):: xyr_HeatFlux (0:imax-1, 1:jmax, 0:kmax)
! 熱フラックス.
! Heat flux
real(DP), intent(in):: xyrf_QMixFlux (0:imax-1, 1:jmax, 0:kmax, 1:ncmax)
! 質量フラックス.
! Mass flux of constituents
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_VirTemp (0:imax-1, 1:jmax, 0:kmax)
! $ \hat{T}_v $ . 仮温度 (半整数レベル).
! Virtual temperature (half level)
real(DP), intent(in):: xyz_Height (0:imax-1, 1:jmax, 1:kmax)
! 高度 (整数レベル).
! Height (full level)
real(DP), intent(in):: xyr_VelDiffCoef (0:imax-1, 1:jmax, 0:kmax)
! 拡散係数:運動量.
! Diffusion coefficient: velocity
real(DP), intent(in):: xyr_TempDiffCoef (0:imax-1, 1:jmax, 0:kmax)
! 拡散係数:温度.
! Transfer coefficient: temperature
real(DP), intent(in):: xyr_QMixDiffCoef (0:imax-1, 1:jmax, 0:kmax)
! 拡散係数:比湿.
! Diffusion 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):: xyzf_DQMixDt(0:imax-1, 1:jmax, 1:kmax, 1:ncmax)
! $ \DP{q}{t} $ . 比湿変化.
! Temperature tendency
real(DP), intent(out):: xy_DSurfTempDt (0:imax-1, 1:jmax)
! 地表面温度変化率.
! Surface temperature tendency
! 作業変数
! Work variables
!
real(DP) :: xyr_VelTransCoef (0:imax-1, 1:jmax, 0:kmax)
! 輸送係数:運動量.
! Transfer coefficient: velocity
real(DP) :: xyr_TempTransCoef (0:imax-1, 1:jmax, 0:kmax)
! 輸送係数:温度.
! Transfer coefficient: temperature
real(DP) :: xyr_QMixTransCoef(0:imax-1, 1:jmax, 0:kmax)
! 輸送係数:質量.
! Transfer coefficient: mass of constituents
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):: xyza_QMixMtx (0:imax-1, 1:jmax, 1:kmax, -1:1)
! 質量混合比陰解行列.
! Implicit matrix about mass mixing ratio
real(DP):: xyaa_SurfMtx (0:imax-1, 1:jmax, 0:0, -1:1)
! 惑星表面エネルギー収支用陰解行列
! Implicit matrix for surface energy balance
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):: xyza_QMixLUMtx(0:imax-1, 1:jmax, 1:kmax,-1:1)
! LU 行列.
! LU matrix
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
integer:: n ! 組成方向に回る DO ループ用作業変数
! Work variables for DO loop in dimension of constituents
! 実行文 ; Executable statement
!
! 初期化確認
! Initialization check
!
if ( .not. phy_implicit_inited ) then
call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
end if
! 計算時間計測開始
! Start measurement of computation time
!
call TimesetClockStart( module_name )
! バケツモデルの扱いの確認
! Check use of bucket model
!
! Bucket model with this routine is not supported fully, although some
! parts of calculation is implemented.
!
if ( FlagBucketModel ) then
call MessageNotify( 'E', module_name, 'Bucket model cannot be used with this routine.' )
end if
! 雪の扱いの確認
! Check about treatment of snow
!
if ( FlagSnow ) then
call MessageNotify( 'E', module_name, 'Snow melt is not treated in this routine.' )
end if
! 輸送係数の計算
! Calculate transfer coefficient
!
xyr_VelTransCoef (:,:,0) = 0.0_DP
xyr_VelTransCoef (:,:,kmax) = 0.0_DP
xyr_TempTransCoef(:,:,0) = 0.0_DP
xyr_TempTransCoef(:,:,kmax) = 0.0_DP
xyr_QMixTransCoef(:,:,0) = 0.0_DP
xyr_QMixTransCoef(:,:,kmax) = 0.0_DP
do k = 1, kmax-1
xyr_VelTransCoef(:,:,k) = xyr_VelDiffCoef(:,:,k) * xyr_Press(:,:,k) / ( GasRDry * xyr_VirTemp(:,:,k) ) / ( xyz_Height(:,:,k+1) - xyz_Height(:,:,k) )
xyr_TempTransCoef(:,:,k) = xyr_TempDiffCoef(:,:,k) * xyr_Press(:,:,k) / ( GasRDry * xyr_VirTemp(:,:,k) ) / ( xyz_Height(:,:,k+1) - xyz_Height(:,:,k) )
xyr_QMixTransCoef(:,:,k) = xyr_QMixDiffCoef(:,:,k) * xyr_Press(:,:,k) / ( GasRDry * xyr_VirTemp(:,:,k) ) / ( xyz_Height(:,:,k+1) - xyz_Height(:,:,k) )
end do
! 陰解法のための行列作成
! Create matrices for implicit scheme
!
! 鉛直拡散スキームの輸送係数から陰解行列の計算 (速度)
! Calculate implicit matrices from transfer coefficient of vertical diffusion scheme (velocity)
!
k = 1
xyza_UVMtx (:,:,k,-1) = 0.0_DP
xyza_UVMtx (:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2.0_DP * DelTime ) + xy_SurfVelTransCoef(:,:) + xyr_VelTransCoef(:,:,k )
xyza_UVMtx (:,:,k, 1) = - xyr_VelTransCoef(:,:,k)
do k = 2, kmax-1
xyza_UVMtx (:,:,k,-1) = - xyr_VelTransCoef(:,:,k-1)
xyza_UVMtx (:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2.0_DP * DelTime ) + xyr_VelTransCoef(:,:,k-1) + xyr_VelTransCoef(:,:,k )
xyza_UVMtx (:,:,k, 1) = - xyr_VelTransCoef(:,:,k)
end do
k = kmax
xyza_UVMtx (:,:,k,-1) = - xyr_VelTransCoef(:,:,k-1)
xyza_UVMtx (:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2.0_DP * DelTime ) + xyr_VelTransCoef(:,:,k-1)
xyza_UVMtx (:,:,k, 1) = 0.0_DP
! 鉛直拡散スキームの輸送係数から陰解行列の計算 (温度)
! Calculate implicit matrices from transfer coefficient of vertical diffusion scheme (temperature)
!
k = 1
xyra_TempMtx(:,:,k,-1) = - CpDry * xy_SurfTempTransCoef(:,:)
xyra_TempMtx(:,:,k, 0) = - CpDry * ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2.0_DP * DelTime ) + CpDry * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k ) * xy_SurfTempTransCoef(:,:) + CpDry * xyr_Exner(:,:,k ) / xyz_Exner(:,:,k ) * xyr_TempTransCoef(:,:,k )
xyra_TempMtx(:,:,k, 1) = - CpDry * xyr_Exner(:,:,k ) / xyz_Exner(:,:,k+1) * xyr_TempTransCoef(:,:,k )
do k = 2, kmax-1
xyra_TempMtx(:,:,k,-1) = - CpDry * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k-1) * xyr_TempTransCoef(:,:,k-1)
xyra_TempMtx(:,:,k, 0) = - CpDry * ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2.0_DP * DelTime ) + CpDry * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k ) * xyr_TempTransCoef(:,:,k-1) + CpDry * xyr_Exner(:,:,k ) / xyz_Exner(:,:,k ) * xyr_TempTransCoef(:,:,k )
xyra_TempMtx(:,:,k, 1) = - CpDry * xyr_Exner(:,:,k ) / xyz_Exner(:,:,k+1) * xyr_TempTransCoef(:,:,k )
end do
k = kmax
xyra_TempMtx(:,:,k,-1) = - CpDry * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k-1) * xyr_TempTransCoef(:,:,k-1)
xyra_TempMtx(:,:,k, 0) = - CpDry * ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2.0_DP * DelTime ) + CpDry * xyr_Exner(:,:,k-1) / xyz_Exner(:,:,k ) * xyr_TempTransCoef(:,:,k-1)
xyra_TempMtx(:,:,k, 1) = 0.0_DP
! 鉛直拡散スキームの輸送係数から陰解行列の計算 (比湿)
! Calculate implicit matrices from transfer coefficient of vertical diffusion scheme (specific humidity)
!
! 飽和比湿の計算
! Calculate saturated specific humidity
!
xy_SurfQVapSat = xy_CalcQVapSat ( xy_SurfTemp, xyr_Press(:,:,0) )
xy_SurfDQVapSatDTemp = xy_CalcDQVapSatDTemp( xy_SurfTemp, xy_SurfQVapSat )
k = 1
if ( FlagBucketModel ) then
xyza_QVapMtx(:,:,k,-1) = 0.0d0
xyza_QVapMtx(:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2.0_DP * DelTime ) + xyr_QMixTransCoef(:,:,k )
else
xyza_QVapMtx(:,:,k,-1) = - xy_SurfHumidCoef(:,:) * xy_SurfQVapTransCoef(:,:) * xy_SurfDQVapSatDTemp(:,:)
xyza_QVapMtx(:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2.0_DP * DelTime ) + xy_SurfHumidCoef(:,:) * xy_SurfQVapTransCoef(:,:) + xyr_QMixTransCoef(:,:,k )
end if
xyza_QVapMtx(:,:,k, 1) = - xyr_QMixTransCoef(:,:,k )
do k = 2, kmax-1
xyza_QVapMtx(:,:,k,-1) = - xyr_QMixTransCoef(:,:,k-1)
xyza_QVapMtx(:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2.0_DP * DelTime ) + xyr_QMixTransCoef(:,:,k-1) + xyr_QMixTransCoef(:,:,k )
xyza_QVapMtx(:,:,k, 1) = - xyr_QMixTransCoef(:,:,k )
end do
k = kmax
xyza_QVapMtx(:,:,k,-1) = - xyr_QMixTransCoef(:,:,k-1)
xyza_QVapMtx(:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2.0_DP * DelTime ) + xyr_QMixTransCoef(:,:,k-1)
xyza_QVapMtx(:,:,k, 1) = 0.0_DP
! 鉛直拡散スキームの輸送係数から陰解行列の計算 (比湿を除く混合比)
! Calculate implicit matrices from transfer coefficient of vertical diffusion scheme (mixing ratio except for specific humidity)
!
k = 1
xyza_QMixMtx(:,:,k,-1) = 0.0d0
xyza_QMixMtx(:,:,k, 0) = - ( xyr_Press(:,:,k) - xyr_Press(:,:,k-1) ) / Grav / ( 2.0_DP * DelTime ) + xyr_QMixTransCoef(:,:,k )
xyza_QMixMtx(:,:,k, 1) = - xyr_QMixTransCoef(:,:,k )
do l = -1, 1
do k = 2, kmax
xyza_QMixMtx(:,:,k,l) = xyza_QVapMtx(:,:,k,l)
end do
end do
! 地表面過程の輸送係数から陰解行列の計算
! Calculate implicit matrices from transfer coefficient of surface process
!
do i = 0, imax-1
do j = 1, jmax
if ( xy_SurfCond(i,j) == 1 ) then
if ( FlagBucketModel ) then
xyaa_SurfMtx(i,j,0,-1) = 0.0_DP
xyaa_SurfMtx(i,j,0, 0) = xy_SurfHeatCapacity(i,j) / ( 2.0_DP * DelTime ) + CpDry * xy_SurfTempTransCoef(i,j) + xyra_DelRadLFlux(i,j,0,0)
else
xyaa_SurfMtx(i,j,0,-1) = - LatentHeat * xy_SurfHumidCoef(i,j) * xy_SurfQVapTransCoef(i,j)
xyaa_SurfMtx(i,j,0, 0) = xy_SurfHeatCapacity(i,j) / ( 2.0_DP * DelTime ) + CpDry * xy_SurfTempTransCoef(i,j) + xyra_DelRadLFlux(i,j,0,0) + LatentHeat * xy_SurfHumidCoef(i,j) * xy_SurfQVapTransCoef(i,j) * xy_SurfDQVapSatDTemp(i,j)
end if
xyaa_SurfMtx(i,j,0, 1) = - CpDry * xyr_Exner(i,j,0) / xyz_Exner(i,j,1) * xy_SurfTempTransCoef(i,j) + xyra_DelRadLFlux(i,j,0,1)
else
xyaa_SurfMtx(i,j,0,-1) = 0.0_DP
xyaa_SurfMtx(i,j,0, 0) = 1.0_DP
xyaa_SurfMtx(i,j,0, 1) = 0.0_DP
end if
end do
end do
! 東西風速, 南北風速の計算
! Calculate eastward and northward wind
!
xyza_UVLUMtx = xyza_UVMtx
call PhyImplLUDecomp3( xyza_UVLUMtx, imax * jmax, kmax ) ! (in)
do k = 1, kmax
xyz_DUDt(:,:,k) = - ( xyr_MomFluxX(:,:,k) - xyr_MomFluxX(:,:,k-1) )
xyz_DVDt(:,:,k) = - ( xyr_MomFluxY(:,:,k) - xyr_MomFluxY(:,:,k-1) )
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) = xyza_QVapMtx(:,:,k,l)
end do
k = 0
xyza_TempQVapLUMtx(:,:, k, l) = xyaa_SurfMtx(:,:,0,l)
do k = 1, kmax
xyza_TempQVapLUMtx(:,:, k, l) = xyra_TempMtx(:,:,k,l)
end do
end do
call PhyImplLUDecomp3( xyza_TempQVapLUMtx, imax * jmax, (2 * kmax) + 1 ) ! (in)
do k = 1, kmax
xyz_DelTempQVap(:,:,-k) = - ( xyrf_QMixFlux(:,:,k,IndexH2OVap) - xyrf_QMixFlux(:,:,k-1,IndexH2OVap) )
end do
k = 0
do j = 1, jmax
do i = 0, imax-1
if ( xy_SurfCond(i,j) == 1 ) then
xyz_DelTempQVap(i,j,k) = - xyr_RadSFlux(i,j,0) - xyr_RadLFlux(i,j,0) - xyr_HeatFlux(i,j,0) - LatentHeat * xyrf_QMixFlux(i,j,0,IndexH2OVap) + xy_GroundTempFlux(i,j)
else
xyz_DelTempQVap(i,j,k) = 0.0_DP
end if
end do
end do
do k = 1, kmax
xyz_DelTempQVap(:,:,k) = - ( xyr_HeatFlux(:,:, k) - xyr_HeatFlux(:,:, k-1 ) )
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.0_DP * DelTime )
xyz_DVDt(:,:,k) = xyz_DVDt(:,:,k) / ( 2.0_DP * DelTime )
xyz_DTempDt(:,:,k) = xyz_DelTempQVap(:,:, k) / ( 2.0_DP * DelTime )
xyzf_DQMixDt(:,:,k,IndexH2OVap) = xyz_DelTempQVap(:,:,-k) / ( 2.0_DP * DelTime )
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.0_DP * DelTime )
else
xy_DSurfTempDt(i,j) = 0.0_DP
end if
end do
end do
! 比湿を除く質量混合比
! Calculate mass mixing ratio except for specific humidity
!
xyza_QMixLUMtx = xyza_QMixMtx
call PhyImplLUDecomp3( xyza_QMixLUMtx, imax * jmax, kmax )
do n = 1, ncmax
if ( n == IndexH2OVap ) cycle
do k = 1, kmax
xyzf_DQMixDt(:,:,k,n) = - ( xyrf_QMixFlux(:,:,k,n) - xyrf_QMixFlux(:,:,k-1,n) )
end do
call PhyImplLUSolve3( xyzf_DQMixDt(:,:,:,n), xyza_QMixLUMtx, 1, imax * jmax, kmax )
end do
!#########################################################
! code for debug, this will be removed, (Y. O. Takahashi, 2009/04/07)
!!$ i = 1
!!$ j = jmax / 2
!!$ n = IndexH2OVap
!!$ write( 6, * ) &
!!$ & - xyr_RadSFlux(i,j,0), &
!!$ & - ( xyr_RadLFlux(i,j,0) &
!!$ & + xyra_DelRadLFlux(i,j,0,0) * xy_DSurfTempDt(i,j) * ( 2.0d0 * DelTime ) &
!!$ & + xyra_DelRadLFlux(i,j,0,1) * xyz_DTempDt(i,j,1) * ( 2.0d0 * DelTime ) ), &
!!$ & - ( xyr_HeatFlux(i,j,0) &
!!$ & - CpDry * xyr_Exner(i,j,0) * xy_SurfTempTransCoef(i,j) &
!!$ & * ( xyz_DTempDt(i,j,1) / xyz_Exner(i,j,1) &
!!$ & - xy_DSurfTempDt(i,j) / xyr_Exner(i,j,0) ) * ( 2.0d0 * DelTime ) ), &
!!$ & - LatentHeat &
!!$ & * ( xyrf_QMixFlux(i,j,0,n) &
!!$ & - xy_SurfQVapTransCoef(i,j) &
!!$ & * ( xyzf_DQMixDt(i,j,1,n) &
!!$ & - xy_SurfDQVapSatDTemp(i,j) * xy_DSurfTempDt(i,j) ) &
!!$ & * ( 2.0d0 * DelTime ) ) !, &
!!$! & + xy_GroundTempFlux(i,j)
!!$ write( 6, * ) &
!!$ & - xyr_RadSFlux(i,j,0) &
!!$ & - ( xyr_RadLFlux(i,j,0) &
!!$ & + xyra_DelRadLFlux(i,j,0,0) * xy_DSurfTempDt(i,j) * ( 2.0d0 * DelTime ) &
!!$ & + xyra_DelRadLFlux(i,j,0,1) * xyz_DTempDt(i,j,1) * ( 2.0d0 * DelTime ) ) &
!!$ & - ( xyr_HeatFlux(i,j,0) &
!!$ & - CpDry * xyr_Exner(i,j,0) * xy_SurfTempTransCoef(i,j) &
!!$ & * ( xyz_DTempDt(i,j,1) / xyz_Exner(i,j,1) &
!!$ & - xy_DSurfTempDt(i,j) / xyr_Exner(i,j,0) ) * ( 2.0d0 * DelTime ) ) &
!!$ & - LatentHeat &
!!$ & * ( xyrf_QMixFlux(i,j,0,n) &
!!$ & - xy_SurfQVapTransCoef(i,j) &
!!$ & * ( xyzf_DQMixDt(i,j,1,n) &
!!$ & - xy_SurfDQVapSatDTemp(i,j) * xy_DSurfTempDt(i,j) ) &
!!$ & * ( 2.0d0 * DelTime ) ) !, &
!!$! & + xy_GroundTempFlux(i,j)
!!$
!!$ xy_SurfQVapSat(i,j) = &
!!$ & - xyr_RadSFlux(i,j,0) &
!!$ & - ( xyr_RadLFlux(i,j,0) &
!!$ & + xyra_DelRadLFlux(i,j,0,0) * xy_DSurfTempDt(i,j) * ( 2.0d0 * DelTime ) &
!!$ & + xyra_DelRadLFlux(i,j,0,1) * xyz_DTempDt(i,j,1) * ( 2.0d0 * DelTime ) ) &
!!$ & - ( xyr_HeatFlux(i,j,0) &
!!$ & - CpDry * xyr_Exner(i,j,0) * xy_SurfTempTransCoef(i,j) &
!!$ & * ( xyz_DTempDt(i,j,1) / xyz_Exner(i,j,1) &
!!$ & - xy_DSurfTempDt(i,j) / xyr_Exner(i,j,0) ) * ( 2.0d0 * DelTime ) ) &
!!$ & - LatentHeat &
!!$! & * ( xyr_QVapFlux(i,j,0) &
!!$! & - xy_SurfQVapTransCoef(i,j) &
!!$! & * ( xyz_DQVapDt(i,j,1) &
!!$! & - xy_SurfDQVapSatDTemp(i,j) * xy_DSurfTempDt(i,j) ) &
!!$! & * ( 2.0d0 * DelTime ) )
!!$ & * xyrf_QMixFlux(i,j,0,n)
!!$ write( 6, * ) xy_SurfQVapSat(i,j)
!#########################################################
! 計算時間計測一時停止
! Pause measurement of computation time
!
call TimesetClockStop( module_name )
end subroutine PhyImplTendency