Class composition
In: setup/composition.f90

組成に関わる配列の設定

Settings of array for atmospheric composition

Note that Japanese and English are described in parallel.

組成に関わる配列の設定

Settings of array for atmospheric composition

Variables List

ncmax :微量成分の数
———— :————
ncmax :Number of minor comstituents

Procedures List

CompositionInit :組成に関わる配列の設定
————— :————
CompositionInit :Settings of parameters for array for composition

NAMELIST

NAMELIST#composition_nml

Methods

Included Modules

dc_types dc_message namelist_util dc_iounit

Public Instance methods

Subroutine :

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

"composition" module is initialized. NAMELIST#composition_nml is loaded in this procedure.

This procedure input/output NAMELIST#composition_nml .

[Source]

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

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

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

    ! ファイル入出力補助
    ! 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

    character(len=STRING) :: Names      (1:MaxNmlArySize)
    character(len=STRING) :: LongNames  (1:MaxNmlArySize)
    logical               :: FlagMassFix(1:MaxNmlArySize)
    logical               :: FlagAdv    (1:MaxNmlArySize)
    logical               :: FlagVDiff  (1:MaxNmlArySize)

    integer:: n               ! 組成方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in dimension of constituents

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /composition_nml/ ncmax, Names, LongNames, FlagMassFix, FlagAdv, FlagVDiff
          !
          ! デフォルト値については初期化手続 "composition#CompositionInit" 
          ! のソースコードを参照のこと. 
          !
          ! Refer to source codes in the initialization procedure
          ! "composition#CompositionInit" for the default values. 
          !

    ! 実行文 ; Executable statement
    !

    if ( composition_inited ) return
!!$    call InitCheck

    ! デフォルト値の設定
    ! Default values settings
    !
    ncmax        = 1

    Names    (1) = NameH2OVap
    LongNames(1) = 'specific humidity'
    do n = 2, MaxNmlArySize
      Names    (n) = ''
      LongNames(n) = ''
    end do

    FlagMassFix = .true.

    FlagAdv     = .true.

    FlagVDiff   = .true.


    ! NAMELIST の読み込み
    ! NAMELIST is input
    !
    if ( trim(namelist_filename) /= '' ) then
      call FileOpen( unit_nml, namelist_filename, mode = 'r' ) ! (in)

      rewind( unit_nml )
      read( unit_nml, nml = composition_nml, iostat = iostat_nml )
      close( unit_nml )

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

    ! 配列サイズのチェック
    ! Check number of array size
    !
    if ( ncmax < 1 ) then
      call MessageNotify( 'E', module_name, 'number of composition has to be greater than 0. ' // 'ncmax=%d' , i = (/ ncmax /) )
    end if


    ! Set variable name of constituents for ouptut
    !
    allocate( a_QMixName    ( ncmax ) )
    allocate( a_QMixLongName( ncmax ) )
    if ( ncmax > 99 ) then
      call MessageNotify( 'E', module_name, 'number of composition greater than 99 is inappropriate for current version. ' // 'ncmax=%d' , i = (/ ncmax /) )
    end if
    do n = 1, ncmax
      if ( Names(n) == '' ) then
        write( a_QMixName(n), '(a,i3.3)' ) "QMix", n
      else
        a_QMixName(n) = Names(n)
      end if
      if ( LongNames(n) == '' ) then
        a_QMixLongName(n) = a_QMixName(n)
      else
        a_QMixLongName(n) = LongNames(n)
      end if
    end do

    allocate( a_FlagMassFix( ncmax ) )
    do n = 1, ncmax
      a_FlagMassFix( n ) = FlagMassFix(n)
    end do

    allocate( a_FlagAdv( ncmax ) )
    do n = 1, ncmax
      a_FlagAdv( n ) = FlagAdv(n)
    end do

    allocate( a_FlagVDiff( ncmax ) )
    do n = 1, ncmax
      a_FlagVDiff( n ) = FlagVDiff(n)
    end do

    IndexH2OVap = -1
    do n = 1, ncmax
      if ( a_QMixName(n) == NameH2OVap ) then
        IndexH2OVap = n
        exit
      end if
    end do
    !
    ! 水蒸気のインデックスのチェック
    ! Check index for water vapor
    !
    if ( ( IndexH2OVap < 1 ) .or. ( IndexH2OVap > ncmax ) ) then
      call MessageNotify( 'E', module_name, 'IndexH2OVap has to be greater than or equal to 1 and less than or equal to ncmax. ' // 'IndexH2OVap=%d' , i = (/ IndexH2OVap /) )
    end if

    IndexH2OLiq = -1
    do n = 1, ncmax
      if ( a_QMixName(n) == NameH2OLiq ) then
        IndexH2OLiq = n
        exit
      end if
    end do
    IndexH2OSol = -1
    do n = 1, ncmax
      if ( a_QMixName(n) == NameH2OSol ) then
        IndexH2OSol = n
        exit
      end if
    end do
    IndexTKE = -1
    do n = 1, ncmax
      if ( a_QMixName(n) == NameTKE ) then
        IndexTKE = n
        exit
      end if
    end do

    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, '    ncmax              = %d', i = (/ ncmax /) )
    do n = 1, ncmax
      call MessageNotify( 'M', module_name, '    QMixName(%d)     = %c', i = (/ n /), c1 = trim(a_QMixName(n)) )
    end do
    do n = 1, ncmax
      call MessageNotify( 'M', module_name, '    QMixLongName(%d) = %c', i = (/ n /), c1 = trim(a_QMixLongName(n)) )
    end do
    do n = 1, ncmax
      call MessageNotify( 'M', module_name, '    FlagMassFix(%d)  = %b', i = (/ n /), l = (/ a_FlagMassFix(n) /) )
    end do
    do n = 1, ncmax
      call MessageNotify( 'M', module_name, '    FlagAdv(%d)      = %b', i = (/ n /), l = (/ a_FlagAdv(n) /) )
    end do
    do n = 1, ncmax
      call MessageNotify( 'M', module_name, '    FlagVDiff(%d)    = %b', i = (/ n /), l = (/ a_FlagVDiff(n) /) )
    end do
    call MessageNotify( 'M', module_name, '    IndexH2OVap        = %d', i = (/ IndexH2OVap /) )
    call MessageNotify( 'M', module_name, '    IndexH2OLiq        = %d', i = (/ IndexH2OLiq /) )
    call MessageNotify( 'M', module_name, '    IndexH2OSol        = %d', i = (/ IndexH2OSol /) )
    call MessageNotify( 'M', module_name, '    IndexTKE           = %d', i = (/ IndexTKE    /) )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    composition_inited = .true.
  end subroutine CompositionInit
Function :
FlagAdv :logical
Index :integer, intent(in )

Name で与えられた名前の成分の移流フラグを返します.

Inquiry of flag for advection of constituent

[Source]

  function CompositionInqFlagAdv( Index ) result( FlagAdv )
    !
    ! Name で与えられた名前の成分の移流フラグを返します. 
    !
    ! Inquiry of flag for advection of constituent
    !

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

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

    ! 宣言文 ; Declaration statements
    !
    implicit none

    integer, intent(in ) :: Index

    logical :: FlagAdv


    ! 実行文 ; Executable statement
    !

    ! 初期化
    ! Initialization
    !
    if ( .not. composition_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    if ( ( Index < 1 ) .or. ( Index > ncmax ) ) then
      call MessageNotify( 'E', module_name, 'Index is inappropriate, Index = %d', i = (/ Index /) )
    end if

    FlagAdv = a_FlagAdv( Index )


  end function CompositionInqFlagAdv
Function :
FlagMassFix :logical
Index :integer, intent(in )

Name で与えられた名前の成分の移流フラグを返しします.

Inquiry of flag for advection of constituent

[Source]

  function CompositionInqFlagMassFix( Index ) result( FlagMassFix )
    !
    ! Name で与えられた名前の成分の移流フラグを返しします.
    !
    ! Inquiry of flag for advection of constituent
    !

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

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

    ! 宣言文 ; Declaration statements
    !
    implicit none

    integer, intent(in ) :: Index

    logical :: FlagMassFix


    ! 実行文 ; Executable statement
    !

    ! 初期化
    ! Initialization
    !
    if ( .not. composition_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    if ( ( Index < 1 ) .or. ( Index > ncmax ) ) then
      call MessageNotify( 'E', module_name, 'Index is inappropriate, Index = %d', i = (/ Index /) )
    end if

    FlagMassFix = a_FlagMassFix( Index )


  end function CompositionInqFlagMassFix
Function :
FlagVDiff :logical
Index :integer, intent(in )

Name で与えられた名前の成分の移流フラグを返しします.

Inquiry of flag for advection of constituent

[Source]

  function CompositionInqFlagVDiff( Index ) result( FlagVDiff )
    !
    ! Name で与えられた名前の成分の移流フラグを返しします.
    !
    ! Inquiry of flag for advection of constituent
    !

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

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

    ! 宣言文 ; Declaration statements
    !
    implicit none

    integer, intent(in ) :: Index

    logical :: FlagVDiff


    ! 実行文 ; Executable statement
    !

    ! 初期化
    ! Initialization
    !
    if ( .not. composition_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    if ( ( Index < 1 ) .or. ( Index > ncmax ) ) then
      call MessageNotify( 'E', module_name, 'Index is inappropriate, Index = %d', i = (/ Index /) )
    end if

    FlagVDiff = a_FlagVDiff( Index )


  end function CompositionInqFlagVDiff
Function :
Index :integer
Name :character(len=*), intent(in )

Name で与えられた名前の成分のインデクスを返しします.

Inquiry of index of constituent

[Source]

  function CompositionInqIndex( Name ) result( Index )
    !
    ! Name で与えられた名前の成分のインデクスを返しします. 
    !
    ! Inquiry of index of constituent
    !

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

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

    ! 宣言文 ; Declaration statements
    !
    implicit none

    character(len=*), intent(in ) :: Name

    integer:: Index

    integer:: n               ! 組成方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in dimension of constituents

    ! 実行文 ; Executable statement
    !

    ! 初期化
    ! Initialization
    !
    if ( .not. composition_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if


    Index = -1

    do n = 1, ncmax
      if ( a_QMixName(n) == trim( Name ) ) then
        Index = n
        exit
      end if
    end do

!!$    if ( Index == -1 ) then
!!$      call MessageNotify( 'E', module_name, &
!!$        & 'Cannot find %c in a list of constituents.', c1 = trim( Name ) )
!!$    end if


  end function CompositionInqIndex
IndexH2OLiq
Variable :
IndexH2OLiq :integer , save, public
IndexH2OSol
Variable :
IndexH2OSol :integer , save, public
IndexH2OVap
Variable :
IndexH2OVap :integer , save, public
: 水蒸気のインデックス Index for water vapor
IndexTKE
Variable :
IndexTKE :integer , save, public
a_QMixLongName
Variable :
a_QMixLongName(:) :character(STRING), save, public, allocatable
: 成分の長い変数名 Long name of variables for composition
a_QMixName
Variable :
a_QMixName(:) :character(STRING), save, public, allocatable
: 成分の変数名 Name of variables for composition
composition_inited
Variable :
composition_inited = .false. :logical , save, public
: 初期設定フラグ. Initialization flag
ncmax
Variable :
ncmax :integer , save, public
: 成分の数 Number of composition

Private Instance methods

NameH2OLiq
Constant :
NameH2OLiq = "QH2OLiq" :character(STRING), parameter
NameH2OSol
Constant :
NameH2OSol = "QH2OSol" :character(STRING), parameter
NameH2OVap
Constant :
NameH2OVap = "QH2OVap" :character(STRING), parameter
NameTKE
Constant :
NameTKE = "TKE" :character(STRING), parameter
a_FlagAdv
Variable :
a_FlagAdv(:) :logical , save, allocatable
: 成分の移流計算のフラグ Flag for advection calculation
a_FlagMassFix
Variable :
a_FlagMassFix(:) :logical , save, allocatable
: 質量保存補正のフラグ Flag for mass fix
a_FlagVDiff
Variable :
a_FlagVDiff(:) :logical , save, allocatable
: 鉛直混合の移流計算のフラグ Flag for vertical diffusion
module_name
Constant :
module_name = ‘composition :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: $’ // ’$Id: composition.f90,v 1.6 2012/09/08 15:25:53 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version