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_DensCloudBl(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
xz_PotTempNl(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
xz_ExnerNl(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
xz_KhBl(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
xz_DensCloud_NoFillNeg(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(out)
|
subroutine DensityCloud( xz_DensCloud_in, DelTime, pz_VelXNl, xr_VelZNl, xz_DensCloudNl, xz_DensCloudBl, xz_PotTempNl, xz_ExnerNl, xz_KhBl, xz_DensCloud_NoFillNeg) !(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
use basicset, only: Grav, PressSfc, GasRDry, CpDry, xz_ExnerBasicZ, xz_PotTempBasicZ, xz_DensBasicZ
use cloudset, only: DensIce, NumAerosol, RadiAerosol ! 凝結核半径
use Turbulence, only: xz_TurbScalar
use NumDiffusion, only: xz_NumDiffScalar
use StoreDensCloud, only: StoreDensCloudAdv, StoreDensCloudTurb, StoreDensCloudDiff, StoreDensCloudFall
!=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_ExnerNl(DimXMin:DimXMax, DimZMin:DimZMax)
! エクスナー関数の擾乱成分 [1]
real(8), intent(in) :: xz_PotTempNl(DimXMin:DimXMax, DimZMin:DimZMax)
! 温位の擾乱成分 [K]
real(8), intent(in) :: xz_DensCloudNl(DimXMin:DimXMax, DimZMin:DimZMax)
! 雲の密度 [kg/m^3]
real(8), intent(in) :: xz_DensCloudBl(DimXMin:DimXMax, DimZMin:DimZMax)
! 雲の密度 [kg/m^3]
real(8), intent(in) :: xz_KhBl(DimXMin:DimXMax, DimZMin:DimZMax)
! 乱流拡散係数 [kg/m^3]
!== Output
real(8), intent(out) :: xz_DensCloud_NoFillNeg(DimXMin:DimXMax, DimZMin:DimZMax)
! 雲の密度(負の密度調整前) [kg/m^3]
!=end
!== Work
real(8) :: xz_FluxDensCloud(DimXMin:DimXMax, DimZMin:DimZMax)
! フラックスによる雲密度の時間変化率 [kg/m^3 s]
real(8) :: xz_FallDensCloud(DimXMin:DimXMax, DimZMin:DimZMax)
! 雲粒落下による雲密度の時間変化率 [kg/m^3 s]
real(8) :: xz_TurbDensCloud(DimXMin:DimXMax, DimZMin:DimZMax)
! 乱流拡散による雲密度の時間変化率 [kg/m^3 s]
real(8) :: xz_NumDensCloud(DimXMin:DimXMax, DimZMin:DimZMax)
! 数値粘性による雲密度の時間変化率 [kg/m^3 s]
real(8) :: xz_RadiCloudNl(DimXMin:DimXMax, DimZMin:DimZMax)
! 雲粒の半径 [m]
real(8) :: xz_VelZTerm(DimXMin:DimXMax, DimZMin:DimZMax)
! 雲粒の終端速度
real(8) :: xz_VisCoffCO2(DimXMin:DimXMax, DimZMin:DimZMax)
! CO2 の粘性係数
real(8) :: xz_MeanFreePath(DimXMin:DimXMax, DimZMin:DimZMax)
! CO2 の平均自由行程
real(8) :: xz_KnudsenNum(DimXMin:DimXMax, DimZMin:DimZMax)
! 雲粒に対する Knudsen 数
real, parameter :: Pi = 3.1415926535897932385d0
! 円周率
real, parameter :: BoltzmannConst = 1.38065d-23
! ボルツマン定数
real, parameter :: CO2Diam = 3.30d-10
! CO2 分子の kinetic diameter
real, parameter :: SutherlandConst = 2.40d2
! CO2 分子のサザランド定数
real, parameter :: VisCoffRef = 1.47d-5
! 粘性係数の参照値
real, parameter :: TempRef = 2.93d2
! 粘性係数の参照温度
call BeginSub("DensityCloud", fmt="%c", c1="Time integration of density of cloud.")
! write(*,*) xz_DensCloudNl(1,:)
!=== 粘性係数の計算
xz_VisCoffCO2 = VisCoffRef * (TempRef + SutherlandConst) / ( (xz_PotTempBasicZ + xz_PotTempNl)*(xz_ExnerBasicZ + xz_ExnerNl) + SutherlandConst ) * ( (xz_PotTempBasicZ + xz_PotTempNl)*(xz_ExnerBasicZ + xz_ExnerNl) / TempRef )**1.5d0
!=== 平均自由行程の計算
xz_MeanFreePath = BoltzMannConst * (xz_PotTempBasicZ + xz_PotTempNl)*(xz_ExnerBasicZ + xz_ExnerNl) / ( 2.0**0.5d0 * Pi * (CO2Diam)**2.0d0 * PressSfc * (xz_ExnerBasicZ + xz_ExnerNl)**(CpDry/GasRDry) )
!=== 雲粒半径の計算
xz_RadiCloudNl = ( 3.0d0 * xz_DensCloudNl / (4.0d0 * Pi * DensIce * NumAerosol * xz_DensBasicZ) + RadiAerosol**3.0d0 )**(1.0d0 / 3.0d0)
!=== Knudsen 数の計算
xz_KnudsenNum = xz_MeanFreePath / xz_RadiCloudNl
!=== 終端速度の計算
xz_VelZTerm = ( 1.0d0 + xz_KnudsenNum * (1.257d0 + 0.400d0 * exp(- 1.10d0 / xz_KnudsenNum) ) ) * 2.0d0 * xz_RadiCloudNl * xz_RadiCloudNl * Grav * DensIce / (9.0d0 * xz_VisCoffCO2)
! xz_VelZTerm = 1.0d0 ! テスト用
!=== 落下項の計算. 4 次精度中心差分を用いて計算
xz_FallDensCloud = xz_dz_xr(xr_avr_xz(xz_VelZTerm) * xr_avr_xz(xz_DensCloudNl))
call StoreDensCloudFall( xz_FallDensCloud )
!=== フラックス項の計算. 4 次精度中心差分を用いて計算
xz_FluxDensCloud = - xz_dx_pz(pz_VelXNl * pz_avr_xz(xz_DensCloudNl)) - xz_dz_xr(xr_VelZNl * xr_avr_xz(xz_DensCloudNl))
call StoreDensCloudAdv( xz_FluxDensCloud )
call DbgMessage( fmt="*** Debug Message [DensityCloud] *** xz_DensCloud_out ")
!=== 乱流拡散項の計算
xz_TurbDensCloud = xz_TurbScalar(xz_DensCloudBl, xz_KhBl)
call StoreDensCloudTurb( xz_TurbDensCloud )
!=== 数値粘性項の計算
xz_NumDensCloud = xz_NumDiffScalar(xz_DensCloudBl)
call StoreDensCloudDiff( xz_NumDensCloud )
xz_DensCloud_NoFillNeg = xz_DensCloud_in + DelTime * ( xz_FallDensCloud + xz_FluxDensCloud + xz_TurbDensCloud + xz_NumDensCloud )
! !=== 境界条件
! call boundary(xz_BC, xz_DensCloud_out)
!
call EndSub("DensityCloud")
end subroutine DensityCloud