Class | constants |
In: |
../src/setup/constants.f90
|
Subroutine : |
namelist の設定に基づいて物性を決める. 順序は以下の通り.
This procedure input/output NAMELIST#constants_nml .
subroutine constants_Init ! != 初期化ルーチン ! ! namelist の設定に基づいて物性を決める. 順序は以下の通り. ! ! * namelist から定圧比熱 (CpDry) と平均分子量 (MolWtDry) が与えられた場合は, ! それらを元に気体定数(GasRDry), 定積比熱 (CvDry), 定圧モル比熱 (CpDryMol) を決める. ! * namelist から定圧比熱 (CpDry) と気体定数 (GasRDry) が与えられた場合は, ! それらを元に平均分子量(MolWtDry), 定積比熱 (CvDry), 定圧モル比熱 (CpDryMol) を決める. ! * namelist から乾燥成分のモル比 (SpcDryMolFr) が与えられた場合には, ! それと熱力学テーブルより, 定圧比熱 (CpDry), 気体定数 (GasRDry), ! 平均分子量(MolWtDry), 定積比熱 (CvDry), 定圧モル比熱 (CpDryMol) を決める. ! !暗黙の型宣言禁止 implicit none !変数定義 integer :: SpcDryNum !乾燥成分の化学種の数 character(20) :: SpcDrySymbol(5) !乾燥成分の化学種名 real(DP) :: SpcDryMolFr(5) !乾燥成分の化学種の存在度 integer, allocatable :: SpcDryID(:) !乾燥成分の化学種のID real(DP), allocatable :: PropertyDry(:) !作業配列 integer :: s !作業変数 integer :: unit !装置番号 logical :: flag = .false. !NAMELIST の定義 NAMELIST /constants_nml/ Grav, PressBasis, TempSfc, PressSfc, TempTop, PressTop, SpcDrySymbol, SpcDryMolFr, DayTime, CpDry, MolWtDry, GasRDry SpcDrySymbol = '' SpcDryMolFr = 0.0d0 !ファイルオープン. 情報取得. call FileOpen(unit, file=namelist_filename, mode='r') read(unit, NML=constants_nml) close(unit) if (CpDry /= 0.0d0 .AND. MolWtDry /= 0.0d0) then ! namelist から CpDry の値が入力された場合 !定圧モル比熱 CpDryMol = CpDry * MolWtDry !気体定数 GasRDry = GasRUniv / MolWtDry !定積比熱 CvDry = CpDry - GasRDry elseif (CpDry /= 0.0d0 .AND. GasRDry /= 0.0d0) then ! namelist から CpDry と GasRDry の値が入力された場合 !乾燥成分の分子量 MolWtDry = GasRUniv / GasRDry !定圧モル比熱 CpDryMol = CpDry * MolWtDry !定積比熱 CvDry = CpDry - GasRDry elseif (SpcDryMolFr(1) /= 0.0d0) then ! namelist から モル比が入力された場合 flag = .true. !---------------------------------------------------------- ! 乾燥成分の物性値の初期化 ! !乾燥成分の個数を数える SpcDryNum = count(SpcDrySymbol /= "") !化学種の ID を取得 allocate(SpcDryID(SpcDryNum)) do s = 1, SpcDryNum SpcDryID(s) = ChemData_OneSpcID( SpcDrySymbol(s) ) end do !作業配列の準備 allocate(PropertyDry(SpcDryNum)) !分子量 do s = 1, SpcDryNum PropertyDry(s) = ChemData_MolWt(SpcDryID(s)) end do MolWtDry = dot_product(PropertyDry, SpcDryMolFr(1:SpcDryNum)) !定圧比熱(モル当量) do s = 1, SpcDryNum PropertyDry(s) = ChemData_CpPerMolRef(SpcDryID(s)) end do CpDryMol = dot_product(PropertyDry, SpcDryMolFr(1:SpcDryNum)) !定圧比熱 CpDry = CpDryMol / MolWtDry !気体定数 GasRDry = GasRUniv / MolWtDry !定積比熱 CvDry = CpDry - GasRDry else call MessageNotify( "E", "constants_init", "Enough Variables are not set" ) end if !---------------------------------------------------------- ! 確認 !---------------------------------------------------------- if (myrank == 0) then call MessageNotify( "M", "constants_init", "Grav = %f", d=(/Grav/) ) call MessageNotify( "M", "constants_init", "PressBasis = %f", d=(/PressBasis/)) if (TempSfc /= 0) then ! 地表面温度が与えられた場合 ! call MessageNotify( "M", "constants_init", "TempSfc = %f", d=(/TempSfc/) ) call MessageNotify( "M", "constants_init", "PressSfc = %f", d=(/PressSfc/) ) end if if (TempTop /= 0) then ! 上部境界の温度が与えられた場合 ! call MessageNotify( "M", "constants_init", "TempTop = %f", d=(/TempTop/) ) call MessageNotify( "M", "constants_init", "PressTop = %f", d=(/PressTop/) ) end if if (flag) then ! モル比が与えられた場合にはプロットする. ! do s = 1, SpcDryNum call MessageNotify( "M", "constants_init", "SpcDryID = %d", i=(/SpcDryID(s)/)) call MessageNotify( "M", "constants_init", "SpcDrySymbol = %c", c1=trim(SpcDrySymbol(s))) call MessageNotify( "M", "constants_init", "SpcDryMolFr = %f", d=(/SpcDryMolFr(s)/)) end do end if call MessageNotify( "M", "constants_init", "CpDry = %f", d=(/CpDry/) ) call MessageNotify( "M", "constants_init", "CpDryMol = %f", d=(/CpDryMol/) ) call MessageNotify( "M", "constants_init", "CvDry = %f", d=(/CvDry/) ) call MessageNotify( "M", "constants_init", "GasRDry = %f", d=(/GasRDry/) ) call MessageNotify( "M", "constants_init", "MolWtDry = %f", d=(/MolWtDry/) ) call MessageNotify( "M", "constants_init", "DayTime = %f", d=(/DayTime/) ) end if end subroutine Constants_Init