Subroutine : |
|
DelTimeShort : | real(8), intent(in)
|
xz_SatRatioNs(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
xz_PotTempNs(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
xz_ExnerNs(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
xz_DensCloudAs(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(inout)
|
xz_MassCondNs(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(out)
|
subroutine MassCondense( DelTimeShort, xz_SatRatioNs, xz_PotTempNs, xz_ExnerNs, xz_DensCloudAs, xz_MassCondNs) !(out) 凝結量
!=begin
!==Dependency
use dc_trace, only: BeginSub, EndSub, DbgMessage
use gridset, only: DimXMin, DimXMax, DimZMin, DimZMax
use cloudset, only: DensIce, NumAerosol, RadiAerosol, Kd, SatRatioCr ! 臨界飽和比
use basicset, only: xz_PotTempBasicZ, xz_ExnerBasicZ, xz_DensBasicZ, GasRDry, PressSfc, CpDry, LatHeatPerMass ! CO2 の単位質量あたりの潜熱
use StoreDensCloud, only: StoreDensCloudCond
!=end
!==暗黙の型宣言を禁止
implicit none
!==Input
real(8), intent(in) :: DelTimeShort ! 時間ステップ [s]
real(8), intent(in) :: xz_SatRatioNs(DimXMin:DimXMax, DimZMin:DimZMax)
! 飽和比 [1]
real(8), intent(in) :: xz_PotTempNs(DimXMin:DimXMax, DimZMin:DimZMax)
! 温位 [K]
real(8), intent(in) :: xz_ExnerNs(DimXMin:DimXMax, DimZMin:DimZMax)
! 無次元圧力 [1]
real(8), intent(inout) :: xz_DensCloudAs(DimXMin:DimXMax, DimZMin:DimZMax)
! 雲の密度 [kg/m^3]
!==Output
real(8), intent(out) :: xz_MassCondNs(DimXMin:DimXMax, DimZMin:DimZMax)
! 凝結量 [kg/m^3 s]
!==Work
real(8) :: xz_Rh(DimXMin:DimXMax, DimZMin:DimZMax)
! 熱拡散による抵抗係数
real(8) :: xz_RadiCloud(DimXMin:DimXMax, DimZMin:DimZMax)
! 雲粒の半径 [m]
real, parameter :: Pi = 3.1415926535897932385d0
! 円周率
! real(8) :: DensCloudThreshold = 1.0d-5
! real(8) :: DensCloudThreshold = 5.0d-6
real(8) :: DensCloudThreshold = 5.0d-5
integer :: i,k ! ループ変数
call BeginSub("MassCondense", fmt="%c", c1="Calculate latent heat per unit mass.")
call DbgMessage( fmt="*** Debug Message [MassCondense] *** xz_Rh ")
!write(*,*) xz_PotTempNs(1,:)
xz_Rh = LatHeatPerMass ** 2.0d0 / ( Kd * GasRDry * (xz_ExnerBasicZ + xz_ExnerNs)** 2.0d0 * (xz_PotTempBasicZ + xz_PotTempNs)**2.0d0 )
! 山下, 20080410,
! L**2 / (k * R * Exner**2 * theta**2) という値を計算している.
! 参考文献 : 北守修論の付録 A に導出過程が記述されている.
! 但し, 北守修論本文の Rh の式 (3.4) には誤植があるので注意されたい.
! 教科書として 1 次参照元はどれなのか調べる必要がありそうだ.
call DbgMessage( fmt="*** Debug Message [MassCondense] *** xz_RadiCloud ")
xz_RadiCloud = ( 3.0d0 * xz_DensCloudAs / (4.0d0 * Pi * DensIce * NumAerosol * xz_DensBasicZ) + RadiAerosol**3.0d0 )**(1.0d0 / 3.0d0)
do k= DimZMin, DimZMax
do i = DimXMin, DimXMax
if (xz_SatRatioNs(i,k) - SatRatioCr > epsilon(0.0d0)) then
! 飽和比が定めた臨界飽和比よりも大きい場合
! CO2 の飽和比の式が地球大気で用いられる飽和比の式と等価なものなのか確認する必要がある.
xz_MassCondNs(i,k) = max( - xz_DensCloudAs(i,k) / DelTimeShort, 4.0d0 * Pi * xz_RadiCloud(i,k) * NumAerosol * xz_DensBasicZ(i,k) / xz_Rh(i,k) * (xz_SatRatioNs(i,k) - 1.0d0) )
else if (xz_DensCloudAs(i,k) > DensCloudThreshold ) then
! 臨界飽和比を超えていないが, 雲密度が閾値以上である場合
! 飽和比が 1 以上のときに凝結, 1 以下のときに蒸発.
xz_MassCondNs(i,k) = max( - xz_DensCloudAs(i,k) / DelTimeShort, 4.0d0 * Pi * xz_RadiCloud(i,k) * NumAerosol * xz_DensBasicZ(i,k) / xz_Rh(i,k) * (xz_SatRatioNs(i,k) - 1.0d0) )
else if (xz_DensCloudAs(i,k) /= 0.0d0 .and. xz_SatRatioNs(i,k) < 1.0d0 ) then
! 雲密度が閾値未満で, 飽和比 1.0 未満である場合に蒸発.
xz_MassCondNs(i,k) = max( - xz_DensCloudAs(i,k) / DelTimeShort, 4.0d0 * Pi * xz_RadiCloud(i,k) * NumAerosol * xz_DensBasicZ(i,k) / xz_Rh(i,k) * (xz_SatRatioNs(i,k) - 1.0d0) )
else
! それ以外の場合は何も生じない.
xz_MassCondNs(i,k) = 0.0d0
end if
! call StoreDensCloudCond( xz_MassCondNs(i,k) )
xz_DensCloudAs(i,k) = xz_DensCloudAs(i,k) + DelTimeShort * xz_MassCondNs(i,k)
end do
end do
call StoreDensCloudCond( xz_MassCondNs )
call EndSub("MassCondense")
end subroutine MassCondense