Subroutine : |
|
xz_DensCloud_in(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
DelTime : | real(8), intent(in)
|
pz_VelXNl(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
xr_VelZNl(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
xz_DensCloudNl(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
xz_MasscondNl(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(inout)
: | 凝結量 [kg/m^3 s] =end == Work
|
|
xz_DensCloud_out(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(out)
|
subroutine DensityCloud( xz_DensCloud_in, DelTime, pz_VelXNl, xr_VelZNl, xz_DensCloudNl, xz_MasscondNl, xz_DensCloud_out) !(out) 雲の密度
!=begin
!== Dependancy
use dc_trace, only: BeginSub, EndSub, DbgMessage
use gridset, only: DimXMin, DimXMax, DimZMin, DimZMax
! use bcset, only: xz_BC
use average, only: pz_avr_xz, xr_avr_xz
use differentiate_center4, only: xz_dx_pz, xz_dz_xr
!=end
!== 暗黙の型宣言禁止
implicit none
!=begin
!== Input
real(8), intent(in) :: DelTime ! leapfrog スキームの時間間隔 [s]
real(8), intent(in) :: xz_DensCloud_in(DimXMin:DimXMax, DimZMin:DimZMax)
! 雲の密度 [kg/m^3]
real(8), intent(in) :: pz_VelXNl(DimXMin:DimXMax, DimZMin:DimZMax)
! 水平速度 [m/s]
real(8), intent(in) :: xr_VelZNl(DimXMin:DimXMax, DimZMin:DimZMax)
! 鉛直速度 [m/s]
real(8), intent(in) :: xz_DensCloudNl(DimXMin:DimXMax, DimZMin:DimZMax)
! 雲の密度 [kg/m^3]
!== Output
real(8), intent(out) :: xz_DensCloud_out(DimXMin:DimXMax, DimZMin:DimZMax)
! 雲の密度 [kg/m^3]
!== In/Out
real(8), intent(inout):: xz_MasscondNl(DimXMin:DimXMax, DimZMin:DimZMax)
! 凝結量 [kg/m^3 s]
!=end
!== Work
real(8) :: xz_FluxDensCloud(DimXMin:DimXMax, DimZMin:DimZMax)
! フラックスによる雲密度の時間変化率 [kg/m^3 s]
real(8) :: xz_SrcDensCloud(DimXMin:DimXMax, DimZMin:DimZMax)
! 生成消滅による雲密度の時間変化率 [kg/m^3 s]
call BeginSub("DensityCloud", fmt="%c", c1="Time integration of density of cloud.")
! write(*,*) xz_DensCloudNl(1,:)
!=== フラックス項の計算. 4 次精度中心差分を用いて計算
xz_FluxDensCloud = - xz_dx_pz(pz_VelXNl * pz_avr_xz(xz_DensCloudNl)) - xz_dz_xr(xr_VelZNl * xr_avr_xz(xz_DensCloudNl))
! xz_FluxDensCloud = 0.0d0
call DbgMessage( fmt="*** Debug Message [DensityCloud] *** xz_SrcDensCloud ")
!=== 生成項の計算
xz_SrcDensCloud = xz_MasscondNl
call DbgMessage( fmt="*** Debug Message [DensityCloud] *** xz_DensCloud_out ")
xz_DensCloud_out = xz_DensCloud_in + DelTime * (xz_FluxDensCloud + xz_SrcDensCloud)
call DbgMessage( fmt="*** Debug Message [DensityCloud] *** xz_DensCloud_out(adjustment) ")
! 雲の量が負になる場合は 0 にする.
! 差し引き分だけ蒸発量を減らす
xz_DensCloud_out = max(xz_DensCloud_out, 0.0d0)
xz_MasscondNl = max(xz_MasscondNl, - (xz_DensCloudNl / DelTime + xz_FluxDensCloud))
! !=== 境界条件
! call boundary(xz_BC, xz_DensCloud_out)
!
call EndSub("DensityCloud")
end subroutine DensityCloud