Class | ChemCalc |
In: |
moist/chemcalc.f90
|
化学関連の諸量を計算するためのモジュール. AMP と Antoine の飽和蒸気圧式を用いて以下を求める. デフォルトでは AMP 式を使うようにしてある.
* 飽和蒸気圧 * 飽和蒸気圧の温度微分 * 潜熱
Subroutine : |
初期化ルーチン
subroutine ChemCalc_Init( ) ! !初期化ルーチン ! !暗黙の型宣言禁止 implicit none !入出力変数 character(20) :: Name integer :: id !----------------------------------------------------------- ! 初期化 !----------------------------------------------------------- !Antoine の飽和蒸気圧式の係数 a_antA = ChemData_SvapPress_AntoineA a_antB = ChemData_SvapPress_AntoineB a_antC = ChemData_SvapPress_AntoineC a_antU = ChemData_SvapPress_AntoineUnit !AMP 式の飽和蒸気圧式の係数 a_ampA = ChemData_SvapPress_AMPA a_ampB = ChemData_SvapPress_AMPB a_ampC = ChemData_SvapPress_AMPC a_ampD = ChemData_SvapPress_AMPD a_ampE = ChemData_SvapPress_AMPE !分子量 a_MolWt = ChemData_MolWt !NH4SH の反応熱の初期化 ! NH4SH 1kg に対する反応熱にする. Name = 'NH4SH-s' id = ChemData_OneSpcID( Name ) ReactHeatNH4SHPerMol = GasRUniv * 10834.0d0 ReactHeatNH4SH = GasRUniv * 10834.0d0 / MolWt( id ) !----------------------------------------------------------- ! 確認 !----------------------------------------------------------- if (cpurank == 0) then call MessageNotify( "M", "ChemCalcSpc_Init", "SpcNum = %d", i=(/SpcNum/) ) call MessageNotify( "M", "ChemCalc_Init", "ReactHeatNH4SH = %f", d=(/ReactHeatNH4SH/) ) call MessageNotify( "M", "ChemCalc_Init", "NH4SH MolWt = %f", d=(/MolWt(id)/) ) end if end subroutine ChemCalc_Init
Function : | |||
CpPerMolRef : | real(8)
| ||
ID : | integer, intent(in)
|
引数で与えられた化学種に対して, 標準状態での単位モル当たりの定圧比熱を計算
function CpPerMolRef(ID) ! !引数で与えられた化学種に対して, 標準状態での単位モル当たりの定圧比熱を計算 ! !暗黙の型宣言禁止 implicit none !入出力変数 real(8) :: CpPerMolRef !標準状態での単位モル当たりの比熱 integer, intent(in) :: ID !化学種の ID !データベースから情報取得 CpPerMolRef = ChemData_CpPerMolRef(ID) end function CpPerMolRef
Function : | |||
CpRef : | real(8)
| ||
ID : | integer, intent(in)
|
引数で与えられた化学種に対して, 標準状態での単位質量当たりの定圧比熱を計算
function CpRef(ID) ! !引数で与えられた化学種に対して, 標準状態での単位質量当たりの定圧比熱を計算 ! !暗黙の型宣言禁止 implicit none !入出力変数 real(8) :: CpRef !標準状態での単位質量当たりの比熱 integer, intent(in) :: ID !化学種の ID !データベースから情報取得 CpRef = ChemData_CpRef(ID) end function CpRef
Function : | |||
CvRef : | real(8)
| ||
ID : | integer, intent(in)
|
引数で与えられた化学種に対して, 標準状態での単位質量当たりの定圧比熱を計算
function CvRef(ID) ! !引数で与えられた化学種に対して, 標準状態での単位質量当たりの定圧比熱を計算 ! !暗黙の型宣言禁止 implicit none !入出力変数 real(8) :: CvRef !標準状態での単位質量当たりの比熱 integer, intent(in) :: ID !化学種の ID !データベースから情報取得 CvRef = ChemData_CvRef(ID) end function CvRef
Function : | |||
LatentHeatPerMol : | real(8)
| ||
ID : | integer, intent(in)
| ||
Temp : | real(8), intent(in)
|
引数で与えられた化学種と温度に対して, 潜熱を計算
function LatentHeatPerMol(ID, Temp) ! !引数で与えられた化学種と温度に対して, 潜熱を計算 ! !暗黙の型宣言禁止 implicit none !入出力変数 real(8) :: LatentHeatPerMol !潜熱 real(8), intent(in) :: Temp !温度 integer, intent(in) :: ID !化学種名 !内部変数 real(8) :: DLogSvapPressDTemp real(8), parameter :: GasRUniv = 8.314d0 !普遍気体定数 !飽和蒸気圧の温度微分 DLogSvapPressDTemp = - a_ampA(ID) / (Temp ** 2.0d0) + a_ampC(ID) / Temp + a_ampD(ID) + a_ampE(ID) * 2.0d0 * Temp !潜熱の計算 LatentHeatPerMol = DLogSvapPressDTemp * GasRUniv * (Temp ** 2.0d0) end function LatentHeatPerMol
Function : | |||
SvapPress : | real(8)
| ||
ID : | integer, intent(in)
| ||
Temp : | real(8), intent(in)
|
引数で与えられた化学種と温度に対して, 飽和蒸気圧を計算. AMP 式を利用
function SvapPress(ID, Temp) ! !引数で与えられた化学種と温度に対して, 飽和蒸気圧を計算. AMP 式を利用 ! !暗黙の型宣言禁止 implicit none !入出力変数 real(8) :: SvapPress !飽和蒸気圧 real(8), intent(in) :: Temp !温度 [K] integer, intent(in) :: ID !化学種の ID !内部変数 real(8) :: LogSvapPress !飽和蒸気圧の対数を計算 !対数が大きくなりすぎないようにする. ! Fujitsu Fortran Compiler では 700 より大きい数の exp を取ると警告が出る. LogSvapPress = min( ( a_ampA(ID) / Temp + a_ampB(ID) + a_ampC(ID) * dlog( Temp ) + a_ampD(ID) * Temp + a_ampE(ID) * ( temp ** 2 ) + dlog(1.0d-1) ), 7.0d2 ) !飽和蒸気圧を計算 SvapPress = dexp( LogSvapPress ) end function SvapPress
Function : | |||
xz_DSvapPressDTemp(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8)
| ||
ID : | integer, intent(in)
| ||
xz_Temp(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
引数で与えられた化学種と温度に対して, 飽和蒸気圧の温度微分を計算
function xz_DSvapPressDTemp(ID, xz_Temp) ! !引数で与えられた化学種と温度に対して, 飽和蒸気圧の温度微分を計算 ! !暗黙の型宣言禁止 implicit none !入出力変数 real(8) :: xz_DSvapPressDTemp(DimXMin:DimXMax, DimZMin:DimZMax) !飽和蒸気圧の温度微分 [Pa/K] real(8), intent(in) :: xz_Temp(DimXMin:DimXMax, DimZMin:DimZMax) !温度 [K] integer, intent(in) :: ID !化学種名 !内部変数 real(8) :: xz_LogSvapPress(DimXMin:DimXMax, DimZMin:DimZMax) real(8) :: xz_DLogSvapPressDTemp(DimXMin:DimXMax, DimZMin:DimZMax) !飽和蒸気圧の対数を計算 !対数が大きくなりすぎないようにする. ! Fujitsu Fortran Compiler では 700 より大きい数の exp を取ると警告が出る. xz_LogSvapPress = min( ( a_ampA(ID) / xz_Temp + a_ampB(ID) + a_ampC(ID) * dlog( xz_Temp ) + a_ampD(ID) * xz_Temp + a_ampE(ID) * ( xz_temp ** 2 ) + dlog(1.0d-1) ), 7.0d2 ) !飽和蒸気圧の温度微分 xz_DLogSvapPressDTemp = - a_ampA(ID) / (xz_Temp ** 2.0d0) + a_ampC(ID) / xz_Temp + a_ampD(ID) + a_ampE(ID) * 2.0d0 * xz_Temp !飽和蒸気圧の温度微分 xz_DSvapPressDTemp = xz_DLogSvapPressDTemp * dexp( xz_LogSvapPress ) end function xz_DSvapPressDTemp
Function : | |||
xz_LatentHeat(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8)
| ||
ID : | integer, intent(in)
| ||
xz_Temp(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
引数で与えられた化学種と温度に対して, 潜熱を計算
function xz_LatentHeat(ID, xz_Temp) ! !引数で与えられた化学種と温度に対して, 潜熱を計算 ! !暗黙の型宣言禁止 implicit none !入出力変数 real(8) :: xz_LatentHeat(DimXMin:DimXMax, DimZMin:DimZMax) !潜熱 real(8), intent(in) :: xz_Temp(DimXMin:DimXMax, DimZMin:DimZMax) !温度 integer, intent(in) :: ID !化学種名 !内部変数 real(8) :: xz_DLogSvapPressDTemp(DimXMin:DimXMax, DimZMin:DimZMax) real(8), parameter :: GasRUniv = 8.314d0 !普遍気体定数 !飽和蒸気圧の温度微分 xz_DLogSvapPressDTemp = - a_ampA(ID) / (xz_Temp ** 2.0d0) + a_ampC(ID) / xz_Temp + a_ampD(ID) + a_ampE(ID) * 2.0d0 * xz_Temp !潜熱の計算 xz_LatentHeat = xz_DLogSvapPressDTemp * GasRUniv * (xz_Temp ** 2.0d0) / a_MolWt(ID) end function xz_LatentHeat
Function : | |||
xz_SvapPress(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8)
| ||
ID : | integer, intent(in)
| ||
xz_Temp(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
|
引数で与えられた化学種と温度に対して, 飽和蒸気圧を計算. AMP 式を利用
function xz_SvapPress(ID, xz_Temp) ! !引数で与えられた化学種と温度に対して, 飽和蒸気圧を計算. AMP 式を利用 ! !暗黙の型宣言禁止 implicit none !入出力変数 real(8) :: xz_SvapPress(DimXMin:DimXMax, DimZMin:DimZMax) !飽和蒸気圧 real(8), intent(in) :: xz_Temp(DimXMin:DimXMax, DimZMin:DimZMax) !温度 [K] integer, intent(in) :: ID !化学種の ID !内部変数 real(8) :: xz_LogSvapPress(DimXMin:DimXMax, DimZMin:DimZMax) !飽和蒸気圧の対数を計算 !対数が大きくなりすぎないようにする. ! Fujitsu Fortran Compiler では 700 より大きい数の exp を取ると警告が出る. xz_LogSvapPress = min( ( a_ampA(ID) / xz_Temp + a_ampB(ID) + a_ampC(ID) * dlog( xz_Temp ) + a_ampD(ID) * xz_Temp + a_ampE(ID) * ( xz_temp ** 2 ) + dlog(1.0d-1) ), 7.0d2 ) !飽和蒸気圧を計算 xz_SvapPress = dexp( xz_LogSvapPress ) end function xz_SvapPress