Class | ECCM |
In: |
moist/eccm.f90
|
Copyright (C) GFD Dennou Club, 2005. All rights reserved.
* Developer: SUGIYAMA Ko-ichiro * Version: $Id: eccm.f90,v 1.1.1.1 2006/04/25 03:43:58 deepconv Exp $ * Tag Name: $Name: $ * Change History:
断熱的に上昇する気塊の温度減率を計算し, 静水圧平衡から圧力を決める
* 比熱は定数, 平均分子量は配列 * オイラースキームだと精度が足りないので, ルンゲクッタスキームを利用.
subroutine ECCM_DTempDZ ( Temp, Press, MolFr, MolWtMean, DTempDZ ) !暗黙の型宣言禁止 implicit none !変数定義 real(8), intent(in) :: Temp real(8), intent(in) :: Press real(8), intent(inout) :: MolFr(0:SpcNum) !モル分率 real(8), intent(out):: MolWtMean real(8), intent(out):: DTempDZ real(8) :: MolWt(0:SpcNum) real(8) :: ReactHeat real(8) :: Heat(SpcNum) real(8) :: DelMolFr real(8) :: SatPress real(8) :: VapPress real(8) :: A, B integer :: s !初期化 DTempDZ = 0.0d0 ReactHeat = 0.0d0 Heat = 0.0d0 DelMolFr = 0.0d0 SatPress = 0.0d0 VapPress = 0.0d0 MolWt(0) = MolWtDry MolWt(1:SpcNum) = MolWtWet(1:SpcNum) !------------------------------------------------------------ !NH4SH 以外の化学種の平衡条件 !------------------------------------------------------------ do s = 1, LoopNum !飽和蒸気圧 SatPress = SvapPress( SpcWetID(CloudNum(s)), Temp ) !潜熱. Heat(GasNum(s)) = & & LatentHeatPerMol( SpcWetID(CloudNum(s)), Temp ) !元々のモル分率を用いて現在の蒸気圧を計算 VapPress = MolFr(GasNum(s)) * Press !飽和蒸気圧から凝結の有無を決める if ( VapPress < SatPress ) then !モル比は変化しない MolFr(GasNum(s)) = MolFr(GasNum(s)) !凝結していないので潜熱なし. Heat(GasNum(s)) = 0.0d0 else !飽和蒸気圧と圧力から現在のモル比を計算 MolFr(GasNum(s)) = max(SatPress / Press, 1.0d-16) end if end do !------------------------------------------------------------ !NH4SH の平衡条件 !------------------------------------------------------------ if ( NH3Num /= 0 ) then DelMolFr = DelMolFrNH4SH(Temp, Press, & & MolFr(NH3Num), MolFr(H2SNum) ) MolFr(NH3Num) = MolFr(NH3Num) - DelMolFr MolFr(H2SNum) = MolFr(H2SNum) - DelMolFr ReactHeat = ReactHeatNH4SHPerMol * DelMolFr end if !------------------------------------------------------------ !温度勾配を計算 !------------------------------------------------------------ !モル比 MolFr(0) = 1.0d0 - sum( MolFr(1:SpcNum) ) !分子量 MolWtMean = dot_product( MolWt, MolFr ) !係数. 温度 Temp(i) で評価 A = dot_product( Heat(1:SpcNum), MolFr(1:SpcNum)) & & / ( GasRUniv * Temp ) B = dot_product(( Heat(1:SpcNum) ** 2.0d0), MolFr(1:SpcNum)) & & / ( CpDryMol * GasRUniv * ( Temp ** 2.0d0 ) ) !断熱温度減率 DTempDZ = - Grav * MolWtMean * ( 1.0d0 + A ) & & / ( CpDryMol * ( 1.0d0 + B ) ) & & + ReactHeat / ( CpDryMol * DelZ ) end subroutine
subroutine ECCM_Init ( ) !暗黙の型宣言禁止 implicit none !変数定義 integer :: s integer :: n1 !----------------------------------------------------------- ! 雲粒と気体の ID の組を作る !----------------------------------------------------------- !初期化 LoopNum = 0 !化学種の中から雲粒を作るものを選び, その配列添え字と分子量を保管. SelectCloud: do s = 1, SpcNum ! NH4SH については無視 if ( trim(SpcWetSymbol(s)) == 'NH4SH-s-Cloud' ) then cycle SelectCloud end if !'Cloud' という文字列が含まれるものの個数を数える n1 = index(SpcWetSymbol(s), '-Cloud' ) if (n1 /= 0) then LoopNum = LoopNum + 1 CloudNum(LoopNum)= s GasNum(LoopNum) = minloc(SpcWetID, 1, SpcWetID == ChemData_OneSpcID(SpcWetSymbol(s)(1:n1-3) // '-g')) end if end do SelectCloud !----------------------------------------------------------- ! 硫化アンモニウム, およびアンモニアと硫化水素の ID を取得 !----------------------------------------------------------- NH3Num = minloc(SpcWetID, 1, SpcWetID == ChemData_OneSpcID('NH3-g')) H2SNum = minloc(SpcWetID, 1, SpcWetID == ChemData_OneSpcID('H2S-g')) !----------------------------------------------------------- ! 確認 !----------------------------------------------------------- if ( LoopNum == 0 ) then write(*,*) "MoistAdjust: CloudNum = 0, please comment out of MoistAdjust" stop end if write(*,*) "ECCM_Init, LoopNum: ", LoopNum write(*,*) "ECCM_Init, CloudNum: ", CloudNum write(*,*) "ECCM_Init, GasNum: ", GasNum write(*,*) "ECCM_Init, NH3Num: ", NH3Num write(*,*) "ECCM_Init, H2SNum: ", H2SNum end subroutine
subroutine ECCM_MolFr ( MolFrIni, Humidity, z_Temp, z_Press, z_MolFr ) !暗黙の型宣言禁止 implicit none real(8), intent(in):: MolFrIni(1:SpcNum) real(8), intent(in):: Humidity real(8), intent(in) :: z_Temp(DimZMin:DimZMax) real(8), intent(in) :: z_Press(DimZMin:DimZMax) real(8), intent(out):: z_MolFr(DimZMin:DimZMax, 1:SpcNum) real(8) :: SatPress real(8) :: DelMolFr integer :: k, s !----------------------------------------------------------- ! 配列の初期化 !----------------------------------------------------------- !初期モル比 z_MolFr = 0.0d0 z_MolFr(RegZMin, 1:SpcNum) = MolFrIni(1:SpcNum) z_MolFr(RegZMin-1, 1:SpcNum) = MolFrIni(1:SpcNum) !----------------------------------------------------------- ! 断熱減率 dT/dz の計算. !----------------------------------------------------------- do k = RegZMin, DimZMax z_MolFr(k,:) = z_MolFr(k-1,:) !------------------------------------------------------------ !NH4SH 以外の化学種の平衡条件 !------------------------------------------------------------ do s = 1, LoopNum !飽和蒸気圧 SatPress = SvapPress( SpcWetID(CloudNum(s)), z_Temp(k) ) !飽和蒸気圧と圧力からモル比を計算 ! 前のステップでの値よりも大きくならないようにする z_MolFr(k,GasNum(s)) = & & min( SatPress * Humidity / z_Press(k), z_MolFr(k-1,GasNum(s)) ) end do !------------------------------------------------------------ !NH4SH の平衡条件 !------------------------------------------------------------ if ( NH3Num /= 0 ) then DelMolFr = DelMolFrNH4SH(z_Temp(k), z_Press(k), & & z_MolFr(k,NH3Num), z_MolFr(k,H2SNum) ) z_MolFr(k,NH3Num) = z_MolFr(k,NH3Num) - DelMolFr z_MolFr(k,H2SNum) = z_MolFr(k,H2SNum) - DelMolFr end if end do end subroutine
subroutine ECCM_Temp_Press ( MolFrIni, z_Temp, z_Press ) !暗黙の型宣言禁止 implicit none real(8), intent(in) :: MolFrIni(1:SpcNum) real(8), intent(out):: z_Temp(DimZMin:DimZMax) real(8), intent(out):: z_Press(DimZMin:DimZMax) real(8) :: z_MolFr(DimZMin:DimZMax, 0:SpcNum) !モル分率 real(8) :: z_DTempDZ(DimZMin:DimZMax) real(8) :: MolFr(0:SpcNum) real(8) :: z_MolWtMean(DimZMin:DimZMax) integer :: k real(8) :: Temp1, Press1, MolWtMean1, DTempDZ1 real(8) :: Temp2, Press2, MolWtMean2, DTempDZ2 real(8) :: Temp3, Press3, MolWtMean3, DTempDZ3 real(8) :: Temp4, Press4, MolWtMean4, DTempDZ4 !------------------------------------------------------------- ! 配列の初期化 !------------------------------------------------------------- !初期モル比 z_MolFr = 0.0d0 z_MolFr(RegZMin, 1:SpcNum) = MolFrIni(1:SpcNum) z_MolFr(RegZMin, 0) = 1.0d0 - sum(MolFrIni) z_MolFr(RegZMin-1, 1:SpcNum) = MolFrIni(1:SpcNum) z_MolFr(RegZMin-1, 0) = 1.0d0 - sum(MolFrIni) !地表面での温度(RegZMin は, 高度 DelZ / 2 に相当) z_Temp = 0.0d0 z_Temp(RegZMin) = TempSfc - Grav * MolWtDry / CpDryMol * DelZ * 5.0d-1 !地表面での圧力(RegZMin は, 高度 DelZ / 2 に相当) z_Press = 0.0d0 z_Press(RegZMin) = PressSfc *((TempSfc / z_Temp(RegZMin)) & & ** ( - CpDryMol / GasRUniv )) !分子量の初期化 z_MolWtMean = 0.0d0 !----------------------------------------------------------- ! 断熱減率 dT/dz の計算. !----------------------------------------------------------- do k = RegZMin, DimZMax-1 !初期化 z_MolFr(k,:) = z_MolFr(k-1,:) !============================================================== !湿潤断熱減率をルンゲクッタ法を用いて計算 Temp1 = z_Temp(k) Press1 = z_Press(k) MolFr = z_MolFr(k-1,:) call ECCM_DTempDZ( Temp1, Press1, z_MolFr(k,:), MolWtMean1, DTempDZ1 ) Temp2 = Temp1 + DTempDZ1 * DelZ * 5.0d-1 Press2 = & & Press1 * ( ( Temp1 / Temp2 ) & & ** (Grav * MolWtMean1 / (GasRUniv * DTempDZ1) ) ) call ECCM_DTempDZ( Temp2, Press2, MolFr, MolWtMean2, DTempDZ2 ) Temp3 = Temp1 + DTempDZ2 * DelZ * 5.0d-1 Press3 = & & Press1 * ( ( Temp1 / Temp3 ) & & ** (Grav * MolWtMean2 / (GasRUniv * DTempDZ2) ) ) call ECCM_DTempDZ( Temp3, Press3, MolFr, MolWtMean3, DTempDZ3 ) Temp4 = Temp1 + DTempDZ3 * DelZ Press4 = & & Press1 * ( ( Temp1 / Temp4 ) & & ** (Grav * MolWtMean3 / (GasRUniv * DTempDZ3) ) ) call ECCM_DTempDZ( Temp4, Press4, MolFr, MolWtMean4, DTempDZ4 ) z_DTempDZ(k) = ( + DTempDZ1 + DTempDZ2 * 2.0d0 & & + DTempDZ3 * 2.0d0 + DTempDZ4 ) / 6.0d0 z_MolWtMean(k) = ( + MolWtMean1 + MolWtMean2 * 2.0d0 & & + MolWtMean3 * 2.0d0 + MolWtMean4 ) / 6.0d0 !============================================================== !温度を計算 z_Temp(k+1) = z_Temp(k) + z_DTempDz(k) * DelZ !圧力を静水圧平衡から計算 z_Press(k+1) = & & z_Press(k) * ( ( z_Temp(k) / z_Temp(k+1)) & & ** (Grav * z_MolWtMean(k) / (GasRUniv * z_DTempDZ(k)) ) ) ! write(*,*) 'ECCM, MolFr: ', k, z_MolFr(k,1:3) ! write(*,*) 'ECCM, Temp: ', k, z_Temp(k) ! write(*,*) 'ECCM, Press: ', k, z_Press(k) ! write(*,*) 'ECCM, DTempDZ: ', k, z_DTempDZ(k) ! write(*,*) 'ECCM, DTempDZ: ', k, DTempDZ1, DTempDZ2, DTempDZ3, DTempDZ4 !モル比 (DimZMax でも値を持つようにするため) z_MolFr(k+1,:) = z_MolFr(k,:) z_MolFr(k+1,0) = 1.0d0 - sum( z_MolFr(k,1:SpcNum) ) end do end subroutine