Class constants
In: setup/constants.f90

物理定数設定

Physical constants settings

Note that Japanese and English are described in parallel.

物理定数の設定および保管を行います. デフォルト値は地球大気を想定した値が設定されています. これらの値は NAMELIST 変数群 constants_nml によって変更することが 可能です.

Physical constants are set and stored. By default, values on atmosphere of earth are set. These values can be changed by NAMELIST group name "constants_nml".

Procedures List

ConstantsInit :物理定数の設定
———— :————
ConstantsInit :Settings of physical constants

NAMELIST

NAMELIST#constants_nml

Methods

Included Modules

dc_types namelist_util dc_iounit dc_message

Public Instance methods

Subroutine :

constants モジュールの初期化を行います. NAMELIST#constants_nml の読み込みはこの手続きで行われます.

"constants" module is initialized. NAMELIST#constants_nml is loaded in this procedure.

This procedure input/output NAMELIST#constants_nml .

[Source]

  subroutine ConstantsInit
    !
    ! constants モジュールの初期化を行います. 
    ! NAMELIST#constants_nml の読み込みはこの手続きで行われます. 
    !
    ! "constants" module is initialized. 
    ! NAMELIST#constants_nml is loaded in this procedure. 
    !

    ! モジュール引用 ; USE statements
    !

    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_filename, NmlutilMsg

    ! ファイル入出力補助
    ! File I/O support
    !
    use dc_iounit, only: FileOpen

    ! 種別型パラメタ
    ! Kind type parameter
    !
    use dc_types, only: STDOUT ! 標準出力の装置番号. Unit number of standard output

    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify

    ! 宣言文 ; Declaration statements
    !
    implicit none

    integer:: unit_nml        ! NAMELIST ファイルオープン用装置番号. 
                              ! Unit number for NAMELIST file open
    integer:: iostat_nml      ! NAMELIST 読み込み時の IOSTAT. 
                              ! IOSTAT of NAMELIST read

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /constants_nml/ RPlanet, Omega, Grav, MolWtDry, CpDry, GasRDry, MolWtWet, CpWet, GasRWet, LatentHeat, EpsV
          !
          ! デフォルト値については初期化手続 "constants#ConstantsInit" 
          ! のソースコードを参照のこと. 
          !
          ! Refer to source codes in the initialization procedure
          ! "constants#ConstantsInit" for the default values. 
          !

    ! 実行文 ; Executable statement
    !

    if ( constants_inited ) return
    call InitCheck

    ! デフォルト値の設定
    ! Default values settings
    !
    RPlanet          = 6.371e6_DP
    Omega            = 2.0_DP * PI / ( 60.0_DP * 60.0_DP * 23.9345_DP )
    Grav             = 9.8_DP

    CpDry            = 1004.6_DP
    MolWtDry         = 28.964e-3_DP
    GasRDry          = GasRUniv / MolWtDry

    CpWet            = 1810.0_DP
    MolWtWet         = 18.01528e-3_DP
    GasRWet          = GasRUniv / MolWtWet
    LatentHeat       = 2.5e6_DP 
    LatentHeatFusion = 334.0d3

    EpsV             = MolWtWet / MolWtDry

    ! NAMELIST からの入力
    ! Input from NAMELIST
    !
    if ( trim(namelist_filename) /= '' ) then
      call FileOpen( unit_nml, namelist_filename, mode = 'r' ) ! (in)

      rewind( unit_nml )
      read( unit_nml, nml = constants_nml, iostat = iostat_nml )   ! (out)
      close( unit_nml )

      call NmlutilMsg( iostat_nml, module_name ) ! (in)
      if ( iostat_nml == 0 ) write( STDOUT, nml = constants_nml )
    end if

    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, '  PI               = %f', d = (/ PI               /) )
    call MessageNotify( 'M', module_name, '  GasRUniv         = %f', d = (/ GasRUniv         /) )
    call MessageNotify( 'M', module_name, '  StB              = %f', d = (/ StB              /) )
    call MessageNotify( 'M', module_name, '  FKarm            = %f', d = (/ FKarm            /) )
    
    call MessageNotify( 'M', module_name, '  RPlanet          = %f', d = (/ RPlanet          /) )
    call MessageNotify( 'M', module_name, '  Omega            = %f', d = (/ Omega            /) )
    call MessageNotify( 'M', module_name, '  Grav             = %f', d = (/ Grav             /) )
    
    call MessageNotify( 'M', module_name, '  CpDry            = %f', d = (/ CpDry            /) )
    call MessageNotify( 'M', module_name, '  MolWtDry         = %f', d = (/ MolWtDry         /) )
    call MessageNotify( 'M', module_name, '  GasRDry          = %f', d = (/ GasRDry          /) )
    
    call MessageNotify( 'M', module_name, '  CpWet            = %f', d = (/ CpWet            /) )
    call MessageNotify( 'M', module_name, '  MolWtWet         = %f', d = (/ MolWtWet         /) )
    call MessageNotify( 'M', module_name, '  GasRWet          = %f', d = (/ GasRWet          /) )
    call MessageNotify( 'M', module_name, '  LatentHeat       = %f', d = (/ LatentHeat       /) )
    call MessageNotify( 'M', module_name, '  LatentHeatFusion = %f', d = (/ LatentHeatFusion /) )

    call MessageNotify( 'M', module_name, '  EpsV       = %f', d = (/ EpsV       /) )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    constants_inited = .true.
  end subroutine ConstantsInit
CpDry
Variable :
CpDry :real(DP), save, public
: $ C_p $ [J kg-1 K-1]. 乾燥大気の定圧比熱. Specific heat of air at constant pressure
CpWet
Variable :
CpWet :real(DP), save, public
: $ C_v $ [J kg-1 K-1] . 凝結成分の定圧比熱. Specific heat of condensible elements at constant pressure
EpsV
Variable :
EpsV :real(DP), save, public
: $ epsilon_v $ . 水蒸気分子量比. Molecular weight of water vapor
FKarm
Constant :
FKarm = 0.4_DP :real(DP), parameter, public
: $ k $ . カルマン定数. Karman constant
GasRDry
Variable :
GasRDry :real(DP), save, public
: $ R $ [J kg-1 K-1]. 乾燥大気の気体定数. Gas constant of air
GasRUniv
Constant :
GasRUniv = 8.314_DP :real(DP), parameter, public
: $ R^{*} $ [J K-1 mol-1]. 普遍気体定数. Universal gas constant
GasRWet
Variable :
GasRWet :real(DP), save, public
: $ R_v $ [J kg-1 K-1]. 凝結成分の気体定数. Gas constant of condensible elements
Grav
Variable :
Grav :real(DP), save, public
: $ g $ [m s-2]. 重力加速度. Gravitational acceleration
LatentHeat
Variable :
LatentHeat :real(DP), save, public
: $ L $ [J kg-1] . 凝結の潜熱. Latent heat of condensation
LatentHeatFusion
Variable :
LatentHeatFusion :real(DP), save, public
: $ L $ [J kg-1] . 融解の潜熱. Latent heat of fusion
MolWtDry
Variable :
MolWtDry :real(DP), save, public
: $ M $ [kg mol-1]. 乾燥大気の平均分子量. Mean molecular weight of dry air
MolWtWet
Variable :
MolWtWet :real(DP), save, public
: $ M_v $ [kg mol-1]. 凝結成分の平均分子量. Mean molecular weight of condensible elements
Omega
Variable :
Omega :real(DP), save, public
: $ Omega $ [s-1]. 回転角速度. Angular velocity
PI
Constant :
PI = 3.1415926535897932_DP :real(DP), parameter, public
: $ pi $ . 円周率. Circular constant
RPlanet
Variable :
RPlanet :real(DP), save, public
: $ a $ [m]. 惑星半径. Radius of planet
StB
Constant :
StB = 5.67e-8_DP :real(DP), parameter, public
: $ sigma_{SB} $ . ステファンボルツマン定数. Stefan-Boltzmann constant
constants_inited
Variable :
constants_inited = .false. :logical, save, public
: 初期設定フラグ. Initialization flag

Private Instance methods

Subroutine :

依存モジュールの初期化チェック

Check initialization of dependency modules

[Source]

  subroutine InitCheck
    !
    ! 依存モジュールの初期化チェック
    !
    ! Check initialization of dependency modules

    ! モジュール引用 ; USE statements
    !

    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_util_inited

    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify

    ! 実行文 ; Executable statement
    !

    if ( .not. namelist_util_inited ) call MessageNotify( 'E', module_name, '"namelist_util" module is not initialized.' )

  end subroutine InitCheck
module_name
Constant :
module_name = ‘constants :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: dcpam5-20100224 $’ // ’$Id: constants.f90,v 1.4 2009-05-24 02:36:07 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version

[Validate]